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

Advection for the adjoint form of the linearised Navier-Stokes equations. More...

#include <AdjointAdvection.h>

Inheritance diagram for Nektar::AdjointAdvection:
[legend]

Static Public Member Functions

static SolverUtils::AdvectionSharedPtr create (std::string)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::LinearisedAdvection
static SolverUtils::AdvectionSharedPtr create (std::string)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::LinearisedAdvection
static std::string className
 Name of class. More...
 

Protected Member Functions

 AdjointAdvection ()
 
virtual ~AdjointAdvection ()
 
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, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Advects a vector field. More...
 
- Protected Member Functions inherited from Nektar::LinearisedAdvection
DNekBlkMatSharedPtr GetFloquetBlockMatrix (FloquetMatType mattype, bool UseContCoeffs=false) const
 
DNekBlkMatSharedPtr GenFloquetBlockMatrix (FloquetMatType mattype, bool UseContCoeffs=false) const
 
 LinearisedAdvection ()
 
virtual ~LinearisedAdvection ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialises the advection object. More...
 
virtual void v_SetBaseFlow (const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 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 UpdateGradBase (const int var, const MultiRegions::ExpListSharedPtr &field)
 
void DFT (const std::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...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Advection
virtual SOLVER_UTILS_EXPORT void v_AdvectVolumeFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &advVel, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &time)
 Advects Volume Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_AdvectTraceFlux (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 >> &pTraceFlux, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Advects Trace Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_AdvectCoeffs (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, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 
virtual SOLVER_UTILS_EXPORT void v_AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 

Friends

class MemoryManager< AdjointAdvection >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Advection
virtual SOLVER_UTILS_EXPORT ~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, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Interface function to advect the vector field. More...
 
SOLVER_UTILS_EXPORT void AdvectVolumeFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &pTime)
 Interface function to advect the Volume field. More...
 
SOLVER_UTILS_EXPORT void AdvectTraceFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, Array< OneD, Array< OneD, NekDouble >> &pTraceFlux, const NekDouble &pTime, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Interface function to advect the Trace field. More...
 
SOLVER_UTILS_EXPORT void AdvectCoeffs (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, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Similar with Advection::Advect(): calculate the advection flux The difference is in the outarray: it is the coefficients of basis for AdvectCoeffs() it is the physics on quadrature points for Advect() 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, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 Set the base flow used for linearised advection objects. More...
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
SOLVER_UTILS_EXPORT void AddTraceJacToMat (const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType >> &TracePntJacGradSign)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
void CalcJacobTraceInteg (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int m, const int n, const Array< OneD, const TypeNekBlkMatSharedPtr > &PntJac, const Array< OneD, const Array< OneD, DataType >> &PntJacSign, Array< OneD, DNekMatSharedPtr > &TraceJacFwd, Array< OneD, DNekMatSharedPtr > &TraceJacBwd)
 
SOLVER_UTILS_EXPORT void AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
void AddTraceJacToMat (const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr >> &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType >> &TracePntJacGradSign)
 
- Protected Attributes inherited from Nektar::LinearisedAdvection
LibUtilities::SessionReaderSharedPtr m_session
 
int m_spacedim
 
int m_expdim
 
Array< OneD, Array< OneD, NekDouble > > m_baseflow
 Storage for base flow. More...
 
Array< OneD, Array< OneD, NekDouble > > m_gradBase
 
int m_start
 number of slices More...
 
int m_skip
 
int m_slices
 
NekDouble m_period
 period length More...
 
int m_isperiodic
 
int m_interporder
 
Array< OneD, Array< OneD, NekDouble > > m_interp
 interpolation vector More...
 
LibUtilities::NektarFFTSharedPtr m_FFT
 auxiliary variables More...
 
Array< OneD, NekDoublem_tmpIN
 
Array< OneD, NekDoublem_tmpOUT
 
bool m_useFFTW
 
bool m_singleMode
 flag to determine if use single mode or not More...
 
bool m_halfMode
 flag to determine if use half mode or not More...
 
bool m_multipleModes
 flag to determine if use multiple mode or not More...
 
FloquetBlockMatrixMapShPtr m_FloquetBlockMat
 
- Protected Attributes inherited from Nektar::SolverUtils::Advection
AdvectionFluxVecCB m_fluxVector
 Callback function to the flux vector (set when advection is in conservative form). More...
 
RiemannSolverSharedPtr m_riemann
 Riemann solver for DG-type schemes. More...
 
int m_spaceDim
 Storage for space dimension. Used for homogeneous extension. More...
 

Detailed Description

Advection for the adjoint form of the linearised Navier-Stokes equations.

Definition at line 44 of file AdjointAdvection.h.

Constructor & Destructor Documentation

◆ AdjointAdvection()

Nektar::AdjointAdvection::AdjointAdvection ( )
protected

Definition at line 49 of file AdjointAdvection.cpp.

◆ ~AdjointAdvection()

Nektar::AdjointAdvection::~AdjointAdvection ( )
protectedvirtual

Definition at line 53 of file AdjointAdvection.cpp.

54 {
55 }

Member Function Documentation

◆ create()

static SolverUtils::AdvectionSharedPtr Nektar::AdjointAdvection::create ( std::string  )
inlinestatic

Creates an instance of this class.

Definition at line 50 of file AdjointAdvection.h.

51  {
53  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ v_Advect()

void Nektar::AdjointAdvection::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,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd = NullNekDoubleArrayOfArray,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd = NullNekDoubleArrayOfArray 
)
protectedvirtual

Advects a vector field.

Reimplemented from Nektar::LinearisedAdvection.

Definition at line 57 of file AdjointAdvection.cpp.

65 {
66  ASSERTL1(nConvectiveFields == inarray.size(),
67  "Number of convective fields and Inarray are not compatible");
68 
69  int nPointsTot = fields[0]->GetNpoints();
70  int ndim = advVel.size();
71  int nBaseDerivs = (m_halfMode || m_singleMode) ? 2 : m_spacedim;
72  int nDerivs = (m_halfMode) ? 2 : m_spacedim;
73 
74  Array<OneD, Array<OneD, NekDouble>> velocity(ndim);
75  int nScalar = nConvectiveFields - ndim;
76  Array<OneD, Array<OneD, NekDouble>> scalar(nScalar);
77 
78  for (int i = 0; i < ndim; ++i)
79  {
80  if (fields[i]->GetWaveSpace() && !m_singleMode && !m_halfMode)
81  {
82  velocity[i] = Array<OneD, NekDouble>(nPointsTot, 0.0);
83  fields[i]->HomogeneousBwdTrans(advVel[i], velocity[i]);
84  }
85  else
86  {
87  velocity[i] = advVel[i];
88  }
89  }
90  if (nScalar > 0) // add for temperature field
91  {
92  for (int jj = ndim; jj < nConvectiveFields; ++jj)
93  {
94  scalar[jj - ndim] = inarray[jj];
95  }
96  }
97 
98  Array<OneD, Array<OneD, NekDouble>> grad(nDerivs);
99  for (int i = 0; i < nDerivs; ++i)
100  {
101  grad[i] = Array<OneD, NekDouble>(nPointsTot);
102  }
103 
104  // Evaluation of the base flow for periodic cases
105  if (m_slices > 1)
106  {
107  for (int i = 0; i < ndim; ++i)
108  {
110  m_period);
111  UpdateGradBase(i, fields[i]);
112  }
113  }
114 
115  // Evaluate the linearised advection term
116  for (int i = 0; i < ndim; ++i)
117  {
118  // Calculate gradient
119  switch (nDerivs)
120  {
121  case 1:
122  {
123  fields[i]->PhysDeriv(inarray[i], grad[0]);
124  }
125  break;
126  case 2:
127  {
128  fields[i]->PhysDeriv(inarray[i], grad[0], grad[1]);
129  }
130  break;
131  case 3:
132  {
133  fields[i]->PhysDeriv(inarray[i], grad[0], grad[1], grad[2]);
134  if (m_multipleModes)
135  {
136  // transform gradients into physical Fourier space
137  fields[i]->HomogeneousBwdTrans(grad[0], grad[0]);
138  fields[i]->HomogeneousBwdTrans(grad[1], grad[1]);
139  fields[i]->HomogeneousBwdTrans(grad[2], grad[2]);
140  }
141  }
142  break;
143  }
144 
145  // Calculate -U_j du'_i/dx_j
146  Vmath::Vmul(nPointsTot, grad[0], 1, m_baseflow[0], 1, outarray[i], 1);
147  for (int j = 1; j < nDerivs; ++j)
148  {
149  Vmath::Vvtvp(nPointsTot, grad[j], 1, m_baseflow[j], 1, outarray[i],
150  1, outarray[i], 1);
151  }
152  Vmath::Neg(nPointsTot, outarray[i], 1);
153 
154  // Add u'_j U_j/ dx_i
155  int lim = (m_halfMode) ? 2 : ndim;
156  if ((m_halfMode || m_singleMode) && i == 2)
157  {
158  lim = 0;
159  }
160  for (int j = 0; j < lim; ++j)
161  {
162  Vmath::Vvtvp(nPointsTot, m_gradBase[j * nBaseDerivs + i], 1,
163  velocity[j], 1, outarray[i], 1, outarray[i], 1);
164  }
165  // Add Tprime*Grad_Tbase in u, v equations
166  if (nScalar > 0 && i < ndim)
167  {
168  for (int s = 0; s < nScalar; ++s)
169  {
170  Vmath::Vvtvp(nPointsTot,
171  m_gradBase[(ndim + s) * nBaseDerivs + i], 1,
172  scalar[s], 1, outarray[i], 1, outarray[i], 1);
173  }
174  }
175 
176  if (m_multipleModes)
177  {
178  fields[i]->HomogeneousFwdTrans(outarray[i], outarray[i]);
179  }
180  Vmath::Neg(nPointsTot, outarray[i], 1);
181  }
182 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
NekDouble m_period
period length
void UpdateGradBase(const int var, const MultiRegions::ExpListSharedPtr &field)
bool m_singleMode
flag to determine if use single mode or not
Array< OneD, Array< OneD, NekDouble > > m_baseflow
Storage for base flow.
Array< OneD, Array< OneD, NekDouble > > m_gradBase
void UpdateBase(const NekDouble m_slices, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble m_time, const NekDouble m_period)
bool m_multipleModes
flag to determine if use multiple mode or not
bool m_halfMode
flag to determine if use half mode or not
Array< OneD, Array< OneD, NekDouble > > m_interp
interpolation vector
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.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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.cpp:574

References ASSERTL1, Nektar::LinearisedAdvection::m_baseflow, Nektar::LinearisedAdvection::m_gradBase, Nektar::LinearisedAdvection::m_halfMode, Nektar::LinearisedAdvection::m_interp, Nektar::LinearisedAdvection::m_multipleModes, Nektar::LinearisedAdvection::m_period, Nektar::LinearisedAdvection::m_singleMode, Nektar::LinearisedAdvection::m_slices, Nektar::LinearisedAdvection::m_spacedim, Vmath::Neg(), Nektar::LinearisedAdvection::UpdateBase(), Nektar::LinearisedAdvection::UpdateGradBase(), Vmath::Vmul(), and Vmath::Vvtvp().

Friends And Related Function Documentation

◆ MemoryManager< AdjointAdvection >

friend class MemoryManager< AdjointAdvection >
friend

Definition at line 1 of file AdjointAdvection.h.

Member Data Documentation

◆ className

string Nektar::AdjointAdvection::className
static
Initial value:
=
static SolverUtils::AdvectionSharedPtr create(std::string)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47

Name of class.

Definition at line 56 of file AdjointAdvection.h.