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

#include <UnsteadyDiffusion.h>

Inheritance diagram for Nektar::UnsteadyDiffusion:
Inheritance graph
[legend]
Collaboration diagram for Nektar::UnsteadyDiffusion:
Collaboration graph
[legend]

Public Member Functions

virtual ~UnsteadyDiffusion ()
 Destructor.
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object.
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
template<class T >
boost::shared_ptr< T > as ()
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name.
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available.
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow.
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $.
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection.
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection.
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields.
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point.
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time.
SOLVER_UTILS_EXPORT int GetNcoeffs ()
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
SOLVER_UTILS_EXPORT int GetNumExpModes ()
SOLVER_UTILS_EXPORT const
Array< OneD, int > 
GetNumExpModesPerExp ()
SOLVER_UTILS_EXPORT int GetNvariables ()
SOLVER_UTILS_EXPORT const
std::string 
GetVariable (unsigned int i)
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
SOLVER_UTILS_EXPORT int GetExpSize ()
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
SOLVER_UTILS_EXPORT int GetTotPoints ()
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
SOLVER_UTILS_EXPORT int GetNpoints ()
SOLVER_UTILS_EXPORT int GetNumElmVelocity ()
SOLVER_UTILS_EXPORT int GetSteps ()
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
SOLVER_UTILS_EXPORT void SetStepsToOne ()
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
SOLVER_UTILS_EXPORT void FwdTransFields ()
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison.

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 Creates an instance of this class.

Static Public Attributes

static std::string className = GetEquationSystemFactory().RegisterCreatorFunction("UnsteadyDiffusion", UnsteadyDiffusion::create)
 Name of class.

Protected Member Functions

virtual void v_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters.
 UnsteadyDiffusion (const LibUtilities::SessionReaderSharedPtr &pSession)
virtual void v_InitObject ()
 Initialisation object for the unsteady diffusion problem.
void GetFluxVector (const int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
 Return the flux vector for the unsteady diffusion problem.
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the projection for the unsteady diffusion problem.
virtual void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
 Implicit solution of the unsteady diffusion problem.
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control.
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop.
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members.
int nocase_cmp (const string &s1, const string &s2)
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time.
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution.
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution.
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space.
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
virtual SOLVER_UTILS_EXPORT void v_Output (void)
virtual SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
v_GetPressure (void)
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)

Protected Attributes

bool m_useSpecVanVisc
NekDouble m_sVVCutoffRatio
NekDouble m_sVVDiffCoeff
SolverUtils::DiffusionSharedPtr m_diffusion
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme.
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
NekDouble m_epsilon
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used.
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used.
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used.
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state.
std::vector< int > m_intVariables
std::vector< FilterSharedPtrm_filters
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
std::string m_sessionName
 Name of the session.
NekDouble m_time
 Current time of simulation.
NekDouble m_fintime
 Finish time of the simulation.
NekDouble m_timestep
 Time step size.
NekDouble m_lambda
 Lambda constant in real system if one required.
NekDouble m_checktime
 Time between checkpoints.
int m_steps
 Number of steps to take.
int m_checksteps
 Number of steps between checkpoints.
int m_spacedim
 Spatial dimension (>= expansion dim).
int m_expdim
 Expansion dimension.
bool m_SingleMode
 Flag to determine if single homogeneous mode is used.
bool m_HalfMode
 Flag to determine if half homogeneous mode is used.
bool m_MultipleModes
 Flag to determine if use multiple homogenenous modes are used.
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
enum HomogeneousType m_HomogeneousType
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
int m_npointsX
 number of points in X direction (if homogeneous)
int m_npointsY
 number of points in Y direction (if homogeneous)
int m_npointsZ
 number of points in Z direction (if homogeneous)
int m_HomoDirec
 number of homogenous directions
int m_NumMode
 Mode to use in case of single mode analysis.

Private Attributes

NekDouble m_waveFreq
NekDouble m_epsilon
StdRegions::VarCoeffMap m_varcoeff

Friends

class MemoryManager< UnsteadyDiffusion >

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1).
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...

Detailed Description

Definition at line 46 of file UnsteadyDiffusion.h.

Constructor & Destructor Documentation

Nektar::UnsteadyDiffusion::~UnsteadyDiffusion ( )
virtual

Destructor.

Unsteady diffusion problem destructor.

Definition at line 133 of file UnsteadyDiffusion.cpp.

{
}
Nektar::UnsteadyDiffusion::UnsteadyDiffusion ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 44 of file UnsteadyDiffusion.cpp.

: UnsteadySystem(pSession)
{
}

Member Function Documentation

static EquationSystemSharedPtr Nektar::UnsteadyDiffusion::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file UnsteadyDiffusion.h.

void Nektar::UnsteadyDiffusion::DoImplicitSolve ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
NekDouble  time,
NekDouble  lambda 
)
protectedvirtual

Implicit solution of the unsteady diffusion problem.

Definition at line 223 of file UnsteadyDiffusion.cpp.

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorTau, m_epsilon, Nektar::SolverUtils::EquationSystem::m_fields, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, m_varcoeff, Nektar::NullFlagList, and Vmath::Smul().

Referenced by v_InitObject().

{
int nvariables = inarray.num_elements();
int npoints = m_fields[0]->GetNpoints();
factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
factors[StdRegions::eFactorTau] = 1.0;
{
}
// We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
// inarray = input: \hat{rhs} -> output: \hat{Y}
// outarray = output: nabla^2 \hat{Y}
// where \hat = modal coeffs
for (int i = 0; i < nvariables; ++i)
{
// Multiply 1.0/timestep/lambda
Vmath::Smul(npoints,
inarray[i], 1,
m_fields[i]->UpdatePhys(), 1);
// Solve a system of equations with Helmholtz solver
m_fields[i]->HelmSolve(m_fields[i]->GetPhys(),
m_fields[i]->UpdateCoeffs(),
factors,
m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
m_fields[i]->UpdatePhys());
m_fields[i]->SetPhysState(false);
// The solution is Y[i]
outarray[i] = m_fields[i]->GetPhys();
}
}
void Nektar::UnsteadyDiffusion::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Compute the projection for the unsteady diffusion problem.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 178 of file UnsteadyDiffusion.cpp.

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Vcopy().

Referenced by v_InitObject().

{
int i;
int nvariables = inarray.num_elements();
{
{
// Just copy over array
int npoints = GetNpoints();
for(i = 0; i < nvariables; ++i)
{
Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
}
break;
}
{
Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
for(i = 0; i < nvariables; ++i)
{
m_fields[i]->FwdTrans(inarray[i], coeffs);
m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
}
break;
}
default:
{
ASSERTL0(false, "Unknown projection scheme");
break;
}
}
}
void Nektar::UnsteadyDiffusion::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 156 of file UnsteadyDiffusion.cpp.

References m_diffusion, and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by v_InitObject().

{
// Number of fields (variables of the problem)
int nVariables = inarray.num_elements();
// RHS computation using the new advection base class
m_diffusion->Diffuse(nVariables,
inarray,
outarray);
}
void Nektar::UnsteadyDiffusion::GetFluxVector ( const int  i,
const int  j,
const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  derivatives,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protected

Return the flux vector for the unsteady diffusion problem.

Definition at line 274 of file UnsteadyDiffusion.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), Vmath::Vcopy(), and Vmath::Zero().

Referenced by v_InitObject().

{
for(int k = 0; k < flux.num_elements(); ++k)
{
Vmath::Zero(GetNpoints(), flux[k], 1);
}
Vmath::Vcopy(GetNpoints(), physfield[i], 1, flux[j], 1);
}
void Nektar::UnsteadyDiffusion::v_GenerateSummary ( SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Prints a summary with some information regards the time-stepping.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 137 of file UnsteadyDiffusion.cpp.

References Nektar::SolverUtils::AddSummaryItem(), m_sVVCutoffRatio, m_sVVDiffCoeff, and m_useSpecVanVisc.

{
{
stringstream ss;
ss << "SVV (cut off = " << m_sVVCutoffRatio
<< ", coeff = " << m_sVVDiffCoeff << ")";
AddSummaryItem(s, "Smoothing", ss.str());
}
}
void Nektar::UnsteadyDiffusion::v_InitObject ( )
protectedvirtual

Initialisation object for the unsteady diffusion problem.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 53 of file UnsteadyDiffusion.cpp.

References ASSERTL0, Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoImplicitSolve(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::SolverUtils::GetDiffusionFactory(), GetFluxVector(), m_diffusion, m_epsilon, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_session, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, m_varcoeff, and m_waveFreq.

{
m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
m_session->LoadParameter("epsilon", m_epsilon, 0.0);
m_session->MatchSolverInfo(
"SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
{
m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
}
int npoints = m_fields[0]->GetNpoints();
if(m_session->DefinesParameter("d00"))
{
= Array<OneD, NekDouble>(npoints, m_session->GetParameter("d00"));
}
if(m_session->DefinesParameter("d11"))
{
= Array<OneD, NekDouble>(npoints, m_session->GetParameter("d11"));
}
if(m_session->DefinesParameter("d22"))
{
= Array<OneD, NekDouble>(npoints, m_session->GetParameter("d22"));
}
{
{
std::string diffName;
// Do not forwards transform initial condition
m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
CreateInstance(diffName, diffName);
m_diffusion->SetFluxVector(&UnsteadyDiffusion::
GetFluxVector, this);
break;
}
{
// In case of Galerkin explicit diffusion gives an error
{
ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
}
// In case of Galerkin implicit diffusion: do nothing
}
}
{
}
else
{
}
}

Friends And Related Function Documentation

friend class MemoryManager< UnsteadyDiffusion >
friend

Definition at line 49 of file UnsteadyDiffusion.h.

Member Data Documentation

string Nektar::UnsteadyDiffusion::className = GetEquationSystemFactory().RegisterCreatorFunction("UnsteadyDiffusion", UnsteadyDiffusion::create)
static

Name of class.

Definition at line 59 of file UnsteadyDiffusion.h.

SolverUtils::DiffusionSharedPtr Nektar::UnsteadyDiffusion::m_diffusion
protected

Definition at line 68 of file UnsteadyDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

NekDouble Nektar::UnsteadyDiffusion::m_epsilon
private

Definition at line 100 of file UnsteadyDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

SolverUtils::RiemannSolverSharedPtr Nektar::UnsteadyDiffusion::m_riemannSolver
protected

Definition at line 69 of file UnsteadyDiffusion.h.

NekDouble Nektar::UnsteadyDiffusion::m_sVVCutoffRatio
protected

Definition at line 66 of file UnsteadyDiffusion.h.

Referenced by DoImplicitSolve(), v_GenerateSummary(), and v_InitObject().

NekDouble Nektar::UnsteadyDiffusion::m_sVVDiffCoeff
protected

Definition at line 67 of file UnsteadyDiffusion.h.

Referenced by DoImplicitSolve(), v_GenerateSummary(), and v_InitObject().

bool Nektar::UnsteadyDiffusion::m_useSpecVanVisc
protected

Definition at line 65 of file UnsteadyDiffusion.h.

Referenced by DoImplicitSolve(), v_GenerateSummary(), and v_InitObject().

StdRegions::VarCoeffMap Nektar::UnsteadyDiffusion::m_varcoeff
private

Definition at line 101 of file UnsteadyDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

NekDouble Nektar::UnsteadyDiffusion::m_waveFreq
private

Definition at line 99 of file UnsteadyDiffusion.h.

Referenced by v_InitObject().