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

Contents of /trunk/finley/src/Assemble_PDE_System2_1D.c

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.26