/[escript]/branches/domexper/dudley/src/Assemble_LumpedSystem.c
ViewVC logotype

Contents of /branches/domexper/dudley/src/Assemble_LumpedSystem.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3221 - (show annotations)
Wed Sep 29 01:00:21 2010 UTC (8 years, 6 months ago) by jfenwick
File MIME type: text/plain
File size: 13896 byte(s)
Comment stripping

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 #include "ShapeTable.h"
32
33 /* Disabled until the tests pass */
34 /* #define NEW_LUMPING */
35
36 /**************************************************************/
37
38 void Dudley_Assemble_LumpedSystem(Dudley_NodeFile* nodes,Dudley_ElementFile* elements, escriptDataC* lumpedMat, escriptDataC* D)
39 {
40
41 bool_t reducedIntegrationOrder=FALSE, expandedD;
42 char error_msg[LenErrorMsg_MAX];
43 Assemble_Parameters p;
44 dim_t dimensions[ESCRIPT_MAX_DATA_RANK], k, e, len_EM_lumpedMat, q, s;
45 type_t funcspace;
46 index_t color,*row_index=NULL;
47 __const double *D_p=NULL;
48 const double* S=NULL;
49 double *EM_lumpedMat=NULL, *lumpedMat_p=NULL;
50 register double rtmp;
51 size_t len_EM_lumpedMat_size;
52 #ifdef NEW_LUMPING
53 register double m_t=0., diagS=0.;
54 #endif
55
56 Dudley_resetError();
57
58 if (nodes==NULL || elements==NULL) return;
59 if (isEmpty(lumpedMat) || isEmpty(D)) return;
60 if (isEmpty(lumpedMat) && !isEmpty(D)) {
61 Dudley_setError(TYPE_ERROR,"Assemble_LumpedSystem: coefficients are non-zero but no lumped matrix is given.");
62 return;
63 }
64 funcspace=getFunctionSpaceType(D);
65 /* check if all function spaces are the same */
66 if (funcspace==DUDLEY_ELEMENTS) {
67 reducedIntegrationOrder=FALSE;
68 } else if (funcspace==DUDLEY_FACE_ELEMENTS) {
69 reducedIntegrationOrder=FALSE;
70 } else if (funcspace==DUDLEY_REDUCED_ELEMENTS) {
71 reducedIntegrationOrder=TRUE;
72 } else if (funcspace==DUDLEY_REDUCED_FACE_ELEMENTS) {
73 reducedIntegrationOrder=TRUE;
74 } else {
75 Dudley_setError(TYPE_ERROR,"Assemble_LumpedSystem: assemblage failed because of illegal function space.");
76 }
77 if (! Dudley_noError()) return;
78
79 /* set all parameters in p*/
80 Assemble_getAssembleParameters(nodes,elements,NULL,lumpedMat, reducedIntegrationOrder, &p);
81 if (! Dudley_noError()) return;
82
83 /* check if all function spaces are the same */
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 Dudley_setError(TYPE_ERROR,error_msg);
87 }
88
89 /* check the dimensions: */
90 if (p.numEqu==1) {
91 if (!isEmpty(D)) {
92 if (!isDataPointShapeEqual(D,0,dimensions)) {
93 Dudley_setError(TYPE_ERROR,"Assemble_LumpedSystem: coefficient D, rank 0 expected.");
94 }
95
96 }
97 } else {
98 if (!isEmpty(D)) {
99 dimensions[0]=p.numEqu;
100 if (!isDataPointShapeEqual(D,1,dimensions)) {
101 sprintf(error_msg,"Assemble_LumpedSystem: coefficient D, expected shape (%d,)",dimensions[0]);
102 Dudley_setError(TYPE_ERROR,error_msg);
103 }
104 }
105 }
106 if (Dudley_noError()) {
107 requireWrite(lumpedMat);
108 lumpedMat_p=getSampleDataRW(lumpedMat,0);
109 len_EM_lumpedMat=p.numShapes*p.numEqu;
110 len_EM_lumpedMat_size=len_EM_lumpedMat*sizeof(double);
111
112 expandedD=isExpanded(D);
113 if (!getQuadShape(elements->numDim, reducedIntegrationOrder, &S))
114 {
115 Dudley_setError(TYPE_ERROR, "Assemble_LumpedSystem: Unable to locate shape function.");
116 }
117
118 #ifdef NEW_LUMPING
119 #pragma omp parallel private(color, EM_lumpedMat, row_index, Vol, D_p, s, q, k, rtmp, diagS, m_t)
120 #else
121 #pragma omp parallel private(color, EM_lumpedMat, row_index, Vol, D_p, s, q, k, rtmp)
122 #endif
123 {
124 EM_lumpedMat=THREAD_MEMALLOC(len_EM_lumpedMat,double);
125 row_index=THREAD_MEMALLOC(p.numShapes,index_t);
126 if ( !Dudley_checkPtr(EM_lumpedMat) && !Dudley_checkPtr(row_index) ) {
127 if (p.numEqu == 1) { /* single equation */
128 if (expandedD) { /* with expanded D */
129 for (color=elements->minColor;color<=elements->maxColor;color++) {
130 /* open loop over all elements: */
131 #pragma omp for private(e) schedule(static)
132 for(e=0;e<elements->numElements;e++){
133 if (elements->Color[e]==color) {
134 double vol=p.row_jac->absD[e]*p.row_jac->quadweight;
135 D_p=getSampleDataRO(D,e);
136 #ifdef NEW_LUMPING /* HRZ lumping */
137 m_t=0; /* mass of the element: m_t */
138 for (q=0;q<p.numQuad;q++) m_t+=vol*D_p[INDEX2(q, 0,p.numQuad) ];
139 diagS=0; /* diagonal sum: S */
140 for (s=0;s<p.numShapes;s++) {
141 rtmp=0;
142 for (q=0;q<p.numQuad;q++) rtmp+=vol*D_p[INDEX2(q, 0, p.numQuad)]*S[INDEX2(s,q,p.numShapes)]*S[INDEX2(s,q,p.numShapes)];
143 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp;
144 diagS+=rtmp;
145 }
146 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
147 rtmp=m_t/diagS;
148 for (s=0;s<p.numShapes;s++) EM_lumpedMat[INDEX2(0,s,p.numEqu)]*=rtmp;
149
150 #else /* row-sum lumping */
151 for (s=0;s<p.numShapes;s++)
152 {
153 rtmp=0;
154 for (q=0;q<p.numQuad;q++) rtmp+=vol*S[INDEX2(s,q,p.numShapes)]*D_p[INDEX2(q, 0, p.numQuad)];
155 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp;
156 }
157 #endif
158 for (q=0;q<p.numShapes;q++)
159 {
160 row_index[q]=p.row_DOF[elements->Nodes[INDEX2(q,e,p.NN)]];
161 }
162 Dudley_Util_AddScatter(p.numShapes,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
163 } /* end color check */
164 } /* end element loop */
165 } /* end color loop */
166 } else { /* with constant D */
167
168 for (color=elements->minColor;color<=elements->maxColor;color++) {
169 /* open loop over all elements: */
170 #pragma omp for private(e) schedule(static)
171 for(e=0;e<elements->numElements;e++){
172 if (elements->Color[e]==color) {
173 double vol=p.row_jac->absD[e]*p.row_jac->quadweight;
174 D_p=getSampleDataRO(D,e);
175 #ifdef NEW_LUMPING /* HRZ lumping */
176 m_t=0; /* mass of the element: m_t */
177 for (q=0;q<p.numQuad;q++) m_t+=vol;
178 diagS=0; /* diagonal sum: S */
179 for (s=0;s<p.numShapes;s++) {
180 rtmp=0;
181 for (q=0;q<p.numQuad;q++)
182 {
183 rtmp+=vol*S[INDEX2(s,q,p.numShapes)]*S[INDEX2(s,q,p.numShapes)];
184 }
185 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp;
186 diagS+=rtmp;
187 }
188 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
189 rtmp=m_t/diagS*D_p[0];
190 for (s=0;s<p.numShapes;s++) EM_lumpedMat[INDEX2(0,s,p.numEqu)]*=rtmp;
191 #else /* row-sum lumping */
192 for (s=0;s<p.numShapes;s++) {
193 rtmp=0;
194 for (q=0;q<p.numQuad;q++) rtmp+=vol*S[INDEX2(s,q,p.numShapes)];
195 EM_lumpedMat[INDEX2(0,s,p.numEqu)]=rtmp*D_p[0];
196 }
197 #endif
198 for (q=0;q<p.numShapes;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(q,e,p.NN)]];
199 Dudley_Util_AddScatter(p.numShapes,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
200 } /* end color check */
201 } /* end element loop */
202 } /* end color loop */
203
204 }
205 } else { /* system of equation */
206 if (expandedD) { /* with expanded D */
207 for (color=elements->minColor;color<=elements->maxColor;color++) {
208 /* open loop over all elements: */
209 #pragma omp for private(e) schedule(static)
210 for(e=0;e<elements->numElements;e++){
211 if (elements->Color[e]==color) {
212 double vol=p.row_jac->absD[e]*p.row_jac->quadweight;
213 D_p=getSampleDataRO(D,e);
214
215 #ifdef NEW_LUMPING /* HRZ lumping */
216 for (k=0;k<p.numEqu;k++) {
217 m_t=0; /* mass of the element: m_t */
218 for (q=0;q<p.numQuad;q++) m_t+=vol*D_p[INDEX3(k,q,0,p.numEqu,p.numQuad)];
219
220 diagS=0; /* diagonal sum: S */
221 for (s=0;s<p.numShapes;s++) {
222 rtmp=0;
223 for (q=0;q<p.numQuad;q++) rtmp+=vol*D_p[INDEX3(k,q,0,p.numEqu,p.numQuad)]*S[INDEX2(s,q,p.numShapes)]*S[INDEX2(s,q,p.numShapes)];
224 EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp;
225 diagS+=rtmp;
226 }
227 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
228 rtmp=m_t/diagS;
229 for (s=0;s<p.numShapes;s++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]*=rtmp;
230 }
231 #else /* row-sum lumping */
232 for (s=0;s<p.numShapes;s++) {
233 for (k=0;k<p.numEqu;k++) {
234 rtmp=0.;
235 for (q=0;q<p.numQuad;q++) rtmp+=vol*S[INDEX2(s,q,p.numShapes)]*D_p[INDEX3(k,q,0,p.numEqu,p.numQuad)];
236 EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp;
237 }
238 }
239 #endif
240 for (q=0;q<p.numShapes;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(q,e,p.NN)]];
241 Dudley_Util_AddScatter(p.numShapes,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
242 } /* end color check */
243 } /* end element loop */
244 } /* end color loop */
245 } else { /* with constant D */
246 for (color=elements->minColor;color<=elements->maxColor;color++) {
247 /* open loop over all elements: */
248 #pragma omp for private(e) schedule(static)
249 for(e=0;e<elements->numElements;e++){
250 if (elements->Color[e]==color) {
251 double vol=p.row_jac->absD[e]*p.row_jac->quadweight; D_p=getSampleDataRO(D,e);
252
253 #ifdef NEW_LUMPING /* HRZ lumping */
254 m_t=0; /* mass of the element: m_t */
255 for (q=0;q<p.numQuad;q++) m_t+=vol;
256 diagS=0; /* diagonal sum: S */
257 for (s=0;s<p.numShapes;s++) {
258 rtmp=0;
259 for (q=0;q<p.numQuad;q++) rtmp+=vol*S[INDEX2(s,q,p.numShapes)]*S[INDEX2(s,q,p.numShapes)];
260 for (k=0;k<p.numEqu;k++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp;
261 diagS+=rtmp;
262 }
263 /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
264 rtmp=m_t/diagS;
265 for (s=0;s<p.numShapes;s++) {
266 for (k=0;k<p.numEqu;k++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]*=rtmp*D_p[k];
267 }
268 #else /* row-sum lumping */
269 for (s=0;s<p.numShapes;s++) {
270 for (k=0;k<p.numEqu;k++) {
271 rtmp=0.;
272 for (q=0;q<p.numQuad;q++) rtmp+=vol*S[INDEX2(s,q,p.numShapes)];
273 EM_lumpedMat[INDEX2(k,s,p.numEqu)]=rtmp*D_p[k];
274 }
275 }
276 #endif
277 for (q=0;q<p.numShapes;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(q,e,p.NN)]];
278 Dudley_Util_AddScatter(p.numShapes,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
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