Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Field definition for 3D domain with boundary
33 // conditions and a 2 homogeneous directions
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
39 
40 namespace Nektar
41 {
42  namespace MultiRegions
43  {
44 
47  {
48  }
49 
52  {
53 
54  ContField1DSharedPtr zero_line = boost::dynamic_pointer_cast<ContField1D> (In.m_lines[0]);
55 
56  for(int n = 0; n < m_lines.num_elements(); ++n)
57  {
59  }
60 
61  SetCoeffPhys();
62  }
63 
65  {
66  }
67 
69  const LibUtilities::BasisKey &HomoBasis_y,
70  const LibUtilities::BasisKey &HomoBasis_z,
71  const NekDouble lhom_y,
72  const NekDouble lhom_z,
73  const bool useFFT,
74  const bool dealiasing,
76  const std::string &variable):
77  DisContField3DHomogeneous2D(pSession,HomoBasis_y,HomoBasis_z,lhom_y,lhom_z,useFFT,dealiasing)
78  {
79  int i,n,nel;
80  ContField1DSharedPtr line_zero;
81  SpatialDomains::BoundaryConditions bcs(pSession, graph1D);
82 
83  m_lines[0] = line_zero = MemoryManager<ContField1D>::AllocateSharedPtr(pSession,graph1D,variable);
84 
86  nel = m_lines[0]->GetExpSize();
87 
88  for(i = 0; i < nel; ++i)
89  {
90  (*m_exp).push_back(m_lines[0]->GetExp(i));
91  }
92 
93  int nylines = m_homogeneousBasis_y->GetNumPoints();
94  int nzlines = m_homogeneousBasis_z->GetNumPoints();
95 
96  for(n = 1; n < nylines*nzlines; ++n)
97  {
98  m_lines[n] = MemoryManager<ContField1D>::AllocateSharedPtr(pSession,graph1D,variable);
99 
100  for(i = 0; i < nel; ++i)
101  {
102  (*m_exp).push_back((*m_exp)[i]);
103  }
104  }
105 
106  // Setup Default optimisation information.
107  nel = GetExpSize();
108 
111 
112  SetCoeffPhys();
113 
114  SetupBoundaryConditions(HomoBasis_y,HomoBasis_z,lhom_y,lhom_z,bcs);
115  }
116 
117 
119  {
121  int ncoeffs = m_lines[0]->GetNcoeffs();
122 
123  for(int n = 0; n < m_lines.num_elements(); ++n)
124  {
125  m_lines[n]->ImposeDirichletConditions(tmp = outarray +
126  n*ncoeffs);
127  }
128  }
129 
130 
131  /**
132  *
133  */
135  {
136  for(int n = 0; n < m_lines.num_elements(); ++n)
137  {
138  m_lines[n]->LocalToGlobal();
139  }
140  };
141 
142 
143  /**
144  *
145  */
147  {
148  for(int n = 0; n < m_lines.num_elements(); ++n)
149  {
150  m_lines[n]->GlobalToLocal();
151  }
152  };
153 
154 
156  const Array<OneD, const NekDouble> &inarray,
157  Array<OneD, NekDouble> &outarray,
158  const FlagList &flags,
159  const StdRegions::ConstFactorMap &factors,
160  const StdRegions::VarCoeffMap &varcoeff,
161  const Array<OneD, const NekDouble> &dirForcing)
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.num_elements());
175 
176  if(m_WaveSpace)
177  {
178  fce = inarray;
179  }
180  else
181  {
182  // Fourier transform forcing function
183  HomogeneousFwdTrans(inarray,fce,(flags.isSet(eUseGlobal))?eGlobal:eLocal);
184  }
185 
186  int l =0;
187  for(n = 0; n < nhom_modes_z; ++n)
188  {
189  for(m = 0; m < nhom_modes_y; ++m, l++)
190  {
191  beta_z = 2*M_PI*(n/2)/m_lhom_z;
192  beta_y = 2*M_PI*(m/2)/m_lhom_y;
193  beta = beta_y*beta_y + beta_z*beta_z;
194  new_factors = factors;
195  new_factors[StdRegions::eFactorLambda] += beta;
196 
197  m_lines[l]->HelmSolve(fce + cnt,
198  e_out = outarray + cnt1,
199  flags, new_factors, varcoeff, dirForcing);
200 
201  cnt += m_lines[l]->GetTotPoints();
202 
203  cnt1 += m_lines[l]->GetNcoeffs();
204 
205  }
206  }
207  }
208 
209  /**
210  * Reset the GlobalLinSys Manager
211  */
213  {
214  for(int n = 0; n < m_lines.num_elements(); ++n)
215  {
216  m_lines[n]->ClearGlobalLinSysManager();
217  }
218  }
219 
220  } // end of namespace
221 } //end of namespace
void HomogeneousFwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
boost::shared_ptr< ContField1D > ContField1DSharedPtr
Definition: ContField1D.h:237
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 Array< OneD, const NekDouble > &dirForcing)
Solves the three-dimensional Helmholtz equation, subject to the boundary conditions specified...
Local coefficients.
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1001
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1917
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis_y, const LibUtilities::BasisKey &HomoBasis_z, const NekDouble lhom_y, const NekDouble lhom_z, SpatialDomains::BoundaryConditions &bcs)
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
NekDouble m_lhom_z
Width of homogeneous direction z.
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1896
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
Global coefficients.
Array< OneD, ExpListSharedPtr > m_lines
Vector of ExpList, will be filled with ExpList1D.
virtual void v_GlobalToLocal(void)
Template method virtual forwarded for GlobalToLocal()
bool isSet(const FlagType &key) const
LibUtilities::BasisSharedPtr m_homogeneousBasis_z
Base expansion in z direction.
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
Abstraction of a global continuous one-dimensional spectral/hp element expansion which approximates t...
Definition: ContField1D.h:56
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
NekDouble m_lhom_y
Width of homogeneous direction y.
virtual void v_LocalToGlobal(void)
Template method virtual forwarded for LocalToGlobal()
double NekDouble
Defines a list of flags.
LibUtilities::BasisSharedPtr m_homogeneousBasis_y
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
Describes the specification for a Basis.
Definition: Basis.h:50