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

#include <DiffusionLDG.h>

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

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 DiffusionLDG ()
 
void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 
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, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
 
void v_DiffuseCoeffs (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, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
 
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) override
 Diffusion Flux, calculate the physical derivatives. More...
 
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) override
 Diffusion Volume Flux. More...
 
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) override
 Diffusion term Trace Flux. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual SOLVER_UTILS_EXPORT void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)=0
 
virtual SOLVER_UTILS_EXPORT void v_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, const Array< OneD, Array< OneD, NekDouble > > &pBwd)=0
 
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_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 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 ()
 

Private Member Functions

void NumFluxforScalar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &uflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 
void ApplyScalarBCs (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const Array< OneD, const NekDouble > &ufield, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &penaltyflux)
 
void NumFluxforVector (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to penalty term. More...
 
void ApplyVectorBCs (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
 

Private Attributes

std::string m_shockCaptureType
 
NekDouble m_C11
 Coefficient of penalty term. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 
LibUtilities::SessionReaderSharedPtr m_session
 

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 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=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, 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 SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
 
SOLVER_UTILS_EXPORT TensorOfArray3D< NekDouble > & GetFluxTensor ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetTraceNormal ()
 Get trace normal. More...
 
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 SetSpecialBndTreat (FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetDiffusionSymmFluxCons (FuncPointerT func, ObjectPointerT obj)
 
SOLVER_UTILS_EXPORT void SetGridVelocityTrace (Array< OneD, Array< OneD, NekDouble > > &gridVelocityTrace)
 
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
Array< OneD, NekDoublem_divVel
 Params for Ducros sensor. More...
 
Array< OneD, NekDoublem_divVelSquare
 
Array< OneD, NekDoublem_curlVelSquare
 
DiffusionFluxVecCB m_fluxVector
 
DiffusionFluxVecCBNS m_fluxVectorNS
 
DiffusionFluxPenaltyNS m_fluxPenaltyNS
 
DiffusionFluxCons m_FunctorDiffusionfluxCons
 
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
 
SpecialBndTreat m_SpecialBndTreat
 
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
NekDouble m_time = 0.0
 

Detailed Description

Definition at line 42 of file DiffusionLDG.h.

Constructor & Destructor Documentation

◆ DiffusionLDG()

Nektar::SolverUtils::DiffusionLDG::DiffusionLDG ( )
protected

Definition at line 47 of file DiffusionLDG.cpp.

48{
49}

Referenced by create().

Member Function Documentation

◆ ApplyScalarBCs()

void Nektar::SolverUtils::DiffusionLDG::ApplyScalarBCs ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const std::size_t  var,
const Array< OneD, const NekDouble > &  ufield,
const Array< OneD, const NekDouble > &  Fwd,
const Array< OneD, const NekDouble > &  Bwd,
Array< OneD, NekDouble > &  penaltyflux 
)
private

Definition at line 284 of file DiffusionLDG.cpp.

291{
292 // Number of boundary regions
293 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
294 std::size_t cnt = 0;
295 for (std::size_t i = 0; i < nBndRegions; ++i)
296 {
297 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
299 {
300 continue;
301 }
302
303 // Number of boundary expansion related to that region
304 std::size_t nBndEdges =
305 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
306
307 // Weakly impose boundary conditions by modifying flux values
308 for (std::size_t e = 0; e < nBndEdges; ++e)
309 {
310 std::size_t nBndEdgePts = fields[var]
311 ->GetBndCondExpansions()[i]
312 ->GetExp(e)
313 ->GetTotPoints();
314
315 std::size_t id1 =
316 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
317
318 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
319 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
320
321 // AV boundary conditions
322 if (boost::iequals(
323 fields[var]->GetBndConditions()[i]->GetUserDefined(),
324 "Wall") ||
325 boost::iequals(
326 fields[var]->GetBndConditions()[i]->GetUserDefined(),
327 "Symmetry") ||
328 boost::iequals(
329 fields[var]->GetBndConditions()[i]->GetUserDefined(),
330 "WallViscous") ||
331 boost::iequals(
332 fields[var]->GetBndConditions()[i]->GetUserDefined(),
333 "WallAdiabatic") ||
334 boost::iequals(
335 fields[var]->GetBndConditions()[i]->GetUserDefined(),
336 "WallRotational"))
337 {
338 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
339 }
340 // For Dirichlet boundary condition: uflux = g_D
341 else if (fields[var]
342 ->GetBndConditions()[i]
343 ->GetBoundaryConditionType() ==
345 {
347 nBndEdgePts,
348 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
349 1, &penaltyflux[id2], 1);
350 }
351 // For Neumann boundary condition: uflux = u+
352 else if ((fields[var]->GetBndConditions()[i])
353 ->GetBoundaryConditionType() ==
355 {
356 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
357 }
358 }
359 }
360}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, and Vmath::Vcopy().

Referenced by NumFluxforScalar().

◆ ApplyVectorBCs()

void Nektar::SolverUtils::DiffusionLDG::ApplyVectorBCs ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const std::size_t  var,
const std::size_t  dir,
const Array< OneD, const NekDouble > &  qfield,
const Array< OneD, const NekDouble > &  qFwd,
const Array< OneD, const NekDouble > &  qBwd,
Array< OneD, NekDouble > &  penaltyflux 
)
private

Diffusion: Imposing weak boundary condition for q with flux uflux = g_D on Dirichlet boundary condition uflux = u_Fwd on Neumann boundary condition

Definition at line 427 of file DiffusionLDG.cpp.

434{
435 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
436 std::size_t cnt = 0;
437
438 for (std::size_t i = 0; i < nBndRegions; ++i)
439 {
440 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
442 {
443 continue;
444 }
445 std::size_t nBndEdges =
446 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
447
448 // Weakly impose boundary conditions by modifying flux values
449 for (std::size_t e = 0; e < nBndEdges; ++e)
450 {
451 std::size_t nBndEdgePts = fields[var]
452 ->GetBndCondExpansions()[i]
453 ->GetExp(e)
454 ->GetTotPoints();
455
456 std::size_t id1 =
457 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
458
459 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
460 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
461
462 // AV boundary conditions
463 if (boost::iequals(
464 fields[var]->GetBndConditions()[i]->GetUserDefined(),
465 "Wall") ||
466 boost::iequals(
467 fields[var]->GetBndConditions()[i]->GetUserDefined(),
468 "Symmetry") ||
469 boost::iequals(
470 fields[var]->GetBndConditions()[i]->GetUserDefined(),
471 "WallViscous") ||
472 boost::iequals(
473 fields[var]->GetBndConditions()[i]->GetUserDefined(),
474 "WallAdiabatic") ||
475 boost::iequals(
476 fields[var]->GetBndConditions()[i]->GetUserDefined(),
477 "WallRotational"))
478 {
479 Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
480 }
481 // For Dirichlet boundary condition:
482 // qflux = q+ - C_11 (u+ - g_D) (nx, ny)
483 else if (fields[var]
484 ->GetBndConditions()[i]
485 ->GetBoundaryConditionType() ==
487 {
488 Vmath::Vmul(nBndEdgePts, &m_traceNormals[dir][id2], 1,
489 &qFwd[id2], 1, &penaltyflux[id2], 1);
490 }
491 // For Neumann boundary condition: qflux = g_N
492 else if ((fields[var]->GetBndConditions()[i])
493 ->GetBoundaryConditionType() ==
495 {
497 nBndEdgePts, &m_traceNormals[dir][id2], 1,
498 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
499 1, &penaltyflux[id2], 1);
500 }
501 }
502 }
503}
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionLDG.h:104
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, m_traceNormals, Vmath::Vmul(), and Vmath::Zero().

Referenced by NumFluxforVector().

◆ create()

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

Definition at line 45 of file DiffusionLDG.h.

46 {
47 return DiffusionSharedPtr(new DiffusionLDG());
48 }
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:55

References DiffusionLDG().

◆ NumFluxforScalar()

void Nektar::SolverUtils::DiffusionLDG::NumFluxforScalar ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
TensorOfArray3D< NekDouble > &  uflux,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd 
)
private

Definition at line 235 of file DiffusionLDG.cpp.

241{
242 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
243 std::size_t nvariables = fields.size();
244 std::size_t nDim = fields[0]->GetCoordim(0);
245
246 Array<OneD, NekDouble> Fwd{nTracePts};
247 Array<OneD, NekDouble> Bwd{nTracePts};
248 Array<OneD, NekDouble> fluxtemp{nTracePts, 0.0};
249
250 // Get the sign of (v \cdot n), v = an arbitrary vector
251 // Evaluate upwind flux:
252 // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
253 for (std::size_t i = 0; i < nvariables; ++i)
254 {
255 // Compute Fwd and Bwd value of ufield of i direction
256 if (pFwd == NullNekDoubleArrayOfArray ||
258 {
259 fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
260 }
261 else
262 {
263 Fwd = pFwd[i];
264 Bwd = pBwd[i];
265 }
266
267 // Upwind
268 Vmath::Vcopy(nTracePts, Fwd, 1, fluxtemp, 1);
269
270 // Imposing weak boundary condition with flux
271 if (fields[0]->GetBndCondExpansions().size())
272 {
273 ApplyScalarBCs(fields, i, ufield[i], Fwd, Bwd, fluxtemp);
274 }
275
276 for (std::size_t j = 0; j < nDim; ++j)
277 {
278 Vmath::Vmul(nTracePts, m_traceNormals[j], 1, fluxtemp, 1,
279 uflux[j][i], 1);
280 }
281 }
282}
void ApplyScalarBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const Array< OneD, const NekDouble > &ufield, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &penaltyflux)
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray

References ApplyScalarBCs(), m_traceNormals, Nektar::NullNekDoubleArrayOfArray, Vmath::Vcopy(), and Vmath::Vmul().

Referenced by v_DiffuseCalcDerivative().

◆ NumFluxforVector()

void Nektar::SolverUtils::DiffusionLDG::NumFluxforVector ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
TensorOfArray3D< NekDouble > &  qfield,
Array< OneD, Array< OneD, NekDouble > > &  qflux 
)
private

Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to penalty term.

Definition at line 366 of file DiffusionLDG.cpp.

371{
372 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
373 std::size_t nvariables = fields.size();
374 std::size_t nDim = qfield.size();
375
376 Array<OneD, NekDouble> Fwd{nTracePts};
377 Array<OneD, NekDouble> Bwd{nTracePts};
378 Array<OneD, NekDouble> qFwd{nTracePts};
379 Array<OneD, NekDouble> qBwd{nTracePts};
380 Array<OneD, NekDouble> qfluxtemp{nTracePts, 0.0};
381 Array<OneD, NekDouble> uterm{nTracePts};
382
383 // Evaulate upwind flux:
384 // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
385 for (std::size_t i = 0; i < nvariables; ++i)
386 {
387 // Generate Stability term = - C11 ( u- - u+ )
388 fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
389 Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, uterm, 1);
390 Vmath::Smul(nTracePts, -m_C11, uterm, 1, uterm, 1);
391
392 qflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
393 for (std::size_t j = 0; j < nDim; ++j)
394 {
395 // Compute Fwd and Bwd value of ufield of jth direction
396 fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
397
398 // Downwind
399 Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
400
401 Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
402 qfluxtemp, 1);
403
404 // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
405 Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
406
407 // Imposing weak boundary condition with flux
408 if (fields[0]->GetBndCondExpansions().size())
409 {
410 ApplyVectorBCs(fields, i, j, qfield[j][i], qFwd, qBwd,
411 qfluxtemp);
412 }
413
414 // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
415 // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
416 // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
417 Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
418 }
419 }
420}
void ApplyVectorBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
NekDouble m_C11
Coefficient of penalty term.
Definition: DiffusionLDG.h:102
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.hpp:180
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 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 ApplyVectorBCs(), m_C11, m_traceNormals, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by v_DiffuseTraceFlux().

◆ v_Diffuse()

void Nektar::SolverUtils::DiffusionLDG::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,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Diffusion.

Definition at line 74 of file DiffusionLDG.cpp.

81{
82 std::size_t nCoeffs = fields[0]->GetNcoeffs();
83
84 Array<OneD, Array<OneD, NekDouble>> tmp{nConvectiveFields};
85 for (std::size_t i = 0; i < nConvectiveFields; ++i)
86 {
87 tmp[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
88 }
89
90 DiffusionLDG::v_DiffuseCoeffs(nConvectiveFields, fields, inarray, tmp, pFwd,
91 pBwd);
92
93 for (std::size_t i = 0; i < nConvectiveFields; ++i)
94 {
95 fields[i]->BwdTrans(tmp[i], outarray[i]);
96 }
97}
void v_DiffuseCoeffs(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, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override

References v_DiffuseCoeffs().

◆ v_DiffuseCalcDerivative()

void Nektar::SolverUtils::DiffusionLDG::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 
)
overrideprotectedvirtual

Diffusion Flux, calculate the physical derivatives.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 174 of file DiffusionLDG.cpp.

180{
181 std::size_t nConvectiveFields = fields.size();
182 std::size_t nDim = fields[0]->GetCoordim(0);
183 std::size_t nCoeffs = fields[0]->GetNcoeffs();
184 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
185
186 Array<OneD, NekDouble> tmp{nCoeffs};
187 TensorOfArray3D<NekDouble> flux{nDim};
188 for (std::size_t j = 0; j < nDim; ++j)
189 {
190 flux[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
191 for (std::size_t i = 0; i < nConvectiveFields; ++i)
192 {
193 flux[j][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
194 }
195 }
196
197 NumFluxforScalar(fields, inarray, flux, pFwd, pBwd);
198
199 for (std::size_t j = 0; j < nDim; ++j)
200 {
201 for (std::size_t i = 0; i < nConvectiveFields; ++i)
202 {
203 fields[i]->IProductWRTDerivBase(j, inarray[i], tmp);
204 Vmath::Neg(nCoeffs, tmp, 1);
205 fields[i]->AddTraceIntegral(flux[j][i], tmp);
206 fields[i]->SetPhysState(false);
207 fields[i]->MultiplyByElmtInvMass(tmp, tmp);
208 fields[i]->BwdTrans(tmp, qfield[j][i]);
209 }
210 }
211}
void NumFluxforScalar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &uflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292

References Vmath::Neg(), and NumFluxforScalar().

◆ v_DiffuseCoeffs()

void Nektar::SolverUtils::DiffusionLDG::v_DiffuseCoeffs ( 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,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 99 of file DiffusionLDG.cpp.

106{
107 if (fields[0]->GetGraph()->GetMovement()->GetMoveFlag()) // i.e. if
108 // m_ALESolver
109 {
110 fields[0]->GetTrace()->GetNormals(m_traceNormals);
111 }
112
113 std::size_t nDim = fields[0]->GetCoordim(0);
114 std::size_t nPts = fields[0]->GetTotPoints();
115 std::size_t nCoeffs = fields[0]->GetNcoeffs();
116 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
117
118 TensorOfArray3D<NekDouble> qfield{nDim};
119 for (std::size_t j = 0; j < nDim; ++j)
120 {
121 qfield[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
122 for (std::size_t i = 0; i < nConvectiveFields; ++i)
123 {
124 qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
125 }
126 }
127
128 Array<OneD, Array<OneD, NekDouble>> traceflux{nConvectiveFields};
129 for (std::size_t i = 0; i < nConvectiveFields; ++i)
130 {
131 traceflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
132 }
133
134 DiffuseCalcDerivative(fields, inarray, qfield, pFwd, pBwd);
135
136 // Initialize viscous tensor
137 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> viscTensor{nDim};
138 for (std::size_t j = 0; j < nDim; ++j)
139 {
140 viscTensor[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
141 for (std::size_t i = 0; i < nConvectiveFields; ++i)
142 {
143 viscTensor[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
144 }
145 }
146 DiffuseVolumeFlux(fields, inarray, qfield, viscTensor);
147 DiffuseTraceFlux(fields, inarray, qfield, viscTensor, traceflux, pFwd,
148 pBwd);
149
150 Array<OneD, Array<OneD, NekDouble>> qdbase{nDim};
151 for (std::size_t i = 0; i < nConvectiveFields; ++i)
152 {
153 for (std::size_t j = 0; j < nDim; ++j)
154 {
155 qdbase[j] = viscTensor[j][i];
156 }
157 fields[i]->IProductWRTDerivBase(qdbase, outarray[i]);
158
159 Vmath::Neg(nCoeffs, outarray[i], 1);
160 fields[i]->AddTraceIntegral(traceflux[i], outarray[i]);
161 fields[i]->SetPhysState(false);
162 }
163
164 if (!fields[0]->GetGraph()->GetMovement()->GetMoveFlag()) // i.e. if
165 // m_ALESolver
166 {
167 for (std::size_t i = 0; i < nConvectiveFields; ++i)
168 {
169 fields[i]->MultiplyByElmtInvMass(outarray[i], outarray[i]);
170 }
171 }
172}
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)
Definition: Diffusion.h:202
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.
Definition: Diffusion.h:215
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.
Definition: Diffusion.h:226

References Nektar::SolverUtils::Diffusion::DiffuseCalcDerivative(), Nektar::SolverUtils::Diffusion::DiffuseTraceFlux(), Nektar::SolverUtils::Diffusion::DiffuseVolumeFlux(), m_traceNormals, and Vmath::Neg().

Referenced by v_Diffuse().

◆ v_DiffuseTraceFlux()

void Nektar::SolverUtils::DiffusionLDG::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 
)
overrideprotectedvirtual

Diffusion term Trace Flux.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 222 of file DiffusionLDG.cpp.

231{
232 NumFluxforVector(fields, inarray, viscTensor, TraceFlux);
233}
void NumFluxforVector(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to p...

References NumFluxforVector().

◆ v_DiffuseVolumeFlux()

void Nektar::SolverUtils::DiffusionLDG::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 
)
overrideprotectedvirtual

Diffusion Volume Flux.

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 213 of file DiffusionLDG.cpp.

218{
219 m_fluxVector(inarray, qfield, viscTensor);
220}
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:353

References Nektar::SolverUtils::Diffusion::m_fluxVector.

◆ v_InitObject()

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

Implements Nektar::SolverUtils::Diffusion.

Definition at line 51 of file DiffusionLDG.cpp.

54{
55 m_session = pSession;
56
57 m_session->LoadSolverInfo("ShockCaptureType", m_shockCaptureType, "Off");
58
59 // Set up penalty term for LDG
60 m_session->LoadParameter("LDGc11", m_C11, 1.0);
61
62 // Setting up the normals
63 std::size_t nDim = pFields[0]->GetCoordim(0);
64 std::size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
65
66 m_traceNormals = Array<OneD, Array<OneD, NekDouble>>{nDim};
67 for (std::size_t i = 0; i < nDim; ++i)
68 {
69 m_traceNormals[i] = Array<OneD, NekDouble>{nTracePts};
70 }
71 pFields[0]->GetTrace()->GetNormals(m_traceNormals);
72}
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionLDG.h:105

References m_C11, m_session, m_shockCaptureType, and m_traceNormals.

Member Data Documentation

◆ m_C11

NekDouble Nektar::SolverUtils::DiffusionLDG::m_C11
private

Coefficient of penalty term.

Definition at line 102 of file DiffusionLDG.h.

Referenced by NumFluxforVector(), and v_InitObject().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::DiffusionLDG::m_session
private

Definition at line 105 of file DiffusionLDG.h.

Referenced by v_InitObject().

◆ m_shockCaptureType

std::string Nektar::SolverUtils::DiffusionLDG::m_shockCaptureType
private

Definition at line 99 of file DiffusionLDG.h.

Referenced by v_InitObject().

◆ m_traceNormals

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionLDG::m_traceNormals
private

◆ type

std::string Nektar::SolverUtils::DiffusionLDG::type
static
Initial value:
"LDG", DiffusionLDG::create, "Local Discontinuous Galkerin")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionLDG.h:45
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:39

Definition at line 50 of file DiffusionLDG.h.