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

Contents of /trunk-mpi-branch/finley/src/Assemble_PDE_System2_3D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1253 - (show annotations)
Fri Aug 17 04:09:29 2007 UTC (12 years, 9 months ago) by gross
File MIME type: text/plain
File size: 21097 byte(s)
finally everthing for MPI is in the code. now testing is needed
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 for (color=elements->minColor;color<=elements->maxColor;color++) {
85 /* open loop over all elements: */
86 #pragma omp for private(e) schedule(static)
87 for(e=0;e<elements->numElements;e++){
88 if (elements->Color[e]==color) {
89 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
90 DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]);
91 for (q=0;q<len_EM_S;++q) EM_S[q]=0;
92 for (q=0;q<len_EM_F;++q) EM_F[q]=0;
93 add_EM_F=FALSE;
94 add_EM_S=FALSE;
95 /**************************************************************/
96 /* process A: */
97 /**************************************************************/
98 A_p=getSampleData(A,e);
99 if (NULL!=A_p) {
100 add_EM_S=TRUE;
101 if (extendedA) {
102 for (s=0;s<p.row_NS;s++) {
103 for (r=0;r<p.col_NS;r++) {
104 for (k=0;k<p.numEqu;k++) {
105 for (m=0;m<p.numComp;m++) {
106 rtmp=0;
107 for (q=0;q<p.numQuad;q++) {
108 rtmp+=Vol[q]*(
109 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)]
110 +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)]
111 +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)]
112 +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)]
113 +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)]
114 +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)]
115 +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)]
116 +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)]
117 +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)]);
118
119 }
120 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
121 }
122 }
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 rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
139 rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
140 rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
141 rtmp02+=rtmp0*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
142
143 rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
144 rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
145 rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
146 rtmp12+=rtmp1*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
147
148 rtmp2=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
149 rtmp20+=rtmp2*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
150 rtmp21+=rtmp2*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
151 rtmp22+=rtmp2*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
152 }
153 for (k=0;k<p.numEqu;k++) {
154 for (m=0;m<p.numComp;m++) {
155 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)]
156 +rtmp01*A_p[INDEX4(k,0,m,1,p.numEqu,DIM,p.numComp)]
157 +rtmp02*A_p[INDEX4(k,0,m,2,p.numEqu,DIM,p.numComp)]
158 +rtmp10*A_p[INDEX4(k,1,m,0,p.numEqu,DIM,p.numComp)]
159 +rtmp11*A_p[INDEX4(k,1,m,1,p.numEqu,DIM,p.numComp)]
160 +rtmp12*A_p[INDEX4(k,1,m,2,p.numEqu,DIM,p.numComp)]
161 +rtmp20*A_p[INDEX4(k,2,m,0,p.numEqu,DIM,p.numComp)]
162 +rtmp21*A_p[INDEX4(k,2,m,1,p.numEqu,DIM,p.numComp)]
163 +rtmp22*A_p[INDEX4(k,2,m,2,p.numEqu,DIM,p.numComp)];
164 }
165 }
166 }
167 }
168 }
169 }
170 /**************************************************************/
171 /* process B: */
172 /**************************************************************/
173 B_p=getSampleData(B,e);
174 if (NULL!=B_p) {
175 add_EM_S=TRUE;
176 if (extendedB) {
177 for (s=0;s<p.row_NS;s++) {
178 for (r=0;r<p.col_NS;r++) {
179 for (k=0;k<p.numEqu;k++) {
180 for (m=0;m<p.numComp;m++) {
181 rtmp=0;
182 for (q=0;q<p.numQuad;q++) {
183 rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]*(
184 DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX4(k,0,m,q,p.numEqu,DIM,p.numComp)]
185 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX4(k,1,m,q,p.numEqu,DIM,p.numComp)]
186 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*B_p[INDEX4(k,2,m,q,p.numEqu,DIM,p.numComp)] );
187 }
188 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
189 }
190 }
191 }
192 }
193 } else {
194 for (s=0;s<p.row_NS;s++) {
195 for (r=0;r<p.col_NS;r++) {
196 rtmp0=0;
197 rtmp1=0;
198 rtmp2=0;
199 for (q=0;q<p.numQuad;q++) {
200 rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)];
201 rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
202 rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
203 rtmp2+=rtmp*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
204 }
205 for (k=0;k<p.numEqu;k++) {
206 for (m=0;m<p.numComp;m++) {
207 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp0*B_p[INDEX3(k,0,m,p.numEqu,DIM)]
208 +rtmp1*B_p[INDEX3(k,1,m,p.numEqu,DIM)]
209 +rtmp2*B_p[INDEX3(k,2,m,p.numEqu,DIM)];
210 }
211 }
212 }
213 }
214 }
215 }
216 /**************************************************************/
217 /* process C: */
218 /**************************************************************/
219 C_p=getSampleData(C,e);
220 if (NULL!=C_p) {
221 add_EM_S=TRUE;
222 if (extendedC) {
223 for (s=0;s<p.row_NS;s++) {
224 for (r=0;r<p.col_NS;r++) {
225 for (k=0;k<p.numEqu;k++) {
226 for (m=0;m<p.numComp;m++) {
227 rtmp=0;
228 for (q=0;q<p.numQuad;q++) {
229 rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*(
230 C_p[INDEX4(k,m,0,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
231 +C_p[INDEX4(k,m,1,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
232 +C_p[INDEX4(k,m,2,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]);
233 }
234 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
235 }
236 }
237 }
238 }
239 } else {
240 for (s=0;s<p.row_NS;s++) {
241 for (r=0;r<p.col_NS;r++) {
242 rtmp0=0;
243 rtmp1=0;
244 rtmp2=0;
245 for (q=0;q<p.numQuad;q++) {
246 rtmp=Vol[q]*S[INDEX2(s,q,p.row_NS)];
247 rtmp0+=rtmp*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
248 rtmp1+=rtmp*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
249 rtmp2+=rtmp*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
250 }
251 for (k=0;k<p.numEqu;k++) {
252 for (m=0;m<p.numComp;m++) {
253 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)]
254 +rtmp1*C_p[INDEX3(k,m,1,p.numEqu,p.numComp)]
255 +rtmp2*C_p[INDEX3(k,m,2,p.numEqu,p.numComp)];
256 }
257 }
258 }
259 }
260 }
261 }
262 /************************************************************* */
263 /* process D */
264 /**************************************************************/
265 D_p=getSampleData(D,e);
266 if (NULL!=D_p) {
267 add_EM_S=TRUE;
268 if (extendedD) {
269 for (s=0;s<p.row_NS;s++) {
270 for (r=0;r<p.col_NS;r++) {
271 for (k=0;k<p.numEqu;k++) {
272 for (m=0;m<p.numComp;m++) {
273 rtmp=0;
274 for (q=0;q<p.numQuad;q++) {
275 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)];
276 }
277 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
278 }
279 }
280 }
281 }
282 } else {
283 for (s=0;s<p.row_NS;s++) {
284 for (r=0;r<p.col_NS;r++) {
285 rtmp=0;
286 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
287 for (k=0;k<p.numEqu;k++) {
288 for (m=0;m<p.numComp;m++) {
289 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[INDEX2(k,m,p.numEqu)];
290 }
291 }
292 }
293 }
294 }
295 }
296 /**************************************************************/
297 /* process X: */
298 /**************************************************************/
299 X_p=getSampleData(X,e);
300 if (NULL!=X_p) {
301 add_EM_F=TRUE;
302 if (extendedX) {
303 for (s=0;s<p.row_NS;s++) {
304 for (k=0;k<p.numEqu;k++) {
305 rtmp=0;
306 for (q=0;q<p.numQuad;q++) {
307 rtmp+=Vol[q]* ( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX3(k,0,q,p.numEqu,DIM)]
308 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX3(k,1,q,p.numEqu,DIM)]
309 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*X_p[INDEX3(k,2,q,p.numEqu,DIM)]);
310 }
311 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
312 }
313 }
314 } else {
315 for (s=0;s<p.row_NS;s++) {
316 rtmp0=0;
317 rtmp1=0;
318 rtmp2=0;
319 for (q=0;q<p.numQuad;q++) {
320 rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
321 rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
322 rtmp2+=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
323 }
324 for (k=0;k<p.numEqu;k++) {
325 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp0*X_p[INDEX2(k,0,p.numEqu)]
326 +rtmp1*X_p[INDEX2(k,1,p.numEqu)]
327 +rtmp2*X_p[INDEX2(k,2,p.numEqu)];
328 }
329 }
330 }
331 }
332 /**************************************************************/
333 /* process Y: */
334 /**************************************************************/
335 Y_p=getSampleData(Y,e);
336 if (NULL!=Y_p) {
337 add_EM_F=TRUE;
338 if (extendedY) {
339 for (s=0;s<p.row_NS;s++) {
340 for (k=0;k<p.numEqu;k++) {
341 rtmp=0.;
342 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[INDEX2(k,q,p.numEqu)];
343 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
344 }
345 }
346 } else {
347 for (s=0;s<p.row_NS;s++) {
348 rtmp=0;
349 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
350 for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*Y_p[k];
351 }
352 }
353 }
354 /***********************************************************************************************/
355 /* add the element matrices onto the matrix and right hand side */
356 /***********************************************************************************************/
357 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
358 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
359 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
360
361 } /* end color check */
362 } /* end element loop */
363 } /* end color loop */
364
365 THREAD_MEMFREE(EM_S);
366 THREAD_MEMFREE(EM_F);
367 THREAD_MEMFREE(row_index);
368
369 } /* end of pointer check */
370 } /* end parallel region */
371 }
372 /*
373 * $Log$
374 */

  ViewVC Help
Powered by ViewVC 1.1.26