Nektar++
Loading...
Searching...
No Matches
ProcessEquiSpacedOutput.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: ProcessEquiSpacedOutput.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Set up fields as interpolation to equispaced output
32//
33////////////////////////////////////////////////////////////////////////////////
34#include <iostream>
35#include <string>
36using namespace std;
37
42
44
45namespace Nektar::FieldUtils
46{
47
50 ModuleKey(eProcessModule, "equispacedoutput"),
52 "Write data as equi-spaced output using simplices to represent the "
53 "data for connecting points");
54
56 : ProcessModule(f)
57{
58 m_config["tetonly"] =
59 ConfigOption(true, "0", "Only process tetrahedral elements");
60
61 m_config["modalenergy"] =
62 ConfigOption(true, "0", "Write output as modal energy");
63}
64
68
69void ProcessEquiSpacedOutput::v_Process(po::variables_map &vm)
70{
71 m_f->SetUpExp(vm);
72
73 int nel = m_f->m_exp[0]->GetExpSize();
74 if (!nel)
75 {
76 // Create empty PtsField
77 int nfields = m_f->m_variables.size();
78 int coordim = 3;
79
80 Array<OneD, Array<OneD, NekDouble>> pts(nfields + coordim);
81 for (int i = 0; i < nfields + coordim; ++i)
82 {
83 pts[i] = Array<OneD, NekDouble>(0);
84 }
85 vector<Array<OneD, int>> ptsConn;
86
87 m_f->m_fieldPts =
89 coordim, m_f->m_variables, pts);
90 m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsTetBlock);
91 m_f->m_fieldPts->SetConnectivity(ptsConn);
92 return;
93 }
94
95 int coordim = m_f->m_exp[0]->GetCoordim(0);
96 int shapedim = m_f->m_exp[0]->GetExp(0)->GetShapeDimension();
97 int npts = m_f->m_exp[0]->GetTotPoints();
99
100 // Check if we have a homogeneous expansion
101 bool homogeneous1D = false;
102 if (m_f->m_numHomogeneousDir == 1)
103 {
104 coordim++;
105 shapedim++;
106 homogeneous1D = true;
107 }
108 else if (m_f->m_numHomogeneousDir == 2)
109 {
110 ASSERTL0(false, "Homegeneous2D case not supported");
111 }
112
113 // set up the number of points in each element
114 int newpoints = 0;
115 int newtotpoints = 0;
116
117 Array<OneD, int> conn;
118 int prevNcoeffs = 0;
119 int prevNpoints = 0;
120 int cnt = 0;
121
122 // identify face 1 connectivity for prisms
123 map<int, StdRegions::Orientation> face0orient;
124 set<int> prismorient;
126
127 // prepare PtsField
128 vector<int> ppe;
129 vector<Array<OneD, int>> ptsConn;
130 int nfields;
131
132 for (int i = 0; i < nel; ++i)
133 {
134 e = m_f->m_exp[0]->GetExp(i);
135 if (e->DetShapeType() == LibUtilities::ePrism)
136 {
137 StdRegions::Orientation forient = e->GetTraceOrient(0);
138 int fid = e->GetGeom()->GetFid(0);
139 if (face0orient.count(fid))
140 { // face 1 meeting face 1 so reverse this id
141 prismorient.insert(i);
142 }
143 else
144 {
145 // just store if Dir 1 is fwd or bwd
146 if ((forient == StdRegions::eDir1BwdDir1_Dir2FwdDir2) ||
150 {
151 face0orient[fid] = StdRegions::eBwd;
152 }
153 else
154 {
155 face0orient[fid] = StdRegions::eFwd;
156 }
157 }
158 }
159 }
160
161 for (int i = 0; i < nel; ++i)
162 {
163 e = m_f->m_exp[0]->GetExp(i);
164 if (e->DetShapeType() == LibUtilities::ePrism)
165 {
166 int fid = e->GetGeom()->GetFid(2);
167 // check to see if face 2 meets face 1
168 if (face0orient.count(fid))
169 {
170 // check to see how face 2 is orientated
171 StdRegions::Orientation forient2 = e->GetTraceOrient(2);
172 StdRegions::Orientation forient0 = face0orient[fid];
173
174 // If dir 1 or forient2 is bwd then check agains
175 // face 1 value
176 if ((forient2 == StdRegions::eDir1BwdDir1_Dir2FwdDir2) ||
180 {
181 if (forient0 == StdRegions::eFwd)
182 {
183 prismorient.insert(i);
184 }
185 }
186 else
187 {
188 if (forient0 == StdRegions::eBwd)
189 {
190 prismorient.insert(i);
191 }
192 }
193 }
194 }
195 }
196
197 for (int i = 0; i < nel; ++i)
198 {
199 e = m_f->m_exp[0]->GetExp(i);
200 if (m_config["tetonly"].m_beenSet && m_config["tetonly"].as<bool>())
201 {
202 if (m_f->m_exp[0]->GetExp(i)->DetShapeType() !=
204 {
205 continue;
206 }
207 }
208
209 switch (e->DetShapeType())
210 {
212 {
213 int npoints0 = e->GetBasis(0)->GetNumPoints();
214
215 newpoints =
217 }
218 break;
220 {
221 int np0 = e->GetBasis(0)->GetNumPoints();
222 int np1 = e->GetBasis(1)->GetNumPoints();
223 int np = max(np0, np1);
224 newpoints =
226 }
227 break;
229 {
230 int np0 = e->GetBasis(0)->GetNumPoints();
231 int np1 = e->GetBasis(1)->GetNumPoints();
232 int np = max(np0, np1);
233
234 newpoints =
236 }
237 break;
239 {
240 int np0 = e->GetBasis(0)->GetNumPoints();
241 int np1 = e->GetBasis(1)->GetNumPoints();
242 int np2 = e->GetBasis(2)->GetNumPoints();
243 int np = max(np0, max(np1, np2));
244
246 np, np, np);
247 }
248 break;
250 {
251 int np0 = e->GetBasis(0)->GetNumPoints();
252 int np1 = e->GetBasis(1)->GetNumPoints();
253 int np2 = e->GetBasis(2)->GetNumPoints();
254 int np = max(np0, max(np1, np2));
255
257 np, np, np);
258 }
259 break;
261 {
262 int np0 = e->GetBasis(0)->GetNumPoints();
263 int np1 = e->GetBasis(1)->GetNumPoints();
264 int np2 = e->GetBasis(2)->GetNumPoints();
265 int np = max(np0, max(np1, np2));
266
268 np, np, np);
269 }
270 break;
272 {
273 int np0 = e->GetBasis(0)->GetNumPoints();
274 int np1 = e->GetBasis(1)->GetNumPoints();
275 int np2 = e->GetBasis(2)->GetNumPoints();
276 int np = max(np0, max(np1, np2));
277
279 np, np, np);
280 }
281 break;
282 default:
283 {
284 ASSERTL0(false, "Points not known");
285 }
286 }
287
288 ppe.push_back(newpoints);
289 newtotpoints += newpoints;
290
291 if (e->DetShapeType() == LibUtilities::ePrism)
292 {
293 bool standard = true;
294
295 if (prismorient.count(i))
296 {
297 standard = false; // reverse direction
298 }
299
300 e->GetSimplexEquiSpacedConnectivity(conn, standard);
301 }
302 else
303 {
304
305 if ((prevNcoeffs != e->GetNcoeffs()) ||
306 (prevNpoints != e->GetTotPoints()))
307 {
308 prevNcoeffs = e->GetNcoeffs();
309 prevNpoints = e->GetTotPoints();
310
311 e->GetSimplexEquiSpacedConnectivity(conn);
312 }
313 }
314 Array<OneD, int> newconn(conn.size());
315 for (int j = 0; j < conn.size(); ++j)
316 {
317 newconn[j] = conn[j] + cnt;
318 }
319
320 ptsConn.push_back(newconn);
321 cnt += newpoints;
322 }
323
324 nfields = m_f->m_variables.size();
325
326 Array<OneD, Array<OneD, NekDouble>> pts(nfields + coordim);
327
328 for (int i = 0; i < nfields + coordim; ++i)
329 {
330 pts[i] = Array<OneD, NekDouble>(newtotpoints);
331 }
332
333 // Interpolate coordinates
334 for (int i = 0; i < coordim; ++i)
335 {
336 coords[i] = Array<OneD, NekDouble>(npts);
337 }
338
339 for (int i = coordim; i < 3; ++i)
340 {
341 coords[i] = NullNekDouble1DArray;
342 }
343
344 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
345
347
348 for (int n = 0; n < coordim; ++n)
349 {
350 cnt = 0;
351 int cnt1 = 0;
352 for (int i = 0; i < nel; ++i)
353 {
354 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
355 coords[n] + cnt, tmp = pts[n] + cnt1);
356 cnt1 += ppe[i];
357 cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
358 }
359 }
360
361 for (int n = 0; n < m_f->m_variables.size(); ++n)
362 {
363 cnt = 0;
364 int cnt1 = 0;
365
366 if (m_config["modalenergy"].m_beenSet &&
367 m_config["modalenergy"].as<bool>())
368 {
369 Array<OneD, const NekDouble> phys = m_f->m_exp[n]->GetPhys();
370 for (int i = 0; i < nel; ++i)
371 {
372 GenOrthoModes(i, phys + cnt, tmp = pts[coordim + n] + cnt1);
373 cnt1 += ppe[i];
374 cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
375 }
376 }
377 else
378 {
379 Array<OneD, const NekDouble> phys = m_f->m_exp[n]->GetPhys();
380 for (int i = 0; i < nel; ++i)
381 {
382 m_f->m_exp[0]->GetExp(i)->PhysInterpToSimplexEquiSpaced(
383 phys + cnt, tmp = pts[coordim + n] + cnt1);
384 cnt1 += ppe[i];
385 cnt += m_f->m_exp[0]->GetExp(i)->GetTotPoints();
386 }
387 }
388 }
389
391 coordim, m_f->m_variables, pts);
392 if (shapedim == 1)
393 {
394 m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsSegBlock);
395 }
396
397 if (shapedim == 2)
398 {
399 m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsTriBlock);
400 }
401 else if (shapedim == 3)
402 {
403 m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsTetBlock);
404 }
405 m_f->m_fieldPts->SetConnectivity(ptsConn);
406 if (homogeneous1D)
407 {
409 }
410
411 // Save points per element in the point field
412 m_f->m_fieldPts->SetPointsPerElement(ppe);
413
414 // Clear m_exp
415 m_f->m_exp = vector<MultiRegions::ExpListSharedPtr>();
416}
417
419{
421 int nel = m_f->m_exp[0]->GetPlane(0)->GetExpSize();
422 int nPlanes = m_f->m_exp[0]->GetZIDs().size();
423 int npts = m_f->m_fieldPts->GetNpoints();
424 int nptsPerPlane = npts / nPlanes;
425 int coordim = 3;
426
427 if (m_f->m_exp[0]->GetHomogeneousBasis()->GetBasisType() ==
429 m_f->m_exp[0]->GetHomogeneousBasis()->GetPointsType() ==
431 {
432 // Write points with extra plane
434 m_f->m_fieldPts->GetPts(pts);
435 Array<OneD, Array<OneD, NekDouble>> newPts(pts.size());
436 for (int i = 0; i < pts.size(); i++)
437 {
438 newPts[i] = Array<OneD, NekDouble>(npts + nptsPerPlane);
439 // Copy old points
440 Vmath::Vcopy(npts, pts[i], 1, newPts[i], 1);
441 // Get new plane
442 Array<OneD, NekDouble> extraPlane(nptsPerPlane, newPts[i] + npts);
443 if (m_f->m_comm->GetColumnComm()->GetSize() == 1)
444 {
445 if (i == (coordim - 1))
446 {
447 // z-coordinate
448 Vmath::Sadd(nptsPerPlane, m_f->m_exp[0]->GetHomoLen(),
449 pts[i], 1, extraPlane, 1);
450 }
451 else
452 {
453 Vmath::Vcopy(nptsPerPlane, pts[i], 1, extraPlane, 1);
454 }
455 }
456 else
457 {
458 // Determine to and from rank for communication
459 int size = m_f->m_comm->GetColumnComm()->GetSize();
460 int rank = m_f->m_comm->GetColumnComm()->GetRank();
461 int fromRank = (rank + 1) % size;
462 int toRank = (rank == 0) ? size - 1 : rank - 1;
463 // Communicate using SendRecv
464 Array<OneD, NekDouble> send(nptsPerPlane, pts[i]);
465 m_f->m_comm->GetColumnComm()->SendRecv(toRank, send, fromRank,
466 extraPlane);
467 // Adjust z-coordinate of last process
468 if (i == (coordim - 1) && (rank == size - 1))
469 {
470 Vmath::Sadd(nptsPerPlane, m_f->m_exp[0]->GetHomoLen(),
471 extraPlane, 1, extraPlane, 1);
472 }
473 }
474 }
475 m_f->m_fieldPts->SetPts(newPts);
476 // Now extend plane connectivity
477 vector<Array<OneD, int>> oldConn;
478 Array<OneD, int> conn;
479 m_f->m_fieldPts->GetConnectivity(oldConn);
480 vector<Array<OneD, int>> newConn = oldConn;
481 int connPerPlane = oldConn.size() / nPlanes;
482 for (int i = 0; i < connPerPlane; i++)
483 {
484 conn = Array<OneD, int>(oldConn[i].size());
485 for (int j = 0; j < conn.size(); j++)
486 {
487 conn[j] = oldConn[i][j] + npts;
488 }
489 newConn.push_back(conn);
490 }
491 m_f->m_fieldPts->SetConnectivity(newConn);
492
493 nPlanes++;
494 npts += nptsPerPlane;
495 }
496
497 vector<Array<OneD, int>> oldConn;
498 vector<Array<OneD, int>> newConn;
499 Array<OneD, int> conn;
500 m_f->m_fieldPts->GetConnectivity(oldConn);
501
502 // 2DH1D case
503 if (m_f->m_fieldPts->GetPtsType() == LibUtilities::ePtsTriBlock)
504 {
505 for (int i = 0; i < nel; ++i)
506 {
507 int nLines = oldConn[i].size() / 2;
508 // Create array for new connectivity
509 // (2 triangles between each plane for each line)
510 conn = Array<OneD, int>(2 * 3 * nLines * (nPlanes - 1));
511 int cnt = 0;
512 for (int n = 0; n < nLines; n++)
513 {
514 // Define new connectivity with triangles
515 for (int p = 0; p < nPlanes - 1; p++)
516 {
517 conn[cnt++] = oldConn[i + p * nel][2 * n + 0];
518 conn[cnt++] = oldConn[i + p * nel][2 * n + 1];
519 conn[cnt++] = oldConn[i + (p + 1) * nel][2 * n + 0];
520
521 conn[cnt++] = oldConn[i + p * nel][2 * n + 1];
522 conn[cnt++] = oldConn[i + (p + 1) * nel][2 * n + 0];
523 conn[cnt++] = oldConn[i + (p + 1) * nel][2 * n + 1];
524 }
525 }
526 newConn.push_back(conn);
527 }
528 }
529 else if (m_f->m_fieldPts->GetPtsType() == LibUtilities::ePtsTetBlock)
530 {
531 // Get maximum number of points per direction in the mesh
532 int maxN = 0;
533 for (int i = 0; i < nel; ++i)
534 {
535 e = m_f->m_exp[0]->GetPlane(0)->GetExp(i);
536
537 int np0 = e->GetBasis(0)->GetNumPoints();
538 int np1 = e->GetBasis(1)->GetNumPoints();
539 int np = max(np0, np1);
540 maxN = max(np, maxN);
541 }
542
543 // Create a global numbering for points in plane 0
544 Array<OneD, int> vId(nptsPerPlane);
545 int cnt1 = 0;
546 int cnt2 = 0;
547 for (int i = 0; i < nel; ++i)
548 {
549 e = m_f->m_exp[0]->GetPlane(0)->GetExp(i);
550 int np0 = e->GetBasis(0)->GetNumPoints();
551 int np1 = e->GetBasis(1)->GetNumPoints();
552 int np = max(np0, np1);
553 switch (e->DetShapeType())
554 {
556 {
557 int newpoints =
559 np);
560
561 // Vertex numbering
562 vId[cnt1] = e->GetGeom()->GetVid(0);
563 vId[cnt1 + np - 1] = e->GetGeom()->GetVid(1);
564 vId[cnt1 + newpoints - 1] = e->GetGeom()->GetVid(2);
565
566 // Edge numbering
567 StdRegions::Orientation edgeOrient;
568 int edge1 = 0;
569 int edge2 = 0;
570 for (int n = 1; n < np - 1; n++)
571 {
572 // Edge 0
573 edgeOrient = e->GetGeom()->GetEorient(0);
574 if (edgeOrient == StdRegions::eForwards)
575 {
576 vId[cnt1 + n] =
577 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
578 }
579 else
580 {
581 vId[cnt1 + np - 1 - n] =
582 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
583 }
584
585 // Edge 1
586 edgeOrient = e->GetGeom()->GetEorient(1);
587 if (edgeOrient == StdRegions::eForwards)
588 {
589 edge1 += np - n;
590 vId[cnt1 + np - 1 + edge1] =
591 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
592 }
593 else
594 {
595 edge1 += n;
596 vId[cnt1 + newpoints - 1 - edge1] =
597 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
598 }
599
600 // Edge 2
601 edgeOrient = e->GetGeom()->GetEorient(2);
602 if (edgeOrient == StdRegions::eForwards)
603 {
604 edge2 += n + 1;
605 vId[cnt1 + newpoints - 1 - edge2] =
606 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
607 }
608 else
609 {
610 edge2 += np + 1 - n;
611 vId[cnt1 + edge2] =
612 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
613 }
614 }
615
616 // Interior numbering
617 edge2 = 0;
618 for (int n = 1; n < np - 1; n++)
619 {
620 edge2 += np + 1 - n;
621 for (int m = 1; m < np - n - 1; m++)
622 {
623 vId[cnt1 + edge2 + m] =
624 4 * nel + maxN * 4 * nel + cnt2;
625 cnt2++;
626 }
627 }
628 cnt1 += newpoints;
629 }
630 break;
632 {
633 int newpoints =
635 np);
636 // Vertex numbering
637 vId[cnt1] = e->GetGeom()->GetVid(0);
638 vId[cnt1 + np - 1] = e->GetGeom()->GetVid(1);
639 vId[cnt1 + np * np - 1] = e->GetGeom()->GetVid(2);
640 vId[cnt1 + np * (np - 1)] = e->GetGeom()->GetVid(3);
641
642 // Edge numbering
643 StdRegions::Orientation edgeOrient;
644 for (int n = 1; n < np - 1; n++)
645 {
646 // Edge 0
647 edgeOrient = e->GetGeom()->GetEorient(0);
648 if (edgeOrient == StdRegions::eForwards)
649 {
650 vId[cnt1 + n] =
651 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
652 }
653 else
654 {
655 vId[cnt1 + np - 1 - n] =
656 4 * nel + maxN * e->GetGeom()->GetEid(0) + n;
657 }
658
659 // Edge 1
660 edgeOrient = e->GetGeom()->GetEorient(1);
661 if (edgeOrient == StdRegions::eForwards)
662 {
663 vId[cnt1 + np - 1 + n * np] =
664 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
665 }
666 else
667 {
668 vId[cnt1 + np * np - 1 - n * np] =
669 4 * nel + maxN * e->GetGeom()->GetEid(1) + n;
670 }
671
672 // Edge 2
673 edgeOrient = e->GetGeom()->GetEorient(2);
674 if (edgeOrient == StdRegions::eForwards)
675 {
676 vId[cnt1 + np * np - 1 - n] =
677 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
678 }
679 else
680 {
681 vId[cnt1 + np * (np - 1) + n] =
682 4 * nel + maxN * e->GetGeom()->GetEid(2) + n;
683 }
684
685 // Edge 3
686 edgeOrient = e->GetGeom()->GetEorient(3);
687 if (edgeOrient == StdRegions::eForwards)
688 {
689 vId[cnt1 + np * (np - 1) - n * np] =
690 4 * nel + maxN * e->GetGeom()->GetEid(3) + n;
691 }
692 else
693 {
694 vId[cnt1 + n * np] =
695 4 * nel + maxN * e->GetGeom()->GetEid(3) + n;
696 }
697 }
698
699 // Interior numbering
700 for (int n = 1; n < np - 1; n++)
701 {
702 for (int m = 1; m < np - 1; m++)
703 {
704 vId[cnt1 + m * np + n] =
705 4 * nel + maxN * 4 * nel + cnt2;
706 cnt2++;
707 }
708 }
709 cnt1 += newpoints;
710 }
711 break;
712 default:
713 {
714 ASSERTL0(false, "Points not known");
715 }
716 }
717 }
718
719 // Crete new connectivity using homogeneous information
720 for (int i = 0; i < nel; ++i)
721 {
722 int nTris = oldConn[i].size() / 3;
723 // Create array for new connectivity
724 // (3 tetrahedra between each plane for each triangle)
725 conn = Array<OneD, int>(4 * 3 * nTris * (nPlanes - 1));
726 cnt2 = 0;
727 for (int n = 0; n < nTris; n++)
728 {
729 // Get id of vertices in this triangle (on plane 0)
730 int vId0 = vId[oldConn[i][3 * n + 0]];
731 int vId1 = vId[oldConn[i][3 * n + 1]];
732 int vId2 = vId[oldConn[i][3 * n + 2]];
733
734 // Determine ordering based on global Id of pts
735 Array<OneD, int> rot(3);
736 if ((vId0 < vId1) && (vId0 < vId2))
737 {
738 rot[0] = 0;
739 if (vId1 < vId2)
740 {
741 rot[1] = 1;
742 rot[2] = 2;
743 }
744 else
745 {
746 rot[1] = 2;
747 rot[2] = 1;
748 }
749 }
750 else if ((vId1 < vId0) && (vId1 < vId2))
751 {
752 rot[0] = 1;
753 if (vId0 < vId2)
754 {
755 rot[1] = 0;
756 rot[2] = 2;
757 }
758 else
759 {
760 rot[1] = 2;
761 rot[2] = 0;
762 }
763 }
764 else
765 {
766 rot[0] = 2;
767 if (vId0 < vId1)
768 {
769 rot[1] = 0;
770 rot[2] = 1;
771 }
772 else
773 {
774 rot[1] = 1;
775 rot[2] = 0;
776 }
777 }
778
779 // Define new connectivity with tetrahedra
780 for (int p = 0; p < nPlanes - 1; p++)
781 {
782 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[0]];
783 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[1]];
784 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[2]];
785 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[2]];
786
787 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[0]];
788 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[1]];
789 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[2]];
790 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[1]];
791
792 conn[cnt2++] = oldConn[i + p * nel][3 * n + rot[0]];
793 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[1]];
794 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[2]];
795 conn[cnt2++] = oldConn[i + (p + 1) * nel][3 * n + rot[0]];
796 }
797 }
798 newConn.push_back(conn);
799 }
800 }
801 else
802 {
803 ASSERTL0(false, "Points type not supported.");
804 }
805
806 m_f->m_fieldPts->SetConnectivity(newConn);
807}
808
810 int n, const Array<OneD, const NekDouble> &phys,
812{
814 e = m_f->m_exp[0]->GetExp(n);
815
816 switch (e->DetShapeType())
817 {
819 {
820 int np0 = e->GetBasis(0)->GetNumPoints();
821 int np1 = e->GetBasis(1)->GetNumPoints();
822 int np = max(np0, np1);
823
824 // to ensure points are correctly projected to np need
825 // to increase the order slightly of coordinates
826 LibUtilities::PointsKey pa(np + 1, e->GetPointsType(0));
827 LibUtilities::PointsKey pb(np, e->GetPointsType(1));
828 Array<OneD, NekDouble> tophys(np * (np + 1));
829
832 StdRegions::StdTriExp OrthoExp(Ba, Bb);
833
834 // interpolate points to new phys points!
835 LibUtilities::Interp2D(e->GetBasis(0)->GetBasisKey(),
836 e->GetBasis(1)->GetBasisKey(), phys, Ba, Bb,
837 tophys);
838
839 OrthoExp.FwdTrans(tophys, coeffs);
840 break;
841 }
843 {
844 int np0 = e->GetBasis(0)->GetNumPoints();
845 int np1 = e->GetBasis(1)->GetNumPoints();
846 int np = max(np0, np1);
847
848 LibUtilities::PointsKey pa(np + 1, e->GetPointsType(0));
849 LibUtilities::PointsKey pb(np + 1, e->GetPointsType(1));
850 Array<OneD, NekDouble> tophys((np + 1) * (np + 1));
851
854 StdRegions::StdQuadExp OrthoExp(Ba, Bb);
855
856 // interpolate points to new phys points!
857 LibUtilities::Interp2D(e->GetBasis(0)->GetBasisKey(),
858 e->GetBasis(1)->GetBasisKey(), phys, Ba, Bb,
859 tophys);
860
861 OrthoExp.FwdTrans(phys, coeffs);
862 break;
863 }
864 default:
865 ASSERTL0(false, "Shape needs setting up");
866 break;
867 }
868}
869} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
FieldSharedPtr m_f
Field object.
Definition Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition Module.h:272
void v_Process(po::variables_map &vm) override
Write mesh to output file.
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
void GenOrthoModes(int n, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &coeffs)
Abstract base class for processing modules.
Definition Module.h:301
Describes the specification for a Basis.
Definition Basis.h:45
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Defines a specification for a set of points.
Definition Points.h:50
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
std::shared_ptr< Field > FieldSharedPtr
Definition Field.hpp:1026
std::pair< ModuleType, std::string > ModuleKey
Definition Module.h:180
ModuleFactory & GetModuleFactory()
Definition Module.cpp:47
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb)
int getNumberOfCoefficients(int Na, int Nb, int Nc)
int getNumberOfCoefficients(int Na, int Nb)
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis,...
Definition Interp.cpp:101
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition PointsType.h:74
@ eOrtho_A
Principle Orthogonal Functions .
Definition BasisType.h:42
@ eOrtho_B
Principle Orthogonal Functions .
Definition BasisType.h:44
@ eFourier
Fourier Expansion .
Definition BasisType.h:55
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition Expansion.h:66
static Array< OneD, NekDouble > NullNekDouble1DArray
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition Vmath.hpp:194
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825
STL namespace.
scalarT< T > max(scalarT< T > lhs, scalarT< T > rhs)
Definition scalar.hpp:305
Represents a command-line configuration option.
Definition Module.h:129