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

Base class for unsteady solvers. More...

#include <UnsteadySystem.h>

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

Public Member Functions

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.

Public Attributes

NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1).

Protected Member Functions

SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members.
virtual SOLVER_UTILS_EXPORT void v_InitObject ()
 Init object for UnsteadySystem class.
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_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters.
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

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 Member Functions

void WeakPenaltyforScalar (const int var, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble time=0.0)
void WeakPenaltyforVector (const int var, const int dir, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble C11, NekDouble time=0.0)

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...

Detailed Description

Base class for unsteady solvers.

Provides the underlying timestepping framework for unsteady solvers including the general timestepping routines. This class is not intended to be directly instantiated, but rather is a base class on which to define unsteady solvers.

For details on implementing unsteady solvers see sectionADRSolverModuleImplementation here

Definition at line 48 of file UnsteadySystem.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::UnsteadySystem::~UnsteadySystem ( )
virtual

Destructor.

Destructor for the class UnsteadyAdvection.

Definition at line 122 of file UnsteadySystem.cpp.

{
}
Nektar::SolverUtils::UnsteadySystem::UnsteadySystem ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Initialises UnsteadySystem class members.

Processes SolverInfo parameters from the session file and sets up timestepping-specific code.

Parameters
pSessionSession object to read parameters from.

Definition at line 65 of file UnsteadySystem.cpp.

: EquationSystem(pSession),
{
}

Member Function Documentation

void Nektar::SolverUtils::UnsteadySystem::CheckForRestartTime ( NekDouble time)
protected

Definition at line 658 of file UnsteadySystem.cpp.

References Nektar::LibUtilities::eFunctionTypeFile, Nektar::iterator, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_fld, Nektar::SolverUtils::EquationSystem::m_session, Nektar::LibUtilities::NullFieldMetaDataMap, and Nektar::LibUtilities::PortablePath().

Referenced by v_DoInitialise().

{
if (m_session->DefinesFunction("InitialConditions"))
{
for (int i = 0; i < m_fields.num_elements(); ++i)
{
vType = m_session->GetFunctionType(
"InitialConditions", m_session->GetVariable(i));
{
std::string filename
= m_session->GetFunctionFilename(
"InitialConditions", m_session->GetVariable(i));
fs::path pfilename(filename);
// redefine path for parallel file which is in directory
if(fs::is_directory(pfilename))
{
fs::path metafile("Info.xml");
fs::path fullpath = pfilename / metafile;
filename = LibUtilities::PortablePath(fullpath);
}
m_fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
// check to see if time defined
{
iter = m_fieldMetaDataMap.find("Time");
if (iter != m_fieldMetaDataMap.end())
{
time = boost::lexical_cast<NekDouble>(
iter->second);
}
}
break;
}
}
}
}
NekDouble Nektar::SolverUtils::UnsteadySystem::GetTimeStep ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray)

Calculate the larger time-step mantaining the problem stable.

Return the timestep to be used for the next step in the time-marching loop.

This function can be overloaded to facilitate solver which utilise a CFL (or other) parameter to determine a maximum timestep under which the problem remains stable.

Definition at line 865 of file UnsteadySystem.cpp.

References v_GetTimeStep().

{
return v_GetTimeStep(inarray);
}
NekDouble Nektar::SolverUtils::UnsteadySystem::MaxTimeStepEstimator ( )
protected

Get the maximum timestep estimator for cfl control.

Returns the maximum time estimator for CFL control.

Definition at line 129 of file UnsteadySystem.cpp.

References ASSERTL0, Nektar::LibUtilities::eAdamsBashforthOrder1, Nektar::LibUtilities::eAdamsBashforthOrder2, Nektar::LibUtilities::eClassicalRungeKutta4, Nektar::LibUtilities::eForwardEuler, Nektar::LibUtilities::eRungeKutta2_ImprovedEuler, Nektar::LibUtilities::eRungeKutta2_ModifiedEuler, and m_intScheme.

{
NekDouble TimeStability = 0.0;
switch(m_intScheme->GetIntegrationMethod())
{
{
TimeStability = 2.784;
break;
}
{
TimeStability = 2.0;
break;
}
{
TimeStability = 1.0;
break;
}
default:
{
false,
"No CFL control implementation for this time"
"integration scheme");
}
}
return TimeStability;
}
void Nektar::SolverUtils::UnsteadySystem::v_AppendOutput1D ( Array< OneD, Array< OneD, NekDouble > > &  solution1D)
protectedvirtual

Print the solution at each solution point in a txt file.

Stores the solution in a file for 1D problems only. This method has been implemented to facilitate the post-processing for 1D problems.

Definition at line 462 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by v_DoSolve().

{
// Coordinates of the quadrature points in the real physical space
Array<OneD,NekDouble> x(GetNpoints());
Array<OneD,NekDouble> y(GetNpoints());
Array<OneD,NekDouble> z(GetNpoints());
m_fields[0]->GetCoords(x, y, z);
// Print out the solution in a txt file
ofstream outfile;
outfile.open("solution1D.txt");
for(int i = 0; i < GetNpoints(); i++)
{
outfile << scientific << setw (17) << setprecision(16) << x[i]
<< " " << solution1D[0][i] << endl;
}
outfile << endl << endl;
outfile.close();
}
void Nektar::SolverUtils::UnsteadySystem::v_DoInitialise ( void  )
protectedvirtual
void Nektar::SolverUtils::UnsteadySystem::v_DoSolve ( void  )
protectedvirtual

Solves an unsteady problem.

Initialises the time integration scheme (as specified in the session file), and perform the time integration.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::PulseWaveSystem.

Definition at line 167 of file UnsteadySystem.cpp.

References ASSERTL0, Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::GetTimeStep(), Nektar::iterator, Nektar::NekConstants::kNekZeroTol, m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_checksteps, Nektar::SolverUtils::EquationSystem::m_checktime, Nektar::SolverUtils::EquationSystem::m_fields, m_filters, Nektar::SolverUtils::EquationSystem::m_fintime, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_homoInitialFwd, m_infosteps, m_intScheme, m_intSoln, m_intVariables, m_ode, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, Nektar::Timer::Start(), Nektar::Timer::Stop(), Nektar::Timer::TimePerTest(), v_AppendOutput1D(), v_PostIntegrate(), and v_PreIntegrate().

{
ASSERTL0(m_intScheme != 0, "No time integration scheme.");
int i, nchk = 1;
int nvariables = 0;
int nfields = m_fields.num_elements();
if (m_intVariables.empty())
{
for (i = 0; i < nfields; ++i)
{
m_intVariables.push_back(i);
}
nvariables = nfields;
}
else
{
nvariables = m_intVariables.size();
}
// Integrate in wave-space if using homogeneous1D
{
for(i = 0; i < nfields; ++i)
{
m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetPhys(),
m_fields[i]->UpdatePhys());
m_fields[i]->SetWaveSpace(true);
m_fields[i]->SetPhysState(false);
}
}
// Set up wrapper to fields data storage.
Array<OneD, Array<OneD, NekDouble> > fields(nvariables);
Array<OneD, Array<OneD, NekDouble> > tmp (nvariables);
// Order storage to list time-integrated fields first.
for(i = 0; i < nvariables; ++i)
{
fields[i] = m_fields[m_intVariables[i]]->GetPhys();
m_fields[m_intVariables[i]]->SetPhysState(false);
}
// Initialise time integration scheme
m_intSoln = m_intScheme->InitializeScheme(
m_timestep, fields, m_time, m_ode);
// Initialise filters
for (x = m_filters.begin(); x != m_filters.end(); ++x)
{
(*x)->Initialise(m_fields, m_time);
}
// Ensure that there is no conflict of parameters
{
// Check final condition
ASSERTL0(m_fintime == 0.0 || m_steps == 0,
"Final condition not unique: "
"fintime > 0.0 and Nsteps > 0");
// Check timestep condition
"Timestep not unique: timestep > 0.0 & CFL > 0.0");
}
// Check uniqueness of checkpoint output
ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
(m_checktime > 0.0 && m_checksteps == 0) ||
(m_checktime == 0.0 && m_checksteps > 0),
"Only one of IO_CheckTime and IO_CheckSteps "
"should be set!");
Timer timer;
bool doCheckTime = false;
int step = 0;
NekDouble intTime = 0.0;
NekDouble lastCheckTime = 0.0;
NekDouble cpuTime = 0.0;
NekDouble elapsed = 0.0;
while (step < m_steps ||
{
{
// Ensure that the final timestep finishes at the final
// time, or at a prescribed IO_CheckTime.
{
}
else if (m_checktime &&
m_time + m_timestep - lastCheckTime >= m_checktime)
{
lastCheckTime += m_checktime;
m_timestep = lastCheckTime - m_time;
doCheckTime = true;
}
}
// Perform any solver-specific pre-integration steps
if (v_PreIntegrate(step))
{
break;
}
timer.Start();
fields = m_intScheme->TimeIntegrate(
timer.Stop();
elapsed = timer.TimePerTest(1);
intTime += elapsed;
cpuTime += elapsed;
// Write out status information
if (m_session->GetComm()->GetRank() == 0 &&
!((step+1) % m_infosteps))
{
cout << "Steps: " << setw(8) << left << step+1 << " "
<< "Time: " << setw(12) << left << m_time;
{
cout << " Time-step: " << setw(12)
<< left << m_timestep;
}
stringstream ss;
ss << cpuTime << "s";
cout << " CPU Time: " << setw(8) << left
<< ss.str() << endl;
cpuTime = 0.0;
}
// Perform any solver-specific post-integration steps
if (v_PostIntegrate(step))
{
break;
}
// Transform data into coefficient space
for (i = 0; i < nvariables; ++i)
{
m_fields[m_intVariables[i]]->SetPhys(fields[i]);
m_fields[m_intVariables[i]]->FwdTrans_IterPerExp(
fields[i],
m_fields[m_intVariables[i]]->UpdateCoeffs());
m_fields[m_intVariables[i]]->SetPhysState(false);
}
// Update filters
for (x = m_filters.begin(); x != m_filters.end(); ++x)
{
(*x)->Update(m_fields, m_time);
}
// Write out checkpoint files
if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
doCheckTime)
{
{
vector<bool> transformed(nfields, false);
for(i = 0; i < nfields; i++)
{
if (m_fields[i]->GetWaveSpace())
{
m_fields[i]->SetWaveSpace(false);
m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
m_fields[i]->UpdatePhys());
m_fields[i]->SetPhysState(true);
transformed[i] = true;
}
}
for(i = 0; i < nfields; i++)
{
if (transformed[i])
{
m_fields[i]->SetWaveSpace(true);
m_fields[i]->HomogeneousFwdTrans(
m_fields[i]->GetPhys(),
m_fields[i]->UpdatePhys());
m_fields[i]->SetPhysState(false);
}
}
}
else
{
}
doCheckTime = false;
}
// Step advance
++step;
}
// Print out summary statistics
if (m_session->GetComm()->GetRank() == 0)
{
if (m_cflSafetyFactor > 0.0)
{
cout << "CFL safety factor : " << m_cflSafetyFactor << endl
<< "CFL time-step : " << m_timestep << endl;
}
cout << "Time-integration : " << intTime << "s" << endl;
}
// If homogeneous, transform back into physical space if necessary.
{
for(i = 0; i < nfields; i++)
{
if (m_fields[i]->GetWaveSpace())
{
m_fields[i]->SetWaveSpace(false);
m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
m_fields[i]->UpdatePhys());
m_fields[i]->SetPhysState(true);
}
}
}
else
{
for(i = 0; i < nvariables; ++i)
{
m_fields[m_intVariables[i]]->SetPhys(fields[i]);
m_fields[m_intVariables[i]]->SetPhysState(true);
}
}
// Finalise filters
for (x = m_filters.begin(); x != m_filters.end(); ++x)
{
(*x)->Finalise(m_fields, m_time);
}
// Print for 1D problems
if(m_spacedim == 1)
{
v_AppendOutput1D(fields);
}
}
void Nektar::SolverUtils::UnsteadySystem::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::EquationSystem.

Reimplemented in Nektar::CompressibleFlowSystem, Nektar::UnsteadyAdvectionDiffusion, Nektar::CFLtester, Nektar::UnsteadyAdvection, Nektar::PulseWavePropagation, Nektar::Monodomain, Nektar::Bidomain, Nektar::ShallowWaterSystem, Nektar::LinearSWE, Nektar::EulerADCFE, Nektar::EulerCFE, Nektar::NonlinearSWE, Nektar::ImageWarpingSystem, Nektar::NavierStokesCFE, and Nektar::UnsteadyDiffusion.

Definition at line 435 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::AddSummaryItem(), Nektar::SolverUtils::EquationSystem::m_checksteps, m_explicitAdvection, m_explicitDiffusion, m_explicitReaction, m_intScheme, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_timestep, and Nektar::LibUtilities::TimeIntegrationMethodMap.

{
AddSummaryItem(s, "Advection",
m_explicitAdvection ? "explicit" : "implicit");
AddSummaryItem(s, "Diffusion",
m_explicitDiffusion ? "explicit" : "implicit");
if (m_session->GetSolverInfo("EQTYPE")
== "SteadyAdvectionDiffusionReaction")
{
AddSummaryItem(s, "Reaction",
m_explicitReaction ? "explicit" : "implicit");
}
AddSummaryItem(s, "Time Step", m_timestep);
AddSummaryItem(s, "No. of Steps", m_steps);
AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
AddSummaryItem(s, "Integration Type",
m_intScheme->GetIntegrationMethod()]);
}
NekDouble Nektar::SolverUtils::UnsteadySystem::v_GetTimeStep ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray)
protectedvirtual

Return the timestep to be used for the next step in the time-marching loop.

See Also
UnsteadySystem::GetTimeStep

Reimplemented in Nektar::CompressibleFlowSystem.

Definition at line 877 of file UnsteadySystem.cpp.

References ASSERTL0.

Referenced by GetTimeStep().

{
ASSERTL0(false, "Not defined for this class");
return 0.0;
}
void Nektar::SolverUtils::UnsteadySystem::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::CoupledLinearNS, Nektar::PulseWaveSystem, Nektar::CompressibleFlowSystem, Nektar::UnsteadyAdvectionDiffusion, Nektar::IncNavierStokes, Nektar::CFLtester, Nektar::UnsteadyAdvection, Nektar::UnsteadyInviscidBurger, Nektar::ShallowWaterSystem, Nektar::EulerADCFE, Nektar::APE, Nektar::EulerCFE, Nektar::NavierStokesCFE, Nektar::PulseWavePropagation, Nektar::ImageWarpingSystem, Nektar::UnsteadyDiffusion, Nektar::Monodomain, Nektar::Bidomain, Nektar::LinearSWE, Nektar::NonlinearSWE, Nektar::PulseWaveSystemOutput, and Nektar::VelocityCorrectionScheme.

Definition at line 76 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::GetFilterFactory(), Nektar::LibUtilities::GetTimeIntegrationWrapperFactory(), m_cflSafetyFactor, m_explicitAdvection, m_explicitDiffusion, m_explicitReaction, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, m_filters, m_homoInitialFwd, m_infosteps, m_intScheme, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::SolverUtils::EquationSystem::m_time.

{
// Load SolverInfo parameters
m_session->MatchSolverInfo("DIFFUSIONADVANCEMENT","Explicit",
m_session->MatchSolverInfo("ADVECTIONADVANCEMENT","Explicit",
m_session->MatchSolverInfo("REACTIONADVANCEMENT", "Explicit",
// For steady problems, we do not initialise the time integration
if (m_session->DefinesSolverInfo("TIMEINTEGRATIONMETHOD"))
{
CreateInstance(m_session->GetSolverInfo(
"TIMEINTEGRATIONMETHOD"));
// Load generic input parameters
m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
// Set up time to be dumped in field information
boost::lexical_cast<std::string>(m_time);
}
// By default attempt to forward transform initial condition.
// Set up filters
LibUtilities::FilterMap::const_iterator x;
for (x = f.begin(); x != f.end(); ++x)
{
m_filters.push_back(GetFilterFactory().CreateInstance(
x->first,
x->second));
}
}
void Nektar::SolverUtils::UnsteadySystem::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::IncNavierStokes, Nektar::CFLtester, Nektar::APE, Nektar::PulseWavePropagation, and Nektar::ImageWarpingSystem.

Definition at line 483 of file UnsteadySystem.cpp.

References ASSERTL0.

{
ASSERTL0(false,
"This function is not implemented for this equation.");
}
void Nektar::SolverUtils::UnsteadySystem::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX,
Array< OneD, Array< OneD, NekDouble > > &  numfluxY 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::APE.

Definition at line 491 of file UnsteadySystem.cpp.

References ASSERTL0.

{
ASSERTL0(false,
"This function is not implemented for this equation.");
}
void Nektar::SolverUtils::UnsteadySystem::v_NumFluxforScalar ( const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  uflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 500 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Vmul(), and WeakPenaltyforScalar().

{
int i, j;
int nTraceNumPoints = GetTraceNpoints();
int nvariables = m_fields.num_elements();
int nqvar = uflux.num_elements();
Array<OneD, NekDouble > Fwd (nTraceNumPoints);
Array<OneD, NekDouble > Bwd (nTraceNumPoints);
Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
Array<OneD, NekDouble > fluxtemp(nTraceNumPoints, 0.0);
// Get the sign of (v \cdot n), v = an arbitrary vector
// Evaulate upwind flux:
// uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
for (j = 0; j < nqvar; ++j)
{
for (i = 0; i < nvariables ; ++i)
{
// Compute Fwd and Bwd value of ufield of i direction
m_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
m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
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(m_fields[0]->GetBndCondExpansions().num_elements())
{
WeakPenaltyforScalar(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(nTraceNumPoints,
fluxtemp, 1,
uflux[j][i], 1);
}
}
}
void Nektar::SolverUtils::UnsteadySystem::v_NumFluxforVector ( const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  qfield,
Array< OneD, Array< OneD, NekDouble > > &  qflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 570 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Smul(), Vmath::Vadd(), Vmath::Vmul(), Vmath::Vsub(), and WeakPenaltyforVector().

{
int nTraceNumPoints = GetTraceNpoints();
int nvariables = m_fields.num_elements();
int nqvar = qfield.num_elements();
NekDouble C11 = 1.0;
Array<OneD, NekDouble > Fwd(nTraceNumPoints);
Array<OneD, NekDouble > Bwd(nTraceNumPoints);
Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
Array<OneD, NekDouble > qFwd (nTraceNumPoints);
Array<OneD, NekDouble > qBwd (nTraceNumPoints);
Array<OneD, NekDouble > qfluxtemp(nTraceNumPoints, 0.0);
Array<OneD, NekDouble > uterm(nTraceNumPoints);
// Evaulate upwind flux:
// qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
for (int i = 0; i < nvariables; ++i)
{
qflux[i] = Array<OneD, NekDouble> (nTraceNumPoints, 0.0);
for (int j = 0; j < nqvar; ++j)
{
// Compute Fwd and Bwd value of ufield of jth direction
m_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+
m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
qBwd, qFwd,
qfluxtemp);
Vmath::Vmul(nTraceNumPoints,
qfluxtemp, 1,
qfluxtemp, 1);
// Generate Stability term = - C11 ( u- - u+ )
m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
Vmath::Vsub(nTraceNumPoints,
Fwd, 1, Bwd, 1,
uterm, 1);
Vmath::Smul(nTraceNumPoints,
-1.0 * C11, uterm, 1,
uterm, 1);
// Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
Vmath::Vadd(nTraceNumPoints,
uterm, 1,
qfluxtemp, 1,
qfluxtemp, 1);
// Imposing weak boundary condition with flux
if (m_fields[0]->GetBndCondExpansions().num_elements())
{
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(nTraceNumPoints,
qfluxtemp, 1,
qflux[i], 1,
qflux[i], 1);
}
}
}
bool Nektar::SolverUtils::UnsteadySystem::v_PostIntegrate ( int  step)
protectedvirtual

Reimplemented in Nektar::IncNavierStokes.

Definition at line 889 of file UnsteadySystem.cpp.

Referenced by v_DoSolve().

{
return false;
}
bool Nektar::SolverUtils::UnsteadySystem::v_PreIntegrate ( int  step)
protectedvirtual

Reimplemented in Nektar::IncNavierStokes.

Definition at line 884 of file UnsteadySystem.cpp.

Referenced by v_DoSolve().

{
return false;
}
void Nektar::SolverUtils::UnsteadySystem::WeakPenaltyforScalar ( const int  var,
const Array< OneD, const NekDouble > &  physfield,
Array< OneD, NekDouble > &  penaltyflux,
NekDouble  time = 0.0 
)
private

Definition at line 706 of file UnsteadySystem.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SolverUtils::EquationSystem::GetPhys_Offset(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, and Vmath::Vcopy().

Referenced by v_NumFluxforScalar().

{
int i, e, npoints, id1, id2;
// Number of boundary regions
int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
int Nfps, numBDEdge;
int nTraceNumPoints = GetTraceNpoints();
int cnt = 0;
Array<OneD, NekDouble > uplus(nTraceNumPoints);
m_fields[var]->ExtractTracePhys(physfield, uplus);
for (i = 0; i < nbnd; ++i)
{
// Number of boundary expansion related to that region
numBDEdge = m_fields[var]->
GetBndCondExpansions()[i]->GetExpSize();
// Evaluate boundary values g_D or g_N from input files
m_session->GetFunction("InitialConditions", 0);
npoints = m_fields[var]->
GetBndCondExpansions()[i]->GetNpoints();
Array<OneD,NekDouble> BDphysics(npoints);
Array<OneD,NekDouble> x0(npoints,0.0);
Array<OneD,NekDouble> x1(npoints,0.0);
Array<OneD,NekDouble> x2(npoints,0.0);
m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
ifunc->Evaluate(x0,x1,x2,time,BDphysics);
// Weakly impose boundary conditions by modifying flux values
for (e = 0; e < numBDEdge ; ++e)
{
// Number of points on the expansion
Nfps = m_fields[var]->
GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
id1 = m_fields[var]->
GetBndCondExpansions()[i]->GetPhys_Offset(e);
id2 = m_fields[0]->GetTrace()->
GetPhys_Offset(m_fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
// For Dirichlet boundary condition: uflux = g_D
if (m_fields[var]->GetBndConditions()[i]->
GetBoundaryConditionType() == SpatialDomains::eDirichlet)
{
Vmath::Vcopy(Nfps,
&BDphysics[id1], 1,
&penaltyflux[id2], 1);
}
// For Neumann boundary condition: uflux = u+
else if ((m_fields[var]->GetBndConditions()[i])->
GetBoundaryConditionType() == SpatialDomains::eNeumann)
{
Vmath::Vcopy(Nfps,
&uplus[id2], 1,
&penaltyflux[id2], 1);
}
}
}
}
void Nektar::SolverUtils::UnsteadySystem::WeakPenaltyforVector ( const int  var,
const int  dir,
const Array< OneD, const NekDouble > &  physfield,
Array< OneD, NekDouble > &  penaltyflux,
NekDouble  C11,
NekDouble  time = 0.0 
)
private

Diffusion: Imposing weak boundary condition for q with flux uflux = g_D on Dirichlet boundary condition uflux = u_Fwd on Neumann boundary condition

Definition at line 783 of file UnsteadySystem.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SolverUtils::EquationSystem::GetPhys_Offset(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_traceNormals, and Vmath::Vmul().

Referenced by v_NumFluxforVector().

{
int i, e, npoints, id1, id2;
int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
int numBDEdge, Nfps;
int nTraceNumPoints = GetTraceNpoints();
Array<OneD, NekDouble > uterm(nTraceNumPoints);
Array<OneD, NekDouble > qtemp(nTraceNumPoints);
int cnt = 0;
m_fields[var]->ExtractTracePhys(physfield,qtemp);
for (i = 0; i < nbnd; ++i)
{
numBDEdge = m_fields[var]->
GetBndCondExpansions()[i]->GetExpSize();
// Evaluate boundary values g_D or g_N from input files
m_session->GetFunction("InitialConditions", 0);
npoints = m_fields[var]->
GetBndCondExpansions()[i]->GetNpoints();
Array<OneD,NekDouble> BDphysics(npoints);
Array<OneD,NekDouble> x0(npoints,0.0);
Array<OneD,NekDouble> x1(npoints,0.0);
Array<OneD,NekDouble> x2(npoints,0.0);
m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
ifunc->Evaluate(x0,x1,x2,time,BDphysics);
// Weakly impose boundary conditions by modifying flux values
for (e = 0; e < numBDEdge ; ++e)
{
Nfps = m_fields[var]->
GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
id1 = m_fields[var]->
GetBndCondExpansions()[i]->GetPhys_Offset(e);
id2 = m_fields[0]->GetTrace()->
GetPhys_Offset(m_fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
// For Dirichlet boundary condition:
//qflux = q+ - C_11 (u+ - g_D) (nx, ny)
if(m_fields[var]->GetBndConditions()[i]->
GetBoundaryConditionType() == SpatialDomains::eDirichlet)
{
Vmath::Vmul(Nfps,
&m_traceNormals[dir][id2], 1,
&qtemp[id2], 1,
&penaltyflux[id2], 1);
}
// For Neumann boundary condition: qflux = g_N
else if((m_fields[var]->GetBndConditions()[i])->
GetBoundaryConditionType() == SpatialDomains::eNeumann)
{
&m_traceNormals[dir][id2], 1,
&BDphysics[id1], 1,
&penaltyflux[id2], 1);
}
}
}
}

Member Data Documentation

NekDouble Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor

CFL safety factor (comprise between 0 to 1).

Definition at line 59 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

NekDouble Nektar::SolverUtils::UnsteadySystem::m_epsilon
protected

Definition at line 71 of file UnsteadySystem.h.

bool Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection
protected
bool Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion
protected
bool Nektar::SolverUtils::UnsteadySystem::m_explicitReaction
protected

Indicates if explicit or implicit treatment of reaction is used.

Definition at line 77 of file UnsteadySystem.h.

Referenced by v_GenerateSummary(), and v_InitObject().

std::vector<FilterSharedPtr> Nektar::SolverUtils::UnsteadySystem::m_filters
protected
bool Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd
protected
int Nektar::SolverUtils::UnsteadySystem::m_infosteps
protected
LibUtilities::TimeIntegrationWrapperSharedPtr Nektar::SolverUtils::UnsteadySystem::m_intScheme
protected
LibUtilities::TimeIntegrationSolutionSharedPtr Nektar::SolverUtils::UnsteadySystem::m_intSoln
protected
std::vector<int> Nektar::SolverUtils::UnsteadySystem::m_intVariables
protected
LibUtilities::TimeIntegrationSchemeOperators Nektar::SolverUtils::UnsteadySystem::m_ode
protected