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)
 
- 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
 
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 273 of file DiffusionLDG.cpp.

280{
281 // Number of boundary regions
282 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
283 std::size_t cnt = 0;
284 for (std::size_t i = 0; i < nBndRegions; ++i)
285 {
286 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
288 {
289 continue;
290 }
291
292 // Number of boundary expansion related to that region
293 std::size_t nBndEdges =
294 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
295
296 // Weakly impose boundary conditions by modifying flux values
297 for (std::size_t e = 0; e < nBndEdges; ++e)
298 {
299 std::size_t nBndEdgePts = fields[var]
300 ->GetBndCondExpansions()[i]
301 ->GetExp(e)
302 ->GetTotPoints();
303
304 std::size_t id1 =
305 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
306
307 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
308 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
309
310 // AV boundary conditions
311 if (boost::iequals(
312 fields[var]->GetBndConditions()[i]->GetUserDefined(),
313 "Wall") ||
314 boost::iequals(
315 fields[var]->GetBndConditions()[i]->GetUserDefined(),
316 "Symmetry") ||
317 boost::iequals(
318 fields[var]->GetBndConditions()[i]->GetUserDefined(),
319 "WallViscous") ||
320 boost::iequals(
321 fields[var]->GetBndConditions()[i]->GetUserDefined(),
322 "WallAdiabatic"))
323 {
324 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
325 }
326 // For Dirichlet boundary condition: uflux = g_D
327 else if (fields[var]
328 ->GetBndConditions()[i]
329 ->GetBoundaryConditionType() ==
331 {
333 nBndEdgePts,
334 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
335 1, &penaltyflux[id2], 1);
336 }
337 // For Neumann boundary condition: uflux = u+
338 else if ((fields[var]->GetBndConditions()[i])
339 ->GetBoundaryConditionType() ==
341 {
342 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
343 }
344 }
345 }
346}
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 413 of file DiffusionLDG.cpp.

420{
421 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
422 std::size_t cnt = 0;
423
424 for (std::size_t i = 0; i < nBndRegions; ++i)
425 {
426 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
428 {
429 continue;
430 }
431 std::size_t nBndEdges =
432 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
433
434 // Weakly impose boundary conditions by modifying flux values
435 for (std::size_t e = 0; e < nBndEdges; ++e)
436 {
437 std::size_t nBndEdgePts = fields[var]
438 ->GetBndCondExpansions()[i]
439 ->GetExp(e)
440 ->GetTotPoints();
441
442 std::size_t id1 =
443 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
444
445 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
446 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
447
448 // AV boundary conditions
449 if (boost::iequals(
450 fields[var]->GetBndConditions()[i]->GetUserDefined(),
451 "Wall") ||
452 boost::iequals(
453 fields[var]->GetBndConditions()[i]->GetUserDefined(),
454 "Symmetry") ||
455 boost::iequals(
456 fields[var]->GetBndConditions()[i]->GetUserDefined(),
457 "WallViscous") ||
458 boost::iequals(
459 fields[var]->GetBndConditions()[i]->GetUserDefined(),
460 "WallAdiabatic"))
461 {
462 Vmath::Zero(nBndEdgePts, &penaltyflux[id2], 1);
463 }
464 // For Dirichlet boundary condition:
465 // qflux = q+ - C_11 (u+ - g_D) (nx, ny)
466 else if (fields[var]
467 ->GetBndConditions()[i]
468 ->GetBoundaryConditionType() ==
470 {
471 Vmath::Vmul(nBndEdgePts, &m_traceNormals[dir][id2], 1,
472 &qFwd[id2], 1, &penaltyflux[id2], 1);
473 }
474 // For Neumann boundary condition: qflux = g_N
475 else if ((fields[var]->GetBndConditions()[i])
476 ->GetBoundaryConditionType() ==
478 {
480 nBndEdgePts, &m_traceNormals[dir][id2], 1,
481 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
482 1, &penaltyflux[id2], 1);
483 }
484 }
485 }
486}
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:54

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 224 of file DiffusionLDG.cpp.

230{
231 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
232 std::size_t nvariables = fields.size();
233 std::size_t nDim = fields[0]->GetCoordim(0);
234
235 Array<OneD, NekDouble> Fwd{nTracePts};
236 Array<OneD, NekDouble> Bwd{nTracePts};
237 Array<OneD, NekDouble> fluxtemp{nTracePts, 0.0};
238
239 // Get the sign of (v \cdot n), v = an arbitrary vector
240 // Evaluate upwind flux:
241 // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
242 for (std::size_t i = 0; i < nvariables; ++i)
243 {
244 // Compute Fwd and Bwd value of ufield of i direction
245 if (pFwd == NullNekDoubleArrayOfArray ||
247 {
248 fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
249 }
250 else
251 {
252 Fwd = pFwd[i];
253 Bwd = pBwd[i];
254 }
255
256 // Upwind
257 Vmath::Vcopy(nTracePts, Fwd, 1, fluxtemp, 1);
258
259 // Imposing weak boundary condition with flux
260 if (fields[0]->GetBndCondExpansions().size())
261 {
262 ApplyScalarBCs(fields, i, ufield[i], Fwd, Bwd, fluxtemp);
263 }
264
265 for (std::size_t j = 0; j < nDim; ++j)
266 {
267 Vmath::Vmul(nTracePts, m_traceNormals[j], 1, fluxtemp, 1,
268 uflux[j][i], 1);
269 }
270 }
271}
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 352 of file DiffusionLDG.cpp.

357{
358 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
359 std::size_t nvariables = fields.size();
360 std::size_t nDim = qfield.size();
361
362 Array<OneD, NekDouble> Fwd{nTracePts};
363 Array<OneD, NekDouble> Bwd{nTracePts};
364 Array<OneD, NekDouble> qFwd{nTracePts};
365 Array<OneD, NekDouble> qBwd{nTracePts};
366 Array<OneD, NekDouble> qfluxtemp{nTracePts, 0.0};
367 Array<OneD, NekDouble> uterm{nTracePts};
368
369 // Evaulate upwind flux:
370 // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
371 for (std::size_t i = 0; i < nvariables; ++i)
372 {
373 // Generate Stability term = - C11 ( u- - u+ )
374 fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
375 Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, uterm, 1);
376 Vmath::Smul(nTracePts, -m_C11, uterm, 1, uterm, 1);
377
378 qflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
379 for (std::size_t j = 0; j < nDim; ++j)
380 {
381 // Compute Fwd and Bwd value of ufield of jth direction
382 fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
383
384 // Downwind
385 Vmath::Vcopy(nTracePts, qBwd, 1, qfluxtemp, 1);
386
387 Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
388 qfluxtemp, 1);
389
390 // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
391 Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
392
393 // Imposing weak boundary condition with flux
394 if (fields[0]->GetBndCondExpansions().size())
395 {
396 ApplyVectorBCs(fields, i, j, qfield[j][i], qFwd, qBwd,
397 qfluxtemp);
398 }
399
400 // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
401 // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
402 // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
403 Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
404 }
405 }
406}
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 163 of file DiffusionLDG.cpp.

169{
170 std::size_t nConvectiveFields = fields.size();
171 std::size_t nDim = fields[0]->GetCoordim(0);
172 std::size_t nCoeffs = fields[0]->GetNcoeffs();
173 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
174
175 Array<OneD, NekDouble> tmp{nCoeffs};
176 TensorOfArray3D<NekDouble> flux{nDim};
177 for (std::size_t j = 0; j < nDim; ++j)
178 {
179 flux[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
180 for (std::size_t i = 0; i < nConvectiveFields; ++i)
181 {
182 flux[j][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
183 }
184 }
185
186 NumFluxforScalar(fields, inarray, flux, pFwd, pBwd);
187
188 for (std::size_t j = 0; j < nDim; ++j)
189 {
190 for (std::size_t i = 0; i < nConvectiveFields; ++i)
191 {
192 fields[i]->IProductWRTDerivBase(j, inarray[i], tmp);
193 Vmath::Neg(nCoeffs, tmp, 1);
194 fields[i]->AddTraceIntegral(flux[j][i], tmp);
195 fields[i]->SetPhysState(false);
196 fields[i]->MultiplyByElmtInvMass(tmp, tmp);
197 fields[i]->BwdTrans(tmp, qfield[j][i]);
198 }
199 }
200}
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 std::size_t nDim = fields[0]->GetCoordim(0);
108 std::size_t nPts = fields[0]->GetTotPoints();
109 std::size_t nCoeffs = fields[0]->GetNcoeffs();
110 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
111
112 Array<OneD, NekDouble> tmp{nCoeffs};
113
114 TensorOfArray3D<NekDouble> qfield{nDim};
115 for (std::size_t j = 0; j < nDim; ++j)
116 {
117 qfield[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
118 for (std::size_t i = 0; i < nConvectiveFields; ++i)
119 {
120 qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
121 }
122 }
123
124 Array<OneD, Array<OneD, NekDouble>> traceflux{nConvectiveFields};
125 for (std::size_t i = 0; i < nConvectiveFields; ++i)
126 {
127 traceflux[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
128 }
129
130 DiffuseCalcDerivative(fields, inarray, qfield, pFwd, pBwd);
131
132 // Initialize viscous tensor
133 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> viscTensor{nDim};
134 for (std::size_t j = 0; j < nDim; ++j)
135 {
136 viscTensor[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
137 for (std::size_t i = 0; i < nConvectiveFields; ++i)
138 {
139 viscTensor[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
140 }
141 }
142 DiffuseVolumeFlux(fields, inarray, qfield, viscTensor);
143 DiffuseTraceFlux(fields, inarray, qfield, viscTensor, traceflux, pFwd,
144 pBwd);
145
146 Array<OneD, Array<OneD, NekDouble>> qdbase{nDim};
147
148 for (std::size_t i = 0; i < nConvectiveFields; ++i)
149 {
150 for (std::size_t j = 0; j < nDim; ++j)
151 {
152 qdbase[j] = viscTensor[j][i];
153 }
154 fields[i]->IProductWRTDerivBase(qdbase, tmp);
155
156 Vmath::Neg(nCoeffs, tmp, 1);
157 fields[i]->AddTraceIntegral(traceflux[i], tmp);
158 fields[i]->SetPhysState(false);
159 fields[i]->MultiplyByElmtInvMass(tmp, outarray[i]);
160 }
161}
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:201
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:214
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:225

References Nektar::SolverUtils::Diffusion::DiffuseCalcDerivative(), Nektar::SolverUtils::Diffusion::DiffuseTraceFlux(), Nektar::SolverUtils::Diffusion::DiffuseVolumeFlux(), 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 211 of file DiffusionLDG.cpp.

220{
221 NumFluxforVector(fields, inarray, viscTensor, TraceFlux);
222}
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 202 of file DiffusionLDG.cpp.

207{
208 m_fluxVector(inarray, qfield, viscTensor);
209}
DiffusionFluxVecCB m_fluxVector
Definition: Diffusion.h:346

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

Definition at line 104 of file DiffusionLDG.h.

Referenced by ApplyVectorBCs(), NumFluxforScalar(), NumFluxforVector(), and v_InitObject().

◆ 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.
Definition: NekFactory.hpp:197
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionLDG.h:45
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:39

Definition at line 50 of file DiffusionLDG.h.