Nektar++
DisContField3DHomogeneous1D.cpp
Go to the documentation of this file.
1//////////////////////////////////////////////////////////////////////////////
2//
3// File: DisContField3DHomogeneous1D.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: Field definition for 3D domain with boundary
32// conditions using LDG flux and a 1D homogeneous direction
33//
34///////////////////////////////////////////////////////////////////////////////
35
40
42{
43
45 : ExpList3DHomogeneous1D(), m_bndCondExpansions(), m_bndCondBndWeight(),
46 m_bndConditions()
47{
48}
49
52 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
53 const bool useFFT, const bool dealiasing)
54 : ExpList3DHomogeneous1D(pSession, HomoBasis, lhom, useFFT, dealiasing),
55 m_bndCondExpansions(), m_bndCondBndWeight(), m_bndConditions()
56{
57}
58
60 const DisContField3DHomogeneous1D &In, const bool DeclarePlanesSetCoeffPhys)
61 : ExpList3DHomogeneous1D(In, false),
62 m_bndCondExpansions(In.m_bndCondExpansions),
63 m_bndCondBndWeight(In.m_bndCondBndWeight),
64 m_bndConditions(In.m_bndConditions)
65{
66 if (DeclarePlanesSetCoeffPhys)
67 {
68 DisContFieldSharedPtr zero_plane =
69 std::dynamic_pointer_cast<DisContField>(In.m_planes[0]);
70
71 for (int n = 0; n < m_planes.size(); ++n)
72 {
74 *zero_plane, false);
75 }
76
78 }
79}
80
83 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
84 const bool useFFT, const bool dealiasing,
86 const std::string &variable, const Collections::ImplementationType ImpType)
87 : ExpList3DHomogeneous1D(pSession, HomoBasis, lhom, useFFT, dealiasing),
88 m_bndCondExpansions(), m_bndCondBndWeight(), m_bndConditions()
89{
90 int i, n, nel;
91 DisContFieldSharedPtr plane_zero;
93
94 // note that nzplanes can be larger than nzmodes
96 pSession, graph2D, variable, true, false, ImpType);
97
99
100 nel = m_planes[0]->GetExpSize();
101
102 for (i = 0; i < nel; ++i)
103 {
104 (*m_exp).push_back(m_planes[0]->GetExp(i));
105 }
106
107 for (n = 1; n < m_planes.size(); ++n)
108 {
110 *plane_zero, graph2D, variable, true, false);
111
112 for (i = 0; i < nel; ++i)
113 {
114 (*m_exp).push_back((*m_exp)[i]);
115 }
116 }
117
118 // Set up trace object.
120 for (n = 0; n < m_planes.size(); ++n)
121 {
122 trace[n] = m_planes[n]->GetTrace();
123 }
124
126 pSession, HomoBasis, lhom, useFFT, dealiasing, trace);
127
128 // Setup default optimisation information
129 nel = GetExpSize();
130
131 SetCoeffPhys();
132
133 // Do not set up BCs if default variable
134 if (variable.compare("DefaultVar") != 0)
135 {
136 SetupBoundaryConditions(HomoBasis, lhom, bcs, variable);
137 }
138
139 SetUpDG();
140}
141
142/**
143 * @brief Default destructor
144 */
146{
147}
148
150 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
151 SpatialDomains::BoundaryConditions &bcs, const std::string variable)
152{
153 int n, cnt = 0;
154
155 // Setup an ExpList2DHomogeneous1D expansion for boundary
156 // conditions and link to class declared in m_planes
158 bcs.GetBoundaryRegions();
161
162 int nplanes = m_planes.size();
163
166 m_bndConditions = m_planes[0]->UpdateBndConditions();
167
168 m_bndCondBndWeight = Array<OneD, NekDouble>{bregions.size(), 0.0};
169
170 Array<OneD, MultiRegions::ExpListSharedPtr> PlanesBndCondExp(nplanes);
171
172 for (auto &it : bregions)
173 {
175 GetBoundaryCondition(bconditions, it.first, variable);
176
177 for (n = 0; n < nplanes; ++n)
178 {
179 PlanesBndCondExp[n] = m_planes[n]->UpdateBndCondExpansion(cnt);
180 }
181
182 // Initialise comm for the boundary regions
183 auto comm = boundaryCondition->GetComm();
184 int size = boundaryCondition->GetComm()->GetSize();
185
186 if (size > 1)
187 {
188 // It seems to work either way
189 // comm->SplitComm(1,size);
190 comm->SplitComm(m_StripZcomm->GetSize(),
191 size / m_StripZcomm->GetSize());
192 }
193
194 m_bndCondExpansions[cnt++] =
196 AllocateSharedPtr(m_session, HomoBasis, lhom, m_useFFT, false,
197 PlanesBndCondExp, comm);
198 }
199 v_EvaluateBoundaryConditions(0.0, variable);
200}
201
202/**
203 * @brief Set up all DG member variables and maps.
204 */
206{
207 const int nPlanes = m_planes.size();
208 const int nTracePlane = m_planes[0]->GetTrace()->GetExpSize();
209
210 // Get trace map from first plane.
211 AssemblyMapDGSharedPtr traceMap = m_planes[0]->GetTraceMap();
212 const Array<OneD, const int> &traceBndMap =
213 traceMap->GetBndCondIDToGlobalTraceID();
214 int mapSize = traceBndMap.size();
215
216 // Set up trace boundary map
217 m_traceBndMap = Array<OneD, int>(nPlanes * mapSize);
218
219 int i, n, e, cnt = 0, cnt1 = 0;
220
221 for (i = 0; i < m_bndCondExpansions.size(); ++i)
222 {
223 int nExp = m_bndCondExpansions[i]->GetExpSize();
224 int nPlaneExp = nExp / nPlanes;
225
226 if (m_bndConditions[i]->GetBoundaryConditionType() ==
228 {
229 continue;
230 }
231
232 for (n = 0; n < nPlanes; ++n)
233 {
234 const int offset = n * nTracePlane;
235 for (e = 0; e < nPlaneExp; ++e)
236 {
237 m_traceBndMap[cnt++] = offset + traceBndMap[cnt1 + e];
238 }
239 }
240
241 cnt1 += nPlaneExp;
242 }
243}
244
246 const Array<OneD, const NekDouble> &inarray,
248 const StdRegions::VarCoeffMap &varcoeff,
249 const MultiRegions::VarFactorsMap &varfactors,
250 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
251{
252 int n;
253 int cnt = 0;
254 int cnt1 = 0;
256 StdRegions::ConstFactorMap new_factors;
257
258 int npts_fce = PhysSpaceForcing ? m_npoints : m_ncoeffs;
260 Array<OneD, NekDouble> fce(npts_fce);
262
263 // Transform forcing function in half-physical space if required
264 if (m_WaveSpace)
265 {
266 fce = inarray;
267 }
268 else
269 {
270 HomogeneousFwdTrans(npts_fce, inarray, fce);
271 }
272
273 for (n = 0; n < m_planes.size(); ++n)
274 {
275 if (n != 1 || m_transposition->GetK(n) != 0)
276 {
277 beta = 2 * M_PI * (m_transposition->GetK(n)) / m_lhom;
278 new_factors = factors;
279 // add in Homogeneous Fourier direction and SVV if turned on
280 new_factors[StdRegions::eFactorLambda] +=
281 beta * beta * (1 + GetSpecVanVisc(n));
282
283 wfce = (PhysSpaceForcing) ? fce + cnt : fce + cnt1;
284 m_planes[n]->HelmSolve(wfce, e_out = outarray + cnt1, new_factors,
285 varcoeff, varfactors, dirForcing,
286 PhysSpaceForcing);
287 }
288
289 cnt += m_planes[n]->GetTotPoints();
290 cnt1 += m_planes[n]->GetNcoeffs();
291 }
292 return NullGlobalLinSysKey;
293}
294
295/// @todo Fix in another way considering all the planes
297{
298 return m_trace;
299}
300
301/// @todo Fix in another way considering all the planes
303{
304 return m_planes[0]->GetTraceMap();
305}
306
308 Array<OneD, NekDouble> &outarray)
309{
310 ASSERTL1(m_physState == true,
311 "Field must be in physical state to extract trace space.");
312
313 v_ExtractTracePhys(m_phys, outarray);
314}
315
316/**
317 * @brief This method extracts the trace (edges in 2D) for each plane
318 * from the field @a inarray and puts the values in @a outarray.
319 *
320 * It assumes the field is C0 continuous so that it can overwrite the
321 * edge data when visited by the two adjacent elements.
322 *
323 * @param inarray An array containing the 2D data from which we wish
324 * to extract the edge data.
325 * @param outarray The resulting edge information.
326 */
328 const Array<OneD, const NekDouble> &inarray,
329 Array<OneD, NekDouble> &outarray)
330{
331 int nPoints_plane = m_planes[0]->GetTotPoints();
332 int nTracePts = m_planes[0]->GetTrace()->GetTotPoints();
333
334 for (int i = 0; i < m_planes.size(); ++i)
335 {
336 Array<OneD, NekDouble> inarray_plane(nPoints_plane, 0.0);
337 Array<OneD, NekDouble> outarray_plane(nPoints_plane, 0.0);
338
339 Vmath::Vcopy(nPoints_plane, &inarray[i * nPoints_plane], 1,
340 &inarray_plane[0], 1);
341
342 m_planes[i]->ExtractTracePhys(inarray_plane, outarray_plane);
343
344 Vmath::Vcopy(nTracePts, &outarray_plane[0], 1, &outarray[i * nTracePts],
345 1);
346 }
347}
348
350{
351 return m_traceBndMap;
352}
353
355 int i, std::shared_ptr<ExpList> &result, const bool DeclareCoeffPhysArrays)
356{
357 int n, cnt, nq;
358 int offsetOld, offsetNew;
359
360 std::vector<unsigned int> eIDs;
361 Array<OneD, int> ElmtID, EdgeID;
362 GetBoundaryToElmtMap(ElmtID, EdgeID);
363
364 // Skip other boundary regions
365 for (cnt = n = 0; n < i; ++n)
366 {
367 cnt += m_bndCondExpansions[n]->GetExpSize();
368 }
369
370 // Populate eIDs with information from BoundaryToElmtMap
371 for (n = 0; n < m_bndCondExpansions[i]->GetExpSize(); ++n)
372 {
373 eIDs.push_back(ElmtID[cnt + n]);
374 }
375
376 // Create expansion list
377 // Note: third arguemnt declares phys coeffs that are not
378 // required but currently it is needed to declare the
379 // planes because bool is liked
381 *this, eIDs, true, Collections::eNoCollection);
382
383 // Copy phys and coeffs to new explist
384 if (DeclareCoeffPhysArrays)
385 {
386 Array<OneD, NekDouble> tmp1, tmp2;
387 for (n = 0; n < result->GetExpSize(); ++n)
388 {
389 nq = GetExp(ElmtID[cnt + n])->GetTotPoints();
390 offsetOld = GetPhys_Offset(ElmtID[cnt + n]);
391 offsetNew = result->GetPhys_Offset(n);
392 Vmath::Vcopy(nq, tmp1 = GetPhys() + offsetOld, 1,
393 tmp2 = result->UpdatePhys() + offsetNew, 1);
394
395 nq = GetExp(ElmtID[cnt + n])->GetNcoeffs();
396 offsetOld = GetCoeff_Offset(ElmtID[cnt + n]);
397 offsetNew = result->GetCoeff_Offset(n);
398 Vmath::Vcopy(nq, tmp1 = GetCoeffs() + offsetOld, 1,
399 tmp2 = result->UpdateCoeffs() + offsetNew, 1);
400 }
401 }
402
403 // Set wavespace value
404 result->SetWaveSpace(GetWaveSpace());
405}
406
408 Array<OneD, int> &ElmtID, Array<OneD, int> &EdgeID)
409{
410
411 if (m_BCtoElmMap.size() == 0)
412 {
413 Array<OneD, int> ElmtID_tmp;
414 Array<OneD, int> EdgeID_tmp;
415
416 m_planes[0]->GetBoundaryToElmtMap(ElmtID_tmp, EdgeID_tmp);
417 int nel_per_plane = m_planes[0]->GetExpSize();
418 int nplanes = m_planes.size();
419
420 int MapSize = ElmtID_tmp.size();
421
422 m_BCtoElmMap = Array<OneD, int>(nplanes * MapSize);
423 m_BCtoEdgMap = Array<OneD, int>(nplanes * MapSize);
424
425 // If this mesh (or partition) has no BCs, skip this step
426 if (MapSize > 0)
427 {
428 int i, j, n, cnt;
429 int cntPlane = 0;
430 for (cnt = n = 0; n < m_bndCondExpansions.size(); ++n)
431 {
432 int planeExpSize =
433 m_planes[0]->GetBndCondExpansions()[n]->GetExpSize();
434 for (i = 0; i < planeExpSize; ++i, ++cntPlane)
435 {
436 for (j = 0; j < nplanes; j++)
437 {
438 m_BCtoElmMap[cnt + i + j * planeExpSize] =
439 ElmtID_tmp[cntPlane] + j * nel_per_plane;
440 m_BCtoEdgMap[cnt + i + j * planeExpSize] =
441 EdgeID_tmp[cntPlane];
442 }
443 }
444 cnt += m_bndCondExpansions[n]->GetExpSize();
445 }
446 }
447 }
448 ElmtID = m_BCtoElmMap;
449 EdgeID = m_BCtoEdgMap;
450}
451
453 Array<OneD, NekDouble> &BndVals, const Array<OneD, NekDouble> &TotField,
454 int BndID)
455{
458
461
462 int cnt = 0;
463 int pos = 0;
464 int exp_size, exp_size_per_plane, elmtID, boundaryID;
465 int offset, exp_dim;
466
467 for (int k = 0; k < m_planes.size(); k++)
468 {
469 for (int n = 0; n < m_bndConditions.size(); ++n)
470 {
471 exp_size = m_bndCondExpansions[n]->GetExpSize();
472 exp_size_per_plane = exp_size / m_planes.size();
473
474 for (int i = 0; i < exp_size_per_plane; i++)
475 {
476 if (n == BndID)
477 {
478 elmtID = m_BCtoElmMap[cnt];
479 boundaryID = m_BCtoEdgMap[cnt];
480 exp_dim = m_bndCondExpansions[n]
481 ->GetExp(i + k * exp_size_per_plane)
482 ->GetTotPoints();
483 offset = GetPhys_Offset(elmtID);
484 elmt = GetExp(elmtID);
485 temp_BC_exp =
486 std::dynamic_pointer_cast<LocalRegions::Expansion1D>(
488 i + k * exp_size_per_plane));
489
490 elmt->GetTracePhysVals(boundaryID, temp_BC_exp,
491 tmp_Tot = TotField + offset,
492 tmp_BC = BndVals + pos);
493 pos += exp_dim;
494 }
495 cnt++;
496 }
497 }
498 }
499}
500
503 Array<OneD, NekDouble> &outarray, int BndID)
504{
507
510 Array<OneD, NekDouble> tmp_outarray;
511
512 int cnt = 0;
513 int exp_size, exp_size_per_plane, elmtID, Phys_offset, Coef_offset;
514
515 for (int k = 0; k < m_planes.size(); k++)
516 {
517 for (int n = 0; n < m_bndConditions.size(); ++n)
518 {
519 exp_size = m_bndCondExpansions[n]->GetExpSize();
520 exp_size_per_plane = exp_size / m_planes.size();
521
522 for (int i = 0; i < exp_size_per_plane; i++)
523 {
524 if (n == BndID)
525 {
526 elmtID = m_BCtoElmMap[cnt];
527
528 Phys_offset = m_bndCondExpansions[n]->GetPhys_Offset(
529 i + k * exp_size_per_plane);
530 Coef_offset = m_bndCondExpansions[n]->GetCoeff_Offset(
531 i + k * exp_size_per_plane);
532
533 elmt = GetExp(elmtID);
534 temp_BC_exp =
535 std::dynamic_pointer_cast<LocalRegions::Expansion1D>(
537 i + k * exp_size_per_plane));
538
539 temp_BC_exp->NormVectorIProductWRTBase(
540 tmp_V1 = V1 + Phys_offset, tmp_V2 = V2 + Phys_offset,
541 tmp_outarray = outarray + Coef_offset);
542 }
543 cnt++;
544 }
545 }
546 }
547}
548
549/**
550 */
552 int i, Array<OneD, Array<OneD, NekDouble>> &normals)
553{
554 int expdim = GetCoordim(0);
555 int coordim = 3;
558
559 Array<OneD, int> ElmtID, EdgeID;
560 GetBoundaryToElmtMap(ElmtID, EdgeID);
561
562 // Initialise result
563 normals = Array<OneD, Array<OneD, NekDouble>>(coordim);
564 for (int j = 0; j < coordim; ++j)
565 {
566 normals[j] = Array<OneD, NekDouble>(
568 }
569
570 // Skip other boundary regions
571 int cnt = 0;
572 for (int n = 0; n < i; ++n)
573 {
574 cnt += GetBndCondExpansions()[n]->GetExpSize();
575 }
576
577 int offset;
578 for (int n = 0; n < GetBndCondExpansions()[i]->GetExpSize(); ++n)
579 {
580 offset = GetBndCondExpansions()[i]->GetPhys_Offset(n);
581 int nq = GetBndCondExpansions()[i]->GetExp(n)->GetTotPoints();
582
583 elmt = GetExp(ElmtID[cnt + n]);
585 elmt->GetTraceNormal(EdgeID[cnt + n]);
586 // Copy to result
587 for (int j = 0; j < expdim; ++j)
588 {
589 Vmath::Vcopy(nq, normalsElmt[j], 1, tmp = normals[j] + offset, 1);
590 }
591 }
592}
593
594std::map<int, RobinBCInfoSharedPtr> DisContField3DHomogeneous1D::
596{
597 return std::map<int, RobinBCInfoSharedPtr>();
598}
599
601 const NekDouble time, const std::string varName,
602 [[maybe_unused]] const NekDouble x2_in,
603 [[maybe_unused]] const NekDouble x3_in)
604{
605 int i;
606 int npoints;
607 int nbnd = m_bndCondExpansions.size();
609
610 for (i = 0; i < nbnd; ++i)
611 {
612 if (time == 0.0 || m_bndConditions[i]->IsTimeDependent())
613 {
614 locExpList = m_bndCondExpansions[i];
615 npoints = locExpList->GetNpoints();
616
617 Array<OneD, NekDouble> x0(npoints, 0.0);
618 Array<OneD, NekDouble> x1(npoints, 0.0);
619 Array<OneD, NekDouble> x2(npoints, 0.0);
620 Array<OneD, NekDouble> valuesFile(npoints, 1.0),
621 valuesExp(npoints, 1.0);
622
623 locExpList->GetCoords(x0, x1, x2);
624
625 if (m_bndConditions[i]->GetBoundaryConditionType() ==
627 {
629 std::static_pointer_cast<
631 m_bndConditions[i]);
632 std::string bcfilename = bcPtr->m_filename;
633 std::string exprbcs = bcPtr->m_expr;
634
635 if (bcfilename != "")
636 {
637 locExpList->ExtractCoeffsFromFile(
638 bcfilename, bcPtr->GetComm(), varName,
639 locExpList->UpdateCoeffs());
640 locExpList->BwdTrans(locExpList->GetCoeffs(),
641 locExpList->UpdatePhys());
642 valuesFile = locExpList->GetPhys();
643 }
644
645 if (exprbcs != "")
646 {
647 LibUtilities::Equation condition =
648 std::static_pointer_cast<
651 ->m_dirichletCondition;
652
653 condition.Evaluate(x0, x1, x2, time, valuesExp);
654 }
655
656 Vmath::Vmul(npoints, valuesExp, 1, valuesFile, 1,
657 locExpList->UpdatePhys(), 1);
658 // set wave space to false since have set up phys values
659 locExpList->SetWaveSpace(false);
660 locExpList->FwdTransBndConstrained(locExpList->GetPhys(),
661 locExpList->UpdateCoeffs());
662 }
663 else if (m_bndConditions[i]->GetBoundaryConditionType() ==
665 {
666 SpatialDomains::NeumannBCShPtr bcPtr = std::static_pointer_cast<
668 m_bndConditions[i]);
669
670 std::string bcfilename = bcPtr->m_filename;
671
672 if (bcfilename != "")
673 {
674 locExpList->ExtractCoeffsFromFile(
675 bcfilename, bcPtr->GetComm(), varName,
676 locExpList->UpdateCoeffs());
677 locExpList->BwdTrans(locExpList->GetCoeffs(),
678 locExpList->UpdatePhys());
679 }
680 else
681 {
682 LibUtilities::Equation condition =
683 std::static_pointer_cast<
686 ->m_neumannCondition;
687
688 condition.Evaluate(x0, x1, x2, time,
689 locExpList->UpdatePhys());
690 }
691
692 locExpList->IProductWRTBase(locExpList->GetPhys(),
693 locExpList->UpdateCoeffs());
694 }
695 else if (m_bndConditions[i]->GetBoundaryConditionType() ==
697 {
698 SpatialDomains::RobinBCShPtr bcPtr = std::static_pointer_cast<
700
701 std::string bcfilename = bcPtr->m_filename;
702
703 if (bcfilename != "")
704 {
705 locExpList->ExtractCoeffsFromFile(
706 bcfilename, bcPtr->GetComm(), varName,
707 locExpList->UpdateCoeffs());
708 locExpList->BwdTrans(locExpList->GetCoeffs(),
709 locExpList->UpdatePhys());
710 }
711 else
712 {
713 LibUtilities::Equation condition =
714 std::static_pointer_cast<
717 ->m_robinFunction;
718
719 condition.Evaluate(x0, x1, x2, time,
720 locExpList->UpdatePhys());
721 }
722
723 locExpList->IProductWRTBase(locExpList->GetPhys(),
724 locExpList->UpdateCoeffs());
725 }
726 else if (m_bndConditions[i]->GetBoundaryConditionType() ==
728 {
729 continue;
730 }
731 else
732 {
733 ASSERTL0(false, "This type of BC not implemented yet");
734 }
735 }
736 }
737}
738
741{
742 return m_bndCondExpansions;
743}
744
747{
748 return m_bndConditions;
749}
750
752 int i)
753{
754 return m_bndCondExpansions[i];
755}
756
759{
760 return m_bndConditions;
761}
762
764 const NekDouble value)
765{
766 m_bndCondBndWeight[index] = value;
767}
768} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Describes the specification for a Basis.
Definition: Basis.h:45
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID) override
GlobalLinSysKey v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing) override
void SetUpDG()
Set up all DG member variables and maps.
void v_ExtractTracePhys(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
This method extracts the trace (edges in 2D) for each plane from the field inarray and puts the value...
void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays) override
std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i) override
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
Array< OneD, int > m_BCtoElmMap
Storage space for the boundary to element and boundary to trace map. This member variable is really a...
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, SpatialDomains::BoundaryConditions &bcs, const std::string variable)
void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID) override
void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID) override
void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals) override
void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble) override
This function evaluates the boundary conditions at a certain time-level.
const Array< OneD, const MultiRegions::ExpListSharedPtr > & v_GetBndCondExpansions(void) override
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions() override
const Array< OneD, const int > & v_GetTraceBndMap() override
std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo() override
void v_SetBndCondBwdWeight(const int index, const NekDouble value) override
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions() override
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
Array< OneD, ExpListSharedPtr > m_planes
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:2121
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2266
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:5429
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1952
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2038
int GetCoeff_Offset(int n) const
Get the start offset position for a local contiguous list of coeffs correspoinding to element n.
Definition: ExpList.h:2087
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2046
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1104
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1115
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1060
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2079
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1055
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1604
int GetPhys_Offset(int n) const
Get the start offset position for a local contiguous list of quadrature points in a full array corres...
Definition: ExpList.h:2094
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1096
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1848
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1557
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1911
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:234
const BoundaryConditionCollection & GetBoundaryConditions(void) const
Definition: Conditions.h:244
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
Definition: Expansion1D.h:50
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:46
std::shared_ptr< DisContField > DisContFieldSharedPtr
Definition: DisContField.h:345
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static GlobalLinSysKey NullGlobalLinSysKey(StdRegions::eNoMatrixType)
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::shared_ptr< DirichletBoundaryCondition > DirichletBCShPtr
Definition: Conditions.h:213
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:219
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< NeumannBoundaryCondition > NeumannBCShPtr
Definition: Conditions.h:214
std::shared_ptr< RobinBoundaryCondition > RobinBCShPtr
Definition: Conditions.h:215
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:346
StdRegions::ConstFactorMap factors
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825