/[escript]/branches/domexper/dudley/src/Assemble_PDE_Single2_1D.c
ViewVC logotype

Contents of /branches/domexper/dudley/src/Assemble_PDE_Single2_1D.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2271 - (show annotations)
Mon Feb 16 05:08:29 2009 UTC (10 years, 2 months ago) by jfenwick
Original Path: trunk/finley/src/Assemble_PDE_Single2_1D.c
File MIME type: text/plain
File size: 12315 byte(s)
Merging version 2269 to trunk

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 1D 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 = 1 x 1 */
29 /* B = 1 */
30 /* C = 1 */
31 /* D = scalar */
32 /* X = 1 */
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
48 void Finley_Assemble_PDE_Single2_1D(Assemble_Parameters p, Finley_ElementFile* elements,
49 Paso_SystemMatrix* Mat, escriptDataC* F,
50 escriptDataC* A, escriptDataC* B, escriptDataC* C, escriptDataC* D, escriptDataC* X, escriptDataC* Y) {
51
52 #define DIM 1
53 index_t color;
54 dim_t e;
55 __const double *A_p, *B_p, *C_p, *D_p, *X_p, *Y_p;
56 double *EM_S, *EM_F, *Vol, *DSDX;
57 index_t *row_index;
58 register dim_t q, s,r;
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=(requireWrite(F), getSampleDataRW(F,0)); /* use comma, to get around the mixed code and declarations thing */
69 double *S=p.row_jac->ReferenceElement->S;
70 dim_t len_EM_S=p.row_NN*p.col_NN;
71 dim_t len_EM_F=p.row_NN;
72
73 void* ABuff=allocSampleBuffer(A);
74 void* BBuff=allocSampleBuffer(B);
75 void* CBuff=allocSampleBuffer(C);
76 void* DBuff=allocSampleBuffer(D);
77 void* XBuff=allocSampleBuffer(X);
78 void* YBuff=allocSampleBuffer(Y);
79 #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)
80 {
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
86 if (!Finley_checkPtr(EM_S) && !Finley_checkPtr(EM_F) && !Finley_checkPtr(row_index) ) {
87
88 for (color=elements->minColor;color<=elements->maxColor;color++) {
89 /* open loop over all elements: */
90 #pragma omp for private(e) schedule(static)
91 for(e=0;e<elements->numElements;e++){
92 if (elements->Color[e]==color) {
93 Vol=&(p.row_jac->volume[INDEX2(0,e,p.numQuad)]);
94 DSDX=&(p.row_jac->DSDX[INDEX4(0,0,0,e,p.row_NN,DIM,p.numQuad)]);
95 for (q=0;q<len_EM_S;++q) EM_S[q]=0;
96 for (q=0;q<len_EM_F;++q) EM_F[q]=0;
97 add_EM_F=FALSE;
98 add_EM_S=FALSE;
99 /**************************************************************/
100 /* process A: */
101 /**************************************************************/
102 A_p=getSampleDataRO(A,e,ABuff);
103 if (NULL!=A_p) {
104 add_EM_S=TRUE;
105 if (extendedA) {
106 for (s=0;s<p.row_NS;s++) {
107 for (r=0;r<p.col_NS;r++) {
108 rtmp=0;
109 for (q=0;q<p.numQuad;q++) {
110 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)];
111 }
112 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
113 }
114 }
115 } else {
116 for (s=0;s<p.row_NS;s++) {
117 for (r=0;r<p.col_NS;r++) {
118 rtmp=0;
119 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)];
120 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*A_p[INDEX2(0,0,DIM)];
121 }
122 }
123 }
124 }
125 /**************************************************************/
126 /* process B: */
127 /**************************************************************/
128 B_p=getSampleDataRO(B,e,BBuff);
129 if (NULL!=B_p) {
130 add_EM_S=TRUE;
131 if (extendedB) {
132 for (s=0;s<p.row_NS;s++) {
133 for (r=0;r<p.col_NS;r++) {
134 rtmp=0;
135 for (q=0;q<p.numQuad;q++) {
136 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)];
137 }
138 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
139 }
140 }
141 } else {
142 for (s=0;s<p.row_NS;s++) {
143 for (r=0;r<p.col_NS;r++) {
144 rtmp=0;
145 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)];
146 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*B_p[0];
147 }
148 }
149 }
150 }
151 /**************************************************************/
152 /* process C: */
153 /**************************************************************/
154 C_p=getSampleDataRO(C,e,CBuff);
155 if (NULL!=C_p) {
156 add_EM_S=TRUE;
157 if (extendedC) {
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++) {
162 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)];
163 }
164 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
165 }
166 }
167 } else {
168 for (s=0;s<p.row_NS;s++) {
169 for (r=0;r<p.col_NS;r++) {
170 rtmp=0;
171 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)];
172 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*C_p[0];
173 }
174 }
175 }
176 }
177 /************************************************************* */
178 /* process D */
179 /**************************************************************/
180 D_p=getSampleDataRO(D,e,DBuff);
181 if (NULL!=D_p) {
182 add_EM_S=TRUE;
183 if (extendedD) {
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++) {
188 rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*D_p[q]*S[INDEX2(r,q,p.row_NS)];
189 }
190 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp;
191 }
192 }
193 } else {
194 for (s=0;s<p.row_NS;s++) {
195 for (r=0;r<p.col_NS;r++) {
196 rtmp=0;
197 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*S[INDEX2(r,q,p.row_NS)];
198 EM_S[INDEX4(0,0,s,r,p.numEqu,p.numComp,p.row_NN)]+=rtmp*D_p[0];
199 }
200 }
201 }
202 }
203 /**************************************************************/
204 /* process X: */
205 /**************************************************************/
206 X_p=getSampleDataRO(X,e,XBuff);
207 if (NULL!=X_p) {
208 add_EM_F=TRUE;
209 if (extendedX) {
210 for (s=0;s<p.row_NS;s++) {
211 rtmp=0;
212 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)];
213 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
214 }
215 } else {
216 for (s=0;s<p.row_NS;s++) {
217 rtmp=0;
218 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*DSDX[INDEX3(s,0,q,p.row_NN,DIM)];
219 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*X_p[0];
220 }
221 }
222 }
223 /**************************************************************/
224 /* process Y: */
225 /**************************************************************/
226 Y_p=getSampleDataRO(Y,e,YBuff);
227 if (NULL!=Y_p) {
228 add_EM_F=TRUE;
229 if (extendedY) {
230 for (s=0;s<p.row_NS;s++) {
231 rtmp=0;
232 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)]*Y_p[q];
233 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp;
234 }
235 } else {
236 for (s=0;s<p.row_NS;s++) {
237 rtmp=0;
238 for (q=0;q<p.numQuad;q++) rtmp+=Vol[q]*S[INDEX2(s,q,p.row_NS)];
239 EM_F[INDEX2(0,s,p.numEqu)]+=rtmp*Y_p[0];
240 }
241 }
242 }
243 /***********************************************************************************************/
244 /* add the element matrices onto the matrix and right hand side */
245 /***********************************************************************************************/
246 for (q=0;q<p.row_NN;q++) row_index[q]=p.row_DOF[elements->Nodes[INDEX2(p.row_node[q],e,p.NN)]];
247 if (add_EM_F) Finley_Util_AddScatter(p.row_NN,row_index,p.numEqu,EM_F,F_p, p.row_DOF_UpperBound);
248 if (add_EM_S) Finley_Assemble_addToSystemMatrix(Mat,p.row_NN,row_index,p.numEqu,p.col_NN,row_index,p.numComp,EM_S);
249
250 } /* end color check */
251 } /* end element loop */
252 } /* end color loop */
253
254 THREAD_MEMFREE(EM_S); /* these FREEs appear to be inside the if because if any of the allocs */
255 THREAD_MEMFREE(EM_F); /* failed it means an out of memory (which is not recoverable anyway) */
256 THREAD_MEMFREE(row_index);
257
258 } /* end of pointer check */
259 } /* end parallel region */
260 freeSampleBuffer(ABuff);
261 freeSampleBuffer(BBuff);
262 freeSampleBuffer(CBuff);
263 freeSampleBuffer(DBuff);
264 freeSampleBuffer(XBuff);
265 freeSampleBuffer(YBuff);
266 }
267 /*
268 * $Log$
269 */

  ViewVC Help
Powered by ViewVC 1.1.26