Nektar++
SubSteppingExtrapolate.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: SubSteppingExtrapolate.h
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: Abstract base class for SubSteppingExtrapolate.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_SOLVERS_SUBSTEPPINGEXTRAPOLATE_H
36 #define NEKTAR_SOLVERS_SUBSTEPPINGEXTRAPOLATE_H
37 
41 #include <MultiRegions/ExpList.h>
46 
47 namespace Nektar
48 {
49  //--------
50  // Substepping
51  // --------
52 
54 
55  typedef std::shared_ptr<SubSteppingExtrapolate> SubSteppingExtrapolateSharedPtr;
56 
58  {
59  public:
60 
61  /// Creates an instance of this class
66  const Array<OneD, int> &pVel,
67  const SolverUtils::AdvectionSharedPtr &advObject)
68  {
70  ::AllocateSharedPtr(pSession,pFields,pPressure,pVel,advObject);
71  return p;
72  }
73 
74  /// Name of class
75  static std::string className;
76 
81  const Array<OneD, int> pVel,
82  const SolverUtils::AdvectionSharedPtr advObject);
83 
84  virtual ~SubSteppingExtrapolate();
85 
86  protected:
87  virtual void v_EvaluatePressureBCs(const Array<OneD, const Array<OneD, NekDouble> > &fields,
88  const Array<OneD, const Array<OneD, NekDouble> > &N,
89  NekDouble kinvis);
90 
91  virtual void v_SubSteppingTimeIntegration(
92  int intMethod,
94  &IntegrationScheme);
95 
96  virtual void v_SubStepSaveFields(
97  int nstep);
98 
99  virtual void v_SubStepSetPressureBCs(
100  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
101  NekDouble Aii_Dt,
102  NekDouble kinvis);
103 
104  virtual void v_SubStepAdvance(
105  const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln,
106  int nstep,
107  NekDouble time);
108 
109  virtual void v_MountHOPBCs(
110  int HBCdata,
111  NekDouble kinvis,
114 
116 
117  void SubStepAdvection(
118  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
119  Array<OneD, Array<OneD,NekDouble> > &outarray,
120  const NekDouble time);
121 
122  void SubStepProjection(
123  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
124  Array<OneD, Array<OneD, NekDouble> > &outarray,
125  const NekDouble time);
126 
128  NekDouble toff,
129  Array< OneD, Array<OneD, NekDouble> > &ExtVel);
130 
132  const Array<OneD, const Array<OneD, NekDouble> > &velfield,
133  const Array<OneD, const Array<OneD, NekDouble> > &physfield,Array<OneD,
134  Array<OneD, NekDouble> > &outarray);
135 
137 
140 
142 
147  };
148 }
149 
150 #endif
151 
virtual void v_EvaluatePressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void SubStepExtrapolateField(NekDouble toff, Array< OneD, Array< OneD, NekDouble > > &ExtVel)
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
virtual void v_SubStepSetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, NekDouble Aii_Dt, NekDouble kinvis)
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:170
std::shared_ptr< Extrapolate > ExtrapolateSharedPtr
Definition: Extrapolate.h:59
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble > > &velfield, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &outarray)
virtual void v_SubStepSaveFields(int nstep)
Array< OneD, Array< OneD, NekDouble > > m_previousVelFields
virtual void v_SubSteppingTimeIntegration(int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)
std::shared_ptr< TimeIntegrationWrapper > TimeIntegrationWrapperSharedPtr
static std::string className
Name of class.
std::shared_ptr< SubSteppingExtrapolate > SubSteppingExtrapolateSharedPtr
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual void v_SubStepAdvance(const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln, int nstep, NekDouble time)
double NekDouble
SubSteppingExtrapolate(const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const SolverUtils::AdvectionSharedPtr advObject)
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
virtual void v_MountHOPBCs(int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection)
static ExtrapolateSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, MultiRegions::ExpListSharedPtr &pPressure, const Array< OneD, int > &pVel, const SolverUtils::AdvectionSharedPtr &advObject)
Creates an instance of this class.
virtual LibUtilities::TimeIntegrationMethod v_GetSubStepIntegrationMethod(void)
std::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
An abstract base class encapsulating the concept of advection of a vector field.
Definition: Advection.h:69
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps