Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  const ContField3DHomogeneous1D &In,
71  const std::string &variable):
72  DisContField3DHomogeneous1D (In, false)
73  {
74  ContField2DSharedPtr zero_plane_old =
75  boost::dynamic_pointer_cast<ContField2D> (In.m_planes[0]);
76 
77  ContField2DSharedPtr zero_plane =
79  AllocateSharedPtr(*zero_plane_old,graph2D,
80  variable);
81 
82  for(int n = 0; n < m_planes.num_elements(); ++n)
83  {
85  AllocateSharedPtr(*zero_plane,graph2D,
86  variable);
87  }
88 
89  SetCoeffPhys();
90 
91  if(variable.compare("DefaultVar") != 0)
92  {
95  m_lhom, bcs,variable);
96  }
97  }
98 
100  {
101  }
102 
104  const LibUtilities::SessionReaderSharedPtr &pSession,
105  const LibUtilities::BasisKey &HomoBasis,
106  const NekDouble lhom,
107  const bool useFFT,
108  const bool dealiasing,
109  const SpatialDomains::MeshGraphSharedPtr &graph2D,
110  const std::string &variable,
111  const bool CheckIfSingularSystem):
112  DisContField3DHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing)
113  {
114  int i,n,nel;
115  ContField2DSharedPtr plane_zero;
116  ContField2DSharedPtr plane_two;
117 
119  m_graph = graph2D;
120 
121  // Plane zero (k=0 - cos) - singularaty check required for Poisson
122  // problems
124  pSession, graph2D, variable, false,
125  CheckIfSingularSystem);
126 
128  pSession, graph2D, variable, false,
129  false);
130 
133 
134  for(n = 0; n < m_planes.num_elements(); ++n)
135  {
136  // Plane zero and one (k=0 - cos and sin) - singularaty check
137  // required for Poisson problems
138  if(m_transposition->GetK(n) == 0)
139  {
141  ::AllocateSharedPtr(*plane_zero, graph2D, variable,
142  false, CheckIfSingularSystem);
143  }
144  else
145  {
146  // For k > 0 singularty check not required anymore -
147  // creating another ContField2D to avoid Assembly Map copy
148  // TODO: We may want to deal with it in a more efficient
149  // way in the future.
151  ::AllocateSharedPtr(*plane_two, graph2D, variable,
152  false, false);
153  }
154 
155  nel = m_planes[n]->GetExpSize();
156 
157  for(i = 0; i < nel; ++i)
158  {
159  (*m_exp).push_back(m_planes[n]->GetExp(i));
160  }
161  }
162 
163  nel = GetExpSize();
164 
166  AllocateSharedPtr(nel);
167 
168  SetCoeffPhys();
169 
170  // Do not set up BCs if default variable
171  if(variable.compare("DefaultVar") != 0)
172  {
173  SetupBoundaryConditions(HomoBasis,lhom,bcs,variable);
174  }
175  }
176 
177 
179  Array<OneD,NekDouble>& outarray)
180  {
182  int ncoeffs = m_planes[0]->GetNcoeffs();
183 
184  for(int n = 0; n < m_planes.num_elements(); ++n)
185  {
186  m_planes[n]->ImposeDirichletConditions(tmp = outarray +
187  n*ncoeffs);
188  }
189  }
190 
192  {
193  for(int n = 0; n < m_planes.num_elements(); ++n)
194  {
195  m_planes[n]->FillBndCondFromField();
196  }
197  }
198 
200  {
201  for(int n = 0; n < m_planes.num_elements(); ++n)
202  {
203  m_planes[n]->FillBndCondFromField(nreg);
204  }
205  }
206 
207  /**
208  *
209  */
211  {
212  for(int n = 0; n < m_planes.num_elements(); ++n)
213  {
214  m_planes[n]->LocalToGlobal(useComm);
215  }
216  };
217 
218 
219  /**
220  *
221  */
223  {
224  for(int n = 0; n < m_planes.num_elements(); ++n)
225  {
226  m_planes[n]->GlobalToLocal();
227  }
228  };
229 
230 
231  /**
232  *
233  */
235  Array<OneD,NekDouble> &field)
236  {
237  int cnt = 0;
239 
240  for(int n = 0; n < m_planes.num_elements(); ++n)
241  {
242  m_planes[n]->SmoothField(tmp = field + cnt);
243 
244  cnt += m_planes[n]->GetTotPoints();
245  }
246  }
247 
248 
250  const Array<OneD, const NekDouble> &inarray,
251  Array<OneD, NekDouble> &outarray,
252  const FlagList &flags,
253  const StdRegions::ConstFactorMap &factors,
254  const StdRegions::VarCoeffMap &varcoeff,
255  const Array<OneD, const NekDouble> &dirForcing,
256  const bool PhysSpaceForcing)
257  {
258 
259  int n;
260  int cnt = 0;
261  int cnt1 = 0;
262  NekDouble beta;
263  StdRegions::ConstFactorMap new_factors;
264 
266  Array<OneD, NekDouble> fce(inarray.num_elements());
268 
269  // Fourier transform forcing function
270  if(m_WaveSpace)
271  {
272  fce = inarray;
273  }
274  else
275  {
276  HomogeneousFwdTrans(inarray, fce,
277  (flags.isSet(eUseGlobal))?eGlobal:eLocal);
278  }
279 
280  bool smode = false;
281 
282  if (m_homogeneousBasis->GetBasisType() ==
284  m_homogeneousBasis->GetBasisType() ==
286  {
287  smode = true;
288  }
289 
290  for(n = 0; n < m_planes.num_elements(); ++n)
291  {
292  if(n != 1 || m_transposition->GetK(n) != 0 || smode)
293  {
294 
295  beta = 2*M_PI*(m_transposition->GetK(n))/m_lhom;
296  new_factors = factors;
297  // add in Homogeneous Fourier direction and SVV if turned on
298  new_factors[StdRegions::eFactorLambda] +=
299  beta*beta*(1+GetSpecVanVisc(n));
300 
301  wfce = (PhysSpaceForcing)? fce+cnt:fce+cnt1;
302  m_planes[n]->HelmSolve(wfce,
303  e_out = outarray + cnt1,
304  flags, new_factors, varcoeff,
305  dirForcing,
306  PhysSpaceForcing);
307  }
308 
309  cnt += m_planes[n]->GetTotPoints();
310  cnt1 += m_planes[n]->GetNcoeffs();
311  }
312  }
313 
314  /**
315  * Reset the GlobalLinSys Manager
316  */
318  {
319  for(int n = 0; n < m_planes.num_elements(); ++n)
320  {
321  m_planes[n]->ClearGlobalLinSysManager();
322  }
323  }
324 
325  } // end of namespace
326 } //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:1060
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:2075
boost::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:287
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
NekDouble m_lhom
Width of homogeneous direction.
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2054
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
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:1036
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:227
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:972
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:969
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: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()
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, const bool PhysSpaceForcing)
Solves the three-dimensional Helmholtz equation, subject to the boundary conditions specified...
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
Describes the specification for a Basis.
Definition: Basis.h:50
virtual void v_LocalToGlobal(bool useComm)
Template method virtual forwarded for LocalToGlobal()