Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Nektar::SubSteppingExtrapolate Class Reference

#include <SubSteppingExtrapolate.h>

Inheritance diagram for Nektar::SubSteppingExtrapolate:
[legend]

Public Member Functions

 SubSteppingExtrapolate (const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const SolverUtils::AdvectionSharedPtr advObject)
 
 ~SubSteppingExtrapolate () override
 
- Public Member Functions inherited from Nektar::Extrapolate
 Extrapolate (const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const SolverUtils::AdvectionSharedPtr advObject)
 
virtual ~Extrapolate ()
 
void GenerateHOPBCMap (const LibUtilities::SessionReaderSharedPtr &pSsession)
 
void UpdateRobinPrimCoeff (void)
 
void SubSteppingTimeIntegration (const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme)
 
void SubStepSaveFields (const int nstep)
 
void SubStepSetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const NekDouble Aii_DT, NekDouble kinvis)
 
void SubStepAdvance (const int nstep, NekDouble time)
 
void MountHOPBCs (int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection)
 
void EvaluatePressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
 
void SetForcing (const std::vector< SolverUtils::ForcingSharedPtr > &forcing)
 
void AddDuDt (void)
 
void AddVelBC (void)
 
void ExtrapolatePressureHBCs (void)
 
void CopyPressureHBCsToPbndExp (void)
 
Array< OneD, NekDoubleGetMaxStdVelocity (const Array< OneD, Array< OneD, NekDouble > > inarray)
 
void CorrectPressureBCs (const Array< OneD, NekDouble > &pressure)
 
void IProductNormVelocityOnHBC (const Array< OneD, const Array< OneD, NekDouble > > &Vel, Array< OneD, NekDouble > &IprodVn)
 
void IProductNormVelocityBCOnHBC (Array< OneD, NekDouble > &IprodVn)
 
std::string GetSubStepName (void)
 
void ExtrapolateArray (Array< OneD, Array< OneD, NekDouble > > &array)
 
void EvaluateBDFArray (Array< OneD, Array< OneD, NekDouble > > &array)
 
void ExtrapolateArray (Array< OneD, Array< OneD, NekDouble > > &oldarrays, Array< OneD, NekDouble > &newarray, Array< OneD, NekDouble > &outarray)
 
void AddNormVelOnOBC (const int nbcoeffs, const int nreg, Array< OneD, Array< OneD, NekDouble > > &u)
 
void AddPressureToOutflowBCs (NekDouble kinvis)
 
void GenerateBndElmtExpansion (void)
 

Static Public Member Functions

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. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 

Protected Member Functions

void v_EvaluatePressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis) override
 
void v_SubSteppingTimeIntegration (const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme) override
 
void v_SubStepSaveFields (int nstep) override
 
void v_SubStepSetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, NekDouble Aii_Dt, NekDouble kinvis) override
 
void v_AccelerationBDF (Array< OneD, Array< OneD, NekDouble > > &array) override
 
void v_SubStepAdvance (int nstep, NekDouble time) override
 
void v_MountHOPBCs (int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection) override
 
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)
 
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, Array< OneD, NekDouble > &ExtVn)
 
void AddAdvectionPenaltyFlux (const Array< OneD, NekDouble > &Vn, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
NekDouble GetSubstepTimeStep ()
 
- Protected Member Functions inherited from Nektar::Extrapolate
virtual void v_EvaluatePressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)=0
 
virtual void v_SubSteppingTimeIntegration (const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme)=0
 
virtual void v_SubStepSaveFields (int nstep)=0
 
virtual void v_SubStepSetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, NekDouble Aii_DT, NekDouble kinvis)=0
 
virtual void v_SubStepAdvance (int nstep, NekDouble time)=0
 
virtual void v_MountHOPBCs (int HBCdata, NekDouble kinvis, Array< OneD, NekDouble > &Q, Array< OneD, const NekDouble > &Advection)=0
 
virtual std::string v_GetSubStepName (void)
 
virtual void v_AccelerationBDF (Array< OneD, Array< OneD, NekDouble > > &array)
 
void CalcNeumannPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
 
virtual void v_CalcNeumannPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
 
virtual void v_CorrectPressureBCs (const Array< OneD, NekDouble > &pressure)
 
virtual void v_AddNormVelOnOBC (const int nbcoeffs, const int nreg, Array< OneD, Array< OneD, NekDouble > > &u)
 
void CalcOutflowBCs (const Array< OneD, const Array< OneD, NekDouble > > &fields, NekDouble kinvis)
 
void RollOver (Array< OneD, Array< OneD, NekDouble > > &input)
 

Protected Attributes

LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_subStepIntegrationScheme
 
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
 
Array< OneD, Array< OneD, NekDouble > > m_previousVelFields
 
Array< OneD, Array< OneD, NekDouble > > m_previousVnFields
 
NekDouble m_cflSafetyFactor
 
int m_infosteps
 
int m_minsubsteps
 
int m_maxsubsteps
 
- Protected Attributes inherited from Nektar::Extrapolate
LibUtilities::SessionReaderSharedPtr m_session
 
LibUtilities::CommSharedPtr m_comm
 
Array< OneD, HBCTypem_hbcType
 Array of type of high order BCs for splitting shemes. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Velocity fields. More...
 
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field. More...
 
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w); More...
 
SolverUtils::AdvectionSharedPtr m_advObject
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 
Array< OneD, Array< OneD, NekDouble > > m_previousVelFields
 
int m_curl_dim
 Curl-curl dimensionality. More...
 
int m_bnd_dim
 bounday dimensionality More...
 
Array< OneD, const SpatialDomains::BoundaryConditionShPtrm_PBndConds
 pressure boundary conditions container More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_PBndExp
 pressure boundary conditions expansion container More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_bndElmtExps
 Boundary expansions on each domain boundary. More...
 
int m_pressureCalls
 number of times the high-order pressure BCs have been called More...
 
int m_numHBCDof
 
int m_HBCnumber
 
int m_intSteps
 Maximum points used in pressure BC evaluation. More...
 
NekDouble m_timestep
 
Array< OneD, Array< OneD, NekDouble > > m_pressureHBCs
 Storage for current and previous levels of high order pressure boundary conditions. More...
 
Array< OneD, Array< OneD, NekDouble > > m_iprodnormvel
 Storage for current and previous levels of the inner product of normal velocity. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 
HighOrderOutflowSharedPtr m_houtflow
 

Additional Inherited Members

- Static Protected Attributes inherited from Nektar::Extrapolate
static NekDouble StifflyStable_Betaq_Coeffs [3][3]
 
static NekDouble StifflyStable_Alpha_Coeffs [3][3]
 
static NekDouble StifflyStable_Gamma0_Coeffs [3] = {1.0, 1.5, 11.0 / 6.0}
 

Detailed Description

Definition at line 57 of file SubSteppingExtrapolate.h.

Constructor & Destructor Documentation

◆ SubSteppingExtrapolate()

Nektar::SubSteppingExtrapolate::SubSteppingExtrapolate ( const LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields,
MultiRegions::ExpListSharedPtr  pPressure,
const Array< OneD, int >  pVel,
const SolverUtils::AdvectionSharedPtr  advObject 
)

Definition at line 50 of file SubSteppingExtrapolate.cpp.

55 : Extrapolate(pSession, pFields, pPressure, pVel, advObject)
56{
57 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
58 m_session->LoadParameter("SubStepCFL", m_cflSafetyFactor, 0.5);
59 m_session->LoadParameter("MinSubSteps", m_minsubsteps, 1);
60 m_session->LoadParameter("MaxSubSteps", m_maxsubsteps, 100);
61
62 size_t dim = m_fields[0]->GetCoordim(0);
63 m_traceNormals = Array<OneD, Array<OneD, NekDouble>>(dim);
64
65 size_t nTracePts = m_fields[0]->GetTrace()->GetNpoints();
66 for (size_t i = 0; i < dim; ++i)
67 {
68 m_traceNormals[i] = Array<OneD, NekDouble>(nTracePts);
69 }
70 m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
71}
Extrapolate(const LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields, MultiRegions::ExpListSharedPtr pPressure, const Array< OneD, int > pVel, const SolverUtils::AdvectionSharedPtr advObject)
Definition: Extrapolate.cpp:54
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Velocity fields.
Definition: Extrapolate.h:201
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: Extrapolate.h:253
LibUtilities::SessionReaderSharedPtr m_session
Definition: Extrapolate.h:193

References m_cflSafetyFactor, Nektar::Extrapolate::m_fields, m_infosteps, m_maxsubsteps, m_minsubsteps, Nektar::Extrapolate::m_session, and Nektar::Extrapolate::m_traceNormals.

◆ ~SubSteppingExtrapolate()

Nektar::SubSteppingExtrapolate::~SubSteppingExtrapolate ( )
override

Definition at line 73 of file SubSteppingExtrapolate.cpp.

74{
75}

Member Function Documentation

◆ AddAdvectionPenaltyFlux()

void Nektar::SubSteppingExtrapolate::AddAdvectionPenaltyFlux ( const Array< OneD, NekDouble > &  Vn,
const Array< OneD, const Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  Outarray 
)
protected

Add the advection penalty term \( (\hat{u} - u^e)V_n \) given the normal velocity Vn at this time level and the physfield values containing the velocity field at this time level

Number of trace points

Forward state array

Backward state array

upwind numerical flux state array

Extract forwards/backwards trace spaces Note it is important to use the zeroth field but with the specialised method to use boudnary conditions from other fields since trace spaces may not be the same if there are mixed boundary conditions

Upwind between elements

Construct difference between numflux and Fwd,Bwd

Calculate the numerical fluxes multipling Fwd, Bwd and numflux by the normal advection velocity

Definition at line 522 of file SubSteppingExtrapolate.cpp.

526{
527 ASSERTL1(physfield.size() == Outarray.size(),
528 "Physfield and outarray are of different dimensions");
529
530 size_t i;
532 timer.Start();
533
534 /// Number of trace points
535 size_t nTracePts = m_fields[0]->GetTrace()->GetNpoints();
536
537 /// Forward state array
538 Array<OneD, NekDouble> Fwd(3 * nTracePts);
539
540 /// Backward state array
541 Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
542
543 /// upwind numerical flux state array
544 Array<OneD, NekDouble> numflux = Bwd + nTracePts;
545
546 for (i = 0; i < physfield.size(); ++i)
547 {
548 /// Extract forwards/backwards trace spaces
549 /// Note it is important to use the zeroth field but with the
550 /// specialised method to use boudnary conditions from other
551 /// fields since trace spaces may not be the same if there are
552 /// mixed boundary conditions
553 std::dynamic_pointer_cast<MultiRegions::DisContField>(m_fields[0])
554 ->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd,
555 m_fields[i]->GetBndConditions(),
556 m_fields[i]->GetBndCondExpansions());
557
558 /// Upwind between elements
559 m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
560
561 /// Construct difference between numflux and Fwd,Bwd
562 Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
563 Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
564
565 /// Calculate the numerical fluxes multipling Fwd, Bwd and
566 /// numflux by the normal advection velocity
567 Vmath::Vmul(nTracePts, Fwd, 1, Vn, 1, Fwd, 1);
568 Vmath::Vmul(nTracePts, Bwd, 1, Vn, 1, Bwd, 1);
569
570 m_fields[0]->AddFwdBwdTraceIntegral(Fwd, Bwd, Outarray[i]);
571 }
572 timer.Stop();
573 timer.AccumulateRegion("SubSteppingExtrapolate:AddAdvectionPenaltyFlux",
574 10);
575}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References Nektar::LibUtilities::Timer::AccumulateRegion(), ASSERTL1, Nektar::Extrapolate::m_fields, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by SubStepAdvection().

◆ create()

static ExtrapolateSharedPtr Nektar::SubSteppingExtrapolate::create ( const LibUtilities::SessionReaderSharedPtr pSession,
Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
MultiRegions::ExpListSharedPtr pPressure,
const Array< OneD, int > &  pVel,
const SolverUtils::AdvectionSharedPtr advObject 
)
inlinestatic

Creates an instance of this class.

Definition at line 61 of file SubSteppingExtrapolate.h.

66 {
69 pSession, pFields, pPressure, pVel, advObject);
70 return p;
71 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Extrapolate > ExtrapolateSharedPtr
Definition: Extrapolate.h:60

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ GetSubstepTimeStep()

NekDouble Nektar::SubSteppingExtrapolate::GetSubstepTimeStep ( )
protected

Definition at line 481 of file SubSteppingExtrapolate.cpp.

482{
484 timer.Start();
485 size_t n_element = m_fields[0]->GetExpSize();
486
487 const Array<OneD, int> ExpOrder = m_fields[0]->EvalBasisNumModesMaxPerExp();
488
489 const NekDouble cLambda = 0.2; // Spencer book pag. 317
490
491 Array<OneD, NekDouble> tstep(n_element, 0.0);
492 Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
493 Array<OneD, Array<OneD, NekDouble>> velfields(m_velocity.size());
494
495 for (size_t i = 0; i < m_velocity.size(); ++i)
496 {
497 velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
498 }
499 stdVelocity = GetMaxStdVelocity(velfields);
500
501 for (size_t el = 0; el < n_element; ++el)
502 {
503 tstep[el] =
504 m_cflSafetyFactor / (stdVelocity[el] * cLambda *
505 (ExpOrder[el] - 1) * (ExpOrder[el] - 1));
506 }
507
508 NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
509 m_comm->AllReduce(TimeStep, LibUtilities::ReduceMin);
510 timer.Stop();
511 timer.AccumulateRegion("SubSteppingExtrapolate:GetSubStepTimeStep", 10);
512
513 return TimeStep;
514}
Array< OneD, NekDouble > GetMaxStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
Definition: Extrapolate.h:208
LibUtilities::CommSharedPtr m_comm
Definition: Extrapolate.h:195
double NekDouble
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.hpp:725

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::Extrapolate::GetMaxStdVelocity(), m_cflSafetyFactor, Nektar::Extrapolate::m_comm, Nektar::Extrapolate::m_fields, Nektar::Extrapolate::m_velocity, Nektar::LibUtilities::ReduceMin, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Vmath::Vmin().

Referenced by v_SubStepAdvance().

◆ SubStepAdvection()

void Nektar::SubSteppingExtrapolate::SubStepAdvection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Explicit Advection terms used by SubStepAdvance time integration

Get the number of coefficients

Define an auxiliary variable to compute the RHS

Operations to compute the RHS

Multiply the flux by the inverse of the mass matrix

Store in outarray the physical values of the RHS

Definition at line 161 of file SubSteppingExtrapolate.cpp.

164{
165 Nektar::LibUtilities::Timer timer, timer1;
166 size_t i;
167 size_t nVariables = inarray.size();
168 size_t nQuadraturePts = inarray[0].size();
169
170 timer.Start();
171 /// Get the number of coefficients
172 size_t ncoeffs = m_fields[0]->GetNcoeffs();
173
174 /// Define an auxiliary variable to compute the RHS
175 Array<OneD, Array<OneD, NekDouble>> WeakAdv(nVariables);
176 WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs * nVariables);
177 for (i = 1; i < nVariables; ++i)
178 {
179 WeakAdv[i] = WeakAdv[i - 1] + ncoeffs;
180 }
181
182 Array<OneD, Array<OneD, NekDouble>> Velfields(m_velocity.size());
183
184 Velfields[0] = Array<OneD, NekDouble>(nQuadraturePts * m_velocity.size());
185
186 for (i = 1; i < m_velocity.size(); ++i)
187 {
188 Velfields[i] = Velfields[i - 1] + nQuadraturePts;
189 }
190
191 Array<OneD, NekDouble> Vn(m_fields[0]->GetTrace()->GetTotPoints());
192
193 SubStepExtrapolateField(fmod(time, m_timestep), Velfields, Vn);
194
195 for (auto &x : m_forcing)
196 {
197 x->PreApply(m_fields, Velfields, Velfields, time);
198 }
199 timer1.Start();
200 m_advObject->Advect(m_velocity.size(), m_fields, Velfields, inarray,
201 outarray, time);
202 timer1.Stop();
203 timer1.AccumulateRegion("SubStepAdvection:Advect", 10);
204
205 for (auto &x : m_forcing)
206 {
207 x->Apply(m_fields, outarray, outarray, time);
208 }
209
210 for (i = 0; i < nVariables; ++i)
211 {
212 m_fields[i]->IProductWRTBase(outarray[i], WeakAdv[i]);
213 // negation requried due to sign of DoAdvection term to be consistent
214 Vmath::Neg(ncoeffs, WeakAdv[i], 1);
215 }
216
217 AddAdvectionPenaltyFlux(Vn, inarray, WeakAdv);
218
219 /// Operations to compute the RHS
220 for (i = 0; i < nVariables; ++i)
221 {
222 // Negate the RHS
223 Vmath::Neg(ncoeffs, WeakAdv[i], 1);
224
225 /// Multiply the flux by the inverse of the mass matrix
226 m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
227
228 /// Store in outarray the physical values of the RHS
229 m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
230 }
231 timer.Stop();
232 timer.AccumulateRegion("SubSteppingExtrapolate:SubStepAdvection", 10);
233}
NekDouble m_timestep
Definition: Extrapolate.h:243
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Definition: Extrapolate.h:212
SolverUtils::AdvectionSharedPtr m_advObject
Definition: Extrapolate.h:210
void AddAdvectionPenaltyFlux(const Array< OneD, NekDouble > &Vn, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &outarray)
void SubStepExtrapolateField(NekDouble toff, Array< OneD, Array< OneD, NekDouble > > &ExtVel, Array< OneD, NekDouble > &ExtVn)
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292

References Nektar::LibUtilities::Timer::AccumulateRegion(), AddAdvectionPenaltyFlux(), Nektar::Extrapolate::m_advObject, Nektar::Extrapolate::m_fields, Nektar::Extrapolate::m_forcing, Nektar::Extrapolate::m_timestep, Nektar::Extrapolate::m_velocity, Vmath::Neg(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and SubStepExtrapolateField().

Referenced by v_SubSteppingTimeIntegration().

◆ SubStepExtrapolateField()

void Nektar::SubSteppingExtrapolate::SubStepExtrapolateField ( NekDouble  toff,
Array< OneD, Array< OneD, NekDouble > > &  ExtVel,
Array< OneD, NekDouble > &  ExtVn 
)
protected

Extrapolate field using equally time spaced field un,un-1,un-2, (at dt intervals) to time for substep n+t at order m_intSteps. Also extrapolate the normal velocity along the trace at the same order

Definition at line 582 of file SubSteppingExtrapolate.cpp.

585{
586 size_t npts = m_fields[0]->GetTotPoints();
587 size_t ntracepts = m_fields[0]->GetTrace()->GetTotPoints();
588 size_t nvel = m_velocity.size();
589 size_t i, j;
590 Array<OneD, NekDouble> l(4);
592 timer.Start();
593
594 size_t ord = m_intSteps;
595
596 // calculate Lagrange interpolants
597 Vmath::Fill(4, 1.0, l, 1);
598
599 for (i = 0; i <= ord; ++i)
600 {
601 for (j = 0; j <= ord; ++j)
602 {
603 if (i != j)
604 {
605 l[i] *= (j * m_timestep + toff);
606 l[i] /= (j * m_timestep - i * m_timestep);
607 }
608 }
609 }
610
611 for (i = 0; i < nvel; ++i)
612 {
613 Vmath::Smul(npts, l[0], m_previousVelFields[i], 1, ExtVel[i], 1);
614
615 for (j = 1; j <= ord; ++j)
616 {
617 Blas::Daxpy(npts, l[j], m_previousVelFields[j * nvel + i], 1,
618 ExtVel[i], 1);
619 }
620 }
621
622 Vmath::Smul(ntracepts, l[0], m_previousVnFields[0], 1, ExtVn, 1);
623 for (j = 1; j <= ord; ++j)
624 {
625 Blas::Daxpy(ntracepts, l[j], m_previousVnFields[j], 1, ExtVn, 1);
626 }
627 timer.Stop();
628 timer.AccumulateRegion("SubSteppingExtrapolate:SubStepExtrapolateFields",
629 10);
630}
int m_intSteps
Maximum points used in pressure BC evaluation.
Definition: Extrapolate.h:241
Array< OneD, Array< OneD, NekDouble > > m_previousVnFields
Array< OneD, Array< OneD, NekDouble > > m_previousVelFields
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54

References Nektar::LibUtilities::Timer::AccumulateRegion(), Blas::Daxpy(), Vmath::Fill(), Nektar::Extrapolate::m_fields, Nektar::Extrapolate::m_intSteps, m_previousVelFields, m_previousVnFields, Nektar::Extrapolate::m_timestep, Nektar::Extrapolate::m_velocity, Vmath::Smul(), Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by SubStepAdvection().

◆ SubStepProjection()

void Nektar::SubSteppingExtrapolate::SubStepProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Projection used by SubStepAdvance time integration

Definition at line 238 of file SubSteppingExtrapolate.cpp.

242{
244 ASSERTL1(inarray.size() == outarray.size(),
245 "Inarray and outarray of different sizes ");
246
247 timer.Start();
248 for (size_t i = 0; i < inarray.size(); ++i)
249 {
250 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
251 }
252 timer.Stop();
253 timer.AccumulateRegion("SubSteppingExtrapolate:SubStepProjection", 10);
254}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::LibUtilities::Timer::AccumulateRegion(), ASSERTL1, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Vmath::Vcopy().

Referenced by v_SubSteppingTimeIntegration().

◆ v_AccelerationBDF()

void Nektar::SubSteppingExtrapolate::v_AccelerationBDF ( Array< OneD, Array< OneD, NekDouble > > &  array)
overrideprotectedvirtual

At the start, the newest value is stored in array[nlevels-1] and the previous values in the first positions At the end, the acceleration from BDF is stored in array[nlevels-1] and the storage has been updated to included the new value

Reimplemented from Nektar::Extrapolate.

Definition at line 295 of file SubSteppingExtrapolate.cpp.

297{
299 int nlevels = array.size();
300 int nPts = array[0].size();
301
302 timer.Start();
303 if (nPts)
304 {
305 // Update array
306 RollOver(array);
307
308 // Calculate acceleration using Backward Differentiation Formula
309 Array<OneD, NekDouble> accelerationTerm(nPts, 0.0);
310
311 int acc_order = std::min(m_pressureCalls, m_intSteps);
312 Vmath::Smul(nPts, StifflyStable_Gamma0_Coeffs[acc_order - 1], array[0],
313 1, accelerationTerm, 1);
314
315 for (int i = 0; i < acc_order; i++)
316 {
318 nPts, -1 * StifflyStable_Alpha_Coeffs[acc_order - 1][i],
319 array[i + 1], 1, accelerationTerm, 1, accelerationTerm, 1);
320 }
321 array[nlevels - 1] = accelerationTerm;
322 }
323 timer.Stop();
324 timer.AccumulateRegion("SubSteppingExtrapolate:v_AccelerationBDF", 10);
325}
static NekDouble StifflyStable_Alpha_Coeffs[3][3]
Definition: Extrapolate.h:257
void RollOver(Array< OneD, Array< OneD, NekDouble > > &input)
int m_pressureCalls
number of times the high-order pressure BCs have been called
Definition: Extrapolate.h:232
static NekDouble StifflyStable_Gamma0_Coeffs[3]
Definition: Extrapolate.h:258
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::Extrapolate::m_intSteps, Nektar::Extrapolate::m_pressureCalls, Nektar::Extrapolate::RollOver(), Vmath::Smul(), Nektar::LibUtilities::Timer::Start(), Nektar::Extrapolate::StifflyStable_Alpha_Coeffs, Nektar::Extrapolate::StifflyStable_Gamma0_Coeffs, Nektar::LibUtilities::Timer::Stop(), and Vmath::Svtvp().

◆ v_EvaluatePressureBCs()

void Nektar::SubSteppingExtrapolate::v_EvaluatePressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
const Array< OneD, const Array< OneD, NekDouble > > &  N,
NekDouble  kinvis 
)
overrideprotectedvirtual

Implements Nektar::Extrapolate.

Definition at line 77 of file SubSteppingExtrapolate.cpp.

81{
82 ASSERTL0(false, "This method should not be called by Substepping routine");
83}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

References ASSERTL0.

◆ v_GetSubStepName()

std::string Nektar::SubSteppingExtrapolate::v_GetSubStepName ( void  )
overrideprotectedvirtual

Reimplemented from Nektar::Extrapolate.

Definition at line 642 of file SubSteppingExtrapolate.cpp.

643{
644 return m_subStepIntegrationScheme->GetName();
645}
LibUtilities::TimeIntegrationSchemeSharedPtr m_subStepIntegrationScheme

References m_subStepIntegrationScheme.

◆ v_MountHOPBCs()

void Nektar::SubSteppingExtrapolate::v_MountHOPBCs ( int  HBCdata,
NekDouble  kinvis,
Array< OneD, NekDouble > &  Q,
Array< OneD, const NekDouble > &  Advection 
)
overrideprotectedvirtual

Implements Nektar::Extrapolate.

Definition at line 635 of file SubSteppingExtrapolate.cpp.

638{
639 Vmath::Smul(HBCdata, -kinvis, Q, 1, Q, 1);
640}

References Vmath::Smul().

◆ v_SubStepAdvance()

void Nektar::SubSteppingExtrapolate::v_SubStepAdvance ( int  nstep,
NekDouble  time 
)
overrideprotectedvirtual

Implements Nektar::Extrapolate.

Definition at line 415 of file SubSteppingExtrapolate.cpp.

416{
417 int n;
418 int nsubsteps;
420
421 NekDouble dt;
422
423 Array<OneD, Array<OneD, NekDouble>> fields;
424
425 timer.Start();
426 static int ncalls = 1;
427 size_t nint = std::min(ncalls++, m_intSteps);
428
429 // this needs to change
430 m_comm = m_fields[0]->GetComm()->GetRowComm();
431
432 // Get the proper time step with CFL control
433 dt = GetSubstepTimeStep();
434
435 nsubsteps = (m_timestep > dt) ? ((int)(m_timestep / dt) + 1) : 1;
436 nsubsteps = std::max(m_minsubsteps, nsubsteps);
437
438 ASSERTL0(nsubsteps < m_maxsubsteps,
439 "Number of substeps has exceeded maximum");
440
441 dt = m_timestep / nsubsteps;
442
443 if (m_infosteps && !((nstep + 1) % m_infosteps) && m_comm->GetRank() == 0)
444 {
445 std::cout << "Sub-integrating using " << nsubsteps
446 << " steps over Dt = " << m_timestep
447 << " (SubStep CFL=" << m_cflSafetyFactor << ")" << std::endl;
448 }
449
450 const TripleArray &solutionVector = m_intScheme->GetSolutionVector();
451
452 for (size_t m = 0; m < nint; ++m)
453 {
454 // We need to update the fields held by the m_intScheme
455 fields = solutionVector[m];
456
457 // Initialise NS solver which is set up to use a GLM method
458 // with calls to EvaluateAdvection_SetPressureBCs and
459 // SolveUnsteadyStokesSystem
460 m_subStepIntegrationScheme->InitializeScheme(dt, fields, time,
462
463 for (n = 0; n < nsubsteps; ++n)
464 {
465 fields = m_subStepIntegrationScheme->TimeIntegrate(n, dt);
466 }
467
468 // set up HBC m_acceleration field for Pressure BCs
470
471 // Reset time integrated solution in m_intScheme
472 m_intScheme->SetSolutionVector(m, fields);
473 }
474 timer.Stop();
475 timer.AccumulateRegion("SubSteppingExtrapolate:v_SubStepAdvance", 10);
476}
Array< OneD, Array< OneD, NekDouble > > m_iprodnormvel
Storage for current and previous levels of the inner product of normal velocity.
Definition: Extrapolate.h:251
void IProductNormVelocityOnHBC(const Array< OneD, const Array< OneD, NekDouble > > &Vel, Array< OneD, NekDouble > &IprodVn)
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
AT< AT< AT< NekDouble > > > TripleArray

References Nektar::LibUtilities::Timer::AccumulateRegion(), ASSERTL0, GetSubstepTimeStep(), Nektar::Extrapolate::IProductNormVelocityOnHBC(), m_cflSafetyFactor, Nektar::Extrapolate::m_comm, Nektar::Extrapolate::m_fields, m_infosteps, m_intScheme, Nektar::Extrapolate::m_intSteps, Nektar::Extrapolate::m_iprodnormvel, m_maxsubsteps, m_minsubsteps, m_subStepIntegrationOps, m_subStepIntegrationScheme, Nektar::Extrapolate::m_timestep, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

◆ v_SubSteppingTimeIntegration()

void Nektar::SubSteppingExtrapolate::v_SubSteppingTimeIntegration ( const LibUtilities::TimeIntegrationSchemeSharedPtr IntegrationScheme)
overrideprotectedvirtual

Implements Nektar::Extrapolate.

Definition at line 85 of file SubSteppingExtrapolate.cpp.

87{
89 m_intScheme = IntegrationScheme;
90
91 timer.Start();
92 size_t order = IntegrationScheme->GetOrder();
93
94 // Set to 1 for first step and it will then be increased in
95 // time advance routines
96 if ((IntegrationScheme->GetName() == "Euler" &&
97 IntegrationScheme->GetVariant() == "Backward") ||
98
99 (IntegrationScheme->GetName() == "BDFImplicit" &&
100 (order == 1 || order == 2)))
101 {
102 // Note RK first order SSP is just Forward Euler.
103 std::string vSubStepIntScheme = "RungeKutta";
104 std::string vSubStepIntSchemeVariant = "SSP";
105 int vSubStepIntSchemeOrder = order;
106
107 if (m_session->DefinesSolverInfo("SubStepIntScheme"))
108 {
109 vSubStepIntScheme = m_session->GetSolverInfo("SubStepIntScheme");
110 vSubStepIntSchemeVariant = "";
111 vSubStepIntSchemeOrder = order;
112 }
113
116 vSubStepIntScheme, vSubStepIntSchemeVariant,
117 vSubStepIntSchemeOrder, std::vector<NekDouble>());
118
119 size_t nvel = m_velocity.size();
120 size_t ndim = order + 1;
121
122 // Fields for linear/quadratic interpolation
123 m_previousVelFields = Array<OneD, Array<OneD, NekDouble>>(ndim * nvel);
124 int ntotpts = m_fields[0]->GetTotPoints();
125 m_previousVelFields[0] = Array<OneD, NekDouble>(ndim * nvel * ntotpts);
126
127 for (size_t i = 1; i < ndim * nvel; ++i)
128 {
129 m_previousVelFields[i] = m_previousVelFields[i - 1] + ntotpts;
130 }
131 // Vn fields
132 m_previousVnFields = Array<OneD, Array<OneD, NekDouble>>(ndim);
133 ntotpts = m_fields[0]->GetTrace()->GetTotPoints();
134 m_previousVnFields[0] = Array<OneD, NekDouble>(ndim * ntotpts);
135 for (size_t i = 1; i < ndim; ++i)
136 {
137 m_previousVnFields[i] = m_previousVnFields[i - 1] + ntotpts;
138 }
139 }
140 else
141 {
142 ASSERTL0(0, "Integration method not suitable: Options include "
143 "BackwardEuler or BDFImplicitOrder{1,2}");
144 }
145
146 m_intSteps = IntegrationScheme->GetNumIntegrationPhases();
147
148 // set explicit time-integration class operators
153 timer.Stop();
154 timer.AccumulateRegion(
155 "SubSteppingExtrapolate:v_SubSteppingTimeIntegration", 10);
156}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
TimeIntegrationSchemeFactory & GetTimeIntegrationSchemeFactory()

References Nektar::LibUtilities::Timer::AccumulateRegion(), ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), Nektar::LibUtilities::GetTimeIntegrationSchemeFactory(), Nektar::Extrapolate::m_fields, m_intScheme, Nektar::Extrapolate::m_intSteps, m_previousVelFields, m_previousVnFields, Nektar::Extrapolate::m_session, m_subStepIntegrationOps, m_subStepIntegrationScheme, Nektar::Extrapolate::m_velocity, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), SubStepAdvection(), and SubStepProjection().

◆ v_SubStepSaveFields()

void Nektar::SubSteppingExtrapolate::v_SubStepSaveFields ( int  nstep)
overrideprotectedvirtual

Save the current fields in m_previousVelFields and cycle out older previous fields so that it can be extrapolated to the new substeps of the next time step. Also extract the normal trace of the velocity field into m_previousVnFields along the trace so that this can also be extrapolated.

Implements Nektar::Extrapolate.

Definition at line 334 of file SubSteppingExtrapolate.cpp.

335{
337 size_t i, n;
338 timer.Start();
339 size_t nvel = m_velocity.size();
340 size_t npts = m_fields[0]->GetTotPoints();
341 size_t ntracepts = m_fields[0]->GetTrace()->GetTotPoints();
342
343 // rotate fields
344 size_t nblocks = m_previousVelFields.size() / nvel;
345 Array<OneD, NekDouble> save;
346
347 // rotate storage space
348 for (n = 0; n < nvel; ++n)
349 {
350 save = m_previousVelFields[(nblocks - 1) * nvel + n];
351
352 for (i = nblocks - 1; i > 0; --i)
353 {
354 m_previousVelFields[i * nvel + n] =
355 m_previousVelFields[(i - 1) * nvel + n];
356 }
357
358 m_previousVelFields[n] = save;
359 }
360
361 save = m_previousVnFields[nblocks - 1];
362 for (i = nblocks - 1; i > 0; --i)
363 {
365 }
366 m_previousVnFields[0] = save;
367
368 // Put previous field
369 for (i = 0; i < nvel; ++i)
370 {
371 m_fields[m_velocity[i]]->BwdTrans(
372 m_fields[m_velocity[i]]->GetCoeffs(),
373 m_fields[m_velocity[i]]->UpdatePhys());
374 Vmath::Vcopy(npts, m_fields[m_velocity[i]]->GetPhys(), 1,
375 m_previousVelFields[i], 1);
376 }
377
378 Array<OneD, NekDouble> Fwd(ntracepts);
379
380 // calculate Vn
381 m_fields[0]->ExtractTracePhys(m_previousVelFields[0], Fwd);
382 Vmath::Vmul(ntracepts, m_traceNormals[0], 1, Fwd, 1, m_previousVnFields[0],
383 1);
384 for (i = 1; i < m_bnd_dim; ++i)
385 {
386 m_fields[0]->ExtractTracePhys(m_previousVelFields[i], Fwd);
387 Vmath::Vvtvp(ntracepts, m_traceNormals[i], 1, Fwd, 1,
389 }
390
391 if (nstep == 0) // initialise all levels with first field
392 {
393 for (n = 0; n < nvel; ++n)
394 {
395 for (i = 1; i < nblocks; ++i)
396 {
397 Vmath::Vcopy(npts, m_fields[m_velocity[n]]->GetPhys(), 1,
398 m_previousVelFields[i * nvel + n], 1);
399 }
400 }
401
402 for (i = 1; i < nblocks; ++i)
403 {
404 Vmath::Vcopy(ntracepts, m_previousVnFields[0], 1,
405 m_previousVnFields[i], 1);
406 }
407 }
408 timer.Stop();
409 timer.AccumulateRegion("SubSteppingExtrapolate:v_SubStepSaveFields", 10);
410}
int m_bnd_dim
bounday dimensionality
Definition: Extrapolate.h:220
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::Extrapolate::m_bnd_dim, Nektar::Extrapolate::m_fields, m_previousVelFields, m_previousVnFields, Nektar::Extrapolate::m_traceNormals, Nektar::Extrapolate::m_velocity, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_SubStepSetPressureBCs()

void Nektar::SubSteppingExtrapolate::v_SubStepSetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
NekDouble  Aii_Dt,
NekDouble  kinvis 
)
overrideprotectedvirtual

Implements Nektar::Extrapolate.

Reimplemented in Nektar::SubSteppingExtrapolateWeakPressure.

Definition at line 259 of file SubSteppingExtrapolate.cpp.

262{
264 // int nConvectiveFields =m_fields.size()-1;
265 Array<OneD, Array<OneD, NekDouble>> nullvelfields;
266
267 timer.Start();
269
270 // Calculate viscous BCs at current level and
271 // put in m_pressureHBCs[0]
272 CalcNeumannPressureBCs(m_previousVelFields, nullvelfields, kinvis);
273
274 // Extrapolate to m_pressureHBCs to n+1
276
277 // Add (phi,Du/Dt) term to m_presureHBC
278 AddDuDt();
279
280 // Copy m_pressureHBCs to m_PbndExp
282
283 // Evaluate High order outflow conditiosn if required.
284 CalcOutflowBCs(inarray, kinvis);
285 timer.Stop();
286 timer.AccumulateRegion("SubSteppingExtrapolate:SubStepSetPressureBCs", 10);
287}
Array< OneD, Array< OneD, NekDouble > > m_pressureHBCs
Storage for current and previous levels of high order pressure boundary conditions.
Definition: Extrapolate.h:247
void CopyPressureHBCsToPbndExp(void)
void CalcNeumannPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &fields, const Array< OneD, const Array< OneD, NekDouble > > &N, NekDouble kinvis)
Definition: Extrapolate.h:172
void ExtrapolateArray(Array< OneD, Array< OneD, NekDouble > > &array)
void CalcOutflowBCs(const Array< OneD, const Array< OneD, NekDouble > > &fields, NekDouble kinvis)

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::Extrapolate::AddDuDt(), Nektar::Extrapolate::CalcNeumannPressureBCs(), Nektar::Extrapolate::CalcOutflowBCs(), Nektar::Extrapolate::CopyPressureHBCsToPbndExp(), Nektar::Extrapolate::ExtrapolateArray(), Nektar::Extrapolate::m_pressureCalls, Nektar::Extrapolate::m_pressureHBCs, m_previousVelFields, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Member Data Documentation

◆ className

std::string Nektar::SubSteppingExtrapolate::className
static
Initial value:
=
"SubStepping", SubSteppingExtrapolate::create, "SubStepping")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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.
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:48

Name of class.

Registers the class with the Factory.

Definition at line 74 of file SubSteppingExtrapolate.h.

◆ m_cflSafetyFactor

NekDouble Nektar::SubSteppingExtrapolate::m_cflSafetyFactor
protected

◆ m_infosteps

int Nektar::SubSteppingExtrapolate::m_infosteps
protected

Definition at line 136 of file SubSteppingExtrapolate.h.

Referenced by SubSteppingExtrapolate(), and v_SubStepAdvance().

◆ m_intScheme

LibUtilities::TimeIntegrationSchemeSharedPtr Nektar::SubSteppingExtrapolate::m_intScheme
protected

Definition at line 128 of file SubSteppingExtrapolate.h.

Referenced by v_SubStepAdvance(), and v_SubSteppingTimeIntegration().

◆ m_maxsubsteps

int Nektar::SubSteppingExtrapolate::m_maxsubsteps
protected

Definition at line 138 of file SubSteppingExtrapolate.h.

Referenced by SubSteppingExtrapolate(), and v_SubStepAdvance().

◆ m_minsubsteps

int Nektar::SubSteppingExtrapolate::m_minsubsteps
protected

Definition at line 137 of file SubSteppingExtrapolate.h.

Referenced by SubSteppingExtrapolate(), and v_SubStepAdvance().

◆ m_previousVelFields

Array<OneD, Array<OneD, NekDouble> > Nektar::SubSteppingExtrapolate::m_previousVelFields
protected

◆ m_previousVnFields

Array<OneD, Array<OneD, NekDouble> > Nektar::SubSteppingExtrapolate::m_previousVnFields
protected

◆ m_subStepIntegrationOps

LibUtilities::TimeIntegrationSchemeOperators Nektar::SubSteppingExtrapolate::m_subStepIntegrationOps
protected

Definition at line 130 of file SubSteppingExtrapolate.h.

Referenced by v_SubStepAdvance(), and v_SubSteppingTimeIntegration().

◆ m_subStepIntegrationScheme

LibUtilities::TimeIntegrationSchemeSharedPtr Nektar::SubSteppingExtrapolate::m_subStepIntegrationScheme
protected