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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.26