Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DisContField3DHomogeneous1D.cpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 //
3 // File DisContField3DHomogeneous1D.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 using LDG flux and a 1D homogeneous direction
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
40 
41 
42 namespace Nektar
43 {
44  namespace MultiRegions
45  {
46 
49  m_bndCondExpansions(),
50  m_bndConditions()
51  {
52  }
53 
56  const LibUtilities::BasisKey &HomoBasis,
57  const NekDouble lhom,
58  const bool useFFT,
59  const bool dealiasing)
60  : ExpList3DHomogeneous1D(pSession,HomoBasis,lhom,useFFT,dealiasing),
61  m_bndCondExpansions(),
62  m_bndConditions()
63  {
64  }
65 
68  const bool DeclarePlanesSetCoeffPhys)
69  : ExpList3DHomogeneous1D (In,false),
70  m_bndCondExpansions (In.m_bndCondExpansions),
71  m_bndConditions (In.m_bndConditions)
72  {
73  if (DeclarePlanesSetCoeffPhys)
74  {
75  DisContField2DSharedPtr zero_plane =
76  boost::dynamic_pointer_cast<DisContField2D> (In.m_planes[0]);
77 
78  for(int n = 0; n < m_planes.num_elements(); ++n)
79  {
80  m_planes[n] =
82  AllocateSharedPtr(*zero_plane, false);
83  }
84 
85  SetCoeffPhys();
86  }
87  }
88 
91  const LibUtilities::BasisKey &HomoBasis,
92  const NekDouble lhom,
93  const bool useFFT,
94  const bool dealiasing,
96  const std::string &variable)
97  : ExpList3DHomogeneous1D(pSession, HomoBasis, lhom, useFFT,
98  dealiasing),
99  m_bndCondExpansions(),
100  m_bndConditions()
101  {
102  int i, n, nel;
103  DisContField2DSharedPtr plane_zero;
105 
106  // note that nzplanes can be larger than nzmodes
107  m_planes[0] = plane_zero = MemoryManager<DisContField2D>::
108  AllocateSharedPtr(pSession, graph2D, variable, true, false);
109 
112 
113  nel = m_planes[0]->GetExpSize();
114 
115  for (i = 0; i < nel; ++i)
116  {
117  (*m_exp).push_back(m_planes[0]->GetExp(i));
118  }
119 
120  for (n = 1; n < m_planes.num_elements(); ++n)
121  {
123  AllocateSharedPtr(*plane_zero, graph2D,
124  variable, true, false);
125  for(i = 0; i < nel; ++i)
126  {
127  (*m_exp).push_back((*m_exp)[i]);
128  }
129  }
130 
131  // Set up trace object.
132  Array<OneD, ExpListSharedPtr> trace(m_planes.num_elements());
133  for (n = 0; n < m_planes.num_elements(); ++n)
134  {
135  trace[n] = m_planes[n]->GetTrace();
136  }
137 
139  pSession, HomoBasis, lhom, useFFT, dealiasing, trace);
140 
141  // Setup default optimisation information
142  nel = GetExpSize();
143 
145  AllocateSharedPtr(nel);
146 
147  SetCoeffPhys();
148 
149  // Do not set up BCs if default variable
150  if(variable.compare("DefaultVar") != 0)
151  {
152  SetupBoundaryConditions(HomoBasis, lhom, bcs, variable);
153  }
154 
155  SetUpDG();
156  }
157 
158  /**
159  * @brief Default destructor
160  */
162  {
163  }
164 
166  const LibUtilities::BasisKey &HomoBasis,
167  const NekDouble lhom,
169  const std::string variable)
170  {
171  int n;
172 
173  // Setup an ExpList2DHomogeneous1D expansion for boundary
174  // conditions and link to class declared in m_planes
176  bcs.GetBoundaryRegions();
178  bcs.GetBoundaryConditions();
179  SpatialDomains::BoundaryRegionCollection::const_iterator it;
180 
181  // count the number of non-periodic boundary regions
182  int cnt = 0;
183  for (it = bregions.begin(); it != bregions.end(); ++it)
184  {
185  SpatialDomains::BoundaryConditionShPtr boundaryCondition =
186  GetBoundaryCondition(bconditions, it->first, variable);
187  if (boundaryCondition->GetBoundaryConditionType()
189  {
190  cnt++;
191  }
192  }
193 
195  ExpListSharedPtr>(cnt);
196  m_bndConditions = m_planes[0]->UpdateBndConditions();
197 
198  int nplanes = m_planes.num_elements();
200  PlanesBndCondExp(nplanes);
201 
202  cnt = 0;
203  for (it = bregions.begin(); it != bregions.end(); ++it)
204  {
205  SpatialDomains::BoundaryConditionShPtr boundaryCondition =
206  GetBoundaryCondition(bconditions, it->first, variable);
207  if(boundaryCondition->GetBoundaryConditionType() !=
209  {
210  for (n = 0; n < nplanes; ++n)
211  {
212  PlanesBndCondExp[n] = m_planes[n]->
214  }
215 
216  m_bndCondExpansions[cnt++] =
218  AllocateSharedPtr(m_session, HomoBasis, lhom,
219  m_useFFT, false,
220  PlanesBndCondExp);
221  }
222  }
223  EvaluateBoundaryConditions(0.0, variable);
224  }
225 
227  const NekDouble time,
228  const std::string varName)
229  {
230  int n;
232  Array<OneD, NekDouble> local_z(m_planes.num_elements());
233 
234  for (n = 0; n < m_planes.num_elements(); n++)
235  {
236  local_z[n] = z[m_transposition->GetPlaneID(n)];
237  }
238 
239  for (n = 0; n < m_planes.num_elements(); ++n)
240  {
241  m_planes[n]->EvaluateBoundaryConditions(
242  time, varName, 0.5*m_lhom*(1.0+local_z[n]));
243  }
244 
245  // Fourier transform coefficient space boundary values
246  // This will only be undertaken for time dependent
247  // boundary conditions unless time == 0.0 which is the
248  // case when the method is called from the constructor.
249  for (n = 0; n < m_bndCondExpansions.num_elements(); ++n)
250  {
251  if (time == 0.0 ||
252  m_bndConditions[n]->IsTimeDependent() ||
253  boost::iequals(m_bndConditions[n]->GetUserDefined(),
254  "MovingBody"))
255  {
256  m_bndCondExpansions[n]->HomogeneousFwdTrans(
259  }
260  }
261  }
262 
264  const Array<OneD, const NekDouble> &inarray,
265  Array<OneD, NekDouble> &outarray,
266  const FlagList &flags,
267  const StdRegions::ConstFactorMap &factors,
268  const StdRegions::VarCoeffMap &varcoeff,
269  const Array<OneD, const NekDouble> &dirForcing)
270  {
271  int n;
272  int cnt = 0;
273  int cnt1 = 0;
274  NekDouble beta;
275  StdRegions::ConstFactorMap new_factors;
276 
278  Array<OneD, NekDouble> fce(inarray.num_elements());
279 
280  // Transform forcing function in half-physical space if required
281  if (m_WaveSpace)
282  {
283  fce = inarray;
284  }
285  else
286  {
287  HomogeneousFwdTrans(inarray,fce);
288  }
289 
290  for (n = 0; n < m_planes.num_elements(); ++n)
291  {
292  if (n != 1 || m_transposition->GetK(n) != 0)
293  {
294  beta = 2*M_PI*(m_transposition->GetK(n))/m_lhom;
295  new_factors = factors;
296  // add in Homogeneous Fourier direction and SVV if turned on
297  new_factors[StdRegions::eFactorLambda] +=
298  beta*beta*(1+GetSpecVanVisc(n));
299  m_planes[n]->HelmSolve(
300  fce + cnt,
301  e_out = outarray + cnt1,
302  flags, new_factors, varcoeff, dirForcing);
303  }
304 
305  cnt += m_planes[n]->GetTotPoints();
306  cnt1 += m_planes[n]->GetNcoeffs();
307  }
308  }
309 
311  const NekDouble time,
312  const std::string varName,
313  const NekDouble x2_in,
314  const NekDouble x3_in)
315  {
316  EvaluateBoundaryConditions(time, varName);
317  }
318 
319  boost::shared_ptr<ExpList> &DisContField3DHomogeneous1D::
321  {
322  return UpdateBndCondExpansion(i);
323  }
324 
327  {
328  return UpdateBndConditions();
329  }
330 
332  Array<OneD, int> &ElmtID,
333  Array<OneD,int> &EdgeID)
334  {
335 
336  if(m_BCtoElmMap.num_elements() == 0)
337  {
338  Array<OneD, int> ElmtID_tmp;
339  Array<OneD, int> EdgeID_tmp;
340 
341  m_planes[0]->GetBoundaryToElmtMap(ElmtID_tmp, EdgeID_tmp);
342  int nel_per_plane = m_planes[0]->GetExpSize();
343  int nplanes = m_planes.num_elements();
344 
345  int MapSize = ElmtID_tmp.num_elements();
346 
347  ElmtID = Array<OneD, int>(nplanes*MapSize);
348  EdgeID = Array<OneD, int>(nplanes*MapSize);
349 
350  // If this mesh (or partition) has no BCs, skip this step
351  if (MapSize > 0)
352  {
353  int i ,j, n, cnt;
354  int cntPlane = 0;
355  for (cnt=n=0; n < m_bndCondExpansions.num_elements(); ++n)
356  {
357  int planeExpSize = m_planes[0]
358  ->GetBndCondExpansions()[n]
359  ->GetExpSize();
360  for (i = 0; i < planeExpSize ; ++i, ++cntPlane)
361  {
362  for(j = 0; j < nplanes; j++)
363  {
364  ElmtID[cnt+i+j*planeExpSize] =
365  ElmtID_tmp[cntPlane]+j*nel_per_plane;
366  EdgeID[cnt+i+j*planeExpSize] =
367  EdgeID_tmp[cntPlane];
368  }
369  }
370  cnt += m_bndCondExpansions[n]->GetExpSize();
371  }
372 
373  m_BCtoElmMap = Array<OneD, int>(nplanes*MapSize);
374  m_BCtoEdgMap = Array<OneD, int>(nplanes*MapSize);
375 
376  Vmath::Vcopy(nplanes*MapSize,ElmtID,1,m_BCtoElmMap,1);
377  Vmath::Vcopy(nplanes*MapSize,EdgeID,1,m_BCtoEdgMap,1);
378  }
379  }
380  else
381  {
382  int MapSize = m_BCtoElmMap.num_elements();
383 
384  ElmtID = Array<OneD, int>(MapSize);
385  EdgeID = Array<OneD, int>(MapSize);
386 
387  Vmath::Vcopy(MapSize, m_BCtoElmMap, 1, ElmtID, 1);
388  Vmath::Vcopy(MapSize, m_BCtoEdgMap, 1, EdgeID, 1);
389  }
390  }
391 
393  boost::shared_ptr<ExpList> &result)
394  {
395  int n, cnt, nq;
396  int offsetOld, offsetNew;
397  Array<OneD, NekDouble> tmp1, tmp2;
398  std::vector<unsigned int> eIDs;
399 
400  Array<OneD, int> ElmtID,EdgeID;
401  GetBoundaryToElmtMap(ElmtID,EdgeID);
402 
403  // Skip other boundary regions
404  for (cnt = n = 0; n < i; ++n)
405  {
406  cnt += m_bndCondExpansions[n]->GetExpSize();
407  }
408 
409  // Populate eIDs with information from BoundaryToElmtMap
410  for (n = 0; n < m_bndCondExpansions[i]->GetExpSize(); ++n)
411  {
412  eIDs.push_back(ElmtID[cnt+n]);
413  }
414 
415  // Create expansion list
416  result =
418 
419  // Copy phys and coeffs to new explist
420  for (n = 0; n < result->GetExpSize(); ++n)
421  {
422  nq = GetExp(ElmtID[cnt+n])->GetTotPoints();
423  offsetOld = GetPhys_Offset(ElmtID[cnt+n]);
424  offsetNew = result->GetPhys_Offset(n);
425  Vmath::Vcopy(nq, tmp1 = GetPhys()+ offsetOld, 1,
426  tmp2 = result->UpdatePhys()+ offsetNew, 1);
427 
428  nq = GetExp(ElmtID[cnt+n])->GetNcoeffs();
429  offsetOld = GetCoeff_Offset(ElmtID[cnt+n]);
430  offsetNew = result->GetCoeff_Offset(n);
431  Vmath::Vcopy(nq, tmp1 = GetCoeffs()+ offsetOld, 1,
432  tmp2 = result->UpdateCoeffs()+ offsetNew, 1);
433  }
434 
435  // Set wavespace value
436  result->SetWaveSpace(GetWaveSpace());
437  }
438 
440  Array<OneD, NekDouble> &BndVals,
441  const Array<OneD, NekDouble> &TotField,
442  int BndID)
443  {
446 
448  Array<OneD, NekDouble> tmp_BC;
449 
450  int cnt = 0;
451  int pos = 0;
452  int exp_size, exp_size_per_plane, elmtID, boundaryID;
453  int offset, exp_dim;
454 
455  for (int k = 0; k < m_planes.num_elements(); k++)
456  {
457  for (int n = 0; n < m_bndConditions.num_elements(); ++n)
458  {
459  exp_size = m_bndCondExpansions[n]->GetExpSize();
460  exp_size_per_plane = exp_size/m_planes.num_elements();
461 
462  for (int i = 0; i < exp_size_per_plane; i++)
463  {
464  if(n == BndID)
465  {
466  elmtID = m_BCtoElmMap[cnt];
467  boundaryID = m_BCtoEdgMap[cnt];
468  exp_dim = m_bndCondExpansions[n]->
469  GetExp(i+k*exp_size_per_plane)->GetTotPoints();
470  offset = GetPhys_Offset(elmtID);
471  elmt = GetExp(elmtID);
472  temp_BC_exp = boost::dynamic_pointer_cast<
474  m_bndCondExpansions[n]->GetExp(
475  i+k*exp_size_per_plane));
476 
477  elmt->GetEdgePhysVals(boundaryID, temp_BC_exp,
478  tmp_Tot = TotField + offset,
479  tmp_BC = BndVals + pos);
480  pos += exp_dim;
481  }
482  cnt++;
483  }
484  }
485  }
486  }
487 
491  Array<OneD, NekDouble> &outarray,
492  int BndID)
493  {
496 
497  Array<OneD, NekDouble> tmp_V1;
498  Array<OneD, NekDouble> tmp_V2;
499  Array<OneD, NekDouble> tmp_outarray;
500 
501  int cnt = 0;
502  int exp_size, exp_size_per_plane, elmtID, Phys_offset, Coef_offset;
503 
504  for(int k = 0; k < m_planes.num_elements(); k++)
505  {
506  for(int n = 0; n < m_bndConditions.num_elements(); ++n)
507  {
508  exp_size = m_bndCondExpansions[n]->GetExpSize();
509  exp_size_per_plane = exp_size/m_planes.num_elements();
510 
511  for(int i = 0; i < exp_size_per_plane; i++)
512  {
513  if(n == BndID)
514  {
515  elmtID = m_BCtoElmMap[cnt];
516 
517  Phys_offset = m_bndCondExpansions[n]->
518  GetPhys_Offset(i+k*exp_size_per_plane);
519  Coef_offset = m_bndCondExpansions[n]->
520  GetCoeff_Offset(i+k*exp_size_per_plane);
521 
522  elmt = GetExp(elmtID);
523  temp_BC_exp = boost::dynamic_pointer_cast<
525  m_bndCondExpansions[n]->GetExp(
526  i+k*exp_size_per_plane));
527 
528  temp_BC_exp->NormVectorIProductWRTBase(
529  tmp_V1 = V1 + Phys_offset,
530  tmp_V2 = V2 + Phys_offset,
531  tmp_outarray = outarray + Coef_offset);
532  }
533  cnt++;
534  }
535  }
536  }
537  }
538 
540  Array<OneD, NekDouble> &outarray)
541  {
542  ASSERTL1(m_physState == true,
543  "Field must be in physical state to extract trace space.");
544 
545  v_ExtractTracePhys(m_phys, outarray);
546  }
547 
548  /**
549  * @brief This method extracts the trace (edges in 2D) for each plane
550  * from the field @a inarray and puts the values in @a outarray.
551  *
552  * It assumes the field is C0 continuous so that it can overwrite the
553  * edge data when visited by the two adjacent elements.
554  *
555  * @param inarray An array containing the 2D data from which we wish
556  * to extract the edge data.
557  * @param outarray The resulting edge information.
558  */
560  const Array<OneD, const NekDouble> &inarray,
561  Array<OneD, NekDouble> &outarray)
562  {
563  int nPoints_plane = m_planes[0]->GetTotPoints();
564  int nTracePts = m_planes[0]->GetTrace()->GetTotPoints();
565 
566  for (int i = 0; i < m_planes.num_elements(); ++i)
567  {
568  Array<OneD, NekDouble> inarray_plane(nPoints_plane, 0.0);
569  Array<OneD, NekDouble> outarray_plane(nPoints_plane, 0.0);
570 
571  Vmath::Vcopy(nPoints_plane,
572  &inarray[i*nPoints_plane], 1,
573  &inarray_plane[0], 1);
574 
575  m_planes[i]->ExtractTracePhys(inarray_plane, outarray_plane);
576 
577  Vmath::Vcopy(nTracePts,
578  &outarray_plane[0], 1,
579  &outarray[i*nTracePts], 1);
580  }
581  }
582 
583  /**
584  */
586  Array<OneD, Array<OneD, NekDouble> > &normals)
587  {
588  int j, n, cnt, nq;
589  int expdim = GetCoordim(0);
590  int coordim = 3;
593 
594  Array<OneD, int> ElmtID,EdgeID;
595  GetBoundaryToElmtMap(ElmtID,EdgeID);
596 
597  // Initialise result
598  normals = Array<OneD, Array<OneD, NekDouble> > (coordim);
599  for (j = 0; j < coordim; ++j)
600  {
601  normals[j] = Array<OneD, NekDouble> (
602  GetBndCondExpansions()[i]->GetTotPoints(), 0.0);
603  }
604 
605  // Skip other boundary regions
606  for (cnt = n = 0; n < i; ++n)
607  {
608  cnt += GetBndCondExpansions()[n]->GetExpSize();
609  }
610 
611  int offset;
612  for (n = 0; n < GetBndCondExpansions()[i]->GetExpSize(); ++n)
613  {
614  offset = GetBndCondExpansions()[i]->GetPhys_Offset(n);
615  nq = GetBndCondExpansions()[i]->GetExp(n)->GetTotPoints();
616 
617  elmt = GetExp(ElmtID[cnt+n]);
618  const Array<OneD, const Array<OneD, NekDouble> > normalsElmt
619  = elmt->GetSurfaceNormal(EdgeID[cnt+n]);
620  // Copy to result
621  for (j = 0; j < expdim; ++j)
622  {
623  Vmath::Vcopy(nq, normalsElmt[j], 1,
624  tmp = normals[j] + offset, 1);
625  }
626  }
627  }
628 
629  /**
630  * @brief Set up all DG member variables and maps.
631  */
633  {
634  const int nPlanes = m_planes.num_elements();
635  const int nTracePlane = m_planes[0]->GetTrace()->GetExpSize();
636 
637  // Get trace map from first plane.
638  AssemblyMapDGSharedPtr traceMap = m_planes[0]->GetTraceMap();
639  const Array<OneD, const int> &traceBndMap
640  = traceMap->GetBndCondTraceToGlobalTraceMap();
641  int mapSize = traceBndMap.num_elements();
642 
643  // Set up trace boundary map
644  m_traceBndMap = Array<OneD, int>(nPlanes * mapSize);
645 
646  int i, n, e, cnt = 0, cnt1 = 0;
647 
648  for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
649  {
650  int nExp = m_bndCondExpansions[i]->GetExpSize();
651  int nPlaneExp = nExp / nPlanes;
652 
653  for (n = 0; n < nPlanes; ++n)
654  {
655  const int offset = n * nTracePlane;
656  for (e = 0; e < nPlaneExp; ++e)
657  {
658  m_traceBndMap[cnt++] = offset + traceBndMap[cnt1+e];
659  }
660  }
661 
662  cnt1 += nPlaneExp;
663  }
664  }
665  } // end of namespace
666 } //end of namespace
boost::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:49
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1834
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n...
Definition: ExpList.h:1926
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1001
LibUtilities::TranspositionSharedPtr m_transposition
int GetPhys_Offset(int n) const
Get the start offset position for a global list of m_phys correspoinding to element n...
Definition: ExpList.h:1934
const BoundaryConditionCollection & GetBoundaryConditions(void) const
Definition: Conditions.h:237
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1953
void SetUpDG()
Set up all DG member variables and maps.
boost::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
const boost::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:1917
void NormVectorIProductWRTBase(const Array< OneD, const NekDouble > &Fx, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:727
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)
Array< OneD, int > m_BCtoElmMap
Storage space for the boundary to element and boundary to trace map. This member variable is really a...
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
NekDouble m_lhom
Width of homogeneous direction.
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:956
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
virtual void v_ExtractTracePhys(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This method extracts the trace (edges in 2D) for each plane from the field inarray and puts the value...
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
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
Definition of the total number of degrees of freedom and quadrature points. Sets up the storage for m...
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
virtual boost::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
boost::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:977
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:965
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:206
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
boost::shared_ptr< DisContField2D > DisContField2DSharedPtr
boost::shared_ptr< StdExpansion1D > StdExpansion1DSharedPtr
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:910
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)
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:215
void SetupBoundaryConditions(const LibUtilities::BasisKey &HomoBasis, const NekDouble lhom, SpatialDomains::BoundaryConditions &bcs, const std::string variable)
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:1887
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
const Array< OneD, const MultiRegions::ExpListSharedPtr > & GetBndCondExpansions()
void GetEdgePhysVals(const int edge, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Extract the physical values along edge edge from inarray into outarray following the local edge orien...
Definition: StdExpansion.h:871
virtual void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
This funtion extract form a vector containing a full 3D-homogenous-1D field the value associated with...
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
This function calculate the inner product of two vectors (V1 and V2) respect to the basis along a bou...
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Set up a list of element ids and edge ids the link to the boundary conditions.
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:2860
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:227
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1502
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1794
boost::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:208
virtual void v_GetBndElmtExpansion(int i, boost::shared_ptr< ExpList > &result)
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="")
This function evaluates the boundary conditions at a certaintime-level.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
Abstraction of a two-dimensional multi-elemental expansion which is merely a collection of local expa...
Describes the specification for a Basis.
Definition: Basis.h:50
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...