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

Contents of /trunk-mpi-branch/finley/src/Assemble_PDE_Single2_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: 17548 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_{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 #ifdef _OPENMP
45 #include <omp.h>
46 #endif
47
48
49 /**************************************************************/
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 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 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 #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 {
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 rtmp=0;
105 for (q=0;q<p.numQuad;q++) {
106 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)]
107 + 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)]
108 + 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)]
109 + 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)]
110 + 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)]
111 + 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)]
112 + 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)]
113 + 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)]
114 + 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)]);
115 }
116 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
117 }
118 }
119 } else {
120 for (s=0;s<p.row_NS;s++) {
121 for (r=0;r<p.col_NS;r++) {
122 rtmp00=0;
123 rtmp01=0;
124 rtmp02=0;
125 rtmp10=0;
126 rtmp11=0;
127 rtmp12=0;
128 rtmp20=0;
129 rtmp21=0;
130 rtmp22=0;
131 for (q=0;q<p.numQuad;q++) {
132
133 rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
134 rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
135 rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
136 rtmp02+=rtmp0*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
137
138 rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
139 rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
140 rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
141 rtmp12+=rtmp1*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
142
143 rtmp2=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
144 rtmp20+=rtmp2*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
145 rtmp21+=rtmp2*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
146 rtmp22+=rtmp2*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
147 }
148 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp00*A_p[INDEX2(0,0,DIM)]
149 +rtmp01*A_p[INDEX2(0,1,DIM)]
150 +rtmp02*A_p[INDEX2(0,2,DIM)]
151 +rtmp10*A_p[INDEX2(1,0,DIM)]
152 +rtmp11*A_p[INDEX2(1,1,DIM)]
153 +rtmp12*A_p[INDEX2(1,2,DIM)]
154 +rtmp20*A_p[INDEX2(2,0,DIM)]
155 +rtmp21*A_p[INDEX2(2,1,DIM)]
156 +rtmp22*A_p[INDEX2(2,2,DIM)];
157 }
158 }
159 }
160 }
161 /**************************************************************/
162 /* process B: */
163 /**************************************************************/
164 B_p=getSampleData(B,e);
165 if (NULL!=B_p) {
166 add_EM_S=TRUE;
167 if (extendedB) {
168 for (s=0;s<p.row_NS;s++) {
169 for (r=0;r<p.col_NS;r++) {
170 rtmp=0;
171 for (q=0;q<p.numQuad;q++) {
172 rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]*
173 ( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX2(0,q,DIM)]
174 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX2(1,q,DIM)]
175 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*B_p[INDEX2(2,q,DIM)]);
176 }
177 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
178 }
179 }
180 } else {
181 for (s=0;s<p.row_NS;s++) {
182 for (r=0;r<p.col_NS;r++) {
183 rtmp0=0;
184 rtmp1=0;
185 rtmp2=0;
186 for (q=0;q<p.numQuad;q++) {
187 rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)];
188 rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
189 rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
190 rtmp2+=rtmp*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
191 }
192 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];
193 }
194 }
195 }
196 }
197 /**************************************************************/
198 /* process C: */
199 /**************************************************************/
200 C_p=getSampleData(C,e);
201 if (NULL!=C_p) {
202 add_EM_S=TRUE;
203 if (extendedC) {
204 for (s=0;s<p.row_NS;s++) {
205 for (r=0;r<p.col_NS;r++) {
206 rtmp=0;
207 for (q=0;q<p.numQuad;q++) {
208 rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*
209 ( C_p[INDEX2(0,q,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
210 + C_p[INDEX2(1,q,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
211 + C_p[INDEX2(2,q,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]);
212 }
213 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
214 }
215 }
216 } else {
217 for (s=0;s<p.row_NS;s++) {
218 for (r=0;r<p.col_NS;r++) {
219 rtmp0=0;
220 rtmp1=0;
221 rtmp2=0;
222 for (q=0;q<p.numQuad;q++) {
223 rtmp=Vol[q]*S[INDEX2(s,q,p.row_NS)];
224 rtmp0+=rtmp*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
225 rtmp1+=rtmp*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
226 rtmp2+=rtmp*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
227 }
228 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];
229 }
230 }
231 }
232 }
233 /************************************************************* */
234 /* process D */
235 /**************************************************************/
236 D_p=getSampleData(D,e);
237 if (NULL!=D_p) {
238 add_EM_S=TRUE;
239 if (extendedD) {
240 for (s=0;s<p.row_NS;s++) {
241 for (r=0;r<p.col_NS;r++) {
242 rtmp=0;
243 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)];
244 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
245 }
246 }
247 } else {
248 for (s=0;s<p.row_NS;s++) {
249 for (r=0;r<p.col_NS;r++) {
250 rtmp=0;
251 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
252 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[0];
253 }
254 }
255 }
256 }
257 /**************************************************************/
258 /* process X: */
259 /**************************************************************/
260 X_p=getSampleData(X,e);
261 if (NULL!=X_p) {
262 add_EM_F=TRUE;
263 if (extendedX) {
264 for (s=0;s<p.row_NS;s++) {
265 rtmp=0;
266 for (q=0;q<p.numQuad;q++) {
267 rtmp+=Vol[q]*( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX2(0,q,DIM)]
268 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX2(1,q,DIM)]
269 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*X_p[INDEX2(2,q,DIM)]);
270 }
271 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
272 }
273 } else {
274 for (s=0;s<p.row_NS;s++) {
275 rtmp0=0;
276 rtmp1=0;
277 rtmp2=0;
278 for (q=0;q<p.numQuad;q++) {
279 rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
280 rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
281 rtmp2+=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
282 }
283 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp0*X_p[0]+rtmp1*X_p[1]+rtmp2*X_p[2];
284 }
285 }
286 }
287 /**************************************************************/
288 /* process Y: */
289 /**************************************************************/
290 Y_p=getSampleData(Y,e);
291 if (NULL!=Y_p) {
292 add_EM_F=TRUE;
293 if (extendedY) {
294 for (s=0;s<p.row_NS;s++) {
295 rtmp=0;
296 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[q];
297 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
298 }
299 } else {
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)];
303 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*Y_p[0];
304 }
305 }
306 }
307 /***********************************************************************************************/
308 /* add the element matrices onto the matrix and right hand side */
309 /***********************************************************************************************/
310 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
311 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
312 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
313
314 } /* end color check */
315 } /* end element loop */
316 } /* end color loop */
317
318 THREAD_MEMFREE(EM_S);
319 THREAD_MEMFREE(EM_F);
320 THREAD_MEMFREE(row_index);
321
322 } /* end of pointer check */
323 } /* end parallel region */
324 }
325 /*
326 * $Log$
327 */

  ViewVC Help
Powered by ViewVC 1.1.26