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::DiffusionLDGNS Class Reference

#include <DiffusionLDGNS.h>

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

Static Public Member Functions

static DiffusionSharedPtr create (std::string diffType)

Static Public Attributes

static std::string type

Protected Member Functions

 DiffusionLDGNS ()
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)
 Calculate weak DG Diffusion in the LDG form for the Navier-Stokes (NS) equations:
virtual void v_NumericalFluxO1 (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &numericalFluxO1)
 Builds the numerical flux for the 1st order derivatives.
virtual void v_WeakPenaltyO1 (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &penaltyfluxO1)
 Imposes appropriate bcs for the 1st order derivatives.
virtual void v_NumericalFluxO2 (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)
 Build the numerical flux for the 2nd order derivatives.
virtual void v_WeakPenaltyO2 (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const int var, const int dir, const Array< OneD, const NekDouble > &qfield, Array< OneD, NekDouble > &penaltyflux)
 Imposes appropriate bcs for the 2nd order derivatives.
virtual void v_SetHomoDerivs (Array< OneD, Array< OneD, NekDouble > > &deriv)
virtual Array< OneD, Array
< OneD, Array< OneD, NekDouble > > > & 
v_GetFluxTensor ()

Protected Attributes

Array< OneD, Array< OneD,
NekDouble > > 
m_traceVel
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
LibUtilities::SessionReaderSharedPtr m_session
NekDouble m_gamma
NekDouble m_gasConstant
NekDouble m_Twall
std::string m_ViscosityType
NekDouble m_mu
NekDouble m_thermalConductivity
NekDouble m_rhoInf
NekDouble m_pInf
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_viscTensor
Array< OneD, Array< OneD,
NekDouble > > 
m_homoDerivs
int m_spaceDim
int m_diffDim
- 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 DiffusionLDGNS.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::DiffusionLDGNS::DiffusionLDGNS ( )
protected

Definition at line 49 of file DiffusionLDGNS.cpp.

Referenced by create().

{
}

Member Function Documentation

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

Definition at line 48 of file DiffusionLDGNS.h.

References DiffusionLDGNS().

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

Calculate weak DG Diffusion in the LDG form for the Navier-Stokes (NS) equations:

$ \langle\psi, \hat{u}\cdot n\rangle - \langle\nabla\psi \cdot u\rangle \langle\phi, \hat{q}\cdot n\rangle - (\nabla \phi \cdot q) \rangle $

The equations that need a diffusion operator are those related with the velocities and with the energy.

Implements Nektar::SolverUtils::Diffusion.

Definition at line 105 of file DiffusionLDGNS.cpp.

References m_diffDim, Nektar::SolverUtils::Diffusion::m_fluxVectorNS, m_homoDerivs, m_spaceDim, m_viscTensor, Vmath::Neg(), v_NumericalFluxO1(), v_NumericalFluxO2(), and Vmath::Vadd().

{
int i, j;
int nDim = fields[0]->GetCoordim(0);
int nScalars = inarray.num_elements();
int nPts = fields[0]->GetTotPoints();
int nCoeffs = fields[0]->GetNcoeffs();
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
Array<OneD, NekDouble> tmp1(nCoeffs);
Array<OneD, Array<OneD, NekDouble> > tmp2(nConvectiveFields);
Array<OneD, Array<OneD, Array<OneD, NekDouble> > >
numericalFluxO1(m_spaceDim);
Array<OneD, Array<OneD, Array<OneD, NekDouble> > >
derivativesO1(m_spaceDim);
Array<OneD, Array<OneD, NekDouble> > fluxvector(m_spaceDim);
for (j = 0; j < m_spaceDim; ++j)
{
numericalFluxO1[j] = Array<OneD, Array<OneD, NekDouble> >(
nScalars);
derivativesO1[j] = Array<OneD, Array<OneD, NekDouble> >(
nScalars);
for (i = 0; i < nScalars; ++i)
{
numericalFluxO1[j][i] = Array<OneD, NekDouble>(
nTracePts, 0.0);
derivativesO1[j][i] = Array<OneD, NekDouble>(nPts, 0.0);
}
}
// Compute the numerical fluxes for the first order derivatives
v_NumericalFluxO1(fields, inarray, numericalFluxO1);
for (j = 0; j < nDim; ++j)
{
for (i = 0; i < nScalars; ++i)
{
fields[i]->IProductWRTDerivBase (j, inarray[i], tmp1);
Vmath::Neg (nCoeffs, tmp1, 1);
fields[i]->AddTraceIntegral (numericalFluxO1[j][i],
tmp1);
fields[i]->SetPhysState (false);
fields[i]->MultiplyByElmtInvMass(tmp1, tmp1);
fields[i]->BwdTrans (tmp1, derivativesO1[j][i]);
}
}
// For 3D Homogeneous 1D only take derivatives in 3rd direction
if (m_diffDim == 1)
{
for (i = 0; i < nScalars; ++i)
{
derivativesO1[2][i] = m_homoDerivs[i];
}
}
// Initialisation viscous tensor
m_viscTensor = Array<OneD, Array<OneD, Array<OneD, NekDouble> > >
Array<OneD, Array<OneD, NekDouble> > viscousFlux(nConvectiveFields);
for (j = 0; j < m_spaceDim; ++j)
{
m_viscTensor[j] = Array<OneD, Array<OneD, NekDouble> >(
nScalars+1);
for (i = 0; i < nScalars+1; ++i)
{
m_viscTensor[j][i] = Array<OneD, NekDouble>(nPts, 0.0);
}
}
for (i = 0; i < nConvectiveFields; ++i)
{
viscousFlux[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
}
m_fluxVectorNS(inarray, derivativesO1, m_viscTensor);
// Compute u from q_{\eta} and q_{\xi}
// Obtain numerical fluxes
v_NumericalFluxO2(fields, inarray, m_viscTensor, viscousFlux);
for (i = 0; i < nConvectiveFields; ++i)
{
tmp2[i] = Array<OneD, NekDouble>(nCoeffs, 0.0);
for (j = 0; j < nDim; ++j)
{
fields[i]->IProductWRTDerivBase(j, m_viscTensor[j][i], tmp1);
Vmath::Vadd(nCoeffs, tmp1, 1, tmp2[i], 1, tmp2[i], 1);
}
// Evaulate <\phi, \hat{F}\cdot n> - outarray[i]
Vmath::Neg (nCoeffs, tmp2[i], 1);
fields[i]->AddTraceIntegral (viscousFlux[i], tmp2[i]);
fields[i]->SetPhysState (false);
fields[i]->MultiplyByElmtInvMass(tmp2[i], tmp2[i]);
fields[i]->BwdTrans (tmp2[i], outarray[i]);
}
}
virtual Array<OneD, Array<OneD, Array<OneD, NekDouble> > >& Nektar::SolverUtils::DiffusionLDGNS::v_GetFluxTensor ( )
inlineprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 117 of file DiffusionLDGNS.h.

References m_viscTensor.

{
return m_viscTensor;
}
void Nektar::SolverUtils::DiffusionLDGNS::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 53 of file DiffusionLDGNS.cpp.

References m_diffDim, m_gamma, m_gasConstant, m_mu, m_pInf, m_rhoInf, m_session, m_spaceDim, m_thermalConductivity, m_traceNormals, m_traceVel, m_Twall, and m_ViscosityType.

{
m_session = pSession;
m_session->LoadParameter ("Gamma", m_gamma, 1.4);
m_session->LoadParameter ("GasConstant", m_gasConstant, 287.058);
m_session->LoadParameter ("Twall", m_Twall, 300.15);
m_session->LoadSolverInfo("ViscosityType", m_ViscosityType,
"Constant");
m_session->LoadParameter ("mu", m_mu, 1.78e-05);
m_session->LoadParameter ("thermalConductivity",
m_session->LoadParameter ("rhoInf", m_rhoInf, 1.225);
m_session->LoadParameter ("pInf", m_pInf, 101325);
// Setting up the normals
int i;
int nDim = pFields[0]->GetCoordim(0);
int nTracePts = pFields[0]->GetTrace()->GetTotPoints();
m_spaceDim = nDim;
if (pSession->DefinesSolverInfo("HOMOGENEOUS"))
{
}
m_traceVel = Array<OneD, Array<OneD, NekDouble> >(m_spaceDim);
m_traceNormals = Array<OneD, Array<OneD, NekDouble> >(m_spaceDim);
for(i = 0; i < m_spaceDim; ++i)
{
m_traceVel[i] = Array<OneD, NekDouble> (nTracePts, 0.0);
m_traceNormals[i] = Array<OneD, NekDouble> (nTracePts);
}
pFields[0]->GetTrace()->GetNormals(m_traceNormals);
}
void Nektar::SolverUtils::DiffusionLDGNS::v_NumericalFluxO1 ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  numericalFluxO1 
)
protectedvirtual

Builds the numerical flux for the 1st order derivatives.

Definition at line 218 of file DiffusionLDGNS.cpp.

References m_spaceDim, m_traceNormals, Vmath::Svtvp(), v_WeakPenaltyO1(), and Vmath::Vmul().

Referenced by v_Diffuse().

{
int i, j;
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
int nScalars = inarray.num_elements();
int nDim = fields[0]->GetCoordim(0);
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);
}
// Store forwards/backwards space along trace space
Array<OneD, Array<OneD, NekDouble> > Fwd (nScalars);
Array<OneD, Array<OneD, NekDouble> > Bwd (nScalars);
Array<OneD, Array<OneD, NekDouble> > numflux(nScalars);
for (i = 0; i < nScalars; ++i)
{
Fwd[i] = Array<OneD, NekDouble>(nTracePts);
Bwd[i] = Array<OneD, NekDouble>(nTracePts);
numflux[i] = Array<OneD, NekDouble>(nTracePts);
fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
fields[0]->GetTrace()->Upwind(Vn, Fwd[i], Bwd[i], numflux[i]);
}
// Modify the values in case of boundary interfaces
if (fields[0]->GetBndCondExpansions().num_elements())
{
v_WeakPenaltyO1(fields, inarray, numflux);
}
// Splitting the numerical flux into the dimensions
for (j = 0; j < m_spaceDim; ++j)
{
for (i = 0; i < nScalars; ++i)
{
Vmath::Vmul(nTracePts, m_traceNormals[j], 1,
numflux[i], 1, numericalFluxO1[j][i], 1);
}
}
}
void Nektar::SolverUtils::DiffusionLDGNS::v_NumericalFluxO2 ( 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

Build the numerical flux for the 2nd order derivatives.

Definition at line 472 of file DiffusionLDGNS.cpp.

References m_traceNormals, Vmath::Svtvp(), v_WeakPenaltyO2(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by v_Diffuse().

{
int i, j;
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
int nVariables = fields.num_elements();
int nDim = fields[0]->GetCoordim(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);
// 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 Riemann flux
// qflux = \hat{q} \cdot u = q \cdot n
// Notice: i = 1 (first row of the viscous tensor is zero)
for (i = 1; i < nVariables; ++i)
{
qflux[i] = Array<OneD, NekDouble> (nTracePts, 0.0);
for (j = 0; j < nDim; ++j)
{
// Compute qFwd and qBwd value of qfield in position 'ji'
fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
// Get Riemann flux of qflux --> LDG implies upwind
fields[i]->GetTrace()->Upwind(Vn, qBwd, qFwd, qfluxtemp);
// Multiply the Riemann flux by the trace normals
Vmath::Vmul(nTracePts, m_traceNormals[j], 1, qfluxtemp, 1,
qfluxtemp, 1);
// Impose weak boundary condition with flux
if (fields[0]->GetBndCondExpansions().num_elements())
{
v_WeakPenaltyO2(fields, i, j, qfield[j][i], qfluxtemp);
}
// Store the final flux into qflux
Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1,
qflux[i], 1);
}
}
}
virtual void Nektar::SolverUtils::DiffusionLDGNS::v_SetHomoDerivs ( Array< OneD, Array< OneD, NekDouble > > &  deriv)
inlineprotectedvirtual

Reimplemented from Nektar::SolverUtils::Diffusion.

Definition at line 111 of file DiffusionLDGNS.h.

References m_homoDerivs.

{
m_homoDerivs = deriv;
}
void Nektar::SolverUtils::DiffusionLDGNS::v_WeakPenaltyO1 ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  penaltyfluxO1 
)
protectedvirtual

Imposes appropriate bcs for the 1st order derivatives.

Definition at line 274 of file DiffusionLDGNS.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::eWallViscous, m_gamma, m_gasConstant, m_Twall, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vsub(), and Vmath::Zero().

Referenced by v_NumericalFluxO1().

{
int cnt;
int i, j, e;
int id1, id2;
int nBndEdgePts, nBndEdges, nBndRegions;
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
int nScalars = inarray.num_elements();
Array<OneD, NekDouble> tmp1(nTracePts, 0.0);
Array<OneD, NekDouble> tmp2(nTracePts, 0.0);
Array<OneD, NekDouble> Tw(nTracePts, m_Twall);
Array< OneD, Array<OneD, NekDouble > > scalarVariables(nScalars);
Array< OneD, Array<OneD, NekDouble > > uplus(nScalars);
// Extract internal values of the scalar variables for Neumann bcs
for (i = 0; i < nScalars; ++i)
{
scalarVariables[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
uplus[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
fields[i]->ExtractTracePhys(inarray[i], uplus[i]);
}
// Compute boundary conditions for velocities
for (i = 0; i < nScalars-1; ++i)
{
// Note that cnt has to loop on nBndRegions and nBndEdges
// and has to be reset to zero per each equation
cnt = 0;
nBndRegions = fields[i+1]->
GetBndCondExpansions().num_elements();
for (j = 0; j < nBndRegions; ++j)
{
nBndEdges = fields[i+1]->
GetBndCondExpansions()[j]->GetExpSize();
for (e = 0; e < nBndEdges; ++e)
{
nBndEdgePts = fields[i+1]->
GetBndCondExpansions()[j]->GetExp(e)->GetTotPoints();
id1 = fields[i+1]->
GetBndCondExpansions()[j]->GetPhys_Offset(e);
id2 = fields[0]->GetTrace()->
GetPhys_Offset(fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
// Reinforcing bcs for velocity in case of Wall bcs
if (fields[i]->GetBndConditions()[j]->
GetUserDefined() ==
{
Vmath::Zero(nBndEdgePts,
&scalarVariables[i][id2], 1);
}
// Imposing velocity bcs if not Wall
else if (fields[i]->GetBndConditions()[j]->
GetBoundaryConditionType() ==
{
Vmath::Vdiv(nBndEdgePts,
&(fields[i+1]->GetBndCondExpansions()[j]->
UpdatePhys())[id1], 1,
&(fields[0]->GetBndCondExpansions()[j]->
UpdatePhys())[id1], 1,
&scalarVariables[i][id2], 1);
}
// For Dirichlet boundary condition: uflux = u_bcs
if (fields[i]->GetBndConditions()[j]->
GetBoundaryConditionType() ==
{
Vmath::Vcopy(nBndEdgePts,
&scalarVariables[i][id2], 1,
&penaltyfluxO1[i][id2], 1);
}
// For Neumann boundary condition: uflux = u_+
else if ((fields[i]->GetBndConditions()[j])->
GetBoundaryConditionType() ==
{
Vmath::Vcopy(nBndEdgePts,
&uplus[i][id2], 1,
&penaltyfluxO1[i][id2], 1);
}
// Building kinetic energy to be used for T bcs
Vmath::Vmul(nBndEdgePts,
&scalarVariables[i][id2], 1,
&scalarVariables[i][id2], 1,
&tmp1[id2], 1);
Vmath::Smul(nBndEdgePts, 0.5,
&tmp1[id2], 1,
&tmp1[id2], 1);
Vmath::Vadd(nBndEdgePts,
&tmp2[id2], 1,
&tmp1[id2], 1,
&tmp2[id2], 1);
}
}
}
// Compute boundary conditions for temperature
cnt = 0;
nBndRegions = fields[nScalars]->
GetBndCondExpansions().num_elements();
for (j = 0; j < nBndRegions; ++j)
{
nBndEdges = fields[nScalars]->
GetBndCondExpansions()[j]->GetExpSize();
for (e = 0; e < nBndEdges; ++e)
{
nBndEdgePts = fields[nScalars]->
GetBndCondExpansions()[j]->GetExp(e)->GetTotPoints();
id1 = fields[nScalars]->
GetBndCondExpansions()[j]->GetPhys_Offset(e);
id2 = fields[0]->GetTrace()->
GetPhys_Offset(fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
// Imposing Temperature Twall at the wall
if (fields[i]->GetBndConditions()[j]->
GetUserDefined() ==
{
Vmath::Vcopy(nBndEdgePts,
&Tw[0], 1,
&scalarVariables[nScalars-1][id2], 1);
}
// Imposing Temperature through condition on the Energy
// for no wall boundaries (e.g. farfield)
else if (fields[i]->GetBndConditions()[j]->
GetBoundaryConditionType() ==
{
// Divide E by rho
Vmath::Vdiv(nBndEdgePts,
&(fields[nScalars]->
GetBndCondExpansions()[j]->
GetPhys())[id1], 1,
&(fields[0]->
GetBndCondExpansions()[j]->
GetPhys())[id1], 1,
&scalarVariables[nScalars-1][id2], 1);
// Subtract kinetic energy to E/rho
Vmath::Vsub(nBndEdgePts,
&scalarVariables[nScalars-1][id2], 1,
&tmp2[id2], 1,
&scalarVariables[nScalars-1][id2], 1);
// Multiply by constant factor (gamma-1)/R
Vmath::Smul(nBndEdgePts, (m_gamma - 1)/m_gasConstant,
&scalarVariables[nScalars-1][id2], 1,
&scalarVariables[nScalars-1][id2], 1);
}
// For Dirichlet boundary condition: uflux = u_bcs
if (fields[nScalars]->GetBndConditions()[j]->
GetBoundaryConditionType() ==
{
Vmath::Vcopy(nBndEdgePts,
&scalarVariables[nScalars-1][id2], 1,
&penaltyfluxO1[nScalars-1][id2], 1);
}
// For Neumann boundary condition: uflux = u_+
else if ((fields[nScalars]->GetBndConditions()[j])->
GetBoundaryConditionType() ==
{
Vmath::Vcopy(nBndEdgePts,
&uplus[nScalars-1][id2], 1,
&penaltyfluxO1[nScalars-1][id2], 1);
}
}
}
}
void Nektar::SolverUtils::DiffusionLDGNS::v_WeakPenaltyO2 ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const int  var,
const int  dir,
const Array< OneD, const NekDouble > &  qfield,
Array< OneD, NekDouble > &  penaltyflux 
)
protectedvirtual

Imposes appropriate bcs for the 2nd order derivatives.

Definition at line 534 of file DiffusionLDGNS.cpp.

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

Referenced by v_NumericalFluxO2().

{
int cnt = 0;
int nBndEdges, nBndEdgePts;
int i, e;
int id2;
int nTracePts = fields[0]->GetTrace()->GetTotPoints();
int nBndRegions = fields[var]->GetBndCondExpansions().num_elements();
Array<OneD, NekDouble > uterm(nTracePts);
Array<OneD, NekDouble > qtemp(nTracePts);
// Extract the physical values of the solution at the boundaries
fields[var]->ExtractTracePhys(qfield, qtemp);
// Loop on the boundary regions to apply appropriate bcs
for (i = 0; i < nBndRegions; ++i)
{
// Number of boundary regions related to region 'i'
nBndEdges = fields[var]->
GetBndCondExpansions()[i]->GetExpSize();
// Weakly impose bcs by modifying flux values
for (e = 0; e < nBndEdges; ++e)
{
nBndEdgePts = fields[var]->
GetBndCondExpansions()[i]->GetExp(e)->GetTotPoints();
id2 = fields[0]->GetTrace()->
GetPhys_Offset(fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
// In case of Dirichlet bcs:
// uflux = gD
if(fields[var]->GetBndConditions()[i]->
GetBoundaryConditionType() == SpatialDomains::eDirichlet)
{
Vmath::Vmul(nBndEdgePts,
&m_traceNormals[dir][id2], 1,
&qtemp[id2], 1,
&penaltyflux[id2], 1);
}
// 3.4) In case of Neumann bcs:
// uflux = u+
else if((fields[var]->GetBndConditions()[i])->
GetBoundaryConditionType() == SpatialDomains::eNeumann)
{
ASSERTL0(false,
"Neumann bcs not implemented for LDGNS");
/*
Vmath::Vmul(nBndEdgePts,
&m_traceNormals[dir][id2], 1,
&(fields[var]->
GetBndCondExpansions()[i]->
UpdatePhys())[id1], 1,
&penaltyflux[id2], 1);
*/
}
}
}
}

Member Data Documentation

int Nektar::SolverUtils::DiffusionLDGNS::m_diffDim
protected

Definition at line 75 of file DiffusionLDGNS.h.

Referenced by v_Diffuse(), and v_InitObject().

NekDouble Nektar::SolverUtils::DiffusionLDGNS::m_gamma
protected

Definition at line 61 of file DiffusionLDGNS.h.

Referenced by v_InitObject(), and v_WeakPenaltyO1().

NekDouble Nektar::SolverUtils::DiffusionLDGNS::m_gasConstant
protected

Definition at line 62 of file DiffusionLDGNS.h.

Referenced by v_InitObject(), and v_WeakPenaltyO1().

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionLDGNS::m_homoDerivs
protected

Definition at line 72 of file DiffusionLDGNS.h.

Referenced by v_Diffuse(), and v_SetHomoDerivs().

NekDouble Nektar::SolverUtils::DiffusionLDGNS::m_mu
protected

Definition at line 65 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

NekDouble Nektar::SolverUtils::DiffusionLDGNS::m_pInf
protected

Definition at line 68 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

NekDouble Nektar::SolverUtils::DiffusionLDGNS::m_rhoInf
protected

Definition at line 67 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::DiffusionLDGNS::m_session
protected

Definition at line 60 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

int Nektar::SolverUtils::DiffusionLDGNS::m_spaceDim
protected

Definition at line 74 of file DiffusionLDGNS.h.

Referenced by v_Diffuse(), v_InitObject(), and v_NumericalFluxO1().

NekDouble Nektar::SolverUtils::DiffusionLDGNS::m_thermalConductivity
protected

Definition at line 66 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionLDGNS::m_traceNormals
protected
Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::DiffusionLDGNS::m_traceVel
protected

Definition at line 58 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

NekDouble Nektar::SolverUtils::DiffusionLDGNS::m_Twall
protected

Definition at line 63 of file DiffusionLDGNS.h.

Referenced by v_InitObject(), and v_WeakPenaltyO1().

std::string Nektar::SolverUtils::DiffusionLDGNS::m_ViscosityType
protected

Definition at line 64 of file DiffusionLDGNS.h.

Referenced by v_InitObject().

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::SolverUtils::DiffusionLDGNS::m_viscTensor
protected

Definition at line 70 of file DiffusionLDGNS.h.

Referenced by v_Diffuse(), and v_GetFluxTensor().

std::string Nektar::SolverUtils::DiffusionLDGNS::type
static
Initial value:
RegisterCreatorFunction("LDGNS", DiffusionLDGNS::create)

Definition at line 53 of file DiffusionLDGNS.h.