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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1645 - (hide annotations)
Tue Jul 15 05:26:38 2008 UTC (10 years, 8 months ago) by phornby
File MIME type: text/plain
File size: 14865 byte(s)
unused vars.
1 jgs 82
2 ksteube 1312 /* $Id$ */
3 jgs 150
4 ksteube 1312 /*******************************************************
5     *
6     * Copyright 2003-2007 by ACceSS MNRF
7     * Copyright 2007 by University of Queensland
8     *
9     * http://esscc.uq.edu.au
10     * Primary Business: Queensland, Australia
11     * Licensed under the Open Software License version 3.0
12     * http://www.opensource.org/licenses/osl-3.0.php
13     *
14     *******************************************************/
15    
16 jgs 82 /**************************************************************/
17    
18 gross 1204 /* assembles the mass matrix in lumped form */
19 jgs 82
20 gross 1204 /* The coefficient D has to be defined on the integration points or not present. */
21 jgs 82
22 gross 1204 /* lumpedMat has to be initialized before the routine is called. */
23 jgs 82
24     /**************************************************************/
25    
26     #include "Assemble.h"
27     #include "Util.h"
28     #ifdef _OPENMP
29     #include <omp.h>
30     #endif
31    
32    
33 ksteube 1221 /* Disabled until the tests pass */
34     /* #define NEW_LUMPING */ /* */
35 btully 1220
36 jgs 82 /**************************************************************/
37    
38 gross 1204 void Finley_Assemble_LumpedSystem(Finley_NodeFile* nodes,Finley_ElementFile* elements, escriptDataC* lumpedMat, escriptDataC* D)
39     {
40 jgs 82
41 gross 1204 bool_t reducedIntegrationOrder=FALSE, expandedD;
42 jgs 150 char error_msg[LenErrorMsg_MAX];
43 gross 798 Assemble_Parameters p;
44 gross 1204 dim_t dimensions[ESCRIPT_MAX_DATA_RANK], k, e, len_EM_lumpedMat, q, s;
45 gross 1028 type_t funcspace;
46 gross 1204 index_t color,*row_index=NULL;
47     double *S=NULL, *EM_lumpedMat=NULL, *Vol=NULL, *D_p=NULL, *lumpedMat_p=NULL;
48 phornby 1628 register double rtmp;
49 gross 1204 size_t len_EM_lumpedMat_size;
50 phornby 1645
51     #ifdef NEW_LUMPING
52 phornby 1628 register double m_t, diagS;
53 phornby 1645 #endif
54 gross 798
55 jgs 150 Finley_resetError();
56 jgs 82
57     if (nodes==NULL || elements==NULL) return;
58 gross 1204 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 gross 798 }
63 gross 1204 funcspace=getFunctionSpaceType(D);
64 jgs 82 /* check if all function spaces are the same */
65 gross 798 if (funcspace==FINLEY_ELEMENTS) {
66     reducedIntegrationOrder=FALSE;
67     } else if (funcspace==FINLEY_FACE_ELEMENTS) {
68     reducedIntegrationOrder=FALSE;
69 gross 1072 } else if (funcspace==FINLEY_REDUCED_ELEMENTS) {
70     reducedIntegrationOrder=TRUE;
71     } else if (funcspace==FINLEY_REDUCED_FACE_ELEMENTS) {
72     reducedIntegrationOrder=TRUE;
73 gross 798 } else {
74 gross 1204 Finley_setError(TYPE_ERROR,"Assemble_LumpedSystem: assemblage failed because of illegal function space.");
75 gross 798 }
76     if (! Finley_noError()) return;
77 jgs 82
78 gross 798 /* set all parameters in p*/
79 gross 1204 Assemble_getAssembleParameters(nodes,elements,NULL,lumpedMat, reducedIntegrationOrder, &p);
80 gross 798 if (! Finley_noError()) return;
81    
82     /* check if all function spaces are the same */
83    
84 jgs 82 if (! numSamplesEqual(D,p.numQuad,elements->numElements) ) {
85 gross 1204 sprintf(error_msg,"Assemble_LumpedSystem: sample points of coefficient D don't match (%d,%d)",p.numQuad,elements->numElements);
86 jgs 150 Finley_setError(TYPE_ERROR,error_msg);
87 jgs 82 }
88    
89     /* check the dimensions: */
90    
91 gross 1204 if (p.numEqu==1) {
92 jgs 82 if (!isEmpty(D)) {
93     if (!isDataPointShapeEqual(D,0,dimensions)) {
94 gross 1204 Finley_setError(TYPE_ERROR,"Assemble_LumpedSystem: coefficient D, rank 0 expected.");
95 jgs 82 }
96 gross 1204
97 jgs 82 }
98     } else {
99     if (!isEmpty(D)) {
100     dimensions[0]=p.numEqu;
101 gross 1204 if (!isDataPointShapeEqual(D,1,dimensions)) {
102     sprintf(error_msg,"Assemble_LumpedSystem: coefficient D, expected shape (%d,)",dimensions[0]);
103 jgs 150 Finley_setError(TYPE_ERROR,error_msg);
104 jgs 82 }
105     }
106     }
107 gross 1204
108 jgs 150 if (Finley_noError()) {
109 gross 1204 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 btully 1220 #pragma omp parallel private(color, EM_lumpedMat, row_index, Vol, D_p, s, q, k, rtmp, m_t, diagS)
116 gross 1204 {
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 btully 1220 #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 gross 1204 for (s=0;s<p.row_NS;s++) {
141     rtmp=0;
142 btully 1220 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 gross 1204 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 btully 1220 #endif
159 gross 1204 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 btully 1220 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 gross 1204 } /* end color check */
205     } /* end element loop */
206     } /* end color loop */
207     }
208 gross 798 } else {
209 gross 1204 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 btully 1220 #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 gross 1204 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 btully 1220 #endif
251 gross 1204 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 ksteube 1312 /* open loop over all elements: */
258 gross 1204 for (color=elements->minColor;color<=elements->maxColor;color++) {
259     #pragma omp for private(e) schedule(static)
260     for(e=0;e<elements->numElements;e++){
261     if (elements->Color[e]==color) {
262     Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
263     memset(EM_lumpedMat,0,len_EM_lumpedMat_size);
264     D_p=getSampleData(D,e);
265 btully 1220 #ifdef NEW_LUMPING /* HRZ lumping */
266     /*
267     * Number of PDEs: Multiple
268     * D_p varies over element: False
269     */
270     for (k=0;k<p.numEqu;k++) {
271     m_t=0; /* mass of the element: m_t */
272     for (q=0;q<p.numQuad;q++) {
273     m_t+=Vol[q]*D_p[k];
274     }
275     diagS=0; /* diagonal sum: S */
276     for (s=0;s<p.row_NS;s++) {
277     rtmp=0;
278     for (q=0;q<p.numQuad;q++) {
279     rtmp+=Vol[q]*D_p[k]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(s,q,p.row_NS)];
280     }
281     EM_lumpedMat[INDEX2(k,s,p.numEqu)]+=rtmp;
282     diagS+=EM_lumpedMat[INDEX2(k,s,p.numEqu)];
283     }
284     /* rescale diagonals by m_t/diagS to ensure consistent mass over element */
285     for (s=0;s<p.row_NS;s++) {
286     EM_lumpedMat[INDEX2(k,s,p.numEqu)]*=m_t/diagS;
287     }
288     }
289     #else /* row-sum lumping */
290 gross 1204 for (s=0;s<p.row_NS;s++) {
291     rtmp=0;
292     for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
293     for (k=0;k<p.numEqu;k++) EM_lumpedMat[INDEX2(k,s,p.numEqu)]+=rtmp*D_p[k];
294     }
295 btully 1220 #endif
296 gross 1204 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
297     Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_lumpedMat,lumpedMat_p, p.row_DOF_UpperBound);
298     } /* end color check */
299     } /* end element loop */
300     } /* end color loop */
301     }
302 gross 798 }
303 gross 1204 } /* end of pointer check */
304     THREAD_MEMFREE(EM_lumpedMat);
305     THREAD_MEMFREE(row_index);
306     } /* end parallel region */
307 jgs 82 }
308     }

Properties

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

  ViewVC Help
Powered by ViewVC 1.1.26