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

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...
 
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...
 
- 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_slices
 number of slices More...
 
NekDouble m_period
 period length More...
 
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 45 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 54 of file AdjointAdvection.cpp.

55 {
56 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 51 of file AdjointAdvection.h.

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

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

◆ 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 58 of file AdjointAdvection.cpp.

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().

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

Friends And Related Function Documentation

◆ MemoryManager< AdjointAdvection >

friend class MemoryManager< AdjointAdvection >
friend

Definition at line 48 of file AdjointAdvection.h.

Member Data Documentation

◆ className

string Nektar::AdjointAdvection::className
static
Initial value:

Name of class.

Definition at line 57 of file AdjointAdvection.h.