36 #include <boost/core/ignore_unused.hpp>
46 namespace MultiRegions
51 m_bndCondExpansions(),
62 const bool dealiasing):
64 m_bndCondExpansions(),
72 const bool DeclarePlanesSetCoeffPhys)
74 m_bndCondExpansions (In.m_bndCondExpansions),
75 m_bndCondBndWeight (In.m_bndCondBndWeight),
76 m_bndConditions (In.m_bndConditions)
78 if (DeclarePlanesSetCoeffPhys)
81 std::dynamic_pointer_cast<DisContField> (In.
m_planes[0]);
83 for(
int n = 0; n <
m_planes.size(); ++n)
99 const bool dealiasing,
101 const std::string &variable,
105 m_bndCondExpansions(),
106 m_bndCondBndWeight(),
123 for (i = 0; i < nel; ++i)
128 for (n = 1; n <
m_planes.size(); ++n)
132 variable,
true,
false);
133 for(i = 0; i < nel; ++i)
135 (*m_exp).push_back((*
m_exp)[i]);
141 for (n = 0; n <
m_planes.size(); ++n)
147 pSession, HomoBasis, lhom, useFFT,
156 if(variable.compare(
"DefaultVar") != 0)
175 const std::string variable)
194 PlanesBndCondExp(nplanes);
196 for (
auto &it : bregions)
200 for (n = 0; n < nplanes; ++n)
207 auto comm = boundaryCondition->GetComm();
208 int size = boundaryCondition->GetComm()->GetSize();
222 PlanesBndCondExp, comm);
234 const bool PhysSpaceForcing)
256 for (n = 0; n <
m_planes.size(); ++n)
261 new_factors = factors;
266 wfce = (PhysSpaceForcing)? fce+cnt:fce+cnt1;
269 e_out = outarray + cnt1,
270 new_factors, varcoeff, varfactors,
282 const std::string varName,
286 boost::ignore_unused(x2_in,x3_in);
292 for (i = 0; i < nbnd; ++i)
297 npoints = locExpList->GetNpoints();
304 locExpList->GetCoords(x0, x1, x2);
313 std::string exprbcs = bcPtr->m_expr;
318 valuesFile = locExpList->GetPhys();
324 std::static_pointer_cast<SpatialDomains::
325 DirichletBoundaryCondition >(
328 condition.
Evaluate(x0, x1, x2, time, valuesExp);
332 locExpList->UpdatePhys(), 1);
335 locExpList->SetWaveSpace(
false);
337 locExpList->FwdTrans_BndConstrained(
338 locExpList->GetPhys(),
339 locExpList->UpdateCoeffs());
357 static_pointer_cast<SpatialDomains::
358 NeumannBoundaryCondition>(
361 condition.
Evaluate(x0, x1, x2, time,
362 locExpList->UpdatePhys());
364 locExpList->IProductWRTBase(locExpList->GetPhys(),
365 locExpList->UpdateCoeffs());
383 static_pointer_cast<SpatialDomains::
384 RobinBoundaryCondition>(
387 condition.
Evaluate(x0, x1, x2, time,
388 locExpList->UpdatePhys());
392 locExpList->IProductWRTBase(locExpList->GetPhys(),
393 locExpList->UpdateCoeffs());
403 ASSERTL0(
false,
"This type of BC not implemented yet");
431 m_planes[0]->GetBoundaryToElmtMap(ElmtID_tmp, EdgeID_tmp);
432 int nel_per_plane =
m_planes[0]->GetExpSize();
435 int MapSize = ElmtID_tmp.size();
448 ->GetBndCondExpansions()[n]
450 for (i = 0; i < planeExpSize ; ++i, ++cntPlane)
452 for(j = 0; j < nplanes; j++)
455 ElmtID_tmp[cntPlane]+j*nel_per_plane;
457 EdgeID_tmp[cntPlane];
469 std::shared_ptr<ExpList> &result,
470 const bool DeclareCoeffPhysArrays)
473 int offsetOld, offsetNew;
475 std::vector<unsigned int> eIDs;
480 for (cnt = n = 0; n < i; ++n)
488 eIDs.push_back(ElmtID[cnt+n]);
497 if ( DeclareCoeffPhysArrays)
500 for (n = 0; n < result->GetExpSize(); ++n)
502 nq =
GetExp(ElmtID[cnt+n])->GetTotPoints();
504 offsetNew = result->GetPhys_Offset(n);
506 tmp2 = result->UpdatePhys()+ offsetNew, 1);
508 nq =
GetExp(ElmtID[cnt+n])->GetNcoeffs();
510 offsetNew = result->GetCoeff_Offset(n);
512 tmp2 = result->UpdateCoeffs()+ offsetNew, 1);
532 int exp_size, exp_size_per_plane, elmtID, boundaryID;
535 for (
int k = 0; k <
m_planes.size(); k++)
540 exp_size_per_plane = exp_size/
m_planes.size();
542 for (
int i = 0; i < exp_size_per_plane; i++)
549 GetExp(i+k*exp_size_per_plane)->GetTotPoints();
552 temp_BC_exp = std::dynamic_pointer_cast<
555 i + k * exp_size_per_plane )
559 tmp_Tot = TotField + offset,
560 tmp_BC = BndVals + pos);
583 int exp_size, exp_size_per_plane, elmtID, Phys_offset, Coef_offset;
585 for(
int k = 0; k <
m_planes.size(); k++)
590 exp_size_per_plane = exp_size/
m_planes.size();
592 for(
int i = 0; i < exp_size_per_plane; i++)
604 temp_BC_exp = std::dynamic_pointer_cast<
607 i + k * exp_size_per_plane )
611 tmp_V1 = V1 + Phys_offset,
612 tmp_V2 = V2 + Phys_offset,
613 tmp_outarray = outarray + Coef_offset);
625 "Field must be in physical state to extract trace space.");
645 int nPoints_plane =
m_planes[0]->GetTotPoints();
646 int nTracePts =
m_planes[0]->GetTrace()->GetTotPoints();
648 for (
int i = 0; i <
m_planes.size(); ++i)
654 &inarray[i*nPoints_plane], 1,
655 &inarray_plane[0], 1);
657 m_planes[i]->ExtractTracePhys(inarray_plane, outarray_plane);
660 &outarray_plane[0], 1,
661 &outarray[i*nTracePts], 1);
680 for (
int j = 0; j < coordim; ++j)
687 for(
int n = 0; n < i; ++n )
698 elmt =
GetExp(ElmtID[cnt+n]);
700 = elmt->GetTraceNormal(EdgeID[cnt+n]);
702 for (
int j = 0; j < expdim; ++j)
705 tmp = normals[j] + offset, 1);
715 const int nPlanes =
m_planes.size();
716 const int nTracePlane =
m_planes[0]->GetTrace()->GetExpSize();
721 = traceMap->GetBndCondIDToGlobalTraceID();
722 int mapSize = traceBndMap.size();
727 int i, n, e, cnt = 0, cnt1 = 0;
732 int nPlaneExp = nExp / nPlanes;
740 for (n = 0; n < nPlanes; ++n)
742 const int offset = n * nTracePlane;
743 for (e = 0; e < nPlaneExp; ++e)
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Describes the specification for a Basis.
NekDouble Evaluate() const
void GetTracePhysVals(const int trace, const StdRegions::StdExpansionSharedPtr &TraceExp, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, StdRegions::Orientation orient=StdRegions::eNoOrientation)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void SetUpDG()
Set up all DG member variables and maps.
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
This funtion extract form a vector containing a full 3D-homogenous-1D field the value associated with...
Array< OneD, int > m_traceBndMap
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
This function calculate the inner product of two vectors (V1 and V2) respect to the basis along a bou...
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
virtual void v_ExtractTracePhys(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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)
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...
Array< OneD, int > m_BCtoEdgMap
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, SpatialDomains::BoundaryConditions &bcs, const std::string variable)
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
virtual ~DisContField3DHomogeneous1D()
Destructor.
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
virtual void 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)
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)
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 MultiRegions::ExpListSharedPtr > & GetBndCondExpansions()
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Set up a list of element ids and edge ids the link to the boundary conditions.
Array< OneD, NekDouble > m_bndCondBndWeight
DisContField3DHomogeneous1D()
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
NekDouble GetSpecVanVisc(const int k)
LibUtilities::CommSharedPtr m_StripZcomm
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
bool m_useFFT
FFT variables.
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n.
int GetExpSize(void)
This function returns the number of elements in the expansion.
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
bool m_physState
The state of the array m_phys.
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
LibUtilities::SessionReaderSharedPtr m_session
Session.
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
int GetPhys_Offset(int n) const
Get the start offset position for a global list of m_phys correspoinding to element n.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
const BoundaryConditionCollection & GetBoundaryConditions(void) const
const BoundaryRegionCollection & GetBoundaryRegions(void) const
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< Expansion1D > Expansion1DSharedPtr
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
std::shared_ptr< DisContField > DisContFieldSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
std::shared_ptr< DirichletBoundaryCondition > DirichletBCShPtr
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::shared_ptr< NeumannBoundaryCondition > NeumannBCShPtr
std::shared_ptr< RobinBoundaryCondition > RobinBCShPtr
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
std::map< ConstFactorType, NekDouble > ConstFactorMap
The above copyright notice and this permission notice shall be included.
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.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)