/[escript]/branches/arrexp_2137_win_merge/finley/src/Assemble_PDE_Single2_3D.c
ViewVC logotype

Contents of /branches/arrexp_2137_win_merge/finley/src/Assemble_PDE_Single2_3D.c

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.26