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

Contents of /trunk/finley/src/Assemble_PDE_System2_3D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 853 - (show annotations)
Wed Sep 20 05:56:36 2006 UTC (13 years, 1 month ago) by gross
File MIME type: text/plain
File size: 21239 byte(s)
some performance problems wit openmp fixed.
1 /*
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_{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 */
20
21 /* u has p.numComp components 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 = p.numEqu x 3 x p.numComp x 3 */
27 /* B = 3 x p.numEqu x p.numComp */
28 /* C = p.numEqu x 3 x p.numComp */
29 /* D = p.numEqu x p.numComp */
30 /* X = p.numEqu x 3 */
31 /* Y = p.numEqu */
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 #ifdef _OPENMP
45 #include <omp.h>
46 #endif
47
48
49 /**************************************************************/
50
51 void Finley_Assemble_PDE_System2_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 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,k,m;
61 register double rtmp, rtmp0, rtmp1, rtmp2, rtmp00, rtmp01, rtmp02, rtmp10, rtmp11, rtmp12, rtmp20, rtmp21, rtmp22;
62 bool_t add_EM_F, add_EM_S;
63
64 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*p.numEqu*p.numComp;
73 dim_t len_EM_F=p.row_NN*p.numEqu;
74
75
76 #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)
77 {
78 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
84 #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 for (k=0;k<p.numEqu;k++) {
111 for (m=0;m<p.numComp;m++) {
112 rtmp=0;
113 for (q=0;q<p.numQuad;q++) {
114 rtmp+=Vol[q]*(
115 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)]
116 +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)]
117 +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)]
118 +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)]
119 +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)]
120 +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)]
121 +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)]
122 +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)]
123 +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)]);
124
125 }
126 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
127 }
128 }
129 }
130 }
131 } else {
132 for (s=0;s<p.row_NS;s++) {
133 for (r=0;r<p.col_NS;r++) {
134 rtmp00=0;
135 rtmp01=0;
136 rtmp02=0;
137 rtmp10=0;
138 rtmp11=0;
139 rtmp12=0;
140 rtmp20=0;
141 rtmp21=0;
142 rtmp22=0;
143 for (q=0;q<p.numQuad;q++) {
144 rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
145 rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
146 rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
147 rtmp02+=rtmp0*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
148
149 rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
150 rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
151 rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
152 rtmp12+=rtmp1*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
153
154 rtmp2=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
155 rtmp20+=rtmp2*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
156 rtmp21+=rtmp2*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
157 rtmp22+=rtmp2*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
158 }
159 for (k=0;k<p.numEqu;k++) {
160 for (m=0;m<p.numComp;m++) {
161 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)]
162 +rtmp01*A_p[INDEX4(k,0,m,1,p.numEqu,DIM,p.numComp)]
163 +rtmp02*A_p[INDEX4(k,0,m,2,p.numEqu,DIM,p.numComp)]
164 +rtmp10*A_p[INDEX4(k,1,m,0,p.numEqu,DIM,p.numComp)]
165 +rtmp11*A_p[INDEX4(k,1,m,1,p.numEqu,DIM,p.numComp)]
166 +rtmp12*A_p[INDEX4(k,1,m,2,p.numEqu,DIM,p.numComp)]
167 +rtmp20*A_p[INDEX4(k,2,m,0,p.numEqu,DIM,p.numComp)]
168 +rtmp21*A_p[INDEX4(k,2,m,1,p.numEqu,DIM,p.numComp)]
169 +rtmp22*A_p[INDEX4(k,2,m,2,p.numEqu,DIM,p.numComp)];
170 }
171 }
172 }
173 }
174 }
175 }
176 /**************************************************************/
177 /* process B: */
178 /**************************************************************/
179 B_p=getSampleData(B,e);
180 if (NULL!=B_p) {
181 add_EM_S=TRUE;
182 if (extendedB) {
183 for (s=0;s<p.row_NS;s++) {
184 for (r=0;r<p.col_NS;r++) {
185 for (k=0;k<p.numEqu;k++) {
186 for (m=0;m<p.numComp;m++) {
187 rtmp=0;
188 for (q=0;q<p.numQuad;q++) {
189 rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]*(
190 DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX4(k,0,m,q,p.numEqu,DIM,p.numComp)]
191 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX4(k,1,m,q,p.numEqu,DIM,p.numComp)]
192 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*B_p[INDEX4(k,2,m,q,p.numEqu,DIM,p.numComp)] );
193 }
194 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
195 }
196 }
197 }
198 }
199 } else {
200 for (s=0;s<p.row_NS;s++) {
201 for (r=0;r<p.col_NS;r++) {
202 rtmp0=0;
203 rtmp1=0;
204 rtmp2=0;
205 for (q=0;q<p.numQuad;q++) {
206 rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)];
207 rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
208 rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
209 rtmp2+=rtmp*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
210 }
211 for (k=0;k<p.numEqu;k++) {
212 for (m=0;m<p.numComp;m++) {
213 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp0*B_p[INDEX3(k,0,m,p.numEqu,DIM)]
214 +rtmp1*B_p[INDEX3(k,1,m,p.numEqu,DIM)]
215 +rtmp2*B_p[INDEX3(k,2,m,p.numEqu,DIM)];
216 }
217 }
218 }
219 }
220 }
221 }
222 /**************************************************************/
223 /* process C: */
224 /**************************************************************/
225 C_p=getSampleData(C,e);
226 if (NULL!=C_p) {
227 add_EM_S=TRUE;
228 if (extendedC) {
229 for (s=0;s<p.row_NS;s++) {
230 for (r=0;r<p.col_NS;r++) {
231 for (k=0;k<p.numEqu;k++) {
232 for (m=0;m<p.numComp;m++) {
233 rtmp=0;
234 for (q=0;q<p.numQuad;q++) {
235 rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*(
236 C_p[INDEX4(k,m,0,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
237 +C_p[INDEX4(k,m,1,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
238 +C_p[INDEX4(k,m,2,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]);
239 }
240 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
241 }
242 }
243 }
244 }
245 } else {
246 for (s=0;s<p.row_NS;s++) {
247 for (r=0;r<p.col_NS;r++) {
248 rtmp0=0;
249 rtmp1=0;
250 rtmp2=0;
251 for (q=0;q<p.numQuad;q++) {
252 rtmp=Vol[q]*S[INDEX2(s,q,p.row_NS)];
253 rtmp0+=rtmp*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
254 rtmp1+=rtmp*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
255 rtmp2+=rtmp*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
256 }
257 for (k=0;k<p.numEqu;k++) {
258 for (m=0;m<p.numComp;m++) {
259 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)]
260 +rtmp1*C_p[INDEX3(k,m,1,p.numEqu,p.numComp)]
261 +rtmp2*C_p[INDEX3(k,m,2,p.numEqu,p.numComp)];
262 }
263 }
264 }
265 }
266 }
267 }
268 /************************************************************* */
269 /* process D */
270 /**************************************************************/
271 D_p=getSampleData(D,e);
272 if (NULL!=D_p) {
273 add_EM_S=TRUE;
274 if (extendedD) {
275 for (s=0;s<p.row_NS;s++) {
276 for (r=0;r<p.col_NS;r++) {
277 for (k=0;k<p.numEqu;k++) {
278 for (m=0;m<p.numComp;m++) {
279 rtmp=0;
280 for (q=0;q<p.numQuad;q++) {
281 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)];
282 }
283 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
284 }
285 }
286 }
287 }
288 } else {
289 for (s=0;s<p.row_NS;s++) {
290 for (r=0;r<p.col_NS;r++) {
291 rtmp=0;
292 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
293 for (k=0;k<p.numEqu;k++) {
294 for (m=0;m<p.numComp;m++) {
295 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[INDEX2(k,m,p.numEqu)];
296 }
297 }
298 }
299 }
300 }
301 }
302 /**************************************************************/
303 /* process X: */
304 /**************************************************************/
305 X_p=getSampleData(X,e);
306 if (NULL!=X_p) {
307 add_EM_F=TRUE;
308 if (extendedX) {
309 for (s=0;s<p.row_NS;s++) {
310 for (k=0;k<p.numEqu;k++) {
311 rtmp=0;
312 for (q=0;q<p.numQuad;q++) {
313 rtmp+=Vol[q]* ( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX3(k,0,q,p.numEqu,DIM)]
314 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX3(k,1,q,p.numEqu,DIM)]
315 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*X_p[INDEX3(k,2,q,p.numEqu,DIM)]);
316 }
317 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
318 }
319 }
320 } else {
321 for (s=0;s<p.row_NS;s++) {
322 rtmp0=0;
323 rtmp1=0;
324 rtmp2=0;
325 for (q=0;q<p.numQuad;q++) {
326 rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
327 rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
328 rtmp2+=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
329 }
330 for (k=0;k<p.numEqu;k++) {
331 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp0*X_p[INDEX2(k,0,p.numEqu)]
332 +rtmp1*X_p[INDEX2(k,1,p.numEqu)]
333 +rtmp2*X_p[INDEX2(k,2,p.numEqu)];
334 }
335 }
336 }
337 }
338 /**************************************************************/
339 /* process Y: */
340 /**************************************************************/
341 Y_p=getSampleData(Y,e);
342 if (NULL!=Y_p) {
343 add_EM_F=TRUE;
344 if (extendedY) {
345 for (s=0;s<p.row_NS;s++) {
346 for (k=0;k<p.numEqu;k++) {
347 rtmp=0.;
348 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[INDEX2(k,q,p.numEqu)];
349 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
350 }
351 }
352 } else {
353 for (s=0;s<p.row_NS;s++) {
354 rtmp=0;
355 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
356 for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*Y_p[k];
357 }
358 }
359 }
360 /***********************************************************************************************/
361 /* add the element matrices onto the matrix and right hand side */
362 /***********************************************************************************************/
363 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
364 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
365 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
366
367 } /* end color check */
368 } /* end element loop */
369 } /* end color loop */
370
371 THREAD_MEMFREE(EM_S);
372 THREAD_MEMFREE(EM_F);
373 THREAD_MEMFREE(row_index);
374
375 } /* end of pointer check */
376 } /* end parallel region */
377 }
378 /*
379 * $Log$
380 */

  ViewVC Help
Powered by ViewVC 1.1.26