/[escript]/trunk-mpi-branch/finley/src/Assemble_LumpedSystem.c
ViewVC logotype

Contents of /trunk-mpi-branch/finley/src/Assemble_LumpedSystem.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1295 - (show annotations)
Mon Sep 10 06:07:09 2007 UTC (12 years, 8 months ago) by ksteube
File MIME type: text/plain
File size: 15081 byte(s)
Have now merged latest trunk features into MPI branch in preparation for
ending the MPI branch.
Compiles but has run time problems in bandwith reduction.

1 /*
2 ************************************************************
3 * Copyright 2006,2007 by ACcENULLNULL MNRF *
4 * *
5 * http://www.access.edu.au *
6 * Primary Business: Queensland, Australia *
7 * Licensed under the Open NULLoftware License version 3.0 *
8 * http://www.opensource.org/licenses/osl-3.0.php *
9 * *
10 ************************************************************
11 */
12
13
14 /**************************************************************/
15
16 /* assembles the mass matrix in lumped form */
17
18 /* The coefficient D has to be defined on the integration points or not present. */
19
20 /* lumpedMat has to be initialized before the routine is called. */
21
22 /**************************************************************/
23
24 /* Author: gross@access.edu.au */
25 /* Version: $Id$ */
26
27 /**************************************************************/
28
29 #include "Assemble.h"
30 #include "Util.h"
31 #ifdef _OPENMP
32 #include <omp.h>
33 #endif
34
35
36 /* Disabled until the tests pass */
37 /* #define NEW_LUMPING */ /* */
38
39 /**************************************************************/
40
41 void Finley_Assemble_LumpedSystem(Finley_NodeFile* nodes,Finley_ElementFile* elements, escriptDataC* lumpedMat, escriptDataC* D)
42 {
43
44 bool_t reducedIntegrationOrder=FALSE, expandedD;
45 char error_msg[LenErrorMsg_MAX];
46 Assemble_Parameters p;
47 double time0;
48 dim_t dimensions[ESCRIPT_MAX_DATA_RANK], k, e, len_EM_lumpedMat, q, s;
49 type_t funcspace;
50 index_t color,*row_index=NULL;
51 double *S=NULL, *EM_lumpedMat=NULL, *Vol=NULL, *D_p=NULL, *lumpedMat_p=NULL;
52 register double rtmp, m_t, diagS;
53 size_t len_EM_lumpedMat_size;
54
55 Finley_resetError();
56
57 if (nodes==NULL || elements==NULL) return;
58 if (isEmpty(lumpedMat) || isEmpty(D)) return;
59 if (isEmpty(lumpedMat) && !isEmpty(D)) {
60 Finley_setError(TYPE_ERROR,"Assemble_LumpedSystem: coefficients are non-zero but no lumped matrix is given.");
61 return;
62 }
63 funcspace=getFunctionSpaceType(D);
64 /* check if all function spaces are the same */
65 if (funcspace==FINLEY_ELEMENTS) {
66 reducedIntegrationOrder=FALSE;
67 } else if (funcspace==FINLEY_FACE_ELEMENTS) {
68 reducedIntegrationOrder=FALSE;
69 } else if (funcspace==FINLEY_REDUCED_ELEMENTS) {
70 reducedIntegrationOrder=TRUE;
71 } else if (funcspace==FINLEY_REDUCED_FACE_ELEMENTS) {
72 reducedIntegrationOrder=TRUE;
73 } else {
74 Finley_setError(TYPE_ERROR,"Assemble_LumpedSystem: assemblage failed because of illegal function space.");
75 }
76 if (! Finley_noError()) return;
77
78 /* set all parameters in p*/
79 Assemble_getAssembleParameters(nodes,elements,NULL,lumpedMat, reducedIntegrationOrder, &p);
80 if (! Finley_noError()) return;
81
82 /* check if all function spaces are the same */
83
84 if (! numSamplesEqual(D,p.numQuad,elements->numElements) ) {
85 sprintf(error_msg,"Assemble_LumpedSystem: sample points of coefficient D don't match (%d,%d)",p.numQuad,elements->numElements);
86 Finley_setError(TYPE_ERROR,error_msg);
87 }
88
89 /* check the dimensions: */
90
91 if (p.numEqu==1) {
92 if (!isEmpty(D)) {
93 if (!isDataPointShapeEqual(D,0,dimensions)) {
94 Finley_setError(TYPE_ERROR,"Assemble_LumpedSystem: coefficient D, rank 0 expected.");
95 }
96
97 }
98 } else {
99 if (!isEmpty(D)) {
100 dimensions[0]=p.numEqu;
101 if (!isDataPointShapeEqual(D,1,dimensions)) {
102 sprintf(error_msg,"Assemble_LumpedSystem: coefficient D, expected shape (%d,)",dimensions[0]);
103 Finley_setError(TYPE_ERROR,error_msg);
104 }
105 }
106 }
107
108 if (Finley_noError()) {
109 lumpedMat_p=getSampleData(lumpedMat,0);
110 len_EM_lumpedMat=p.row_NN*p.numEqu;
111 len_EM_lumpedMat_size=len_EM_lumpedMat*sizeof(double);
112 expandedD=isExpanded(D);
113 S=p.row_jac->ReferenceElement->S;
114
115 #pragma omp parallel private(color, EM_lumpedMat, row_index, Vol, D_p, s, q, k, rtmp, m_t, diagS)
116 {
117 EM_lumpedMat=THREAD_MEMALLOC(len_EM_lumpedMat,double);
118 row_index=THREAD_MEMALLOC(p.row_NN,index_t);
119 if ( !Finley_checkPtr(EM_lumpedMat) && !Finley_checkPtr(row_index) ) {
120 if (p.numEqu == 1) {
121 if (expandedD) {
122 for (color=elements->minColor;color<=elements->maxColor;color++) {
123 /* open loop over all elements: */
124 #pragma omp for private(e) schedule(static)
125 for(e=0;e<elements->numElements;e++){
126 if (elements->Color[e]==color) {
127 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
128 memset(EM_lumpedMat,0,len_EM_lumpedMat_size);
129 D_p=getSampleData(D,e);
130 #ifdef NEW_LUMPING /* HRZ lumping */
131 /*
132 * Number of PDEs: 1
133 * D_p varies over element: True
134 */
135 m_t=0; /* mass of the element: m_t */
136 for (q=0;q<p.numQuad;q++) {
137 m_t+=Vol[q]*D_p[q];
138 }
139 diagS=0; /* diagonal sum: S */
140 for (s=0;s<p.row_NS;s++) {
141 rtmp=0;
142 for (q=0;q<p.numQuad;q++) {
143 rtmp+=Vol[q]*D_p[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(s,q,p.row_NS)];
144 }
145 EM_lumpedMat[INDEX2(0,s,p.numEqu)]+=rtmp;
146 diagS+=EM_lumpedMat[INDEX2(0,s,p.numEqu)];
147 }
148 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
149 for (s=0;s<p.row_NS;s++) {
150 EM_lumpedMat[INDEX2(0,s,p.numEqu)]*=m_t/diagS;
151 }
152 #else /* row-sum lumping */
153 for (s=0;s<p.row_NS;s++) {
154 rtmp=0;
155 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*D_p[q];
156 EM_lumpedMat[INDEX2(0,s,p.numEqu)]+=rtmp;
157 }
158 #endif
159 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
160 Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
161 } /* end color check */
162 } /* end element loop */
163 } /* end color loop */
164 } else {
165 for (color=elements->minColor;color<=elements->maxColor;color++) {
166 /* open loop over all elements: */
167 #pragma omp for private(e) schedule(static)
168 for(e=0;e<elements->numElements;e++){
169 if (elements->Color[e]==color) {
170 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
171 memset(EM_lumpedMat,0,len_EM_lumpedMat_size);
172 D_p=getSampleData(D,e);
173 #ifdef NEW_LUMPING /* HRZ lumping */
174 /*
175 * Number of PDEs: 1
176 * D_p varies over element: False
177 */
178 m_t=0; /* mass of the element: m_t */
179 for (q=0;q<p.numQuad;q++) {
180 m_t+=Vol[q]*D_p[0];
181 }
182 diagS=0; /* diagonal sum: S */
183 for (s=0;s<p.row_NS;s++) {
184 rtmp=0;
185 for (q=0;q<p.numQuad;q++) {
186 rtmp+=Vol[q]*D_p[0]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(s,q,p.row_NS)];
187 }
188 EM_lumpedMat[INDEX2(0,s,p.numEqu)]+=rtmp;
189 diagS+=EM_lumpedMat[INDEX2(0,s,p.numEqu)];
190 }
191 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
192 for (s=0;s<p.row_NS;s++) {
193 EM_lumpedMat[INDEX2(0,s,p.numEqu)]*=m_t/diagS;
194 }
195 #else /* row-sum lumping */
196 for (s=0;s<p.row_NS;s++) {
197 rtmp=0;
198 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
199 EM_lumpedMat[INDEX2(0,s,p.numEqu)]+=rtmp*D_p[0];
200 }
201 #endif
202 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
203 Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
204 } /* end color check */
205 } /* end element loop */
206 } /* end color loop */
207 }
208 } else {
209 if (expandedD) {
210 for (color=elements->minColor;color<=elements->maxColor;color++) {
211 /* open loop over all elements: */
212 #pragma omp for private(e) schedule(static)
213 for(e=0;e<elements->numElements;e++){
214 if (elements->Color[e]==color) {
215 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
216 memset(EM_lumpedMat,0,len_EM_lumpedMat_size);
217 D_p=getSampleData(D,e);
218 #ifdef NEW_LUMPING /* HRZ lumping */
219 /*
220 * Number of PDEs: Multiple
221 * D_p varies over element: True
222 */
223 for (k=0;k<p.numEqu;k++) {
224 m_t=0; /* mass of the element: m_t */
225 for (q=0;q<p.numQuad;q++) {
226 m_t+=Vol[q]*D_p[INDEX2(k,q,p.numEqu)];
227 }
228 diagS=0; /* diagonal sum: S */
229 for (s=0;s<p.row_NS;s++) {
230 rtmp=0;
231 for (q=0;q<p.numQuad;q++) {
232 rtmp+=Vol[q]*D_p[INDEX2(k,q,p.numEqu)]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(s,q,p.row_NS)];
233 }
234 EM_lumpedMat[INDEX2(k,s,p.numEqu)]+=rtmp;
235 diagS+=EM_lumpedMat[INDEX2(k,s,p.numEqu)];
236 }
237 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
238 for (s=0;s<p.row_NS;s++) {
239 EM_lumpedMat[INDEX2(k,s,p.numEqu)]*=m_t/diagS;
240 }
241 }
242 #else /* row-sum lumping */
243 for (s=0;s<p.row_NS;s++) {
244 for (k=0;k<p.numEqu;k++) {
245 rtmp=0.;
246 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*D_p[INDEX2(k,q,p.numEqu)];
247 EM_lumpedMat[INDEX2(k,s,p.numEqu)]+=rtmp;
248 }
249 }
250 #endif
251 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
252 Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
253 } /* end color check */
254 } /* end element loop */
255 } /* end color loop */
256 } else {
257 /* open loop over all elements: */
258 #pragma omp for private(e) schedule(static)
259 for(e=0;e<elements->numElements;e++){
260 if (elements->Color[e]==color) {
261 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
262 memset(EM_lumpedMat,0,len_EM_lumpedMat_size);
263 D_p=getSampleData(D,e);
264 #ifdef NEW_LUMPING /* HRZ lumping */
265 /*
266 * Number of PDEs: Multiple
267 * D_p varies over element: False
268 */
269 for (k=0;k<p.numEqu;k++) {
270 m_t=0; /* mass of the element: m_t */
271 for (q=0;q<p.numQuad;q++) {
272 m_t+=Vol[q]*D_p[k];
273 }
274 diagS=0; /* diagonal sum: S */
275 for (s=0;s<p.row_NS;s++) {
276 rtmp=0;
277 for (q=0;q<p.numQuad;q++) {
278 rtmp+=Vol[q]*D_p[k]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(s,q,p.row_NS)];
279 }
280 EM_lumpedMat[INDEX2(k,s,p.numEqu)]+=rtmp;
281 diagS+=EM_lumpedMat[INDEX2(k,s,p.numEqu)];
282 }
283 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
284 for (s=0;s<p.row_NS;s++) {
285 EM_lumpedMat[INDEX2(k,s,p.numEqu)]*=m_t/diagS;
286 }
287 }
288 #else /* row-sum lumping */
289 for (s=0;s<p.row_NS;s++) {
290 rtmp=0;
291 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
292 for (k=0;k<p.numEqu;k++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]+=rtmp*D_p[k];
293 }
294 #endif
295 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
296 Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
297 } /* end color check */
298 } /* end element loop */
299 }
300 }
301 } /* end of pointer check */
302 THREAD_MEMFREE(EM_lumpedMat);
303 THREAD_MEMFREE(row_index);
304 } /* end parallel region */
305 }
306 #ifdef Finley_TRACE
307 printf("timing: assemblage lumped PDE: %.4e sec\n",Finley_timer()-time0);
308 #endif
309 }

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26