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 
44 #include <MultiRegions/ExpList.h>
46 
47 namespace Nektar
48 {
49 //--------
50 // Substepping
51 // --------
52 
53 class SubSteppingExtrapolate;
54 
55 typedef std::shared_ptr<SubSteppingExtrapolate> SubSteppingExtrapolateSharedPtr;
56 
58 {
59 public:
60  /// Creates an instance of this class
64  MultiRegions::ExpListSharedPtr &pPressure, const Array<OneD, int> &pVel,
65  const SolverUtils::AdvectionSharedPtr &advObject)
66  {
69  pSession, pFields, pPressure, pVel, advObject);
70  return p;
71  }
72 
73  /// Name of class
74  static std::string className;
75 
79  const Array<OneD, int> pVel,
80  const SolverUtils::AdvectionSharedPtr advObject);
81 
82  virtual ~SubSteppingExtrapolate();
83 
84 protected:
85  virtual void v_EvaluatePressureBCs(
86  const Array<OneD, const Array<OneD, NekDouble>> &fields,
87  const Array<OneD, const Array<OneD, NekDouble>> &N,
88  NekDouble kinvis) override;
89 
90  virtual void v_SubSteppingTimeIntegration(
91  const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme)
92  override;
93 
94  virtual void v_SubStepSaveFields(int nstep) override;
95 
96  virtual void v_SubStepSetPressureBCs(
97  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
98  NekDouble Aii_Dt, NekDouble kinvis) override;
99 
100  virtual void v_SubStepAdvance(int nstep, NekDouble time) override;
101 
102  virtual void v_MountHOPBCs(
103  int HBCdata, NekDouble kinvis, Array<OneD, NekDouble> &Q,
105 
106  virtual std::string v_GetSubStepName(void) override;
107 
108  void SubStepAdvection(
109  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
110  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time);
111 
112  void SubStepProjection(
113  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
114  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time);
115 
117  Array<OneD, Array<OneD, NekDouble>> &ExtVel);
118 
120  const Array<OneD, const Array<OneD, NekDouble>> &velfield,
121  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
122  Array<OneD, Array<OneD, NekDouble>> &outarray);
123 
125 
129 
131 
136 };
137 } // namespace Nektar
138 
139 #endif
Binds a set of functions for use by time integration schemes.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
An abstract base class encapsulating the concept of advection of a vector field.
Definition: Advection.h:70
void SubStepExtrapolateField(NekDouble toff, Array< OneD, Array< OneD, NekDouble >> &ExtVel)
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
virtual void v_MountHOPBCs(int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection) override
virtual void v_EvaluatePressureBCs(const Array< OneD, const Array< OneD, NekDouble >> &fields, const Array< OneD, const Array< OneD, NekDouble >> &N, NekDouble kinvis) override
LibUtilities::TimeIntegrationSchemeSharedPtr m_subStepIntegrationScheme
static std::string className
Name of class.
virtual void v_SubStepSaveFields(int nstep) override
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.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
virtual void v_SubSteppingTimeIntegration(const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme) override
virtual void v_SubStepAdvance(int nstep, NekDouble time) override
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 std::string v_GetSubStepName(void) override
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
SubSteppingExtrapolate(const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const SolverUtils::AdvectionSharedPtr advObject)
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
virtual void v_SubStepSetPressureBCs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, NekDouble Aii_Dt, NekDouble kinvis) override
Array< OneD, Array< OneD, NekDouble > > m_previousVelFields
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:278
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< SubSteppingExtrapolate > SubSteppingExtrapolateSharedPtr
std::shared_ptr< Extrapolate > ExtrapolateSharedPtr
Definition: Extrapolate.h:59
double NekDouble