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