Nektar++
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:
[legend]

Public Member Functions

virtual ~UnsteadyDiffusion ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
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. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
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. More...
 
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. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
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 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 SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

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

Protected Member Functions

virtual void v_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters. More...
 
 UnsteadyDiffusion (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual void v_InitObject ()
 Initialisation object for the unsteady diffusion problem. More...
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
 Return the flux vector for the unsteady diffusion problem. More...
 
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. More...
 
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. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. More...
 
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. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck ()
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
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. More...
 
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. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
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)
 
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. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda =0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm =-1.0
 
int m_TotLinItePerStep =0
 
int m_StagesPerStep =1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts =0
 
int m_TotLinIts =0
 
int m_TotImpStages =0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
bool m_root
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_timestepMax = -1.0
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

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). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 45 of file UnsteadyDiffusion.h.

Constructor & Destructor Documentation

◆ ~UnsteadyDiffusion()

Nektar::UnsteadyDiffusion::~UnsteadyDiffusion ( )
virtual

Destructor.

Unsteady diffusion problem destructor.

Definition at line 139 of file UnsteadyDiffusion.cpp.

140  {
141  }

◆ UnsteadyDiffusion()

Nektar::UnsteadyDiffusion::UnsteadyDiffusion ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
protected

Definition at line 48 of file UnsteadyDiffusion.cpp.

51  : UnsteadySystem(pSession, pGraph)
52  {
53  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ create()

static EquationSystemSharedPtr Nektar::UnsteadyDiffusion::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 51 of file UnsteadyDiffusion.h.

54  {
56  ::AllocateSharedPtr(pSession, pGraph);
57  p->InitObject();
58  return p;
59  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ DoImplicitSolve()

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 231 of file UnsteadyDiffusion.cpp.

236  {
237  boost::ignore_unused(time);
238 
240 
241  int nvariables = inarray.size();
242  int npoints = m_fields[0]->GetNpoints();
243  factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
244  factors[StdRegions::eFactorTau] = 1.0;
245 
246  if(m_useSpecVanVisc)
247  {
250  }
251 
252  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
253  // inarray = input: \hat{rhs} -> output: \hat{Y}
254  // outarray = output: nabla^2 \hat{Y}
255  // where \hat = modal coeffs
256  for (int i = 0; i < nvariables; ++i)
257  {
258  // Multiply 1.0/timestep/lambda
259  Vmath::Smul(npoints,
260  -factors[StdRegions::eFactorLambda],
261  inarray[i], 1,
262  outarray[i], 1);
263 
264  // Solve a system of equations with Helmholtz solver
265  m_fields[i]->HelmSolve(outarray[i],
266  m_fields[i]->UpdateCoeffs(),
267  factors,
268  m_varcoeff);
269 
270  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
271  outarray[i]);
272 
273  m_fields[i]->SetPhysState(false);
274  }
275  }
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
StdRegions::VarCoeffMap m_varcoeff
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225

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, and Vmath::Smul().

Referenced by v_InitObject().

◆ DoOdeProjection()

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 186 of file UnsteadyDiffusion.cpp.

190  {
191  int i;
192  int nvariables = inarray.size();
193  SetBoundaryConditions(time);
194 
195  switch(m_projectionType)
196  {
198  {
199  // Just copy over array
200  int npoints = GetNpoints();
201 
202  for(i = 0; i < nvariables; ++i)
203  {
204  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
205  }
206  break;
207  }
210  {
211  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
212 
213  for(i = 0; i < nvariables; ++i)
214  {
215  m_fields[i]->FwdTrans(inarray[i], coeffs);
216  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
217  }
218  break;
219  }
220  default:
221  {
222  ASSERTL0(false, "Unknown projection scheme");
223  break;
224  }
225  }
226  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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().

◆ DoOdeRhs()

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 162 of file UnsteadyDiffusion.cpp.

166  {
167  boost::ignore_unused(time);
168 
169  // Number of fields (variables of the problem)
170  int nVariables = inarray.size();
171 
172  // RHS computation using the new advection base class
173  m_diffusion->Diffuse(nVariables,
174  m_fields,
175  inarray,
176  outarray);
177  }
SolverUtils::DiffusionSharedPtr m_diffusion

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

Referenced by v_InitObject().

◆ GetFluxVector()

void Nektar::UnsteadyDiffusion::GetFluxVector ( const Array< OneD, Array< OneD, NekDouble > > &  inarray,
const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  qfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  viscousTensor 
)
protected

Return the flux vector for the unsteady diffusion problem.

Definition at line 280 of file UnsteadyDiffusion.cpp.

284  {
285  boost::ignore_unused(inarray);
286 
287  unsigned int nDim = qfield.size();
288  unsigned int nConvectiveFields = qfield[0].size();
289  unsigned int nPts = qfield[0][0].size();
290 
291  for (unsigned int j = 0; j < nDim; ++j)
292  {
293  for (unsigned int i = 0; i < nConvectiveFields; ++i)
294  {
295  Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1,
296  viscousTensor[j][i], 1 );
297  }
298  }
299  }

References m_epsilon, and Vmath::Smul().

Referenced by v_InitObject().

◆ v_GenerateSummary()

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 143 of file UnsteadyDiffusion.cpp.

144  {
146  if(m_useSpecVanVisc)
147  {
148  stringstream ss;
149  ss << "SVV (cut off = " << m_sVVCutoffRatio
150  << ", coeff = " << m_sVVDiffCoeff << ")";
151  AddSummaryItem(s, "Smoothing", ss.str());
152  }
153  }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47

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

◆ v_InitObject()

void Nektar::UnsteadyDiffusion::v_InitObject ( )
protectedvirtual

Initialisation object for the unsteady diffusion problem.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 58 of file UnsteadyDiffusion.cpp.

59  {
61 
62  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
63  m_session->LoadParameter("epsilon", m_epsilon, 1.0);
64 
65  m_session->MatchSolverInfo(
66  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
67 
69  {
70  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
71  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
72  }
73 
74  int npoints = m_fields[0]->GetNpoints();
75 
76  if(m_session->DefinesParameter("d00"))
77  {
79  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d00"));
80  }
81  if(m_session->DefinesParameter("d11"))
82  {
84  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d11"));
85  }
86  if(m_session->DefinesParameter("d22"))
87  {
89  = Array<OneD, NekDouble>(npoints, m_session->GetParameter("d22"));
90  }
91 
92  switch (m_projectionType)
93  {
95  {
96  std::string diffName;
97 
98  // Do not forwards transform initial condition
99  m_homoInitialFwd = false;
100 
101  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
103  CreateInstance(diffName, diffName);
104  m_diffusion->SetFluxVector(&UnsteadyDiffusion::
105  GetFluxVector, this);
106  m_diffusion->InitObject(m_session, m_fields);
107  break;
108  }
109 
112  {
113  // In case of Galerkin explicit diffusion gives an error
115  {
116  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
117  }
118  // In case of Galerkin implicit diffusion: do nothing
119  }
120  }
121 
123  {
126  }
127  else
128  {
133  }
134  }
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
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.
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
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)
UnsteadyDiffusion(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
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.
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41

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, m_waveFreq, and Nektar::SolverUtils::UnsteadySystem::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< UnsteadyDiffusion >

friend class MemoryManager< UnsteadyDiffusion >
friend

Definition at line 1 of file UnsteadyDiffusion.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 61 of file UnsteadyDiffusion.h.

◆ m_diffusion

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

Definition at line 70 of file UnsteadyDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_epsilon

NekDouble Nektar::UnsteadyDiffusion::m_epsilon
private

Definition at line 101 of file UnsteadyDiffusion.h.

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

◆ m_riemannSolver

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

Definition at line 71 of file UnsteadyDiffusion.h.

◆ m_sVVCutoffRatio

NekDouble Nektar::UnsteadyDiffusion::m_sVVCutoffRatio
protected

Definition at line 68 of file UnsteadyDiffusion.h.

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

◆ m_sVVDiffCoeff

NekDouble Nektar::UnsteadyDiffusion::m_sVVDiffCoeff
protected

Definition at line 69 of file UnsteadyDiffusion.h.

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

◆ m_useSpecVanVisc

bool Nektar::UnsteadyDiffusion::m_useSpecVanVisc
protected

Definition at line 67 of file UnsteadyDiffusion.h.

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

◆ m_varcoeff

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

Definition at line 102 of file UnsteadyDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_waveFreq

NekDouble Nektar::UnsteadyDiffusion::m_waveFreq
private

Definition at line 100 of file UnsteadyDiffusion.h.

Referenced by v_InitObject().