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_{i,j} u_,j)_i-(B_{i} u)_i+C_{j} u_,j-D u_m and -(X_,i)_i + Y */ |
22 |
|
23 |
/* 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 = 1 x 1 */ |
29 |
/* B = 1 */ |
30 |
/* C = 1 */ |
31 |
/* D = scalar */ |
32 |
/* X = 1 */ |
33 |
/* Y = scalar */ |
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_Single2_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 |
double *EM_S, *EM_F, *Vol, *DSDX, *A_p, *B_p, *C_p, *D_p, *X_p, *Y_p; |
56 |
index_t *row_index; |
57 |
register dim_t q, s,r; |
58 |
register double rtmp; |
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=getSampleData(F,0); |
68 |
double *S=p.row_jac->ReferenceElement->S; |
69 |
dim_t len_EM_S=p.row_NN*p.col_NN; |
70 |
dim_t len_EM_F=p.row_NN; |
71 |
|
72 |
|
73 |
#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,add_EM_F, add_EM_S) |
74 |
{ |
75 |
EM_S=THREAD_MEMALLOC(len_EM_S,double); |
76 |
EM_F=THREAD_MEMALLOC(len_EM_F,double); |
77 |
row_index=THREAD_MEMALLOC(p.row_NN,index_t); |
78 |
|
79 |
if (!Finley_checkPtr(EM_S) && !Finley_checkPtr(EM_F) && !Finley_checkPtr(row_index) ) { |
80 |
|
81 |
for (color=elements->minColor;color<=elements->maxColor;color++) { |
82 |
/* open loop over all elements: */ |
83 |
#pragma omp for private(e) schedule(static) |
84 |
for(e=0;e<elements->numElements;e++){ |
85 |
if (elements->Color[e]==color) { |
86 |
Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]); |
87 |
DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]); |
88 |
for (q=0;q<len_EM_S;++q) EM_S[q]=0; |
89 |
for (q=0;q<len_EM_F;++q) EM_F[q]=0; |
90 |
add_EM_F=FALSE; |
91 |
add_EM_S=FALSE; |
92 |
/**************************************************************/ |
93 |
/* process A: */ |
94 |
/**************************************************************/ |
95 |
A_p=getSampleData(A,e); |
96 |
if (NULL!=A_p) { |
97 |
add_EM_S=TRUE; |
98 |
if (extendedA) { |
99 |
for (s=0;s<p.row_NS;s++) { |
100 |
for (r=0;r<p.col_NS;r++) { |
101 |
rtmp=0; |
102 |
for (q=0;q<p.numQuad;q++) { |
103 |
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)]; |
104 |
} |
105 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
106 |
} |
107 |
} |
108 |
} else { |
109 |
for (s=0;s<p.row_NS;s++) { |
110 |
for (r=0;r<p.col_NS;r++) { |
111 |
rtmp=0; |
112 |
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)]; |
113 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*A_p[INDEX2(0,0,DIM)]; |
114 |
} |
115 |
} |
116 |
} |
117 |
} |
118 |
/**************************************************************/ |
119 |
/* process B: */ |
120 |
/**************************************************************/ |
121 |
B_p=getSampleData(B,e); |
122 |
if (NULL!=B_p) { |
123 |
add_EM_S=TRUE; |
124 |
if (extendedB) { |
125 |
for (s=0;s<p.row_NS;s++) { |
126 |
for (r=0;r<p.col_NS;r++) { |
127 |
rtmp=0; |
128 |
for (q=0;q<p.numQuad;q++) { |
129 |
rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX2(0,q,DIM)]*S[INDEX2(r,q,p.row_NS)]; |
130 |
} |
131 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
132 |
} |
133 |
} |
134 |
} else { |
135 |
for (s=0;s<p.row_NS;s++) { |
136 |
for (r=0;r<p.col_NS;r++) { |
137 |
rtmp=0; |
138 |
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)]; |
139 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*B_p[0]; |
140 |
} |
141 |
} |
142 |
} |
143 |
} |
144 |
/**************************************************************/ |
145 |
/* process C: */ |
146 |
/**************************************************************/ |
147 |
C_p=getSampleData(C,e); |
148 |
if (NULL!=C_p) { |
149 |
add_EM_S=TRUE; |
150 |
if (extendedC) { |
151 |
for (s=0;s<p.row_NS;s++) { |
152 |
for (r=0;r<p.col_NS;r++) { |
153 |
rtmp=0; |
154 |
for (q=0;q<p.numQuad;q++) { |
155 |
rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*C_p[INDEX2(0,q,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)]; |
156 |
} |
157 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
158 |
} |
159 |
} |
160 |
} else { |
161 |
for (s=0;s<p.row_NS;s++) { |
162 |
for (r=0;r<p.col_NS;r++) { |
163 |
rtmp=0; |
164 |
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)]; |
165 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*C_p[0]; |
166 |
} |
167 |
} |
168 |
} |
169 |
} |
170 |
/************************************************************* */ |
171 |
/* process D */ |
172 |
/**************************************************************/ |
173 |
D_p=getSampleData(D,e); |
174 |
if (NULL!=D_p) { |
175 |
add_EM_S=TRUE; |
176 |
if (extendedD) { |
177 |
for (s=0;s<p.row_NS;s++) { |
178 |
for (r=0;r<p.col_NS;r++) { |
179 |
rtmp=0; |
180 |
for (q=0;q<p.numQuad;q++) { |
181 |
rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*D_p[q]*S[INDEX2(r,q,p.row_NS)]; |
182 |
} |
183 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp; |
184 |
} |
185 |
} |
186 |
} else { |
187 |
for (s=0;s<p.row_NS;s++) { |
188 |
for (r=0;r<p.col_NS;r++) { |
189 |
rtmp=0; |
190 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)]; |
191 |
EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[0]; |
192 |
} |
193 |
} |
194 |
} |
195 |
} |
196 |
/**************************************************************/ |
197 |
/* process X: */ |
198 |
/**************************************************************/ |
199 |
X_p=getSampleData(X,e); |
200 |
if (NULL!=X_p) { |
201 |
add_EM_F=TRUE; |
202 |
if (extendedX) { |
203 |
for (s=0;s<p.row_NS;s++) { |
204 |
rtmp=0; |
205 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX2(0,q,DIM)]; |
206 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp; |
207 |
} |
208 |
} else { |
209 |
for (s=0;s<p.row_NS;s++) { |
210 |
rtmp=0; |
211 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]; |
212 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*X_p[0]; |
213 |
} |
214 |
} |
215 |
} |
216 |
/**************************************************************/ |
217 |
/* process Y: */ |
218 |
/**************************************************************/ |
219 |
Y_p=getSampleData(Y,e); |
220 |
if (NULL!=Y_p) { |
221 |
add_EM_F=TRUE; |
222 |
if (extendedY) { |
223 |
for (s=0;s<p.row_NS;s++) { |
224 |
rtmp=0; |
225 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[q]; |
226 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp; |
227 |
} |
228 |
} else { |
229 |
for (s=0;s<p.row_NS;s++) { |
230 |
rtmp=0; |
231 |
for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]; |
232 |
EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*Y_p[0]; |
233 |
} |
234 |
} |
235 |
} |
236 |
/***********************************************************************************************/ |
237 |
/* add the element matrices onto the matrix and right hand side */ |
238 |
/***********************************************************************************************/ |
239 |
for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]]; |
240 |
if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound); |
241 |
if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S); |
242 |
|
243 |
} /* end color check */ |
244 |
} /* end element loop */ |
245 |
} /* end color loop */ |
246 |
|
247 |
THREAD_MEMFREE(EM_S); |
248 |
THREAD_MEMFREE(EM_F); |
249 |
THREAD_MEMFREE(row_index); |
250 |
|
251 |
} /* end of pointer check */ |
252 |
} /* end parallel region */ |
253 |
} |
254 |
/* |
255 |
* $Log$ |
256 |
*/ |