Nektar++
ContField3DHomogeneous2D.cpp
Go to the documentation of this file.
1//////////////////////////////////////////////////////////////////////////////
2//
3// File: ContField3DHomogeneous2D.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 and a 2 homogeneous directions
33//
34///////////////////////////////////////////////////////////////////////////////
35
38
40{
41
44{
45}
46
50{
51
52 ContFieldSharedPtr zero_line =
53 std::dynamic_pointer_cast<ContField>(In.m_lines[0]);
54
55 for (int n = 0; n < m_lines.size(); ++n)
56 {
58 }
59
61}
62
64{
65}
66
69 const LibUtilities::BasisKey &HomoBasis_y,
70 const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y,
71 const NekDouble lhom_z, const bool useFFT, const bool dealiasing,
73 const std::string &variable, const Collections::ImplementationType ImpType)
74 : DisContField3DHomogeneous2D(pSession, HomoBasis_y, HomoBasis_z, lhom_y,
75 lhom_z, useFFT, dealiasing, ImpType)
76{
77 int i, n, nel;
78 ContFieldSharedPtr line_zero;
79 SpatialDomains::BoundaryConditions bcs(pSession, graph1D);
80
82 pSession, graph1D, variable, false, false, ImpType);
83
85 nel = m_lines[0]->GetExpSize();
86
87 for (i = 0; i < nel; ++i)
88 {
89 (*m_exp).push_back(m_lines[0]->GetExp(i));
90 }
91
92 int nylines = m_homogeneousBasis_y->GetNumPoints();
93 int nzlines = m_homogeneousBasis_z->GetNumPoints();
94
95 for (n = 1; n < nylines * nzlines; ++n)
96 {
98 pSession, graph1D, variable, false, false, ImpType);
99
100 for (i = 0; i < nel; ++i)
101 {
102 (*m_exp).push_back((*m_exp)[i]);
103 }
104 }
105
106 // Setup Default optimisation information.
107 nel = GetExpSize();
108
109 SetCoeffPhys();
110
111 SetupBoundaryConditions(HomoBasis_y, HomoBasis_z, lhom_y, lhom_z, bcs,
112 variable);
113}
114
116 Array<OneD, NekDouble> &outarray)
117{
119 int ncoeffs = m_lines[0]->GetNcoeffs();
120
121 for (int n = 0; n < m_lines.size(); ++n)
122 {
123 m_lines[n]->ImposeDirichletConditions(tmp = outarray + n * ncoeffs);
124 }
125}
126
127/**
128 *
129 */
131{
132 for (int n = 0; n < m_lines.size(); ++n)
133 {
134 m_lines[n]->LocalToGlobal(useComm);
135 }
136}
137
138/**
139 *
140 */
142{
143 for (int n = 0; n < m_lines.size(); ++n)
144 {
145 m_lines[n]->GlobalToLocal();
146 }
147}
148
150 const Array<OneD, const NekDouble> &inarray,
152 const StdRegions::VarCoeffMap &varcoeff,
153 const MultiRegions::VarFactorsMap &varfactors,
154 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
155{
156 int n, m;
157 int cnt = 0;
158 int cnt1 = 0;
159 int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
160 int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
161 NekDouble beta_y;
162 NekDouble beta_z;
164 StdRegions::ConstFactorMap new_factors;
165
166 int npts_fce = PhysSpaceForcing ? m_npoints : m_ncoeffs;
168 Array<OneD, NekDouble> fce(npts_fce);
170
171 GlobalLinSysKey gkey(NullGlobalLinSysKey); // Default: return Null
172
173 if (m_WaveSpace)
174 {
175 fce = inarray;
176 }
177 else
178 {
179 // Fourier transform forcing function
180 HomogeneousFwdTrans(npts_fce, inarray, fce);
181 }
182
183 int l = 0;
184 for (n = 0; n < nhom_modes_z; ++n)
185 {
186 for (m = 0; m < nhom_modes_y; ++m, l++)
187 {
188 beta_z = 2 * M_PI * (n / 2) / m_lhom_z;
189 beta_y = 2 * M_PI * (m / 2) / m_lhom_y;
190 beta = beta_y * beta_y + beta_z * beta_z;
191 new_factors = factors;
192 new_factors[StdRegions::eFactorLambda] += beta;
193
194 wfce = (PhysSpaceForcing) ? fce + cnt : fce + cnt1;
195 auto gkey = m_lines[l]->HelmSolve(wfce, e_out = outarray + cnt1,
196 new_factors, varcoeff, varfactors,
197 dirForcing, PhysSpaceForcing);
198
199 cnt += m_lines[l]->GetTotPoints();
200 cnt1 += m_lines[l]->GetNcoeffs();
201 }
202 }
203 return gkey;
204}
205
206/**
207 * Reset the GlobalLinSys Manager
208 */
210{
211 for (int n = 0; n < m_lines.size(); ++n)
212 {
213 m_lines[n]->ClearGlobalLinSysManager();
214 }
215}
216
217} // namespace Nektar::MultiRegions
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_GlobalToLocal(void) override
Template method virtual forwarded for GlobalToLocal()
void v_LocalToGlobal(bool useComm) override
Template method virtual forwarded for LocalToGlobal()
void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray) 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
Solves the three-dimensional Helmholtz equation, subject to the boundary conditions specified.
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, SpatialDomains::BoundaryConditions &bcs, const std::string variable)
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
NekDouble m_lhom_z
Width of homogeneous direction z.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
NekDouble m_lhom_y
Width of homogeneous direction y.
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2038
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1118
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1063
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2070
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1840
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
static GlobalLinSysKey NullGlobalLinSysKey(StdRegions::eNoMatrixType)
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:268
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
double NekDouble