/[escript]/branches/arrexp_2137_win_merge/finley/src/Assemble_PDE_Single2_2D.c
ViewVC logotype

Contents of /branches/arrexp_2137_win_merge/finley/src/Assemble_PDE_Single2_2D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2222 - (show annotations)
Tue Jan 20 04:52:39 2009 UTC (11 years, 1 month ago) by jfenwick
File MIME type: text/plain
File size: 14901 byte(s)
Saving work
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 2D 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 = 2 x 2 */
29 /* B = 2 */
30 /* C = 2 */
31 /* D = scalar */
32 /* X = 2 */
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 void Finley_Assemble_PDE_Single2_2D(Assemble_Parameters p, Finley_ElementFile* elements,
48 Paso_SystemMatrix* Mat, escriptDataC* F,
49 escriptDataC* A, escriptDataC* B, escriptDataC* C, escriptDataC* D, escriptDataC* X, escriptDataC* Y) {
50
51 #define DIM 2
52 index_t color;
53 dim_t e;
54 __const double *A_p, *B_p, *C_p, *D_p, *X_p, *Y_p;
55 double *EM_S, *EM_F, *Vol, *DSDX;
56 index_t *row_index;
57 register dim_t q, s,r;
58 register double rtmp00, rtmp01, rtmp10, rtmp11, rtmp, rtmp0, rtmp1;
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=(requireWrite(F), getSampleDataRW(F,0)); /* use comma, to get around the mixed code and declarations thing */
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,rtmp00, rtmp01, rtmp10, rtmp11, rtmp, rtmp0, rtmp1,add_EM_F, add_EM_S)
74 {
75 void* ABuff=allocSampleBuffer(A);
76 void* BBuff=allocSampleBuffer(B);
77 void* CBuff=allocSampleBuffer(C);
78 void* DBuff=allocSampleBuffer(D);
79 void* XBuff=allocSampleBuffer(X);
80 void* YBuff=allocSampleBuffer(Y);
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 rtmp=0;
108 for (q=0;q<p.numQuad;q++) {
109 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)]
110 + 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)]
111 + 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)]
112 + 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)] );
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 rtmp10=0;
123 rtmp11=0;
124 for (q=0;q<p.numQuad;q++) {
125 rtmp0=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
126 rtmp1=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
127 rtmp00+=rtmp0*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
128 rtmp01+=rtmp0*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
129 rtmp10+=rtmp1*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
130 rtmp11+=rtmp1*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
131 }
132 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+= rtmp00*A_p[INDEX2(0,0,DIM)]
133 + rtmp01*A_p[INDEX2(0,1,DIM)]
134 + rtmp10*A_p[INDEX2(1,0,DIM)]
135 + rtmp11*A_p[INDEX2(1,1,DIM)];
136 }
137 }
138 }
139 }
140 /**************************************************************/
141 /* process B: */
142 /**************************************************************/
143 B_p=getSampleDataRO(B,e,BBuff);
144 if (NULL!=B_p) {
145 add_EM_S=TRUE;
146 if (extendedB) {
147 for (s=0;s<p.row_NS;s++) {
148 for (r=0;r<p.col_NS;r++) {
149 rtmp=0.;
150 for (q=0;q<p.numQuad;q++) {
151 rtmp+=Vol[q]*S[INDEX2(r,q,p.row_NS)]*(DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX2(0,q,DIM)]
152 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*B_p[INDEX2(1,q,DIM)]);
153 }
154 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
155 }
156 }
157 } else {
158 for (s=0;s<p.row_NS;s++) {
159 for (r=0;r<p.col_NS;r++) {
160 rtmp0=0;
161 rtmp1=0;
162 for (q=0;q<p.numQuad;q++) {
163 rtmp=Vol[q]*S[INDEX2(r,q,p.row_NS)];
164 rtmp0+=rtmp*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
165 rtmp1+=rtmp*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
166 }
167 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp0*B_p[0]+rtmp1*B_p[1];
168 }
169 }
170 }
171 }
172 /**************************************************************/
173 /* process C: */
174 /**************************************************************/
175 C_p=getSampleDataRO(C,e,CBuff);
176 if (NULL!=C_p) {
177 add_EM_S=TRUE;
178 if (extendedC) {
179 for (s=0;s<p.row_NS;s++) {
180 for (r=0;r<p.col_NS;r++) {
181 rtmp=0;
182 for (q=0;q<p.numQuad;q++) {
183 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)]
184 + C_p[INDEX2(1,q,DIM)]*DSDX[INDEX3(r,1,q,p.row_NN,DIM)]);
185 }
186 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
187 }
188 }
189 } else {
190 for (s=0;s<p.row_NS;s++) {
191 for (r=0;r<p.col_NS;r++) {
192 rtmp0=0;
193 rtmp1=0;
194 for (q=0;q<p.numQuad;q++) {
195 rtmp=Vol[q]*S[INDEX2(s,q,p.row_NS)];
196 rtmp0+=rtmp*DSDX[INDEX3(r,0,q,p.row_NN,DIM)];
197 rtmp1+=rtmp*DSDX[INDEX3(r,1,q,p.row_NN,DIM)];
198 }
199 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp0*C_p[0]+rtmp1*C_p[1];
200 }
201 }
202 }
203 }
204 /************************************************************* */
205 /* process D */
206 /**************************************************************/
207 D_p=getSampleDataRO(D,e,DBuff);
208 if (NULL!=D_p) {
209 add_EM_S=TRUE;
210 if (extendedD) {
211 for (s=0;s<p.row_NS;s++) {
212 for (r=0;r<p.col_NS;r++) {
213 rtmp=0;
214 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)];
215 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
216 }
217 }
218 } else {
219 for (s=0;s<p.row_NS;s++) {
220 for (r=0;r<p.col_NS;r++) {
221 rtmp=0;
222 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
223 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[0];
224 }
225 }
226 }
227 }
228 /**************************************************************/
229 /* process X: */
230 /**************************************************************/
231 X_p=getSampleDataRO(X,e,XBuff);
232 if (NULL!=X_p) {
233 add_EM_F=TRUE;
234 if (extendedX) {
235 for (s=0;s<p.row_NS;s++) {
236 rtmp=0.;
237 for (q=0;q<p.numQuad;q++) {
238 rtmp+=Vol[q]*( DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX2(0,q,DIM)]
239 + DSDX[INDEX3(s,1,q,p.row_NN,DIM)]*X_p[INDEX2(1,q,DIM)]);
240 }
241 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
242 }
243 } else {
244 for (s=0;s<p.row_NS;s++) {
245 rtmp0=0.;
246 rtmp1=0.;
247 for (q=0;q<p.numQuad;q++) {
248 rtmp0+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
249 rtmp1+=Vol[q]*DSDX[INDEX3(s,1,q,p.row_NN,DIM)];
250 }
251 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp0*X_p[0]+rtmp1*X_p[1];
252 }
253 }
254 }
255 /**************************************************************/
256 /* process Y: */
257 /**************************************************************/
258 Y_p=getSampleDataRO(Y,e,YBuff);
259 if (NULL!=Y_p) {
260 add_EM_F=TRUE;
261 if (extendedY) {
262 for (s=0;s<p.row_NS;s++) {
263 rtmp=0;
264 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[q];
265 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
266 }
267 } else {
268 for (s=0;s<p.row_NS;s++) {
269 rtmp=0;
270 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
271 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*Y_p[0];
272 }
273 }
274 }
275 /***********************************************************************************************/
276 /* add the element matrices onto the matrix and right hand side */
277 /***********************************************************************************************/
278 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
279 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
280 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
281
282 } /* end color check */
283 } /* end element loop */
284 } /* end color loop */
285
286 THREAD_MEMFREE(EM_S);
287 THREAD_MEMFREE(EM_F);
288 THREAD_MEMFREE(row_index);
289
290 } /* end of pointer check */
291 freeSampleBuffer(ABuff);
292 freeSampleBuffer(BBuff);
293 freeSampleBuffer(CBuff);
294 freeSampleBuffer(DBuff);
295 freeSampleBuffer(XBuff);
296 freeSampleBuffer(YBuff);
297 } /* end parallel region */
298 }
299 /*
300 * $Log$
301 */

  ViewVC Help
Powered by ViewVC 1.1.26