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

Contents of /trunk-mpi-branch/finley/src/Assemble_PDE_Single2_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: 11897 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 1.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_{i,j} u_,j)_i-(B_{i} u)_i+C_{j} u_,j-D u_m and -(X_,i)_i + Y */
20
21 /* 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 = 1 x 1 */
27 /* B = 1 */
28 /* C = 1 */
29 /* D = scalar */
30 /* X = 1 */
31 /* Y = scalar */
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_Single2_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;
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;
73 dim_t len_EM_F=p.row_NN;
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,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 rtmp=0;
105 for (q=0;q<p.numQuad;q++) {
106 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)];
107 }
108 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
109 }
110 }
111 } else {
112 for (s=0;s<p.row_NS;s++) {
113 for (r=0;r<p.col_NS;r++) {
114 rtmp=0;
115 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)];
116 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*A_p[INDEX2(0,0,DIM)];
117 }
118 }
119 }
120 }
121 /**************************************************************/
122 /* process B: */
123 /**************************************************************/
124 B_p=getSampleData(B,e);
125 if (NULL!=B_p) {
126 add_EM_S=TRUE;
127 if (extendedB) {
128 for (s=0;s<p.row_NS;s++) {
129 for (r=0;r<p.col_NS;r++) {
130 rtmp=0;
131 for (q=0;q<p.numQuad;q++) {
132 rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*B_p[INDEX2(0,q,DIM)]*S[INDEX2(r,q,p.row_NS)];
133 }
134 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
135 }
136 }
137 } else {
138 for (s=0;s<p.row_NS;s++) {
139 for (r=0;r<p.col_NS;r++) {
140 rtmp=0;
141 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)];
142 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*B_p[0];
143 }
144 }
145 }
146 }
147 /**************************************************************/
148 /* process C: */
149 /**************************************************************/
150 C_p=getSampleData(C,e);
151 if (NULL!=C_p) {
152 add_EM_S=TRUE;
153 if (extendedC) {
154 for (s=0;s<p.row_NS;s++) {
155 for (r=0;r<p.col_NS;r++) {
156 rtmp=0;
157 for (q=0;q<p.numQuad;q++) {
158 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)];
159 }
160 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
161 }
162 }
163 } else {
164 for (s=0;s<p.row_NS;s++) {
165 for (r=0;r<p.col_NS;r++) {
166 rtmp=0;
167 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)];
168 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*C_p[0];
169 }
170 }
171 }
172 }
173 /************************************************************* */
174 /* process D */
175 /**************************************************************/
176 D_p=getSampleData(D,e);
177 if (NULL!=D_p) {
178 add_EM_S=TRUE;
179 if (extendedD) {
180 for (s=0;s<p.row_NS;s++) {
181 for (r=0;r<p.col_NS;r++) {
182 rtmp=0;
183 for (q=0;q<p.numQuad;q++) {
184 rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*D_p[q]*S[INDEX2(r,q,p.row_NS)];
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 rtmp=0;
193 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
194 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[0];
195 }
196 }
197 }
198 }
199 /**************************************************************/
200 /* process X: */
201 /**************************************************************/
202 X_p=getSampleData(X,e);
203 if (NULL!=X_p) {
204 add_EM_F=TRUE;
205 if (extendedX) {
206 for (s=0;s<p.row_NS;s++) {
207 rtmp=0;
208 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)]*X_p[INDEX2(0,q,DIM)];
209 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
210 }
211 } else {
212 for (s=0;s<p.row_NS;s++) {
213 rtmp=0;
214 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
215 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*X_p[0];
216 }
217 }
218 }
219 /**************************************************************/
220 /* process Y: */
221 /**************************************************************/
222 Y_p=getSampleData(Y,e);
223 if (NULL!=Y_p) {
224 add_EM_F=TRUE;
225 if (extendedY) {
226 for (s=0;s<p.row_NS;s++) {
227 rtmp=0;
228 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[q];
229 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
230 }
231 } else {
232 for (s=0;s<p.row_NS;s++) {
233 rtmp=0;
234 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
235 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*Y_p[0];
236 }
237 }
238 }
239 /***********************************************************************************************/
240 /* add the element matrices onto the matrix and right hand side */
241 /***********************************************************************************************/
242 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
243 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
244 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
245
246 } /* end color check */
247 } /* end element loop */
248 } /* end color loop */
249
250 THREAD_MEMFREE(EM_S);
251 THREAD_MEMFREE(EM_F);
252 THREAD_MEMFREE(row_index);
253
254 } /* end of pointer check */
255 } /* end parallel region */
256 }
257 /*
258 * $Log$
259 */

  ViewVC Help
Powered by ViewVC 1.1.26