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

Annotation of /branches/domexper/dudley/src/Assemble_PDE_System2_3D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1388 - (hide annotations)
Fri Jan 11 07:45:58 2008 UTC (11 years, 5 months ago) by trankine
Original Path: trunk/finley/src/Assemble_PDE_System2_3D.c
File MIME type: text/plain
File size: 20848 byte(s)
And get the *(&(*&(* name right
1 gross 798
2 ksteube 1312 /* $Id$ */
3    
4     /*******************************************************
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 gross 798 /**************************************************************/
17    
18     /* assembles the system of numEq PDEs into the stiffness matrix S right hand side F */
19     /* the shape functions for test and solution must be identical */
20    
21    
22     /* -(A_{k,i,m,j} u_m,j)_i-(B_{k,i,m} u_m)_i+C_{k,m,j} u_m,j-D_{k,m} u_m and -(X_{k,i})_i + Y_k */
23    
24     /* u has p.numComp components in a 3D domain. The shape functions for test and solution must be identical */
25     /* and row_NS == row_NN */
26    
27     /* Shape of the coefficients: */
28    
29     /* A = p.numEqu x 3 x p.numComp x 3 */
30     /* B = 3 x p.numEqu x p.numComp */
31     /* C = p.numEqu x 3 x p.numComp */
32     /* D = p.numEqu x p.numComp */
33     /* X = p.numEqu x 3 */
34     /* Y = p.numEqu */
35    
36    
37     /**************************************************************/
38    
39    
40     #include "Assemble.h"
41     #include "Util.h"
42 gross 853 #ifdef _OPENMP
43     #include <omp.h>
44     #endif
45 gross 798
46 gross 853
47 gross 798 /**************************************************************/
48    
49     void Finley_Assemble_PDE_System2_3D(Assemble_Parameters p, Finley_ElementFile* elements,
50     Paso_SystemMatrix* Mat, escriptDataC* F,
51     escriptDataC* A, escriptDataC* B, escriptDataC* C, escriptDataC* D, escriptDataC* X, escriptDataC* Y) {
52    
53     #define DIM 3
54     index_t color;
55     dim_t e;
56 gross 853 double *EM_S, *EM_F, *Vol, *DSDX, *A_p, *B_p, *C_p, *D_p, *X_p, *Y_p;
57     index_t *row_index;
58     register dim_t q, s,r,k,m;
59     register double rtmp, rtmp0, rtmp1, rtmp2, rtmp00, rtmp01, rtmp02, rtmp10, rtmp11, rtmp12, rtmp20, rtmp21, rtmp22;
60     bool_t add_EM_F, add_EM_S;
61    
62 gross 798 bool_t extendedA=isExpanded(A);
63     bool_t extendedB=isExpanded(B);
64     bool_t extendedC=isExpanded(C);
65     bool_t extendedD=isExpanded(D);
66     bool_t extendedX=isExpanded(X);
67     bool_t extendedY=isExpanded(Y);
68     double *F_p=getSampleData(F,0);
69     double *S=p.row_jac->ReferenceElement->S;
70     dim_t len_EM_S=p.row_NN*p.col_NN*p.numEqu*p.numComp;
71     dim_t len_EM_F=p.row_NN*p.numEqu;
72    
73    
74 gross 853 #pragma omp parallel private(color,EM_S, EM_F, Vol, DSDX, A_p, B_p, C_p, D_p, X_p, Y_p,row_index,q, s,r,k,m,rtmp, rtmp0, rtmp1, rtmp2, rtmp00, rtmp01, rtmp02, rtmp10, rtmp11, rtmp12, rtmp20, rtmp21, rtmp22,add_EM_F, add_EM_S)
75 gross 798 {
76 gross 853 EM_S=THREAD_MEMALLOC(len_EM_S,double);
77     EM_F=THREAD_MEMALLOC(len_EM_F,double);
78     row_index=THREAD_MEMALLOC(p.row_NN,index_t);
79    
80     if (!Finley_checkPtr(EM_S) && !Finley_checkPtr(EM_F) && !Finley_checkPtr(row_index) ) {
81 gross 798
82 gross 853 for (color=elements->minColor;color<=elements->maxColor;color++) {
83     /* open loop over all elements: */
84     #pragma omp for private(e) schedule(static)
85     for(e=0;e<elements->numElements;e++){
86     if (elements->Color[e]==color) {
87     Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
88     DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]);
89     for (q=0;q<len_EM_S;++q) EM_S[q]=0;
90     for (q=0;q<len_EM_F;++q) EM_F[q]=0;
91     add_EM_F=FALSE;
92     add_EM_S=FALSE;
93     /**************************************************************/
94     /* process A: */
95     /**************************************************************/
96     A_p=getSampleData(A,e);
97     if (NULL!=A_p) {
98     add_EM_S=TRUE;
99     if (extendedA) {
100     for (s=0;s<p.row_NS;s++) {
101     for (r=0;r<p.col_NS;r++) {
102     for (k=0;k<p.numEqu;k++) {
103 gross 798 for (m=0;m<p.numComp;m++) {
104 gross 853 rtmp=0;
105     for (q=0;q<p.numQuad;q++) {
106     rtmp+=Vol[q]*(
107     DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*A_p[INDEX5(k,0,m,0,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
108     +DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*A_p[INDEX5(k,0,m,1,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
109     +DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*A_p[INDEX5(k,0,m,2,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]
110     +DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*A_p[INDEX5(k,1,m,0,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
111     +DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*A_p[INDEX5(k,1,m,1,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
112     +DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*A_p[INDEX5(k,1,m,2,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]
113     +DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*A_p[INDEX5(k,2,m,0,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
114     +DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*A_p[INDEX5(k,2,m,1,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
115     +DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*A_p[INDEX5(k,2,m,2,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]);
116    
117     }
118     EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
119 gross 798 }
120 gross 853 }
121     }
122     }
123     } else {
124     for (s=0;s<p.row_NS;s++) {
125     for (r=0;r<p.col_NS;r++) {
126     rtmp00=0;
127     rtmp01=0;
128     rtmp02=0;
129     rtmp10=0;
130     rtmp11=0;
131     rtmp12=0;
132     rtmp20=0;
133     rtmp21=0;
134     rtmp22=0;
135 gross 798 for (q=0;q<p.numQuad;q++) {
136 gross 853 rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
137     rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
138     rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
139     rtmp02+=rtmp0*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
140    
141     rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
142     rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
143     rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
144     rtmp12+=rtmp1*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
145    
146     rtmp2=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
147     rtmp20+=rtmp2*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
148     rtmp21+=rtmp2*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
149     rtmp22+=rtmp2*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
150 gross 798 }
151     for (k=0;k<p.numEqu;k++) {
152     for (m=0;m<p.numComp;m++) {
153 gross 853 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+= rtmp00*A_p[INDEX4(k,0,m,0,p.numEqu,DIM,p.numComp)]
154     +rtmp01*A_p[INDEX4(k,0,m,1,p.numEqu,DIM,p.numComp)]
155     +rtmp02*A_p[INDEX4(k,0,m,2,p.numEqu,DIM,p.numComp)]
156     +rtmp10*A_p[INDEX4(k,1,m,0,p.numEqu,DIM,p.numComp)]
157     +rtmp11*A_p[INDEX4(k,1,m,1,p.numEqu,DIM,p.numComp)]
158     +rtmp12*A_p[INDEX4(k,1,m,2,p.numEqu,DIM,p.numComp)]
159     +rtmp20*A_p[INDEX4(k,2,m,0,p.numEqu,DIM,p.numComp)]
160     +rtmp21*A_p[INDEX4(k,2,m,1,p.numEqu,DIM,p.numComp)]
161     +rtmp22*A_p[INDEX4(k,2,m,2,p.numEqu,DIM,p.numComp)];
162 gross 798 }
163     }
164 gross 853 }
165     }
166     }
167     }
168     /**************************************************************/
169     /* process B: */
170     /**************************************************************/
171     B_p=getSampleData(B,e);
172     if (NULL!=B_p) {
173     add_EM_S=TRUE;
174     if (extendedB) {
175     for (s=0;s<p.row_NS;s++) {
176     for (r=0;r<p.col_NS;r++) {
177     for (k=0;k<p.numEqu;k++) {
178     for (m=0;m<p.numComp;m++) {
179     rtmp=0;
180     for (q=0;q<p.numQuad;q++) {
181     rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]*(
182     DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX4(k,0,m,q,p.numEqu,DIM,p.numComp)]
183     + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX4(k,1,m,q,p.numEqu,DIM,p.numComp)]
184     + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*B_p[INDEX4(k,2,m,q,p.numEqu,DIM,p.numComp)] );
185     }
186     EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
187     }
188     }
189     }
190     }
191     } else {
192     for (s=0;s<p.row_NS;s++) {
193     for (r=0;r<p.col_NS;r++) {
194     rtmp0=0;
195     rtmp1=0;
196     rtmp2=0;
197     for (q=0;q<p.numQuad;q++) {
198     rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)];
199     rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
200     rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
201     rtmp2+=rtmp*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
202     }
203     for (k=0;k<p.numEqu;k++) {
204     for (m=0;m<p.numComp;m++) {
205     EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp0*B_p[INDEX3(k,0,m,p.numEqu,DIM)]
206     +rtmp1*B_p[INDEX3(k,1,m,p.numEqu,DIM)]
207     +rtmp2*B_p[INDEX3(k,2,m,p.numEqu,DIM)];
208     }
209     }
210     }
211     }
212 gross 798 }
213     }
214 gross 853 /**************************************************************/
215     /* process C: */
216     /**************************************************************/
217     C_p=getSampleData(C,e);
218     if (NULL!=C_p) {
219     add_EM_S=TRUE;
220     if (extendedC) {
221     for (s=0;s<p.row_NS;s++) {
222     for (r=0;r<p.col_NS;r++) {
223     for (k=0;k<p.numEqu;k++) {
224     for (m=0;m<p.numComp;m++) {
225     rtmp=0;
226     for (q=0;q<p.numQuad;q++) {
227     rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*(
228     C_p[INDEX4(k,m,0,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
229     +C_p[INDEX4(k,m,1,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
230     +C_p[INDEX4(k,m,2,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]);
231     }
232     EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
233 gross 798 }
234 gross 853 }
235     }
236     }
237     } else {
238     for (s=0;s<p.row_NS;s++) {
239     for (r=0;r<p.col_NS;r++) {
240     rtmp0=0;
241     rtmp1=0;
242     rtmp2=0;
243     for (q=0;q<p.numQuad;q++) {
244     rtmp=Vol[q]*S[INDEX2(s,q,p.row_NS)];
245     rtmp0+=rtmp*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
246     rtmp1+=rtmp*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
247     rtmp2+=rtmp*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
248     }
249     for (k=0;k<p.numEqu;k++) {
250     for (m=0;m<p.numComp;m++) {
251     EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp0*C_p[INDEX3(k,m,0,p.numEqu,p.numComp)]
252     +rtmp1*C_p[INDEX3(k,m,1,p.numEqu,p.numComp)]
253     +rtmp2*C_p[INDEX3(k,m,2,p.numEqu,p.numComp)];
254     }
255     }
256     }
257     }
258     }
259     }
260     /************************************************************* */
261     /* process D */
262     /**************************************************************/
263     D_p=getSampleData(D,e);
264     if (NULL!=D_p) {
265     add_EM_S=TRUE;
266     if (extendedD) {
267     for (s=0;s<p.row_NS;s++) {
268     for (r=0;r<p.col_NS;r++) {
269     for (k=0;k<p.numEqu;k++) {
270     for (m=0;m<p.numComp;m++) {
271     rtmp=0;
272     for (q=0;q<p.numQuad;q++) {
273     rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*D_p[INDEX3(k,m,q,p.numEqu,p.numComp)]*S[INDEX2(r,q,p.row_NS)];
274     }
275     EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
276     }
277     }
278     }
279     }
280     } else {
281     for (s=0;s<p.row_NS;s++) {
282     for (r=0;r<p.col_NS;r++) {
283     rtmp=0;
284     for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
285 gross 798 for (k=0;k<p.numEqu;k++) {
286 gross 853 for (m=0;m<p.numComp;m++) {
287     EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[INDEX2(k,m,p.numEqu)];
288     }
289 gross 798 }
290 gross 853 }
291     }
292     }
293     }
294     /**************************************************************/
295     /* process X: */
296     /**************************************************************/
297     X_p=getSampleData(X,e);
298     if (NULL!=X_p) {
299     add_EM_F=TRUE;
300     if (extendedX) {
301     for (s=0;s<p.row_NS;s++) {
302     for (k=0;k<p.numEqu;k++) {
303 gross 798 rtmp=0;
304     for (q=0;q<p.numQuad;q++) {
305 gross 853 rtmp+=Vol[q]* ( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX3(k,0,q,p.numEqu,DIM)]
306     + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX3(k,1,q,p.numEqu,DIM)]
307     + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*X_p[INDEX3(k,2,q,p.numEqu,DIM)]);
308 gross 798 }
309 gross 853 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
310     }
311     }
312     } else {
313     for (s=0;s<p.row_NS;s++) {
314     rtmp0=0;
315     rtmp1=0;
316     rtmp2=0;
317     for (q=0;q<p.numQuad;q++) {
318     rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
319     rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
320     rtmp2+=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
321 gross 798 }
322 gross 853 for (k=0;k<p.numEqu;k++) {
323     EM_F[INDEX2(k,s,p.numEqu)]+=rtmp0*X_p[INDEX2(k,0,p.numEqu)]
324     +rtmp1*X_p[INDEX2(k,1,p.numEqu)]
325     +rtmp2*X_p[INDEX2(k,2,p.numEqu)];
326     }
327 gross 798 }
328 gross 853 }
329 gross 798 }
330 gross 853 /**************************************************************/
331     /* process Y: */
332     /**************************************************************/
333     Y_p=getSampleData(Y,e);
334     if (NULL!=Y_p) {
335     add_EM_F=TRUE;
336     if (extendedY) {
337     for (s=0;s<p.row_NS;s++) {
338     for (k=0;k<p.numEqu;k++) {
339     rtmp=0.;
340     for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[INDEX2(k,q,p.numEqu)];
341     EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
342 gross 798 }
343     }
344 gross 853 } else {
345     for (s=0;s<p.row_NS;s++) {
346     rtmp=0;
347     for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
348     for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*Y_p[k];
349 gross 798 }
350 gross 853 }
351     }
352     /***********************************************************************************************/
353     /* add the element matrices onto the matrix and right hand side */
354     /***********************************************************************************************/
355     for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
356     if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
357     if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
358    
359     } /* end color check */
360     } /* end element loop */
361     } /* end color loop */
362    
363     THREAD_MEMFREE(EM_S);
364     THREAD_MEMFREE(EM_F);
365     THREAD_MEMFREE(row_index);
366 gross 798
367 gross 853 } /* end of pointer check */
368 gross 798 } /* end parallel region */
369     }
370     /*
371     * $Log$
372     */

  ViewVC Help
Powered by ViewVC 1.1.26