/[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 2990 - (show annotations)
Fri Mar 19 05:57:22 2010 UTC (8 years, 10 months ago) by gross
File MIME type: text/plain
File size: 14712 byte(s)
some prints removed.
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 requireWrite(lumpedMat);
106 lumpedMat_p=getSampleDataRW(lumpedMat,0);
107 len_EM_lumpedMat=p.row_numShapesTotal*p.numEqu;
108 len_EM_lumpedMat_size=len_EM_lumpedMat*sizeof(double);
109
110 expandedD=isExpanded(D);
111 S=p.row_jac->BasisFunctions->S;
112
113 #ifdef NEW_LUMPING
114 #pragma omp parallel private(color, EM_lumpedMat, row_index, Vol, D_p, s, q, k, rtmp, diagS, m_t)
115 #else
116 #pragma omp parallel private(color, EM_lumpedMat, row_index, Vol, D_p, s, q, k, rtmp)
117 #endif
118 {
119 EM_lumpedMat=THREAD_MEMALLOC(len_EM_lumpedMat,double);
120 row_index=THREAD_MEMALLOC(p.row_numShapesTotal,index_t);
121 if ( !Finley_checkPtr(EM_lumpedMat) && !Finley_checkPtr(row_index) ) {
122 if (p.numEqu == 1) { /* single equation */
123 if (expandedD) { /* with expanded D */
124 for (color=elements->minColor;color<=elements->maxColor;color++) {
125 /* open loop over all elements: */
126 #pragma omp for private(e) schedule(static)
127 for(e=0;e<elements->numElements;e++){
128 if (elements->Color[e]==color) {
129 for (isub=0; isub<p.numSub; isub++) {
130 Vol=&(p.row_jac->volume[INDEX3(0,isub,e, p.numQuadSub,p.numSub)]);
131 D_p=getSampleDataRO(D,e);
132 #ifdef NEW_LUMPING /* HRZ lumping */
133
134 m_t=0; /* mass of the element: m_t */
135 for (q=0;q<p.numQuadSub;q++) m_t+=Vol[q]*D_p[INDEX2(q, isub,p.numQuadSub) ];
136
137 diagS=0; /* diagonal sum: S */
138 for (s=0;s<p.row_numShapes;s++) {
139 rtmp=0;
140 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)];
141 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp;
142 diagS+=rtmp;
143 }
144 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
145 rtmp=m_t/diagS;
146 for (s=0;s<p.row_numShapes;s++) EM_lumpedMat[INDEX2(0,s,p.numEqu)]*=rtmp;
147
148 #else /* row-sum lumping */
149 for (s=0;s<p.row_numShapes;s++) {
150 rtmp=0;
151 for (q=0;q<p.numQuadSub;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_numShapes)]*D_p[INDEX2(q, isub,p.numQuadSub)];
152 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp;
153 }
154
155 #endif
156 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)]];
157 Finley_Util_AddScatter(p.row_numShapesTotal,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
158 } /* end of isub loop */
159 } /* end color check */
160 } /* end element loop */
161 } /* end color loop */
162 } else { /* with constant D */
163
164 for (color=elements->minColor;color<=elements->maxColor;color++) {
165 /* open loop over all elements: */
166 #pragma omp for private(e) schedule(static)
167 for(e=0;e<elements->numElements;e++){
168 if (elements->Color[e]==color) {
169 for (isub=0; isub<p.numSub; isub++) {
170 Vol=&(p.row_jac->volume[INDEX3(0,isub,e, p.numQuadSub,p.numSub)]);
171 D_p=getSampleDataRO(D,e);
172 #ifdef NEW_LUMPING /* HRZ lumping */
173 m_t=0; /* mass of the element: m_t */
174 for (q=0;q<p.numQuadSub;q++) m_t+=Vol[q];
175 diagS=0; /* diagonal sum: S */
176 for (s=0;s<p.row_numShapes;s++) {
177 rtmp=0;
178 for (q=0;q<p.numQuadSub;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_numShapes)]*S[INDEX2(s,q,p.row_numShapes)];
179 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp;
180 diagS+=rtmp;
181 }
182 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
183 rtmp=m_t/diagS*D_p[0];
184 for (s=0;s<p.row_numShapes;s++) EM_lumpedMat[INDEX2(0,s,p.numEqu)]*=rtmp;
185 #else /* row-sum lumping */
186 for (s=0;s<p.row_numShapes;s++) {
187 rtmp=0;
188 for (q=0;q<p.numQuadSub;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_numShapes)];
189 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp*D_p[0];
190 }
191 #endif
192 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)]];
193 Finley_Util_AddScatter(p.row_numShapesTotal,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
194 } /* end of isub loop */
195 } /* end color check */
196 } /* end element loop */
197 } /* end color loop */
198
199 }
200 } else { /* system of equation */
201 if (expandedD) { /* with expanded D */
202 for (color=elements->minColor;color<=elements->maxColor;color++) {
203 /* open loop over all elements: */
204 #pragma omp for private(e) schedule(static)
205 for(e=0;e<elements->numElements;e++){
206 if (elements->Color[e]==color) {
207 for (isub=0; isub<p.numSub; isub++) {
208 Vol=&(p.row_jac->volume[INDEX3(0,isub,e, p.numQuadSub,p.numSub)]);
209 D_p=getSampleDataRO(D,e);
210
211 #ifdef NEW_LUMPING /* HRZ lumping */
212 for (k=0;k<p.numEqu;k++) {
213 m_t=0; /* mass of the element: m_t */
214 for (q=0;q<p.numQuadSub;q++) m_t+=Vol[q]*D_p[INDEX3(k,q,isub,p.numEqu,p.numQuadSub)];
215
216 diagS=0; /* diagonal sum: S */
217 for (s=0;s<p.row_numShapes;s++) {
218 rtmp=0;
219 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)];
220 EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp;
221 diagS+=rtmp;
222 }
223 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
224 rtmp=m_t/diagS;
225 for (s=0;s<p.row_numShapes;s++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]*=rtmp;
226 }
227 #else /* row-sum lumping */
228 for (s=0;s<p.row_numShapes;s++) {
229 for (k=0;k<p.numEqu;k++) {
230 rtmp=0.;
231 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)];
232 EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp;
233 }
234 }
235 #endif
236 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)]];
237 Finley_Util_AddScatter(p.row_numShapesTotal,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
238 } /* end of isub loop */
239 } /* end color check */
240 } /* end element loop */
241 } /* end color loop */
242 } else { /* with constant D */
243 for (color=elements->minColor;color<=elements->maxColor;color++) {
244 /* open loop over all elements: */
245 #pragma omp for private(e) schedule(static)
246 for(e=0;e<elements->numElements;e++){
247 if (elements->Color[e]==color) {
248 for (isub=0; isub<p.numSub; isub++) {
249 Vol=&(p.row_jac->volume[INDEX3(0,isub,e, p.numQuadSub,p.numSub)]);
250 D_p=getSampleDataRO(D,e);
251
252 #ifdef NEW_LUMPING /* HRZ lumping */
253 m_t=0; /* mass of the element: m_t */
254 for (q=0;q<p.numQuadSub;q++) m_t+=Vol[q];
255 diagS=0; /* diagonal sum: S */
256 for (s=0;s<p.row_numShapes;s++) {
257 rtmp=0;
258 for (q=0;q<p.numQuadSub;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_numShapes)]*S[INDEX2(s,q,p.row_numShapes)];
259 for (k=0;k<p.numEqu;k++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp;
260 diagS+=rtmp;
261 }
262 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
263 rtmp=m_t/diagS;
264 for (s=0;s<p.row_numShapes;s++) {
265 for (k=0;k<p.numEqu;k++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]*=rtmp*D_p[k];
266 }
267 #else /* row-sum lumping */
268 for (s=0;s<p.row_numShapes;s++) {
269 for (k=0;k<p.numEqu;k++) {
270 rtmp=0.;
271 for (q=0;q<p.numQuadSub;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_numShapes)];
272 EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp*D_p[k];
273 }
274 }
275 #endif
276 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)]];
277 Finley_Util_AddScatter(p.row_numShapesTotal,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
278 } /* end of isub loop */
279 } /* end color check */
280 } /* end element loop */
281 } /* end color loop */
282 }
283 }
284 } /* end of pointer check */
285 THREAD_MEMFREE(EM_lumpedMat);
286 THREAD_MEMFREE(row_index);
287 } /* end parallel region */
288 }
289 }

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26