/[escript]/branches/doubleplusgood/dudley/src/Assemble_PDE_Single2_3D.cpp
ViewVC logotype

Annotation of /branches/doubleplusgood/dudley/src/Assemble_PDE_Single2_3D.cpp

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.26