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
39namespace Nektar
40{
41namespace MultiRegions
42{
43
46{
47}
48
52{
53
54 bool False = false;
55 ContFieldSharedPtr zero_plane =
56 std::dynamic_pointer_cast<ContField>(In.m_planes[0]);
57
58 for (int n = 0; n < m_planes.size(); ++n)
59 {
60 m_planes[n] =
62 }
63
65}
66
70 const std::string &variable)
72{
73 ContFieldSharedPtr zero_plane_old =
74 std::dynamic_pointer_cast<ContField>(In.m_planes[0]);
75
77 *zero_plane_old, graph2D, variable);
78
79 for (int n = 0; n < m_planes.size(); ++n)
80 {
82 *zero_plane, graph2D, variable);
83 }
84
86
87 if (variable.compare("DefaultVar") != 0)
88 {
91 variable);
92 }
93}
94
96{
97}
98
101 const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
102 const bool useFFT, const bool dealiasing,
104 const std::string &variable, const bool CheckIfSingularSystem,
106 : DisContField3DHomogeneous1D(pSession, HomoBasis, lhom, useFFT, dealiasing)
107{
108 int i, n, nel;
109 ContFieldSharedPtr plane_zero;
110 ContFieldSharedPtr plane_two;
111
113 m_graph = graph2D;
114
115 // Plane zero (k=0 - cos) - singularaty check required for Poisson
116 // problems
118 pSession, graph2D, variable, false, CheckIfSingularSystem, ImpType);
119
121 pSession, graph2D, variable, false, false, ImpType);
122
124
125 for (n = 0; n < m_planes.size(); ++n)
126 {
127 // Plane zero and one (k=0 - cos and sin) - singularaty check
128 // required for Poisson problems
129 if (m_transposition->GetK(n) == 0)
130 {
132 *plane_zero, graph2D, variable, false, CheckIfSingularSystem);
133 }
134 else
135 {
136 // For k > 0 singularty check not required anymore -
137 // creating another ContField to avoid Assembly Map copy
138 // TODO: We may want to deal with it in a more efficient
139 // way in the future.
141 *plane_two, graph2D, variable, false, false);
142 }
143
144 nel = m_planes[n]->GetExpSize();
145
146 for (i = 0; i < nel; ++i)
147 {
148 (*m_exp).push_back(m_planes[n]->GetExp(i));
149 }
150 }
151
152 nel = GetExpSize();
153
154 SetCoeffPhys();
155
156 // Do not set up BCs if default variable
157 if (variable.compare("DefaultVar") != 0)
158 {
159 SetupBoundaryConditions(HomoBasis, lhom, bcs, variable);
160 }
161}
162
164 Array<OneD, NekDouble> &outarray)
165{
167 int ncoeffs = m_planes[0]->GetNcoeffs();
168
169 for (int n = 0; n < m_planes.size(); ++n)
170 {
171 m_planes[n]->ImposeDirichletConditions(tmp = outarray + n * ncoeffs);
172 }
173}
174
176 const Array<OneD, NekDouble> coeffs)
177{
178 int numcoeffs_per_plane = m_planes[0]->GetNcoeffs();
179 for (int n = 0; n < m_planes.size(); ++n)
180 {
181 m_planes[n]->FillBndCondFromField(coeffs + n * numcoeffs_per_plane);
182 }
183}
184
186 const int nreg, const Array<OneD, NekDouble> coeffs)
187{
188 int numcoeffs_per_plane = m_planes[0]->GetNcoeffs();
189 for (int n = 0; n < m_planes.size(); ++n)
190 {
191 m_planes[n]->FillBndCondFromField(nreg,
192 coeffs + n * numcoeffs_per_plane);
193 }
194}
195
196/**
197 *
198 */
200{
201 for (int n = 0; n < m_planes.size(); ++n)
202 {
203 m_planes[n]->LocalToGlobal(useComm);
204 }
205}
206
207/**
208 *
209 */
211{
212 for (int n = 0; n < m_planes.size(); ++n)
213 {
214 m_planes[n]->GlobalToLocal();
215 }
216}
217
218/**
219 *
220 */
222{
223 int cnt = 0;
225
226 for (int n = 0; n < m_planes.size(); ++n)
227 {
228 m_planes[n]->SmoothField(tmp = field + cnt);
229 cnt += m_planes[n]->GetTotPoints();
230 }
231}
232
234 const Array<OneD, const NekDouble> &inarray,
236 const StdRegions::VarCoeffMap &varcoeff,
237 const MultiRegions::VarFactorsMap &varfactors,
238 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
239{
240
241 int n;
242 int cnt = 0;
243 int cnt1 = 0;
245 StdRegions::ConstFactorMap new_factors;
246
247 int npts_fce = PhysSpaceForcing ? m_npoints : m_ncoeffs;
249 Array<OneD, NekDouble> fce(npts_fce);
251
252 GlobalLinSysKey gkey(NullGlobalLinSysKey); // Default: return Null
253
254 // Fourier transform forcing function
255 if (m_WaveSpace)
256 {
257 fce = inarray;
258 }
259 else
260 {
261 HomogeneousFwdTrans(npts_fce, inarray, fce);
262 }
263
264 bool smode = false;
265
266 if (m_homogeneousBasis->GetBasisType() ==
269 {
270 smode = true;
271 }
272
274 {
275 ContFieldSharedPtr zero_plane =
276 std::dynamic_pointer_cast<ContField>(m_planes[0]);
277 for (n = 1; n < m_planes.size(); ++n)
278 {
279 std::dynamic_pointer_cast<ContField>(m_planes[n])
280 ->SetGJPForcing(zero_plane->GetGJPForcing());
281 }
282 }
283
284 for (n = 0; n < m_planes.size(); ++n)
285 {
286 if (n != 1 || m_transposition->GetK(n) != 0 || smode)
287 {
288 beta = 2 * M_PI * (m_transposition->GetK(n)) / m_lhom;
289 new_factors = factors;
290 // add in Homogeneous Fourier direction and SVV if turned on
291 new_factors[StdRegions::eFactorLambda] +=
292 beta * beta * (1 + GetSpecVanVisc(n));
293
294 wfce = (PhysSpaceForcing) ? fce + cnt : fce + cnt1;
295 auto gkey = m_planes[n]->HelmSolve(
296 wfce, e_out = outarray + cnt1, new_factors, varcoeff,
297 varfactors, dirForcing, PhysSpaceForcing);
298 }
299
300 cnt += m_planes[n]->GetTotPoints();
301 cnt1 += m_planes[n]->GetNcoeffs();
302 }
303 return gkey;
304}
305
306/**
307 * Reset the GlobalLinSys Manager
308 */
310{
311 for (int n = 0; n < m_planes.size(); ++n)
312 {
313 m_planes[n]->ClearGlobalLinSysManager();
314 }
315}
316
317} // namespace MultiRegions
318} // namespace Nektar
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_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray) override
virtual void v_LocalToGlobal(bool useComm) override
Template method virtual forwarded for LocalToGlobal()
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
Solves the three-dimensional Helmholtz equation, subject to the boundary conditions specified.
virtual void v_GlobalToLocal(void) override
Template method virtual forwarded for GlobalToLocal()
virtual void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs) override
virtual 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:2061
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
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1069
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1067
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
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:70
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:68
static GlobalLinSysKey NullGlobalLinSysKey(StdRegions::eNoMatrixType)
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:270
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
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