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 
36 #include <MultiRegions/ContField.h>
38 
39 namespace Nektar
40 {
41 namespace MultiRegions
42 {
43 
46 {
47 }
48 
50  const ContField3DHomogeneous1D &In)
51  : DisContField3DHomogeneous1D(In, false)
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 
64  SetCoeffPhys();
65 }
66 
68  const ContField3DHomogeneous1D &In,
70  const std::string &variable)
71  : DisContField3DHomogeneous1D(In, false)
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 
85  SetCoeffPhys();
86 
87  if (variable.compare("DefaultVar") != 0)
88  {
91  variable);
92  }
93 }
94 
96 {
97 }
98 
100  const LibUtilities::SessionReaderSharedPtr &pSession,
101  const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom,
102  const bool useFFT, const bool dealiasing,
103  const SpatialDomains::MeshGraphSharedPtr &graph2D,
104  const std::string &variable, const bool CheckIfSingularSystem,
105  const Collections::ImplementationType ImpType)
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 {
177  for (int n = 0; n < m_planes.size(); ++n)
178  {
179  m_planes[n]->FillBndCondFromField();
180  }
181 }
182 
184 {
185  for (int n = 0; n < m_planes.size(); ++n)
186  {
187  m_planes[n]->FillBndCondFromField(nreg);
188  }
189 }
190 
191 /**
192  *
193  */
195 {
196  for (int n = 0; n < m_planes.size(); ++n)
197  {
198  m_planes[n]->LocalToGlobal(useComm);
199  }
200 }
201 
202 /**
203  *
204  */
206 {
207  for (int n = 0; n < m_planes.size(); ++n)
208  {
209  m_planes[n]->GlobalToLocal();
210  }
211 }
212 
213 /**
214  *
215  */
217 {
218  int cnt = 0;
220 
221  for (int n = 0; n < m_planes.size(); ++n)
222  {
223  m_planes[n]->SmoothField(tmp = field + cnt);
224 
225  cnt += m_planes[n]->GetTotPoints();
226  }
227 }
228 
230  const Array<OneD, const NekDouble> &inarray,
231  Array<OneD, NekDouble> &outarray, const StdRegions::ConstFactorMap &factors,
232  const StdRegions::VarCoeffMap &varcoeff,
233  const MultiRegions::VarFactorsMap &varfactors,
234  const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
235 {
236 
237  int n;
238  int cnt = 0;
239  int cnt1 = 0;
240  NekDouble beta;
241  StdRegions::ConstFactorMap new_factors;
242 
244  Array<OneD, NekDouble> fce(inarray.size());
246 
247  // Fourier transform forcing function
248  if (m_WaveSpace)
249  {
250  fce = inarray;
251  }
252  else
253  {
254  HomogeneousFwdTrans(inarray, fce);
255  }
256 
257  bool smode = false;
258 
259  if (m_homogeneousBasis->GetBasisType() ==
262  {
263  smode = true;
264  }
265 
266  if (factors.count(StdRegions::eFactorGJP))
267  {
268  ContFieldSharedPtr zero_plane =
269  std::dynamic_pointer_cast<ContField>(m_planes[0]);
270  for (n = 1; n < m_planes.size(); ++n)
271  {
272  std::dynamic_pointer_cast<ContField>(m_planes[n])
273  ->SetGJPForcing(zero_plane->GetGJPForcing());
274  }
275  }
276 
277  for (n = 0; n < m_planes.size(); ++n)
278  {
279  if (n != 1 || m_transposition->GetK(n) != 0 || smode)
280  {
281 
282  beta = 2 * M_PI * (m_transposition->GetK(n)) / m_lhom;
283  new_factors = factors;
284  // add in Homogeneous Fourier direction and SVV if turned on
285  new_factors[StdRegions::eFactorLambda] +=
286  beta * beta * (1 + GetSpecVanVisc(n));
287 
288  wfce = (PhysSpaceForcing) ? fce + cnt : fce + cnt1;
289  m_planes[n]->HelmSolve(wfce, e_out = outarray + cnt1, new_factors,
290  varcoeff, varfactors, dirForcing,
291  PhysSpaceForcing);
292  }
293 
294  cnt += m_planes[n]->GetTotPoints();
295  cnt1 += m_planes[n]->GetNcoeffs();
296  }
297 }
298 
299 /**
300  * Reset the GlobalLinSys Manager
301  */
303 {
304  for (int n = 0; n < m_planes.size(); ++n)
305  {
306  m_planes[n]->ClearGlobalLinSysManager();
307  }
308 }
309 
310 } // namespace MultiRegions
311 } // namespace Nektar
Describes the specification for a Basis.
Definition: Basis.h:50
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)
virtual void v_GlobalToLocal(void)
Template method virtual forwarded for GlobalToLocal()
virtual void v_LocalToGlobal(bool useComm)
Template method virtual forwarded for LocalToGlobal()
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
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)
Solves the three-dimensional Helmholtz equation, subject to the boundary conditions specified.
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...
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2204
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1196
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2223
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1132
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1129
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
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:277
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:240
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:282
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble