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

Annotation of /branches/domexper/dudley/src/Assemble_PDE_System2_3D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3187 - (hide annotations)
Thu Sep 16 02:57:17 2010 UTC (8 years, 9 months ago) by jfenwick
File MIME type: text/plain
File size: 17592 byte(s)
Enforcing indent -linux -i4 -bl -bli0 -l120

1 gross 798
2 ksteube 1312 /*******************************************************
3 ksteube 1811 *
4 jfenwick 2881 * Copyright (c) 2003-2010 by University of Queensland
5 ksteube 1811 * 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 ksteube 1312
14 gross 798 /**************************************************************/
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_{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 */
20    
21     /* u has p.numComp components in a 3D 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 = p.numEqu x 3 x p.numComp x 3 */
27     /* B = 3 x p.numEqu x p.numComp */
28     /* C = p.numEqu x 3 x p.numComp */
29     /* D = p.numEqu x p.numComp */
30     /* X = p.numEqu x 3 */
31     /* Y = p.numEqu */
32    
33     /**************************************************************/
34    
35     #include "Assemble.h"
36     #include "Util.h"
37 gross 853 #ifdef _OPENMP
38     #include <omp.h>
39     #endif
40 gross 798
41     /**************************************************************/
42    
43 jfenwick 3187 void Dudley_Assemble_PDE_System2_3D(Assemble_Parameters p, Dudley_ElementFile * elements,
44     Paso_SystemMatrix * Mat, escriptDataC * F,
45     escriptDataC * A, escriptDataC * B, escriptDataC * C, escriptDataC * D,
46     escriptDataC * X, escriptDataC * Y)
47     {
48 gross 798
49 jfenwick 3187 #define DIM 3
50 gross 798 index_t color;
51 jfenwick 3136 dim_t e;
52 jfenwick 3187 __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;
53 jfenwick 3184 double *EM_S, *EM_F, *DSDX;
54 gross 853 index_t *row_index;
55 jfenwick 3187 register dim_t q, s, r, k, m;
56 gross 853 register double rtmp, rtmp0, rtmp1, rtmp2, rtmp00, rtmp01, rtmp02, rtmp10, rtmp11, rtmp12, rtmp20, rtmp21, rtmp22;
57     bool_t add_EM_F, add_EM_S;
58    
59 jfenwick 3187 bool_t extendedA = isExpanded(A);
60     bool_t extendedB = isExpanded(B);
61     bool_t extendedC = isExpanded(C);
62     bool_t extendedD = isExpanded(D);
63     bool_t extendedX = isExpanded(X);
64     bool_t extendedY = isExpanded(Y);
65     double *F_p = (requireWrite(F), getSampleDataRW(F, 0)); /* use comma, to get around the mixed code and declarations thing */
66     double *S = p.row_jac->BasisFunctions->S;
67     dim_t len_EM_S = p.row_numShapesTotal * p.col_numShapesTotal * p.numEqu * p.numComp;
68     dim_t len_EM_F = p.row_numShapesTotal * p.numEqu;
69 gross 798
70 jfenwick 3187 #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,k,m,rtmp, rtmp0, rtmp1, rtmp2, rtmp00, rtmp01, rtmp02, rtmp10, rtmp11, rtmp12, rtmp20, rtmp21, rtmp22,add_EM_F, add_EM_S)
71 gross 798 {
72 jfenwick 3187 EM_S = THREAD_MEMALLOC(len_EM_S, double);
73     EM_F = THREAD_MEMALLOC(len_EM_F, double);
74     row_index = THREAD_MEMALLOC(p.row_numShapesTotal, index_t);
75 gross 798
76 jfenwick 3187 if (!Dudley_checkPtr(EM_S) && !Dudley_checkPtr(EM_F) && !Dudley_checkPtr(row_index))
77     {
78 gross 2748
79 jfenwick 3187 for (color = elements->minColor; color <= elements->maxColor; color++)
80     {
81     /* open loop over all elements: */
82     #pragma omp for private(e) schedule(static)
83     for (e = 0; e < elements->numElements; e++)
84     {
85     if (elements->Color[e] == color)
86     {
87 gross 2748
88 jfenwick 3187 A_p = getSampleDataRO(A, e);
89     B_p = getSampleDataRO(B, e);
90     C_p = getSampleDataRO(C, e);
91     D_p = getSampleDataRO(D, e);
92     X_p = getSampleDataRO(X, e);
93     Y_p = getSampleDataRO(Y, e);
94     double vol = p.row_jac->absD[e] * p.row_jac->quadweight;
95 gross 2748
96 jfenwick 3187 DSDX = &(p.row_jac->DSDX[INDEX5(0, 0, 0, 0, e, p.row_numShapesTotal, DIM, p.numQuadTotal, 1)]);
97     for (q = 0; q < len_EM_S; ++q)
98     EM_S[q] = 0;
99     for (q = 0; q < len_EM_F; ++q)
100     EM_F[q] = 0;
101     add_EM_F = FALSE;
102     add_EM_S = FALSE;
103 gross 798
104 jfenwick 3187 /**************************************************************/
105     /* process A: */
106     /**************************************************************/
107     if (NULL != A_p)
108     {
109     add_EM_S = TRUE;
110     if (extendedA)
111     {
112     A_q = &(A_p[INDEX6(0, 0, 0, 0, 0, 0, p.numEqu, DIM, p.numComp, DIM, p.numQuadTotal)]);
113     for (s = 0; s < p.row_numShapes; s++)
114     {
115     for (r = 0; r < p.col_numShapes; r++)
116     {
117     for (k = 0; k < p.numEqu; k++)
118     {
119     for (m = 0; m < p.numComp; m++)
120     {
121     rtmp = 0;
122     for (q = 0; q < p.numQuadTotal; q++)
123     {
124     rtmp +=
125     vol * (DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)] *
126     A_q[INDEX5(k, 0, m, 0, q, p.numEqu, DIM, p.numComp, DIM)]
127     * DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)] +
128     DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)] *
129     A_q[INDEX5(k, 0, m, 1, q, p.numEqu, DIM, p.numComp, DIM)]
130     * DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)] +
131     DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)] *
132     A_q[INDEX5(k, 0, m, 2, q, p.numEqu, DIM, p.numComp, DIM)]
133     * DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)] +
134     DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)] *
135     A_q[INDEX5(k, 1, m, 0, q, p.numEqu, DIM, p.numComp, DIM)]
136     * DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)] +
137     DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)] *
138     A_q[INDEX5(k, 1, m, 1, q, p.numEqu, DIM, p.numComp, DIM)]
139     * DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)] +
140     DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)] *
141     A_q[INDEX5(k, 1, m, 2, q, p.numEqu, DIM, p.numComp, DIM)]
142     * DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)] +
143     DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)] *
144     A_q[INDEX5(k, 2, m, 0, q, p.numEqu, DIM, p.numComp, DIM)]
145     * DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)] +
146     DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)] *
147     A_q[INDEX5(k, 2, m, 1, q, p.numEqu, DIM, p.numComp, DIM)]
148     * DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)] +
149     DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)] *
150     A_q[INDEX5(k, 2, m, 2, q, p.numEqu, DIM, p.numComp, DIM)]
151     * DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)]);
152    
153     }
154     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
155     rtmp;
156     }
157     }
158     }
159     }
160     } else
161     {
162     for (s = 0; s < p.row_numShapes; s++)
163     {
164     for (r = 0; r < p.col_numShapes; r++)
165     {
166     rtmp00 = 0;
167     rtmp01 = 0;
168     rtmp02 = 0;
169     rtmp10 = 0;
170     rtmp11 = 0;
171     rtmp12 = 0;
172     rtmp20 = 0;
173     rtmp21 = 0;
174     rtmp22 = 0;
175     for (q = 0; q < p.numQuadTotal; q++)
176     {
177     rtmp0 = vol * DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)];
178     rtmp00 += rtmp0 * DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)];
179     rtmp01 += rtmp0 * DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)];
180     rtmp02 += rtmp0 * DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)];
181    
182     rtmp1 = vol * DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)];
183     rtmp10 += rtmp1 * DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)];
184     rtmp11 += rtmp1 * DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)];
185     rtmp12 += rtmp1 * DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)];
186    
187     rtmp2 = vol * DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)];
188     rtmp20 += rtmp2 * DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)];
189     rtmp21 += rtmp2 * DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)];
190     rtmp22 += rtmp2 * DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)];
191     }
192     for (k = 0; k < p.numEqu; k++)
193     {
194     for (m = 0; m < p.numComp; m++)
195     {
196     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
197     rtmp00 * A_p[INDEX4(k, 0, m, 0, p.numEqu, DIM, p.numComp)] +
198     rtmp01 * A_p[INDEX4(k, 0, m, 1, p.numEqu, DIM, p.numComp)] +
199     rtmp02 * A_p[INDEX4(k, 0, m, 2, p.numEqu, DIM, p.numComp)] +
200     rtmp10 * A_p[INDEX4(k, 1, m, 0, p.numEqu, DIM, p.numComp)] +
201     rtmp11 * A_p[INDEX4(k, 1, m, 1, p.numEqu, DIM, p.numComp)] +
202     rtmp12 * A_p[INDEX4(k, 1, m, 2, p.numEqu, DIM, p.numComp)] +
203     rtmp20 * A_p[INDEX4(k, 2, m, 0, p.numEqu, DIM, p.numComp)] +
204     rtmp21 * A_p[INDEX4(k, 2, m, 1, p.numEqu, DIM, p.numComp)] +
205     rtmp22 * A_p[INDEX4(k, 2, m, 2, p.numEqu, DIM, p.numComp)];
206     }
207     }
208     }
209     }
210     }
211     }
212     /**************************************************************/
213     /* process B: */
214     /**************************************************************/
215     if (NULL != B_p)
216     {
217     add_EM_S = TRUE;
218     if (extendedB)
219     {
220     B_q = &(B_p[INDEX5(0, 0, 0, 0, 0, p.numEqu, DIM, p.numComp, p.numQuadTotal)]);
221     for (s = 0; s < p.row_numShapes; s++)
222     {
223     for (r = 0; r < p.col_numShapes; r++)
224     {
225     for (k = 0; k < p.numEqu; k++)
226     {
227     for (m = 0; m < p.numComp; m++)
228     {
229     rtmp = 0;
230     for (q = 0; q < p.numQuadTotal; q++)
231     {
232     rtmp +=
233     vol * S[INDEX2(r, q, p.row_numShapes)] *
234     (DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)] *
235     B_q[INDEX4(k, 0, m, q, p.numEqu, DIM, p.numComp)] +
236     DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)] *
237     B_q[INDEX4(k, 1, m, q, p.numEqu, DIM, p.numComp)] +
238     DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)] *
239     B_q[INDEX4(k, 2, m, q, p.numEqu, DIM, p.numComp)]);
240     }
241     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
242     rtmp;
243     }
244     }
245     }
246     }
247     } else
248     {
249     for (s = 0; s < p.row_numShapes; s++)
250     {
251     for (r = 0; r < p.col_numShapes; r++)
252     {
253     rtmp0 = 0;
254     rtmp1 = 0;
255     rtmp2 = 0;
256     for (q = 0; q < p.numQuadTotal; q++)
257     {
258     rtmp = vol * S[INDEX2(r, q, p.row_numShapes)];
259     rtmp0 += rtmp * DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)];
260     rtmp1 += rtmp * DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)];
261     rtmp2 += rtmp * DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)];
262     }
263     for (k = 0; k < p.numEqu; k++)
264     {
265     for (m = 0; m < p.numComp; m++)
266     {
267     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
268     rtmp0 * B_p[INDEX3(k, 0, m, p.numEqu, DIM)] +
269     rtmp1 * B_p[INDEX3(k, 1, m, p.numEqu, DIM)] +
270     rtmp2 * B_p[INDEX3(k, 2, m, p.numEqu, DIM)];
271     }
272     }
273     }
274     }
275     }
276     }
277     /**************************************************************/
278     /* process C: */
279     /**************************************************************/
280     if (NULL != C_p)
281     {
282     add_EM_S = TRUE;
283     if (extendedC)
284     {
285     C_q = &(C_p[INDEX5(0, 0, 0, 0, 0, p.numEqu, p.numComp, DIM, p.numQuadTotal)]);
286     for (s = 0; s < p.row_numShapes; s++)
287     {
288     for (r = 0; r < p.col_numShapes; r++)
289     {
290     for (k = 0; k < p.numEqu; k++)
291     {
292     for (m = 0; m < p.numComp; m++)
293     {
294     rtmp = 0;
295     for (q = 0; q < p.numQuadTotal; q++)
296     {
297     rtmp +=
298     vol * S[INDEX2(s, q, p.row_numShapes)] *
299     (C_q[INDEX4(k, m, 0, q, p.numEqu, p.numComp, DIM)] *
300     DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)] +
301     C_q[INDEX4(k, m, 1, q, p.numEqu, p.numComp, DIM)] *
302     DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)] +
303     C_q[INDEX4(k, m, 2, q, p.numEqu, p.numComp, DIM)] *
304     DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)]);
305     }
306     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
307     rtmp;
308     }
309     }
310     }
311     }
312     } else
313     {
314     for (s = 0; s < p.row_numShapes; s++)
315     {
316     for (r = 0; r < p.col_numShapes; r++)
317     {
318     rtmp0 = 0;
319     rtmp1 = 0;
320     rtmp2 = 0;
321     for (q = 0; q < p.numQuadTotal; q++)
322     {
323     rtmp = vol * S[INDEX2(s, q, p.row_numShapes)];
324     rtmp0 += rtmp * DSDX[INDEX3(r, 0, q, p.row_numShapesTotal, DIM)];
325     rtmp1 += rtmp * DSDX[INDEX3(r, 1, q, p.row_numShapesTotal, DIM)];
326     rtmp2 += rtmp * DSDX[INDEX3(r, 2, q, p.row_numShapesTotal, DIM)];
327     }
328     for (k = 0; k < p.numEqu; k++)
329     {
330     for (m = 0; m < p.numComp; m++)
331     {
332     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
333     rtmp0 * C_p[INDEX3(k, m, 0, p.numEqu, p.numComp)] +
334     rtmp1 * C_p[INDEX3(k, m, 1, p.numEqu, p.numComp)] +
335     rtmp2 * C_p[INDEX3(k, m, 2, p.numEqu, p.numComp)];
336     }
337     }
338     }
339     }
340     }
341     }
342     /************************************************************* */
343     /* process D */
344     /**************************************************************/
345     if (NULL != D_p)
346     {
347     add_EM_S = TRUE;
348     if (extendedD)
349     {
350     D_q = &(D_p[INDEX4(0, 0, 0, 0, p.numEqu, p.numComp, p.numQuadTotal)]);
351     for (s = 0; s < p.row_numShapes; s++)
352     {
353     for (r = 0; r < p.col_numShapes; r++)
354     {
355     for (k = 0; k < p.numEqu; k++)
356     {
357     for (m = 0; m < p.numComp; m++)
358     {
359     rtmp = 0;
360     for (q = 0; q < p.numQuadTotal; q++)
361     {
362     rtmp +=
363     vol * S[INDEX2(s, q, p.row_numShapes)] *
364     D_q[INDEX3(k, m, q, p.numEqu, p.numComp)] *
365     S[INDEX2(r, q, p.row_numShapes)];
366     }
367     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
368     rtmp;
369     }
370     }
371     }
372     }
373     } else
374     {
375     for (s = 0; s < p.row_numShapes; s++)
376     {
377     for (r = 0; r < p.col_numShapes; r++)
378     {
379     rtmp = 0;
380     for (q = 0; q < p.numQuadTotal; q++)
381     rtmp +=
382     vol * S[INDEX2(s, q, p.row_numShapes)] *
383     S[INDEX2(r, q, p.row_numShapes)];
384     for (k = 0; k < p.numEqu; k++)
385     {
386     for (m = 0; m < p.numComp; m++)
387     {
388     EM_S[INDEX4(k, m, s, r, p.numEqu, p.numComp, p.row_numShapesTotal)] +=
389     rtmp * D_p[INDEX2(k, m, p.numEqu)];
390     }
391     }
392     }
393     }
394     }
395     }
396     /**************************************************************/
397     /* process X: */
398     /**************************************************************/
399     if (NULL != X_p)
400     {
401     add_EM_F = TRUE;
402     if (extendedX)
403     {
404     X_q = &(X_p[INDEX4(0, 0, 0, 0, p.numEqu, DIM, p.numQuadTotal)]);
405     for (s = 0; s < p.row_numShapes; s++)
406     {
407     for (k = 0; k < p.numEqu; k++)
408     {
409     rtmp = 0;
410     for (q = 0; q < p.numQuadTotal; q++)
411     {
412     rtmp +=
413     vol * (DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)] *
414     X_q[INDEX3(k, 0, q, p.numEqu, DIM)] +
415     DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)] *
416     X_q[INDEX3(k, 1, q, p.numEqu, DIM)] +
417     DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)] *
418     X_q[INDEX3(k, 2, q, p.numEqu, DIM)]);
419     }
420     EM_F[INDEX2(k, s, p.numEqu)] += rtmp;
421     }
422     }
423     } else
424     {
425     for (s = 0; s < p.row_numShapes; s++)
426     {
427     rtmp0 = 0;
428     rtmp1 = 0;
429     rtmp2 = 0;
430     for (q = 0; q < p.numQuadTotal; q++)
431     {
432     rtmp0 += vol * DSDX[INDEX3(s, 0, q, p.row_numShapesTotal, DIM)];
433     rtmp1 += vol * DSDX[INDEX3(s, 1, q, p.row_numShapesTotal, DIM)];
434     rtmp2 += vol * DSDX[INDEX3(s, 2, q, p.row_numShapesTotal, DIM)];
435     }
436     for (k = 0; k < p.numEqu; k++)
437     {
438     EM_F[INDEX2(k, s, p.numEqu)] += rtmp0 * X_p[INDEX2(k, 0, p.numEqu)]
439     + rtmp1 * X_p[INDEX2(k, 1, p.numEqu)] + rtmp2 * X_p[INDEX2(k, 2, p.numEqu)];
440     }
441     }
442     }
443     }
444     /**************************************************************/
445     /* process Y: */
446     /**************************************************************/
447     if (NULL != Y_p)
448     {
449     add_EM_F = TRUE;
450     if (extendedY)
451     {
452     Y_q = &(Y_p[INDEX3(0, 0, 0, p.numEqu, p.numQuadTotal)]);
453     for (s = 0; s < p.row_numShapes; s++)
454     {
455     for (k = 0; k < p.numEqu; k++)
456     {
457     rtmp = 0.;
458     for (q = 0; q < p.numQuadTotal; q++)
459     rtmp +=
460     vol * S[INDEX2(s, q, p.row_numShapes)] * Y_q[INDEX2(k, q, p.numEqu)];
461     EM_F[INDEX2(k, s, p.numEqu)] += rtmp;
462     }
463     }
464     } else
465     {
466     for (s = 0; s < p.row_numShapes; s++)
467     {
468     rtmp = 0;
469     for (q = 0; q < p.numQuadTotal; q++)
470     rtmp += vol * S[INDEX2(s, q, p.row_numShapes)];
471     for (k = 0; k < p.numEqu; k++)
472     EM_F[INDEX2(k, s, p.numEqu)] += rtmp * Y_p[k];
473     }
474     }
475     }
476    
477     /***********************************************************************************************/
478     /* add the element matrices onto the matrix and right hand side */
479     /***********************************************************************************************/
480     for (q = 0; q < p.row_numShapesTotal; q++)
481     row_index[q] = p.row_DOF[elements->Nodes[INDEX2(q, e, p.NN)]];
482    
483     if (add_EM_F)
484     Dudley_Util_AddScatter(p.row_numShapesTotal, row_index, p.numEqu, EM_F, F_p,
485     p.row_DOF_UpperBound);
486     if (add_EM_S)
487     Dudley_Assemble_addToSystemMatrix(Mat, p.row_numShapesTotal, row_index, p.numEqu,
488     p.col_numShapesTotal, row_index, p.numComp, EM_S);
489     } /* end color check */
490     } /* end element loop */
491     } /* end color loop */
492    
493     THREAD_MEMFREE(EM_S);
494     THREAD_MEMFREE(EM_F);
495     THREAD_MEMFREE(row_index);
496    
497     } /* end of pointer check */
498     } /* end parallel region */
499 gross 798 }

  ViewVC Help
Powered by ViewVC 1.1.26