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

Contents of /trunk/finley/src/Assemble_LumpedSystem.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2989 - (show annotations)
Thu Mar 18 01:45:55 2010 UTC (8 years, 10 months ago) by gross
File MIME type: text/plain
File size: 14867 byte(s)
- bug in error handling in the case of lumped mass matrix fixed.
- switched back to rowsum lumping as there are problems with HRZ lumping for order 2 


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

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26