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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1312 - (show annotations)
Mon Sep 24 06:18:44 2007 UTC (11 years, 7 months ago) by ksteube
File MIME type: text/plain
File size: 11648 byte(s)
The MPI branch is hereby closed. All future work should be in trunk.

Previously in revision 1295 I merged the latest changes to trunk into trunk-mpi-branch.
In this revision I copied all files from trunk-mpi-branch over the corresponding
trunk files. I did not use 'svn merge', it was a copy.

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

  ViewVC Help
Powered by ViewVC 1.1.26