1 |
|
2 |
/******************************************************* |
3 |
* |
4 |
* Copyright (c) 2003-2008 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 1D 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 1 x p.numComp x 1 */ |
29 |
/* B = 1 x numEqu x p.numComp */ |
30 |
/* C = p.numEqu x 1 x p.numComp */ |
31 |
/* D = p.numEqu x p.numComp */ |
32 |
/* X = p.numEqu x 1 */ |
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_1D(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 1 |
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; |
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,add_EM_F, add_EM_S) |
80 |
{ |
81 |
EM_S=THREAD_MEMALLOC(len_EM_S,double); |
82 |
EM_F=THREAD_MEMALLOC(len_EM_F,double); |
83 |
row_index=THREAD_MEMALLOC(p.row_NN,index_t); |
84 |
|
85 |
if (!Finley_checkPtr(EM_S) && !Finley_checkPtr(EM_F) && !Finley_checkPtr(row_index) ) { |
86 |
|
87 |
for (color=elements->minColor;color<=elements->maxColor;color++) { |
88 |
/* open loop over all elements: */ |
89 |
#pragma omp for private(e) schedule(static) |
90 |
for(e=0;e<elements->numElements;e++){ |
91 |
if (elements->Color[e]==color) { |
92 |
Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]); |
93 |
DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]); |
94 |
for (q=0;q<len_EM_S;++q) EM_S[q]=0; |
95 |
for (q=0;q<len_EM_F;++q) EM_F[q]=0; |
96 |
add_EM_F=FALSE; |
97 |
add_EM_S=FALSE; |
98 |
/**************************************************************/ |
99 |
/* process A: */ |
100 |
/**************************************************************/ |
101 |
A_p=getSampleDataRO(A,e,ABuff); |
102 |
if (NULL!=A_p) { |
103 |
add_EM_S=TRUE; |
104 |
if (extendedA) { |
105 |
for (s=0;s<p.row_NS;s++) { |
106 |
for (r=0;r<p.col_NS;r++) { |
107 |
for (k=0;k<p.numEqu;k++) { |
108 |
for (m=0;m<p.numComp;m++) { |
109 |
rtmp=0.; |
110 |
for (q=0;q<p.numQuad;q++) { |
111 |
rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]* |
112 |
A_p[INDEX5(k,0,m,0,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
113 |
} |
114 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
115 |
} |
116 |
} |
117 |
} |
118 |
} |
119 |
} else { |
120 |
for (s=0;s<p.row_NS;s++) { |
121 |
for (r=0;r<p.col_NS;r++) { |
122 |
rtmp=0; |
123 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
124 |
for (k=0;k<p.numEqu;k++) { |
125 |
for (m=0;m<p.numComp;m++) { |
126 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*A_p[INDEX4(k,0,m,0,p.numEqu,DIM,p.numComp)]; |
127 |
} |
128 |
} |
129 |
} |
130 |
} |
131 |
} |
132 |
} |
133 |
/**************************************************************/ |
134 |
/* process B: */ |
135 |
/**************************************************************/ |
136 |
B_p=getSampleDataRO(B,e,BBuff); |
137 |
if (NULL!=B_p) { |
138 |
add_EM_S=TRUE; |
139 |
if (extendedB) { |
140 |
for (s=0;s<p.row_NS;s++) { |
141 |
for (r=0;r<p.col_NS;r++) { |
142 |
for (k=0;k<p.numEqu;k++) { |
143 |
for (m=0;m<p.numComp;m++) { |
144 |
rtmp=0.; |
145 |
for (q=0;q<p.numQuad;q++) { |
146 |
rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]* |
147 |
B_p[INDEX4(k,0,m,q,p.numEqu,DIM,p.numComp)]*S[INDEX2(r,q,p.row_NS)]; |
148 |
} |
149 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
150 |
} |
151 |
} |
152 |
} |
153 |
} |
154 |
} else { |
155 |
for (s=0;s<p.row_NS;s++) { |
156 |
for (r=0;r<p.col_NS;r++) { |
157 |
rtmp=0; |
158 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*S[INDEX2(r,q,p.row_NS)]; |
159 |
for (k=0;k<p.numEqu;k++) { |
160 |
for (m=0;m<p.numComp;m++) { |
161 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*B_p[INDEX3(k,0,m,p.numEqu,DIM)]; |
162 |
} |
163 |
} |
164 |
} |
165 |
} |
166 |
} |
167 |
} |
168 |
/**************************************************************/ |
169 |
/* process C: */ |
170 |
/**************************************************************/ |
171 |
C_p=getSampleDataRO(C,e,CBuff); |
172 |
if (NULL!=C_p) { |
173 |
add_EM_S=TRUE; |
174 |
if (extendedC) { |
175 |
for (s=0;s<p.row_NS;s++) { |
176 |
for (r=0;r<p.col_NS;r++) { |
177 |
for (k=0;k<p.numEqu;k++) { |
178 |
for (m=0;m<p.numComp;m++) { |
179 |
rtmp=0; |
180 |
for (q=0;q<p.numQuad;q++) { |
181 |
rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*C_p[INDEX4(k,m,0,q,p.numEqu,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
182 |
} |
183 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
184 |
} |
185 |
} |
186 |
} |
187 |
} |
188 |
} else { |
189 |
for (s=0;s<p.row_NS;s++) { |
190 |
for (r=0;r<p.col_NS;r++) { |
191 |
rtmp=0; |
192 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
193 |
for (k=0;k<p.numEqu;k++) { |
194 |
for (m=0;m<p.numComp;m++) { |
195 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*C_p[INDEX3(k,m,0,p.numEqu,p.numComp)]; |
196 |
} |
197 |
} |
198 |
} |
199 |
} |
200 |
} |
201 |
} |
202 |
/************************************************************* */ |
203 |
/* process D */ |
204 |
/**************************************************************/ |
205 |
D_p=getSampleDataRO(D,e,DBuff); |
206 |
if (NULL!=D_p) { |
207 |
add_EM_S=TRUE; |
208 |
if (extendedD) { |
209 |
for (s=0;s<p.row_NS;s++) { |
210 |
for (r=0;r<p.col_NS;r++) { |
211 |
for (k=0;k<p.numEqu;k++) { |
212 |
for (m=0;m<p.numComp;m++) { |
213 |
rtmp=0; |
214 |
for (q=0;q<p.numQuad;q++) { |
215 |
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)]; |
216 |
} |
217 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
218 |
|
219 |
} |
220 |
} |
221 |
} |
222 |
} |
223 |
} else { |
224 |
for (s=0;s<p.row_NS;s++) { |
225 |
for (r=0;r<p.col_NS;r++) { |
226 |
rtmp=0; |
227 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)]; |
228 |
for (k=0;k<p.numEqu;k++) { |
229 |
for (m=0;m<p.numComp;m++) { |
230 |
EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[INDEX2(k,m,p.numEqu)]; |
231 |
} |
232 |
} |
233 |
} |
234 |
} |
235 |
} |
236 |
} |
237 |
/**************************************************************/ |
238 |
/* process X: */ |
239 |
/**************************************************************/ |
240 |
X_p=getSampleDataRO(X,e,XBuff); |
241 |
if (NULL!=X_p) { |
242 |
add_EM_F=TRUE; |
243 |
if (extendedX) { |
244 |
for (s=0;s<p.row_NS;s++) { |
245 |
for (k=0;k<p.numEqu;k++) { |
246 |
rtmp=0; |
247 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX3(k,0,q,p.numEqu,DIM)]; |
248 |
EM_F[INDEX2(k,s,p.numEqu)]+=rtmp; |
249 |
} |
250 |
} |
251 |
} else { |
252 |
for (s=0;s<p.row_NS;s++) { |
253 |
rtmp=0; |
254 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]; |
255 |
for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*X_p[INDEX2(k,0,p.numEqu)]; |
256 |
} |
257 |
} |
258 |
} |
259 |
/**************************************************************/ |
260 |
/* process Y: */ |
261 |
/**************************************************************/ |
262 |
Y_p=getSampleDataRO(Y,e,YBuff); |
263 |
if (NULL!=Y_p) { |
264 |
add_EM_F=TRUE; |
265 |
if (extendedY) { |
266 |
for (s=0;s<p.row_NS;s++) { |
267 |
for (k=0;k<p.numEqu;k++) { |
268 |
rtmp=0; |
269 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[INDEX2(k,q,p.numEqu)]; |
270 |
EM_F[INDEX2(k,s,p.numEqu)]+=rtmp; |
271 |
} |
272 |
} |
273 |
} else { |
274 |
for (s=0;s<p.row_NS;s++) { |
275 |
rtmp=0; |
276 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]; |
277 |
for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*Y_p[k]; |
278 |
} |
279 |
} |
280 |
} |
281 |
/***********************************************************************************************/ |
282 |
/* add the element matrices onto the matrix and right hand side */ |
283 |
/***********************************************************************************************/ |
284 |
for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]]; |
285 |
if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound); |
286 |
if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S); |
287 |
|
288 |
} /* end color check */ |
289 |
} /* end element loop */ |
290 |
} /* end color loop */ |
291 |
|
292 |
THREAD_MEMFREE(EM_S); |
293 |
THREAD_MEMFREE(EM_F); |
294 |
THREAD_MEMFREE(row_index); |
295 |
|
296 |
} /* end of pointer check */ |
297 |
} /* end parallel region */ |
298 |
freeSampleBuffer(ABuff); |
299 |
freeSampleBuffer(BBuff); |
300 |
freeSampleBuffer(CBuff); |
301 |
freeSampleBuffer(DBuff); |
302 |
freeSampleBuffer(XBuff); |
303 |
freeSampleBuffer(YBuff); |
304 |
} |
305 |
/* |
306 |
* $Log$ |
307 |
*/ |