/[escript]/branches/domexper/dudley/src/Assemble_PDE_Single2_2D.c
ViewVC logotype

Contents of /branches/domexper/dudley/src/Assemble_PDE_Single2_2D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3205 - (show annotations)
Fri Sep 24 00:30:43 2010 UTC (8 years, 10 months ago) by jfenwick
File MIME type: text/plain
File size: 12084 byte(s)
Removing -total from field names

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(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, Vol, 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
90 A_p = getSampleDataRO(A, e);
91 B_p = getSampleDataRO(B, e);
92 C_p = getSampleDataRO(C, e);
93 D_p = getSampleDataRO(D, e);
94 X_p = getSampleDataRO(X, e);
95 Y_p = getSampleDataRO(Y, e);
96
97 double vol = p.row_jac->absD[e] * p.row_jac->quadweight;
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 } else
141 {
142 for (s = 0; s < p.numShapes; s++)
143 {
144 for (r = 0; r < p.numShapes; r++)
145 {
146 rtmp00 = 0;
147 rtmp01 = 0;
148 rtmp10 = 0;
149 rtmp11 = 0;
150 for (q = 0; q < p.numQuad; q++)
151 {
152 rtmp0 = vol * DSDX[INDEX3(s, 0, q, p.numShapes, DIM)];
153 rtmp1 = vol * DSDX[INDEX3(s, 1, q, p.numShapes, DIM)];
154 rtmp00 += rtmp0 * DSDX[INDEX3(r, 0, q, p.numShapes, DIM)];
155 rtmp01 += rtmp0 * DSDX[INDEX3(r, 1, q, p.numShapes, DIM)];
156 rtmp10 += rtmp1 * DSDX[INDEX3(r, 0, q, p.numShapes, DIM)];
157 rtmp11 += rtmp1 * DSDX[INDEX3(r, 1, q, p.numShapes, DIM)];
158 }
159 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] +=
160 rtmp00 * A_p[INDEX2(0, 0, DIM)] + rtmp01 * A_p[INDEX2(0, 1, DIM)] +
161 rtmp10 * A_p[INDEX2(1, 0, DIM)] + rtmp11 * A_p[INDEX2(1, 1, DIM)];
162 }
163 }
164 }
165 }
166 /**************************************************************/
167 /* process B: */
168 /**************************************************************/
169 if (NULL != B_p)
170 {
171 add_EM_S = TRUE;
172 if (extendedB)
173 {
174 B_q = &(B_p[INDEX3(0, 0, 0, DIM, p.numQuad)]);
175 for (s = 0; s < p.numShapes; s++)
176 {
177 for (r = 0; r < p.numShapes; r++)
178 {
179 rtmp = 0.;
180 for (q = 0; q < p.numQuad; q++)
181 {
182 rtmp +=
183 vol * S[INDEX2(r, q, p.numShapes)] *
184 (DSDX[INDEX3(s, 0, q, p.numShapes, DIM)] *
185 B_q[INDEX2(0, q, DIM)] +
186 DSDX[INDEX3(s, 1, q, p.numShapes, DIM)] *
187 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 } else
193 {
194 for (s = 0; s < p.numShapes; s++)
195 {
196 for (r = 0; r < p.numShapes; r++)
197 {
198 rtmp0 = 0;
199 rtmp1 = 0;
200 for (q = 0; q < p.numQuad; q++)
201 {
202 rtmp = vol * S[INDEX2(r, q, p.numShapes)];
203 rtmp0 += rtmp * DSDX[INDEX3(s, 0, q, p.numShapes, DIM)];
204 rtmp1 += rtmp * DSDX[INDEX3(s, 1, q, p.numShapes, DIM)];
205 }
206 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] +=
207 rtmp0 * B_p[0] + rtmp1 * B_p[1];
208 }
209 }
210 }
211 }
212 /**************************************************************/
213 /* process C: */
214 /**************************************************************/
215 if (NULL != C_p)
216 {
217 add_EM_S = TRUE;
218 if (extendedC)
219 {
220 C_q = &(C_p[INDEX3(0, 0, 0, DIM, p.numQuad)]);
221 for (s = 0; s < p.numShapes; s++)
222 {
223 for (r = 0; r < p.numShapes; r++)
224 {
225 rtmp = 0;
226 for (q = 0; q < p.numQuad; q++)
227 {
228 rtmp +=
229 vol * S[INDEX2(s, q, p.numShapes)] * (C_q[INDEX2(0, q, DIM)] *
230 DSDX[INDEX3
231 (r, 0, q,
232 p.numShapes,
233 DIM)] + C_q[INDEX2(1, q,
234 DIM)]
235 *
236 DSDX[INDEX3
237 (r, 1, q,
238 p.numShapes,
239 DIM)]);
240 }
241 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] += rtmp;
242 }
243 }
244 } else
245 {
246 for (s = 0; s < p.numShapes; s++)
247 {
248 for (r = 0; r < p.numShapes; r++)
249 {
250 rtmp0 = 0;
251 rtmp1 = 0;
252 for (q = 0; q < p.numQuad; q++)
253 {
254 rtmp = vol * S[INDEX2(s, q, p.numShapes)];
255 rtmp0 += rtmp * DSDX[INDEX3(r, 0, q, p.numShapes, DIM)];
256 rtmp1 += rtmp * DSDX[INDEX3(r, 1, q, p.numShapes, DIM)];
257 }
258 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] +=
259 rtmp0 * C_p[0] + rtmp1 * C_p[1];
260 }
261 }
262 }
263 }
264 /************************************************************* */
265 /* process D */
266 /**************************************************************/
267 if (NULL != D_p)
268 {
269 add_EM_S = TRUE;
270 if (extendedD)
271 {
272 D_q = &(D_p[INDEX2(0, 0, p.numQuad)]);
273 for (s = 0; s < p.numShapes; s++)
274 {
275 for (r = 0; r < p.numShapes; r++)
276 {
277 rtmp = 0;
278 for (q = 0; q < p.numQuad; q++)
279 rtmp +=
280 vol * S[INDEX2(s, q, p.numShapes)] * D_q[q] *
281 S[INDEX2(r, q, p.numShapes)];
282 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] += rtmp;
283 }
284 }
285 } else
286 {
287 for (s = 0; s < p.numShapes; s++)
288 {
289 for (r = 0; r < p.numShapes; r++)
290 {
291 rtmp = 0;
292 for (q = 0; q < p.numQuad; q++)
293 rtmp +=
294 vol * S[INDEX2(s, q, p.numShapes)] *
295 S[INDEX2(r, q, p.numShapes)];
296 EM_S[INDEX4(0, 0, s, r, p.numEqu, p.numComp, p.numShapes)] +=
297 rtmp * D_p[0];
298 }
299 }
300 }
301 }
302 /**************************************************************/
303 /* process X: */
304 /**************************************************************/
305 if (NULL != X_p)
306 {
307 add_EM_F = TRUE;
308 if (extendedX)
309 {
310 X_q = &(X_p[INDEX3(0, 0, 0, DIM, p.numQuad)]);
311 for (s = 0; s < p.numShapes; s++)
312 {
313 rtmp = 0.;
314 for (q = 0; q < p.numQuad; q++)
315 {
316 rtmp +=
317 vol * (DSDX[INDEX3(s, 0, q, p.numShapes, DIM)] *
318 X_q[INDEX2(0, q, DIM)] +
319 DSDX[INDEX3(s, 1, q, p.numShapes, DIM)] *
320 X_q[INDEX2(1, q, DIM)]);
321 }
322 EM_F[INDEX2(0, s, p.numEqu)] += rtmp;
323 }
324 } else
325 {
326 for (s = 0; s < p.numShapes; s++)
327 {
328 rtmp0 = 0.;
329 rtmp1 = 0.;
330 for (q = 0; q < p.numQuad; q++)
331 {
332 rtmp0 += vol * DSDX[INDEX3(s, 0, q, p.numShapes, DIM)];
333 rtmp1 += vol * DSDX[INDEX3(s, 1, q, p.numShapes, DIM)];
334 }
335 EM_F[INDEX2(0, s, p.numEqu)] += rtmp0 * X_p[0] + rtmp1 * X_p[1];
336 }
337 }
338 }
339 /**************************************************************/
340 /* process Y: */
341 /**************************************************************/
342 if (NULL != Y_p)
343 {
344 add_EM_F = TRUE;
345 if (extendedY)
346 {
347 Y_q = &(Y_p[INDEX2(0, 0, p.numQuad)]);
348 for (s = 0; s < p.numShapes; s++)
349 {
350 rtmp = 0;
351 for (q = 0; q < p.numQuad; q++)
352 rtmp += vol * S[INDEX2(s, q, p.numShapes)] * Y_q[q];
353 EM_F[INDEX2(0, s, p.numEqu)] += rtmp;
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,
374 p.row_DOF_UpperBound);
375 if (add_EM_S)
376 Dudley_Assemble_addToSystemMatrix(Mat, p.numShapes, row_index, p.numEqu,
377 p.numShapes, row_index, p.numComp, EM_S);
378 } /* end color check */
379 } /* end element loop */
380 } /* end color loop */
381
382 THREAD_MEMFREE(EM_S);
383 THREAD_MEMFREE(EM_F);
384 THREAD_MEMFREE(row_index);
385
386 } /* end of pointer check */
387 } /* end parallel region */
388 }

  ViewVC Help
Powered by ViewVC 1.1.26