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

  ViewVC Help
Powered by ViewVC 1.1.26