/[escript]/trunk/dudley/src/Assemble_PDE_Single2_3D.c
ViewVC logotype

Contents of /trunk/dudley/src/Assemble_PDE_Single2_3D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2548 - (show annotations)
Mon Jul 20 06:20:06 2009 UTC (10 years ago) by jfenwick
Original Path: trunk/finley/src/Assemble_PDE_Single2_3D.c
File MIME type: text/plain
File size: 17817 byte(s)
Updating copyright notices
1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2009 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 void* ABuff=allocSampleBuffer(A);
74 void* BBuff=allocSampleBuffer(B);
75 void* CBuff=allocSampleBuffer(C);
76 void* DBuff=allocSampleBuffer(D);
77 void* XBuff=allocSampleBuffer(X);
78 void* YBuff=allocSampleBuffer(Y);
79 #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)
80 {
81 EM_S=THREAD_MEMALLOC(len_EM_S,double);
82 EM_F=THREAD_MEMALLOC(len_EM_F,double);
83 row_index=THREAD_MEMALLOC(p.row_NN,index_t);
84
85 if (!Finley_checkPtr(EM_S) && !Finley_checkPtr(EM_F) && !Finley_checkPtr(row_index) ) {
86
87 for (color=elements->minColor;color<=elements->maxColor;color++) {
88 /* open loop over all elements: */
89 #pragma omp for private(e) schedule(static)
90 for(e=0;e<elements->numElements;e++){
91 if (elements->Color[e]==color) {
92 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
93 DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]);
94 for (q=0;q<len_EM_S;++q) EM_S[q]=0;
95 for (q=0;q<len_EM_F;++q) EM_F[q]=0;
96 add_EM_F=FALSE;
97 add_EM_S=FALSE;
98 /**************************************************************/
99 /* process A: */
100 /**************************************************************/
101 A_p=getSampleDataRO(A,e,ABuff);
102 if (NULL!=A_p) {
103 add_EM_S=TRUE;
104 if (extendedA) {
105 for (s=0;s<p.row_NS;s++) {
106 for (r=0;r<p.col_NS;r++) {
107 rtmp=0;
108 for (q=0;q<p.numQuad;q++) {
109 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)]
110 + 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)]
111 + 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)]
112 + 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)]
113 + 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)]
114 + 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)]
115 + 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)]
116 + 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)]
117 + 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)]);
118 }
119 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
120 }
121 }
122 } else {
123 for (s=0;s<p.row_NS;s++) {
124 for (r=0;r<p.col_NS;r++) {
125 rtmp00=0;
126 rtmp01=0;
127 rtmp02=0;
128 rtmp10=0;
129 rtmp11=0;
130 rtmp12=0;
131 rtmp20=0;
132 rtmp21=0;
133 rtmp22=0;
134 for (q=0;q<p.numQuad;q++) {
135
136 rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
137 rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
138 rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
139 rtmp02+=rtmp0*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
140
141 rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
142 rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
143 rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
144 rtmp12+=rtmp1*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
145
146 rtmp2=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
147 rtmp20+=rtmp2*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
148 rtmp21+=rtmp2*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
149 rtmp22+=rtmp2*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
150 }
151 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp00*A_p[INDEX2(0,0,DIM)]
152 +rtmp01*A_p[INDEX2(0,1,DIM)]
153 +rtmp02*A_p[INDEX2(0,2,DIM)]
154 +rtmp10*A_p[INDEX2(1,0,DIM)]
155 +rtmp11*A_p[INDEX2(1,1,DIM)]
156 +rtmp12*A_p[INDEX2(1,2,DIM)]
157 +rtmp20*A_p[INDEX2(2,0,DIM)]
158 +rtmp21*A_p[INDEX2(2,1,DIM)]
159 +rtmp22*A_p[INDEX2(2,2,DIM)];
160 }
161 }
162 }
163 }
164 /**************************************************************/
165 /* process B: */
166 /**************************************************************/
167 B_p=getSampleDataRO(B,e,BBuff);
168 if (NULL!=B_p) {
169 add_EM_S=TRUE;
170 if (extendedB) {
171 for (s=0;s<p.row_NS;s++) {
172 for (r=0;r<p.col_NS;r++) {
173 rtmp=0;
174 for (q=0;q<p.numQuad;q++) {
175 rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]*
176 ( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX2(0,q,DIM)]
177 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX2(1,q,DIM)]
178 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*B_p[INDEX2(2,q,DIM)]);
179 }
180 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
181 }
182 }
183 } else {
184 for (s=0;s<p.row_NS;s++) {
185 for (r=0;r<p.col_NS;r++) {
186 rtmp0=0;
187 rtmp1=0;
188 rtmp2=0;
189 for (q=0;q<p.numQuad;q++) {
190 rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)];
191 rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
192 rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
193 rtmp2+=rtmp*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
194 }
195 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];
196 }
197 }
198 }
199 }
200 /**************************************************************/
201 /* process C: */
202 /**************************************************************/
203 C_p=getSampleDataRO(C,e,CBuff);
204 if (NULL!=C_p) {
205 add_EM_S=TRUE;
206 if (extendedC) {
207 for (s=0;s<p.row_NS;s++) {
208 for (r=0;r<p.col_NS;r++) {
209 rtmp=0;
210 for (q=0;q<p.numQuad;q++) {
211 rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*
212 ( C_p[INDEX2(0,q,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
213 + C_p[INDEX2(1,q,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]
214 + C_p[INDEX2(2,q,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]);
215 }
216 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
217 }
218 }
219 } else {
220 for (s=0;s<p.row_NS;s++) {
221 for (r=0;r<p.col_NS;r++) {
222 rtmp0=0;
223 rtmp1=0;
224 rtmp2=0;
225 for (q=0;q<p.numQuad;q++) {
226 rtmp=Vol[q]*S[INDEX2(s,q,p.row_NS)];
227 rtmp0+=rtmp*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
228 rtmp1+=rtmp*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
229 rtmp2+=rtmp*DSDX[INDEX3(r,2,q,p.row_NN,DIM)];
230 }
231 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];
232 }
233 }
234 }
235 }
236 /************************************************************* */
237 /* process D */
238 /**************************************************************/
239 D_p=getSampleDataRO(D,e,DBuff);
240 if (NULL!=D_p) {
241 add_EM_S=TRUE;
242 if (extendedD) {
243 for (s=0;s<p.row_NS;s++) {
244 for (r=0;r<p.col_NS;r++) {
245 rtmp=0;
246 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)];
247 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
248 }
249 }
250 } else {
251 for (s=0;s<p.row_NS;s++) {
252 for (r=0;r<p.col_NS;r++) {
253 rtmp=0;
254 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
255 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[0];
256 }
257 }
258 }
259 }
260 /**************************************************************/
261 /* process X: */
262 /**************************************************************/
263 X_p=getSampleDataRO(X,e,XBuff);
264 if (NULL!=X_p) {
265 add_EM_F=TRUE;
266 if (extendedX) {
267 for (s=0;s<p.row_NS;s++) {
268 rtmp=0;
269 for (q=0;q<p.numQuad;q++) {
270 rtmp+=Vol[q]*( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX2(0,q,DIM)]
271 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX2(1,q,DIM)]
272 + DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*X_p[INDEX2(2,q,DIM)]);
273 }
274 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
275 }
276 } else {
277 for (s=0;s<p.row_NS;s++) {
278 rtmp0=0;
279 rtmp1=0;
280 rtmp2=0;
281 for (q=0;q<p.numQuad;q++) {
282 rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
283 rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
284 rtmp2+=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)];
285 }
286 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp0*X_p[0]+rtmp1*X_p[1]+rtmp2*X_p[2];
287 }
288 }
289 }
290 /**************************************************************/
291 /* process Y: */
292 /**************************************************************/
293 Y_p=getSampleDataRO(Y,e,YBuff);
294 if (NULL!=Y_p) {
295 add_EM_F=TRUE;
296 if (extendedY) {
297 for (s=0;s<p.row_NS;s++) {
298 rtmp=0;
299 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[q];
300 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
301 }
302 } else {
303 for (s=0;s<p.row_NS;s++) {
304 rtmp=0;
305 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
306 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*Y_p[0];
307 }
308 }
309 }
310 /***********************************************************************************************/
311 /* add the element matrices onto the matrix and right hand side */
312 /***********************************************************************************************/
313 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
314 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
315 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
316
317 } /* end color check */
318 } /* end element loop */
319 } /* end color loop */
320
321 THREAD_MEMFREE(EM_S);
322 THREAD_MEMFREE(EM_F);
323 THREAD_MEMFREE(row_index);
324
325 } /* end of pointer check */
326 } /* end parallel region */
327 freeSampleBuffer(ABuff);
328 freeSampleBuffer(BBuff);
329 freeSampleBuffer(CBuff);
330 freeSampleBuffer(DBuff);
331 freeSampleBuffer(XBuff);
332 freeSampleBuffer(YBuff);
333 }
334 /*
335 * $Log$
336 */

  ViewVC Help
Powered by ViewVC 1.1.26