Nektar++
|
Advection for the adjoint form of the linearised Navier-Stokes equations. More...
#include <AdjointAdvection.h>
Static Public Member Functions | |
static SolverUtils::AdvectionSharedPtr | create (std::string) |
Creates an instance of this class. More... | |
Static Public Attributes | |
static std::string | className |
Name of class. More... | |
Protected Member Functions | |
DNekBlkMatSharedPtr | GetFloquetBlockMatrix (FloquetMatType mattype, bool UseContCoeffs=false) const |
DNekBlkMatSharedPtr | GenFloquetBlockMatrix (FloquetMatType mattype, bool UseContCoeffs=false) const |
AdjointAdvection () | |
virtual | ~AdjointAdvection () |
virtual void | v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) |
Initialises the advection object. More... | |
virtual void | v_Advect (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) |
Advects a vector field. More... | |
virtual void | v_SetBaseFlow (const Array< OneD, Array< OneD, NekDouble > > &inarray) |
Overrides the base flow used during linearised advection. More... | |
void | UpdateBase (const NekDouble m_slices, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble m_time, const NekDouble m_period) |
void | DFT (const string file, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const NekDouble m_slices) |
void | ImportFldBase (std::string pInfile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, int slice) |
Import Base flow. More... | |
Private Types | |
enum | FloquetMatType { eForwardsCoeff, eForwardsPhys } |
enum | HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous } |
Parameter for homogeneous expansions. More... | |
typedef map< FloquetMatType, DNekBlkMatSharedPtr > | FloquetBlockMatrixMap |
A map between matrix keys and their associated block matrices. More... | |
typedef boost::shared_ptr < FloquetBlockMatrixMap > | FloquetBlockMatrixMapShPtr |
A shared pointer to a BlockMatrixMap. More... | |
Private Attributes | |
bool | m_useFFT |
flag to determine if use or not the FFT for transformations More... | |
enum HomogeneousType | m_HomogeneousType |
NekDouble | m_LhomX |
physical length in X direction (if homogeneous) More... | |
NekDouble | m_LhomY |
physical length in Y direction (if homogeneous) More... | |
NekDouble | m_LhomZ |
physical length in Z direction (if homogeneous) More... | |
int | m_npointsX |
number of points in X direction (if homogeneous) More... | |
int | m_npointsY |
number of points in Y direction (if homogeneous) More... | |
int | m_npointsZ |
number of points in Z direction (if homogeneous) More... | |
int | m_HomoDirec |
number of homogenous directions More... | |
int | m_NumMode |
Mode to use in case of single mode analysis. More... | |
SpatialDomains::BoundaryConditionsSharedPtr | m_boundaryConditions |
Friends | |
class | MemoryManager< AdjointAdvection > |
Additional Inherited Members | |
Public Member Functions inherited from Nektar::SolverUtils::Advection | |
SOLVER_UTILS_EXPORT void | InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) |
Interface function to initialise the advection object. More... | |
SOLVER_UTILS_EXPORT void | Advect (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) |
Interface function to advect the vector field. More... | |
template<typename FuncPointerT , typename ObjectPointerT > | |
void | SetFluxVector (FuncPointerT func, ObjectPointerT obj) |
Set the flux vector callback function. More... | |
void | SetRiemannSolver (RiemannSolverSharedPtr riemann) |
Set a Riemann solver object for this advection object. More... | |
void | SetFluxVector (AdvectionFluxVecCB fluxVector) |
Set the flux vector callback function. More... | |
void | SetBaseFlow (const Array< OneD, Array< OneD, NekDouble > > &inarray) |
Set the base flow used for linearised advection objects. More... | |
Advection for the adjoint form of the linearised Navier-Stokes equations.
Definition at line 47 of file AdjointAdvection.h.
|
private |
A map between matrix keys and their associated block matrices.
Definition at line 56 of file AdjointAdvection.h.
|
private |
A shared pointer to a BlockMatrixMap.
Definition at line 58 of file AdjointAdvection.h.
|
private |
|
private |
Parameter for homogeneous expansions.
Enumerator | |
---|---|
eHomogeneous1D | |
eHomogeneous2D | |
eHomogeneous3D | |
eNotHomogeneous |
Definition at line 153 of file AdjointAdvection.h.
|
protected |
Definition at line 60 of file AdjointAdvection.cpp.
|
protectedvirtual |
Definition at line 258 of file AdjointAdvection.cpp.
|
inlinestatic |
Creates an instance of this class.
Definition at line 64 of file AdjointAdvection.h.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr().
|
protected |
Definition at line 666 of file AdjointAdvection.cpp.
References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), eForwardsPhys, Nektar::eWrapper, GetFloquetBlockMatrix(), Nektar::LibUtilities::GetNektarFFTFactory(), ImportFldBase(), m_baseflow, m_FFT, m_interp, m_period, m_session, m_slices, m_tmpIN, m_tmpOUT, Vmath::Smul(), Vmath::Vcopy(), and Vmath::Zero().
Referenced by v_InitObject().
|
protected |
|
protected |
Definition at line 627 of file AdjointAdvection.cpp.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::eDIAGONAL, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::StdRegions::eFwdTrans, Nektar::StdRegions::StdExpansion::GetStdMatrix(), m_baseflow, and m_slices.
Referenced by DFT().
|
protected |
Import Base flow.
Import field from infile and load into m_fields. This routine will also perform a BwdTrans to ensure data is in both the physical and coefficient storage.
infile | Filename to read. |
Definition at line 497 of file AdjointAdvection.cpp.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, m_baseflow, m_HalfMode, m_interp, m_npointsZ, m_session, m_SingleMode, Vmath::Vcopy(), and Vmath::Zero().
Referenced by DFT(), and v_InitObject().
|
protected |
Definition at line 600 of file AdjointAdvection.cpp.
References m_baseflow, m_period, m_slices, Vmath::Svtvp(), and Vmath::Vcopy().
Referenced by v_Advect().
|
protectedvirtual |
Advects a vector field.
Implements Nektar::SolverUtils::Advection.
Definition at line 263 of file AdjointAdvection.cpp.
References ASSERTL0, ASSERTL1, Nektar::LibUtilities::eFunctionTypeFile, m_baseflow, m_HalfMode, m_interp, m_period, m_session, m_SingleMode, m_slices, Vmath::Neg(), UpdateBase(), Vmath::Vmul(), and Vmath::Vvtvp().
|
protectedvirtual |
Initialises the advection object.
This function should be overridden in derived classes to initialise the specific advection data members. However, this base class function should be called as the first statement of the overridden function to ensure the base class is correctly initialised in order.
pSession | Session information. |
pFields | Expansion lists for scalar fields. |
Reimplemented from Nektar::SolverUtils::Advection.
Definition at line 66 of file AdjointAdvection.cpp.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, DFT(), Nektar::MultiRegions::eDiscontinuous, Nektar::LibUtilities::eFunctionTypeFile, Nektar::MultiRegions::eGalerkin, eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, Nektar::MultiRegions::eLocal, eNotHomogeneous, ImportFldBase(), m_baseflow, m_boundaryConditions, m_CoeffState, m_expdim, m_HalfMode, m_HomoDirec, m_homogen_dealiasing, m_HomogeneousType, m_LhomX, m_LhomY, m_LhomZ, m_MultipleModes, m_npointsX, m_npointsY, m_npointsZ, m_projectionType, m_session, m_SingleMode, m_slices, m_spacedim, and m_useFFT.
|
protectedvirtual |
Overrides the base flow used during linearised advection.
Reimplemented from Nektar::SolverUtils::Advection.
Definition at line 474 of file AdjointAdvection.cpp.
References ASSERTL1, m_baseflow, npts, and Vmath::Vcopy().
|
friend |
Definition at line 61 of file AdjointAdvection.h.
|
static |
Name of class.
Definition at line 69 of file AdjointAdvection.h.
Storage for base flow.
Definition at line 79 of file AdjointAdvection.h.
Referenced by DFT(), GetFloquetBlockMatrix(), ImportFldBase(), UpdateBase(), v_Advect(), v_InitObject(), and v_SetBaseFlow().
|
private |
Definition at line 178 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
protected |
Definition at line 100 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
protected |
Definition at line 76 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
protected |
Definition at line 88 of file AdjointAdvection.h.
Referenced by DFT().
|
protected |
Definition at line 110 of file AdjointAdvection.h.
|
protected |
flag to determine if use half mode or not
Definition at line 96 of file AdjointAdvection.h.
Referenced by ImportFldBase(), v_Advect(), and v_InitObject().
|
private |
number of homogenous directions
Definition at line 174 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
protected |
Definition at line 99 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
private |
Definition at line 164 of file AdjointAdvection.h.
Referenced by v_InitObject().
Definition at line 86 of file AdjointAdvection.h.
Referenced by DFT(), ImportFldBase(), and v_Advect().
|
private |
physical length in X direction (if homogeneous)
Definition at line 166 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
private |
physical length in Y direction (if homogeneous)
Definition at line 167 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
private |
physical length in Z direction (if homogeneous)
Definition at line 168 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
protected |
flag to determine if use multiple mode or not
Definition at line 98 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
private |
number of points in X direction (if homogeneous)
Definition at line 170 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
private |
number of points in Y direction (if homogeneous)
Definition at line 171 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
private |
number of points in Z direction (if homogeneous)
Definition at line 172 of file AdjointAdvection.h.
Referenced by ImportFldBase(), and v_InitObject().
|
private |
Mode to use in case of single mode analysis.
Definition at line 176 of file AdjointAdvection.h.
|
protected |
Definition at line 84 of file AdjointAdvection.h.
Referenced by DFT(), UpdateBase(), and v_Advect().
|
protected |
Definition at line 74 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
protected |
Definition at line 72 of file AdjointAdvection.h.
Referenced by DFT(), ImportFldBase(), v_Advect(), and v_InitObject().
|
protected |
flag to determine if use single mode or not
Definition at line 94 of file AdjointAdvection.h.
Referenced by ImportFldBase(), v_Advect(), and v_InitObject().
|
protected |
Definition at line 82 of file AdjointAdvection.h.
Referenced by DFT(), GetFloquetBlockMatrix(), UpdateBase(), v_Advect(), and v_InitObject().
|
protected |
Definition at line 75 of file AdjointAdvection.h.
Referenced by v_InitObject().
Definition at line 89 of file AdjointAdvection.h.
Referenced by DFT().
Definition at line 90 of file AdjointAdvection.h.
Referenced by DFT().
|
private |
flag to determine if use or not the FFT for transformations
Definition at line 162 of file AdjointAdvection.h.
Referenced by v_InitObject().
|
protected |
Definition at line 91 of file AdjointAdvection.h.