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 
37 #include <MultiRegions/ContField.h>
38 
39 namespace Nektar
40 {
41  namespace MultiRegions
42  {
43 
46  {
47  }
48 
50  const ContField3DHomogeneous1D &In):
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  {
61  AllocateSharedPtr(*zero_plane,False);
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 
76  ContFieldSharedPtr zero_plane =
78  AllocateSharedPtr(*zero_plane_old,graph2D,
79  variable);
80 
81  for(int n = 0; n < m_planes.size(); ++n)
82  {
84  AllocateSharedPtr(*zero_plane,graph2D,
85  variable);
86  }
87 
88  SetCoeffPhys();
89 
90  if(variable.compare("DefaultVar") != 0)
91  {
94  m_lhom, bcs,variable);
95  }
96  }
97 
99  {
100  }
101 
103  const LibUtilities::SessionReaderSharedPtr &pSession,
104  const LibUtilities::BasisKey &HomoBasis,
105  const NekDouble lhom,
106  const bool useFFT,
107  const bool dealiasing,
108  const SpatialDomains::MeshGraphSharedPtr &graph2D,
109  const std::string &variable,
110  const bool CheckIfSingularSystem,
111  const Collections::ImplementationType ImpType):
112  DisContField3DHomogeneous1D(pSession,HomoBasis,lhom,
113  useFFT,dealiasing)
114  {
115  int i,n,nel;
116  ContFieldSharedPtr plane_zero;
117  ContFieldSharedPtr plane_two;
118 
120  m_graph = graph2D;
121 
122  // Plane zero (k=0 - cos) - singularaty check required for Poisson
123  // problems
125  pSession, graph2D, variable, false,
126  CheckIfSingularSystem, ImpType);
127 
129  pSession, graph2D, variable, false,
130  false, ImpType);
131 
134 
135  for(n = 0; n < m_planes.size(); ++n)
136  {
137  // Plane zero and one (k=0 - cos and sin) - singularaty check
138  // required for Poisson problems
139  if(m_transposition->GetK(n) == 0)
140  {
142  ::AllocateSharedPtr(*plane_zero, graph2D, variable,
143  false, CheckIfSingularSystem);
144  }
145  else
146  {
147  // For k > 0 singularty check not required anymore -
148  // creating another ContField to avoid Assembly Map copy
149  // TODO: We may want to deal with it in a more efficient
150  // way in the future.
152  ::AllocateSharedPtr(*plane_two, graph2D, variable,
153  false, false);
154  }
155 
156  nel = m_planes[n]->GetExpSize();
157 
158  for(i = 0; i < nel; ++i)
159  {
160  (*m_exp).push_back(m_planes[n]->GetExp(i));
161  }
162  }
163 
164  nel = GetExpSize();
165 
166  SetCoeffPhys();
167 
168  // Do not set up BCs if default variable
169  if(variable.compare("DefaultVar") != 0)
170  {
171  SetupBoundaryConditions(HomoBasis,lhom,bcs,variable);
172  }
173  }
174 
175 
177  Array<OneD,NekDouble>& outarray)
178  {
180  int ncoeffs = m_planes[0]->GetNcoeffs();
181 
182  for(int n = 0; n < m_planes.size(); ++n)
183  {
184  m_planes[n]->ImposeDirichletConditions(tmp = outarray +
185  n*ncoeffs);
186  }
187  }
188 
190  {
191  for(int n = 0; n < m_planes.size(); ++n)
192  {
193  m_planes[n]->FillBndCondFromField();
194  }
195  }
196 
198  {
199  for(int n = 0; n < m_planes.size(); ++n)
200  {
201  m_planes[n]->FillBndCondFromField(nreg);
202  }
203  }
204 
205  /**
206  *
207  */
209  {
210  for(int n = 0; n < m_planes.size(); ++n)
211  {
212  m_planes[n]->LocalToGlobal(useComm);
213  }
214  }
215 
216 
217  /**
218  *
219  */
221  {
222  for(int n = 0; n < m_planes.size(); ++n)
223  {
224  m_planes[n]->GlobalToLocal();
225  }
226  }
227 
228 
229  /**
230  *
231  */
233  Array<OneD,NekDouble> &field)
234  {
235  int cnt = 0;
237 
238  for(int n = 0; n < m_planes.size(); ++n)
239  {
240  m_planes[n]->SmoothField(tmp = field + cnt);
241 
242  cnt += m_planes[n]->GetTotPoints();
243  }
244  }
245 
246 
248  const Array<OneD, const NekDouble> &inarray,
249  Array<OneD, NekDouble> &outarray,
250  const StdRegions::ConstFactorMap &factors,
251  const StdRegions::VarCoeffMap &varcoeff,
252  const MultiRegions::VarFactorsMap &varfactors,
253  const Array<OneD, const NekDouble> &dirForcing,
254  const bool PhysSpaceForcing)
255  {
256 
257  int n;
258  int cnt = 0;
259  int cnt1 = 0;
260  NekDouble beta;
261  StdRegions::ConstFactorMap new_factors;
262 
264  Array<OneD, NekDouble> fce(inarray.size());
266 
267  // Fourier transform forcing function
268  if(m_WaveSpace)
269  {
270  fce = inarray;
271  }
272  else
273  {
274  HomogeneousFwdTrans(inarray, fce);
275  }
276 
277  bool smode = false;
278 
279  if (m_homogeneousBasis->GetBasisType() ==
281  m_homogeneousBasis->GetBasisType() ==
283  {
284  smode = true;
285  }
286 
287  for(n = 0; n < m_planes.size(); ++n)
288  {
289  if(n != 1 || m_transposition->GetK(n) != 0 || smode)
290  {
291 
292  beta = 2*M_PI*(m_transposition->GetK(n))/m_lhom;
293  new_factors = factors;
294  // add in Homogeneous Fourier direction and SVV if turned on
295  new_factors[StdRegions::eFactorLambda] +=
296  beta*beta*(1+GetSpecVanVisc(n));
297 
298  wfce = (PhysSpaceForcing)? fce+cnt:fce+cnt1;
299  m_planes[n]->HelmSolve(wfce,
300  e_out = outarray + cnt1,
301  new_factors, varcoeff,
302  varfactors, dirForcing,
303  PhysSpaceForcing);
304  }
305 
306  cnt += m_planes[n]->GetTotPoints();
307  cnt1 += m_planes[n]->GetNcoeffs();
308  }
309  }
310 
311  /**
312  * Reset the GlobalLinSys Manager
313  */
315  {
316  for(int n = 0; n < m_planes.size(); ++n)
317  {
318  m_planes[n]->ClearGlobalLinSysManager();
319  }
320  }
321 
322  } // end of namespace
323 } //end of namespace
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:2401
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1290
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2422
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1226
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1223
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode
Definition: BasisType.h:61
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode
Definition: BasisType.h:60
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:292
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble