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 // 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 1D homogeneous direction
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
39 
40 namespace Nektar
41 {
42  namespace MultiRegions
43  {
44 
47  {
48  }
49 
51  const ContField3DHomogeneous1D &In):
53  {
54 
55  bool False = false;
56  ContField2DSharedPtr zero_plane =
57  boost::dynamic_pointer_cast<ContField2D> (In.m_planes[0]);
58 
59  for(int n = 0; n < m_planes.num_elements(); ++n)
60  {
62  AllocateSharedPtr(*zero_plane,False);
63  }
64 
65  SetCoeffPhys();
66  }
67 
69  {
70  }
71 
74  const LibUtilities::BasisKey &HomoBasis,
75  const NekDouble lhom,
76  const bool useFFT,
77  const bool dealiasing,
79  const std::string &variable,
80  const bool CheckIfSingularSystem):
81  DisContField3DHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
82  {
83  int i,n,nel;
84  ContField2DSharedPtr plane_zero;
85  ContField2DSharedPtr plane_two;
86 
88  m_graph = graph2D;
89 
90  // Plane zero (k=0 - cos) - singularaty check required for Poisson
91  // problems
93  pSession, graph2D, variable, false,
94  CheckIfSingularSystem);
95 
97  pSession, graph2D, variable, false,
98  false);
99 
102 
103  for(n = 0; n < m_planes.num_elements(); ++n)
104  {
105  // Plane zero and one (k=0 - cos and sin) - singularaty check
106  // required for Poisson problems
107  if(m_transposition->GetK(n) == 0)
108  {
110  ::AllocateSharedPtr(*plane_zero, graph2D, variable,
111  false, CheckIfSingularSystem);
112  }
113  else
114  {
115  // For k > 0 singularty check not required anymore -
116  // creating another ContField2D to avoid Assembly Map copy
117  // TODO: We may want to deal with it in a more efficient
118  // way in the future.
120  ::AllocateSharedPtr(*plane_two, graph2D, variable,
121  false, false);
122  }
123 
124  nel = m_planes[n]->GetExpSize();
125 
126  for(i = 0; i < nel; ++i)
127  {
128  (*m_exp).push_back(m_planes[n]->GetExp(i));
129  }
130  }
131 
132  nel = GetExpSize();
133 
135  AllocateSharedPtr(nel);
136 
137  SetCoeffPhys();
138 
139  // Do not set up BCs if default variable
140  if(variable.compare("DefaultVar") != 0)
141  {
142  SetupBoundaryConditions(HomoBasis,lhom,bcs,variable);
143  }
144  }
145 
146 
148  Array<OneD,NekDouble>& outarray)
149  {
151  int ncoeffs = m_planes[0]->GetNcoeffs();
152 
153  for(int n = 0; n < m_planes.num_elements(); ++n)
154  {
155  m_planes[n]->ImposeDirichletConditions(tmp = outarray +
156  n*ncoeffs);
157  }
158  }
159 
160  /**
161  *
162  */
164  {
165  for(int n = 0; n < m_planes.num_elements(); ++n)
166  {
167  m_planes[n]->LocalToGlobal();
168  }
169  };
170 
171 
172  /**
173  *
174  */
176  {
177  for(int n = 0; n < m_planes.num_elements(); ++n)
178  {
179  m_planes[n]->GlobalToLocal();
180  }
181  };
182 
183 
184  /**
185  *
186  */
188  Array<OneD,NekDouble> &field)
189  {
190  int cnt = 0;
192 
193  for(int n = 0; n < m_planes.num_elements(); ++n)
194  {
195  m_planes[n]->SmoothField(tmp = field + cnt);
196 
197  cnt += m_planes[n]->GetTotPoints();
198  }
199  }
200 
201 
203  const Array<OneD, const NekDouble> &inarray,
204  Array<OneD, NekDouble> &outarray,
205  const FlagList &flags,
206  const StdRegions::ConstFactorMap &factors,
207  const StdRegions::VarCoeffMap &varcoeff,
208  const Array<OneD, const NekDouble> &dirForcing)
209  {
210 
211  int n;
212  int cnt = 0;
213  int cnt1 = 0;
214  NekDouble beta;
215  StdRegions::ConstFactorMap new_factors;
216 
218  Array<OneD, NekDouble> fce(inarray.num_elements());
219 
220  // Fourier transform forcing function
221  if(m_WaveSpace)
222  {
223  fce = inarray;
224  }
225  else
226  {
227  HomogeneousFwdTrans(inarray, fce,
228  (flags.isSet(eUseGlobal))?eGlobal:eLocal);
229  }
230 
231  bool smode = false;
232 
233  if (m_homogeneousBasis->GetBasisType() ==
235  m_homogeneousBasis->GetBasisType() ==
237  {
238  smode = true;
239  }
240 
241  for(n = 0; n < m_planes.num_elements(); ++n)
242  {
243  if(n != 1 || m_transposition->GetK(n) != 0 || smode)
244  {
245  beta = 2*M_PI*(m_transposition->GetK(n))/m_lhom;
246  new_factors = factors;
247  // add in Homogeneous Fourier direction and SVV if turned on
248  new_factors[StdRegions::eFactorLambda] +=
249  beta*beta*(1+GetSpecVanVisc(n));
250 
251  m_planes[n]->HelmSolve(fce + cnt,
252  e_out = outarray + cnt1,
253  flags, new_factors, varcoeff,
254  dirForcing);
255  }
256 
257  cnt += m_planes[n]->GetTotPoints();
258  cnt1 += m_planes[n]->GetNcoeffs();
259  }
260  }
261  } // end of namespace
262 } //end of namespace
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Local coefficients.
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:971
LibUtilities::TranspositionSharedPtr m_transposition
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1858
boost::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:291
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:248
NekDouble m_lhom
Width of homogeneous direction.
virtual void v_LocalToGlobal(void)
Template method virtual forwarded for LocalToGlobal()
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1837
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
Global coefficients.
void SetCoeffPhys(void)
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
bool isSet(const FlagType &key) const
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:59
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:947
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField2D.h:56
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:225
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:883
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:880
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)
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...
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:60
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()
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:432
Describes the specification for a Basis.
Definition: Basis.h:50