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 
39 namespace Nektar
40 {
41  namespace MultiRegions
42  {
43 
46  {
47  }
48 
50  const ContField3DHomogeneous1D &In):
52  {
53 
54  bool False = false;
55  ContField2DSharedPtr zero_plane =
56  std::dynamic_pointer_cast<ContField2D> (In.m_planes[0]);
57 
58  for(int n = 0; n < m_planes.num_elements(); ++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  ContField2DSharedPtr zero_plane_old =
74  std::dynamic_pointer_cast<ContField2D> (In.m_planes[0]);
75 
76  ContField2DSharedPtr zero_plane =
78  AllocateSharedPtr(*zero_plane_old,graph2D,
79  variable);
80 
81  for(int n = 0; n < m_planes.num_elements(); ++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  ContField2DSharedPtr plane_zero;
117  ContField2DSharedPtr 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.num_elements(); ++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 ContField2D 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 
167  AllocateSharedPtr(nel);
168 
169  SetCoeffPhys();
170 
171  // Do not set up BCs if default variable
172  if(variable.compare("DefaultVar") != 0)
173  {
174  SetupBoundaryConditions(HomoBasis,lhom,bcs,variable);
175  }
176  }
177 
178 
180  Array<OneD,NekDouble>& outarray)
181  {
183  int ncoeffs = m_planes[0]->GetNcoeffs();
184 
185  for(int n = 0; n < m_planes.num_elements(); ++n)
186  {
187  m_planes[n]->ImposeDirichletConditions(tmp = outarray +
188  n*ncoeffs);
189  }
190  }
191 
193  {
194  for(int n = 0; n < m_planes.num_elements(); ++n)
195  {
196  m_planes[n]->FillBndCondFromField();
197  }
198  }
199 
201  {
202  for(int n = 0; n < m_planes.num_elements(); ++n)
203  {
204  m_planes[n]->FillBndCondFromField(nreg);
205  }
206  }
207 
208  /**
209  *
210  */
212  {
213  for(int n = 0; n < m_planes.num_elements(); ++n)
214  {
215  m_planes[n]->LocalToGlobal(useComm);
216  }
217  }
218 
219 
220  /**
221  *
222  */
224  {
225  for(int n = 0; n < m_planes.num_elements(); ++n)
226  {
227  m_planes[n]->GlobalToLocal();
228  }
229  }
230 
231 
232  /**
233  *
234  */
236  Array<OneD,NekDouble> &field)
237  {
238  int cnt = 0;
240 
241  for(int n = 0; n < m_planes.num_elements(); ++n)
242  {
243  m_planes[n]->SmoothField(tmp = field + cnt);
244 
245  cnt += m_planes[n]->GetTotPoints();
246  }
247  }
248 
249 
251  const Array<OneD, const NekDouble> &inarray,
252  Array<OneD, NekDouble> &outarray,
253  const FlagList &flags,
254  const StdRegions::ConstFactorMap &factors,
255  const StdRegions::VarCoeffMap &varcoeff,
256  const MultiRegions::VarFactorsMap &varfactors,
257  const Array<OneD, const NekDouble> &dirForcing,
258  const bool PhysSpaceForcing)
259  {
260 
261  int n;
262  int cnt = 0;
263  int cnt1 = 0;
264  NekDouble beta;
265  StdRegions::ConstFactorMap new_factors;
266 
268  Array<OneD, NekDouble> fce(inarray.num_elements());
270 
271  // Fourier transform forcing function
272  if(m_WaveSpace)
273  {
274  fce = inarray;
275  }
276  else
277  {
278  HomogeneousFwdTrans(inarray, fce,
279  (flags.isSet(eUseGlobal))?eGlobal:eLocal);
280  }
281 
282  bool smode = false;
283 
284  if (m_homogeneousBasis->GetBasisType() ==
286  m_homogeneousBasis->GetBasisType() ==
288  {
289  smode = true;
290  }
291 
292  for(n = 0; n < m_planes.num_elements(); ++n)
293  {
294  if(n != 1 || m_transposition->GetK(n) != 0 || smode)
295  {
296 
297  beta = 2*M_PI*(m_transposition->GetK(n))/m_lhom;
298  new_factors = factors;
299  // add in Homogeneous Fourier direction and SVV if turned on
300  new_factors[StdRegions::eFactorLambda] +=
301  beta*beta*(1+GetSpecVanVisc(n));
302 
303  wfce = (PhysSpaceForcing)? fce+cnt:fce+cnt1;
304  m_planes[n]->HelmSolve(wfce,
305  e_out = outarray + cnt1,
306  flags, new_factors, varcoeff,
307  varfactors, dirForcing,
308  PhysSpaceForcing);
309  }
310 
311  cnt += m_planes[n]->GetTotPoints();
312  cnt1 += m_planes[n]->GetNcoeffs();
313  }
314  }
315 
316  /**
317  * Reset the GlobalLinSys Manager
318  */
320  {
321  for(int n = 0; n < m_planes.num_elements(); ++n)
322  {
323  m_planes[n]->ClearGlobalLinSysManager();
324  }
325  }
326 
327  } // end of namespace
328 } //end of namespace
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Local coefficients.
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1090
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1106
LibUtilities::TranspositionSharedPtr m_transposition
std::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:289
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
NekDouble m_lhom
Width of homogeneous direction.
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2170
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2191
Global coefficients.
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
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:60
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField2D.h:55
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:264
bool isSet(const FlagType &key) const
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1026
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1023
Array< OneD, ExpListSharedPtr > m_planes
double NekDouble
Defines a list of flags.
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:61
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, SpatialDomains::BoundaryConditions &bcs, const std::string variable)
virtual void v_GlobalToLocal(void)
Template method virtual forwarded for GlobalToLocal()
virtual void v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, 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...
Describes the specification for a Basis.
Definition: Basis.h:49
std::shared_ptr< SessionReader > SessionReaderSharedPtr
virtual void v_LocalToGlobal(bool useComm)
Template method virtual forwarded for LocalToGlobal()
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap