Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // 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: Abstract base class for SubSteppingExtrapolate.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_SOLVERS_SUBSTEPPINGEXTRAPOLATE_H
37 #define NEKTAR_SOLVERS_SUBSTEPPINGEXTRAPOLATE_H
38 
42 #include <MultiRegions/ExpList.h>
47 
48 namespace Nektar
49 {
50  //--------
51  // Substepping
52  // --------
53 
55 
56  typedef boost::shared_ptr<SubSteppingExtrapolate> SubSteppingExtrapolateSharedPtr;
57 
59  {
60  public:
61 
62  /// Creates an instance of this class
67  const Array<OneD, int> &pVel,
68  const SolverUtils::AdvectionSharedPtr &advObject)
69  {
71  ::AllocateSharedPtr(pSession,pFields,pPressure,pVel,advObject);
72  return p;
73  }
74 
75  /// Name of class
76  static std::string className;
77 
82  const Array<OneD, int> pVel,
83  const SolverUtils::AdvectionSharedPtr advObject);
84 
85  virtual ~SubSteppingExtrapolate();
86 
87  protected:
88  virtual void v_EvaluatePressureBCs(const Array<OneD, const Array<OneD, NekDouble> > &fields,
89  const Array<OneD, const Array<OneD, NekDouble> > &N,
90  NekDouble kinvis);
91 
92  virtual void v_SubSteppingTimeIntegration(
93  int intMethod,
95  &IntegrationScheme);
96 
97  virtual void v_SubStepSaveFields(
98  int nstep);
99 
100  virtual void v_SubStepSetPressureBCs(
101  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
102  NekDouble Aii_Dt,
103  NekDouble kinvis);
104 
105  virtual void v_SubStepAdvance(
106  const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln,
107  int nstep,
108  NekDouble time);
109 
110  virtual void v_MountHOPBCs(
111  int HBCdata,
112  NekDouble kinvis,
114  Array<OneD, const NekDouble> &Advection);
115 
117 
118  void AddDuDt(void);
119 
120 
121  void SubStepAdvection(
122  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
123  Array<OneD, Array<OneD,NekDouble> > &outarray,
124  const NekDouble time);
125 
126  void SubStepProjection(
127  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
128  Array<OneD, Array<OneD, NekDouble> > &outarray,
129  const NekDouble time);
130 
132  NekDouble toff,
133  Array< OneD, Array<OneD, NekDouble> > &ExtVel);
134 
136  const Array<OneD, const Array<OneD, NekDouble> > &velfield,
137  const Array<OneD, const Array<OneD, NekDouble> > &physfield,Array<OneD,
138  Array<OneD, NekDouble> > &outarray);
139 
141 
144 
146 
151  };
152 }
153 
154 #endif
155 
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
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void v_SubStepSetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, NekDouble Aii_Dt, NekDouble kinvis)
boost::shared_ptr< TimeIntegrationWrapper > TimeIntegrationWrapperSharedPtr
boost::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:158
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble > > &velfield, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &outarray)
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
virtual void v_SubStepSaveFields(int nstep)
Array< OneD, Array< OneD, NekDouble > > m_previousVelFields
virtual void v_SubSteppingTimeIntegration(int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)
static std::string className
Name of class.
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
boost::shared_ptr< Extrapolate > ExtrapolateSharedPtr
Definition: Extrapolate.h:72
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.
boost::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr
virtual LibUtilities::TimeIntegrationMethod v_GetSubStepIntegrationMethod(void)
boost::shared_ptr< SubSteppingExtrapolate > SubSteppingExtrapolateSharedPtr
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps