1 |
|
2 |
/* $Id$ */ |
3 |
|
4 |
/******************************************************* |
5 |
* |
6 |
* Copyright 2003-2007 by ACceSS MNRF |
7 |
* Copyright 2007 by University of Queensland |
8 |
* |
9 |
* http://esscc.uq.edu.au |
10 |
* Primary Business: Queensland, Australia |
11 |
* Licensed under the Open Software License version 3.0 |
12 |
* http://www.opensource.org/licenses/osl-3.0.php |
13 |
* |
14 |
*******************************************************/ |
15 |
|
16 |
/**************************************************************/ |
17 |
|
18 |
/* assembles the system of numEq PDEs into the stiffness matrix S right hand side F */ |
19 |
/* the shape functions for test and solution must be identical */ |
20 |
|
21 |
|
22 |
/* -(A_{i,j} u_,j)_i-(B_{i} u)_i+C_{j} u_,j-D u_m and -(X_,i)_i + Y */ |
23 |
|
24 |
/* in a 3D domain. The shape functions for test and solution must be identical */ |
25 |
/* and row_NS == row_NN */ |
26 |
|
27 |
/* Shape of the coefficients: */ |
28 |
|
29 |
/* A = 3 x 3 */ |
30 |
/* B = 3 */ |
31 |
/* C = 3 */ |
32 |
/* D = scalar */ |
33 |
/* X = 3 */ |
34 |
/* Y = scalar */ |
35 |
|
36 |
|
37 |
/**************************************************************/ |
38 |
|
39 |
|
40 |
#include "Assemble.h" |
41 |
#include "Util.h" |
42 |
#ifdef _OPENMP |
43 |
#include <omp.h> |
44 |
#endif |
45 |
|
46 |
|
47 |
/**************************************************************/ |
48 |
|
49 |
void Finley_Assemble_PDE_Single2_3D(Assemble_Parameters p, Finley_ElementFile* elements, |
50 |
Paso_SystemMatrix* Mat, escriptDataC* F, |
51 |
escriptDataC* A, escriptDataC* B, escriptDataC* C, escriptDataC* D, escriptDataC* X, escriptDataC* Y) { |
52 |
|
53 |
#define DIM 3 |
54 |
index_t color; |
55 |
dim_t e; |
56 |
double *EM_S, *EM_F, *Vol, *DSDX, *A_p, *B_p, *C_p, *D_p, *X_p, *Y_p; |
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=getSampleData(F,0); |
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 |
|
74 |
#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) |
75 |
{ |
76 |
EM_S=THREAD_MEMALLOC(len_EM_S,double); |
77 |
EM_F=THREAD_MEMALLOC(len_EM_F,double); |
78 |
row_index=THREAD_MEMALLOC(p.row_NN,index_t); |
79 |
|
80 |
if (!Finley_checkPtr(EM_S) && !Finley_checkPtr(EM_F) && !Finley_checkPtr(row_index) ) { |
81 |
|
82 |
for (color=elements->minColor;color<=elements->maxColor;color++) { |
83 |
/* open loop over all elements: */ |
84 |
#pragma omp for private(e) schedule(static) |
85 |
for(e=0;e<elements->numElements;e++){ |
86 |
if (elements->Color[e]==color) { |
87 |
Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]); |
88 |
DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]); |
89 |
for (q=0;q<len_EM_S;++q) EM_S[q]=0; |
90 |
for (q=0;q<len_EM_F;++q) EM_F[q]=0; |
91 |
add_EM_F=FALSE; |
92 |
add_EM_S=FALSE; |
93 |
/**************************************************************/ |
94 |
/* process A: */ |
95 |
/**************************************************************/ |
96 |
A_p=getSampleData(A,e); |
97 |
if (NULL!=A_p) { |
98 |
add_EM_S=TRUE; |
99 |
if (extendedA) { |
100 |
for (s=0;s<p.row_NS;s++) { |
101 |
for (r=0;r<p.col_NS;r++) { |
102 |
rtmp=0; |
103 |
for (q=0;q<p.numQuad;q++) { |
104 |
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)] |
105 |
+ 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)] |
106 |
+ 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)] |
107 |
+ 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)] |
108 |
+ 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)] |
109 |
+ 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)] |
110 |
+ 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)] |
111 |
+ 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)] |
112 |
+ 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)]); |
113 |
} |
114 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
115 |
} |
116 |
} |
117 |
} else { |
118 |
for (s=0;s<p.row_NS;s++) { |
119 |
for (r=0;r<p.col_NS;r++) { |
120 |
rtmp00=0; |
121 |
rtmp01=0; |
122 |
rtmp02=0; |
123 |
rtmp10=0; |
124 |
rtmp11=0; |
125 |
rtmp12=0; |
126 |
rtmp20=0; |
127 |
rtmp21=0; |
128 |
rtmp22=0; |
129 |
for (q=0;q<p.numQuad;q++) { |
130 |
|
131 |
rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]; |
132 |
rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
133 |
rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]; |
134 |
rtmp02+=rtmp0*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]; |
135 |
|
136 |
rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)]; |
137 |
rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
138 |
rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]; |
139 |
rtmp12+=rtmp1*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]; |
140 |
|
141 |
rtmp2=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)]; |
142 |
rtmp20+=rtmp2*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
143 |
rtmp21+=rtmp2*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]; |
144 |
rtmp22+=rtmp2*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]; |
145 |
} |
146 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp00*A_p[INDEX2(0,0,DIM)] |
147 |
+rtmp01*A_p[INDEX2(0,1,DIM)] |
148 |
+rtmp02*A_p[INDEX2(0,2,DIM)] |
149 |
+rtmp10*A_p[INDEX2(1,0,DIM)] |
150 |
+rtmp11*A_p[INDEX2(1,1,DIM)] |
151 |
+rtmp12*A_p[INDEX2(1,2,DIM)] |
152 |
+rtmp20*A_p[INDEX2(2,0,DIM)] |
153 |
+rtmp21*A_p[INDEX2(2,1,DIM)] |
154 |
+rtmp22*A_p[INDEX2(2,2,DIM)]; |
155 |
} |
156 |
} |
157 |
} |
158 |
} |
159 |
/**************************************************************/ |
160 |
/* process B: */ |
161 |
/**************************************************************/ |
162 |
B_p=getSampleData(B,e); |
163 |
if (NULL!=B_p) { |
164 |
add_EM_S=TRUE; |
165 |
if (extendedB) { |
166 |
for (s=0;s<p.row_NS;s++) { |
167 |
for (r=0;r<p.col_NS;r++) { |
168 |
rtmp=0; |
169 |
for (q=0;q<p.numQuad;q++) { |
170 |
rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]* |
171 |
( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX2(0,q,DIM)] |
172 |
+ DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX2(1,q,DIM)] |
173 |
+ DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*B_p[INDEX2(2,q,DIM)]); |
174 |
} |
175 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
176 |
} |
177 |
} |
178 |
} else { |
179 |
for (s=0;s<p.row_NS;s++) { |
180 |
for (r=0;r<p.col_NS;r++) { |
181 |
rtmp0=0; |
182 |
rtmp1=0; |
183 |
rtmp2=0; |
184 |
for (q=0;q<p.numQuad;q++) { |
185 |
rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)]; |
186 |
rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]; |
187 |
rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)]; |
188 |
rtmp2+=rtmp*DSDX[INDEX3(s,2,q,p.row_NN,DIM)]; |
189 |
} |
190 |
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]; |
191 |
} |
192 |
} |
193 |
} |
194 |
} |
195 |
/**************************************************************/ |
196 |
/* process C: */ |
197 |
/**************************************************************/ |
198 |
C_p=getSampleData(C,e); |
199 |
if (NULL!=C_p) { |
200 |
add_EM_S=TRUE; |
201 |
if (extendedC) { |
202 |
for (s=0;s<p.row_NS;s++) { |
203 |
for (r=0;r<p.col_NS;r++) { |
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[INDEX2(0,q,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)] |
208 |
+ C_p[INDEX2(1,q,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)] |
209 |
+ C_p[INDEX2(2,q,DIM)]*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]); |
210 |
} |
211 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
212 |
} |
213 |
} |
214 |
} else { |
215 |
for (s=0;s<p.row_NS;s++) { |
216 |
for (r=0;r<p.col_NS;r++) { |
217 |
rtmp0=0; |
218 |
rtmp1=0; |
219 |
rtmp2=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 |
rtmp2+=rtmp*DSDX[INDEX3(r,2,q,p.row_NN,DIM)]; |
225 |
} |
226 |
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]; |
227 |
} |
228 |
} |
229 |
} |
230 |
} |
231 |
/************************************************************* */ |
232 |
/* process D */ |
233 |
/**************************************************************/ |
234 |
D_p=getSampleData(D,e); |
235 |
if (NULL!=D_p) { |
236 |
add_EM_S=TRUE; |
237 |
if (extendedD) { |
238 |
for (s=0;s<p.row_NS;s++) { |
239 |
for (r=0;r<p.col_NS;r++) { |
240 |
rtmp=0; |
241 |
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)]; |
242 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
243 |
} |
244 |
} |
245 |
} else { |
246 |
for (s=0;s<p.row_NS;s++) { |
247 |
for (r=0;r<p.col_NS;r++) { |
248 |
rtmp=0; |
249 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)]; |
250 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[0]; |
251 |
} |
252 |
} |
253 |
} |
254 |
} |
255 |
/**************************************************************/ |
256 |
/* process X: */ |
257 |
/**************************************************************/ |
258 |
X_p=getSampleData(X,e); |
259 |
if (NULL!=X_p) { |
260 |
add_EM_F=TRUE; |
261 |
if (extendedX) { |
262 |
for (s=0;s<p.row_NS;s++) { |
263 |
rtmp=0; |
264 |
for (q=0;q<p.numQuad;q++) { |
265 |
rtmp+=Vol[q]*( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX2(0,q,DIM)] |
266 |
+ DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX2(1,q,DIM)] |
267 |
+ DSDX[INDEX3(s,2,q,p.row_NN,DIM)]*X_p[INDEX2(2,q,DIM)]); |
268 |
} |
269 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp; |
270 |
} |
271 |
} else { |
272 |
for (s=0;s<p.row_NS;s++) { |
273 |
rtmp0=0; |
274 |
rtmp1=0; |
275 |
rtmp2=0; |
276 |
for (q=0;q<p.numQuad;q++) { |
277 |
rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]; |
278 |
rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)]; |
279 |
rtmp2+=Vol[q]*DSDX[INDEX3(s,2,q,p.row_NN,DIM)]; |
280 |
} |
281 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp0*X_p[0]+rtmp1*X_p[1]+rtmp2*X_p[2]; |
282 |
} |
283 |
} |
284 |
} |
285 |
/**************************************************************/ |
286 |
/* process Y: */ |
287 |
/**************************************************************/ |
288 |
Y_p=getSampleData(Y,e); |
289 |
if (NULL!=Y_p) { |
290 |
add_EM_F=TRUE; |
291 |
if (extendedY) { |
292 |
for (s=0;s<p.row_NS;s++) { |
293 |
rtmp=0; |
294 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[q]; |
295 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp; |
296 |
} |
297 |
} else { |
298 |
for (s=0;s<p.row_NS;s++) { |
299 |
rtmp=0; |
300 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]; |
301 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*Y_p[0]; |
302 |
} |
303 |
} |
304 |
} |
305 |
/***********************************************************************************************/ |
306 |
/* add the element matrices onto the matrix and right hand side */ |
307 |
/***********************************************************************************************/ |
308 |
for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]]; |
309 |
if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound); |
310 |
if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S); |
311 |
|
312 |
} /* end color check */ |
313 |
} /* end element loop */ |
314 |
} /* end color loop */ |
315 |
|
316 |
THREAD_MEMFREE(EM_S); |
317 |
THREAD_MEMFREE(EM_F); |
318 |
THREAD_MEMFREE(row_index); |
319 |
|
320 |
} /* end of pointer check */ |
321 |
} /* end parallel region */ |
322 |
} |
323 |
/* |
324 |
* $Log$ |
325 |
*/ |