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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3259 - (show annotations)
Mon Oct 11 01:48:14 2010 UTC (9 years ago) by jfenwick
File MIME type: text/plain
File size: 12083 byte(s)
Merging dudley and scons updates from branches

1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2010 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 /* assembles the system of numEq PDEs into the stiffness matrix S right hand side F */
17 /* the shape functions for test and solution must be identical */
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 2D 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 = 2 x 2 */
27 /* B = 2 */
28 /* C = 2 */
29 /* D = scalar */
30 /* X = 2 */
31 /* Y = scalar */
32
33 /**************************************************************/
34
35 #include "Assemble.h"
36 #include "Util.h"
37 #ifdef _OPENMP
38 #include <omp.h>
39 #endif
40
41 #include "ShapeTable.h"
42
43 /**************************************************************/
44
45 void Dudley_Assemble_PDE_Single2_2D(Dudley_Assemble_Parameters p, Dudley_ElementFile * elements,
46 Paso_SystemMatrix * Mat, escriptDataC * F,
47 escriptDataC * A, escriptDataC * B, escriptDataC * C, escriptDataC * D,
48 escriptDataC * X, escriptDataC * Y)
49 {
50
51 #define DIM 2
52 index_t color;
53 dim_t e;
54 __const double *A_p, *B_p, *C_p, *D_p, *X_p, *Y_p, *A_q, *B_q, *C_q, *D_q, *X_q, *Y_q;
55 double *EM_S, *EM_F, *DSDX;
56 index_t *row_index;
57 register dim_t q, s, r;
58 register double rtmp00, rtmp01, rtmp10, rtmp11, rtmp, rtmp0, rtmp1;
59 bool_t add_EM_F, add_EM_S;
60
61 bool_t extendedA = isExpanded(A);
62 bool_t extendedB = isExpanded(B);
63 bool_t extendedC = isExpanded(C);
64 bool_t extendedD = isExpanded(D);
65 bool_t extendedX = isExpanded(X);
66 bool_t extendedY = isExpanded(Y);
67 double *F_p = (requireWrite(F), getSampleDataRW(F, 0)); /* use comma, to get around the mixed code and declarations thing */
68 const double *S = p.shapeFns;
69 dim_t len_EM_S = p.numShapes * p.numShapes;
70 dim_t len_EM_F = p.numShapes;
71
72 #pragma omp parallel private(color,EM_S, EM_F, DSDX, A_p, B_p, C_p, D_p, X_p, Y_p, A_q, B_q, C_q, D_q, X_q, Y_q,row_index,q, s,r,rtmp00, rtmp01, rtmp10, rtmp11, rtmp, rtmp0, rtmp1,add_EM_F, add_EM_S)
73 {
74 EM_S = THREAD_MEMALLOC(len_EM_S, double);
75 EM_F = THREAD_MEMALLOC(len_EM_F, double);
76 row_index = THREAD_MEMALLOC(p.numShapes, index_t);
77
78 if (!Dudley_checkPtr(EM_S) && !Dudley_checkPtr(EM_F) && !Dudley_checkPtr(row_index))
79 {
80
81 for (color = elements->minColor; color <= elements->maxColor; color++)
82 {
83 /* open loop over all elements: */
84 #pragma omp for private(e) schedule(static)
85 for (e = 0; e < elements->numElements; e++)
86 {
87 if (elements->Color[e] == color)
88 {
89 double vol = p.row_jac->absD[e] * p.row_jac->quadweight;
90
91 A_p = getSampleDataRO(A, e);
92 B_p = getSampleDataRO(B, e);
93 C_p = getSampleDataRO(C, e);
94 D_p = getSampleDataRO(D, e);
95 X_p = getSampleDataRO(X, e);
96 Y_p = getSampleDataRO(Y, e);
97
98
99 DSDX = &(p.row_jac->DSDX[INDEX5(0, 0, 0, 0, e, p.numShapes, DIM, p.numQuad, 1)]);
100 for (q = 0; q < len_EM_S; ++q)
101 EM_S[q] = 0;
102 for (q = 0; q < len_EM_F; ++q)
103 EM_F[q] = 0;
104 add_EM_F = FALSE;
105 add_EM_S = FALSE;
106 /**************************************************************/
107 /* process A: */
108 /**************************************************************/
109 if (NULL != A_p)
110 {
111 add_EM_S = TRUE;
112 if (extendedA)
113 {
114 A_q = &(A_p[INDEX4(0, 0, 0, 0, DIM, DIM, p.numQuad)]);
115 for (s = 0; s < p.numShapes; s++)
116 {
117 for (r = 0; r < p.numShapes; r++)
118 {
119 rtmp = 0;
120 for (q = 0; q < p.numQuad; q++)
121 {
122 rtmp +=
123 vol * (DSDX[INDEX3(s, 0, q, p.numShapes, DIM)] *
124 A_q[INDEX3(0, 0, q, DIM, DIM)] *
125 DSDX[INDEX3(r, 0, q, p.numShapes, DIM)] +
126 DSDX[INDEX3(s, 0, q, p.numShapes, DIM)] *
127 A_q[INDEX3(0, 1, q, DIM, DIM)] *
128 DSDX[INDEX3(r, 1, q, p.numShapes, DIM)] +
129 DSDX[INDEX3(s, 1, q, p.numShapes, DIM)] *
130 A_q[INDEX3(1, 0, q, DIM, DIM)] *
131 DSDX[INDEX3(r, 0, q, p.numShapes, DIM)] +
132 DSDX[INDEX3(s, 1, q, p.numShapes, DIM)] *
133 A_q[INDEX3(1, 1, q, DIM, DIM)] *
134 DSDX[INDEX3(r, 1, q, p.numShapes, DIM)]);
135
136 }
137 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] += rtmp;
138 }
139 }
140 }
141 else
142 {
143 for (s = 0; s < p.numShapes; s++)
144 {
145 for (r = 0; r < p.numShapes; r++)
146 {
147 rtmp00 = 0;
148 rtmp01 = 0;
149 rtmp10 = 0;
150 rtmp11 = 0;
151 for (q = 0; q < p.numQuad; q++)
152 {
153 rtmp0 = vol * DSDX[INDEX3(s, 0, q, p.numShapes, DIM)];
154 rtmp1 = vol * DSDX[INDEX3(s, 1, q, p.numShapes, DIM)];
155 rtmp00 += rtmp0 * DSDX[INDEX3(r, 0, q, p.numShapes, DIM)];
156 rtmp01 += rtmp0 * DSDX[INDEX3(r, 1, q, p.numShapes, DIM)];
157 rtmp10 += rtmp1 * DSDX[INDEX3(r, 0, q, p.numShapes, DIM)];
158 rtmp11 += rtmp1 * DSDX[INDEX3(r, 1, q, p.numShapes, DIM)];
159 }
160 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] +=
161 rtmp00 * A_p[INDEX2(0, 0, DIM)] + rtmp01 * A_p[INDEX2(0, 1, DIM)] +
162 rtmp10 * A_p[INDEX2(1, 0, DIM)] + rtmp11 * A_p[INDEX2(1, 1, DIM)];
163 }
164 }
165 }
166 }
167 /**************************************************************/
168 /* process B: */
169 /**************************************************************/
170 if (NULL != B_p)
171 {
172 add_EM_S = TRUE;
173 if (extendedB)
174 {
175 B_q = &(B_p[INDEX3(0, 0, 0, DIM, p.numQuad)]);
176 for (s = 0; s < p.numShapes; s++)
177 {
178 for (r = 0; r < p.numShapes; r++)
179 {
180 rtmp = 0.;
181 for (q = 0; q < p.numQuad; q++)
182 {
183 rtmp +=
184 vol * S[INDEX2(r, q, p.numShapes)] *
185 (DSDX[INDEX3(s, 0, q, p.numShapes, DIM)] *
186 B_q[INDEX2(0, q, DIM)] +
187 DSDX[INDEX3(s, 1, q, p.numShapes, DIM)] * B_q[INDEX2(1, q, DIM)]);
188 }
189 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] += rtmp;
190 }
191 }
192 }
193 else
194 {
195 for (s = 0; s < p.numShapes; s++)
196 {
197 for (r = 0; r < p.numShapes; r++)
198 {
199 rtmp0 = 0;
200 rtmp1 = 0;
201 for (q = 0; q < p.numQuad; q++)
202 {
203 rtmp = vol * S[INDEX2(r, q, p.numShapes)];
204 rtmp0 += rtmp * DSDX[INDEX3(s, 0, q, p.numShapes, DIM)];
205 rtmp1 += rtmp * DSDX[INDEX3(s, 1, q, p.numShapes, DIM)];
206 }
207 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] +=
208 rtmp0 * B_p[0] + rtmp1 * B_p[1];
209 }
210 }
211 }
212 }
213 /**************************************************************/
214 /* process C: */
215 /**************************************************************/
216 if (NULL != C_p)
217 {
218 add_EM_S = TRUE;
219 if (extendedC)
220 {
221 C_q = &(C_p[INDEX3(0, 0, 0, DIM, p.numQuad)]);
222 for (s = 0; s < p.numShapes; s++)
223 {
224 for (r = 0; r < p.numShapes; r++)
225 {
226 rtmp = 0;
227 for (q = 0; q < p.numQuad; q++)
228 {
229 rtmp +=
230 vol * S[INDEX2(s, q, p.numShapes)] * (C_q[INDEX2(0, q, DIM)] *
231 DSDX[INDEX3
232 (r, 0, q,
233 p.numShapes,
234 DIM)] + C_q[INDEX2(1, q,
235 DIM)]
236 *
237 DSDX[INDEX3
238 (r, 1, q,
239 p.numShapes, DIM)]);
240 }
241 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] += rtmp;
242 }
243 }
244 }
245 else
246 {
247 for (s = 0; s < p.numShapes; s++)
248 {
249 for (r = 0; r < p.numShapes; r++)
250 {
251 rtmp0 = 0;
252 rtmp1 = 0;
253 for (q = 0; q < p.numQuad; q++)
254 {
255 rtmp = vol * S[INDEX2(s, q, p.numShapes)];
256 rtmp0 += rtmp * DSDX[INDEX3(r, 0, q, p.numShapes, DIM)];
257 rtmp1 += rtmp * DSDX[INDEX3(r, 1, q, p.numShapes, DIM)];
258 }
259 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] +=
260 rtmp0 * C_p[0] + rtmp1 * C_p[1];
261 }
262 }
263 }
264 }
265 /************************************************************* */
266 /* process D */
267 /**************************************************************/
268 if (NULL != D_p)
269 {
270 add_EM_S = TRUE;
271 if (extendedD)
272 {
273 D_q = &(D_p[INDEX2(0, 0, p.numQuad)]);
274 for (s = 0; s < p.numShapes; s++)
275 {
276 for (r = 0; r < p.numShapes; r++)
277 {
278 rtmp = 0;
279 for (q = 0; q < p.numQuad; q++)
280 rtmp +=
281 vol * S[INDEX2(s, q, p.numShapes)] * D_q[q] *
282 S[INDEX2(r, q, p.numShapes)];
283 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] += rtmp;
284 }
285 }
286 }
287 else
288 {
289 for (s = 0; s < p.numShapes; s++)
290 {
291 for (r = 0; r < p.numShapes; r++)
292 {
293 rtmp = 0;
294 for (q = 0; q < p.numQuad; q++)
295 rtmp += vol * S[INDEX2(s, q, p.numShapes)] * S[INDEX2(r, q, p.numShapes)];
296 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] += rtmp * D_p[0];
297 }
298 }
299 }
300 }
301 /**************************************************************/
302 /* process X: */
303 /**************************************************************/
304 if (NULL != X_p)
305 {
306 add_EM_F = TRUE;
307 if (extendedX)
308 {
309 X_q = &(X_p[INDEX3(0, 0, 0, DIM, p.numQuad)]);
310 for (s = 0; s < p.numShapes; s++)
311 {
312 rtmp = 0.;
313 for (q = 0; q < p.numQuad; q++)
314 {
315 rtmp +=
316 vol * (DSDX[INDEX3(s, 0, q, p.numShapes, DIM)] *
317 X_q[INDEX2(0, q, DIM)] +
318 DSDX[INDEX3(s, 1, q, p.numShapes, DIM)] * X_q[INDEX2(1, q, DIM)]);
319 }
320 EM_F[INDEX2(0, s, p.numEqu)] += rtmp;
321 }
322 }
323 else
324 {
325 for (s = 0; s < p.numShapes; s++)
326 {
327 rtmp0 = 0.;
328 rtmp1 = 0.;
329 for (q = 0; q < p.numQuad; q++)
330 {
331 rtmp0 += vol * DSDX[INDEX3(s, 0, q, p.numShapes, DIM)];
332 rtmp1 += vol * DSDX[INDEX3(s, 1, q, p.numShapes, DIM)];
333 }
334 EM_F[INDEX2(0, s, p.numEqu)] += rtmp0 * X_p[0] + rtmp1 * X_p[1];
335 }
336 }
337 }
338 /**************************************************************/
339 /* process Y: */
340 /**************************************************************/
341 if (NULL != Y_p)
342 {
343 add_EM_F = TRUE;
344 if (extendedY)
345 {
346 Y_q = &(Y_p[INDEX2(0, 0, p.numQuad)]);
347 for (s = 0; s < p.numShapes; s++)
348 {
349 rtmp = 0;
350 for (q = 0; q < p.numQuad; q++)
351 rtmp += vol * S[INDEX2(s, q, p.numShapes)] * Y_q[q];
352 EM_F[INDEX2(0, s, p.numEqu)] += rtmp;
353 }
354 }
355 else
356 {
357 for (s = 0; s < p.numShapes; s++)
358 {
359 rtmp = 0;
360 for (q = 0; q < p.numQuad; q++)
361 rtmp += vol * S[INDEX2(s, q, p.numShapes)];
362 EM_F[INDEX2(0, s, p.numEqu)] += rtmp * Y_p[0];
363 }
364 }
365 }
366 /***********************************************************************************************/
367 /* add the element matrices onto the matrix and right hand side */
368 /***********************************************************************************************/
369
370 for (q = 0; q < p.numShapes; q++)
371 row_index[q] = p.row_DOF[elements->Nodes[INDEX2(q, e, p.NN)]];
372 if (add_EM_F)
373 Dudley_Util_AddScatter(p.numShapes, row_index, p.numEqu, EM_F, F_p, p.row_DOF_UpperBound);
374 if (add_EM_S)
375 Dudley_Assemble_addToSystemMatrix(Mat, p.numShapes, row_index, p.numEqu,
376 p.numShapes, row_index, p.numComp, EM_S);
377 } /* end color check */
378 } /* end element loop */
379 } /* end color loop */
380
381 THREAD_MEMFREE(EM_S);
382 THREAD_MEMFREE(EM_F);
383 THREAD_MEMFREE(row_index);
384
385 } /* end of pointer check */
386 } /* end parallel region */
387 }

  ViewVC Help
Powered by ViewVC 1.1.26