/[escript]/trunk-mpi-branch/finley/src/Assemble_PDE_System2_1D.c
ViewVC logotype

Contents of /trunk-mpi-branch/finley/src/Assemble_PDE_System2_1D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1253 - (show annotations)
Fri Aug 17 04:09:29 2007 UTC (12 years, 9 months ago) by gross
File MIME type: text/plain
File size: 14469 byte(s)
finally everthing for MPI is in the code. now testing is needed
1 /*
2 ************************************************************
3 * Copyright 2006 by ACcESS MNRF *
4 * *
5 * http://www.access.edu.au *
6 * Primary Business: Queensland, Australia *
7 * Licensed under the Open Software License version 3.0 *
8 * http://www.opensource.org/licenses/osl-3.0.php *
9 * *
10 ************************************************************
11 */
12
13 /**************************************************************/
14
15 /* assembles the system of numEq PDEs into the stiffness matrix S right hand side F */
16 /* the shape functions for test and solution must be identical */
17
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 1D 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 1 x p.numComp x 1 */
27 /* B = 1 x numEqu x p.numComp */
28 /* C = p.numEqu x 1 x p.numComp */
29 /* D = p.numEqu x p.numComp */
30 /* X = p.numEqu x 1 */
31 /* Y = p.numEqu */
32
33
34 /**************************************************************/
35
36 /* Author: gross@access.edu.au */
37 /* Version: $Id:$ */
38
39 /**************************************************************/
40
41
42 #include "Assemble.h"
43 #include "Util.h"
44 #ifdef _OPENMP
45 #include <omp.h>
46 #endif
47
48
49 /**************************************************************/
50
51 void Finley_Assemble_PDE_System2_1D(Assemble_Parameters p, Finley_ElementFile* elements,
52 Paso_SystemMatrix* Mat, escriptDataC* F,
53 escriptDataC* A, escriptDataC* B, escriptDataC* C, escriptDataC* D, escriptDataC* X, escriptDataC* Y) {
54
55 #define DIM 1
56 index_t color;
57 dim_t e;
58 double *EM_S, *EM_F, *Vol, *DSDX, *A_p, *B_p, *C_p, *D_p, *X_p, *Y_p;
59 index_t *row_index;
60 register dim_t q, s,r,k,m;
61 register double rtmp;
62 bool_t add_EM_F, add_EM_S;
63
64 bool_t extendedA=isExpanded(A);
65 bool_t extendedB=isExpanded(B);
66 bool_t extendedC=isExpanded(C);
67 bool_t extendedD=isExpanded(D);
68 bool_t extendedX=isExpanded(X);
69 bool_t extendedY=isExpanded(Y);
70 double *F_p=getSampleData(F,0);
71 double *S=p.row_jac->ReferenceElement->S;
72 dim_t len_EM_S=p.row_NN*p.col_NN*p.numEqu*p.numComp;
73 dim_t len_EM_F=p.row_NN*p.numEqu;
74
75
76 #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)
77 {
78 EM_S=THREAD_MEMALLOC(len_EM_S,double);
79 EM_F=THREAD_MEMALLOC(len_EM_F,double);
80 row_index=THREAD_MEMALLOC(p.row_NN,index_t);
81
82 if (!Finley_checkPtr(EM_S) && !Finley_checkPtr(EM_F) && !Finley_checkPtr(row_index) ) {
83
84 for (color=elements->minColor;color<=elements->maxColor;color++) {
85 /* open loop over all elements: */
86 #pragma omp for private(e) schedule(static)
87 for(e=0;e<elements->numElements;e++){
88 if (elements->Color[e]==color) {
89 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
90 DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]);
91 for (q=0;q<len_EM_S;++q) EM_S[q]=0;
92 for (q=0;q<len_EM_F;++q) EM_F[q]=0;
93 add_EM_F=FALSE;
94 add_EM_S=FALSE;
95 /**************************************************************/
96 /* process A: */
97 /**************************************************************/
98 A_p=getSampleData(A,e);
99 if (NULL!=A_p) {
100 add_EM_S=TRUE;
101 if (extendedA) {
102 for (s=0;s<p.row_NS;s++) {
103 for (r=0;r<p.col_NS;r++) {
104 for (k=0;k<p.numEqu;k++) {
105 for (m=0;m<p.numComp;m++) {
106 rtmp=0.;
107 for (q=0;q<p.numQuad;q++) {
108 rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*
109 A_p[INDEX5(k,0,m,0,q,p.numEqu,DIM,p.numComp,DIM)]*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
110 }
111 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
112 }
113 }
114 }
115 }
116 } else {
117 for (s=0;s<p.row_NS;s++) {
118 for (r=0;r<p.col_NS;r++) {
119 rtmp=0;
120 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)];
121 for (k=0;k<p.numEqu;k++) {
122 for (m=0;m<p.numComp;m++) {
123 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)];
124 }
125 }
126 }
127 }
128 }
129 }
130 /**************************************************************/
131 /* process B: */
132 /**************************************************************/
133 B_p=getSampleData(B,e);
134 if (NULL!=B_p) {
135 add_EM_S=TRUE;
136 if (extendedB) {
137 for (s=0;s<p.row_NS;s++) {
138 for (r=0;r<p.col_NS;r++) {
139 for (k=0;k<p.numEqu;k++) {
140 for (m=0;m<p.numComp;m++) {
141 rtmp=0.;
142 for (q=0;q<p.numQuad;q++) {
143 rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*
144 B_p[INDEX4(k,0,m,q,p.numEqu,DIM,p.numComp)]*S[INDEX2(r,q,p.row_NS)];
145 }
146 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
147 }
148 }
149 }
150 }
151 } else {
152 for (s=0;s<p.row_NS;s++) {
153 for (r=0;r<p.col_NS;r++) {
154 rtmp=0;
155 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)];
156 for (k=0;k<p.numEqu;k++) {
157 for (m=0;m<p.numComp;m++) {
158 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*B_p[INDEX3(k,0,m,p.numEqu,DIM)];
159 }
160 }
161 }
162 }
163 }
164 }
165 /**************************************************************/
166 /* process C: */
167 /**************************************************************/
168 C_p=getSampleData(C,e);
169 if (NULL!=C_p) {
170 add_EM_S=TRUE;
171 if (extendedC) {
172 for (s=0;s<p.row_NS;s++) {
173 for (r=0;r<p.col_NS;r++) {
174 for (k=0;k<p.numEqu;k++) {
175 for (m=0;m<p.numComp;m++) {
176 rtmp=0;
177 for (q=0;q<p.numQuad;q++) {
178 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)];
179 }
180 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
181 }
182 }
183 }
184 }
185 } else {
186 for (s=0;s<p.row_NS;s++) {
187 for (r=0;r<p.col_NS;r++) {
188 rtmp=0;
189 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)];
190 for (k=0;k<p.numEqu;k++) {
191 for (m=0;m<p.numComp;m++) {
192 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)];
193 }
194 }
195 }
196 }
197 }
198 }
199 /************************************************************* */
200 /* process D */
201 /**************************************************************/
202 D_p=getSampleData(D,e);
203 if (NULL!=D_p) {
204 add_EM_S=TRUE;
205 if (extendedD) {
206 for (s=0;s<p.row_NS;s++) {
207 for (r=0;r<p.col_NS;r++) {
208 for (k=0;k<p.numEqu;k++) {
209 for (m=0;m<p.numComp;m++) {
210 rtmp=0;
211 for (q=0;q<p.numQuad;q++) {
212 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)];
213 }
214 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
215
216 }
217 }
218 }
219 }
220 } else {
221 for (s=0;s<p.row_NS;s++) {
222 for (r=0;r<p.col_NS;r++) {
223 rtmp=0;
224 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
225 for (k=0;k<p.numEqu;k++) {
226 for (m=0;m<p.numComp;m++) {
227 EM_S[INDEX4(k,m,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[INDEX2(k,m,p.numEqu)];
228 }
229 }
230 }
231 }
232 }
233 }
234 /**************************************************************/
235 /* process X: */
236 /**************************************************************/
237 X_p=getSampleData(X,e);
238 if (NULL!=X_p) {
239 add_EM_F=TRUE;
240 if (extendedX) {
241 for (s=0;s<p.row_NS;s++) {
242 for (k=0;k<p.numEqu;k++) {
243 rtmp=0;
244 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)];
245 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
246 }
247 }
248 } else {
249 for (s=0;s<p.row_NS;s++) {
250 rtmp=0;
251 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
252 for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*X_p[INDEX2(k,0,p.numEqu)];
253 }
254 }
255 }
256 /**************************************************************/
257 /* process Y: */
258 /**************************************************************/
259 Y_p=getSampleData(Y,e);
260 if (NULL!=Y_p) {
261 add_EM_F=TRUE;
262 if (extendedY) {
263 for (s=0;s<p.row_NS;s++) {
264 for (k=0;k<p.numEqu;k++) {
265 rtmp=0;
266 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[INDEX2(k,q,p.numEqu)];
267 EM_F[INDEX2(k,s,p.numEqu)]+=rtmp;
268 }
269 }
270 } else {
271 for (s=0;s<p.row_NS;s++) {
272 rtmp=0;
273 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
274 for (k=0;k<p.numEqu;k++) EM_F[INDEX2(k,s,p.numEqu)]+=rtmp*Y_p[k];
275 }
276 }
277 }
278 /***********************************************************************************************/
279 /* add the element matrices onto the matrix and right hand side */
280 /***********************************************************************************************/
281 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
282 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
283 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
284
285 } /* end color check */
286 } /* end element loop */
287 } /* end color loop */
288
289 THREAD_MEMFREE(EM_S);
290 THREAD_MEMFREE(EM_F);
291 THREAD_MEMFREE(row_index);
292
293 } /* end of pointer check */
294 } /* end parallel region */
295 }
296 /*
297 * $Log$
298 */

  ViewVC Help
Powered by ViewVC 1.1.26