/[escript]/trunk/ripley/src/Rectangle.cpp
ViewVC logotype

Contents of /trunk/ripley/src/Rectangle.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3795 - (show annotations)
Thu Feb 2 05:54:34 2012 UTC (7 years, 9 months ago) by caltinay
File size: 227421 byte(s)
Had to remove 'nowait' directive from PDE boundary conditions since all
4 (2D) / 8 (3D) values are updated the way it is implemented. This caused
race conditions when updating RHS.

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