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

#include <Diffusion3DHomogeneous1D.h>

Inheritance diagram for Nektar::SolverUtils::Diffusion3DHomogeneous1D:
[legend]

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type []
 

Protected Member Functions

 Diffusion3DHomogeneous1D (std::string diffType)
 Diffusion3DHomogeneous1D uses the 2D WeakDG approach to compute the diffusion term looping on the planes in the z direction and adding the flux in z direction at the end. More...
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initiliase Diffusion3DHomogeneous1D objects and store them before starting the time-stepping. More...
 
virtual void v_Diffuse (const std::size_t nConvective, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 Calculate WeakDG Diffusion for the linear problems using an LDG interface flux and the the flux in the third direction. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT void v_DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &vFwd, const Array< OneD, Array< OneD, NekDouble > > &vBwd, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
 
virtual SOLVER_UTILS_EXPORT void v_ConsVarAveJump (const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble > > &vFwd, const Array< OneD, const Array< OneD, NekDouble > > &vBwd, Array< OneD, Array< OneD, NekDouble > > &aver, Array< OneD, Array< OneD, NekDouble > > &jump)
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseCalcDerivative (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Diffusion Flux, calculate the physical derivatives. More...
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseVolumeFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex)
 Diffusion Volume Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &Qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex)
 Diffusion term Trace Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_DiffuseTraceSymmFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump)
 
virtual SOLVER_UTILS_EXPORT void v_AddDiffusionSymmFluxToCoeff (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 
virtual SOLVER_UTILS_EXPORT void v_AddDiffusionSymmFluxToPhys (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual TensorOfArray3D< NekDouble > & v_GetFluxTensor ()
 
virtual SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & v_GetTraceNormal ()
 
virtual SOLVER_UTILS_EXPORT void v_GetPrimVar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &primVar)
 Compute primary derivatives. More...
 
SOLVER_UTILS_EXPORT void GetDivCurl (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &pVarDer)
 Compute divergence and curl squared. More...
 

Protected Attributes

LibUtilities::TranspositionSharedPtr m_trans
 
std::string m_diffType
 
SolverUtils::DiffusionSharedPtr m_planeDiff
 
NekDouble m_homoLen
 
std::size_t m_numPoints
 
std::size_t m_numPointsPlane
 
std::size_t m_numPlanes
 
std::size_t m_planeCounter
 
Array< OneD, unsigned int > m_planes
 
Array< OneD, unsigned int > m_planePos
 
Array< OneD, Array< OneD, NekDouble > > m_homoDerivStore
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_homoDerivPlane
 
Array< OneD, Array< OneD, NekDouble > > m_inarrayPlane
 
Array< OneD, Array< OneD, NekDouble > > m_outarrayPlane
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsPlane
 
Array< OneD, Array< OneD, NekDouble > > m_advVelPlane
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
DiffusionFluxVecCB m_fluxVector
 
DiffusionFluxVecCBNS m_fluxVectorNS
 
DiffusionFluxPenaltyNS m_fluxPenaltyNS
 
DiffusionArtificialDiffusion m_ArtificialDiffusionVector
 
DiffusionFluxCons m_FunctorDiffusionfluxCons
 
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
 
SpecialBndTreat m_SpecialBndTreat
 
CalcViscosity m_CalcViscosity
 
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
 
NekDouble m_time =0.0
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT ~Diffusion ()
 
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 
SOLVER_UTILS_EXPORT void Diffuse (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 Similar with Diffusion::Diffuse(): calculate diffusion flux The difference is in the outarray: it is the coefficients of basis for DiffuseCoeffs() it is the physics on quadrature points for Diffuse() More...
 
SOLVER_UTILS_EXPORT void Diffuse (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, TensorOfArray3D< NekDouble > &qfield, Array< OneD, int > &nonZeroIndex)
 
SOLVER_UTILS_EXPORT void DiffuseCoeffs (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
 Diffusion Volume FLux. More...
 
SOLVER_UTILS_EXPORT void DiffuseTraceFlux (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
 Diffusion term Trace Flux. More...
 
SOLVER_UTILS_EXPORT void DiffuseTraceSymmFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfield, const TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble > > &solution_Aver=NullNekDoubleArrayOfArray, Array< OneD, Array< OneD, NekDouble > > &solution_jump=NullNekDoubleArrayOfArray)
 
SOLVER_UTILS_EXPORT void AddDiffusionSymmFluxToCoeff (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Add symmetric flux to field in coeff space. More...
 
SOLVER_UTILS_EXPORT void AddDiffusionSymmFluxToPhys (const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Add symmetric flux to field in coeff physical space. More...
 
SOLVER_UTILS_EXPORT void AddSymmFluxIntegralToOffDiag (const int nConvectiveFields, const int nDim, const int nPts, const int nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &Fwdflux, TensorOfArray3D< NekDouble > &Bwdflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void FluxVec (TensorOfArray3D< NekDouble > &fluxvector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVector (FuncPointerT func, ObjectPointerT obj)
 
SOLVER_UTILS_EXPORT void SetFluxVector (DiffusionFluxVecCB fluxVector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVectorNS (FuncPointerT func, ObjectPointerT obj)
 
void SetFluxVectorNS (DiffusionFluxVecCBNS fluxVector)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxPenaltyNS (FuncPointerT func, ObjectPointerT obj)
 
void SetFluxPenaltyNS (DiffusionFluxPenaltyNS flux)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionFluxCons (FuncPointerT func, ObjectPointerT obj)
 
void SetDiffusionFluxCons (DiffusionFluxCons flux)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionFluxConsTrace (FuncPointerT func, ObjectPointerT obj)
 
void SetDiffusionFluxConsTrace (DiffusionFluxCons flux)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetArtificialDiffusionVector (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetSpecialBndTreat (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetCalcViscosity (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionSymmFluxCons (FuncPointerT func, ObjectPointerT obj)
 
void SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
virtual TensorOfArray3D< NekDouble > & GetFluxTensor ()
 
SOLVER_UTILS_EXPORT void GetAVmu (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &muvar, Array< OneD, NekDouble > &MuVarTrace)
 Get the mu of artifical viscosity(AV) More...
 
SOLVER_UTILS_EXPORT void ConsVarAveJump (const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble > > &vFwd, const Array< OneD, const Array< OneD, NekDouble > > &vBwd, Array< OneD, Array< OneD, NekDouble > > &aver, Array< OneD, Array< OneD, NekDouble > > &jump)
 Get the average and jump value of conservative variables on trace. More...
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal ()
 Get trace normal. More...
 
- Public Attributes inherited from Nektar::SolverUtils::Diffusion
Array< OneD, NekDoublem_divVel
 Params for Ducros sensor. More...
 
Array< OneD, NekDoublem_divVelSquare
 
Array< OneD, NekDoublem_curlVelSquare
 

Detailed Description

Definition at line 45 of file Diffusion3DHomogeneous1D.h.

Constructor & Destructor Documentation

◆ Diffusion3DHomogeneous1D()

Nektar::SolverUtils::Diffusion3DHomogeneous1D::Diffusion3DHomogeneous1D ( std::string  diffType)
protected

Diffusion3DHomogeneous1D uses the 2D WeakDG approach to compute the diffusion term looping on the planes in the z direction and adding the flux in z direction at the end.

Definition at line 79 of file Diffusion3DHomogeneous1D.cpp.

80  {
81  // Strip trailing string "3DHomogeneous1D" to determine 2D diffusion
82  // type, and create a diffusion object for the plane.
83  m_diffType = diffType.substr(0, diffType.length()-15);
85  }
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), and Nektar::SolverUtils::GetDiffusionFactory().

Referenced by create().

Member Function Documentation

◆ create()

static DiffusionSharedPtr Nektar::SolverUtils::Diffusion3DHomogeneous1D::create ( std::string  diffType)
inlinestatic

Definition at line 48 of file Diffusion3DHomogeneous1D.h.

49  {
50  return DiffusionSharedPtr(
51  new Diffusion3DHomogeneous1D(diffType));
52  }
Diffusion3DHomogeneous1D(std::string diffType)
Diffusion3DHomogeneous1D uses the 2D WeakDG approach to compute the diffusion term looping on the pla...
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:627

References Diffusion3DHomogeneous1D().

◆ v_Diffuse()

void Nektar::SolverUtils::Diffusion3DHomogeneous1D::v_Diffuse ( const std::size_t  nConvective,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd = NullNekDoubleArrayOfArray,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd = NullNekDoubleArrayOfArray 
)
protectedvirtual

Calculate WeakDG Diffusion for the linear problems using an LDG interface flux and the the flux in the third direction.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 190 of file Diffusion3DHomogeneous1D.cpp.

197  {
198  boost::ignore_unused(pFwd, pBwd);
199 
200  Array<OneD, NekDouble> tmp(m_numPoints), tmp2;
201  Array<OneD, Array<OneD, NekDouble> > viscHComp;
202  const int nPointsTot = fields[0]->GetNpoints();
203  NekDouble beta;
204 
205 
206  if (m_fluxVectorNS)
207  {
208  viscHComp = Array<OneD, Array<OneD, NekDouble> >(nConvectiveFields);
209  for (int i = 0; i < nConvectiveFields - 1; ++i)
210  {
211  fields[0]->PhysDeriv(2, inarray[i], m_homoDerivStore[i]);
212  viscHComp[i] = Array<OneD, NekDouble>(m_numPoints);
213  }
214  }
215 
216 
217  for (int i = 0; i < m_numPlanes; ++i)
218  {
219  // Set up memory references for fields, inarray and outarray for
220  // this plane.
221  for (int j = 0; j < inarray.size(); ++j)
222  {
223  m_inarrayPlane [j] = Array<OneD, NekDouble>(
224  m_numPointsPlane, tmp2 = inarray [j] + m_planePos[i]);
225  }
226 
227  for (int j = 0; j < nConvectiveFields; ++j)
228  {
229  m_fieldsPlane [j] = fields[j]->GetPlane(i);
230  m_outarrayPlane[j] = Array<OneD, NekDouble>(
231  m_numPointsPlane, tmp2 = outarray[j] + m_planePos[i]);
232  }
233 
234 
235  if (m_fluxVectorNS)
236  {
237  m_planeDiff->SetHomoDerivs(m_homoDerivPlane[i]);
238  }
239 
240 
241  if (m_diffType == "LDGNS")
242  {
243  // Store plane Fwd/Bwd traces
244  std::size_t nTracePts = m_fieldsPlane[0]->GetTrace()
245  ->GetTotPoints();
246  std::size_t nScalar = m_inarrayPlane.size();
247  Array<OneD, Array<OneD, NekDouble> > Fwd(nScalar);
248  Array<OneD, Array<OneD, NekDouble> > Bwd(nScalar);
249  {
250  for(std::size_t k = 0; k < nScalar; ++k)
251  {
252  Fwd[k] = Array<OneD, NekDouble>(nTracePts, 0.0);
253  Bwd[k] = Array<OneD, NekDouble>(nTracePts, 0.0);
254  m_fieldsPlane[k]->GetFwdBwdTracePhys(
255  m_inarrayPlane[k], Fwd[k], Bwd[k]);
256  }
257  }
258 
259  m_planeDiff->Diffuse(nConvectiveFields,
263  Fwd,
264  Bwd);
265  }
266  else
267  {
268  m_planeDiff->Diffuse(nConvectiveFields,
272  }
273 
274  if (m_fluxVectorNS)
275  {
276  Array<OneD, Array<OneD, Array<OneD, NekDouble> > >
277  &viscTensor = m_planeDiff->GetFluxTensor();
278 
279  // Extract H (viscTensor[2])
280  for (int j = 0; j < nConvectiveFields - 1; ++j)
281  {
283  viscTensor[2][j+1], 1,
284  tmp2 = viscHComp[j] + m_planePos[i], 1);
285  }
286  }
287  }
288 
289 
290 
291  if (m_fluxVectorNS)
292  {
293  for (int j = 0; j < nConvectiveFields - 1; ++j)
294  {
295  fields[j+1]->PhysDeriv(2, viscHComp[j], tmp);
296  Vmath::Vadd(nPointsTot, outarray[j+1], 1, tmp, 1,
297  outarray[j+1], 1);
298  }
299  }
300  else
301  {
302  for (int j = 0; j < nConvectiveFields; ++j)
303  {
304  fields[j]->HomogeneousFwdTrans(inarray[j], tmp);
305 
306  for (int i = 0; i < m_numPlanes; ++i)
307  {
308  beta = 2*M_PI*m_trans->GetK(i)/m_homoLen;
309  beta *= beta;
310 
312  beta,
313  &tmp[0] + i*m_numPointsPlane, 1,
314  &tmp[0] + i*m_numPointsPlane, 1);
315  }
316 
317  fields[0]->HomogeneousBwdTrans(tmp, tmp);
318 
319  Vmath::Vsub(nPointsTot, outarray[j], 1, tmp, 1,
320  outarray[j], 1);
321  }
322  }
323  }
LibUtilities::TranspositionSharedPtr m_trans
Array< OneD, MultiRegions::ExpListSharedPtr > m_fieldsPlane
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_homoDerivPlane
Array< OneD, Array< OneD, NekDouble > > m_outarrayPlane
Array< OneD, Array< OneD, NekDouble > > m_homoDerivStore
Array< OneD, Array< OneD, NekDouble > > m_inarrayPlane
DiffusionFluxVecCBNS m_fluxVectorNS
Definition: Diffusion.h:467
double NekDouble
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:322
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.cpp:225
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
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.cpp:372

References Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vsub().

◆ v_InitObject()

void Nektar::SolverUtils::Diffusion3DHomogeneous1D::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Initiliase Diffusion3DHomogeneous1D objects and store them before starting the time-stepping.

Parameters
pSessionPointer to session reader.
pFieldsPointer to fields.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 94 of file Diffusion3DHomogeneous1D.cpp.

97  {
98  int nConvectiveFields = pFields.size();
99 
100  Array<OneD, MultiRegions::ExpListSharedPtr> pFields_plane0(
101  nConvectiveFields);
102 
103  // Initialise the plane advection object.
104  for (int i = 0; i < nConvectiveFields; ++i)
105  {
106  pFields_plane0[i] = pFields[i]->GetPlane(0);
107  }
108  m_planeDiff->InitObject(pSession, pFields_plane0);
109 
110  m_numPoints = pFields[0]->GetTotPoints();
111  m_planes = pFields[0]->GetZIDs();
112  m_numPlanes = m_planes.size();
114  m_homoLen = pFields[0]->GetHomoLen();
115  m_trans = pFields[0]->GetTransposition();
116  m_planeCounter = 0;
117 
118  if (m_diffType == "LDG")
119  {
120  // Set viscous flux for LDG
121  m_planeDiff->SetFluxVector(m_fluxVector);
122  }
123  else if (m_diffType == "LDGNS")
124  {
125  // Set viscous flux for LDGNS
126  m_planeDiff->SetFluxVectorNS(m_fluxVectorNS);
127  // Set penalty flux
128  m_planeDiff->SetFluxPenaltyNS(m_fluxPenaltyNS);
129  }
130  else if (m_diffType == "LFRDGNS" ||
131  m_diffType == "LFRHUNS" ||
132  m_diffType == "LFRSDNS" )
133  {
134  // Set viscous flux for FR cases
135  m_planeDiff->SetFluxVectorNS(m_fluxVectorNS);
136  }
137 
138  m_fieldsPlane = Array<OneD, MultiRegions::ExpListSharedPtr>
139  (nConvectiveFields);
140 
141  if (m_fluxVectorNS)
142  {
143  m_inarrayPlane = Array<OneD, Array<OneD, NekDouble> >
144  (nConvectiveFields - 1);
145  }
146  else
147  {
148  m_inarrayPlane = Array<OneD, Array<OneD, NekDouble> >
149  (nConvectiveFields);
150  }
151  m_outarrayPlane = Array<OneD, Array<OneD, NekDouble> >
152  (nConvectiveFields);
153  m_planePos = Array<OneD, unsigned int> (m_numPlanes);
154 
155  for (int i = 0; i < m_numPlanes; ++i)
156  {
157  m_planePos[i] = i * m_numPointsPlane;
158  }
159 
160  if (m_fluxVectorNS)
161  {
162  m_homoDerivStore = Array<OneD, Array<OneD, NekDouble> >(
163  nConvectiveFields);
164  m_homoDerivPlane = Array<OneD, Array<OneD, Array<OneD, NekDouble> > >(
165  m_numPlanes);
166 
167  for (int i = 0; i < nConvectiveFields; ++i)
168  {
169  m_homoDerivStore[i] = Array<OneD, NekDouble>(m_numPoints);
170  }
171 
172  for (int i = 0; i < m_numPlanes; ++i)
173  {
174  m_homoDerivPlane[i] = Array<OneD, Array<OneD, NekDouble> >(nConvectiveFields);
175 
176  for (int j = 0; j < nConvectiveFields; ++j)
177  {
178  m_homoDerivPlane[i][j] = Array<OneD, NekDouble>(
180  m_homoDerivStore[j] + m_planePos[i]);
181  }
182  }
183  }
184  }
DiffusionFluxPenaltyNS m_fluxPenaltyNS
Definition: Diffusion.h:468
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:466

Member Data Documentation

◆ m_advVelPlane

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_advVelPlane
protected

Definition at line 73 of file Diffusion3DHomogeneous1D.h.

◆ m_diffType

std::string Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_diffType
protected

Definition at line 59 of file Diffusion3DHomogeneous1D.h.

◆ m_fieldsPlane

Array<OneD, MultiRegions::ExpListSharedPtr> Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_fieldsPlane
protected

Definition at line 72 of file Diffusion3DHomogeneous1D.h.

◆ m_homoDerivPlane

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_homoDerivPlane
protected

Definition at line 69 of file Diffusion3DHomogeneous1D.h.

◆ m_homoDerivStore

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_homoDerivStore
protected

Definition at line 68 of file Diffusion3DHomogeneous1D.h.

◆ m_homoLen

NekDouble Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_homoLen
protected

Definition at line 61 of file Diffusion3DHomogeneous1D.h.

◆ m_inarrayPlane

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_inarrayPlane
protected

Definition at line 70 of file Diffusion3DHomogeneous1D.h.

◆ m_numPlanes

std::size_t Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_numPlanes
protected

Definition at line 64 of file Diffusion3DHomogeneous1D.h.

◆ m_numPoints

std::size_t Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_numPoints
protected

Definition at line 62 of file Diffusion3DHomogeneous1D.h.

◆ m_numPointsPlane

std::size_t Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_numPointsPlane
protected

Definition at line 63 of file Diffusion3DHomogeneous1D.h.

◆ m_outarrayPlane

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_outarrayPlane
protected

Definition at line 71 of file Diffusion3DHomogeneous1D.h.

◆ m_planeCounter

std::size_t Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_planeCounter
protected

Definition at line 65 of file Diffusion3DHomogeneous1D.h.

◆ m_planeDiff

SolverUtils::DiffusionSharedPtr Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_planeDiff
protected

Definition at line 60 of file Diffusion3DHomogeneous1D.h.

◆ m_planePos

Array<OneD, unsigned int> Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_planePos
protected

Definition at line 67 of file Diffusion3DHomogeneous1D.h.

◆ m_planes

Array<OneD, unsigned int> Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_planes
protected

Definition at line 66 of file Diffusion3DHomogeneous1D.h.

◆ m_trans

LibUtilities::TranspositionSharedPtr Nektar::SolverUtils::Diffusion3DHomogeneous1D::m_trans
protected

Definition at line 58 of file Diffusion3DHomogeneous1D.h.

◆ type

std::string Nektar::SolverUtils::Diffusion3DHomogeneous1D::type
static
Initial value:
= {
"LDG3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRDG3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRSD3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRHU3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRcmin3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRcinf3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LDGNS3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRDGNS3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRSDNS3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRHUNS3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRcminNS3DHomogeneous1D", Diffusion3DHomogeneous1D::create),
"LFRcinfNS3DHomogeneous1D", Diffusion3DHomogeneous1D::create)
}
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
static DiffusionSharedPtr create(std::string diffType)

Definition at line 53 of file Diffusion3DHomogeneous1D.h.