Nektar++
ContField3DHomogeneous2D.cpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ContField3DHomogeneous2D.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 2 homogeneous directions
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <MultiRegions/ContField.h>
38 
39 namespace Nektar
40 {
41 namespace MultiRegions
42 {
43 
46 {
47 }
48 
50  const ContField3DHomogeneous2D &In)
51  : DisContField3DHomogeneous2D(In, false)
52 {
53 
54  ContFieldSharedPtr zero_line =
55  std::dynamic_pointer_cast<ContField>(In.m_lines[0]);
56 
57  for (int n = 0; n < m_lines.size(); ++n)
58  {
60  }
61 
62  SetCoeffPhys();
63 }
64 
66 {
67 }
68 
71  const LibUtilities::BasisKey &HomoBasis_y,
72  const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y,
73  const NekDouble lhom_z, const bool useFFT, const bool dealiasing,
75  const std::string &variable, const Collections::ImplementationType ImpType)
76  : DisContField3DHomogeneous2D(pSession, HomoBasis_y, HomoBasis_z, lhom_y,
77  lhom_z, useFFT, dealiasing, ImpType)
78 {
79  int i, n, nel;
80  ContFieldSharedPtr line_zero;
81  SpatialDomains::BoundaryConditions bcs(pSession, graph1D);
82 
84  pSession, graph1D, variable, false, false, ImpType);
85 
87  nel = m_lines[0]->GetExpSize();
88 
89  for (i = 0; i < nel; ++i)
90  {
91  (*m_exp).push_back(m_lines[0]->GetExp(i));
92  }
93 
94  int nylines = m_homogeneousBasis_y->GetNumPoints();
95  int nzlines = m_homogeneousBasis_z->GetNumPoints();
96 
97  for (n = 1; n < nylines * nzlines; ++n)
98  {
100  pSession, graph1D, variable, false, false, ImpType);
101 
102  for (i = 0; i < nel; ++i)
103  {
104  (*m_exp).push_back((*m_exp)[i]);
105  }
106  }
107 
108  // Setup Default optimisation information.
109  nel = GetExpSize();
110 
111  SetCoeffPhys();
112 
113  SetupBoundaryConditions(HomoBasis_y, HomoBasis_z, lhom_y, lhom_z, bcs,
114  variable);
115 }
116 
118  Array<OneD, NekDouble> &outarray)
119 {
121  int ncoeffs = m_lines[0]->GetNcoeffs();
122 
123  for (int n = 0; n < m_lines.size(); ++n)
124  {
125  m_lines[n]->ImposeDirichletConditions(tmp = outarray + n * ncoeffs);
126  }
127 }
128 
129 /**
130  *
131  */
133 {
134  for (int n = 0; n < m_lines.size(); ++n)
135  {
136  m_lines[n]->LocalToGlobal(useComm);
137  }
138 }
139 
140 /**
141  *
142  */
144 {
145  for (int n = 0; n < m_lines.size(); ++n)
146  {
147  m_lines[n]->GlobalToLocal();
148  }
149 }
150 
152  const Array<OneD, const NekDouble> &inarray,
153  Array<OneD, NekDouble> &outarray, const StdRegions::ConstFactorMap &factors,
154  const StdRegions::VarCoeffMap &varcoeff,
155  const MultiRegions::VarFactorsMap &varfactors,
156  const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
157 {
158  int n, m;
159  int cnt = 0;
160  int cnt1 = 0;
161  int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
162  int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
163  NekDouble beta_y;
164  NekDouble beta_z;
165  NekDouble beta;
166  StdRegions::ConstFactorMap new_factors;
167 
168  int npts_fce = PhysSpaceForcing ? m_npoints : m_ncoeffs;
170  Array<OneD, NekDouble> fce(npts_fce);
172 
173  GlobalLinSysKey gkey(NullGlobalLinSysKey); // Default: return Null
174 
175  if (m_WaveSpace)
176  {
177  fce = inarray;
178  }
179  else
180  {
181  // Fourier transform forcing function
182  HomogeneousFwdTrans(npts_fce, inarray, fce);
183  }
184 
185  int l = 0;
186  for (n = 0; n < nhom_modes_z; ++n)
187  {
188  for (m = 0; m < nhom_modes_y; ++m, l++)
189  {
190  beta_z = 2 * M_PI * (n / 2) / m_lhom_z;
191  beta_y = 2 * M_PI * (m / 2) / m_lhom_y;
192  beta = beta_y * beta_y + beta_z * beta_z;
193  new_factors = factors;
194  new_factors[StdRegions::eFactorLambda] += beta;
195 
196  wfce = (PhysSpaceForcing) ? fce + cnt : fce + cnt1;
197  auto gkey = m_lines[l]->HelmSolve(wfce, e_out = outarray + cnt1,
198  new_factors, varcoeff, varfactors,
199  dirForcing, PhysSpaceForcing);
200 
201  cnt += m_lines[l]->GetTotPoints();
202  cnt1 += m_lines[l]->GetNcoeffs();
203  }
204  }
205  return gkey;
206 }
207 
208 /**
209  * Reset the GlobalLinSys Manager
210  */
212 {
213  for (int n = 0; n < m_lines.size(); ++n)
214  {
215  m_lines[n]->ClearGlobalLinSysManager();
216  }
217 }
218 
219 } // namespace MultiRegions
220 } // 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_GlobalToLocal(void) override
Template method virtual forwarded for GlobalToLocal()
virtual void v_LocalToGlobal(bool useComm) override
Template method virtual forwarded for LocalToGlobal()
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray) override
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.
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, 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...
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
NekDouble m_lhom_z
Width of homogeneous direction z.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
NekDouble m_lhom_y
Width of homogeneous direction y.
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1996
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1111
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1056
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2029
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1798
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61
static GlobalLinSysKey NullGlobalLinSysKey(StdRegions::eNoMatrixType)
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:289
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:343
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble