Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Nektar::SolverUtils::DiffusionLDG Class Reference

#include <DiffusionLDG.h>

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

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)

Static Public Attributes

static std::string type

Protected Member Functions

 DiffusionLDG ()
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
virtual void v_Diffuse (const int nConvective, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
virtual void v_NumFluxforScalar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
virtual void v_WeakPenaltyforScalar (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const int var, const Array< OneD, const NekDouble > &ufield, Array< OneD, NekDouble > &penaltyflux)
virtual void v_NumFluxforVector (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
virtual void v_WeakPenaltyforVector (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const int var, const int dir, const Array< OneD, const NekDouble > &qfield, Array< OneD, NekDouble > &penaltyflux, NekDouble C11)
- Protected Member Functions inherited from Nektar::SolverUtils::Diffusion
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
virtual Array< OneD, Array
< OneD, Array< OneD, NekDouble > > > & 
v_GetFluxTensor ()

Protected Attributes

std::string m_shockCaptureType
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
LibUtilities::SessionReaderSharedPtr m_session
- Protected Attributes inherited from Nektar::SolverUtils::Diffusion
DiffusionFluxVecCB m_fluxVector
DiffusionFluxVecCBNS m_fluxVectorNS
RiemannSolverSharedPtr m_riemann
DiffusionArtificialDiffusion m_ArtificialDiffusionVector

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Diffusion
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
SOLVER_UTILS_EXPORT void Diffuse (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SOLVER_UTILS_EXPORT void FluxVec (Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &fluxvector)
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVector (FuncPointerT func, ObjectPointerT obj)
void SetFluxVectorVec (DiffusionFluxVecCB fluxVector)
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVectorNS (FuncPointerT func, ObjectPointerT obj)
template<typename FuncPointerT , typename ObjectPointerT >
void SetArtificialDiffusionVector (FuncPointerT func, ObjectPointerT obj)
void SetFluxVectorNS (DiffusionFluxVecCBNS fluxVector)
void SetRiemannSolver (RiemannSolverSharedPtr riemann)
void SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
virtual Array< OneD, Array
< OneD, Array< OneD, NekDouble > > > & 
GetFluxTensor ()

Detailed Description

Definition at line 45 of file DiffusionLDG.h.

Constructor & Destructor Documentation

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

Definition at line 47 of file DiffusionLDG.cpp.

Referenced by create().

{
}

Member Function Documentation

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

Definition at line 48 of file DiffusionLDG.h.

References DiffusionLDG().

{
}
void Nektar::SolverUtils::DiffusionLDG::v_Diffuse ( const int  nConvective,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protectedvirtual

Implements Nektar::SolverUtils::Diffusion.

Definition at line 73 of file DiffusionLDG.cpp.

References Nektar::SolverUtils::Diffusion::m_ArtificialDiffusionVector, m_shockCaptureType, Vmath::Neg(), v_NumFluxforScalar(), v_NumFluxforVector(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vmul().

{
int nBndEdgePts, i, j, k;
int nDim = fields[0]->GetCoordim(0);
int nPts = fields[0]->GetTotPoints();
int nCoeffs = fields[0]->GetNcoeffs();
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
Array<OneD, NekDouble> qcoeffs(nCoeffs);
Array<OneD, NekDouble> temp (nCoeffs);
Array<OneD, Array<OneD, NekDouble> > fluxvector(nDim);
Array<OneD, Array<OneD, NekDouble> > tmp(nConvectiveFields);
Array<OneD, Array<OneD, Array<OneD, NekDouble> > > flux (nDim);
Array<OneD, Array<OneD, Array<OneD, NekDouble> > > qfield(nDim);
Array<OneD, Array<OneD, Array<OneD, NekDouble> > > qfieldStd(nDim);
for (j = 0; j < nDim; ++j)
{
qfield[j] =
Array<OneD, Array<OneD, NekDouble> >(nConvectiveFields);
qfieldStd[j] =
Array<OneD, Array<OneD, NekDouble> >(nConvectiveFields);
flux[j] =
Array<OneD, Array<OneD, NekDouble> >(nConvectiveFields);
for (i = 0; i < nConvectiveFields; ++i)
{
qfield[j][i] = Array<OneD, NekDouble>(nPts, 0.0);
qfieldStd[j][i] = Array<OneD, NekDouble>(nPts, 0.0);
flux[j][i] = Array<OneD, NekDouble>(nTracePts, 0.0);
}
}
for (k = 0; k < nDim; ++k)
{
fluxvector[k] = Array<OneD, NekDouble>(nPts, 0.0);
}
// Compute q_{\eta} and q_{\xi}
// Obtain numerical fluxes
v_NumFluxforScalar(fields, inarray, flux);
for (j = 0; j < nDim; ++j)
{
for (i = 0; i < nConvectiveFields; ++i)
{
fields[i]->IProductWRTDerivBase(j, inarray[i], qcoeffs);
Vmath::Neg (nCoeffs, qcoeffs, 1);
fields[i]->AddTraceIntegral (flux[j][i], qcoeffs);
fields[i]->SetPhysState (false);
fields[i]->MultiplyByElmtInvMass(qcoeffs, qcoeffs);
fields[i]->BwdTrans (qcoeffs, qfield[j][i]);
}
}
// Compute u from q_{\eta} and q_{\xi}
// Obtain numerical fluxes
v_NumFluxforVector(fields, inarray, qfield, flux[0]);
{
Array<OneD, NekDouble> muvar(nPts, 0.0);
m_ArtificialDiffusionVector(inarray, muvar);
int numConvFields = nConvectiveFields;
if (m_shockCaptureType == "Smooth")
{
numConvFields = nConvectiveFields - 1;
}
for (j = 0; j < nDim; ++j)
{
for (i = 0; i < numConvFields; ++i)
{
Vmath::Vmul(nPts,qfield[j][i],1,muvar,1,qfield[j][i],1);
}
}
Array<OneD, NekDouble> FwdMuVar(nTracePts, 0.0);
Array<OneD, NekDouble> BwdMuVar(nTracePts, 0.0);
fields[0]->GetFwdBwdTracePhys(muvar,FwdMuVar,BwdMuVar);
int nBndRegions = fields[0]->GetBndCondExpansions().
num_elements();
int cnt = 0;
for (int i = 0; i < nBndRegions; ++i)
{
// Number of boundary expansion related to that region
int nBndEdges = fields[0]->
GetBndCondExpansions()[i]->GetExpSize();
// Weakly impose boundary conditions by modifying flux
// values
for (int e = 0; e < nBndEdges ; ++e)
{
nBndEdgePts = fields[0]->GetBndCondExpansions()[i]
->GetExp(e)->GetTotPoints();
int id2 = fields[0]->GetTrace()->GetPhys_Offset(
fields[0]->GetTraceMap()
->GetBndCondTraceToGlobalTraceMap(cnt++));
for (int k = 0; k < nBndEdgePts; ++k)
{
BwdMuVar[id2+k] = 0.0;
}
}
}
for(i = 0; i < numConvFields; ++i)
{
for(int k = 0; k < nTracePts; ++k)
{
flux[0][i][k] =
0.5 * (FwdMuVar[k] + BwdMuVar[k]) * flux[0][i][k];
}
}
}
for (i = 0; i < nConvectiveFields; ++i)
{
tmp[i] = Array<OneD, NekDouble>(nCoeffs, 0.0);
for (j = 0; j < nDim; ++j)
{
Vmath::Vcopy(nPts, qfield[j][i], 1, fluxvector[j], 1);
fields[i]->IProductWRTDerivBase(j, fluxvector[j], qcoeffs);
Vmath::Vadd(nCoeffs, qcoeffs, 1, tmp[i], 1, tmp[i], 1);
}
// Evaulate <\phi, \hat{F}\cdot n> - outarray[i]
Vmath::Neg (nCoeffs, tmp[i], 1);
fields[i]->AddTraceIntegral (flux[0][i], tmp[i]);
fields[i]->SetPhysState (false);
fields[i]->MultiplyByElmtInvMass(tmp[i], tmp[i]);
fields[i]->BwdTrans (tmp[i], outarray[i]);
}
}
void Nektar::SolverUtils::DiffusionLDG::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 51 of file DiffusionLDG.cpp.

References m_session, m_shockCaptureType, and m_traceNormals.

{
m_session = pSession;
m_session->LoadSolverInfo("ShockCaptureType",
// Setting up the normals
int i;
int nDim = pFields[0]->GetCoordim(0);
int nTracePts = pFields[0]->GetTrace()->GetTotPoints();
m_traceNormals = Array<OneD, Array<OneD, NekDouble> >(nDim);
for(i = 0; i < nDim; ++i)
{
m_traceNormals[i] = Array<OneD, NekDouble> (nTracePts);
}
pFields[0]->GetTrace()->GetNormals(m_traceNormals);
}
void Nektar::SolverUtils::DiffusionLDG::v_NumFluxforScalar ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  uflux 
)
protectedvirtual

Definition at line 222 of file DiffusionLDG.cpp.

References m_traceNormals, Vmath::Svtvp(), v_WeakPenaltyforScalar(), and Vmath::Vmul().

Referenced by v_Diffuse().

{
int i, j;
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
int nvariables = fields.num_elements();
int nDim = uflux.num_elements();
Array<OneD, NekDouble > Fwd (nTracePts);
Array<OneD, NekDouble > Bwd (nTracePts);
Array<OneD, NekDouble > Vn (nTracePts, 0.0);
Array<OneD, NekDouble > fluxtemp(nTracePts, 0.0);
// Get the normal velocity Vn
for(i = 0; i < nDim; ++i)
{
Vmath::Svtvp(nTracePts, 1.0, m_traceNormals[i], 1,
Vn, 1, Vn, 1);
}
// Get the sign of (v \cdot n), v = an arbitrary vector
// Evaluate upwind flux:
// uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
for (j = 0; j < nDim; ++j)
{
for (i = 0; i < nvariables ; ++i)
{
// Compute Fwd and Bwd value of ufield of i direction
fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
// if Vn >= 0, flux = uFwd, i.e.,
// edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uFwd
// edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uFwd
// else if Vn < 0, flux = uBwd, i.e.,
// edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uBwd
// edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uBwd
fields[i]->GetTrace()->Upwind(/*m_traceNormals[j]*/Vn,
Fwd, Bwd, fluxtemp);
// Imposing weak boundary condition with flux
// if Vn >= 0, uflux = uBwd at Neumann, i.e.,
// edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uBwd
// edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uBwd
// if Vn >= 0, uflux = uFwd at Neumann, i.e.,
// edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uFwd
// edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uFwd
if(fields[0]->GetBndCondExpansions().num_elements())
{
v_WeakPenaltyforScalar(fields, i, ufield[i], fluxtemp);
}
// if Vn >= 0, flux = uFwd*(tan_{\xi}^- \cdot \vec{n}),
// i.e,
// edge::eForward, uFwd \(\tan_{\xi}^Fwd \cdot \vec{n})
// edge::eBackward, uFwd \(\tan_{\xi}^Bwd \cdot \vec{n})
// else if Vn < 0, flux = uBwd*(tan_{\xi}^- \cdot \vec{n}),
// i.e,
// edge::eForward, uBwd \(\tan_{\xi}^Fwd \cdot \vec{n})
// edge::eBackward, uBwd \(\tan_{\xi}^Bwd \cdot \vec{n})
Vmath::Vmul(nTracePts,
fluxtemp, 1,
uflux[j][i], 1);
}
}
}
void Nektar::SolverUtils::DiffusionLDG::v_NumFluxforVector ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  qfield,
Array< OneD, Array< OneD, NekDouble > > &  qflux 
)
protectedvirtual

Definition at line 358 of file DiffusionLDG.cpp.

References m_traceNormals, Vmath::Smul(), Vmath::Svtvp(), v_WeakPenaltyforVector(), Vmath::Vadd(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by v_Diffuse().

{
int i, j;
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
int nvariables = fields.num_elements();
int nDim = qfield.num_elements();
NekDouble C11 = 0.0;
Array<OneD, NekDouble > Fwd(nTracePts);
Array<OneD, NekDouble > Bwd(nTracePts);
Array<OneD, NekDouble > Vn (nTracePts, 0.0);
Array<OneD, NekDouble > qFwd (nTracePts);
Array<OneD, NekDouble > qBwd (nTracePts);
Array<OneD, NekDouble > qfluxtemp(nTracePts, 0.0);
Array<OneD, NekDouble > uterm(nTracePts);
/*
// Setting up the normals
m_traceNormals = Array<OneD, Array<OneD, NekDouble> >(nDim);
for(i = 0; i < nDim; ++i)
{
m_traceNormals[i] = Array<OneD, NekDouble> (nTracePts);
}
fields[0]->GetTrace()->GetNormals(m_traceNormals);
*/
// Get the normal velocity Vn
for(i = 0; i < nDim; ++i)
{
Vmath::Svtvp(nTracePts, 1.0, m_traceNormals[i], 1,
Vn, 1, Vn, 1);
}
// Evaulate upwind flux:
// qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
for (i = 0; i < nvariables; ++i)
{
qflux[i] = Array<OneD, NekDouble> (nTracePts, 0.0);
for (j = 0; j < nDim; ++j)
{
// Compute Fwd and Bwd value of ufield of jth direction
fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
// if Vn >= 0, flux = uFwd, i.e.,
// edge::eForward, if V*n>=0 <=> V*n_F>=0, pick
// qflux = qBwd = q+
// edge::eBackward, if V*n>=0 <=> V*n_B<0, pick
// qflux = qBwd = q-
// else if Vn < 0, flux = uBwd, i.e.,
// edge::eForward, if V*n<0 <=> V*n_F<0, pick
// qflux = qFwd = q-
// edge::eBackward, if V*n<0 <=> V*n_B>=0, pick
// qflux = qFwd = q+
fields[i]->GetTrace()->Upwind(/*m_traceNormals[j]*/Vn,
qBwd, qFwd,
qfluxtemp);
Vmath::Vmul(nTracePts,
qfluxtemp, 1,
qfluxtemp, 1);
// Generate Stability term = - C11 ( u- - u+ )
fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
Vmath::Vsub(nTracePts,
Fwd, 1, Bwd, 1,
uterm, 1);
Vmath::Smul(nTracePts,
-1.0 * C11, uterm, 1,
uterm, 1);
// Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
Vmath::Vadd(nTracePts,
uterm, 1,
qfluxtemp, 1,
qfluxtemp, 1);
// Imposing weak boundary condition with flux
if (fields[0]->GetBndCondExpansions().num_elements())
{
v_WeakPenaltyforVector(fields, i, j,
qfield[j][i],
qfluxtemp, C11);
}
// q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
// n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
// n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
Vmath::Vadd(nTracePts,
qfluxtemp, 1,
qflux[i], 1,
qflux[i], 1);
}
}
}
void Nektar::SolverUtils::DiffusionLDG::v_WeakPenaltyforScalar ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const int  var,
const Array< OneD, const NekDouble > &  ufield,
Array< OneD, NekDouble > &  penaltyflux 
)
protectedvirtual

Definition at line 299 of file DiffusionLDG.cpp.

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

Referenced by v_NumFluxforScalar().

{
int i, e, id1, id2;
// Number of boundary regions
int nBndEdgePts, nBndEdges;
int cnt = 0;
int nBndRegions = fields[var]->GetBndCondExpansions().num_elements();
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
Array<OneD, NekDouble > uplus(nTracePts);
fields[var]->ExtractTracePhys(ufield, uplus);
for (i = 0; i < nBndRegions; ++i)
{
// Number of boundary expansion related to that region
nBndEdges = fields[var]->
GetBndCondExpansions()[i]->GetExpSize();
// Weakly impose boundary conditions by modifying flux values
for (e = 0; e < nBndEdges ; ++e)
{
nBndEdgePts = fields[var]->
GetBndCondExpansions()[i]->GetExp(e)->GetTotPoints();
id1 = fields[var]->
GetBndCondExpansions()[i]->GetPhys_Offset(e);
id2 = fields[0]->GetTrace()->
GetPhys_Offset(fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
// For Dirichlet boundary condition: uflux = g_D
if (fields[var]->GetBndConditions()[i]->
GetBoundaryConditionType() == SpatialDomains::eDirichlet)
{
Vmath::Vcopy(nBndEdgePts,
&(fields[var]->
GetBndCondExpansions()[i]->
GetPhys())[id1], 1,
&penaltyflux[id2], 1);
}
// For Neumann boundary condition: uflux = u+
else if ((fields[var]->GetBndConditions()[i])->
GetBoundaryConditionType() == SpatialDomains::eNeumann)
{
Vmath::Vcopy(nBndEdgePts,
&uplus[id2], 1,
&penaltyflux[id2], 1);
}
}
}
}
void Nektar::SolverUtils::DiffusionLDG::v_WeakPenaltyforVector ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const int  var,
const int  dir,
const Array< OneD, const NekDouble > &  qfield,
Array< OneD, NekDouble > &  penaltyflux,
NekDouble  C11 
)
protectedvirtual

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

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

Referenced by v_NumFluxforVector().

{
int i, e, id1, id2;
int nBndEdges, nBndEdgePts;
int nBndRegions = fields[var]->GetBndCondExpansions().num_elements();
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
Array<OneD, NekDouble > uterm(nTracePts);
Array<OneD, NekDouble > qtemp(nTracePts);
int cnt = 0;
/*
// Setting up the normals
m_traceNormals = Array<OneD, Array<OneD, NekDouble> >(nDim);
for(i = 0; i < nDim; ++i)
{
m_traceNormals[i] = Array<OneD, NekDouble> (nTracePts);
}
fields[0]->GetTrace()->GetNormals(m_traceNormals);
*/
fields[var]->ExtractTracePhys(qfield, qtemp);
for (i = 0; i < nBndRegions; ++i)
{
nBndEdges = fields[var]->
GetBndCondExpansions()[i]->GetExpSize();
// Weakly impose boundary conditions by modifying flux values
for (e = 0; e < nBndEdges ; ++e)
{
nBndEdgePts = fields[var]->
GetBndCondExpansions()[i]->GetExp(e)->GetTotPoints();
id1 = fields[var]->
GetBndCondExpansions()[i]->GetPhys_Offset(e);
id2 = fields[0]->GetTrace()->
GetPhys_Offset(fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
// For Dirichlet boundary condition:
//qflux = q+ - C_11 (u+ - g_D) (nx, ny)
if(fields[var]->GetBndConditions()[i]->
GetBoundaryConditionType() == SpatialDomains::eDirichlet)
{
Vmath::Vmul(nBndEdgePts,
&m_traceNormals[dir][id2], 1,
&qtemp[id2], 1,
&penaltyflux[id2], 1);
}
// For Neumann boundary condition: qflux = g_N
else if((fields[var]->GetBndConditions()[i])->
GetBoundaryConditionType() == SpatialDomains::eNeumann)
{
Vmath::Vmul(nBndEdgePts,
&m_traceNormals[dir][id2], 1,
&(fields[var]->
GetBndCondExpansions()[i]->
GetPhys())[id1], 1,
&penaltyflux[id2], 1);
}
}
}
}

Member Data Documentation

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::DiffusionLDG::m_session
protected

Definition at line 61 of file DiffusionLDG.h.

Referenced by v_InitObject().

std::string Nektar::SolverUtils::DiffusionLDG::m_shockCaptureType
protected

Definition at line 58 of file DiffusionLDG.h.

Referenced by v_Diffuse(), and v_InitObject().

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionLDG::m_traceNormals
protected
std::string Nektar::SolverUtils::DiffusionLDG::type
static
Initial value:
RegisterCreatorFunction("LDG", DiffusionLDG::create)

Definition at line 53 of file DiffusionLDG.h.