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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 3913 - (show annotations)
Tue Jun 19 06:21:58 2012 UTC (7 years, 4 months ago) by caltinay
File size: 230292 byte(s)
Fixed misuse of getSampleDataRO with lazy data in ripley.

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=sqrt(hSize*hSize+vSize*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
636 {
637 vector<double> f_00(numComp);
638 vector<double> f_01(numComp);
639 vector<double> f_10(numComp);
640 vector<double> f_11(numComp);
641 #pragma omp for
642 for (index_t k1=0; k1 < m_NE1; ++k1) {
643 for (index_t k0=0; k0 < m_NE0; ++k0) {
644 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,k1, m_N0)), numComp*sizeof(double));
645 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,k1+1, m_N0)), numComp*sizeof(double));
646 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,k1, m_N0)), numComp*sizeof(double));
647 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0)), numComp*sizeof(double));
648 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
649 for (index_t i=0; i < numComp; ++i) {
650 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
651 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
652 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
653 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
654 o[INDEX3(i,0,2,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
655 o[INDEX3(i,1,2,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
656 o[INDEX3(i,0,3,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
657 o[INDEX3(i,1,3,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
658 } // end of component loop i
659 } // end of k0 loop
660 } // end of k1 loop
661 } // end of parallel section
662 } else if (out.getFunctionSpace().getTypeCode() == ReducedElements) {
663 out.requireWrite();
664 #pragma omp parallel
665 {
666 vector<double> f_00(numComp);
667 vector<double> f_01(numComp);
668 vector<double> f_10(numComp);
669 vector<double> f_11(numComp);
670 #pragma omp for
671 for (index_t k1=0; k1 < m_NE1; ++k1) {
672 for (index_t k0=0; k0 < m_NE0; ++k0) {
673 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,k1, m_N0)), numComp*sizeof(double));
674 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,k1+1, m_N0)), numComp*sizeof(double));
675 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,k1, m_N0)), numComp*sizeof(double));
676 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0)), numComp*sizeof(double));
677 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
678 for (index_t i=0; i < numComp; ++i) {
679 o[INDEX3(i,0,0,numComp,2)] = cx5*(f_10[i] + f_11[i]) + cx2*(f_00[i] + f_01[i]);
680 o[INDEX3(i,1,0,numComp,2)] = cy2*(f_00[i] + f_10[i]) + cy5*(f_01[i] + f_11[i]);
681 } // end of component loop i
682 } // end of k0 loop
683 } // end of k1 loop
684 } // end of parallel section
685 } else if (out.getFunctionSpace().getTypeCode() == FaceElements) {
686 out.requireWrite();
687 #pragma omp parallel
688 {
689 vector<double> f_00(numComp);
690 vector<double> f_01(numComp);
691 vector<double> f_10(numComp);
692 vector<double> f_11(numComp);
693 if (m_faceOffset[0] > -1) {
694 #pragma omp for nowait
695 for (index_t k1=0; k1 < m_NE1; ++k1) {
696 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(0,k1, m_N0)), numComp*sizeof(double));
697 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(0,k1+1, m_N0)), numComp*sizeof(double));
698 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(1,k1, m_N0)), numComp*sizeof(double));
699 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(1,k1+1, m_N0)), numComp*sizeof(double));
700 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
701 for (index_t i=0; i < numComp; ++i) {
702 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
703 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy0 + f_01[i]*cy7;
704 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
705 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy0 + f_01[i]*cy7;
706 } // end of component loop i
707 } // end of k1 loop
708 } // end of face 0
709 if (m_faceOffset[1] > -1) {
710 #pragma omp for nowait
711 for (index_t k1=0; k1 < m_NE1; ++k1) {
712 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(m_N0-2,k1, m_N0)), numComp*sizeof(double));
713 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(m_N0-2,k1+1, m_N0)), numComp*sizeof(double));
714 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0)), numComp*sizeof(double));
715 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0)), numComp*sizeof(double));
716 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
717 for (index_t i=0; i < numComp; ++i) {
718 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx1 + f_01[i]*cx3 + f_10[i]*cx6 + f_11[i]*cx4;
719 o[INDEX3(i,1,0,numComp,2)] = f_10[i]*cy0 + f_11[i]*cy7;
720 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx3 + f_01[i]*cx1 + f_10[i]*cx4 + f_11[i]*cx6;
721 o[INDEX3(i,1,1,numComp,2)] = f_10[i]*cy0 + f_11[i]*cy7;
722 } // end of component loop i
723 } // end of k1 loop
724 } // end of face 1
725 if (m_faceOffset[2] > -1) {
726 #pragma omp for nowait
727 for (index_t k0=0; k0 < m_NE0; ++k0) {
728 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,0, m_N0)), numComp*sizeof(double));
729 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,1, m_N0)), numComp*sizeof(double));
730 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,0, m_N0)), numComp*sizeof(double));
731 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,1, m_N0)), numComp*sizeof(double));
732 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
733 for (index_t i=0; i < numComp; ++i) {
734 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx0 + f_10[i]*cx7;
735 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
736 o[INDEX3(i,0,1,numComp,2)] = f_00[i]*cx0 + f_10[i]*cx7;
737 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
738 } // end of component loop i
739 } // end of k0 loop
740 } // end of face 2
741 if (m_faceOffset[3] > -1) {
742 #pragma omp for nowait
743 for (index_t k0=0; k0 < m_NE0; ++k0) {
744 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,m_N1-2, m_N0)), numComp*sizeof(double));
745 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0)), numComp*sizeof(double));
746 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,m_N1-2, m_N0)), numComp*sizeof(double));
747 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0)), numComp*sizeof(double));
748 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
749 for (index_t i=0; i < numComp; ++i) {
750 o[INDEX3(i,0,0,numComp,2)] = f_01[i]*cx0 + f_11[i]*cx7;
751 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy1 + f_01[i]*cy6 + f_10[i]*cy3 + f_11[i]*cy4;
752 o[INDEX3(i,0,1,numComp,2)] = f_01[i]*cx0 + f_11[i]*cx7;
753 o[INDEX3(i,1,1,numComp,2)] = f_00[i]*cy3 + f_01[i]*cy4 + f_10[i]*cy1 + f_11[i]*cy6;
754 } // end of component loop i
755 } // end of k0 loop
756 } // end of face 3
757 } // end of parallel section
758
759 } else if (out.getFunctionSpace().getTypeCode() == ReducedFaceElements) {
760 out.requireWrite();
761 #pragma omp parallel
762 {
763 vector<double> f_00(numComp);
764 vector<double> f_01(numComp);
765 vector<double> f_10(numComp);
766 vector<double> f_11(numComp);
767 if (m_faceOffset[0] > -1) {
768 #pragma omp for nowait
769 for (index_t k1=0; k1 < m_NE1; ++k1) {
770 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(0,k1, m_N0)), numComp*sizeof(double));
771 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(0,k1+1, m_N0)), numComp*sizeof(double));
772 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(1,k1, m_N0)), numComp*sizeof(double));
773 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(1,k1+1, m_N0)), numComp*sizeof(double));
774 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
775 for (index_t i=0; i < numComp; ++i) {
776 o[INDEX3(i,0,0,numComp,2)] = cx5*(f_10[i] + f_11[i]) + cx2*(f_00[i] + f_01[i]);
777 o[INDEX3(i,1,0,numComp,2)] = f_00[i]*cy0 + f_01[i]*cy7;
778 } // end of component loop i
779 } // end of k1 loop
780 } // end of face 0
781 if (m_faceOffset[1] > -1) {
782 #pragma omp for nowait
783 for (index_t k1=0; k1 < m_NE1; ++k1) {
784 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(m_N0-2,k1, m_N0)), numComp*sizeof(double));
785 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(m_N0-2,k1+1, m_N0)), numComp*sizeof(double));
786 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0)), numComp*sizeof(double));
787 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0)), numComp*sizeof(double));
788 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
789 for (index_t i=0; i < numComp; ++i) {
790 o[INDEX3(i,0,0,numComp,2)] = cx5*(f_10[i] + f_11[i]) + cx2*(f_00[i] + f_01[i]);
791 o[INDEX3(i,1,0,numComp,2)] = f_10[i]*cy0 + f_11[i]*cy7;
792 } // end of component loop i
793 } // end of k1 loop
794 } // end of face 1
795 if (m_faceOffset[2] > -1) {
796 #pragma omp for nowait
797 for (index_t k0=0; k0 < m_NE0; ++k0) {
798 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,0, m_N0)), numComp*sizeof(double));
799 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,1, m_N0)), numComp*sizeof(double));
800 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,0, m_N0)), numComp*sizeof(double));
801 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,1, m_N0)), numComp*sizeof(double));
802 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
803 for (index_t i=0; i < numComp; ++i) {
804 o[INDEX3(i,0,0,numComp,2)] = f_00[i]*cx0 + f_10[i]*cx7;
805 o[INDEX3(i,1,0,numComp,2)] = cy2*(f_00[i] + f_10[i]) + cy5*(f_01[i] + f_11[i]);
806 } // end of component loop i
807 } // end of k0 loop
808 } // end of face 2
809 if (m_faceOffset[3] > -1) {
810 #pragma omp for nowait
811 for (index_t k0=0; k0 < m_NE0; ++k0) {
812 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,m_N1-2, m_N0)), numComp*sizeof(double));
813 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0)), numComp*sizeof(double));
814 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,m_N1-2, m_N0)), numComp*sizeof(double));
815 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0)), numComp*sizeof(double));
816 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
817 for (index_t i=0; i < numComp; ++i) {
818 o[INDEX3(i,0,0,numComp,2)] = f_01[i]*cx0 + f_11[i]*cx7;
819 o[INDEX3(i,1,0,numComp,2)] = cy5*(f_01[i] + f_11[i]) + cy2*(f_00[i] + f_10[i]);
820 } // end of component loop i
821 } // end of k0 loop
822 } // end of face 3
823 } // end of parallel section
824 }
825 }
826
827 //protected
828 void Rectangle::assembleIntegrate(vector<double>& integrals, escript::Data& arg) const
829 {
830 const dim_t numComp = arg.getDataPointSize();
831 const double h0 = m_l0/m_gNE0;
832 const double h1 = m_l1/m_gNE1;
833 const index_t left = (m_offset0==0 ? 0 : 1);
834 const index_t bottom = (m_offset1==0 ? 0 : 1);
835 const int fs=arg.getFunctionSpace().getTypeCode();
836 if (fs == Elements && arg.actsExpanded()) {
837 #pragma omp parallel
838 {
839 vector<double> int_local(numComp, 0);
840 const double w = h0*h1/4.;
841 #pragma omp for nowait
842 for (index_t k1 = bottom; k1 < bottom+m_ownNE1; ++k1) {
843 for (index_t k0 = left; k0 < left+m_ownNE0; ++k0) {
844 const double* f = arg.getSampleDataRO(INDEX2(k0, k1, m_NE0));
845 for (index_t i=0; i < numComp; ++i) {
846 const double f0 = f[INDEX2(i,0,numComp)];
847 const double f1 = f[INDEX2(i,1,numComp)];
848 const double f2 = f[INDEX2(i,2,numComp)];
849 const double f3 = f[INDEX2(i,3,numComp)];
850 int_local[i]+=(f0+f1+f2+f3)*w;
851 } // end of component loop i
852 } // end of k0 loop
853 } // end of k1 loop
854 #pragma omp critical
855 for (index_t i=0; i<numComp; i++)
856 integrals[i]+=int_local[i];
857 } // end of parallel section
858
859 } else if (fs==ReducedElements || (fs==Elements && !arg.actsExpanded())) {
860 const double w = h0*h1;
861 #pragma omp parallel
862 {
863 vector<double> int_local(numComp, 0);
864 #pragma omp for nowait
865 for (index_t k1 = bottom; k1 < bottom+m_ownNE1; ++k1) {
866 for (index_t k0 = left; k0 < left+m_ownNE0; ++k0) {
867 const double* f = arg.getSampleDataRO(INDEX2(k0, k1, m_NE0));
868 for (index_t i=0; i < numComp; ++i) {
869 int_local[i]+=f[i]*w;
870 }
871 }
872 }
873 #pragma omp critical
874 for (index_t i=0; i<numComp; i++)
875 integrals[i]+=int_local[i];
876 } // end of parallel section
877
878 } else if (fs == FaceElements && arg.actsExpanded()) {
879 #pragma omp parallel
880 {
881 vector<double> int_local(numComp, 0);
882 const double w0 = h0/2.;
883 const double w1 = h1/2.;
884 if (m_faceOffset[0] > -1) {
885 #pragma omp for nowait
886 for (index_t k1 = bottom; k1 < bottom+m_ownNE1; ++k1) {
887 const double* f = arg.getSampleDataRO(m_faceOffset[0]+k1);
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)*w1;
892 } // end of component loop i
893 } // end of k1 loop
894 }
895
896 if (m_faceOffset[1] > -1) {
897 #pragma omp for nowait
898 for (index_t k1 = bottom; k1 < bottom+m_ownNE1; ++k1) {
899 const double* f = arg.getSampleDataRO(m_faceOffset[1]+k1);
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)*w1;
904 } // end of component loop i
905 } // end of k1 loop
906 }
907
908 if (m_faceOffset[2] > -1) {
909 #pragma omp for nowait
910 for (index_t k0 = left; k0 < left+m_ownNE0; ++k0) {
911 const double* f = arg.getSampleDataRO(m_faceOffset[2]+k0);
912 for (index_t i=0; i < numComp; ++i) {
913 const double f0 = f[INDEX2(i,0,numComp)];
914 const double f1 = f[INDEX2(i,1,numComp)];
915 int_local[i]+=(f0+f1)*w0;
916 } // end of component loop i
917 } // end of k0 loop
918 }
919
920 if (m_faceOffset[3] > -1) {
921 #pragma omp for nowait
922 for (index_t k0 = left; k0 < left+m_ownNE0; ++k0) {
923 const double* f = arg.getSampleDataRO(m_faceOffset[3]+k0);
924 for (index_t i=0; i < numComp; ++i) {
925 const double f0 = f[INDEX2(i,0,numComp)];
926 const double f1 = f[INDEX2(i,1,numComp)];
927 int_local[i]+=(f0+f1)*w0;
928 } // end of component loop i
929 } // end of k0 loop
930 }
931 #pragma omp critical
932 for (index_t i=0; i<numComp; i++)
933 integrals[i]+=int_local[i];
934 } // end of parallel section
935
936 } else if (fs==ReducedFaceElements || (fs==FaceElements && !arg.actsExpanded())) {
937 #pragma omp parallel
938 {
939 vector<double> int_local(numComp, 0);
940 if (m_faceOffset[0] > -1) {
941 #pragma omp for nowait
942 for (index_t k1 = bottom; k1 < bottom+m_ownNE1; ++k1) {
943 const double* f = arg.getSampleDataRO(m_faceOffset[0]+k1);
944 for (index_t i=0; i < numComp; ++i) {
945 int_local[i]+=f[i]*h1;
946 }
947 }
948 }
949
950 if (m_faceOffset[1] > -1) {
951 #pragma omp for nowait
952 for (index_t k1 = bottom; k1 < bottom+m_ownNE1; ++k1) {
953 const double* f = arg.getSampleDataRO(m_faceOffset[1]+k1);
954 for (index_t i=0; i < numComp; ++i) {
955 int_local[i]+=f[i]*h1;
956 }
957 }
958 }
959
960 if (m_faceOffset[2] > -1) {
961 #pragma omp for nowait
962 for (index_t k0 = left; k0 < left+m_ownNE0; ++k0) {
963 const double* f = arg.getSampleDataRO(m_faceOffset[2]+k0);
964 for (index_t i=0; i < numComp; ++i) {
965 int_local[i]+=f[i]*h0;
966 }
967 }
968 }
969
970 if (m_faceOffset[3] > -1) {
971 #pragma omp for nowait
972 for (index_t k0 = left; k0 < left+m_ownNE0; ++k0) {
973 const double* f = arg.getSampleDataRO(m_faceOffset[3]+k0);
974 for (index_t i=0; i < numComp; ++i) {
975 int_local[i]+=f[i]*h0;
976 }
977 }
978 }
979
980 #pragma omp critical
981 for (index_t i=0; i<numComp; i++)
982 integrals[i]+=int_local[i];
983 } // end of parallel section
984 } // function space selector
985 }
986
987 //protected
988 dim_t Rectangle::insertNeighbourNodes(IndexVector& index, index_t node) const
989 {
990 const dim_t nDOF0 = (m_gNE0+1)/m_NX;
991 const dim_t nDOF1 = (m_gNE1+1)/m_NY;
992 const int x=node%nDOF0;
993 const int y=node/nDOF0;
994 dim_t num=0;
995 // loop through potential neighbours and add to index if positions are
996 // within bounds
997 for (int i1=-1; i1<2; i1++) {
998 for (int i0=-1; i0<2; i0++) {
999 // skip node itself
1000 if (i0==0 && i1==0)
1001 continue;
1002 // location of neighbour node
1003 const int nx=x+i0;
1004 const int ny=y+i1;
1005 if (nx>=0 && ny>=0 && nx<nDOF0 && ny<nDOF1) {
1006 index.push_back(ny*nDOF0+nx);
1007 num++;
1008 }
1009 }
1010 }
1011
1012 return num;
1013 }
1014
1015 //protected
1016 void Rectangle::nodesToDOF(escript::Data& out, escript::Data& in) const
1017 {
1018 const dim_t numComp = in.getDataPointSize();
1019 out.requireWrite();
1020
1021 const index_t left = (m_offset0==0 ? 0 : 1);
1022 const index_t bottom = (m_offset1==0 ? 0 : 1);
1023 const dim_t nDOF0 = (m_gNE0+1)/m_NX;
1024 const dim_t nDOF1 = (m_gNE1+1)/m_NY;
1025 #pragma omp parallel for
1026 for (index_t i=0; i<nDOF1; i++) {
1027 for (index_t j=0; j<nDOF0; j++) {
1028 const index_t n=j+left+(i+bottom)*m_N0;
1029 const double* src=in.getSampleDataRO(n);
1030 copy(src, src+numComp, out.getSampleDataRW(j+i*nDOF0));
1031 }
1032 }
1033 }
1034
1035 //protected
1036 void Rectangle::dofToNodes(escript::Data& out, escript::Data& in) const
1037 {
1038 const dim_t numComp = in.getDataPointSize();
1039 Paso_Coupler* coupler = Paso_Coupler_alloc(m_connector, numComp);
1040 in.requireWrite();
1041 Paso_Coupler_startCollect(coupler, in.getSampleDataRW(0));
1042
1043 const dim_t numDOF = getNumDOF();
1044 out.requireWrite();
1045 const double* buffer = Paso_Coupler_finishCollect(coupler);
1046
1047 #pragma omp parallel for
1048 for (index_t i=0; i<getNumNodes(); i++) {
1049 const double* src=(m_dofMap[i]<numDOF ?
1050 in.getSampleDataRO(m_dofMap[i])
1051 : &buffer[(m_dofMap[i]-numDOF)*numComp]);
1052 copy(src, src+numComp, out.getSampleDataRW(i));
1053 }
1054 }
1055
1056 //private
1057 void Rectangle::populateSampleIds()
1058 {
1059 // identifiers are ordered from left to right, bottom to top globablly.
1060
1061 // build node distribution vector first.
1062 // rank i owns m_nodeDistribution[i+1]-nodeDistribution[i] nodes
1063 m_nodeDistribution.assign(m_mpiInfo->size+1, 0);
1064 const dim_t numDOF=getNumDOF();
1065 for (dim_t k=1; k<m_mpiInfo->size; k++) {
1066 m_nodeDistribution[k]=k*numDOF;
1067 }
1068 m_nodeDistribution[m_mpiInfo->size]=getNumDataPointsGlobal();
1069 m_nodeId.resize(getNumNodes());
1070 m_dofId.resize(numDOF);
1071 m_elementId.resize(getNumElements());
1072 m_faceId.resize(getNumFaceElements());
1073
1074 #pragma omp parallel
1075 {
1076 // nodes
1077 #pragma omp for nowait
1078 for (dim_t i1=0; i1<m_N1; i1++) {
1079 for (dim_t i0=0; i0<m_N0; i0++) {
1080 m_nodeId[i0+i1*m_N0] = (m_offset1+i1)*(m_gNE0+1)+m_offset0+i0;
1081 }
1082 }
1083
1084 // degrees of freedom
1085 #pragma omp for nowait
1086 for (dim_t k=0; k<numDOF; k++)
1087 m_dofId[k] = m_nodeDistribution[m_mpiInfo->rank]+k;
1088
1089 // elements
1090 #pragma omp for nowait
1091 for (dim_t i1=0; i1<m_NE1; i1++) {
1092 for (dim_t i0=0; i0<m_NE0; i0++) {
1093 m_elementId[i0+i1*m_NE0]=(m_offset1+i1)*m_gNE0+m_offset0+i0;
1094 }
1095 }
1096
1097 // face elements
1098 #pragma omp for
1099 for (dim_t k=0; k<getNumFaceElements(); k++)
1100 m_faceId[k]=k;
1101 } // end parallel section
1102
1103 m_nodeTags.assign(getNumNodes(), 0);
1104 updateTagsInUse(Nodes);
1105
1106 m_elementTags.assign(getNumElements(), 0);
1107 updateTagsInUse(Elements);
1108
1109 // generate face offset vector and set face tags
1110 const IndexVector facesPerEdge = getNumFacesPerBoundary();
1111 const index_t LEFT=1, RIGHT=2, BOTTOM=10, TOP=20;
1112 const index_t faceTag[] = { LEFT, RIGHT, BOTTOM, TOP };
1113 m_faceOffset.assign(facesPerEdge.size(), -1);
1114 m_faceTags.clear();
1115 index_t offset=0;
1116 for (size_t i=0; i<facesPerEdge.size(); i++) {
1117 if (facesPerEdge[i]>0) {
1118 m_faceOffset[i]=offset;
1119 offset+=facesPerEdge[i];
1120 m_faceTags.insert(m_faceTags.end(), facesPerEdge[i], faceTag[i]);
1121 }
1122 }
1123 setTagMap("left", LEFT);
1124 setTagMap("right", RIGHT);
1125 setTagMap("bottom", BOTTOM);
1126 setTagMap("top", TOP);
1127 updateTagsInUse(FaceElements);
1128 }
1129
1130 //private
1131 void Rectangle::createPattern()
1132 {
1133 const dim_t nDOF0 = (m_gNE0+1)/m_NX;
1134 const dim_t nDOF1 = (m_gNE1+1)/m_NY;
1135 const index_t left = (m_offset0==0 ? 0 : 1);
1136 const index_t bottom = (m_offset1==0 ? 0 : 1);
1137
1138 // populate node->DOF mapping with own degrees of freedom.
1139 // The rest is assigned in the loop further down
1140 m_dofMap.assign(getNumNodes(), 0);
1141 #pragma omp parallel for
1142 for (index_t i=bottom; i<bottom+nDOF1; i++) {
1143 for (index_t j=left; j<left+nDOF0; j++) {
1144 m_dofMap[i*m_N0+j]=(i-bottom)*nDOF0+j-left;
1145 }
1146 }
1147
1148 // build list of shared components and neighbours by looping through
1149 // all potential neighbouring ranks and checking if positions are
1150 // within bounds
1151 const dim_t numDOF=nDOF0*nDOF1;
1152 vector<IndexVector> colIndices(numDOF); // for the couple blocks
1153 RankVector neighbour;
1154 IndexVector offsetInShared(1,0);
1155 IndexVector sendShared, recvShared;
1156 int numShared=0;
1157 const int x=m_mpiInfo->rank%m_NX;
1158 const int y=m_mpiInfo->rank/m_NX;
1159 for (int i1=-1; i1<2; i1++) {
1160 for (int i0=-1; i0<2; i0++) {
1161 // skip this rank
1162 if (i0==0 && i1==0)
1163 continue;
1164 // location of neighbour rank
1165 const int nx=x+i0;
1166 const int ny=y+i1;
1167 if (nx>=0 && ny>=0 && nx<m_NX && ny<m_NY) {
1168 neighbour.push_back(ny*m_NX+nx);
1169 if (i0==0) {
1170 // sharing top or bottom edge
1171 const int firstDOF=(i1==-1 ? 0 : numDOF-nDOF0);
1172 const int firstNode=(i1==-1 ? left : m_N0*(m_N1-1)+left);
1173 offsetInShared.push_back(offsetInShared.back()+nDOF0);
1174 for (dim_t i=0; i<nDOF0; i++, numShared++) {
1175 sendShared.push_back(firstDOF+i);
1176 recvShared.push_back(numDOF+numShared);
1177 if (i>0)
1178 colIndices[firstDOF+i-1].push_back(numShared);
1179 colIndices[firstDOF+i].push_back(numShared);
1180 if (i<nDOF0-1)
1181 colIndices[firstDOF+i+1].push_back(numShared);
1182 m_dofMap[firstNode+i]=numDOF+numShared;
1183 }
1184 } else if (i1==0) {
1185 // sharing left or right edge
1186 const int firstDOF=(i0==-1 ? 0 : nDOF0-1);
1187 const int firstNode=(i0==-1 ? bottom*m_N0 : (bottom+1)*m_N0-1);
1188 offsetInShared.push_back(offsetInShared.back()+nDOF1);
1189 for (dim_t i=0; i<nDOF1; i++, numShared++) {
1190 sendShared.push_back(firstDOF+i*nDOF0);
1191 recvShared.push_back(numDOF+numShared);
1192 if (i>0)
1193 colIndices[firstDOF+(i-1)*nDOF0].push_back(numShared);
1194 colIndices[firstDOF+i*nDOF0].push_back(numShared);
1195 if (i<nDOF1-1)
1196 colIndices[firstDOF+(i+1)*nDOF0].push_back(numShared);
1197 m_dofMap[firstNode+i*m_N0]=numDOF+numShared;
1198 }
1199 } else {
1200 // sharing a node
1201 const int dof=(i0+1)/2*(nDOF0-1)+(i1+1)/2*(numDOF-nDOF0);
1202 const int node=(i0+1)/2*(m_N0-1)+(i1+1)/2*m_N0*(m_N1-1);
1203 offsetInShared.push_back(offsetInShared.back()+1);
1204 sendShared.push_back(dof);
1205 recvShared.push_back(numDOF+numShared);
1206 colIndices[dof].push_back(numShared);
1207 m_dofMap[node]=numDOF+numShared;
1208 ++numShared;
1209 }
1210 }
1211 }
1212 }
1213
1214 // create connector
1215 Paso_SharedComponents *snd_shcomp = Paso_SharedComponents_alloc(
1216 numDOF, neighbour.size(), &neighbour[0], &sendShared[0],
1217 &offsetInShared[0], 1, 0, m_mpiInfo);
1218 Paso_SharedComponents *rcv_shcomp = Paso_SharedComponents_alloc(
1219 numDOF, neighbour.size(), &neighbour[0], &recvShared[0],
1220 &offsetInShared[0], 1, 0, m_mpiInfo);
1221 m_connector = Paso_Connector_alloc(snd_shcomp, rcv_shcomp);
1222 Paso_SharedComponents_free(snd_shcomp);
1223 Paso_SharedComponents_free(rcv_shcomp);
1224
1225 // create main and couple blocks
1226 Paso_Pattern *mainPattern = createMainPattern();
1227 Paso_Pattern *colPattern, *rowPattern;
1228 createCouplePatterns(colIndices, numShared, &colPattern, &rowPattern);
1229
1230 // allocate paso distribution
1231 Paso_Distribution* distribution = Paso_Distribution_alloc(m_mpiInfo,
1232 const_cast<index_t*>(&m_nodeDistribution[0]), 1, 0);
1233
1234 // finally create the system matrix
1235 m_pattern = Paso_SystemMatrixPattern_alloc(MATRIX_FORMAT_DEFAULT,
1236 distribution, distribution, mainPattern, colPattern, rowPattern,
1237 m_connector, m_connector);
1238
1239 Paso_Distribution_free(distribution);
1240
1241 // useful debug output
1242 /*
1243 cout << "--- rcv_shcomp ---" << endl;
1244 cout << "numDOF=" << numDOF << ", numNeighbors=" << neighbour.size() << endl;
1245 for (size_t i=0; i<neighbour.size(); i++) {
1246 cout << "neighbor[" << i << "]=" << neighbour[i]
1247 << " offsetInShared[" << i+1 << "]=" << offsetInShared[i+1] << endl;
1248 }
1249 for (size_t i=0; i<recvShared.size(); i++) {
1250 cout << "shared[" << i << "]=" << recvShared[i] << endl;
1251 }
1252 cout << "--- snd_shcomp ---" << endl;
1253 for (size_t i=0; i<sendShared.size(); i++) {
1254 cout << "shared[" << i << "]=" << sendShared[i] << endl;
1255 }
1256 cout << "--- dofMap ---" << endl;
1257 for (size_t i=0; i<m_dofMap.size(); i++) {
1258 cout << "m_dofMap[" << i << "]=" << m_dofMap[i] << endl;
1259 }
1260 cout << "--- colIndices ---" << endl;
1261 for (size_t i=0; i<colIndices.size(); i++) {
1262 cout << "colIndices[" << i << "].size()=" << colIndices[i].size() << endl;
1263 }
1264 */
1265
1266 /*
1267 cout << "--- main_pattern ---" << endl;
1268 cout << "M=" << mainPattern->numOutput << ", N=" << mainPattern->numInput << endl;
1269 for (size_t i=0; i<mainPattern->numOutput+1; i++) {
1270 cout << "ptr[" << i << "]=" << mainPattern->ptr[i] << endl;
1271 }
1272 for (size_t i=0; i<mainPattern->ptr[mainPattern->numOutput]; i++) {
1273 cout << "index[" << i << "]=" << mainPattern->index[i] << endl;
1274 }
1275 */
1276
1277 /*
1278 cout << "--- colCouple_pattern ---" << endl;
1279 cout << "M=" << colPattern->numOutput << ", N=" << colPattern->numInput << endl;
1280 for (size_t i=0; i<colPattern->numOutput+1; i++) {
1281 cout << "ptr[" << i << "]=" << colPattern->ptr[i] << endl;
1282 }
1283 for (size_t i=0; i<colPattern->ptr[colPattern->numOutput]; i++) {
1284 cout << "index[" << i << "]=" << colPattern->index[i] << endl;
1285 }
1286 */
1287
1288 /*
1289 cout << "--- rowCouple_pattern ---" << endl;
1290 cout << "M=" << rowPattern->numOutput << ", N=" << rowPattern->numInput << endl;
1291 for (size_t i=0; i<rowPattern->numOutput+1; i++) {
1292 cout << "ptr[" << i << "]=" << rowPattern->ptr[i] << endl;
1293 }
1294 for (size_t i=0; i<rowPattern->ptr[rowPattern->numOutput]; i++) {
1295 cout << "index[" << i << "]=" << rowPattern->index[i] << endl;
1296 }
1297 */
1298
1299 Paso_Pattern_free(mainPattern);
1300 Paso_Pattern_free(colPattern);
1301 Paso_Pattern_free(rowPattern);
1302 }
1303
1304 //private
1305 void Rectangle::addToMatrixAndRHS(Paso_SystemMatrix* S, escript::Data& F,
1306 const vector<double>& EM_S, const vector<double>& EM_F, bool addS,
1307 bool addF, index_t firstNode, dim_t nEq, dim_t nComp) const
1308 {
1309 IndexVector rowIndex;
1310 rowIndex.push_back(m_dofMap[firstNode]);
1311 rowIndex.push_back(m_dofMap[firstNode+1]);
1312 rowIndex.push_back(m_dofMap[firstNode+m_N0]);
1313 rowIndex.push_back(m_dofMap[firstNode+m_N0+1]);
1314 if (addF) {
1315 double *F_p=F.getSampleDataRW(0);
1316 for (index_t i=0; i<rowIndex.size(); i++) {
1317 if (rowIndex[i]<getNumDOF()) {
1318 for (index_t eq=0; eq<nEq; eq++) {
1319 F_p[INDEX2(eq, rowIndex[i], nEq)]+=EM_F[INDEX2(eq,i,nEq)];
1320 }
1321 }
1322 }
1323 }
1324 if (addS) {
1325 addToSystemMatrix(S, rowIndex, nEq, rowIndex, nComp, EM_S);
1326 }
1327 }
1328
1329 //protected
1330 void Rectangle::interpolateNodesOnElements(escript::Data& out,
1331 escript::Data& in, bool reduced) const
1332 {
1333 const dim_t numComp = in.getDataPointSize();
1334 if (reduced) {
1335 out.requireWrite();
1336 const double c0 = 0.25;
1337 #pragma omp parallel
1338 {
1339 vector<double> f_00(numComp);
1340 vector<double> f_01(numComp);
1341 vector<double> f_10(numComp);
1342 vector<double> f_11(numComp);
1343 #pragma omp for
1344 for (index_t k1=0; k1 < m_NE1; ++k1) {
1345 for (index_t k0=0; k0 < m_NE0; ++k0) {
1346 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,k1, m_N0)), numComp*sizeof(double));
1347 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,k1+1, m_N0)), numComp*sizeof(double));
1348 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,k1, m_N0)), numComp*sizeof(double));
1349 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0)), numComp*sizeof(double));
1350 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
1351 for (index_t i=0; i < numComp; ++i) {
1352 o[INDEX2(i,numComp,0)] = c0*(f_00[i] + f_01[i] + f_10[i] + f_11[i]);
1353 } /* end of component loop i */
1354 } /* end of k0 loop */
1355 } /* end of k1 loop */
1356 } /* end of parallel section */
1357 } else {
1358 out.requireWrite();
1359 const double c0 = 0.16666666666666666667;
1360 const double c1 = 0.044658198738520451079;
1361 const double c2 = 0.62200846792814621559;
1362 #pragma omp parallel
1363 {
1364 vector<double> f_00(numComp);
1365 vector<double> f_01(numComp);
1366 vector<double> f_10(numComp);
1367 vector<double> f_11(numComp);
1368 #pragma omp for
1369 for (index_t k1=0; k1 < m_NE1; ++k1) {
1370 for (index_t k0=0; k0 < m_NE0; ++k0) {
1371 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,k1, m_N0)), numComp*sizeof(double));
1372 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,k1+1, m_N0)), numComp*sizeof(double));
1373 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,k1, m_N0)), numComp*sizeof(double));
1374 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,k1+1, m_N0)), numComp*sizeof(double));
1375 double* o = out.getSampleDataRW(INDEX2(k0,k1,m_NE0));
1376 for (index_t i=0; i < numComp; ++i) {
1377 o[INDEX2(i,numComp,0)] = c0*(f_01[i] + f_10[i]) + c1*f_11[i] + c2*f_00[i];
1378 o[INDEX2(i,numComp,1)] = c0*(f_00[i] + f_11[i]) + c1*f_01[i] + c2*f_10[i];
1379 o[INDEX2(i,numComp,2)] = c0*(f_00[i] + f_11[i]) + c1*f_10[i] + c2*f_01[i];
1380 o[INDEX2(i,numComp,3)] = c0*(f_01[i] + f_10[i]) + c1*f_00[i] + c2*f_11[i];
1381 } /* end of component loop i */
1382 } /* end of k0 loop */
1383 } /* end of k1 loop */
1384 } /* end of parallel section */
1385 }
1386 }
1387
1388 //protected
1389 void Rectangle::interpolateNodesOnFaces(escript::Data& out, escript::Data& in,
1390 bool reduced) const
1391 {
1392 const dim_t numComp = in.getDataPointSize();
1393 if (reduced) {
1394 out.requireWrite();
1395 const double c0 = 0.5;
1396 #pragma omp parallel
1397 {
1398 vector<double> f_00(numComp);
1399 vector<double> f_01(numComp);
1400 vector<double> f_10(numComp);
1401 vector<double> f_11(numComp);
1402 if (m_faceOffset[0] > -1) {
1403 #pragma omp for nowait
1404 for (index_t k1=0; k1 < m_NE1; ++k1) {
1405 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(0,k1, m_N0)), numComp*sizeof(double));
1406 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(0,k1+1, m_N0)), numComp*sizeof(double));
1407 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
1408 for (index_t i=0; i < numComp; ++i) {
1409 o[INDEX2(i,numComp,0)] = c0*(f_00[i] + f_01[i]);
1410 } /* end of component loop i */
1411 } /* end of k1 loop */
1412 } /* end of face 0 */
1413 if (m_faceOffset[1] > -1) {
1414 #pragma omp for nowait
1415 for (index_t k1=0; k1 < m_NE1; ++k1) {
1416 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0)), numComp*sizeof(double));
1417 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0)), numComp*sizeof(double));
1418 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
1419 for (index_t i=0; i < numComp; ++i) {
1420 o[INDEX2(i,numComp,0)] = c0*(f_10[i] + f_11[i]);
1421 } /* end of component loop i */
1422 } /* end of k1 loop */
1423 } /* end of face 1 */
1424 if (m_faceOffset[2] > -1) {
1425 #pragma omp for nowait
1426 for (index_t k0=0; k0 < m_NE0; ++k0) {
1427 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,0, m_N0)), numComp*sizeof(double));
1428 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,0, m_N0)), numComp*sizeof(double));
1429 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
1430 for (index_t i=0; i < numComp; ++i) {
1431 o[INDEX2(i,numComp,0)] = c0*(f_00[i] + f_10[i]);
1432 } /* end of component loop i */
1433 } /* end of k0 loop */
1434 } /* end of face 2 */
1435 if (m_faceOffset[3] > -1) {
1436 #pragma omp for nowait
1437 for (index_t k0=0; k0 < m_NE0; ++k0) {
1438 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0)), numComp*sizeof(double));
1439 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0)), numComp*sizeof(double));
1440 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
1441 for (index_t i=0; i < numComp; ++i) {
1442 o[INDEX2(i,numComp,0)] = c0*(f_01[i] + f_11[i]);
1443 } /* end of component loop i */
1444 } /* end of k0 loop */
1445 } /* end of face 3 */
1446 } /* end of parallel section */
1447 } else {
1448 out.requireWrite();
1449 const double c0 = 0.21132486540518711775;
1450 const double c1 = 0.78867513459481288225;
1451 #pragma omp parallel
1452 {
1453 vector<double> f_00(numComp);
1454 vector<double> f_01(numComp);
1455 vector<double> f_10(numComp);
1456 vector<double> f_11(numComp);
1457 if (m_faceOffset[0] > -1) {
1458 #pragma omp for nowait
1459 for (index_t k1=0; k1 < m_NE1; ++k1) {
1460 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(0,k1, m_N0)), numComp*sizeof(double));
1461 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(0,k1+1, m_N0)), numComp*sizeof(double));
1462 double* o = out.getSampleDataRW(m_faceOffset[0]+k1);
1463 for (index_t i=0; i < numComp; ++i) {
1464 o[INDEX2(i,numComp,0)] = c0*f_01[i] + c1*f_00[i];
1465 o[INDEX2(i,numComp,1)] = c0*f_00[i] + c1*f_01[i];
1466 } /* end of component loop i */
1467 } /* end of k1 loop */
1468 } /* end of face 0 */
1469 if (m_faceOffset[1] > -1) {
1470 #pragma omp for nowait
1471 for (index_t k1=0; k1 < m_NE1; ++k1) {
1472 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(m_N0-1,k1, m_N0)), numComp*sizeof(double));
1473 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(m_N0-1,k1+1, m_N0)), numComp*sizeof(double));
1474 double* o = out.getSampleDataRW(m_faceOffset[1]+k1);
1475 for (index_t i=0; i < numComp; ++i) {
1476 o[INDEX2(i,numComp,0)] = c1*f_10[i] + c0*f_11[i];
1477 o[INDEX2(i,numComp,1)] = c1*f_11[i] + c0*f_10[i];
1478 } /* end of component loop i */
1479 } /* end of k1 loop */
1480 } /* end of face 1 */
1481 if (m_faceOffset[2] > -1) {
1482 #pragma omp for nowait
1483 for (index_t k0=0; k0 < m_NE0; ++k0) {
1484 memcpy(&f_00[0], in.getSampleDataRO(INDEX2(k0,0, m_N0)), numComp*sizeof(double));
1485 memcpy(&f_10[0], in.getSampleDataRO(INDEX2(k0+1,0, m_N0)), numComp*sizeof(double));
1486 double* o = out.getSampleDataRW(m_faceOffset[2]+k0);
1487 for (index_t i=0; i < numComp; ++i) {
1488 o[INDEX2(i,numComp,0)] = c0*f_10[i] + c1*f_00[i];
1489 o[INDEX2(i,numComp,1)] = c0*f_00[i] + c1*f_10[i];
1490 } /* end of component loop i */
1491 } /* end of k0 loop */
1492 } /* end of face 2 */
1493 if (m_faceOffset[3] > -1) {
1494 #pragma omp for nowait
1495 for (index_t k0=0; k0 < m_NE0; ++k0) {
1496 memcpy(&f_01[0], in.getSampleDataRO(INDEX2(k0,m_N1-1, m_N0)), numComp*sizeof(double));
1497 memcpy(&f_11[0], in.getSampleDataRO(INDEX2(k0+1,m_N1-1, m_N0)), numComp*sizeof(double));
1498 double* o = out.getSampleDataRW(m_faceOffset[3]+k0);
1499 for (index_t i=0; i < numComp; ++i) {
1500 o[INDEX2(i,numComp,0)] = c0*f_11[i] + c1*f_01[i];
1501 o[INDEX2(i,numComp,1)] = c0*f_01[i] + c1*f_11[i];
1502 } /* end of component loop i */
1503 } /* end of k0 loop */
1504 } /* end of face 3 */
1505 } /* end of parallel section */
1506 }
1507 }
1508
1509 //protected
1510 void Rectangle::assemblePDESingle(Paso_SystemMatrix* mat,
1511 escript::Data& rhs, const escript::Data& A, const escript::Data& B,
1512 const escript::Data& C, const escript::Data& D,
1513 const escript::Data& X, const escript::Data& Y) const
1514 {
1515 const double h0 = m_l0/m_gNE0;
1516 const double h1 = m_l1/m_gNE1;
1517 const double w0 = -0.1555021169820365539*h1/h0;
1518 const double w1 = 0.041666666666666666667;
1519 const double w2 = -0.15550211698203655390;
1520 const double w3 = 0.041666666666666666667*h0/h1;
1521 const double w4 = 0.15550211698203655390;
1522 const double w5 = -0.041666666666666666667;
1523 const double w6 = -0.01116454968463011277*h1/h0;
1524 const double w7 = 0.011164549684630112770;
1525 const double w8 = -0.011164549684630112770;
1526 const double w9 = -0.041666666666666666667*h1/h0;
1527 const double w10 = -0.041666666666666666667*h0/h1;
1528 const double w11 = 0.1555021169820365539*h1/h0;
1529 const double w12 = 0.1555021169820365539*h0/h1;
1530 const double w13 = 0.01116454968463011277*h0/h1;
1531 const double w14 = 0.01116454968463011277*h1/h0;
1532 const double w15 = 0.041666666666666666667*h1/h0;
1533 const double w16 = -0.01116454968463011277*h0/h1;
1534 const double w17 = -0.1555021169820365539*h0/h1;
1535 const double w18 = -0.33333333333333333333*h1/h0;
1536 const double w19 = 0.25;
1537 const double w20 = -0.25;
1538 const double w21 = 0.16666666666666666667*h0/h1;
1539 const double w22 = -0.16666666666666666667*h1/h0;
1540 const double w23 = -0.16666666666666666667*h0/h1;
1541 const double w24 = 0.33333333333333333333*h1/h0;
1542 const double w25 = 0.33333333333333333333*h0/h1;
1543 const double w26 = 0.16666666666666666667*h1/h0;
1544 const double w27 = -0.33333333333333333333*h0/h1;
1545 const double w28 = -0.032861463941450536761*h1;
1546 const double w29 = -0.032861463941450536761*h0;
1547 const double w30 = -0.12264065304058601714*h1;
1548 const double w31 = -0.0023593469594139828636*h1;
1549 const double w32 = -0.008805202725216129906*h0;
1550 const double w33 = -0.008805202725216129906*h1;
1551 const double w34 = 0.032861463941450536761*h1;
1552 const double w35 = 0.008805202725216129906*h1;
1553 const double w36 = 0.008805202725216129906*h0;
1554 const double w37 = 0.0023593469594139828636*h1;
1555 const double w38 = 0.12264065304058601714*h1;
1556 const double w39 = 0.032861463941450536761*h0;
1557 const double w40 = -0.12264065304058601714*h0;
1558 const double w41 = -0.0023593469594139828636*h0;
1559 const double w42 = 0.0023593469594139828636*h0;
1560 const double w43 = 0.12264065304058601714*h0;
1561 const double w44 = -0.16666666666666666667*h1;
1562 const double w45 = -0.083333333333333333333*h0;
1563 const double w46 = 0.083333333333333333333*h1;
1564 const double w47 = 0.16666666666666666667*h1;
1565 const double w48 = 0.083333333333333333333*h0;
1566 const double w49 = -0.16666666666666666667*h0;
1567 const double w50 = 0.16666666666666666667*h0;
1568 const double w51 = -0.083333333333333333333*h1;
1569 const double w52 = 0.025917019497006092316*h0*h1;
1570 const double w53 = 0.0018607582807716854616*h0*h1;
1571 const double w54 = 0.0069444444444444444444*h0*h1;
1572 const double w55 = 0.09672363354357992482*h0*h1;
1573 const double w56 = 0.00049858867864229740201*h0*h1;
1574 const double w57 = 0.055555555555555555556*h0*h1;
1575 const double w58 = 0.027777777777777777778*h0*h1;
1576 const double w59 = 0.11111111111111111111*h0*h1;
1577 const double w60 = -0.19716878364870322056*h1;
1578 const double w61 = -0.19716878364870322056*h0;
1579 const double w62 = -0.052831216351296779436*h0;
1580 const double w63 = -0.052831216351296779436*h1;
1581 const double w64 = 0.19716878364870322056*h1;
1582 const double w65 = 0.052831216351296779436*h1;
1583 const double w66 = 0.19716878364870322056*h0;
1584 const double w67 = 0.052831216351296779436*h0;
1585 const double w68 = -0.5*h1;
1586 const double w69 = -0.5*h0;
1587 const double w70 = 0.5*h1;
1588 const double w71 = 0.5*h0;
1589 const double w72 = 0.1555021169820365539*h0*h1;
1590 const double w73 = 0.041666666666666666667*h0*h1;
1591 const double w74 = 0.01116454968463011277*h0*h1;
1592 const double w75 = 0.25*h0*h1;
1593
1594 rhs.requireWrite();
1595 #pragma omp parallel
1596 {
1597 for (index_t k1_0=0; k1_0<2; k1_0++) { // colouring
1598 #pragma omp for
1599 for (index_t k1=k1_0; k1<m_NE1; k1+=2) {
1600 for (index_t k0=0; k0<m_NE0; ++k0) {
1601 bool add_EM_S=false;
1602 bool add_EM_F=false;
1603 vector<double> EM_S(4*4, 0);
1604 vector<double> EM_F(4, 0);
1605 const index_t e = k0 + m_NE0*k1;
1606 ///////////////
1607 // process A //
1608 ///////////////
1609 if (!A.isEmpty()) {
1610 add_EM_S=true;
1611 const double* A_p=const_cast<escript::Data*>(&A)->getSampleDataRO(e);
1612 if (A.actsExpanded()) {
1613 const double A_00_0 = A_p[INDEX3(0,0,0,2,2)];
1614 const double A_10_0 = A_p[INDEX3(1,0,0,2,2)];
1615 const double A_01_0 = A_p[INDEX3(0,1,0,2,2)];
1616 const double A_11_0 = A_p[INDEX3(1,1,0,2,2)];
1617 const double A_00_1 = A_p[INDEX3(0,0,1,2,2)];
1618 const double A_10_1 = A_p[INDEX3(1,0,1,2,2)];
1619 const double A_01_1 = A_p[INDEX3(0,1,1,2,2)];
1620 const double A_11_1 = A_p[INDEX3(1,1,1,2,2)];
1621 const double A_00_2 = A_p[INDEX3(0,0,2,2,2)];
1622 const double A_10_2 = A_p[INDEX3(1,0,2,2,2)];
1623 const double A_01_2 = A_p[INDEX3(0,1,2,2,2)];
1624 const double A_11_2 = A_p[INDEX3(1,1,2,2,2)];
1625 const double A_00_3 = A_p[INDEX3(0,0,3,2,2)];
1626 const double A_10_3 = A_p[INDEX3(1,0,3,2,2)];
1627 const double A_01_3 = A_p[INDEX3(0,1,3,2,2)];
1628 const double A_11_3 = A_p[INDEX3(1,1,3,2,2)];
1629 const double tmp0_0 = A_01_0 + A_01_3;
1630 const double tmp1_0 = A_00_0 + A_00_1;
1631 const double tmp2_0 = A_11_0 + A_11_1 + A_11_2 + A_11_3;
1632 const double tmp3_0 = A_00_2 + A_00_3;
1633 const double tmp4_0 = A_10_1 + A_10_2;
1634 const double tmp5_0 = A_00_0 + A_00_1 + A_00_2 + A_00_3;
1635 const double tmp6_0 = A_01_3 + A_10_0;
1636 const double tmp7_0 = A_01_0 + A_10_3;
1637 const double tmp8_0 = A_01_1 + A_01_2 + A_10_1 + A_10_2;
1638 const double tmp9_0 = A_01_0 + A_10_0;
1639 const double tmp12_0 = A_11_0 + A_11_2;
1640 const double tmp10_0 = A_01_3 + A_10_3;
1641 const double tmp14_0 = A_01_0 + A_01_3 + A_10_0 + A_10_3;
1642 const double tmp13_0 = A_01_2 + A_10_1;
1643 const double tmp11_0 = A_11_1 + A_11_3;
1644 const double tmp18_0 = A_01_1 + A_10_1;
1645 const double tmp15_0 = A_01_1 + A_10_2;
1646 const double tmp16_0 = A_10_0 + A_10_3;
1647 const double tmp17_0 = A_01_1 + A_01_2;
1648 const double tmp19_0 = A_01_2 + A_10_2;
1649 const double tmp0_1 = A_10_3*w8;
1650 const double tmp1_1 = tmp0_0*w1;
1651 const double tmp2_1 = A_01_1*w4;
1652 const double tmp3_1 = tmp1_0*w0;
1653 const double tmp4_1 = A_01_2*w7;
1654 const double tmp5_1 = tmp2_0*w3;
1655 const double tmp6_1 = tmp3_0*w6;
1656 const double tmp7_1 = A_10_0*w2;
1657 const double tmp8_1 = tmp4_0*w5;
1658 const double tmp9_1 = tmp2_0*w10;
1659 const double tmp14_1 = A_10_0*w8;
1660 const double tmp23_1 = tmp3_0*w14;
1661 const double tmp35_1 = A_01_0*w8;
1662 const double tmp54_1 = tmp13_0*w8;
1663 const double tmp20_1 = tmp9_0*w4;
1664 const double tmp25_1 = tmp12_0*w12;
1665 const double tmp44_1 = tmp7_0*w7;
1666 const double tmp26_1 = tmp10_0*w4;
1667 const double tmp52_1 = tmp18_0*w8;
1668 const double tmp48_1 = A_10_1*w7;
1669 const double tmp46_1 = A_01_3*w8;
1670 const double tmp50_1 = A_01_0*w2;
1671 const double tmp56_1 = tmp19_0*w8;
1672 const double tmp19_1 = A_10_3*w2;
1673 const double tmp47_1 = A_10_2*w4;
1674 const double tmp16_1 = tmp3_0*w0;
1675 const double tmp18_1 = tmp1_0*w6;
1676 const double tmp31_1 = tmp11_0*w12;
1677 const double tmp55_1 = tmp15_0*w2;
1678 const double tmp39_1 = A_10_2*w7;
1679 const double tmp11_1 = tmp6_0*w7;
1680 const double tmp40_1 = tmp11_0*w17;
1681 const double tmp34_1 = tmp15_0*w8;
1682 const double tmp33_1 = tmp14_0*w5;
1683 const double tmp24_1 = tmp11_0*w13;
1684 const double tmp43_1 = tmp17_0*w5;
1685 const double tmp15_1 = A_01_2*w4;
1686 const double tmp53_1 = tmp19_0*w2;
1687 const double tmp27_1 = tmp3_0*w11;
1688 const double tmp32_1 = tmp13_0*w2;
1689 const double tmp10_1 = tmp5_0*w9;
1690 const double tmp37_1 = A_10_1*w4;
1691 const double tmp38_1 = tmp5_0*w15;
1692 const double tmp17_1 = A_01_1*w7;
1693 const double tmp12_1 = tmp7_0*w4;
1694 const double tmp22_1 = tmp10_0*w7;
1695 const double tmp57_1 = tmp18_0*w2;
1696 const double tmp28_1 = tmp9_0*w7;
1697 const double tmp29_1 = tmp1_0*w14;
1698 const double tmp51_1 = tmp11_0*w16;
1699 const double tmp42_1 = tmp12_0*w16;
1700 const double tmp49_1 = tmp12_0*w17;
1701 const double tmp21_1 = tmp1_0*w11;
1702 const double tmp45_1 = tmp6_0*w4;
1703 const double tmp13_1 = tmp8_0*w1;
1704 const double tmp36_1 = tmp16_0*w1;
1705 const double tmp41_1 = A_01_3*w2;
1706 const double tmp30_1 = tmp12_0*w13;
1707 EM_S[INDEX2(0,0,4)]+=tmp13_1 + tmp20_1 + tmp21_1 + tmp22_1 + tmp23_1 + tmp24_1 + tmp25_1;
1708 EM_S[INDEX2(1,0,4)]+=tmp36_1 + tmp37_1 + tmp39_1 + tmp3_1 + tmp43_1 + tmp46_1 + tmp50_1 + tmp5_1 + tmp6_1;
1709 EM_S[INDEX2(2,0,4)]+=tmp0_1 + tmp15_1 + tmp17_1 + tmp1_1 + tmp38_1 + tmp49_1 + tmp51_1 + tmp7_1 + tmp8_1;
1710 EM_S[INDEX2(3,0,4)]+=tmp10_1 + tmp32_1 + tmp33_1 + tmp34_1 + tmp9_1;
1711 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1 + tmp4_1 + tmp5_1 + tmp6_1 + tmp7_1 + tmp8_1;
1712 EM_S[INDEX2(1,1,4)]+=tmp21_1 + tmp23_1 + tmp30_1 + tmp31_1 + tmp33_1 + tmp56_1 + tmp57_1;
1713 EM_S[INDEX2(2,1,4)]+=tmp10_1 + tmp13_1 + tmp44_1 + tmp45_1 + tmp9_1;
1714 EM_S[INDEX2(3,1,4)]+=tmp35_1 + tmp36_1 + tmp37_1 + tmp38_1 + tmp39_1 + tmp40_1 + tmp41_1 + tmp42_1 + tmp43_1;
1715 EM_S[INDEX2(0,2,4)]+=tmp36_1 + tmp38_1 + tmp43_1 + tmp46_1 + tmp47_1 + tmp48_1 + tmp49_1 + tmp50_1 + tmp51_1;
1716 EM_S[INDEX2(1,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1 + tmp9_1;
1717 EM_S[INDEX2(2,2,4)]+=tmp24_1 + tmp25_1 + tmp27_1 + tmp29_1 + tmp33_1 + tmp52_1 + tmp53_1;
1718 EM_S[INDEX2(3,2,4)]+=tmp14_1 + tmp15_1 + tmp16_1 + tmp17_1 + tmp18_1 + tmp19_1 + tmp1_1 + tmp5_1 + tmp8_1;
1719 EM_S[INDEX2(0,3,4)]+=tmp10_1 + tmp33_1 + tmp54_1 + tmp55_1 + tmp9_1;
1720 EM_S[INDEX2(1,3,4)]+=tmp14_1 + tmp19_1 + tmp1_1 + tmp2_1 + tmp38_1 + tmp40_1 + tmp42_1 + tmp4_1 + tmp8_1;
1721 EM_S[INDEX2(2,3,4)]+=tmp16_1 + tmp18_1 + tmp35_1 + tmp36_1 + tmp41_1 + tmp43_1 + tmp47_1 + tmp48_1 + tmp5_1;
1722 EM_S[INDEX2(3,3,4)]+=tmp13_1 + tmp26_1 + tmp27_1 + tmp28_1 + tmp29_1 + tmp30_1 + tmp31_1;
1723 } else { // constant data
1724 const double A_00 = A_p[INDEX2(0,0,2)];
1725 const double A_10 = A_p[INDEX2(1,0,2)];
1726 const double A_01 = A_p[INDEX2(0,1,2)];
1727 const double A_11 = A_p[INDEX2(1,1,2)];
1728 const double tmp0_0 = A_01 + A_10;
1729 const double tmp0_1 = A_00*w18;
1730 const double tmp1_1 = A_01*w19;
1731 const double tmp2_1 = A_10*w20;
1732 const double tmp3_1 = A_11*w21;
1733 const double tmp4_1 = A_00*w22;
1734 const double tmp5_1 = tmp0_0*w19;
1735 const double tmp6_1 = A_11*w23;
1736 const double tmp7_1 = A_11*w25;
1737 const double tmp8_1 = A_00*w24;
1738 const double tmp9_1 = tmp0_0*w20;
1739 const double tmp10_1 = A_01*w20;
1740 const double tmp11_1 = A_11*w27;
1741 const double tmp12_1 = A_00*w26;
1742 const double tmp13_1 = A_10*w19;
1743 EM_S[INDEX2(0,0,4)]+=tmp5_1 + tmp7_1 + tmp8_1;
1744 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp10_1 + tmp13_1 + tmp3_1;
1745 EM_S[INDEX2(2,0,4)]+=tmp11_1 + tmp12_1 + tmp1_1 + tmp2_1;
1746 EM_S[INDEX2(3,0,4)]+=tmp4_1 + tmp6_1 + tmp9_1;
1747 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1;
1748 EM_S[INDEX2(1,1,4)]+=tmp7_1 + tmp8_1 + tmp9_1;
1749 EM_S[INDEX2(2,1,4)]+=tmp4_1 + tmp5_1 + tmp6_1;
1750 EM_S[INDEX2(3,1,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1;
1751 EM_S[INDEX2(0,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1;
1752 EM_S[INDEX2(1,2,4)]+=tmp4_1 + tmp5_1 + tmp6_1;
1753 EM_S[INDEX2(2,2,4)]+=tmp7_1 + tmp8_1 + tmp9_1;
1754 EM_S[INDEX2(3,2,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1;
1755 EM_S[INDEX2(0,3,4)]+=tmp4_1 + tmp6_1 + tmp9_1;
1756 EM_S[INDEX2(1,3,4)]+=tmp11_1 + tmp12_1 + tmp1_1 + tmp2_1;
1757 EM_S[INDEX2(2,3,4)]+=tmp0_1 + tmp10_1 + tmp13_1 + tmp3_1;
1758 EM_S[INDEX2(3,3,4)]+=tmp5_1 + tmp7_1 + tmp8_1;
1759 }
1760 }
1761 ///////////////
1762 // process B //
1763 ///////////////
1764 if (!B.isEmpty()) {
1765 add_EM_S=true;
1766 const double* B_p=const_cast<escript::Data*>(&B)->getSampleDataRO(e);
1767 if (B.actsExpanded()) {
1768 const double B_0_0 = B_p[INDEX2(0,0,2)];
1769 const double B_1_0 = B_p[INDEX2(1,0,2)];
1770 const double B_0_1 = B_p[INDEX2(0,1,2)];
1771 const double B_1_1 = B_p[INDEX2(1,1,2)];
1772 const double B_0_2 = B_p[INDEX2(0,2,2)];
1773 const double B_1_2 = B_p[INDEX2(1,2,2)];
1774 const double B_0_3 = B_p[INDEX2(0,3,2)];
1775 const double B_1_3 = B_p[INDEX2(1,3,2)];
1776 const double tmp0_0 = B_1_0 + B_1_1;
1777 const double tmp1_0 = B_1_2 + B_1_3;
1778 const double tmp2_0 = B_0_1 + B_0_3;
1779 const double tmp3_0 = B_0_0 + B_0_2;
1780 const double tmp63_1 = B_1_1*w42;
1781 const double tmp79_1 = B_1_1*w40;
1782 const double tmp37_1 = tmp3_0*w35;
1783 const double tmp8_1 = tmp0_0*w32;
1784 const double tmp71_1 = B_0_1*w34;
1785 const double tmp19_1 = B_0_3*w31;
1786 const double tmp15_1 = B_0_3*w34;
1787 const double tmp9_1 = tmp3_0*w34;
1788 const double tmp35_1 = B_1_0*w36;
1789 const double tmp66_1 = B_0_3*w28;
1790 const double tmp28_1 = B_1_0*w42;
1791 const double tmp22_1 = B_1_0*w40;
1792 const double tmp16_1 = B_1_2*w29;
1793 const double tmp6_1 = tmp2_0*w35;
1794 const double tmp55_1 = B_1_3*w40;
1795 const double tmp50_1 = B_1_3*w42;
1796 const double tmp7_1 = tmp1_0*w29;
1797 const double tmp1_1 = tmp1_0*w32;
1798 const double tmp57_1 = B_0_3*w30;
1799 const double tmp18_1 = B_1_1*w32;
1800 const double tmp53_1 = B_1_0*w41;
1801 const double tmp61_1 = B_1_3*w36;
1802 const double tmp27_1 = B_0_3*w38;
1803 const double tmp64_1 = B_0_2*w30;
1804 const double tmp76_1 = B_0_1*w38;
1805 const double tmp39_1 = tmp2_0*w34;
1806 const double tmp62_1 = B_0_1*w31;
1807 const double tmp56_1 = B_0_0*w31;
1808 const double tmp49_1 = B_1_1*w36;
1809 const double tmp2_1 = B_0_2*w31;
1810 const double tmp23_1 = B_0_2*w33;
1811 const double tmp38_1 = B_1_1*w43;
1812 const double tmp74_1 = B_1_2*w41;
1813 const double tmp43_1 = B_1_1*w41;
1814 const double tmp58_1 = B_0_2*w28;
1815 const double tmp67_1 = B_0_0*w33;
1816 const double tmp33_1 = tmp0_0*w39;
1817 const double tmp4_1 = B_0_0*w28;
1818 const double tmp20_1 = B_0_0*w30;
1819 const double tmp13_1 = B_0_2*w38;
1820 const double tmp65_1 = B_1_2*w43;
1821 const double tmp0_1 = tmp0_0*w29;
1822 const double tmp41_1 = tmp3_0*w33;
1823 const double tmp73_1 = B_0_2*w37;
1824 const double tmp69_1 = B_0_0*w38;
1825 const double tmp48_1 = B_1_2*w39;
1826 const double tmp59_1 = B_0_1*w33;
1827 const double tmp17_1 = B_1_3*w41;
1828 const double tmp5_1 = B_0_3*w33;
1829 const double tmp3_1 = B_0_1*w30;
1830 const double tmp21_1 = B_0_1*w28;
1831 const double tmp42_1 = B_1_0*w29;
1832 const double tmp54_1 = B_1_2*w32;
1833 const double tmp60_1 = B_1_0*w39;
1834 const double tmp32_1 = tmp1_0*w36;
1835 const double tmp10_1 = B_0_1*w37;
1836 const double tmp14_1 = B_0_0*w35;
1837 const double tmp29_1 = B_0_1*w35;
1838 const double tmp26_1 = B_1_2*w36;
1839 const double tmp30_1 = B_1_3*w43;
1840 const double tmp70_1 = B_0_2*w35;
1841 const double tmp34_1 = B_1_3*w39;
1842 const double tmp51_1 = B_1_0*w43;
1843 const double tmp31_1 = B_0_2*w34;
1844 const double tmp45_1 = tmp3_0*w28;
1845 const double tmp11_1 = tmp1_0*w39;
1846 const double tmp52_1 = B_1_1*w29;
1847 const double tmp44_1 = B_1_3*w32;
1848 const double tmp25_1 = B_1_1*w39;
1849 const double tmp47_1 = tmp2_0*w33;
1850 const double tmp72_1 = B_1_3*w29;
1851 const double tmp40_1 = tmp2_0*w28;
1852 const double tmp46_1 = B_1_2*w40;
1853 const double tmp36_1 = B_1_2*w42;
1854 const double tmp24_1 = B_0_0*w37;
1855 const double tmp77_1 = B_0_3*w35;
1856 const double tmp68_1 = B_0_3*w37;
1857 const double tmp78_1 = B_0_0*w34;
1858 const double tmp12_1 = tmp0_0*w36;
1859 const double tmp75_1 = B_1_0*w32;
1860 EM_S[INDEX2(0,0,4)]+=tmp16_1 + tmp17_1 + tmp18_1 + tmp19_1 + tmp20_1 + tmp21_1 + tmp22_1 + tmp23_1;
1861 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp1_1 + tmp68_1 + tmp69_1 + tmp70_1 + tmp71_1;
1862 EM_S[INDEX2(2,0,4)]+=tmp45_1 + tmp47_1 + tmp48_1 + tmp49_1 + tmp50_1 + tmp51_1;
1863 EM_S[INDEX2(3,0,4)]+=tmp32_1 + tmp33_1 + tmp6_1 + tmp9_1;
1864 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1 + tmp4_1 + tmp5_1;
1865 EM_S[INDEX2(1,1,4)]+=tmp72_1 + tmp73_1 + tmp74_1 + tmp75_1 + tmp76_1 + tmp77_1 + tmp78_1 + tmp79_1;
1866 EM_S[INDEX2(2,1,4)]+=tmp32_1 + tmp33_1 + tmp40_1 + tmp41_1;
1867 EM_S[INDEX2(3,1,4)]+=tmp34_1 + tmp35_1 + tmp36_1 + tmp37_1 + tmp38_1 + tmp39_1;
1868 EM_S[INDEX2(0,2,4)]+=tmp42_1 + tmp43_1 + tmp44_1 + tmp45_1 + tmp46_1 + tmp47_1;
1869 EM_S[INDEX2(1,2,4)]+=tmp6_1 + tmp7_1 + tmp8_1 + tmp9_1;
1870 EM_S[INDEX2(2,2,4)]+=tmp60_1 + tmp61_1 + tmp62_1 + tmp63_1 + tmp64_1 + tmp65_1 + tmp66_1 + tmp67_1;
1871 EM_S[INDEX2(3,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1 + tmp14_1 + tmp15_1;
1872 EM_S[INDEX2(0,3,4)]+=tmp40_1 + tmp41_1 + tmp7_1 + tmp8_1;
1873 EM_S[INDEX2(1,3,4)]+=tmp37_1 + tmp39_1 + tmp52_1 + tmp53_1 + tmp54_1 + tmp55_1;
1874 EM_S[INDEX2(2,3,4)]+=tmp11_1 + tmp12_1 + tmp56_1 + tmp57_1 + tmp58_1 + tmp59_1;
1875 EM_S[INDEX2(3,3,4)]+=tmp24_1 + tmp25_1 + tmp26_1 + tmp27_1 + tmp28_1 + tmp29_1 + tmp30_1 + tmp31_1;
1876 } else { // constant data
1877 const double B_0 = B_p[0];
1878 const double B_1 = B_p[1];
1879 const double tmp0_1 = B_0*w44;
1880 const double tmp1_1 = B_1*w45;
1881 const double tmp2_1 = B_0*w46;
1882 const double tmp3_1 = B_0*w47;
1883 const double tmp4_1 = B_1*w48;
1884 const double tmp5_1 = B_1*w49;
1885 const double tmp6_1 = B_1*w50;
1886 const double tmp7_1 = B_0*w51;
1887 EM_S[INDEX2(0,0,4)]+=tmp0_1 + tmp5_1;
1888 EM_S[INDEX2(1,0,4)]+=tmp1_1 + tmp3_1;
1889 EM_S[INDEX2(2,0,4)]+=tmp6_1 + tmp7_1;
1890 EM_S[INDEX2(3,0,4)]+=tmp2_1 + tmp4_1;
1891 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1;
1892 EM_S[INDEX2(1,1,4)]+=tmp3_1 + tmp5_1;
1893 EM_S[INDEX2(2,1,4)]+=tmp4_1 + tmp7_1;
1894 EM_S[INDEX2(3,1,4)]+=tmp2_1 + tmp6_1;
1895 EM_S[INDEX2(0,2,4)]+=tmp5_1 + tmp7_1;
1896 EM_S[INDEX2(1,2,4)]+=tmp1_1 + tmp2_1;
1897 EM_S[INDEX2(2,2,4)]+=tmp0_1 + tmp6_1;
1898 EM_S[INDEX2(3,2,4)]+=tmp3_1 + tmp4_1;
1899 EM_S[INDEX2(0,3,4)]+=tmp1_1 + tmp7_1;
1900 EM_S[INDEX2(1,3,4)]+=tmp2_1 + tmp5_1;
1901 EM_S[INDEX2(2,3,4)]+=tmp0_1 + tmp4_1;
1902 EM_S[INDEX2(3,3,4)]+=tmp3_1 + tmp6_1;
1903 }
1904 }
1905 ///////////////
1906 // process C //
1907 ///////////////
1908 if (!C.isEmpty()) {
1909 add_EM_S=true;
1910 const double* C_p=const_cast<escript::Data*>(&C)->getSampleDataRO(e);
1911 if (C.actsExpanded()) {
1912 const double C_0_0 = C_p[INDEX2(0,0,2)];
1913 const double C_1_0 = C_p[INDEX2(1,0,2)];
1914 const double C_0_1 = C_p[INDEX2(0,1,2)];
1915 const double C_1_1 = C_p[INDEX2(1,1,2)];
1916 const double C_0_2 = C_p[INDEX2(0,2,2)];
1917 const double C_1_2 = C_p[INDEX2(1,2,2)];
1918 const double C_0_3 = C_p[INDEX2(0,3,2)];
1919 const double C_1_3 = C_p[INDEX2(1,3,2)];
1920 const double tmp0_0 = C_1_0 + C_1_1;
1921 const double tmp1_0 = C_1_2 + C_1_3;
1922 const double tmp2_0 = C_0_1 + C_0_3;
1923 const double tmp3_0 = C_0_0 + C_0_2;
1924 const double tmp64_1 = C_0_2*w30;
1925 const double tmp14_1 = C_0_2*w28;
1926 const double tmp19_1 = C_0_3*w31;
1927 const double tmp22_1 = C_1_0*w40;
1928 const double tmp37_1 = tmp3_0*w35;
1929 const double tmp29_1 = C_0_1*w35;
1930 const double tmp73_1 = C_0_2*w37;
1931 const double tmp74_1 = C_1_2*w41;
1932 const double tmp52_1 = C_1_3*w39;
1933 const double tmp25_1 = C_1_1*w39;
1934 const double tmp62_1 = C_0_1*w31;
1935 const double tmp79_1 = C_1_1*w40;
1936 const double tmp43_1 = C_1_1*w36;
1937 const double tmp27_1 = C_0_3*w38;
1938 const double tmp28_1 = C_1_0*w42;
1939 const double tmp63_1 = C_1_1*w42;
1940 const double tmp59_1 = C_0_3*w34;
1941 const double tmp72_1 = C_1_3*w29;
1942 const double tmp40_1 = tmp2_0*w35;
1943 const double tmp13_1 = C_0_3*w30;
1944 const double tmp51_1 = C_1_2*w40;
1945 const double tmp54_1 = C_1_2*w42;
1946 const double tmp12_1 = C_0_0*w31;
1947 const double tmp2_1 = tmp1_0*w32;
1948 const double tmp68_1 = C_0_2*w31;
1949 const double tmp75_1 = C_1_0*w32;
1950 const double tmp49_1 = C_1_1*w41;
1951 const double tmp4_1 = C_0_2*w35;
1952 const double tmp66_1 = C_0_3*w28;
1953 const double tmp56_1 = C_0_1*w37;
1954 const double tmp5_1 = C_0_1*w34;
1955 const double tmp38_1 = tmp2_0*w34;
1956 const double tmp76_1 = C_0_1*w38;
1957 const double tmp21_1 = C_0_1*w28;
1958 const double tmp69_1 = C_0_1*w30;
1959 const double tmp53_1 = C_1_0*w36;
1960 const double tmp42_1 = C_1_2*w39;
1961 const double tmp32_1 = tmp1_0*w29;
1962 const double tmp45_1 = C_1_0*w43;
1963 const double tmp33_1 = tmp0_0*w32;
1964 const double tmp35_1 = C_1_0*w41;
1965 const double tmp26_1 = C_1_2*w36;
1966 const double tmp67_1 = C_0_0*w33;
1967 const double tmp31_1 = C_0_2*w34;
1968 const double tmp20_1 = C_0_0*w30;
1969 const double tmp70_1 = C_0_0*w28;
1970 const double tmp8_1 = tmp0_0*w39;
1971 const double tmp30_1 = C_1_3*w43;
1972 const double tmp0_1 = tmp0_0*w29;
1973 const double tmp17_1 = C_1_3*w41;
1974 const double tmp58_1 = C_0_0*w35;
1975 const double tmp9_1 = tmp3_0*w33;
1976 const double tmp61_1 = C_1_3*w36;
1977 const double tmp41_1 = tmp3_0*w34;
1978 const double tmp50_1 = C_1_3*w32;
1979 const double tmp18_1 = C_1_1*w32;
1980 const double tmp6_1 = tmp1_0*w36;
1981 const double tmp3_1 = C_0_0*w38;
1982 const double tmp34_1 = C_1_1*w29;
1983 const double tmp77_1 = C_0_3*w35;
1984 const double tmp65_1 = C_1_2*w43;
1985 const double tmp71_1 = C_0_3*w33;
1986 const double tmp55_1 = C_1_1*w43;
1987 const double tmp46_1 = tmp3_0*w28;
1988 const double tmp24_1 = C_0_0*w37;
1989 const double tmp10_1 = tmp1_0*w39;
1990 const double tmp48_1 = C_1_0*w29;
1991 const double tmp15_1 = C_0_1*w33;
1992 const double tmp36_1 = C_1_2*w32;
1993 const double tmp60_1 = C_1_0*w39;
1994 const double tmp47_1 = tmp2_0*w33;
1995 const double tmp16_1 = C_1_2*w29;
1996 const double tmp1_1 = C_0_3*w37;
1997 const double tmp7_1 = tmp2_0*w28;
1998 const double tmp39_1 = C_1_3*w40;
1999 const double tmp44_1 = C_1_3*w42;
2000 const double tmp57_1 = C_0_2*w38;
2001 const double tmp78_1 = C_0_0*w34;
2002 const double tmp11_1 = tmp0_0*w36;
2003 const double tmp23_1 = C_0_2*w33;
2004 EM_S[INDEX2(0,0,4)]+=tmp16_1 + tmp17_1 + tmp18_1 + tmp19_1 + tmp20_1 + tmp21_1 + tmp22_1 + tmp23_1;
2005 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp2_1 + tmp68_1 + tmp69_1 + tmp70_1 + tmp71_1;
2006 EM_S[INDEX2(2,0,4)]+=tmp46_1 + tmp47_1 + tmp48_1 + tmp49_1 + tmp50_1 + tmp51_1;
2007 EM_S[INDEX2(3,0,4)]+=tmp32_1 + tmp33_1 + tmp7_1 + tmp9_1;
2008 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1 + tmp4_1 + tmp5_1;
2009 EM_S[INDEX2(1,1,4)]+=tmp72_1 + tmp73_1 + tmp74_1 + tmp75_1 + tmp76_1 + tmp77_1 + tmp78_1 + tmp79_1;
2010 EM_S[INDEX2(2,1,4)]+=tmp32_1 + tmp33_1 + tmp40_1 + tmp41_1;
2011 EM_S[INDEX2(3,1,4)]+=tmp34_1 + tmp35_1 + tmp36_1 + tmp37_1 + tmp38_1 + tmp39_1;
2012 EM_S[INDEX2(0,2,4)]+=tmp42_1 + tmp43_1 + tmp44_1 + tmp45_1 + tmp46_1 + tmp47_1;
2013 EM_S[INDEX2(1,2,4)]+=tmp6_1 + tmp7_1 + tmp8_1 + tmp9_1;
2014 EM_S[INDEX2(2,2,4)]+=tmp60_1 + tmp61_1 + tmp62_1 + tmp63_1 + tmp64_1 + tmp65_1 + tmp66_1 + tmp67_1;
2015 EM_S[INDEX2(3,2,4)]+=tmp10_1 + tmp11_1 + tmp12_1 + tmp13_1 + tmp14_1 + tmp15_1;
2016 EM_S[INDEX2(0,3,4)]+=tmp40_1 + tmp41_1 + tmp6_1 + tmp8_1;
2017 EM_S[INDEX2(1,3,4)]+=tmp37_1 + tmp38_1 + tmp52_1 + tmp53_1 + tmp54_1 + tmp55_1;
2018 EM_S[INDEX2(2,3,4)]+=tmp10_1 + tmp11_1 + tmp56_1 + tmp57_1 + tmp58_1 + tmp59_1;
2019 EM_S[INDEX2(3,3,4)]+=tmp24_1 + tmp25_1 + tmp26_1 + tmp27_1 + tmp28_1 + tmp29_1 + tmp30_1 + tmp31_1;
2020 } else { // constant data
2021 const double C_0 = C_p[0];
2022 const double C_1 = C_p[1];
2023 const double tmp0_1 = C_0*w47;
2024 const double tmp1_1 = C_1*w45;
2025 const double tmp2_1 = C_1*w48;
2026 const double tmp3_1 = C_0*w51;
2027 const double tmp4_1 = C_0*w44;
2028 const double tmp5_1 = C_1*w49;
2029 const double tmp6_1 = C_1*w50;
2030 const double tmp7_1 = C_0*w46;
2031 EM_S[INDEX2(0,0,4)]+=tmp4_1 + tmp5_1;
2032 EM_S[INDEX2(1,0,4)]+=tmp1_1 + tmp4_1;
2033 EM_S[INDEX2(2,0,4)]+=tmp3_1 + tmp5_1;
2034 EM_S[INDEX2(3,0,4)]+=tmp1_1 + tmp3_1;
2035 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1;
2036 EM_S[INDEX2(1,1,4)]+=tmp0_1 + tmp5_1;
2037 EM_S[INDEX2(2,1,4)]+=tmp1_1 + tmp7_1;
2038 EM_S[INDEX2(3,1,4)]+=tmp5_1 + tmp7_1;
2039 EM_S[INDEX2(0,2,4)]+=tmp3_1 + tmp6_1;
2040 EM_S[INDEX2(1,2,4)]+=tmp2_1 + tmp3_1;
2041 EM_S[INDEX2(2,2,4)]+=tmp4_1 + tmp6_1;
2042 EM_S[INDEX2(3,2,4)]+=tmp2_1 + tmp4_1;
2043 EM_S[INDEX2(0,3,4)]+=tmp2_1 + tmp7_1;
2044 EM_S[INDEX2(1,3,4)]+=tmp6_1 + tmp7_1;
2045 EM_S[INDEX2(2,3,4)]+=tmp0_1 + tmp2_1;
2046 EM_S[INDEX2(3,3,4)]+=tmp0_1 + tmp6_1;
2047 }
2048 }
2049 ///////////////
2050 // process D //
2051 ///////////////
2052 if (!D.isEmpty()) {
2053 add_EM_S=true;
2054 const double* D_p=const_cast<escript::Data*>(&D)->getSampleDataRO(e);
2055 if (D.actsExpanded()) {
2056 const double D_0 = D_p[0];
2057 const double D_1 = D_p[1];
2058 const double D_2 = D_p[2];
2059 const double D_3 = D_p[3];
2060 const double tmp0_0 = D_0 + D_1;
2061 const double tmp1_0 = D_2 + D_3;
2062 const double tmp2_0 = D_0 + D_1 + D_2 + D_3;
2063 const double tmp3_0 = D_1 + D_2;
2064 const double tmp4_0 = D_1 + D_3;
2065 const double tmp5_0 = D_0 + D_2;
2066 const double tmp6_0 = D_0 + D_3;
2067 const double tmp0_1 = tmp0_0*w52;
2068 const double tmp1_1 = tmp1_0*w53;
2069 const double tmp2_1 = tmp2_0*w54;
2070 const double tmp3_1 = tmp1_0*w52;
2071 const double tmp4_1 = tmp0_0*w53;
2072 const double tmp5_1 = tmp3_0*w54;
2073 const double tmp6_1 = D_0*w55;
2074 const double tmp7_1 = D_3*w56;
2075 const double tmp8_1 = D_3*w55;
2076 const double tmp9_1 = D_0*w56;
2077 const double tmp10_1 = tmp4_0*w52;
2078 const double tmp11_1 = tmp5_0*w53;
2079 const double tmp12_1 = tmp5_0*w52;
2080 const double tmp13_1 = tmp4_0*w53;
2081 const double tmp14_1 = tmp6_0*w54;
2082 const double tmp15_1 = D_2*w55;
2083 const double tmp16_1 = D_1*w56;
2084 const double tmp17_1 = D_1*w55;
2085 const double tmp18_1 = D_2*w56;
2086 EM_S[INDEX2(0,0,4)]+=tmp5_1 + tmp6_1 + tmp7_1;
2087 EM_S[INDEX2(1,0,4)]+=tmp0_1 + tmp1_1;
2088 EM_S[INDEX2(2,0,4)]+=tmp12_1 + tmp13_1;
2089 EM_S[INDEX2(3,0,4)]+=tmp2_1;
2090 EM_S[INDEX2(0,1,4)]+=tmp0_1 + tmp1_1;
2091 EM_S[INDEX2(1,1,4)]+=tmp14_1 + tmp17_1 + tmp18_1;
2092 EM_S[INDEX2(2,1,4)]+=tmp2_1;
2093 EM_S[INDEX2(3,1,4)]+=tmp10_1 + tmp11_1;
2094 EM_S[INDEX2(0,2,4)]+=tmp12_1 + tmp13_1;
2095 EM_S[INDEX2(1,2,4)]+=tmp2_1;
2096 EM_S[INDEX2(2,2,4)]+=tmp14_1 + tmp15_1 + tmp16_1;
2097 EM_S[INDEX2(3,2,4)]+=tmp3_1 + tmp4_1;
2098 EM_S[INDEX2(0,3,4)]+=tmp2_1;
2099 EM_S[INDEX2(1,3,4)]+=tmp10_1 + tmp11_1;
2100 EM_S[INDEX2(2,3,4)]+=tmp3_1 + tmp4_1;
2101 EM_S[INDEX2(3,3,4)]+=tmp5_1 + tmp8_1 + tmp9_1;
2102 } else { // constant data
2103 const double tmp0_1 = D_p[0]*w57;
2104 const double tmp1_1 = D_p[0]*w58;
2105 const double tmp2_1 = D_p[0]*w59;
2106 EM_S[INDEX2(0,0,4)]+=tmp2_1;
2107 EM_S[INDEX2(1,0,4)]+=tmp0_1;
2108 EM_S[INDEX2(2,0,4)]+=tmp0_1;
2109 EM_S[INDEX2(3,0,4)]+=tmp1_1;
2110 EM_S[INDEX2(0,1,4)]+=tmp0_1;
2111 EM_S[INDEX2(1,1,4)]+=tmp2_1;
2112 EM_S[INDEX2(2,1,4)]+=tmp1_1;
2113 EM_S[INDEX2(3,1,4)]+=tmp0_1;
2114 EM_S[INDEX2(0,2,4)]+=tmp0_1;
2115 EM_S[INDEX2(1,2,4)]+=tmp1_1;
2116 EM_S[INDEX2(2,2,4)]+=tmp2_1;
2117 EM_S[INDEX2(3,2,4)]+=tmp0_1;
2118 EM_S[INDEX2(0,3,4)]+=tmp1_1;
2119 EM_S[INDEX2(1,3,4)]+=tmp0_1;
2120 EM_S[INDEX2(2,3,4)]+=tmp0_1;
2121 EM_S[INDEX2(3,3,4)]+=tmp2_1;
2122 }
2123 }
2124 ///////////////
2125 // process X //
2126 ///////////////
2127 if (!X.isEmpty()) {
2128 add_EM_F=true;
2129 const double* X_p=const_cast<escript::Data*>(&X)->getSampleDataRO(e);
2130 if (X.actsExpanded()) {
2131 const double X_0_0 = X_p[INDEX2(0,0,2)];
2132 const double X_1_0 = X_p[INDEX2(1,0,2)];
2133 const double X_0_1 = X_p[INDEX2(0,1,2)];
2134 const double X_1_1 = X_p[INDEX2(1,1,2)];
2135 const double X_0_2 = X_p[INDEX2(0,2,2)];
2136 const double X_1_2 = X_p[INDEX2(1,2,2)];
2137 const double X_0_3 = X_p[INDEX2(0,3,2)];
2138 const double X_1_3 = X_p[INDEX2(1,3,2)];
2139 const double tmp0_0 = X_0_2 + X_0_3;
2140 const double tmp1_0 = X_1_1 + X_1_3;
2141 const double tmp2_0 = X_1_0 + X_1_2;
2142 const double tmp3_0 = X_0_0 + X_0_1;
2143 const double tmp0_1 = tmp0_0*w63;
2144 const double tmp1_1 = tmp1_0*w62;
2145 const double tmp2_1 = tmp2_0*w61;
2146 const double tmp3_1 = tmp3_0*w60;
2147 const double tmp4_1 = tmp0_0*w65;
2148 const double tmp5_1 = tmp3_0*w64;
2149 const double tmp6_1 = tmp2_0*w62;
2150 const double tmp7_1 = tmp1_0*w61;
2151 const double tmp8_1 = tmp2_0*w66;
2152 const double tmp9_1 = tmp3_0*w63;
2153 const double tmp10_1 = tmp0_0*w60;
2154 const double tmp11_1 = tmp1_0*w67;
2155 const double tmp12_1 = tmp1_0*w66;
2156 const double tmp13_1 = tmp3_0*w65;
2157 const double tmp14_1 = tmp0_0*w64;
2158 const double tmp15_1 = tmp2_0*w67;
2159 EM_F[0]+=tmp0_1 + tmp1_1 + tmp2_1 + tmp3_1;
2160 EM_F[1]+=tmp4_1 + tmp5_1 + tmp6_1 + tmp7_1;
2161 EM_F[2]+=tmp10_1 + tmp11_1 + tmp8_1 + tmp9_1;
2162 EM_F[3]+=tmp12_1 + tmp13_1 + tmp14_1 + tmp15_1;
2163 } else { // constant data
2164 const double tmp0_1 = X_p[1]*w69;
2165 const double tmp1_1 = X_p[0]*w68;
2166 const double tmp2_1 = X_p[0]*w70;
2167 const double tmp3_1 = X_p[1]*w71;
2168 EM_F[0]+=tmp0_1 + tmp1_1;
2169 EM_F[1]+=tmp0_1 + tmp2_1;
2170 EM_F[2]+=tmp1_1 + tmp3_1;
2171 EM_F[3]+=tmp2_1 + tmp3_1;
2172 }
2173 }
2174 ///////////////
2175 // process Y //
2176 ///////////////
2177 if (!Y.isEmpty()) {
2178 add_EM_F=true;
2179 const double* Y_p=const_cast<escript::Data*>(&Y)->getSampleDataRO(e);
2180 if (Y.actsExpanded()) {
2181 const double Y_0 = Y_p[0];
2182 const double Y_1 = Y_p[1];
2183 const double Y_2 = Y_p[2];
2184 const double Y_3 = Y_p[3];
2185 const double tmp0_0 = Y_1 + Y_2;
2186 const double tmp1_0 = Y_0 + Y_3;
2187 const double tmp0_1 = Y_0*w72;
2188 const double tmp1_1 = tmp0_0*w73;
2189 const double tmp2_1 = Y_3*w74;
2190 const double tmp3_1 = Y_1*w72;
2191 const double tmp4_1 = tmp1_0*w73;
2192 const double tmp5_1 = Y_2*w74;
2193 const double tmp6_1 = Y_2*w72;
2194 const double tmp7_1 = Y_1*w74;
2195 const double tmp8_1 = Y_3*w72;
2196 const double tmp9_1 = Y_0*w74;
2197 EM_F[0]+=tmp0_1 + tmp1_1 + tmp2_1;
2198 EM_F[1]+=tmp3_1 + tmp4_1 + tmp5_1;
2199 EM_F[2]+=tmp4_1 + tmp6_1 + tmp7_1;
2200 EM_F[3]+=tmp1_1 + tmp8_1 + tmp9_1;
2201 } else { // constant data
2202 const double tmp0_1 = Y_p[0]*w75;
2203 EM_F[0]+=tmp0_1;
2204 EM_F[1]+=tmp0_1;
2205 EM_F[2]+=tmp0_1;
2206 EM_F[3]+=tmp0_1;
2207 }
2208 }
2209
2210 // add to matrix (if add_EM_S) and RHS (if add_EM_F)
2211 const index_t firstNode=m_N0*k1+k0;
2212 ad