/[escript]/branches/ripleygmg_from_3668/ripley/src/Rectangle.cpp
ViewVC logotype

Contents of /branches/ripleygmg_from_3668/ripley/src/Rectangle.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3753 - (show annotations)
Tue Jan 3 09:01:49 2012 UTC (7 years, 11 months ago) by caltinay
File size: 109555 byte(s)
Took over Node->DOF interpolation to Brick and parallelised some loops.

1
2 /*******************************************************
3 *
4 * Copyright (c) 2003-2011 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 #include <ripley/Rectangle.h>
15 extern "C" {
16 #include "paso/SystemMatrix.h"
17 }
18
19 #if USE_SILO
20 #include <silo.h>
21 #ifdef ESYS_MPI
22 #include <pmpio.h>
23 #endif
24 #endif
25
26 #include <iomanip>
27
28 using namespace std;
29
30 namespace ripley {
31
32 Rectangle::Rectangle(int n0, int n1, double l0, double l1, int d0, int d1) :
33 RipleyDomain(2),
34 m_gNE0(n0),
35 m_gNE1(n1),
36 m_l0(l0),
37 m_l1(l1),
38 m_NX(d0),
39 m_NY(d1)
40 {
41 // ensure number of subdivisions is valid and nodes can be distributed
42 // among number of ranks
43 if (m_NX*m_NY != m_mpiInfo->size)
44 throw RipleyException("Invalid number of spatial subdivisions");
45
46 if ((n0+1)%m_NX > 0 || (n1+1)%m_NY > 0)
47 throw RipleyException("Number of elements+1 must be separable into number of ranks in each dimension");
48
49 if ((m_NX > 1 && (n0+1)/m_NX<2) || (m_NY > 1 && (n1+1)/m_NY<2))
50 throw RipleyException("Too few elements for the number of ranks");
51
52 // local number of elements (including overlap)
53 m_NE0 = (m_NX>1 ? (n0+1)/m_NX : n0);
54 if (m_mpiInfo->rank%m_NX>0 && m_mpiInfo->rank%m_NX<m_NX-1)
55 m_NE0++;
56 m_NE1 = (m_NY>1 ? (n1+1)/m_NY : n1);
57 if (m_mpiInfo->rank/m_NX>0 && m_mpiInfo->rank/m_NX<m_NY-1)
58 m_NE1++;
59
60 // local number of nodes
61 m_N0 = m_NE0+1;
62 m_N1 = m_NE1+1;
63
64 // bottom-left node is at (offset0,offset1) in global mesh
65 m_offset0 = (n0+1)/m_NX*(m_mpiInfo->rank%m_NX);
66 if (m_offset0 > 0)
67 m_offset0--;
68 m_offset1 = (n1+1)/m_NY*(m_mpiInfo->rank/m_NX);
69 if (m_offset1 > 0)
70 m_offset1--;
71
72 populateSampleIds();
73 }
74
75 Rectangle::~Rectangle()
76 {
77 }
78
79 string Rectangle::getDescription() const
80 {
81 return "ripley::Rectangle";
82 }
83
84 bool Rectangle::operator==(const AbstractDomain& other) const
85 {
86 const Rectangle* o=dynamic_cast<const Rectangle*>(&other);
87 if (o) {
88 return (RipleyDomain::operator==(other) &&
89 m_gNE0==o->m_gNE0 && m_gNE1==o->m_gNE1
90 && m_l0==o->m_l0 && m_l1==o->m_l1
91 && m_NX==o->m_NX && m_NY==o->m_NY);
92 }
93
94 return false;
95 }
96
97 void Rectangle::dump(const string& fileName) const
98 {
99 #if USE_SILO
100 string fn(fileName);
101 if (fileName.length() < 6 || fileName.compare(fileName.length()-5, 5, ".silo") != 0) {
102 fn+=".silo";
103 }
104
105 const int NUM_SILO_FILES = 1;
106 const char* blockDirFmt = "/block%04d";
107 int driver=DB_HDF5;
108 string siloPath;
109 DBfile* dbfile = NULL;
110
111 #ifdef ESYS_MPI
112 PMPIO_baton_t* baton = NULL;
113 #endif
114
115 if (m_mpiInfo->size > 1) {
116 #ifdef ESYS_MPI
117 baton = PMPIO_Init(NUM_SILO_FILES, PMPIO_WRITE, m_mpiInfo->comm,
118 0x1337, PMPIO_DefaultCreate, PMPIO_DefaultOpen,
119 PMPIO_DefaultClose, (void*)&driver);
120 // try the fallback driver in case of error
121 if (!baton && driver != DB_PDB) {
122 driver = DB_PDB;
123 baton = PMPIO_Init(NUM_SILO_FILES, PMPIO_WRITE, m_mpiInfo->comm,
124 0x1338, PMPIO_DefaultCreate, PMPIO_DefaultOpen,
125 PMPIO_DefaultClose, (void*)&driver);
126 }
127 if (baton) {
128 char str[64];
129 snprintf(str, 64, blockDirFmt, PMPIO_RankInGroup(baton, m_mpiInfo->rank));
130 siloPath = str;
131 dbfile = (DBfile*) PMPIO_WaitForBaton(baton, fn.c_str(), siloPath.c_str());
132 }
133 #endif
134 } else {
135 dbfile = DBCreate(fn.c_str(), DB_CLOBBER, DB_LOCAL,
136 getDescription().c_str(), driver);
137 // try the fallback driver in case of error
138 if (!dbfile && driver != DB_PDB) {
139 driver = DB_PDB;
140 dbfile = DBCreate(fn.c_str(), DB_CLOBBER, DB_LOCAL,
141 getDescription().c_str(), driver);
142 }
143 }
144
145 if (!dbfile)
146 throw RipleyException("dump: Could not create Silo file");
147
148 /*
149 if (driver==DB_HDF5) {
150 // gzip level 1 already provides good compression with minimal
151 // performance penalty. Some tests showed that gzip levels >3 performed
152 // rather badly on escript data both in terms of time and space
153 DBSetCompression("ERRMODE=FALLBACK METHOD=GZIP LEVEL=1");
154 }
155 */
156
157 boost::scoped_ptr<double> x(new double[m_N0]);
158 boost::scoped_ptr<double> y(new double[m_N1]);
159 double* coords[2] = { x.get(), y.get() };
160 pair<double,double> xdx = getFirstCoordAndSpacing(0);
161 pair<double,double> ydy = getFirstCoordAndSpacing(1);
162 #pragma omp parallel
163 {
164 #pragma omp for nowait
165 for (dim_t i0 = 0; i0 < m_N0; i0++) {
166 coords[0][i0]=xdx.first+i0*xdx.second;
167 }
168 #pragma omp for nowait
169 for (dim_t i1 = 0; i1 < m_N1; i1++) {
170 coords[1][i1]=ydy.first+i1*ydy.second;
171 }
172 }
173 IndexVector dims = getNumNodesPerDim();
174
175 // write mesh
176 DBPutQuadmesh(dbfile, "mesh", NULL, coords, &dims[0], 2, DB_DOUBLE,
177 DB_COLLINEAR, NULL);
178
179 // write node ids
180 DBPutQuadvar1(dbfile, "nodeId", "mesh", (void*)&m_nodeId[0], &dims[0], 2,
181 NULL, 0, DB_INT, DB_NODECENT, NULL);
182
183 // write element ids
184 dims = getNumElementsPerDim();
185 DBPutQuadvar1(dbfile, "elementId", "mesh", (void*)&m_elementId[0],
186 &dims[0], 2, NULL, 0, DB_INT, DB_ZONECENT, NULL);
187
188 // rank 0 writes multimesh and multivar
189 if (m_mpiInfo->rank == 0) {
190 vector<string> tempstrings;
191 vector<char*> names;
192 for (dim_t i=0; i<m_mpiInfo->size; i++) {
193 stringstream path;
194 path << "/block" << setw(4) << setfill('0') << right << i << "/mesh";
195 tempstrings.push_back(path.str());
196 names.push_back((char*)tempstrings.back().c_str());
197 }
198 vector<int> types(m_mpiInfo->size, DB_QUAD_RECT);
199 DBSetDir(dbfile, "/");
200 DBPutMultimesh(dbfile, "multimesh", m_mpiInfo->size, &names[0],
201 &types[0], NULL);
202 tempstrings.clear();
203 names.clear();
204 for (dim_t i=0; i<m_mpiInfo->size; i++) {
205 stringstream path;
206 path << "/block" << setw(4) << setfill('0') << right << i << "/nodeId";
207 tempstrings.push_back(path.str());
208 names.push_back((char*)tempstrings.back().c_str());
209 }
210 types.assign(m_mpiInfo->size, DB_QUADVAR);
211 DBPutMultivar(dbfile, "nodeId", m_mpiInfo->size, &names[0],
212 &types[0], NULL);
213 tempstrings.clear();
214 names.clear();
215 for (dim_t i=0; i<m_mpiInfo->size; i++) {
216 stringstream path;
217 path << "/block" << setw(4) << setfill('0') << right << i << "/elementId";
218 tempstrings.push_back(path.str());
219 names.push_back((char*)tempstrings.back().c_str());
220 }
221 DBPutMultivar(dbfile, "elementId", m_mpiInfo->size, &names[0],
222 &types[0], NULL);
223 }
224
225 if (m_mpiInfo->size > 1) {
226 #ifdef ESYS_MPI
227 PMPIO_HandOffBaton(baton, dbfile);
228 PMPIO_Finish(baton);
229 #endif
230 } else {
231 DBClose(dbfile);
232 }
233
234 #else // USE_SILO
235 throw RipleyException("dump(): no Silo support");
236 #endif
237 }
238
239 const int* Rectangle::borrowSampleReferenceIDs(int fsType) const
240 {
241 switch (fsType) {
242 case Nodes:
243 case ReducedNodes: //FIXME: reduced
244 return &m_nodeId[0];
245 case DegreesOfFreedom:
246 case ReducedDegreesOfFreedom: //FIXME: reduced
247 return &m_dofId[0];
248 case Elements:
249 case ReducedElements:
250 return &m_elementId[0];
251 case FaceElements:
252 case ReducedFaceElements:
253 return &m_faceId[0];
254 default:
255 break;
256 }
257
258 stringstream msg;
259 msg << "borrowSampleReferenceIDs() not implemented for function space type "
260 << functionSpaceTypeAsString(fsType);
261 throw RipleyException(msg.str());
262 }
263
264 bool Rectangle::ownSample(int fsCode, index_t id) const
265 {
266 #ifdef ESYS_MPI
267 if (fsCode == Nodes) {
268 const index_t left = (m_offset0==0 ? 0 : 1);
269 const index_t bottom = (m_offset1==0 ? 0 : 1);
270 const index_t right = (m_mpiInfo->rank%m_NX==m_NX-1 ? m_N0 : m_N0-1);
271 const index_t top = (m_mpiInfo->rank/m_NX==m_NY-1 ? m_N1 : m_N1-1);
272 const index_t x=id%m_N0;
273 const index_t y=id/m_N0;
274 return (x>=left && x<right && y>=bottom && y<top);
275 } else {
276 stringstream msg;
277 msg << "ownSample() not implemented for "
278 << functionSpaceTypeAsString(fsCode);
279 throw RipleyException(msg.str());
280 }
281 #else
282 return true;
283 #endif
284 }
285
286 void Rectangle::setToGradient(escript::Data& out, const escript::Data& cIn) const
287 {
288 escript::Data& in = *const_cast<escript::Data*>(&cIn);
289 const dim_t numComp = in.getDataPointSize();
290 const double h0 = m_l0/m_gNE0;
291 const double h1 = m_l1/m_gNE1;
292 const double cx0 = -1./h0;
293 const double cx1 = -.78867513459481288225/h0;
294 const double cx2 = -.5/h0;
295 const double cx3 = -.21132486540518711775/h0;
296 const double cx4 = .21132486540518711775/h0;
297 const double cx5 = .5/h0;
298 const double cx6 = .78867513459481288225/h0;
299 const double cx7 = 1./h0;
300 const double cy0 = -1./h1;
301 const double cy1 = -.78867513459481288225/h1;
302 const double cy2 = -.5/h1;
303 const double cy3 = -.21132486540518711775/h1;
304 const double cy4 = .21132486540518711775/h1;
305 const double cy5 = .5/h1;
306 const double cy6 = .78867513459481288225/h1;
307 const double cy7 = 1./h1;
308
309 if (out.getFunctionSpace().getTypeCode() == Elements) {
310 /*** GENERATOR SNIP_GRAD_ELEMENTS TOP */
311 #pragma omp parallel for
312 for (index_t k1=0; k1 < m_NE1; ++k1) {
313 for (index_t k0=0; k0 < m_NE0; ++k0) {
314 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,k1, m_N0));
315 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0));
316 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,k1+1, m_N0));
317 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,k1, m_N0));
318 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
319 for (index_t i=0; i < numComp; ++i) {
320 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
321 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
322 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
323 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
324 o[INDEX3(i,0,2,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
325 o[INDEX3(i,1,2,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
326 o[INDEX3(i,0,3,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
327 o[INDEX3(i,1,3,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
328 } /* end of component loop i */
329 } /* end of k0 loop */
330 } /* end of k1 loop */
331 /* GENERATOR SNIP_GRAD_ELEMENTS BOTTOM */
332 } else if (out.getFunctionSpace().getTypeCode() == ReducedElements) {
333 /*** GENERATOR SNIP_GRAD_REDUCED_ELEMENTS TOP */
334 #pragma omp parallel for
335 for (index_t k1=0; k1 < m_NE1; ++k1) {
336 for (index_t k0=0; k0 < m_NE0; ++k0) {
337 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,k1, m_N0));
338 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0));
339 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,k1+1, m_N0));
340 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,k1, m_N0));
341 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
342 for (index_t i=0; i < numComp; ++i) {
343 o[INDEX3(i,0,0,numComp,2)] = cx5*(f_10[i] + f_11[i]) + cx2*(f_00[i] + f_01[i]);
344 o[INDEX3(i,1,0,numComp,2)] = cy2*(f_00[i] + f_10[i]) + cy5*(f_01[i] + f_11[i]);
345 } /* end of component loop i */
346 } /* end of k0 loop */
347 } /* end of k1 loop */
348 /* GENERATOR SNIP_GRAD_REDUCED_ELEMENTS BOTTOM */
349 } else if (out.getFunctionSpace().getTypeCode() == FaceElements) {
350 #pragma omp parallel
351 {
352 /*** GENERATOR SNIP_GRAD_FACES TOP */
353 if (m_faceOffset[0] > -1) {
354 #pragma omp for nowait
355 for (index_t k1=0; k1 < m_NE1; ++k1) {
356 const register double* f_10 = in.getSampleDataRO(INDEX2(1,k1, m_N0));
357 const register double* f_11 = in.getSampleDataRO(INDEX2(1,k1+1, m_N0));
358 const register double* f_01 = in.getSampleDataRO(INDEX2(0,k1+1, m_N0));
359 const register double* f_00 = in.getSampleDataRO(INDEX2(0,k1, m_N0));
360 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
361 for (index_t i=0; i < numComp; ++i) {
362 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
363 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy0 + f_01[i]*cy7;
364 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
365 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy0 + f_01[i]*cy7;
366 } /* end of component loop i */
367 } /* end of k1 loop */
368 } /* end of face 0 */
369 if (m_faceOffset[1] > -1) {
370 #pragma omp for nowait
371 for (index_t k1=0; k1 < m_NE1; ++k1) {
372 const register double* f_10 = in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0));
373 const register double* f_11 = in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0));
374 const register double* f_01 = in.getSampleDataRO(INDEX2(m_N0-2,k1+1, m_N0));
375 const register double* f_00 = in.getSampleDataRO(INDEX2(m_N0-2,k1, m_N0));
376 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
377 for (index_t i=0; i < numComp; ++i) {
378 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
379 o[INDEX3(i,1,0,numComp,2)] = f_10[i]*cy0 + f_11[i]*cy7;
380 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
381 o[INDEX3(i,1,1,numComp,2)] = f_10[i]*cy0 + f_11[i]*cy7;
382 } /* end of component loop i */
383 } /* end of k1 loop */
384 } /* end of face 1 */
385 if (m_faceOffset[2] > -1) {
386 #pragma omp for nowait
387 for (index_t k0=0; k0 < m_NE0; ++k0) {
388 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,0, m_N0));
389 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,0, m_N0));
390 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,1, m_N0));
391 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,1, m_N0));
392 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
393 for (index_t i=0; i < numComp; ++i) {
394 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx0 + f_10[i]*cx7;
395 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
396 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx0 + f_10[i]*cx7;
397 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
398 } /* end of component loop i */
399 } /* end of k0 loop */
400 } /* end of face 2 */
401 if (m_faceOffset[3] > -1) {
402 #pragma omp for nowait
403 for (index_t k0=0; k0 < m_NE0; ++k0) {
404 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0));
405 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0));
406 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,m_N1-2, m_N0));
407 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,m_N1-2, m_N0));
408 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
409 for (index_t i=0; i < numComp; ++i) {
410 o[INDEX3(i,0,0,numComp,2)] = f_01[i]*cx0 + f_11[i]*cx7;
411 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
412 o[INDEX3(i,0,1,numComp,2)] = f_01[i]*cx0 + f_11[i]*cx7;
413 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
414 } /* end of component loop i */
415 } /* end of k0 loop */
416 } /* end of face 3 */
417 /* GENERATOR SNIP_GRAD_FACES BOTTOM */
418 } // end of parallel section
419 } else if (out.getFunctionSpace().getTypeCode() == ReducedFaceElements) {
420 #pragma omp parallel
421 {
422 /*** GENERATOR SNIP_GRAD_REDUCED_FACES TOP */
423 if (m_faceOffset[0] > -1) {
424 #pragma omp for nowait
425 for (index_t k1=0; k1 < m_NE1; ++k1) {
426 const register double* f_10 = in.getSampleDataRO(INDEX2(1,k1, m_N0));
427 const register double* f_11 = in.getSampleDataRO(INDEX2(1,k1+1, m_N0));
428 const register double* f_01 = in.getSampleDataRO(INDEX2(0,k1+1, m_N0));
429 const register double* f_00 = in.getSampleDataRO(INDEX2(0,k1, m_N0));
430 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
431 for (index_t i=0; i < numComp; ++i) {
432 o[INDEX3(i,0,0,numComp,2)] = cx5*(f_10[i] + f_11[i]) + cx2*(f_00[i] + f_01[i]);
433 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy0 + f_01[i]*cy7;
434 } /* end of component loop i */
435 } /* end of k1 loop */
436 } /* end of face 0 */
437 if (m_faceOffset[1] > -1) {
438 #pragma omp for nowait
439 for (index_t k1=0; k1 < m_NE1; ++k1) {
440 const register double* f_10 = in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0));
441 const register double* f_11 = in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0));
442 const register double* f_01 = in.getSampleDataRO(INDEX2(m_N0-2,k1+1, m_N0));
443 const register double* f_00 = in.getSampleDataRO(INDEX2(m_N0-2,k1, m_N0));
444 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
445 for (index_t i=0; i < numComp; ++i) {
446 o[INDEX3(i,0,0,numComp,2)] = cx5*(f_10[i] + f_11[i]) + cx2*(f_00[i] + f_01[i]);
447 o[INDEX3(i,1,0,numComp,2)] = f_10[i]*cy0 + f_11[i]*cy7;
448 } /* end of component loop i */
449 } /* end of k1 loop */
450 } /* end of face 1 */
451 if (m_faceOffset[2] > -1) {
452 #pragma omp for nowait
453 for (index_t k0=0; k0 < m_NE0; ++k0) {
454 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,0, m_N0));
455 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,0, m_N0));
456 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,1, m_N0));
457 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,1, m_N0));
458 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
459 for (index_t i=0; i < numComp; ++i) {
460 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx0 + f_10[i]*cx7;
461 o[INDEX3(i,1,0,numComp,2)] = cy2*(f_00[i] + f_10[i]) + cy5*(f_01[i] + f_11[i]);
462 } /* end of component loop i */
463 } /* end of k0 loop */
464 } /* end of face 2 */
465 if (m_faceOffset[3] > -1) {
466 #pragma omp for nowait
467 for (index_t k0=0; k0 < m_NE0; ++k0) {
468 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0));
469 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0));
470 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,m_N1-2, m_N0));
471 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,m_N1-2, m_N0));
472 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
473 for (index_t i=0; i < numComp; ++i) {
474 o[INDEX3(i,0,0,numComp,2)] = f_01[i]*cx0 + f_11[i]*cx7;
475 o[INDEX3(i,1,0,numComp,2)] = cy5*(f_01[i] + f_11[i]) + cy2*(f_00[i] + f_10[i]);
476 } /* end of component loop i */
477 } /* end of k0 loop */
478 } /* end of face 3 */
479 /* GENERATOR SNIP_GRAD_REDUCED_FACES BOTTOM */
480 } // end of parallel section
481 } else {
482 stringstream msg;
483 msg << "setToGradient() not implemented for "
484 << functionSpaceTypeAsString(out.getFunctionSpace().getTypeCode());
485 throw RipleyException(msg.str());
486 }
487 }
488
489 void Rectangle::setToIntegrals(vector<double>& integrals, const escript::Data& arg) const
490 {
491 escript::Data& in = *const_cast<escript::Data*>(&arg);
492 const dim_t numComp = in.getDataPointSize();
493 const double h0 = m_l0/m_gNE0;
494 const double h1 = m_l1/m_gNE1;
495 if (arg.getFunctionSpace().getTypeCode() == Elements) {
496 const double w_0 = h0*h1/4.;
497 #pragma omp parallel
498 {
499 vector<double> int_local(numComp, 0);
500 #pragma omp for nowait
501 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
502 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
503 const double* f = in.getSampleDataRO(INDEX2(k0, k1, m_NE0));
504 for (index_t i=0; i < numComp; ++i) {
505 const register double f_0 = f[INDEX2(i,0,numComp)];
506 const register double f_1 = f[INDEX2(i,1,numComp)];
507 const register double f_2 = f[INDEX2(i,2,numComp)];
508 const register double f_3 = f[INDEX2(i,3,numComp)];
509 int_local[i]+=(f_0+f_1+f_2+f_3)*w_0;
510 } /* end of component loop i */
511 } /* end of k0 loop */
512 } /* end of k1 loop */
513
514 #pragma omp critical
515 for (index_t i=0; i<numComp; i++)
516 integrals[i]+=int_local[i];
517 } // end of parallel section
518 } else if (arg.getFunctionSpace().getTypeCode() == ReducedElements) {
519 const double w_0 = h0*h1;
520 #pragma omp parallel
521 {
522 vector<double> int_local(numComp, 0);
523 #pragma omp for nowait
524 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
525 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
526 const double* f = in.getSampleDataRO(INDEX2(k0, k1, m_NE0));
527 for (index_t i=0; i < numComp; ++i) {
528 int_local[i]+=f[i]*w_0;
529 } /* end of component loop i */
530 } /* end of k0 loop */
531 } /* end of k1 loop */
532
533 #pragma omp critical
534 for (index_t i=0; i<numComp; i++)
535 integrals[i]+=int_local[i];
536 } // end of parallel section
537 } else if (arg.getFunctionSpace().getTypeCode() == FaceElements) {
538 const double w_0 = h0/2.;
539 const double w_1 = h1/2.;
540 #pragma omp parallel
541 {
542 vector<double> int_local(numComp, 0);
543 if (m_faceOffset[0] > -1) {
544 #pragma omp for nowait
545 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
546 const double* f = in.getSampleDataRO(m_faceOffset[0]+k1);
547 for (index_t i=0; i < numComp; ++i) {
548 const register double f_0 = f[INDEX2(i,0,numComp)];
549 const register double f_1 = f[INDEX2(i,1,numComp)];
550 int_local[i]+=(f_0+f_1)*w_1;
551 } /* end of component loop i */
552 } /* end of k1 loop */
553 }
554
555 if (m_faceOffset[1] > -1) {
556 #pragma omp for nowait
557 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
558 const double* f = in.getSampleDataRO(m_faceOffset[1]+k1);
559 for (index_t i=0; i < numComp; ++i) {
560 const register double f_0 = f[INDEX2(i,0,numComp)];
561 const register double f_1 = f[INDEX2(i,1,numComp)];
562 int_local[i]+=(f_0+f_1)*w_1;
563 } /* end of component loop i */
564 } /* end of k1 loop */
565 }
566
567 if (m_faceOffset[2] > -1) {
568 #pragma omp for nowait
569 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
570 const double* f = in.getSampleDataRO(m_faceOffset[2]+k0);
571 for (index_t i=0; i < numComp; ++i) {
572 const register double f_0 = f[INDEX2(i,0,numComp)];
573 const register double f_1 = f[INDEX2(i,1,numComp)];
574 int_local[i]+=(f_0+f_1)*w_0;
575 } /* end of component loop i */
576 } /* end of k0 loop */
577 }
578
579 if (m_faceOffset[3] > -1) {
580 #pragma omp for nowait
581 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
582 const double* f = in.getSampleDataRO(m_faceOffset[3]+k0);
583 for (index_t i=0; i < numComp; ++i) {
584 const register double f_0 = f[INDEX2(i,0,numComp)];
585 const register double f_1 = f[INDEX2(i,1,numComp)];
586 int_local[i]+=(f_0+f_1)*w_0;
587 } /* end of component loop i */
588 } /* end of k0 loop */
589 }
590
591 #pragma omp critical
592 for (index_t i=0; i<numComp; i++)
593 integrals[i]+=int_local[i];
594 } // end of parallel section
595 } else if (arg.getFunctionSpace().getTypeCode() == ReducedFaceElements) {
596 #pragma omp parallel
597 {
598 vector<double> int_local(numComp, 0);
599 if (m_faceOffset[0] > -1) {
600 #pragma omp for nowait
601 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
602 const double* f = in.getSampleDataRO(m_faceOffset[0]+k1);
603 for (index_t i=0; i < numComp; ++i) {
604 int_local[i]+=f[i]*h1;
605 } /* end of component loop i */
606 } /* end of k1 loop */
607 }
608
609 if (m_faceOffset[1] > -1) {
610 #pragma omp for nowait
611 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
612 const double* f = in.getSampleDataRO(m_faceOffset[1]+k1);
613 for (index_t i=0; i < numComp; ++i) {
614 int_local[i]+=f[i]*h1;
615 } /* end of component loop i */
616 } /* end of k1 loop */
617 }
618
619 if (m_faceOffset[2] > -1) {
620 #pragma omp for nowait
621 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
622 const double* f = in.getSampleDataRO(m_faceOffset[2]+k0);
623 for (index_t i=0; i < numComp; ++i) {
624 int_local[i]+=f[i]*h0;
625 } /* end of component loop i */
626 } /* end of k0 loop */
627 }
628
629 if (m_faceOffset[3] > -1) {
630 #pragma omp for nowait
631 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
632 const double* f = in.getSampleDataRO(m_faceOffset[3]+k0);
633 for (index_t i=0; i < numComp; ++i) {
634 int_local[i]+=f[i]*h0;
635 } /* end of component loop i */
636 } /* end of k0 loop */
637 }
638
639 #pragma omp critical
640 for (index_t i=0; i<numComp; i++)
641 integrals[i]+=int_local[i];
642 } // end of parallel section
643 } else {
644 stringstream msg;
645 msg << "setToIntegrals() not implemented for "
646 << functionSpaceTypeAsString(arg.getFunctionSpace().getTypeCode());
647 throw RipleyException(msg.str());
648 }
649 }
650
651 void Rectangle::setToNormal(escript::Data& out) const
652 {
653 if (out.getFunctionSpace().getTypeCode() == FaceElements) {
654 #pragma omp parallel
655 {
656 if (m_faceOffset[0] > -1) {
657 #pragma omp for nowait
658 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
659 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
660 // set vector at two quadrature points
661 *o++ = -1.;
662 *o++ = 0.;
663 *o++ = -1.;
664 *o = 0.;
665 }
666 }
667
668 if (m_faceOffset[1] > -1) {
669 #pragma omp for nowait
670 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
671 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
672 // set vector at two quadrature points
673 *o++ = 1.;
674 *o++ = 0.;
675 *o++ = 1.;
676 *o = 0.;
677 }
678 }
679
680 if (m_faceOffset[2] > -1) {
681 #pragma omp for nowait
682 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
683 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
684 // set vector at two quadrature points
685 *o++ = 0.;
686 *o++ = -1.;
687 *o++ = 0.;
688 *o = -1.;
689 }
690 }
691
692 if (m_faceOffset[3] > -1) {
693 #pragma omp for nowait
694 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
695 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
696 // set vector at two quadrature points
697 *o++ = 0.;
698 *o++ = 1.;
699 *o++ = 0.;
700 *o = 1.;
701 }
702 }
703 } // end of parallel section
704 } else if (out.getFunctionSpace().getTypeCode() == ReducedFaceElements) {
705 #pragma omp parallel
706 {
707 if (m_faceOffset[0] > -1) {
708 #pragma omp for nowait
709 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
710 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
711 *o++ = -1.;
712 *o = 0.;
713 }
714 }
715
716 if (m_faceOffset[1] > -1) {
717 #pragma omp for nowait
718 for (index_t k1 = 0; k1 < m_NE1; ++k1) {
719 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
720 *o++ = 1.;
721 *o = 0.;
722 }
723 }
724
725 if (m_faceOffset[2] > -1) {
726 #pragma omp for nowait
727 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
728 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
729 *o++ = 0.;
730 *o = -1.;
731 }
732 }
733
734 if (m_faceOffset[3] > -1) {
735 #pragma omp for nowait
736 for (index_t k0 = 0; k0 < m_NE0; ++k0) {
737 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
738 *o++ = 0.;
739 *o = 1.;
740 }
741 }
742 } // end of parallel section
743
744 } else {
745 stringstream msg;
746 msg << "setToNormal() not implemented for "
747 << functionSpaceTypeAsString(out.getFunctionSpace().getTypeCode());
748 throw RipleyException(msg.str());
749 }
750 }
751
752 Paso_SystemMatrixPattern* Rectangle::getPattern(bool reducedRowOrder,
753 bool reducedColOrder) const
754 {
755 if (reducedRowOrder || reducedColOrder)
756 throw RipleyException("getPattern() not implemented for reduced order");
757
758 // connector
759 RankVector neighbour;
760 IndexVector offsetInShared(1,0);
761 IndexVector sendShared, recvShared;
762 const IndexVector faces=getNumFacesPerBoundary();
763 const index_t nDOF0 = (m_gNE0+1)/m_NX;
764 const index_t nDOF1 = (m_gNE1+1)/m_NY;
765 const int numDOF=nDOF0*nDOF1;
766 const index_t left = (m_offset0==0 ? 0 : 1);
767 const index_t bottom = (m_offset1==0 ? 0 : 1);
768 vector<IndexVector> colIndices(numDOF); // for the couple blocks
769 int numShared=0;
770
771 m_dofMap.assign(getNumNodes(), 0);
772 #pragma omp parallel for
773 for (index_t i=bottom; i<m_N1; i++) {
774 for (index_t j=left; j<m_N0; j++) {
775 m_dofMap[i*m_N0+j]=(i-bottom)*nDOF0+j-left;
776 }
777 }
778
779 // corner node from bottom-left
780 if (faces[0] == 0 && faces[2] == 0) {
781 neighbour.push_back(m_mpiInfo->rank-m_NX-1);
782 offsetInShared.push_back(offsetInShared.back()+1);
783 sendShared.push_back(0);
784 recvShared.push_back(numDOF+numShared);
785 colIndices[0].push_back(numShared);
786 m_dofMap[0]=numDOF+numShared;
787 ++numShared;
788 }
789 // bottom edge
790 if (faces[2] == 0) {
791 neighbour.push_back(m_mpiInfo->rank-m_NX);
792 offsetInShared.push_back(offsetInShared.back()+nDOF0);
793 for (dim_t i=0; i<nDOF0; i++, numShared++) {
794 sendShared.push_back(i);
795 recvShared.push_back(numDOF+numShared);
796 if (i>0)
797 colIndices[i-1].push_back(numShared);
798 colIndices[i].push_back(numShared);
799 if (i<nDOF0-1)
800 colIndices[i+1].push_back(numShared);
801 m_dofMap[i+left]=numDOF+numShared;
802 }
803 }
804 // corner node from bottom-right
805 if (faces[1] == 0 && faces[2] == 0) {
806 neighbour.push_back(m_mpiInfo->rank-m_NX+1);
807 offsetInShared.push_back(offsetInShared.back()+1);
808 sendShared.push_back(nDOF0-1);
809 recvShared.push_back(numDOF+numShared);
810 colIndices[nDOF0-2].push_back(numShared);
811 colIndices[nDOF0-1].push_back(numShared);
812 m_dofMap[m_N0-1]=numDOF+numShared;
813 ++numShared;
814 }
815 // right edge
816 if (faces[1] == 0) {
817 neighbour.push_back(m_mpiInfo->rank+1);
818 offsetInShared.push_back(offsetInShared.back()+nDOF1);
819 for (dim_t i=0; i<nDOF1; i++, numShared++) {
820 sendShared.push_back((i+1)*(nDOF0)-1);
821 recvShared.push_back(numDOF+numShared);
822 if (i>0)
823 colIndices[i*(nDOF0)-1].push_back(numShared);
824 colIndices[(i+1)*(nDOF0)-1].push_back(numShared);
825 if (i<nDOF1-1)
826 colIndices[(i+2)*(nDOF0)-1].push_back(numShared);
827 m_dofMap[(i+bottom+1)*m_N0-1]=numDOF+numShared;
828 }
829 }
830 // corner node from top-right
831 if (faces[1] == 0 && faces[3] == 0) {
832 neighbour.push_back(m_mpiInfo->rank+m_NX+1);
833 offsetInShared.push_back(offsetInShared.back()+1);
834 sendShared.push_back(numDOF-1);
835 recvShared.push_back(numDOF+numShared);
836 colIndices[numDOF-1].push_back(numShared);
837 m_dofMap[m_N0*m_N1-1]=numDOF+numShared;
838 ++numShared;
839 }
840 // top edge
841 if (faces[3] == 0) {
842 neighbour.push_back(m_mpiInfo->rank+m_NX);
843 offsetInShared.push_back(offsetInShared.back()+nDOF0);
844 for (dim_t i=0; i<nDOF0; i++, numShared++) {
845 sendShared.push_back(numDOF-nDOF0+i);
846 recvShared.push_back(numDOF+numShared);
847 if (i>0)
848 colIndices[numDOF-nDOF0+i-1].push_back(numShared);
849 colIndices[numDOF-nDOF0+i].push_back(numShared);
850 if (i<nDOF0-1)
851 colIndices[numDOF-nDOF0+i+1].push_back(numShared);
852 m_dofMap[m_N0*(m_N1-1)+i+left]=numDOF+numShared;
853 }
854 }
855 // corner node from top-left
856 if (faces[0] == 0 && faces[3] == 0) {
857 neighbour.push_back(m_mpiInfo->rank+m_NX-1);
858 offsetInShared.push_back(offsetInShared.back()+1);
859 sendShared.push_back(numDOF-nDOF0);
860 recvShared.push_back(numDOF+numShared);
861 colIndices[numDOF-nDOF0].push_back(numShared);
862 m_dofMap[m_N0*(m_N1-1)]=numDOF+numShared;
863 ++numShared;
864 }
865 // left edge
866 if (faces[0] == 0) {
867 neighbour.push_back(m_mpiInfo->rank-1);
868 offsetInShared.push_back(offsetInShared.back()+nDOF1);
869 for (dim_t i=0; i<nDOF1; i++, numShared++) {
870 sendShared.push_back(i*nDOF0);
871 recvShared.push_back(numDOF+numShared);
872 if (i>0)
873 colIndices[(i-1)*nDOF0].push_back(numShared);
874 colIndices[i*nDOF0].push_back(numShared);
875 if (i<nDOF1-1)
876 colIndices[(i+1)*nDOF0].push_back(numShared);
877 m_dofMap[(i+bottom)*m_N0]=numDOF+numShared;
878 }
879 }
880
881 /*
882 cout << "--- rcv_shcomp ---" << endl;
883 cout << "numDOF=" << numDOF << ", numNeighbors=" << neighbour.size() << endl;
884 for (size_t i=0; i<neighbour.size(); i++) {
885 cout << "neighbor[" << i << "]=" << neighbour[i]
886 << " offsetInShared[" << i+1 << "]=" << offsetInShared[i+1] << endl;
887 }
888 for (size_t i=0; i<recvShared.size(); i++) {
889 cout << "shared[" << i << "]=" << recvShared[i] << endl;
890 }
891 cout << "--- snd_shcomp ---" << endl;
892 for (size_t i=0; i<sendShared.size(); i++) {
893 cout << "shared[" << i << "]=" << sendShared[i] << endl;
894 }
895 */
896
897 Paso_SharedComponents *snd_shcomp = Paso_SharedComponents_alloc(
898 numDOF, neighbour.size(), &neighbour[0], &sendShared[0],
899 &offsetInShared[0], 1, 0, m_mpiInfo);
900 Paso_SharedComponents *rcv_shcomp = Paso_SharedComponents_alloc(
901 numDOF, neighbour.size(), &neighbour[0], &recvShared[0],
902 &offsetInShared[0], 1, 0, m_mpiInfo);
903 Paso_Connector* connector = Paso_Connector_alloc(snd_shcomp, rcv_shcomp);
904 Paso_SharedComponents_free(snd_shcomp);
905 Paso_SharedComponents_free(rcv_shcomp);
906
907 // create patterns
908 dim_t M, N;
909 IndexVector ptr(1,0);
910 IndexVector index;
911
912 // main pattern
913 for (index_t i=0; i<numDOF; i++) {
914 // always add the node itself
915 index.push_back(i);
916 const int num=insertNeighbours(index, i);
917 ptr.push_back(ptr.back()+num+1);
918 }
919 M=N=ptr.size()-1;
920 // paso will manage the memory
921 index_t* indexC = MEMALLOC(index.size(),index_t);
922 index_t* ptrC = MEMALLOC(ptr.size(), index_t);
923 copy(index.begin(), index.end(), indexC);
924 copy(ptr.begin(), ptr.end(), ptrC);
925 Paso_Pattern *mainPattern = Paso_Pattern_alloc(MATRIX_FORMAT_DEFAULT,
926 M, N, ptrC, indexC);
927
928 /*
929 cout << "--- main_pattern ---" << endl;
930 cout << "M=" << M << ", N=" << N << endl;
931 for (size_t i=0; i<ptr.size(); i++) {
932 cout << "ptr[" << i << "]=" << ptr[i] << endl;
933 }
934 for (size_t i=0; i<index.size(); i++) {
935 cout << "index[" << i << "]=" << index[i] << endl;
936 }
937 */
938
939 // column & row couple patterns
940 ptr.assign(1, 0);
941 index.clear();
942
943 for (index_t i=0; i<numDOF; i++) {
944 index.insert(index.end(), colIndices[i].begin(), colIndices[i].end());
945 ptr.push_back(ptr.back()+colIndices[i].size());
946 }
947
948 // paso will manage the memory
949 indexC = MEMALLOC(index.size(), index_t);
950 ptrC = MEMALLOC(ptr.size(), index_t);
951 copy(index.begin(), index.end(), indexC);
952 copy(ptr.begin(), ptr.end(), ptrC);
953 M=ptr.size()-1;
954 N=numShared;
955 Paso_Pattern *colCouplePattern=Paso_Pattern_alloc(MATRIX_FORMAT_DEFAULT,
956 M, N, ptrC, indexC);
957
958 /*
959 cout << "--- colCouple_pattern ---" << endl;
960 cout << "M=" << M << ", N=" << N << endl;
961 for (size_t i=0; i<ptr.size(); i++) {
962 cout << "ptr[" << i << "]=" << ptr[i] << endl;
963 }
964 for (size_t i=0; i<index.size(); i++) {
965 cout << "index[" << i << "]=" << index[i] << endl;
966 }
967 */
968
969 // now build the row couple pattern
970 IndexVector ptr2(1,0);
971 IndexVector index2;
972 for (dim_t id=0; id<N; id++) {
973 dim_t n=0;
974 for (dim_t i=0; i<M; i++) {
975 for (dim_t j=ptr[i]; j<ptr[i+1]; j++) {
976 if (index[j]==id) {
977 index2.push_back(i);
978 n++;
979 break;
980 }
981 }
982 }
983 ptr2.push_back(ptr2.back()+n);
984 }
985
986 // paso will manage the memory
987 indexC = MEMALLOC(index2.size(), index_t);
988 ptrC = MEMALLOC(ptr2.size(), index_t);
989 copy(index2.begin(), index2.end(), indexC);
990 copy(ptr2.begin(), ptr2.end(), ptrC);
991 Paso_Pattern *rowCouplePattern=Paso_Pattern_alloc(MATRIX_FORMAT_DEFAULT,
992 N, M, ptrC, indexC);
993
994 /*
995 cout << "--- rowCouple_pattern ---" << endl;
996 cout << "M=" << N << ", N=" << M << endl;
997 for (size_t i=0; i<ptr2.size(); i++) {
998 cout << "ptr[" << i << "]=" << ptr2[i] << endl;
999 }
1000 for (size_t i=0; i<index2.size(); i++) {
1001 cout << "index[" << i << "]=" << index2[i] << endl;
1002 }
1003 */
1004
1005 // allocate paso distribution
1006 Paso_Distribution* distribution = Paso_Distribution_alloc(m_mpiInfo,
1007 const_cast<index_t*>(&m_nodeDistribution[0]), 1, 0);
1008
1009 Paso_SystemMatrixPattern* pattern = Paso_SystemMatrixPattern_alloc(
1010 MATRIX_FORMAT_DEFAULT, distribution, distribution,
1011 mainPattern, colCouplePattern, rowCouplePattern,
1012 connector, connector);
1013 Paso_Pattern_free(mainPattern);
1014 Paso_Pattern_free(colCouplePattern);
1015 Paso_Pattern_free(rowCouplePattern);
1016 Paso_Distribution_free(distribution);
1017 return pattern;
1018 }
1019
1020 void Rectangle::Print_Mesh_Info(const bool full) const
1021 {
1022 RipleyDomain::Print_Mesh_Info(full);
1023 if (full) {
1024 cout << " Id Coordinates" << endl;
1025 cout.precision(15);
1026 cout.setf(ios::scientific, ios::floatfield);
1027 pair<double,double> xdx = getFirstCoordAndSpacing(0);
1028 pair<double,double> ydy = getFirstCoordAndSpacing(1);
1029 for (index_t i=0; i < getNumNodes(); i++) {
1030 cout << " " << setw(5) << m_nodeId[i]
1031 << " " << xdx.first+(i%m_N0)*xdx.second
1032 << " " << ydy.first+(i/m_N0)*ydy.second << endl;
1033 }
1034 }
1035 }
1036
1037 IndexVector Rectangle::getNumNodesPerDim() const
1038 {
1039 IndexVector ret;
1040 ret.push_back(m_N0);
1041 ret.push_back(m_N1);
1042 return ret;
1043 }
1044
1045 IndexVector Rectangle::getNumElementsPerDim() const
1046 {
1047 IndexVector ret;
1048 ret.push_back(m_NE0);
1049 ret.push_back(m_NE1);
1050 return ret;
1051 }
1052
1053 IndexVector Rectangle::getNumFacesPerBoundary() const
1054 {
1055 IndexVector ret(4, 0);
1056 //left
1057 if (m_offset0==0)
1058 ret[0]=m_NE1;
1059 //right
1060 if (m_mpiInfo->rank%m_NX==m_NX-1)
1061 ret[1]=m_NE1;
1062 //bottom
1063 if (m_offset1==0)
1064 ret[2]=m_NE0;
1065 //top
1066 if (m_mpiInfo->rank/m_NX==m_NY-1)
1067 ret[3]=m_NE0;
1068 return ret;
1069 }
1070
1071 pair<double,double> Rectangle::getFirstCoordAndSpacing(dim_t dim) const
1072 {
1073 if (dim==0) {
1074 return pair<double,double>((m_l0*m_offset0)/m_gNE0, m_l0/m_gNE0);
1075 } else if (dim==1) {
1076 return pair<double,double>((m_l1*m_offset1)/m_gNE1, m_l1/m_gNE1);
1077 }
1078 throw RipleyException("getFirstCoordAndSpacing(): invalid argument");
1079 }
1080
1081 //protected
1082 dim_t Rectangle::getNumDOF() const
1083 {
1084 return (m_gNE0+1)/m_NX*(m_gNE1+1)/m_NY;
1085 }
1086
1087 //protected
1088 dim_t Rectangle::getNumFaceElements() const
1089 {
1090 const IndexVector faces = getNumFacesPerBoundary();
1091 dim_t n=0;
1092 for (size_t i=0; i<faces.size(); i++)
1093 n+=faces[i];
1094 return n;
1095 }
1096
1097 //protected
1098 void Rectangle::assembleCoordinates(escript::Data& arg) const
1099 {
1100 escriptDataC x = arg.getDataC();
1101 int numDim = m_numDim;
1102 if (!isDataPointShapeEqual(&x, 1, &numDim))
1103 throw RipleyException("setToX: Invalid Data object shape");
1104 if (!numSamplesEqual(&x, 1, getNumNodes()))
1105 throw RipleyException("setToX: Illegal number of samples in Data object");
1106
1107 pair<double,double> xdx = getFirstCoordAndSpacing(0);
1108 pair<double,double> ydy = getFirstCoordAndSpacing(1);
1109 arg.requireWrite();
1110 #pragma omp parallel for
1111 for (dim_t i1 = 0; i1 < m_N1; i1++) {
1112 for (dim_t i0 = 0; i0 < m_N0; i0++) {
1113 double* point = arg.getSampleDataRW(i0+m_N0*i1);
1114 point[0] = xdx.first+i0*xdx.second;
1115 point[1] = ydy.first+i1*ydy.second;
1116 }
1117 }
1118 }
1119
1120 //private
1121 void Rectangle::populateSampleIds()
1122 {
1123 // identifiers are ordered from left to right, bottom to top globablly.
1124
1125 // build node distribution vector first.
1126 // rank i owns m_nodeDistribution[i+1]-nodeDistribution[i] nodes
1127 m_nodeDistribution.assign(m_mpiInfo->size+1, 0);
1128 const dim_t numDOF=getNumDOF();
1129 for (dim_t k=1; k<m_mpiInfo->size; k++) {
1130 m_nodeDistribution[k]=k*numDOF;
1131 }
1132 m_nodeDistribution[m_mpiInfo->size]=getNumDataPointsGlobal();
1133 m_nodeId.resize(getNumNodes());
1134 m_dofId.resize(numDOF);
1135 m_elementId.resize(getNumElements());
1136 m_faceId.resize(getNumFaceElements());
1137
1138 #pragma omp parallel
1139 {
1140 // nodes
1141 #pragma omp for nowait
1142 for (dim_t i1=0; i1<m_N1; i1++) {
1143 for (dim_t i0=0; i0<m_N0; i0++) {
1144 m_nodeId[i0+i1*m_N0] = (m_offset1+i1)*(m_gNE0+1)+m_offset0+i0;
1145 }
1146 }
1147
1148 // degrees of freedom
1149 #pragma omp for nowait
1150 for (dim_t k=0; k<numDOF; k++)
1151 m_dofId[k] = m_nodeDistribution[m_mpiInfo->rank]+k;
1152
1153 // elements
1154 #pragma omp for nowait
1155 for (dim_t k=0; k<getNumElements(); k++)
1156 m_elementId[k]=k;
1157
1158 // face elements
1159 #pragma omp for
1160 for (dim_t k=0; k<getNumFaceElements(); k++)
1161 m_faceId[k]=k;
1162 } // end parallel section
1163
1164 m_nodeTags.assign(getNumNodes(), 0);
1165 updateTagsInUse(Nodes);
1166
1167 m_elementTags.assign(getNumElements(), 0);
1168 updateTagsInUse(Elements);
1169
1170 // generate face offset vector and set face tags
1171 const IndexVector facesPerEdge = getNumFacesPerBoundary();
1172 const index_t LEFT=1, RIGHT=2, BOTTOM=10, TOP=20;
1173 const index_t faceTag[] = { LEFT, RIGHT, BOTTOM, TOP };
1174 m_faceOffset.assign(facesPerEdge.size(), -1);
1175 m_faceTags.clear();
1176 index_t offset=0;
1177 for (size_t i=0; i<facesPerEdge.size(); i++) {
1178 if (facesPerEdge[i]>0) {
1179 m_faceOffset[i]=offset;
1180 offset+=facesPerEdge[i];
1181 m_faceTags.insert(m_faceTags.end(), facesPerEdge[i], faceTag[i]);
1182 }
1183 }
1184 setTagMap("left", LEFT);
1185 setTagMap("right", RIGHT);
1186 setTagMap("bottom", BOTTOM);
1187 setTagMap("top", TOP);
1188 updateTagsInUse(FaceElements);
1189 }
1190
1191 //private
1192 int Rectangle::insertNeighbours(IndexVector& index, index_t node) const
1193 {
1194 const index_t myN0 = (m_gNE0+1)/m_NX;
1195 const index_t myN1 = (m_gNE1+1)/m_NY;
1196 const int x=node%myN0;
1197 const int y=node/myN0;
1198 int num=0;
1199 if (y>0) {
1200 if (x>0) {
1201 // bottom-left
1202 index.push_back(node-myN0-1);
1203 num++;
1204 }
1205 // bottom
1206 index.push_back(node-myN0);
1207 num++;
1208 if (x<myN0-1) {
1209 // bottom-right
1210 index.push_back(node-myN0+1);
1211 num++;
1212 }
1213 }
1214 if (x<myN0-1) {
1215 // right
1216 index.push_back(node+1);
1217 num++;
1218 if (y<myN1-1) {
1219 // top-right
1220 index.push_back(node+myN0+1);
1221 num++;
1222 }
1223 }
1224 if (y<myN1-1) {
1225 // top
1226 index.push_back(node+myN0);
1227 num++;
1228 if (x>0) {
1229 // top-left
1230 index.push_back(node+myN0-1);
1231 num++;
1232 }
1233 }
1234 if (x>0) {
1235 // left
1236 index.push_back(node-1);
1237 num++;
1238 }
1239
1240 return num;
1241 }
1242
1243 //protected
1244 void Rectangle::interpolateNodesOnElements(escript::Data& out,
1245 escript::Data& in, bool reduced) const
1246 {
1247 const dim_t numComp = in.getDataPointSize();
1248 if (reduced) {
1249 /*** GENERATOR SNIP_INTERPOLATE_REDUCED_ELEMENTS TOP */
1250 const double c0 = .25;
1251 #pragma omp parallel for
1252 for (index_t k1=0; k1 < m_NE1; ++k1) {
1253 for (index_t k0=0; k0 < m_NE0; ++k0) {
1254 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0));
1255 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,k1, m_N0));
1256 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,k1, m_N0));
1257 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,k1+1, m_N0));
1258 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
1259 for (index_t i=0; i < numComp; ++i) {
1260 o[INDEX2(i,numComp,0)] = c0*(f_00[i] + f_01[i] + f_10[i] + f_11[i]);
1261 } /* end of component loop i */
1262 } /* end of k0 loop */
1263 } /* end of k1 loop */
1264 /* GENERATOR SNIP_INTERPOLATE_REDUCED_ELEMENTS BOTTOM */
1265 } else {
1266 /*** GENERATOR SNIP_INTERPOLATE_ELEMENTS TOP */
1267 const double c0 = .16666666666666666667;
1268 const double c1 = .044658198738520451079;
1269 const double c2 = .62200846792814621559;
1270 #pragma omp parallel for
1271 for (index_t k1=0; k1 < m_NE1; ++k1) {
1272 for (index_t k0=0; k0 < m_NE0; ++k0) {
1273 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,k1, m_N0));
1274 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0));
1275 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,k1+1, m_N0));
1276 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,k1, m_N0));
1277 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
1278 for (index_t i=0; i < numComp; ++i) {
1279 o[INDEX2(i,numComp,0)] = f_00[i]*c2 + f_11[i]*c1 + c0*(f_01[i] + f_10[i]);
1280 o[INDEX2(i,numComp,1)] = f_01[i]*c1 + f_10[i]*c2 + c0*(f_00[i] + f_11[i]);
1281 o[INDEX2(i,numComp,2)] = f_01[i]*c2 + f_10[i]*c1 + c0*(f_00[i] + f_11[i]);
1282 o[INDEX2(i,numComp,3)] = f_00[i]*c1 + f_11[i]*c2 + c0*(f_01[i] + f_10[i]);
1283 } /* end of component loop i */
1284 } /* end of k0 loop */
1285 } /* end of k1 loop */
1286 /* GENERATOR SNIP_INTERPOLATE_ELEMENTS BOTTOM */
1287 }
1288 }
1289
1290 //protected
1291 void Rectangle::interpolateNodesOnFaces(escript::Data& out, escript::Data& in,
1292 bool reduced) const
1293 {
1294 const dim_t numComp = in.getDataPointSize();
1295 if (reduced) {
1296 const double c0 = .5;
1297 #pragma omp parallel
1298 {
1299 /*** GENERATOR SNIP_INTERPOLATE_REDUCED_FACES TOP */
1300 if (m_faceOffset[0] > -1) {
1301 #pragma omp for nowait
1302 for (index_t k1=0; k1 < m_NE1; ++k1) {
1303 const register double* f_00 = in.getSampleDataRO(INDEX2(0,k1, m_N0));
1304 const register double* f_01 = in.getSampleDataRO(INDEX2(0,k1+1, m_N0));
1305 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
1306 for (index_t i=0; i < numComp; ++i) {
1307 o[INDEX2(i,numComp,0)] = c0*(f_00[i] + f_01[i]);
1308 } /* end of component loop i */
1309 } /* end of k1 loop */
1310 } /* end of face 0 */
1311 if (m_faceOffset[1] > -1) {
1312 #pragma omp for nowait
1313 for (index_t k1=0; k1 < m_NE1; ++k1) {
1314 const register double* f_10 = in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0));
1315 const register double* f_11 = in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0));
1316 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
1317 for (index_t i=0; i < numComp; ++i) {
1318 o[INDEX2(i,numComp,0)] = c0*(f_10[i] + f_11[i]);
1319 } /* end of component loop i */
1320 } /* end of k1 loop */
1321 } /* end of face 1 */
1322 if (m_faceOffset[2] > -1) {
1323 #pragma omp for nowait
1324 for (index_t k0=0; k0 < m_NE0; ++k0) {
1325 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,0, m_N0));
1326 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,0, m_N0));
1327 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
1328 for (index_t i=0; i < numComp; ++i) {
1329 o[INDEX2(i,numComp,0)] = c0*(f_00[i] + f_10[i]);
1330 } /* end of component loop i */
1331 } /* end of k0 loop */
1332 } /* end of face 2 */
1333 if (m_faceOffset[3] > -1) {
1334 #pragma omp for nowait
1335 for (index_t k0=0; k0 < m_NE0; ++k0) {
1336 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0));
1337 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0));
1338 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
1339 for (index_t i=0; i < numComp; ++i) {
1340 o[INDEX2(i,numComp,0)] = c0*(f_01[i] + f_11[i]);
1341 } /* end of component loop i */
1342 } /* end of k0 loop */
1343 } /* end of face 3 */
1344 /* GENERATOR SNIP_INTERPOLATE_REDUCED_FACES BOTTOM */
1345 } // end of parallel section
1346 } else {
1347 const double c0 = 0.21132486540518711775;
1348 const double c1 = 0.78867513459481288225;
1349 #pragma omp parallel
1350 {
1351 /*** GENERATOR SNIP_INTERPOLATE_FACES TOP */
1352 if (m_faceOffset[0] > -1) {
1353 #pragma omp for nowait
1354 for (index_t k1=0; k1 < m_NE1; ++k1) {
1355 const register double* f_01 = in.getSampleDataRO(INDEX2(0,k1+1, m_N0));
1356 const register double* f_00 = in.getSampleDataRO(INDEX2(0,k1, m_N0));
1357 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
1358 for (index_t i=0; i < numComp; ++i) {
1359 o[INDEX2(i,numComp,0)] = f_00[i]*c1 + f_01[i]*c0;
1360 o[INDEX2(i,numComp,1)] = f_00[i]*c0 + f_01[i]*c1;
1361 } /* end of component loop i */
1362 } /* end of k1 loop */
1363 } /* end of face 0 */
1364 if (m_faceOffset[1] > -1) {
1365 #pragma omp for nowait
1366 for (index_t k1=0; k1 < m_NE1; ++k1) {
1367 const register double* f_10 = in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0));
1368 const register double* f_11 = in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0));
1369 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
1370 for (index_t i=0; i < numComp; ++i) {
1371 o[INDEX2(i,numComp,0)] = f_10[i]*c1 + f_11[i]*c0;
1372 o[INDEX2(i,numComp,1)] = f_10[i]*c0 + f_11[i]*c1;
1373 } /* end of component loop i */
1374 } /* end of k1 loop */
1375 } /* end of face 1 */
1376 if (m_faceOffset[2] > -1) {
1377 #pragma omp for nowait
1378 for (index_t k0=0; k0 < m_NE0; ++k0) {
1379 const register double* f_10 = in.getSampleDataRO(INDEX2(k0+1,0, m_N0));
1380 const register double* f_00 = in.getSampleDataRO(INDEX2(k0,0, m_N0));
1381 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
1382 for (index_t i=0; i < numComp; ++i) {
1383 o[INDEX2(i,numComp,0)] = f_00[i]*c1 + f_10[i]*c0;
1384 o[INDEX2(i,numComp,1)] = f_00[i]*c0 + f_10[i]*c1;
1385 } /* end of component loop i */
1386 } /* end of k0 loop */
1387 } /* end of face 2 */
1388 if (m_faceOffset[3] > -1) {
1389 #pragma omp for nowait
1390 for (index_t k0=0; k0 < m_NE0; ++k0) {
1391 const register double* f_11 = in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0));
1392 const register double* f_01 = in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0));
1393 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
1394 for (index_t i=0; i < numComp; ++i) {
1395 o[INDEX2(i,numComp,0)] = f_01[i]*c1 + f_11[i]*c0;
1396 o[INDEX2(i,numComp,1)] = f_01[i]*c0 + f_11[i]*c1;
1397 } /* end of component loop i */
1398 } /* end of k0 loop */
1399 } /* end of face 3 */
1400 /* GENERATOR SNIP_INTERPOLATE_FACES BOTTOM */
1401 } // end of parallel section
1402 }
1403 }
1404
1405 //protected
1406 void Rectangle::nodesToDOF(escript::Data& out, escript::Data& in) const
1407 {
1408 const dim_t numComp = in.getDataPointSize();
1409 out.requireWrite();
1410
1411 const index_t left = (m_offset0==0 ? 0 : 1);
1412 const index_t bottom = (m_offset1==0 ? 0 : 1);
1413 const index_t right = (m_mpiInfo->rank%m_NX==m_NX-1 ? m_N0 : m_N0-1);
1414 const index_t top = (m_mpiInfo->rank/m_NX==m_NY-1 ? m_N1 : m_N1-1);
1415 index_t n=0;
1416 for (index_t i=bottom; i<top; i++) {
1417 for (index_t j=left; j<right; j++, n++) {
1418 const double* src=in.getSampleDataRO(j+i*m_N0);
1419 copy(src, src+numComp, out.getSampleDataRW(n));
1420 }
1421 }
1422 }
1423
1424 //protected
1425 void Rectangle::assemblePDESingle(Paso_SystemMatrix* mat,
1426 escript::Data& rhs, const escript::Data& A, const escript::Data& B,
1427 const escript::Data& C, const escript::Data& D,
1428 const escript::Data& X, const escript::Data& Y,
1429 const escript::Data& d, const escript::Data& y) const
1430 {
1431 const double h0 = m_l0/m_gNE0;
1432 const double h1 = m_l1/m_gNE1;
1433 /* GENERATOR SNIP_PDE_SINGLE_PRE TOP */
1434 const double w0 = -0.1555021169820365539*h1/h0;
1435 const double w1 = 0.041666666666666666667;
1436 const double w10 = -0.041666666666666666667*h0/h1;
1437 const double w11 = 0.1555021169820365539*h1/h0;
1438 const double w12 = 0.1555021169820365539*h0/h1;
1439 const double w13 = 0.01116454968463011277*h0/h1;
1440 const double w14 = 0.01116454968463011277*h1/h0;
1441 const double w15 = 0.041666666666666666667*h1/h0;
1442 const double w16 = -0.01116454968463011277*h0/h1;
1443 const double w17 = -0.1555021169820365539*h0/h1;
1444 const double w18 = -0.33333333333333333333*h1/h0;
1445 const double w19 = 0.25000000000000000000;
1446 const double w2 = -0.15550211698203655390;
1447 const double w20 = -0.25000000000000000000;
1448 const double w21 = 0.16666666666666666667*h0/h1;
1449 const double w22 = -0.16666666666666666667*h1/h0;
1450 const double w23 = -0.16666666666666666667*h0/h1;
1451 const double w24 = 0.33333333333333333333*h1/h0;
1452 const double w25 = 0.33333333333333333333*h0/h1;
1453 const double w26 = 0.16666666666666666667*h1/h0;
1454 const double w27 = -0.33333333333333333333*h0/h1;
1455 const double w28 = -0.032861463941450536761*h1;
1456 const double w29 = -0.032861463941450536761*h0;
1457 const double w3 = 0.041666666666666666667*h0/h1;
1458 const double w30 = -0.12264065304058601714*h1;
1459 const double w31 = -0.0023593469594139828636*h1;
1460 const double w32 = -0.008805202725216129906*h0;
1461 const double w33 = -0.008805202725216129906*h1;
1462 const double w34 = 0.032861463941450536761*h1;
1463 const double w35 = 0.008805202725216129906*h1;
1464 const double w36 = 0.008805202725216129906*h0;
1465 const double w37 = 0.0023593469594139828636*h1;
1466 const double w38 = 0.12264065304058601714*h1;
1467 const double w39 = 0.032861463941450536761*h0;
1468 const double w4 = 0.15550211698203655390;
1469 const double w40 = -0.12264065304058601714*h0;
1470 const double w41 = -0.0023593469594139828636*h0;
1471 const double w42 = 0.0023593469594139828636*h0;
1472 const double w43 = 0.12264065304058601714*h0;
1473 const double w44 = -0.16666666666666666667*h1;
1474 const double w45 = -0.083333333333333333333*h0;
1475 const double w46 = 0.083333333333333333333*h1;
1476 const double w47 = 0.16666666666666666667*h1;
1477 const double w48 = 0.083333333333333333333*h0;
1478 const double w49 = -0.16666666666666666667*h0;
1479 const double w5 = -0.041666666666666666667;
1480 const double w50 = 0.16666666666666666667*h0;
1481 const double w51 = -0.083333333333333333333*h1;
1482 const double w52 = 0.025917019497006092316*h0*h1;
1483 const double w53 = 0.0018607582807716854616*h0*h1;
1484 const double w54 = 0.0069444444444444444444*h0*h1;
1485 const double w55 = 0.09672363354357992482*h0*h1;
1486 const double w56 = 0.00049858867864229740201*h0*h1;
1487 const double w57 = 0.055555555555555555556*h0*h1;
1488 const double w58 = 0.027777777777777777778*h0*h1;
1489 const double w59 = 0.11111111111111111111*h0*h1;
1490 const double w6 = -0.01116454968463011277*h1/h0;
1491 const double w60 = -0.19716878364870322056*h1;
1492 const double w61 = -0.19716878364870322056*h0;
1493 const double w62 = -0.052831216351296779436*h0;
1494 const double w63 = -0.052831216351296779436*h1;
1495 const double w64 = 0.19716878364870322056*h1;
1496 const double w65 = 0.052831216351296779436*h1;
1497 const double w66 = 0.19716878364870322056*h0;
1498 const double w67 = 0.052831216351296779436*h0;
1499 const double w68 = -0.5*h1;
1500 const double w69 = -0.5*h0;
1501 const double w7 = 0.011164549684630112770;
1502 const double w70 = 0.5*h1;
1503 const double w71 = 0.5*h0;
1504 const double w72 = 0.1555021169820365539*h0*h1;
1505 const double w73 = 0.041666666666666666667*h0*h1;
1506 const double w74 = 0.01116454968463011277*h0*h1;
1507 const double w75 = 0.25*h0*h1;
1508 const double w8 = -0.011164549684630112770;
1509 const double w9 = -0.041666666666666666667*h1/h0;
1510 /* GENERATOR SNIP_PDE_SINGLE_PRE BOTTOM */
1511
1512 rhs.requireWrite();
1513 #pragma omp parallel
1514 {
1515 for (index_t k1_0=0; k1_0<2; k1_0++) { // coloring
1516 #pragma omp for
1517 for (index_t k1=k1_0; k1<m_NE1; k1+=2) {
1518 for (index_t k0=0; k0<m_NE0; ++k0) {
1519 bool add_EM_S=false;
1520 bool add_EM_F=false;
1521 vector<double> EM_S(4*4, 0);
1522 vector<double> EM_F(4, 0);
1523 const index_t e = k0 + m_NE0*k1;
1524 /* GENERATOR SNIP_PDE_SINGLE TOP */
1525 ///////////////
1526 // process A //
1527 ///////////////
1528 if (!A.isEmpty()) {
1529 add_EM_S=true;
1530 const double* A_p=const_cast<escript::Data*>(&A)->getSampleDataRO(e);
1531 if (A.actsExpanded()) {
1532 const register double A_00_0 = A_p[INDEX3(0,0,0,2,2)];
1533 const register double A_01_0 = A_p[INDEX3(0,1,0,2,2)];
1534 const register double A_10_0 = A_p[INDEX3(1,0,0,2,2)];
1535 const register double A_11_0 = A_p[INDEX3(1,1,0,2,2)];
1536 const register double A_00_1 = A_p[INDEX3(0,0,1,2,2)];
1537 const register double A_01_1 = A_p[INDEX3(0,1,1,2,2)];
1538 const register double A_10_1 = A_p[INDEX3(1,0,1,2,2)];
1539 const register double A_11_1 = A_p[INDEX3(1,1,1,2,2)];
1540 const register double A_00_2 = A_p[INDEX3(0,0,2,2,2)];
1541 const register double A_01_2 = A_p[INDEX3(0,1,2,2,2)];
1542 const register double A_10_2 = A_p[INDEX3(1,0,2,2,2)];
1543 const register double A_11_2 = A_p[INDEX3(1,1,2,2,2)];
1544 const register double A_00_3 = A_p[INDEX3(0,0,3,2,2)];
1545 const register double A_01_3 = A_p[INDEX3(0,1,3,2,2)];
1546 const register double A_10_3 = A_p[INDEX3(1,0,3,2,2)];
1547 const register double A_11_3 = A_p[INDEX3(1,1,3,2,2)];
1548 const register double tmp4_0 = A_10_1 + A_10_2;
1549 const register double tmp12_0 = A_11_0 + A_11_2;
1550 const register double tmp2_0 = A_11_0 + A_11_1 + A_11_2 + A_11_3;
1551 const register double tmp10_0 = A_01_3 + A_10_3;
1552 const register double tmp14_0 = A_01_0 + A_01_3 + A_10_0 + A_10_3;
1553 const register double tmp0_0 = A_01_0 + A_01_3;
1554 const register double tmp13_0 = A_01_2 + A_10_1;
1555 const register double tmp3_0 = A_00_2 + A_00_3;
1556 const register double tmp11_0 = A_11_1 + A_11_3;
1557 const register double tmp18_0 = A_01_1 + A_10_1;
1558 const register double tmp1_0 = A_00_0 + A_00_1;
1559 const register double tmp15_0 = A_01_1 + A_10_2;
1560 const register double tmp5_0 = A_00_0 + A_00_1 + A_00_2 + A_00_3;
1561 const register double tmp16_0 = A_10_0 + A_10_3;
1562 const register double tmp6_0 = A_01_3 + A_10_0;
1563 const register double tmp17_0 = A_01_1 + A_01_2;
1564 const register double tmp9_0 = A_01_0 + A_10_0;
1565 const register double tmp7_0 = A_01_0 + A_10_3;
1566 const register double tmp8_0 = A_01_1 + A_01_2 + A_10_1 + A_10_2;
1567 const register double tmp19_0 = A_01_2 + A_10_2;
1568 const register double tmp14_1 = A_10_0*w8;
1569 const register double tmp23_1 = tmp3_0*w14;
1570 const register double tmp35_1 = A_01_0*w8;
1571 const register double tmp54_1 = tmp13_0*w8;
1572 const register double tmp20_1 = tmp9_0*w4;
1573 const register double tmp25_1 = tmp12_0*w12;
1574 const register double tmp2_1 = A_01_1*w4;
1575 const register double tmp44_1 = tmp7_0*w7;
1576 const register double tmp26_1 = tmp10_0*w4;
1577 const register double tmp52_1 = tmp18_0*w8;
1578 const register double tmp48_1 = A_10_1*w7;
1579 const register double tmp46_1 = A_01_3*w8;
1580 const register double tmp50_1 = A_01_0*w2;
1581 const register double tmp8_1 = tmp4_0*w5;
1582 const register double tmp56_1 = tmp19_0*w8;
1583 const register double tmp9_1 = tmp2_0*w10;
1584 const register double tmp19_1 = A_10_3*w2;
1585 const register double tmp47_1 = A_10_2*w4;
1586 const register double tmp16_1 = tmp3_0*w0;
1587 const register double tmp18_1 = tmp1_0*w6;
1588 const register double tmp31_1 = tmp11_0*w12;
1589 const register double tmp55_1 = tmp15_0*w2;
1590 const register double tmp39_1 = A_10_2*w7;
1591 const register double tmp11_1 = tmp6_0*w7;
1592 const register double tmp40_1 = tmp11_0*w17;
1593 const register double tmp34_1 = tmp15_0*w8;
1594 const register double tmp33_1 = tmp14_0*w5;
1595 const register double tmp24_1 = tmp11_0*w13;
1596 const register double tmp3_1 = tmp1_0*w0;
1597 const register double tmp5_1 = tmp2_0*w3;
1598 const register double tmp43_1 = tmp17_0*w5;
1599 const register double tmp15_1 = A_01_2*w4;
1600 const register double tmp53_1 = tmp19_0*w2;
1601 const register double tmp27_1 = tmp3_0*w11;
1602 const register double tmp32_1 = tmp13_0*w2;
1603 const register double tmp10_1 = tmp5_0*w9;
1604 const register double tmp37_1 = A_10_1*w4;
1605 const register double tmp38_1 = tmp5_0*w15;
1606 const register double tmp17_1 = A_01_1*w7;
1607 const register double tmp12_1 = tmp7_0*w4;
1608 const register double tmp22_1 = tmp10_0*w7;
1609 const register double tmp57_1 = tmp18_0*w2;
1610 const register double tmp28_1 = tmp9_0*w7;
1611 const register double tmp29_1 = tmp1_0*w14;
1612 const register double tmp51_1 = tmp11_0*w16;
1613 const register double tmp42_1 = tmp12_0*w16;
1614 const register double tmp49_1 = tmp12_0*w17;
1615 const register double tmp21_1 = tmp1_0*w11;
1616 const register double tmp1_1 = tmp0_0*w1;
1617 const register double tmp45_1 = tmp6_0*w4;
1618 const register double tmp7_1 = A_10_0*w2;
1619 const register double tmp6_1 = tmp3_0*w6;
1620 const register double tmp13_1 = tmp8_0*w1;
1621 const register double tmp36_1 = tmp16_0*w1;
1622 const register double tmp41_1 = A_01_3*w2;
1623 const register double tmp30_1 = tmp12_0*w13;
1624 const register double tmp4_1 = A_01_2*w7;
1625 const register double tmp0_1 = A_10_3*w8;
1626 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1 + tmp4_1 + tmp5_1 + tmp6_1 + tmp7_1 + tmp8_1;
1627 EM_S[INDEX2(1,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1 + tmp9_1;
1628 EM_S[INDEX2(3,2,4)]+=tmp14_1 + tmp15_1 + tmp16_1 + tmp17_1 + tmp18_1 + tmp19_1 + tmp1_1 + tmp5_1 + tmp8_1;
1629 EM_S[INDEX2(0,0,4)]+=tmp13_1 + tmp20_1 + tmp21_1 + tmp22_1 + tmp23_1 + tmp24_1 + tmp25_1;
1630 EM_S[INDEX2(3,3,4)]+=tmp13_1 + tmp26_1 + tmp27_1 + tmp28_1 + tmp29_1 + tmp30_1 + tmp31_1;
1631 EM_S[INDEX2(3,0,4)]+=tmp10_1 + tmp32_1 + tmp33_1 + tmp34_1 + tmp9_1;
1632 EM_S[INDEX2(3,1,4)]+=tmp35_1 + tmp36_1 + tmp37_1 + tmp38_1 + tmp39_1 + tmp40_1 + tmp41_1 + tmp42_1 + tmp43_1;
1633 EM_S[INDEX2(2,1,4)]+=tmp10_1 + tmp13_1 + tmp44_1 + tmp45_1 + tmp9_1;
1634 EM_S[INDEX2(0,2,4)]+=tmp36_1 + tmp38_1 + tmp43_1 + tmp46_1 + tmp47_1 + tmp48_1 + tmp49_1 + tmp50_1 + tmp51_1;
1635 EM_S[INDEX2(2,0,4)]+=tmp0_1 + tmp15_1 + tmp17_1 + tmp1_1 + tmp38_1 + tmp49_1 + tmp51_1 + tmp7_1 + tmp8_1;
1636 EM_S[INDEX2(1,3,4)]+=tmp14_1 + tmp19_1 + tmp1_1 + tmp2_1 + tmp38_1 + tmp40_1 + tmp42_1 + tmp4_1 + tmp8_1;
1637 EM_S[INDEX2(2,3,4)]+=tmp16_1 + tmp18_1 + tmp35_1 + tmp36_1 + tmp41_1 + tmp43_1 + tmp47_1 + tmp48_1 + tmp5_1;
1638 EM_S[INDEX2(2,2,4)]+=tmp24_1 + tmp25_1 + tmp27_1 + tmp29_1 + tmp33_1 + tmp52_1 + tmp53_1;
1639 EM_S[INDEX2(1,0,4)]+=tmp36_1 + tmp37_1 + tmp39_1 + tmp3_1 + tmp43_1 + tmp46_1 + tmp50_1 + tmp5_1 + tmp6_1;
1640 EM_S[INDEX2(0,3,4)]+=tmp10_1 + tmp33_1 + tmp54_1 + tmp55_1 + tmp9_1;
1641 EM_S[INDEX2(1,1,4)]+=tmp21_1 + tmp23_1 + tmp30_1 + tmp31_1 + tmp33_1 + tmp56_1 + tmp57_1;
1642 } else { /* constant data */
1643 const register double A_00 = A_p[INDEX2(0,0,2)];
1644 const register double A_01 = A_p[INDEX2(0,1,2)];
1645 const register double A_10 = A_p[INDEX2(1,0,2)];
1646 const register double A_11 = A_p[INDEX2(1,1,2)];
1647 const register double tmp0_0 = A_01 + A_10;
1648 const register double tmp0_1 = A_00*w18;
1649 const register double tmp10_1 = A_01*w20;
1650 const register double tmp12_1 = A_00*w26;
1651 const register double tmp4_1 = A_00*w22;
1652 const register double tmp8_1 = A_00*w24;
1653 const register double tmp13_1 = A_10*w19;
1654 const register double tmp9_1 = tmp0_0*w20;
1655 const register double tmp3_1 = A_11*w21;
1656 const register double tmp11_1 = A_11*w27;
1657 const register double tmp1_1 = A_01*w19;
1658 const register double tmp6_1 = A_11*w23;
1659 const register double tmp7_1 = A_11*w25;
1660 const register double tmp2_1 = A_10*w20;
1661 const register double tmp5_1 = tmp0_0*w19;
1662 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1;
1663 EM_S[INDEX2(1,2,4)]+=tmp4_1 + tmp5_1 + tmp6_1;
1664 EM_S[INDEX2(3,2,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1;
1665 EM_S[INDEX2(0,0,4)]+=tmp5_1 + tmp7_1 + tmp8_1;
1666 EM_S[INDEX2(3,3,4)]+=tmp5_1 + tmp7_1 + tmp8_1;
1667 EM_S[INDEX2(3,0,4)]+=tmp4_1 + tmp6_1 + tmp9_1;
1668 EM_S[INDEX2(3,1,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1;
1669 EM_S[INDEX2(2,1,4)]+=tmp4_1 + tmp5_1 + tmp6_1;
1670 EM_S[INDEX2(0,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1;
1671 EM_S[INDEX2(2,0,4)]+=tmp11_1 + tmp12_1 + tmp1_1 + tmp2_1;
1672 EM_S[INDEX2(1,3,4)]+=tmp11_1 + tmp12_1 + tmp1_1 + tmp2_1;
1673 EM_S[INDEX2(2,3,4)]+=tmp0_1 + tmp10_1 + tmp13_1 + tmp3_1;
1674 EM_S[INDEX2(2,2,4)]+=tmp7_1 + tmp8_1 + tmp9_1;
1675 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp10_1 + tmp13_1 + tmp3_1;
1676 EM_S[INDEX2(0,3,4)]+=tmp4_1 + tmp6_1 + tmp9_1;
1677 EM_S[INDEX2(1,1,4)]+=tmp7_1 + tmp8_1 + tmp9_1;
1678 }
1679 }
1680 ///////////////
1681 // process B //
1682 ///////////////
1683 if (!B.isEmpty()) {
1684 add_EM_S=true;
1685 const double* B_p=const_cast<escript::Data*>(&B)->getSampleDataRO(e);
1686 if (B.actsExpanded()) {
1687 const register double B_0_0 = B_p[INDEX2(0,0,2)];
1688 const register double B_1_0 = B_p[INDEX2(1,0,2)];
1689 const register double B_0_1 = B_p[INDEX2(0,1,2)];
1690 const register double B_1_1 = B_p[INDEX2(1,1,2)];
1691 const register double B_0_2 = B_p[INDEX2(0,2,2)];
1692 const register double B_1_2 = B_p[INDEX2(1,2,2)];
1693 const register double B_0_3 = B_p[INDEX2(0,3,2)];
1694 const register double B_1_3 = B_p[INDEX2(1,3,2)];
1695 const register double tmp3_0 = B_0_0 + B_0_2;
1696 const register double tmp1_0 = B_1_2 + B_1_3;
1697 const register double tmp2_0 = B_0_1 + B_0_3;
1698 const register double tmp0_0 = B_1_0 + B_1_1;
1699 const register double tmp63_1 = B_1_1*w42;
1700 const register double tmp79_1 = B_1_1*w40;
1701 const register double tmp37_1 = tmp3_0*w35;
1702 const register double tmp8_1 = tmp0_0*w32;
1703 const register double tmp71_1 = B_0_1*w34;
1704 const register double tmp19_1 = B_0_3*w31;
1705 const register double tmp15_1 = B_0_3*w34;
1706 const register double tmp9_1 = tmp3_0*w34;
1707 const register double tmp35_1 = B_1_0*w36;
1708 const register double tmp66_1 = B_0_3*w28;
1709 const register double tmp28_1 = B_1_0*w42;
1710 const register double tmp22_1 = B_1_0*w40;
1711 const register double tmp16_1 = B_1_2*w29;
1712 const register double tmp6_1 = tmp2_0*w35;
1713 const register double tmp55_1 = B_1_3*w40;
1714 const register double tmp50_1 = B_1_3*w42;
1715 const register double tmp7_1 = tmp1_0*w29;
1716 const register double tmp1_1 = tmp1_0*w32;
1717 const register double tmp57_1 = B_0_3*w30;
1718 const register double tmp18_1 = B_1_1*w32;
1719 const register double tmp53_1 = B_1_0*w41;
1720 const register double tmp61_1 = B_1_3*w36;
1721 const register double tmp27_1 = B_0_3*w38;
1722 const register double tmp64_1 = B_0_2*w30;
1723 const register double tmp76_1 = B_0_1*w38;
1724 const register double tmp39_1 = tmp2_0*w34;
1725 const register double tmp62_1 = B_0_1*w31;
1726 const register double tmp56_1 = B_0_0*w31;
1727 const register double tmp49_1 = B_1_1*w36;
1728 const register double tmp2_1 = B_0_2*w31;
1729 const register double tmp23_1 = B_0_2*w33;
1730 const register double tmp38_1 = B_1_1*w43;
1731 const register double tmp74_1 = B_1_2*w41;
1732 const register double tmp43_1 = B_1_1*w41;
1733 const register double tmp58_1 = B_0_2*w28;
1734 const register double tmp67_1 = B_0_0*w33;
1735 const register double tmp33_1 = tmp0_0*w39;
1736 const register double tmp4_1 = B_0_0*w28;
1737 const register double tmp20_1 = B_0_0*w30;
1738 const register double tmp13_1 = B_0_2*w38;
1739 const register double tmp65_1 = B_1_2*w43;
1740 const register double tmp0_1 = tmp0_0*w29;
1741 const register double tmp41_1 = tmp3_0*w33;
1742 const register double tmp73_1 = B_0_2*w37;
1743 const register double tmp69_1 = B_0_0*w38;
1744 const register double tmp48_1 = B_1_2*w39;
1745 const register double tmp59_1 = B_0_1*w33;
1746 const register double tmp17_1 = B_1_3*w41;
1747 const register double tmp5_1 = B_0_3*w33;
1748 const register double tmp3_1 = B_0_1*w30;
1749 const register double tmp21_1 = B_0_1*w28;
1750 const register double tmp42_1 = B_1_0*w29;
1751 const register double tmp54_1 = B_1_2*w32;
1752 const register double tmp60_1 = B_1_0*w39;
1753 const register double tmp32_1 = tmp1_0*w36;
1754 const register double tmp10_1 = B_0_1*w37;
1755 const register double tmp14_1 = B_0_0*w35;
1756 const register double tmp29_1 = B_0_1*w35;
1757 const register double tmp26_1 = B_1_2*w36;
1758 const register double tmp30_1 = B_1_3*w43;
1759 const register double tmp70_1 = B_0_2*w35;
1760 const register double tmp34_1 = B_1_3*w39;
1761 const register double tmp51_1 = B_1_0*w43;
1762 const register double tmp31_1 = B_0_2*w34;
1763 const register double tmp45_1 = tmp3_0*w28;
1764 const register double tmp11_1 = tmp1_0*w39;
1765 const register double tmp52_1 = B_1_1*w29;
1766 const register double tmp44_1 = B_1_3*w32;
1767 const register double tmp25_1 = B_1_1*w39;
1768 const register double tmp47_1 = tmp2_0*w33;
1769 const register double tmp72_1 = B_1_3*w29;
1770 const register double tmp40_1 = tmp2_0*w28;
1771 const register double tmp46_1 = B_1_2*w40;
1772 const register double tmp36_1 = B_1_2*w42;
1773 const register double tmp24_1 = B_0_0*w37;
1774 const register double tmp77_1 = B_0_3*w35;
1775 const register double tmp68_1 = B_0_3*w37;
1776 const register double tmp78_1 = B_0_0*w34;
1777 const register double tmp12_1 = tmp0_0*w36;
1778 const register double tmp75_1 = B_1_0*w32;
1779 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1 + tmp4_1 + tmp5_1;
1780 EM_S[INDEX2(1,2,4)]+=tmp6_1 + tmp7_1 + tmp8_1 + tmp9_1;
1781 EM_S[INDEX2(3,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1 + tmp14_1 + tmp15_1;
1782 EM_S[INDEX2(0,0,4)]+=tmp16_1 + tmp17_1 + tmp18_1 + tmp19_1 + tmp20_1 + tmp21_1 + tmp22_1 + tmp23_1;
1783 EM_S[INDEX2(3,3,4)]+=tmp24_1 + tmp25_1 + tmp26_1 + tmp27_1 + tmp28_1 + tmp29_1 + tmp30_1 + tmp31_1;
1784 EM_S[INDEX2(3,0,4)]+=tmp32_1 + tmp33_1 + tmp6_1 + tmp9_1;
1785 EM_S[INDEX2(3,1,4)]+=tmp34_1 + tmp35_1 + tmp36_1 + tmp37_1 + tmp38_1 + tmp39_1;
1786 EM_S[INDEX2(2,1,4)]+=tmp32_1 + tmp33_1 + tmp40_1 + tmp41_1;
1787 EM_S[INDEX2(0,2,4)]+=tmp42_1 + tmp43_1 + tmp44_1 + tmp45_1 + tmp46_1 + tmp47_1;
1788 EM_S[INDEX2(2,0,4)]+=tmp45_1 + tmp47_1 + tmp48_1 + tmp49_1 + tmp50_1 + tmp51_1;
1789 EM_S[INDEX2(1,3,4)]+=tmp37_1 + tmp39_1 + tmp52_1 + tmp53_1 + tmp54_1 + tmp55_1;
1790 EM_S[INDEX2(2,3,4)]+=tmp11_1 + tmp12_1 + tmp56_1 + tmp57_1 + tmp58_1 + tmp59_1;
1791 EM_S[INDEX2(2,2,4)]+=tmp60_1 + tmp61_1 + tmp62_1 + tmp63_1 + tmp64_1 + tmp65_1 + tmp66_1 + tmp67_1;
1792 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp1_1 + tmp68_1 + tmp69_1 + tmp70_1 + tmp71_1;
1793 EM_S[INDEX2(0,3,4)]+=tmp40_1 + tmp41_1 + tmp7_1 + tmp8_1;
1794 EM_S[INDEX2(1,1,4)]+=tmp72_1 + tmp73_1 + tmp74_1 + tmp75_1 + tmp76_1 + tmp77_1 + tmp78_1 + tmp79_1;
1795 } else { /* constant data */
1796 const register double B_0 = B_p[0];
1797 const register double B_1 = B_p[1];
1798 const register double tmp6_1 = B_1*w50;
1799 const register double tmp1_1 = B_1*w45;
1800 const register double tmp5_1 = B_1*w49;
1801 const register double tmp4_1 = B_1*w48;
1802 const register double tmp0_1 = B_0*w44;
1803 const register double tmp2_1 = B_0*w46;
1804 const register double tmp7_1 = B_0*w51;
1805 const register double tmp3_1 = B_0*w47;
1806 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1;
1807 EM_S[INDEX2(1,2,4)]+=tmp1_1 + tmp2_1;
1808 EM_S[INDEX2(3,2,4)]+=tmp3_1 + tmp4_1;
1809 EM_S[INDEX2(0,0,4)]+=tmp0_1 + tmp5_1;
1810 EM_S[INDEX2(3,3,4)]+=tmp3_1 + tmp6_1;
1811 EM_S[INDEX2(3,0,4)]+=tmp2_1 + tmp4_1;
1812 EM_S[INDEX2(3,1,4)]+=tmp2_1 + tmp6_1;
1813 EM_S[INDEX2(2,1,4)]+=tmp4_1 + tmp7_1;
1814 EM_S[INDEX2(0,2,4)]+=tmp5_1 + tmp7_1;
1815 EM_S[INDEX2(2,0,4)]+=tmp6_1 + tmp7_1;
1816 EM_S[INDEX2(1,3,4)]+=tmp2_1 + tmp5_1;
1817 EM_S[INDEX2(2,3,4)]+=tmp0_1 + tmp4_1;
1818 EM_S[INDEX2(2,2,4)]+=tmp0_1 + tmp6_1;
1819 EM_S[INDEX2(1,0,4)]+=tmp1_1 + tmp3_1;
1820 EM_S[INDEX2(0,3,4)]+=tmp1_1 + tmp7_1;
1821 EM_S[INDEX2(1,1,4)]+=tmp3_1 + tmp5_1;
1822 }
1823 }
1824 ///////////////
1825 // process C //
1826 ///////////////
1827 if (!C.isEmpty()) {
1828 add_EM_S=true;
1829 const double* C_p=const_cast<escript::Data*>(&C)->getSampleDataRO(e);
1830 if (C.actsExpanded()) {
1831 const register double C_0_0 = C_p[INDEX2(0,0,2)];
1832 const register double C_1_0 = C_p[INDEX2(1,0,2)];
1833 const register double C_0_1 = C_p[INDEX2(0,1,2)];
1834 const register double C_1_1 = C_p[INDEX2(1,1,2)];
1835 const register double C_0_2 = C_p[INDEX2(0,2,2)];
1836 const register double C_1_2 = C_p[INDEX2(1,2,2)];
1837 const register double C_0_3 = C_p[INDEX2(0,3,2)];
1838 const register double C_1_3 = C_p[INDEX2(1,3,2)];
1839 const register double tmp2_0 = C_0_1 + C_0_3;
1840 const register double tmp1_0 = C_1_2 + C_1_3;
1841 const register double tmp3_0 = C_0_0 + C_0_2;
1842 const register double tmp0_0 = C_1_0 + C_1_1;
1843 const register double tmp64_1 = C_0_2*w30;
1844 const register double tmp14_1 = C_0_2*w28;
1845 const register double tmp19_1 = C_0_3*w31;
1846 const register double tmp22_1 = C_1_0*w40;
1847 const register double tmp37_1 = tmp3_0*w35;
1848 const register double tmp29_1 = C_0_1*w35;
1849 const register double tmp73_1 = C_0_2*w37;
1850 const register double tmp74_1 = C_1_2*w41;
1851 const register double tmp52_1 = C_1_3*w39;
1852 const register double tmp25_1 = C_1_1*w39;
1853 const register double tmp62_1 = C_0_1*w31;
1854 const register double tmp79_1 = C_1_1*w40;
1855 const register double tmp43_1 = C_1_1*w36;
1856 const register double tmp27_1 = C_0_3*w38;
1857 const register double tmp28_1 = C_1_0*w42;
1858 const register double tmp63_1 = C_1_1*w42;
1859 const register double tmp59_1 = C_0_3*w34;
1860 const register double tmp72_1 = C_1_3*w29;
1861 const register double tmp40_1 = tmp2_0*w35;
1862 const register double tmp13_1 = C_0_3*w30;
1863 const register double tmp51_1 = C_1_2*w40;
1864 const register double tmp54_1 = C_1_2*w42;
1865 const register double tmp12_1 = C_0_0*w31;
1866 const register double tmp2_1 = tmp1_0*w32;
1867 const register double tmp68_1 = C_0_2*w31;
1868 const register double tmp75_1 = C_1_0*w32;
1869 const register double tmp49_1 = C_1_1*w41;
1870 const register double tmp4_1 = C_0_2*w35;
1871 const register double tmp66_1 = C_0_3*w28;
1872 const register double tmp56_1 = C_0_1*w37;
1873 const register double tmp5_1 = C_0_1*w34;
1874 const register double tmp38_1 = tmp2_0*w34;
1875 const register double tmp76_1 = C_0_1*w38;
1876 const register double tmp21_1 = C_0_1*w28;
1877 const register double tmp69_1 = C_0_1*w30;
1878 const register double tmp53_1 = C_1_0*w36;
1879 const register double tmp42_1 = C_1_2*w39;
1880 const register double tmp32_1 = tmp1_0*w29;
1881 const register double tmp45_1 = C_1_0*w43;
1882 const register double tmp33_1 = tmp0_0*w32;
1883 const register double tmp35_1 = C_1_0*w41;
1884 const register double tmp26_1 = C_1_2*w36;
1885 const register double tmp67_1 = C_0_0*w33;
1886 const register double tmp31_1 = C_0_2*w34;
1887 const register double tmp20_1 = C_0_0*w30;
1888 const register double tmp70_1 = C_0_0*w28;
1889 const register double tmp8_1 = tmp0_0*w39;
1890 const register double tmp30_1 = C_1_3*w43;
1891 const register double tmp0_1 = tmp0_0*w29;
1892 const register double tmp17_1 = C_1_3*w41;
1893 const register double tmp58_1 = C_0_0*w35;
1894 const register double tmp9_1 = tmp3_0*w33;
1895 const register double tmp61_1 = C_1_3*w36;
1896 const register double tmp41_1 = tmp3_0*w34;
1897 const register double tmp50_1 = C_1_3*w32;
1898 const register double tmp18_1 = C_1_1*w32;
1899 const register double tmp6_1 = tmp1_0*w36;
1900 const register double tmp3_1 = C_0_0*w38;
1901 const register double tmp34_1 = C_1_1*w29;
1902 const register double tmp77_1 = C_0_3*w35;
1903 const register double tmp65_1 = C_1_2*w43;
1904 const register double tmp71_1 = C_0_3*w33;
1905 const register double tmp55_1 = C_1_1*w43;
1906 const register double tmp46_1 = tmp3_0*w28;
1907 const register double tmp24_1 = C_0_0*w37;
1908 const register double tmp10_1 = tmp1_0*w39;
1909 const register double tmp48_1 = C_1_0*w29;
1910 const register double tmp15_1 = C_0_1*w33;
1911 const register double tmp36_1 = C_1_2*w32;
1912 const register double tmp60_1 = C_1_0*w39;
1913 const register double tmp47_1 = tmp2_0*w33;
1914 const register double tmp16_1 = C_1_2*w29;
1915 const register double tmp1_1 = C_0_3*w37;
1916 const register double tmp7_1 = tmp2_0*w28;
1917 const register double tmp39_1 = C_1_3*w40;
1918 const register double tmp44_1 = C_1_3*w42;
1919 const register double tmp57_1 = C_0_2*w38;
1920 const register double tmp78_1 = C_0_0*w34;
1921 const register double tmp11_1 = tmp0_0*w36;
1922 const register double tmp23_1 = C_0_2*w33;
1923 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1 + tmp4_1 + tmp5_1;
1924 EM_S[INDEX2(1,2,4)]+=tmp6_1 + tmp7_1 + tmp8_1 + tmp9_1;
1925 EM_S[INDEX2(3,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1 + tmp14_1 + tmp15_1;
1926 EM_S[INDEX2(0,0,4)]+=tmp16_1 + tmp17_1 + tmp18_1 + tmp19_1 + tmp20_1 + tmp21_1 + tmp22_1 + tmp23_1;
1927 EM_S[INDEX2(3,3,4)]+=tmp24_1 + tmp25_1 + tmp26_1 + tmp27_1 + tmp28_1 + tmp29_1 + tmp30_1 + tmp31_1;
1928 EM_S[INDEX2(3,0,4)]+=tmp32_1 + tmp33_1 + tmp7_1 + tmp9_1;
1929 EM_S[INDEX2(3,1,4)]+=tmp34_1 + tmp35_1 + tmp36_1 + tmp37_1 + tmp38_1 + tmp39_1;
1930 EM_S[INDEX2(2,1,4)]+=tmp32_1 + tmp33_1 + tmp40_1 + tmp41_1;
1931 EM_S[INDEX2(0,2,4)]+=tmp42_1 + tmp43_1 + tmp44_1 + tmp45_1 + tmp46_1 + tmp47_1;
1932 EM_S[INDEX2(2,0,4)]+=tmp46_1 + tmp47_1 + tmp48_1 + tmp49_1 + tmp50_1 + tmp51_1;
1933 EM_S[INDEX2(1,3,4)]+=tmp37_1 + tmp38_1 + tmp52_1 + tmp53_1 + tmp54_1 + tmp55_1;
1934 EM_S[INDEX2(2,3,4)]+=tmp10_1 + tmp11_1 + tmp56_1 + tmp57_1 + tmp58_1 + tmp59_1;
1935 EM_S[INDEX2(2,2,4)]+=tmp60_1 + tmp61_1 + tmp62_1 + tmp63_1 + tmp64_1 + tmp65_1 + tmp66_1 + tmp67_1;
1936 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp2_1 + tmp68_1 + tmp69_1 + tmp70_1 + tmp71_1;
1937 EM_S[INDEX2(0,3,4)]+=tmp40_1 + tmp41_1 + tmp6_1 + tmp8_1;
1938 EM_S[INDEX2(1,1,4)]+=tmp72_1 + tmp73_1 + tmp74_1 + tmp75_1 + tmp76_1 + tmp77_1 + tmp78_1 + tmp79_1;
1939 } else { /* constant data */
1940 const register double C_0 = C_p[0];
1941 const register double C_1 = C_p[1];
1942 const register double tmp1_1 = C_1*w45;
1943 const register double tmp3_1 = C_0*w51;
1944 const register double tmp4_1 = C_0*w44;
1945 const register double tmp7_1 = C_0*w46;
1946 const register double tmp5_1 = C_1*w49;
1947 const register double tmp2_1 = C_1*w48;
1948 const register double tmp0_1 = C_0*w47;
1949 const register double tmp6_1 = C_1*w50;
1950 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1;
1951 EM_S[INDEX2(1,2,4)]+=tmp2_1 + tmp3_1;
1952 EM_S[INDEX2(3,2,4)]+=tmp2_1 + tmp4_1;
1953 EM_S[INDEX2(0,0,4)]+=tmp4_1 + tmp5_1;
1954 EM_S[INDEX2(3,3,4)]+=tmp0_1 + tmp6_1;
1955 EM_S[INDEX2(3,0,4)]+=tmp1_1 + tmp3_1;
1956 EM_S[INDEX2(3,1,4)]+=tmp5_1 + tmp7_1;
1957 EM_S[INDEX2(2,1,4)]+=tmp1_1 + tmp7_1;
1958 EM_S[INDEX2(0,2,4)]+=tmp3_1 + tmp6_1;
1959 EM_S[INDEX2(2,0,4)]+=tmp3_1 + tmp5_1;
1960 EM_S[INDEX2(1,3,4)]+=tmp6_1 + tmp7_1;
1961 EM_S[INDEX2(2,3,4)]+=tmp0_1 + tmp2_1;
1962 EM_S[INDEX2(2,2,4)]+=tmp4_1 + tmp6_1;
1963 EM_S[INDEX2(1,0,4)]+=tmp1_1 + tmp4_1;
1964 EM_S[INDEX2(0,3,4)]+=tmp2_1 + tmp7_1;
1965 EM_S[INDEX2(1,1,4)]+=tmp0_1 + tmp5_1;
1966 }
1967 }
1968 ///////////////
1969 // process D //
1970 ///////////////
1971 if (!D.isEmpty()) {
1972 add_EM_S=true;
1973 const double* D_p=const_cast<escript::Data*>(&D)->getSampleDataRO(e);
1974 if (D.actsExpanded()) {
1975 const register double D_0 = D_p[0];
1976 const register double D_1 = D_p[1];
1977 const register double D_2 = D_p[2];
1978 const register double D_3 = D_p[3];
1979 const register double tmp4_0 = D_1 + D_3;
1980 const register double tmp2_0 = D_0 + D_1 + D_2 + D_3;
1981 const register double tmp5_0 = D_0 + D_2;
1982 const register double tmp0_0 = D_0 + D_1;
1983 const register double tmp6_0 = D_0 + D_3;
1984 const register double tmp1_0 = D_2 + D_3;
1985 const register double tmp3_0 = D_1 + D_2;
1986 const register double tmp16_1 = D_1*w56;
1987 const register double tmp14_1 = tmp6_0*w54;
1988 const register double tmp8_1 = D_3*w55;
1989 const register double tmp2_1 = tmp2_0*w54;
1990 const register double tmp12_1 = tmp5_0*w52;
1991 const register double tmp4_1 = tmp0_0*w53;
1992 const register double tmp3_1 = tmp1_0*w52;
1993 const register double tmp13_1 = tmp4_0*w53;
1994 const register double tmp10_1 = tmp4_0*w52;
1995 const register double tmp1_1 = tmp1_0*w53;
1996 const register double tmp7_1 = D_3*w56;
1997 const register double tmp0_1 = tmp0_0*w52;
1998 const register double tmp11_1 = tmp5_0*w53;
1999 const register double tmp9_1 = D_0*w56;
2000 const register double tmp5_1 = tmp3_0*w54;
2001 const register double tmp18_1 = D_2*w56;
2002 const register double tmp17_1 = D_1*w55;
2003 const register double tmp6_1 = D_0*w55;
2004 const register double tmp15_1 = D_2*w55;
2005 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1;
2006 EM_S[INDEX2(1,2,4)]+=tmp2_1;
2007 EM_S[INDEX2(3,2,4)]+=tmp3_1 + tmp4_1;
2008 EM_S[INDEX2(0,0,4)]+=tmp5_1 + tmp6_1 + tmp7_1;
2009 EM_S[INDEX2(3,3,4)]+=tmp5_1 + tmp8_1 + tmp9_1;
2010 EM_S[INDEX2(3,0,4)]+=tmp2_1;
2011 EM_S[INDEX2(3,1,4)]+=tmp10_1 + tmp11_1;
2012 EM_S[INDEX2(2,1,4)]+=tmp2_1;
2013 EM_S[INDEX2(0,2,4)]+=tmp12_1 + tmp13_1;
2014 EM_S[INDEX2(2,0,4)]+=tmp12_1 + tmp13_1;
2015 EM_S[INDEX2(1,3,4)]+=tmp10_1 + tmp11_1;
2016 EM_S[INDEX2(2,3,4)]+=tmp3_1 + tmp4_1;
2017 EM_S[INDEX2(2,2,4)]+=tmp14_1 + tmp15_1 + tmp16_1;
2018 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp1_1;
2019 EM_S[INDEX2(0,3,4)]+=tmp2_1;
2020 EM_S[INDEX2(1,1,4)]+=tmp14_1 + tmp17_1 + tmp18_1;
2021 } else { /* constant data */
2022 const register double D_0 = D_p[0];
2023 const register double tmp2_1 = D_0*w59;
2024 const register double tmp1_1 = D_0*w58;
2025 const register double tmp0_1 = D_0*w57;
2026 EM_S[INDEX2(0,1,4)]+=tmp0_1;
2027 EM_S[INDEX2(1,2,4)]+=tmp1_1;
2028 EM_S[INDEX2(3,2,4)]+=tmp0_1;
2029 EM_S[INDEX2(0,0,4)]+=tmp2_1;
2030 EM_S[INDEX2(3,3,4)]+=tmp2_1;
2031 EM_S[INDEX2(3,0,4)]+=tmp1_1;
2032 EM_S[INDEX2(3,1,4)]+=tmp0_1;
2033 EM_S[INDEX2(2,1,4)]+=tmp1_1;
2034 EM_S[INDEX2(0,2,4)]+=tmp0_1;
2035 EM_S[INDEX2(2,0,4)]+=tmp0_1;
2036 EM_S[INDEX2(1,3,4)]+=tmp0_1;
2037 EM_S[INDEX2(2,3,4)]+=tmp0_1;
2038 EM_S[INDEX2(2,2,4)]+=tmp2_1;
2039 EM_S[INDEX2(1,0,4)]+=tmp0_1;
2040 EM_S[INDEX2(0,3,4)]+=tmp1_1;
2041 EM_S[INDEX2(1,1,4)]+=tmp2_1;
2042 }
2043 }
2044 ///////////////
2045 // process X //
2046 ///////////////
2047 if (!X.isEmpty()) {
2048 add_EM_F=true;
2049 const double* X_p=const_cast<escript::Data*>(&X)->getSampleDataRO(e);
2050 if (X.actsExpanded()) {
2051 const register double X_0_0 = X_p[INDEX2(0,0,2)];
2052 const register double X_1_0 = X_p[INDEX2(1,0,2)];
2053 const register double X_0_1 = X_p[INDEX2(0,1,2)];
2054 const register double X_1_1 = X_p[INDEX2(1,1,2)];
2055 const register double X_0_2 = X_p[INDEX2(0,2,2)];
2056 const register double X_1_2 = X_p[INDEX2(1,2,2)];
2057 const register double X_0_3 = X_p[INDEX2(0,3,2)];
2058 const register double X_1_3 = X_p[INDEX2(1,3,2)];
2059 const register double tmp1_0 = X_1_1 + X_1_3;
2060 const register double tmp3_0 = X_0_0 + X_0_1;
2061 const register double tmp2_0 = X_1_0 + X_1_2;
2062 const register double tmp0_0 = X_0_2 + X_0_3;
2063 const register double tmp8_1 = tmp2_0*w66;
2064 const register double tmp5_1 = tmp3_0*w64;
2065 const register double tmp14_1 = tmp0_0*w64;
2066 const register double tmp3_1 = tmp3_0*w60;
2067 const register double tmp9_1 = tmp3_0*w63;
2068 const register double tmp13_1 = tmp3_0*w65;
2069 const register double tmp12_1 = tmp1_0*w66;
2070 const register double tmp10_1 = tmp0_0*w60;
2071 const register double tmp2_1 = tmp2_0*w61;
2072 const register double tmp6_1 = tmp2_0*w62;
2073 const register double tmp4_1 = tmp0_0*w65;
2074 const register double tmp11_1 = tmp1_0*w67;
2075 const register double tmp1_1 = tmp1_0*w62;
2076 const register double tmp7_1 = tmp1_0*w61;
2077 const register double tmp0_1 = tmp0_0*w63;
2078 const register double tmp15_1 = tmp2_0*w67;
2079 EM_F[0]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1;
2080 EM_F[1]+=tmp4_1 + tmp5_1 + tmp6_1 + tmp7_1;
2081 EM_F[2]+=tmp10_1 + tmp11_1 + tmp8_1 + tmp9_1;
2082 EM_F[3]+=tmp12_1 + tmp13_1 + tmp14_1 + tmp15_1;
2083 } else { /* constant data */
2084 const register double X_0 = X_p[0];
2085 const register double X_1 = X_p[1];
2086 const register double tmp3_1 = X_1*w71;
2087 const register double tmp2_1 = X_0*w70;
2088 const register double tmp1_1 = X_0*w68;
2089 const register double tmp0_1 = X_1*w69;
2090 EM_F[0]+=tmp0_1 + tmp1_1;
2091 EM_F[1]+=tmp0_1 + tmp2_1;
2092 EM_F[2]+=tmp1_1 + tmp3_1;
2093 EM_F[3]+=tmp2_1 + tmp3_1;
2094 }
2095 }
2096 ///////////////
2097 // process Y //
2098 ///////////////
2099 if (!Y.isEmpty()) {
2100 add_EM_F=true;
2101 const double* Y_p=const_cast<escript::Data*>(&Y)->getSampleDataRO(e);
2102 if (Y.actsExpanded()) {
2103 const register double Y_0 = Y_p[0];
2104 const register double Y_1 = Y_p[1];
2105 const register double Y_2 = Y_p[2];
2106 const register double Y_3 = Y_p[3];
2107 const register double tmp0_0 = Y_1 + Y_2;
2108 const register double tmp1_0 = Y_0 + Y_3;
2109 const register double tmp9_1 = Y_0*w74;
2110 const register double tmp4_1 = tmp1_0*w73;
2111 const register double tmp5_1 = Y_2*w74;
2112 const register double tmp7_1 = Y_1*w74;
2113 const register double tmp6_1 = Y_2*w72;
2114 const register double tmp2_1 = Y_3*w74;
2115 const register double tmp8_1 = Y_3*w72;
2116 const register double tmp3_1 = Y_1*w72;
2117 const register double tmp0_1 = Y_0*w72;
2118 const register double tmp1_1 = tmp0_0*w73;
2119 EM_F[0]+=tmp0_1 + tmp1_1 + tmp2_1;
2120 EM_F[1]+=tmp3_1 + tmp4_1 + tmp5_1;
2121 EM_F[2]+=tmp4_1 + tmp6_1 + tmp7_1;
2122 EM_F[3]+=tmp1_1 + tmp8_1 + tmp9_1;
2123 } else { /* constant data */
2124 const register double Y_0 = Y_p[0];
2125 const register double tmp0_1 = Y_0*w75;
2126 EM_F[0]+=tmp0_1;
2127 EM_F[1]+=tmp0_1;
2128 EM_F[2]+=tmp0_1;
2129 EM_F[3]+=tmp0_1;
2130 }
2131 }
2132 /* GENERATOR SNIP_PDE_SINGLE BOTTOM */
2133
2134 // add to matrix (if add_EM_S) and RHS (if add_EM_F)
2135 const index_t firstNode=m_N0*k1+k0;
2136 IndexVector rowIndex;
2137 rowIndex.push_back(m_dofMap[firstNode]);
2138 rowIndex.push_back(m_dofMap[firstNode+1]);
2139 rowIndex.push_back(m_dofMap[firstNode+m_N0]);
2140 rowIndex.push_back(m_dofMap[firstNode+m_N0+1]);
2141 if (add_EM_F) {
2142 //cout << "-- AddtoRHS -- " << endl;
2143 double *F_p=rhs.getSampleDataRW(0);
2144 for (index_t i=0; i<4; i++) {
2145 if (rowIndex[i]<getNumDOF()) {
2146 F_p[rowIndex[i]]+=EM_F[i];
2147 //cout << "F[" << rowIndex[i] << "]=" << F_p[rowIndex[i]] << endl;
2148 }
2149 }
2150 //cout << "---"<<endl;
2151 }
2152 if (add_EM_S) {
2153 //cout << "-- AddtoSystem -- " << endl;
2154 addToSystemMatrix(mat, rowIndex, 1, rowIndex, 1, EM_S);
2155 }
2156
2157 } // end k0 loop
2158 } // end k1 loop
2159 } // end of coloring
2160 } // end of parallel region
2161 }
2162
2163 void Rectangle::addToSystemMatrix(Paso_SystemMatrix* mat,
2164 const IndexVector& nodes_Eq, dim_t num_Eq, const IndexVector& nodes_Sol,
2165 dim_t num_Sol, const vector<double>& array) const
2166 {
2167 const dim_t numMyCols = mat->pattern->mainPattern->numInput;
2168 const dim_t numMyRows = mat->pattern->mainPattern->numOutput;
2169
2170 const index_t* mainBlock_ptr = mat->mainBlock->pattern->ptr;
2171 const index_t* mainBlock_index = mat->mainBlock->pattern->index;
2172 double* mainBlock_val = mat->mainBlock->val;
2173 const index_t* col_coupleBlock_ptr = mat->col_coupleBlock->pattern->ptr;
2174 const index_t* col_coupleBlock_index = mat->col_coupleBlock->pattern->index;
2175 double* col_coupleBlock_val = mat->col_coupleBlock->val;
2176 const index_t* row_coupleBlock_ptr = mat->row_coupleBlock->pattern->ptr;
2177 const index_t* row_coupleBlock_index = mat->row_coupleBlock->pattern->index;
2178 double* row_coupleBlock_val = mat->row_coupleBlock->val;
2179
2180 for (dim_t k_Eq = 0; k_Eq < nodes_Eq.size(); ++k_Eq) {
2181 // down columns of array
2182 const dim_t j_Eq = nodes_Eq[k_Eq];
2183 const dim_t i_row = j_Eq;
2184 //printf("row:%d\n", i_row);
2185 // only look at the matrix rows stored on this processor
2186 if (i_row < numMyRows) {
2187 for (dim_t k_Sol = 0; k_Sol < nodes_Sol.size(); ++k_Sol) {
2188 const dim_t i_col = nodes_Sol[k_Sol];
2189 if (i_col < numMyCols) {
2190 for (dim_t k = mainBlock_ptr[i_row]; k < mainBlock_ptr[i_row + 1]; ++k) {
2191 if (mainBlock_index[k] == i_col) {
2192 mainBlock_val[k] += array[INDEX2(k_Eq, k_Sol, nodes_Eq.size())];
2193 break;
2194 }
2195 }