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 
37 #include <MultiRegions/ContField.h>
38 
39 namespace Nektar
40 {
41  namespace MultiRegions
42  {
43 
46  {
47  }
48 
50  const ContField3DHomogeneous2D &In):
52  {
53 
54  ContFieldSharedPtr zero_line = std::dynamic_pointer_cast<ContField> (In.m_lines[0]);
55 
56  for(int n = 0; n < m_lines.size(); ++n)
57  {
59  }
60 
61  SetCoeffPhys();
62  }
63 
65  {
66  }
67 
70  const LibUtilities::BasisKey &HomoBasis_y,
71  const LibUtilities::BasisKey &HomoBasis_z,
72  const NekDouble lhom_y,
73  const NekDouble lhom_z,
74  const bool useFFT,
75  const bool dealiasing,
77  const std::string &variable,
78  const Collections::ImplementationType ImpType):
79  DisContField3DHomogeneous2D(pSession,HomoBasis_y,HomoBasis_z,lhom_y,lhom_z,useFFT,dealiasing,ImpType)
80  {
81  int i,n,nel;
82  ContFieldSharedPtr line_zero;
83  SpatialDomains::BoundaryConditions bcs(pSession, graph1D);
84 
85  m_lines[0] = line_zero = MemoryManager<ContField>::AllocateSharedPtr(pSession,graph1D,variable,ImpType);
86 
88  nel = m_lines[0]->GetExpSize();
89 
90  for(i = 0; i < nel; ++i)
91  {
92  (*m_exp).push_back(m_lines[0]->GetExp(i));
93  }
94 
95  int nylines = m_homogeneousBasis_y->GetNumPoints();
96  int nzlines = m_homogeneousBasis_z->GetNumPoints();
97 
98  for(n = 1; n < nylines*nzlines; ++n)
99  {
100  m_lines[n] = MemoryManager<ContField>::AllocateSharedPtr(pSession,graph1D,variable,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  }
115 
116 
118  {
120  int ncoeffs = m_lines[0]->GetNcoeffs();
121 
122  for(int n = 0; n < m_lines.size(); ++n)
123  {
124  m_lines[n]->ImposeDirichletConditions(tmp = outarray +
125  n*ncoeffs);
126  }
127  }
128 
129 
130  /**
131  *
132  */
134  {
135  for(int n = 0; n < m_lines.size(); ++n)
136  {
137  m_lines[n]->LocalToGlobal(useComm);
138  }
139  }
140 
141 
142  /**
143  *
144  */
146  {
147  for(int n = 0; n < m_lines.size(); ++n)
148  {
149  m_lines[n]->GlobalToLocal();
150  }
151  }
152 
153 
155  const Array<OneD, const NekDouble> &inarray,
156  Array<OneD, NekDouble> &outarray,
157  const StdRegions::ConstFactorMap &factors,
158  const StdRegions::VarCoeffMap &varcoeff,
159  const MultiRegions::VarFactorsMap &varfactors,
160  const Array<OneD, const NekDouble> &dirForcing,
161  const bool PhysSpaceForcing)
162  {
163  int n,m;
164  int cnt = 0;
165  int cnt1 = 0;
166  int nhom_modes_y = m_homogeneousBasis_y->GetNumModes();
167  int nhom_modes_z = m_homogeneousBasis_z->GetNumModes();
168  NekDouble beta_y;
169  NekDouble beta_z;
170  NekDouble beta;
171  StdRegions::ConstFactorMap new_factors;
172 
174  Array<OneD, NekDouble> fce(inarray.size());
176 
177  if(m_WaveSpace)
178  {
179  fce = inarray;
180  }
181  else
182  {
183  // Fourier transform forcing function
184  HomogeneousFwdTrans(inarray,fce);
185  }
186 
187  int l =0;
188  for(n = 0; n < nhom_modes_z; ++n)
189  {
190  for(m = 0; m < nhom_modes_y; ++m, l++)
191  {
192  beta_z = 2*M_PI*(n/2)/m_lhom_z;
193  beta_y = 2*M_PI*(m/2)/m_lhom_y;
194  beta = beta_y*beta_y + beta_z*beta_z;
195  new_factors = factors;
196  new_factors[StdRegions::eFactorLambda] += beta;
197 
198  wfce = (PhysSpaceForcing)? fce+cnt:fce+cnt1;
199  m_lines[l]->HelmSolve(wfce,
200  e_out = outarray + cnt1,
201  new_factors, varcoeff, varfactors,
202  dirForcing,
203  PhysSpaceForcing);
204 
205  cnt += m_lines[l]->GetTotPoints();
206  cnt1 += m_lines[l]->GetNcoeffs();
207  }
208  }
209  }
210 
211  /**
212  * Reset the GlobalLinSys Manager
213  */
215  {
216  for(int n = 0; n < m_lines.size(); ++n)
217  {
218  m_lines[n]->ClearGlobalLinSysManager();
219  }
220  }
221 
222  } // end of namespace
223 } //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_GlobalToLocal(void)
Template method virtual forwarded for GlobalToLocal()
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.
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
virtual void v_LocalToGlobal(bool useComm)
Template method virtual forwarded for LocalToGlobal()
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, SpatialDomains::BoundaryConditions &bcs)
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.
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
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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