Nektar++
ContField3DHomogeneous1D.cpp
Go to the documentation of this file.
1//////////////////////////////////////////////////////////////////////////////
2//
3// File: ContField3DHomogeneous1D.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 1D homogeneous direction
33//
34///////////////////////////////////////////////////////////////////////////////
35
38
40{
41
44{
45}
46
50{
51
52 bool False = false;
53 ContFieldSharedPtr zero_plane =
54 std::dynamic_pointer_cast<ContField>(In.m_planes[0]);
55
56 for (int n = 0; n < m_planes.size(); ++n)
57 {
58 m_planes[n] =
60 }
61
63}
64
68 const std::string &variable)
70{
71 ContFieldSharedPtr zero_plane_old =
72 std::dynamic_pointer_cast<ContField>(In.m_planes[0]);
73
75 *zero_plane_old, graph2D, variable);
76
77 for (int n = 0; n < m_planes.size(); ++n)
78 {
80 *zero_plane, graph2D, variable);
81 }
82
84
85 if (variable.compare("DefaultVar") != 0)
86 {
89 variable);
90 }
91}
92
94{
95}
96
99 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
100 const bool useFFT, const bool dealiasing,
102 const std::string &variable, const bool CheckIfSingularSystem,
104 : DisContField3DHomogeneous1D(pSession, HomoBasis, lhom, useFFT, dealiasing)
105{
106 int i, n, nel;
107 ContFieldSharedPtr plane_zero;
108 ContFieldSharedPtr plane_two;
109
111 m_graph = graph2D;
112
113 // Plane zero (k=0 - cos) - singularaty check required for Poisson
114 // problems
116 pSession, graph2D, variable, false, CheckIfSingularSystem, ImpType);
117
119 pSession, graph2D, variable, false, false, ImpType);
120
122
123 for (n = 0; n < m_planes.size(); ++n)
124 {
125 // Plane zero and one (k=0 - cos and sin) - singularaty check
126 // required for Poisson problems
127 if (m_transposition->GetK(n) == 0)
128 {
130 *plane_zero, graph2D, variable, false, CheckIfSingularSystem);
131 }
132 else
133 {
134 // For k > 0 singularty check not required anymore -
135 // creating another ContField to avoid Assembly Map copy
136 // TODO: We may want to deal with it in a more efficient
137 // way in the future.
139 *plane_two, graph2D, variable, false, false);
140 }
141
142 nel = m_planes[n]->GetExpSize();
143
144 for (i = 0; i < nel; ++i)
145 {
146 (*m_exp).push_back(m_planes[n]->GetExp(i));
147 }
148 }
149
150 nel = GetExpSize();
151
152 SetCoeffPhys();
153
154 // Do not set up BCs if default variable
155 if (variable.compare("DefaultVar") != 0)
156 {
157 SetupBoundaryConditions(HomoBasis, lhom, bcs, variable);
158 }
159}
160
162 Array<OneD, NekDouble> &outarray)
163{
165 int ncoeffs = m_planes[0]->GetNcoeffs();
166
167 for (int n = 0; n < m_planes.size(); ++n)
168 {
169 m_planes[n]->ImposeDirichletConditions(tmp = outarray + n * ncoeffs);
170 }
171}
172
174 const Array<OneD, NekDouble> coeffs)
175{
176 int numcoeffs_per_plane = m_planes[0]->GetNcoeffs();
177 for (int n = 0; n < m_planes.size(); ++n)
178 {
179 m_planes[n]->FillBndCondFromField(coeffs + n * numcoeffs_per_plane);
180 }
181}
182
184 const int nreg, const Array<OneD, NekDouble> coeffs)
185{
186 int numcoeffs_per_plane = m_planes[0]->GetNcoeffs();
187 for (int n = 0; n < m_planes.size(); ++n)
188 {
189 m_planes[n]->FillBndCondFromField(nreg,
190 coeffs + n * numcoeffs_per_plane);
191 }
192}
193
194/**
195 *
196 */
198{
199 for (int n = 0; n < m_planes.size(); ++n)
200 {
201 m_planes[n]->LocalToGlobal(useComm);
202 }
203}
204
205/**
206 *
207 */
209{
210 for (int n = 0; n < m_planes.size(); ++n)
211 {
212 m_planes[n]->GlobalToLocal();
213 }
214}
215
216/**
217 *
218 */
220{
221 int cnt = 0;
223
224 for (int n = 0; n < m_planes.size(); ++n)
225 {
226 m_planes[n]->SmoothField(tmp = field + cnt);
227 cnt += m_planes[n]->GetTotPoints();
228 }
229}
230
232 const Array<OneD, const NekDouble> &inarray,
234 const StdRegions::VarCoeffMap &varcoeff,
235 const MultiRegions::VarFactorsMap &varfactors,
236 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
237{
238
239 int n;
240 int cnt = 0;
241 int cnt1 = 0;
243 StdRegions::ConstFactorMap new_factors;
244
245 int npts_fce = PhysSpaceForcing ? m_npoints : m_ncoeffs;
247 Array<OneD, NekDouble> fce(npts_fce);
249
250 GlobalLinSysKey gkey(NullGlobalLinSysKey); // Default: return Null
251
252 // Fourier transform forcing function
253 if (m_WaveSpace)
254 {
255 fce = inarray;
256 }
257 else
258 {
259 HomogeneousFwdTrans(npts_fce, inarray, fce);
260 }
261
262 bool smode = false;
263
264 if (m_homogeneousBasis->GetBasisType() ==
267 {
268 smode = true;
269 }
270
272 {
273 ContFieldSharedPtr zero_plane =
274 std::dynamic_pointer_cast<ContField>(m_planes[0]);
275 for (n = 1; n < m_planes.size(); ++n)
276 {
277 std::dynamic_pointer_cast<ContField>(m_planes[n])
278 ->SetGJPForcing(zero_plane->GetGJPForcing());
279 }
280 }
281
282 for (n = 0; n < m_planes.size(); ++n)
283 {
284 if (n != 1 || m_transposition->GetK(n) != 0 || smode)
285 {
286 beta = 2 * M_PI * (m_transposition->GetK(n)) / m_lhom;
287 new_factors = factors;
288 // add in Homogeneous Fourier direction and SVV if turned on
289 new_factors[StdRegions::eFactorLambda] +=
290 beta * beta * (1 + GetSpecVanVisc(n));
291
292 wfce = (PhysSpaceForcing) ? fce + cnt : fce + cnt1;
293 auto gkey = m_planes[n]->HelmSolve(
294 wfce, e_out = outarray + cnt1, new_factors, varcoeff,
295 varfactors, dirForcing, PhysSpaceForcing);
296 }
297
298 cnt += m_planes[n]->GetTotPoints();
299 cnt1 += m_planes[n]->GetNcoeffs();
300 }
301 return gkey;
302}
303
304/**
305 * Reset the GlobalLinSys Manager
306 */
308{
309 for (int n = 0; n < m_planes.size(); ++n)
310 {
311 m_planes[n]->ClearGlobalLinSysManager();
312 }
313}
314
315} // 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_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray) override
void v_LocalToGlobal(bool useComm) override
Template method virtual forwarded for LocalToGlobal()
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 v_GlobalToLocal(void) override
Template method virtual forwarded for GlobalToLocal()
void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs) override
void v_SmoothField(Array< OneD, NekDouble > &field) override
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, 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...
NekDouble m_lhom
Width of homogeneous direction.
LibUtilities::TranspositionSharedPtr m_transposition
Array< OneD, ExpListSharedPtr > m_planes
LibUtilities::BasisSharedPtr m_homogeneousBasis
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
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
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1060
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1058
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
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:68
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:66
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