/[escript]/branches/doubleplusgood/dudley/src/Assemble_PDE_System2_2D.cpp
ViewVC logotype

Contents of /branches/doubleplusgood/dudley/src/Assemble_PDE_System2_2D.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2548 - (show annotations)
Mon Jul 20 06:20:06 2009 UTC (9 years, 10 months ago) by jfenwick
Original Path: trunk/finley/src/Assemble_PDE_System2_2D.c
File MIME type: text/plain
File size: 17662 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_{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 */
22
23 /* u has p.numComp components in a 2D 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 = p.numEqu x 2 x p.numComp x 2 */
29 /* B = 2 x p.numEqu x p.numComp */
30 /* C = p.numEqu x 2 x p.numComp */
31 /* D = p.numEqu x p.numComp */
32 /* X = p.numEqu x 2 */
33 /* Y = p.numEqu */
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_System2_2D(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 2
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,k,m;
59 register double rtmp, rtmp0, rtmp1, rtmp00, rtmp10, rtmp01, rtmp11;
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*p.numEqu*p.numComp;
71 dim_t len_EM_F=p.row_NN*p.numEqu;
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,k,m,rtmp, rtmp0, rtmp1, rtmp00, rtmp10, rtmp01, rtmp11,add_EM_F, add_EM_S)
80 {
81
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 for (k=0;k<p.numEqu;k++) {
109 for (m=0;m<p.numComp;m++) {
110 rtmp=0;
111 for (q=0;q<p.numQuad;q++) {
112 rtmp+=Vol[q]* (
113 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)]
114 +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)]
115 +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)]
116 +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)]);
117 }
118 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
119 }
120 }
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 rtmp10=0;
129 rtmp11=0;
130 for (q=0;q<p.numQuad;q++) {
131 rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
132 rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
133 rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
134 rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
135 rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
136 rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
137 }
138 for (k=0;k<p.numEqu;k++) {
139 for (m=0;m<p.numComp;m++) {
140 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=
141 rtmp00*A_p[INDEX4(k,0,m,0,p.numEqu,DIM,p.numComp)]
142 +rtmp01*A_p[INDEX4(k,0,m,1,p.numEqu,DIM,p.numComp)]
143 +rtmp10*A_p[INDEX4(k,1,m,0,p.numEqu,DIM,p.numComp)]
144 +rtmp11*A_p[INDEX4(k,1,m,1,p.numEqu,DIM,p.numComp)];
145 }
146 }
147 }
148 }
149 }
150 }
151 /**************************************************************/
152 /* process B: */
153 /**************************************************************/
154 B_p=getSampleDataRO(B,e,BBuff);
155 if (NULL!=B_p) {
156 add_EM_S=TRUE;
157 if (extendedB) {
158 for (s=0;s<p.row_NS;s++) {
159 for (r=0;r<p.col_NS;r++) {
160 for (k=0;k<p.numEqu;k++) {
161 for (m=0;m<p.numComp;m++) {
162 rtmp=0;
163 for (q=0;q<p.numQuad;q++) {
164 rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]*
165 ( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX4(k,0,m,q,p.numEqu,DIM,p.numComp)]
166 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX4(k,1,m,q,p.numEqu,DIM,p.numComp)]);
167 }
168 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
169 }
170 }
171 }
172 }
173 } else {
174 for (s=0;s<p.row_NS;s++) {
175 for (r=0;r<p.col_NS;r++) {
176 rtmp0=0;
177 rtmp1=0;
178 for (q=0;q<p.numQuad;q++) {
179 rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)];
180 rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
181 rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
182 }
183 for (k=0;k<p.numEqu;k++) {
184 for (m=0;m<p.numComp;m++) {
185 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+= rtmp0*B_p[INDEX3(k,0,m,p.numEqu,DIM)]
186 + rtmp1*B_p[INDEX3(k,1,m,p.numEqu,DIM)];
187 }
188 }
189 }
190 }
191 }
192 }
193 /**************************************************************/
194 /* process C: */
195 /**************************************************************/
196 C_p=getSampleDataRO(C,e,CBuff);
197 if (NULL!=C_p) {
198 add_EM_S=TRUE;
199 if (extendedC) {
200 for (s=0;s<p.row_NS;s++) {
201 for (r=0;r<p.col_NS;r++) {
202 for (k=0;k<p.numEqu;k++) {
203 for (m=0;m<p.numComp;m++) {
204 rtmp=0;
205 for (q=0;q<p.numQuad;q++) {
206 rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*
207 ( C_p[INDEX4(k,m,0,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]
208 + C_p[INDEX4(k,m,1,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]);
209 }
210 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
211 }
212 }
213 }
214 }
215 } else {
216 for (s=0;s<p.row_NS;s++) {
217 for (r=0;r<p.col_NS;r++) {
218 rtmp0=0;
219 rtmp1=0;
220 for (q=0;q<p.numQuad;q++) {
221 rtmp=Vol[q]*S[INDEX2(s,q,p.row_NS)];
222 rtmp0+=rtmp*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
223 rtmp1+=rtmp*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
224 }
225 for (k=0;k<p.numEqu;k++) {
226 for (m=0;m<p.numComp;m++) {
227 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)]
228 +rtmp1*C_p[INDEX3(k,m,1,p.numEqu,p.numComp)];
229 }
230 }
231 }
232 }
233 }
234 }
235 /************************************************************* */
236 /* process D */
237 /**************************************************************/
238 D_p=getSampleDataRO(D,e,DBuff);
239 if (NULL!=D_p) {
240 add_EM_S=TRUE;
241 if (extendedD) {
242 for (s=0;s<p.row_NS;s++) {
243 for (r=0;r<p.col_NS;r++) {
244 for (k=0;k<p.numEqu;k++) {
245 for (m=0;m<p.numComp;m++) {
246 rtmp=0;
247 for (q=0;q<p.numQuad;q++) {
248 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)];
249 }
250 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
251 }
252 }
253 }
254 }
255 } else {
256 for (s=0;s<p.row_NS;s++) {
257 for (r=0;r<p.col_NS;r++) {
258 rtmp=0;
259 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
260 for (k=0;k<p.numEqu;k++) {
261 for (m=0;m<p.numComp;m++) {
262 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[INDEX2(k,m,p.numEqu)];
263 }
264 }
265 }
266 }
267 }
268 }
269 /**************************************************************/
270 /* process X: */
271 /**************************************************************/
272 X_p=getSampleDataRO(X,e,XBuff);
273 if (NULL!=X_p) {
274 add_EM_F=TRUE;
275 if (extendedX) {
276 for (s=0;s<p.row_NS;s++) {
277 for (k=0;k<p.numEqu;k++) {
278 rtmp=0;
279 for (q=0;q<p.numQuad;q++) {
280 rtmp+=Vol[q]*(DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX3(k,0,q,p.numEqu,DIM)]
281 +DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX3(k,1,q,p.numEqu,DIM)]);
282 }
283 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
284 }
285 }
286 } else {
287 for (s=0;s<p.row_NS;s++) {
288 rtmp0=0;
289 rtmp1=0;
290 for (q=0;q<p.numQuad;q++) {
291 rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
292 rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
293 }
294 for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp0*X_p[INDEX2(k,0,p.numEqu)]+rtmp1*X_p[INDEX2(k,1,p.numEqu)];
295 }
296 }
297 }
298 /**************************************************************/
299 /* process Y: */
300 /**************************************************************/
301 Y_p=getSampleDataRO(Y,e,YBuff);
302 if (NULL!=Y_p) {
303 add_EM_F=TRUE;
304 if (extendedY) {
305 for (s=0;s<p.row_NS;s++) {
306 for (k=0;k<p.numEqu;k++) {
307 rtmp=0;
308 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[INDEX2(k,q,p.numEqu)];
309 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
310 }
311 }
312 } else {
313 for (s=0;s<p.row_NS;s++) {
314 rtmp=0;
315 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
316 for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*Y_p[k];
317 }
318 }
319 }
320 /***********************************************************************************************/
321 /* add the element matrices onto the matrix and right hand side */
322 /***********************************************************************************************/
323 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
324 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
325 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
326
327 } /* end color check */
328 } /* end element loop */
329 } /* end color loop */
330
331 THREAD_MEMFREE(EM_S);
332 THREAD_MEMFREE(EM_F);
333 THREAD_MEMFREE(row_index);
334
335 } /* end of pointer check */
336 } /* end parallel region */
337 freeSampleBuffer(ABuff);
338 freeSampleBuffer(BBuff);
339 freeSampleBuffer(CBuff);
340 freeSampleBuffer(DBuff);
341 freeSampleBuffer(XBuff);
342 freeSampleBuffer(YBuff);
343 }
344 /*
345 * $Log$
346 */

  ViewVC Help
Powered by ViewVC 1.1.26