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

#include <UnsteadyAdvectionDiffusion.h>

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

Public Member Functions

virtual ~UnsteadyAdvectionDiffusion ()
 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.
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 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
 Name of class.

Protected Member Functions

 UnsteadyAdvectionDiffusion (const LibUtilities::SessionReaderSharedPtr &pSession)
 Session reader.
void GetFluxVectorAdv (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Evaluate the flux at each solution point for the advection part.
void GetFluxVectorDiff (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)
 Evaluate the flux at each solution point for the diffusion part.
virtual void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the RHS.
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform the projection.
virtual void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
 Solve implicitly the diffusion term.
Array< OneD, NekDouble > & GetNormalVelocity ()
 Get the normal velocity.
virtual void v_InitObject ()
 Initialise the object.
virtual void v_GenerateSummary (SummaryList &s)
 Print Summary.
- 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

SolverUtils::RiemannSolverSharedPtr m_riemannSolver
SolverUtils::AdvectionSharedPtr m_advection
SolverUtils::DiffusionSharedPtr m_diffusion
Array< OneD, Array< OneD,
NekDouble > > 
m_velocity
Array< OneD, NekDoublem_traceVn
int m_planeNumber
- 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_filename
 Filename.
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

Friends

class MemoryManager< UnsteadyAdvectionDiffusion >

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 48 of file UnsteadyAdvectionDiffusion.h.

Constructor & Destructor Documentation

Nektar::UnsteadyAdvectionDiffusion::~UnsteadyAdvectionDiffusion ( )
virtual

Destructor.

Unsteady linear advection diffusion equation destructor.

Definition at line 150 of file UnsteadyAdvectionDiffusion.cpp.

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

Session reader.

Definition at line 46 of file UnsteadyAdvectionDiffusion.cpp.

References m_planeNumber.

: UnsteadySystem(pSession)
{
}

Member Function Documentation

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

Creates an instance of this class.

Definition at line 54 of file UnsteadyAdvectionDiffusion.h.

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

Solve implicitly the diffusion term.

Definition at line 292 of file UnsteadyAdvectionDiffusion.cpp.

References Nektar::StdRegions::eFactorLambda, m_epsilon, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::NullFlagList, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Smul().

Referenced by v_InitObject().

{
int nvariables = inarray.num_elements();
int nq = m_fields[0]->GetNpoints();
factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsilon;
Array<OneD, Array< OneD, NekDouble> > F(nvariables);
F[0] = Array<OneD, NekDouble> (nq*nvariables);
for (int n = 1; n < nvariables; ++n)
{
F[n] = F[n-1] + nq;
}
// 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
inarray[i], 1, F[i], 1);
}
//Setting boundary conditions
for (int i = 0; i < nvariables; ++i)
{
// Solve a system of equations with Helmholtz solver
m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
NullFlagList, factors);
m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
}
}
void Nektar::UnsteadyAdvectionDiffusion::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Perform the projection.

Compute the projection for the unsteady advection diffusion problem.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 243 of file UnsteadyAdvectionDiffusion.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::UnsteadyAdvectionDiffusion::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protectedvirtual

Compute the RHS.

Compute the right-hand side for the unsteady linear advection diffusion problem.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 193 of file UnsteadyAdvectionDiffusion.cpp.

References Nektar::MultiRegions::eDiscontinuous, Nektar::SolverUtils::EquationSystem::GetNpoints(), m_advection, m_diffusion, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, m_velocity, Vmath::Neg(), and Vmath::Vadd().

Referenced by v_InitObject().

{
// Number of fields (variables of the problem)
int nVariables = inarray.num_elements();
// Number of solution points
int nSolutionPts = GetNpoints();
Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
for (int i = 0; i < nVariables; ++i)
{
outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
}
// RHS computation using the new advection base class
m_advection->Advect(nVariables, m_fields, m_velocity,
inarray, outarray);
// Negate the RHS
for (int i = 0; i < nVariables; ++i)
{
Vmath::Neg(nSolutionPts, outarray[i], 1);
}
// No explicit diffusion for CG
{
m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
for (int i = 0; i < nVariables; ++i)
{
Vmath::Vadd(nSolutionPts, &outarray[i][0], 1,
&outarrayDiff[i][0], 1, &outarray[i][0], 1);
}
}
}
void Nektar::UnsteadyAdvectionDiffusion::GetFluxVectorAdv ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Evaluate the flux at each solution point for the advection part.

Return the flux vector for the advection part.

Parameters
physfieldFields.
fluxResulting flux.

Definition at line 342 of file UnsteadyAdvectionDiffusion.cpp.

References ASSERTL1, Nektar::SolverUtils::EquationSystem::m_fields, m_velocity, and Vmath::Vmul().

Referenced by v_InitObject().

{
ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
"Dimension of flux array and velocity array do not match");
const int nq = m_fields[0]->GetNpoints();
for (int i = 0; i < flux.num_elements(); ++i)
{
for (int j = 0; j < flux[0].num_elements(); ++j)
{
Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
flux[i][j], 1);
}
}
}
void Nektar::UnsteadyAdvectionDiffusion::GetFluxVectorDiff ( 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

Evaluate the flux at each solution point for the diffusion part.

Return the flux vector for the diffusion part.

Parameters
iEquation number.
jSpatial direction.
physfieldFields.
derivativesFirst order derivatives.
fluxResulting flux.

Definition at line 370 of file UnsteadyAdvectionDiffusion.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);
}
Array< OneD, NekDouble > & Nektar::UnsteadyAdvectionDiffusion::GetNormalVelocity ( )
protected

Get the normal velocity.

Get the normal velocity for the unsteady linear advection diffusion equation.

Definition at line 158 of file UnsteadyAdvectionDiffusion.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, m_traceVn, m_velocity, Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_InitObject().

{
// Number of trace (interface) points
int i;
int nTracePts = GetTraceNpoints();
// Auxiliary variable to compute the normal velocity
Array<OneD, NekDouble> tmp(nTracePts);
m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
// Reset the normal velocity
Vmath::Zero(nTracePts, m_traceVn, 1);
for (i = 0; i < m_velocity.num_elements(); ++i)
{
m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
Vmath::Vvtvp(nTracePts,
tmp, 1,
m_traceVn, 1);
}
return m_traceVn;
}
void Nektar::UnsteadyAdvectionDiffusion::v_GenerateSummary ( SummaryList s)
protectedvirtual

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 384 of file UnsteadyAdvectionDiffusion.cpp.

void Nektar::UnsteadyAdvectionDiffusion::v_InitObject ( )
protectedvirtual

Initialise the object.

Initialisation object for the unsteady linear advection diffusion equation.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 57 of file UnsteadyAdvectionDiffusion.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::SolverUtils::EquationSystem::EvaluateFunction(), Nektar::SolverUtils::GetAdvectionFactory(), Nektar::SolverUtils::GetDiffusionFactory(), GetFluxVectorAdv(), GetFluxVectorDiff(), GetNormalVelocity(), Nektar::SolverUtils::GetRiemannSolverFactory(), m_advection, 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, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, and m_waveFreq.

{
m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
m_session->LoadParameter("epsilon", m_epsilon, 0.0);
// Define Velocity fields
m_velocity = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
std::vector<std::string> vel;
vel.push_back("Vx");
vel.push_back("Vy");
vel.push_back("Vz");
vel.resize(m_spacedim);
EvaluateFunction(vel, m_velocity, "AdvectionVelocity");
// Type of advection and diffusion classes to be used
{
// Discontinuous field
{
// Do not forwards transform initial condition
// Advection term
string advName;
string riemName;
m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
CreateInstance(advName, advName);
m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
CreateInstance(riemName);
m_advection->SetRiemannSolver(m_riemannSolver);
m_advection->InitObject (m_session, m_fields);
// Diffusion term
std::string diffName;
m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
CreateInstance(diffName, diffName);
break;
}
// Continuous field
{
// Advection term
std::string advName;
m_session->LoadSolverInfo("AdvectionType", advName,
"NonConservative");
CreateInstance(advName, advName);
GetFluxVectorAdv, this);
// 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
break;
}
default:
{
ASSERTL0(false, "Unsupported projection type.");
break;
}
}
{
}
}

Friends And Related Function Documentation

friend class MemoryManager< UnsteadyAdvectionDiffusion >
friend

Definition at line 51 of file UnsteadyAdvectionDiffusion.h.

Member Data Documentation

string Nektar::UnsteadyAdvectionDiffusion::className
static
Initial value:
RegisterCreatorFunction("UnsteadyAdvectionDiffusion",

Name of class.

Definition at line 63 of file UnsteadyAdvectionDiffusion.h.

SolverUtils::AdvectionSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_advection
protected

Definition at line 70 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

SolverUtils::DiffusionSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_diffusion
protected

Definition at line 71 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_epsilon
private

Definition at line 126 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

int Nektar::UnsteadyAdvectionDiffusion::m_planeNumber
protected

Definition at line 77 of file UnsteadyAdvectionDiffusion.h.

Referenced by UnsteadyAdvectionDiffusion().

SolverUtils::RiemannSolverSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_riemannSolver
protected

Definition at line 69 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject().

Array<OneD, NekDouble> Nektar::UnsteadyAdvectionDiffusion::m_traceVn
protected

Definition at line 73 of file UnsteadyAdvectionDiffusion.h.

Referenced by GetNormalVelocity().

Array<OneD, Array<OneD, NekDouble> > Nektar::UnsteadyAdvectionDiffusion::m_velocity
protected
NekDouble Nektar::UnsteadyAdvectionDiffusion::m_waveFreq
private

Definition at line 125 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject().