Nektar++
Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | Friends | List of all members
Nektar::MMFDiffusion Class Reference

A model for cardiac conduction. More...

#include <MMFDiffusion.h>

Inheritance diagram for Nektar::MMFDiffusion:
[legend]

Public Member Functions

 ~MMFDiffusion () override
 Desctructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::MMFSystem
SOLVER_UTILS_EXPORT MMFSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
SOLVER_UTILS_EXPORT ~MMFSystem () override
 
SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Virtual function for generating summary information. More...
 
SOLVER_UTILS_EXPORT void MMFInitObject (const Array< OneD, const Array< OneD, NekDouble > > &Anisotropy, const int TangentXelem=-1)
 
SOLVER_UTILS_EXPORT void CopyBoundaryTrace (const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, const BoundaryCopyType BDCopyType, const int var=0, const std::string btype="NoUserDefined")
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT ~UnsteadySystem () override
 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 NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 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 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 NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. 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 GetTime ()
 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, const Array< OneD, const NekDouble > &input)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
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 int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::ALEHelper
virtual ~ALEHelper ()=default
 
virtual SOLVER_UTILS_EXPORT void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
SOLVER_UTILS_EXPORT void InitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity (const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_ALEPreMultiplyMass (Array< OneD, Array< OneD, NekDouble > > &fields)
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMass (Array< OneD, Array< OneD, NekDouble > > &traceNormals, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble time)
 Update m_fields with u^n by multiplying by inverse mass matrix. That's then used in e.g. checkpoint output and L^2 error calculation. More...
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void MoveMesh (const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
 
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace ()
 
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Static Public Member Functions

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

Public Attributes

TestType m_TestType
 
- Public Attributes inherited from Nektar::SolverUtils::MMFSystem
NekDouble m_pi
 
int m_shapedim
 
SurfaceType m_surfaceType
 
UpwindType m_upwindType
 
TestMaxwellType m_TestMaxwellType
 
PolType m_PolType
 
IncType m_IncType
 
Array< OneD, NekDoublem_MMFfactors
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

 MMFDiffusion (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
 Solve for the diffusion term. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Computes the reaction terms \(f(u,v)\) and \(g(u,v)\). More...
 
void TestPlaneProblem (const NekDouble time, Array< OneD, NekDouble > &outfield)
 
void TestCubeProblem (const NekDouble time, Array< OneD, NekDouble > &outfield)
 
void Morphogenesis (const NekDouble time, unsigned int field, Array< OneD, NekDouble > &outfield)
 
Array< OneD, NekDoublePlanePhiWave ()
 
void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
 Sets a custom initial condition. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Prints a summary of the model parameters. More...
 
void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::MMFSystem
void SetUpMovingFrames (const Array< OneD, const Array< OneD, NekDouble > > &Anisotropy, const int TangentXelem)
 
void CheckMovingFrames (const Array< OneD, const Array< OneD, NekDouble > > &movingframes)
 
SOLVER_UTILS_EXPORT void ComputencdotMF ()
 
SOLVER_UTILS_EXPORT void ComputeDivCurlMF ()
 
SOLVER_UTILS_EXPORT void ComputeMFtrace ()
 
SOLVER_UTILS_EXPORT void VectorDotProd (const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, NekDouble > &v3)
 
SOLVER_UTILS_EXPORT void VectorCrossProd (const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &v3)
 
SOLVER_UTILS_EXPORT void VectorCrossProd (const Array< OneD, NekDouble > &v1, const Array< OneD, NekDouble > &v2, Array< OneD, NekDouble > &v3)
 
SOLVER_UTILS_EXPORT void ComputeCurl (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleCartesianToMovingframes (const Array< OneD, const Array< OneD, NekDouble > > &uvec, unsigned int field)
 
SOLVER_UTILS_EXPORT void DeriveCrossProductMF (Array< OneD, Array< OneD, NekDouble > > &CrossProductMF)
 
SOLVER_UTILS_EXPORT void ComputeNtimesMF ()
 
SOLVER_UTILS_EXPORT void ComputeNtimesFz (const int dir, const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &imFwd, const Array< OneD, const NekDouble > &imBwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void ComputeNtimesF12 (const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &im1Fwd, const Array< OneD, const NekDouble > &im1Bwd, const Array< OneD, const NekDouble > &im2Fwd, const Array< OneD, const NekDouble > &im2Bwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void ComputeNtimestimesdFz (const int dir, const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &imFwd, const Array< OneD, const NekDouble > &imBwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void ComputeNtimestimesdF12 (const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &im1Fwd, const Array< OneD, const NekDouble > &im1Bwd, const Array< OneD, const NekDouble > &im2Fwd, const Array< OneD, const NekDouble > &im2Bwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void CartesianToSpherical (const NekDouble x0j, const NekDouble x1j, const NekDouble x2j, NekDouble &sin_varphi, NekDouble &cos_varphi, NekDouble &sin_theta, NekDouble &cos_theta)
 
SOLVER_UTILS_EXPORT void ComputeZimYim (Array< OneD, Array< OneD, NekDouble > > &epsvec, Array< OneD, Array< OneD, NekDouble > > &muvec)
 
SOLVER_UTILS_EXPORT void AdddedtMaxwell (const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void GetMaxwellFluxVector (const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetMaxwellFlux1D (const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetMaxwellFlux2D (const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void LaxFriedrichMaxwellFlux1D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
 
SOLVER_UTILS_EXPORT void UpwindMaxwellFlux1D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
 
SOLVER_UTILS_EXPORT void AverageMaxwellFlux1D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
 
SOLVER_UTILS_EXPORT void NumericalMaxwellFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time=0.0)
 
SOLVER_UTILS_EXPORT void NumericalMaxwellFluxTM (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time)
 
SOLVER_UTILS_EXPORT void NumericalMaxwellFluxTE (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetIncidentField (const int var, const NekDouble time)
 
SOLVER_UTILS_EXPORT void Computedemdxicdote ()
 
SOLVER_UTILS_EXPORT NekDouble AvgInt (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble AvgAbsInt (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble AbsIntegral (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble RootMeanSquare (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble VectorAvgMagnitude (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 
SOLVER_UTILS_EXPORT void GramSchumitz (const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &outarray, bool KeepTheMagnitude=true)
 
SOLVER_UTILS_EXPORT void BubbleSort (Array< OneD, NekDouble > &refarray, Array< OneD, NekDouble > &sortarray)
 
- 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 void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. 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)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
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...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- 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 void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. 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_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. 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 bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

InitWaveType m_InitWaveType
 
NekDouble m_InitPtx
 
NekDouble m_InitPty
 
NekDouble m_InitPtz
 
- Protected Attributes inherited from Nektar::SolverUtils::MMFSystem
NekDouble m_alpha
 
NekDouble m_Incfreq
 
int m_SmoothFactor
 
NekDouble m_SFinit
 
Array< OneD, Array< OneD, NekDouble > > m_movingframes
 
Array< OneD, Array< OneD, NekDouble > > m_surfaceNormal
 
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFFwd
 
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFBwd
 
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFFwd
 
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFBwd
 
Array< OneD, Array< OneD, NekDouble > > m_DivMF
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_CurlMF
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceFwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceBwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimesMFFwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimesMFBwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimes_ntimesMFFwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimes_ntimesMFBwd
 
Array< OneD, Array< OneD, NekDouble > > m_ZimFwd
 
Array< OneD, Array< OneD, NekDouble > > m_ZimBwd
 
Array< OneD, Array< OneD, NekDouble > > m_YimFwd
 
Array< OneD, Array< OneD, NekDouble > > m_YimBwd
 
Array< OneD, Array< OneD, NekDouble > > m_epsvec
 
Array< OneD, Array< OneD, NekDouble > > m_muvec
 
Array< OneD, Array< OneD, NekDouble > > m_negepsvecminus1
 
Array< OneD, Array< OneD, NekDouble > > m_negmuvecminus1
 
Array< OneD, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > > m_dedxi_cdot_e
 
SpatialDomains::GeomMMF m_MMFdir
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
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...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
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_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_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. 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...
 
Array< OneD, NekDoublem_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z. More...
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData 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...
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_ImplicitALESolver = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 

Private Attributes

StdRegions::VarCoeffMap m_varcoeff
 Variable diffusivity. More...
 
Array< OneD, NekDoublem_epsilon
 
Array< OneD, NekDoublem_epsu
 

Friends

class MemoryManager< MMFDiffusion >
 

Additional Inherited Members

- 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 []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

A model for cardiac conduction.

Definition at line 79 of file MMFDiffusion.h.

Constructor & Destructor Documentation

◆ ~MMFDiffusion()

Nektar::MMFDiffusion::~MMFDiffusion ( )
override

Desctructor.

Definition at line 181 of file MMFDiffusion.cpp.

182{
183}

◆ MMFDiffusion()

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

Constructor.

Definition at line 57 of file MMFDiffusion.cpp.

59 : UnsteadySystem(pSession, pGraph), MMFSystem(pSession, pGraph)
60{
61}
SOLVER_UTILS_EXPORT MMFSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: MMFSystem.cpp:41
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 85 of file MMFDiffusion.h.

88 {
91 p->InitObject();
92 return p;
93 }
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::MMFDiffusion::DoImplicitSolve ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
NekDouble  time,
NekDouble  lambda 
)
protected

Solve for the diffusion term.

OdeRhs

Parameters
inarrayInput array.
outarrayOutput array.
timeCurrent simulation time.
lambdaTimestep.

Definition at line 191 of file MMFDiffusion.cpp.

195{
196 int nvariables = inarray.size();
197 int nq = m_fields[0]->GetNpoints();
198
201
203 factors[StdRegions::eFactorLambda] = 1.0 / lambda;
204 F[0] = Array<OneD, NekDouble>(nq * nvariables);
205
206 for (int n = 1; n < nvariables; ++n)
207 {
208 F[n] = F[n - 1] + nq;
209 }
210
211 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
212 // inarray = input: \hat{rhs} -> output: \hat{Y}
213 // outarray = output: nabla^2 \hat{Y}
214 // where \hat = modal coeffs
216
217 for (int i = 0; i < nvariables; ++i)
218 {
219 factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsu[i];
220
221 // Multiply 1.0/timestep
222 Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
223 F[i], 1);
224
225 // Solve a system of equations with Helmholtz solver and transform
226 // back into physical space.
227 m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(), factors,
228 m_varcoeff);
229
230 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
231 }
232}
StdRegions::VarCoeffMap m_varcoeff
Variable diffusivity.
Definition: MMFDiffusion.h:150
Array< OneD, NekDouble > m_epsu
Definition: MMFDiffusion.h:153
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
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.hpp:100

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorTau, Nektar::VarcoeffHashingTest::factors, m_epsu, Nektar::SolverUtils::EquationSystem::m_fields, m_varcoeff, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Smul().

Referenced by v_InitObject().

◆ DoOdeRhs()

void Nektar::MMFDiffusion::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Computes the reaction terms \(f(u,v)\) and \(g(u,v)\).

Definition at line 237 of file MMFDiffusion.cpp.

240{
241 int nq = GetTotPoints();
242
243 switch (m_TestType)
244 {
245 case eTestPlane:
246 {
247
251
252 m_fields[0]->GetCoords(x, y, z);
253
254 for (int k = 0; k < nq; k++)
255 {
256 outarray[0][k] = (m_epsilon[0] + m_epsilon[1] - 1.0) * m_pi *
257 m_pi * exp(-1.0 * m_pi * m_pi * time) *
258 sin(m_pi * x[k]) * cos(m_pi * y[k]);
259 }
260 }
261 break;
262
263 case eTestCube:
264 {
265
269
270 m_fields[0]->GetCoords(x, y, z);
271
272 for (int k = 0; k < nq; k++)
273 {
274 outarray[0][k] =
275 (m_epsilon[0] + m_epsilon[1] + m_epsilon[2] - 1.0) * m_pi *
276 m_pi * exp(-1.0 * m_pi * m_pi * time) * sin(m_pi * x[k]) *
277 sin(m_pi * y[k]) * sin(m_pi * z[k]);
278 }
279 }
280 break;
281
283 {
284 Array<OneD, NekDouble> temp(nq);
285
286 NekDouble A = 2.0;
287 NekDouble B = 5.0;
288
289 NekDouble m_a, m_b, m_c, m_d;
290 m_a = B - 1.0;
291 m_b = A * A;
292 m_c = -1.0 * B;
293 m_d = -1.0 * A * A;
294
295 temp = Array<OneD, NekDouble>(nq, 0.0);
296 Vmath::Svtvp(nq, m_a, &inarray[0][0], 1, &temp[0], 1, &temp[0], 1);
297 Vmath::Svtvp(nq, m_b, &inarray[1][0], 1, &temp[0], 1,
298 &outarray[0][0], 1);
299
300 temp = Array<OneD, NekDouble>(nq, 0.0);
301 Vmath::Svtvp(nq, m_c, &inarray[0][0], 1, &temp[0], 1, &temp[0], 1);
302 Vmath::Svtvp(nq, m_d, &inarray[1][0], 1, &temp[0], 1,
303 &outarray[1][0], 1);
304 }
305 break;
306
308 {
309 NekDouble A = 2.0;
310 NekDouble B = 5.0;
311
312 Array<OneD, NekDouble> Aonevec(nq, A);
313
314 // cube = phys0*phys0*phy1
315 Array<OneD, NekDouble> cube(nq);
316 Vmath::Vmul(nq, &inarray[0][0], 1, &inarray[0][0], 1, &cube[0], 1);
317 Vmath::Vmul(nq, &inarray[1][0], 1, &cube[0], 1, &cube[0], 1);
318
319 // outarray[0] = A - B*phy0 + phy0*phy0*phy1 - phy0
320 NekDouble coeff = -1.0 * B - 1.0;
322 Vmath::Svtvp(nq, coeff, &inarray[0][0], 1, &cube[0], 1, &tmp[0], 1);
323 Vmath::Vadd(nq, &Aonevec[0], 1, &tmp[0], 1, &outarray[0][0], 1);
324
325 // outarray[1] = B*phys0 - phy0*phy0*phy1
326 Vmath::Svtvm(nq, B, &inarray[0][0], 1, &cube[0], 1, &outarray[1][0],
327 1);
328 }
329 break;
330
331 case eFHNStandard:
332 {
333 // \phi - \phi^3/3 - \psi
334 NekDouble a = 0.12;
335 NekDouble b = 0.011;
336 NekDouble c1 = 0.175;
337 NekDouble c2 = 0.03;
338 NekDouble d = 0.55;
339
341
342 // Reaction for \phi = c1 \phi ( \phi - a)*(1 - \phi) - c2 v
343 Vmath::Smul(nq, -1.0 * c1, inarray[0], 1, outarray[0], 1);
344 Vmath::Sadd(nq, -1.0 * a, inarray[0], 1, tmp, 1);
345 Vmath::Vmul(nq, tmp, 1, inarray[0], 1, outarray[0], 1);
346 Vmath::Sadd(nq, -1.0, inarray[0], 1, tmp, 1);
347 Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
348
349 Vmath::Smul(nq, -1.0 * c2, inarray[1], 1, tmp, 1);
350 Vmath::Vadd(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
351
352 // Reaction for \psi = b (\phi - d \psi )
353 Vmath::Svtvp(nq, -1.0 * d, inarray[1], 1, inarray[0], 1,
354 outarray[1], 1);
355 Vmath::Smul(nq, b, outarray[1], 1, outarray[1], 1);
356 }
357 break;
358
359 case eFHNRogers:
360 {
361 NekDouble a = 0.13;
362 NekDouble b = 0.013;
363 NekDouble c1 = 0.26;
364 NekDouble c2 = 0.1;
365 NekDouble d = 1.0;
366
368
369 // Reaction for \phi = c1 \phi ( \phi - a)*(1 - \phi) - c2 u v
370 Vmath::Smul(nq, -1.0 * c1, inarray[0], 1, outarray[0], 1);
371 Vmath::Sadd(nq, -1.0 * a, inarray[0], 1, tmp, 1);
372 Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
373 Vmath::Sadd(nq, -1.0, inarray[0], 1, tmp, 1);
374 Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
375
376 Vmath::Vmul(nq, inarray[0], 1, inarray[1], 1, tmp, 1);
377 Vmath::Smul(nq, -1.0 * c2, tmp, 1, tmp, 1);
378 Vmath::Vadd(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
379
380 // Reaction for \psi = b (\phi - d \psi )
381 Vmath::Svtvp(nq, -1.0 * d, inarray[1], 1, inarray[0], 1,
382 outarray[1], 1);
383 Vmath::Smul(nq, b, outarray[1], 1, outarray[1], 1);
384 }
385 break;
386
387 case eFHNAlievPanf:
388 {
389
390 NekDouble a = 0.15;
391 NekDouble c1 = 8.0;
392 NekDouble c2 = 1.0;
393 NekDouble c0 = 0.002;
394 NekDouble mu1 = 0.2;
395 NekDouble mu2 = 0.3;
396
398
399 // Reaction for \phi = c1 \phi ( \phi - a)*(1 - \phi) - c2 u v
400 Vmath::Smul(nq, -1.0 * c1, inarray[0], 1, outarray[0], 1);
401 Vmath::Sadd(nq, -1.0 * a, inarray[0], 1, tmp, 1);
402 Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
403 Vmath::Sadd(nq, -1.0, inarray[0], 1, tmp, 1);
404 Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
405
406 Vmath::Vmul(nq, inarray[0], 1, inarray[1], 1, tmp, 1);
407 Vmath::Smul(nq, -1.0 * c2, tmp, 1, tmp, 1);
408 Vmath::Vadd(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
409
410 // Reaction for \psi = (c0 + (\mu1 \psi/(\mu2+\phi) ) )*(-\psi - c1
411 // * \phi*(\phi - a - 1) )
412
413 Vmath::Smul(nq, mu1, inarray[1], 1, outarray[1], 1);
414 Vmath::Sadd(nq, mu2, inarray[0], 1, tmp, 1);
415 Vmath::Vdiv(nq, outarray[1], 1, tmp, 1, outarray[1], 1);
416 Vmath::Sadd(nq, c0, outarray[1], 1, outarray[1], 1);
417
418 Vmath::Sadd(nq, (-a - 1.0), inarray[0], 1, tmp, 1);
419 Vmath::Vmul(nq, inarray[0], 1, tmp, 1, tmp, 1);
420 Vmath::Smul(nq, c1, tmp, 1, tmp, 1);
421 Vmath::Vadd(nq, inarray[1], 1, tmp, 1, tmp, 1);
422 Vmath::Neg(nq, tmp, 1);
423
424 Vmath::Vmul(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
425 }
426 break;
427
428 default:
429 break;
430 }
431}
Array< OneD, NekDouble > m_epsilon
Definition: MMFDiffusion.h:152
SOLVER_UTILS_EXPORT int GetTotPoints()
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
@ eFHNStandard
Definition: MMFDiffusion.h:52
@ eTestLinearSphere
Definition: MMFDiffusion.h:50
@ eTestPlane
Definition: MMFDiffusion.h:48
@ eTestCube
Definition: MMFDiffusion.h:49
@ eTestNonlinearSphere
Definition: MMFDiffusion.h:51
@ eFHNRogers
Definition: MMFDiffusion.h:53
@ eFHNAlievPanf
Definition: MMFDiffusion.h:54
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Svtvm(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvm (scalar times vector minus vector): z = alpha*x - y.
Definition: Vmath.hpp:424
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.hpp:126
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194

References Nektar::UnitTests::d(), Nektar::eFHNAlievPanf, Nektar::eFHNRogers, Nektar::eFHNStandard, Nektar::eTestCube, Nektar::eTestLinearSphere, Nektar::eTestNonlinearSphere, Nektar::eTestPlane, Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_epsilon, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_pi, m_TestType, Vmath::Neg(), Vmath::Sadd(), Vmath::Smul(), Vmath::Svtvm(), Vmath::Svtvp(), Vmath::Vadd(), Vmath::Vdiv(), Vmath::Vmul(), and Nektar::UnitTests::z().

Referenced by v_InitObject().

◆ Morphogenesis()

void Nektar::MMFDiffusion::Morphogenesis ( const NekDouble  time,
unsigned int  field,
Array< OneD, NekDouble > &  outfield 
)
protected

Definition at line 548 of file MMFDiffusion.cpp.

550{
551 int nq = GetTotPoints();
552
553 int i, m, n, ind;
554 NekDouble a_n, d_n, gamma_n;
555 NekDouble A_mn, C_mn, theta, phi, radius;
556
557 std::complex<double> Spericharmonic, delta_n, temp;
558 std::complex<double> varphi0, varphi1;
559 std::complex<double> B_mn, D_mn;
560
561 // Set some parameter values
562 int Maxn = 6;
563 int Maxm = 2 * Maxn - 1;
564
565 NekDouble A = 2.0;
566 NekDouble B = 5.0;
567
568 NekDouble m_mu = 0.001;
569 NekDouble m_nu = 0.002;
570
571 NekDouble m_a, m_b, m_c, m_d;
572
573 m_a = B - 1.0;
574 m_b = A * A;
575 m_c = -1.0 * B;
576 m_d = -1.0 * A * A;
577
580
581 for (i = 0; i < Maxn; ++i)
582 {
583 Ainit[i] = Array<OneD, NekDouble>(Maxm, 0.0);
584 Binit[i] = Array<OneD, NekDouble>(Maxm, 0.0);
585 }
586
587 Ainit[5][0] = -0.5839;
588 Ainit[5][1] = -0.8436;
589 Ainit[5][2] = -0.4764;
590 Ainit[5][3] = 0.6475;
591 Ainit[5][4] = 0.1886;
592 Ainit[5][5] = 0.8709;
593 Ainit[5][6] = -0.8338;
594 Ainit[5][7] = 0.1795;
595 Ainit[5][8] = -0.7873;
596 Ainit[5][9] = 0.8842;
597 Ainit[5][10] = 0.2943;
598
599 Binit[5][0] = -0.6263;
600 Binit[5][1] = 0.9803;
601 Binit[5][2] = 0.7222;
602 Binit[5][3] = 0.5945;
603 Binit[5][4] = 0.6026;
604 Binit[5][5] = -0.2076;
605 Binit[5][6] = 0.4556;
606 Binit[5][7] = 0.6024;
607 Binit[5][8] = 0.9695;
608 Binit[5][9] = -0.4936;
609 Binit[5][10] = 0.1098;
610
616
617 m_fields[0]->GetCoords(x, y, z);
618 for (int i = 0; i < nq; ++i)
619 {
620 radius = sqrt(x[i] * x[i] + y[i] * y[i] + z[i] * z[i]);
621
622 // theta is in [0, pi]
623 theta = asin(z[i] / radius) + 0.5 * m_pi;
624
625 // phi is in [0, 2*pi]
626 phi = atan2(y[i], x[i]) + m_pi;
627
628 varphi0 = 0.0 * varphi0;
629 varphi1 = 0.0 * varphi1;
630 for (n = 0; n < Maxn; ++n)
631 {
632 // Set up parameters
633 a_n = m_a - m_mu * (n * (n + 1) / radius / radius);
634 d_n = m_d - m_nu * (n * (n + 1) / radius / radius);
635
636 gamma_n = 0.5 * (a_n + d_n);
637
638 temp = (a_n + d_n) * (a_n + d_n) - 4.0 * (a_n * d_n - m_b * m_c);
639 delta_n = 0.5 * sqrt(temp);
640
641 for (m = -n; m <= n; ++m)
642 {
643 ind = m + n;
644 A_mn = Ainit[n][ind];
645 C_mn = Binit[n][ind];
646
647 B_mn = ((a_n - gamma_n) * Ainit[n][ind] + m_b * Binit[n][ind]) /
648 delta_n;
649 D_mn = (m_c * Ainit[n][ind] + (d_n - gamma_n) * Binit[n][ind]) /
650 delta_n;
651
652 Spericharmonic =
653 boost::math::spherical_harmonic(n, m, theta, phi);
654 varphi0 += exp(gamma_n * time) *
655 (A_mn * cosh(delta_n * time) +
656 B_mn * sinh(delta_n * time)) *
657 Spericharmonic;
658 varphi1 += exp(gamma_n * time) *
659 (C_mn * cosh(delta_n * time) +
660 D_mn * sinh(delta_n * time)) *
661 Spericharmonic;
662 }
663 }
664
665 u[i] = varphi0.real();
666 v[i] = varphi1.real();
667 }
668
669 switch (field)
670 {
671 case 0:
672 {
673 outfield = u;
674 }
675 break;
676
677 case 1:
678 {
679 outfield = v;
680 }
681 break;
682 }
683}
NekDouble m_mu
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References FilterPython_Function::field, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, m_mu, Nektar::SolverUtils::MMFSystem::m_pi, tinysimd::sqrt(), and Nektar::UnitTests::z().

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

◆ PlanePhiWave()

Array< OneD, NekDouble > Nektar::MMFDiffusion::PlanePhiWave ( )
protected

Definition at line 685 of file MMFDiffusion.cpp.

686{
687 int nq = GetTotPoints();
688 Array<OneD, NekDouble> outarray(nq, 0.0);
689
693
694 m_fields[0]->GetCoords(x, y, z);
695
696 NekDouble xmin, ymin, xmax;
697
698 xmin = Vmath::Vmin(nq, x, 1);
699 xmax = Vmath::Vmax(nq, x, 1);
700 ymin = Vmath::Vmin(nq, y, 1);
701
702 NekDouble xp, yp, xp2;
703 for (int i = 0; i < nq; i++)
704 {
705 switch (m_InitWaveType)
706 {
707 case eLeft:
708 {
709 NekDouble radiusofinit = 4.0;
710 NekDouble frontstiff = 0.1;
711
712 xp = x[i] - xmin;
713 outarray[i] =
714 1.0 / (1.0 + exp((xp - radiusofinit) / frontstiff));
715 }
716 break;
717
718 case eBothEnds:
719 {
720 NekDouble radiusofinit = 3.0;
721 NekDouble frontstiff = 0.1;
722
723 xp = x[i] - xmin;
724 xp2 = x[i] - xmax;
725
726 outarray[i] =
727 1.0 / (1.0 +
728 exp((sqrt(xp * xp) - radiusofinit) / frontstiff)) +
729 1.0 / (1.0 +
730 exp((sqrt(xp2 * xp2) - radiusofinit) / frontstiff));
731 }
732 break;
733
734 case eCenter:
735 {
736 NekDouble radiusofinit = 6.0;
737 NekDouble frontstiff = 0.1;
738
739 xp = x[i] - xmin;
740 outarray[i] =
741 1.0 / (1.0 + exp((xp - radiusofinit) / frontstiff));
742 }
743 break;
744
746 {
747 NekDouble radiusofinit = 6.0;
748 NekDouble frontstiff = 0.1;
749 NekDouble bs = 2.0;
750
751 xp = x[i] - xmin;
752 yp = y[i] - ymin;
753 outarray[i] =
754 1.0 /
755 (1.0 + exp((sqrt(xp * xp + yp * yp) / bs - radiusofinit) /
756 frontstiff));
757 }
758 break;
759
760 case ePoint:
761 {
762 NekDouble xloc, yloc, zloc, rad;
763 NekDouble radiusofinit = 10.0;
764
765 xloc = x[i] - m_InitPtx;
766 yloc = y[i] - m_InitPty;
767 zloc = z[i] - m_InitPtz;
768
769 rad = sqrt(xloc * xloc + yloc * yloc + zloc * zloc);
770
771 xloc = xloc / radiusofinit;
772 yloc = yloc / radiusofinit;
773 zloc = zloc / radiusofinit;
774
775 if (rad < radiusofinit)
776 {
777 outarray[i] =
778 exp(-(1.0 / 2.0) *
779 (xloc * xloc + yloc * yloc + zloc * zloc));
780 }
781
782 else
783 {
784 outarray[i] = 0.0;
785 }
786 }
787 break;
788
789 case eSpiralDock:
790 {
791 NekDouble radiusofinit = 3.0;
792 NekDouble frontstiff = 0.1;
793 xp = x[i] - 4.0;
794 yp = y[i];
795 outarray[i] =
796 (1.0 / (1.0 + exp(2.0 * yp))) *
797 (1.0 / (1.0 + exp(-2.0 * xp))) *
798 (1.0 / (1.0 + exp((xp - radiusofinit) / frontstiff)));
799 }
800 break;
801
802 default:
803 break;
804 }
805 }
806
807 return outarray;
808}
InitWaveType m_InitWaveType
Definition: MMFDiffusion.h:108
static NekDouble rad(NekDouble x, NekDouble y)
@ eLeftBottomCorner
Definition: MMFDiffusion.h:68
@ eBothEnds
Definition: MMFDiffusion.h:66
@ eSpiralDock
Definition: MMFDiffusion.h:70
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.hpp:725
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.hpp:644

References Nektar::eBothEnds, Nektar::eCenter, Nektar::eLeft, Nektar::eLeftBottomCorner, Nektar::ePoint, Nektar::eSpiralDock, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, m_InitPtx, m_InitPty, m_InitPtz, m_InitWaveType, Nektar::LibUtilities::rad(), tinysimd::sqrt(), Vmath::Vmax(), Vmath::Vmin(), and Nektar::UnitTests::z().

Referenced by v_SetInitialConditions().

◆ TestCubeProblem()

void Nektar::MMFDiffusion::TestCubeProblem ( const NekDouble  time,
Array< OneD, NekDouble > &  outfield 
)
protected

Definition at line 528 of file MMFDiffusion.cpp.

531{
532 int nq = GetTotPoints();
533
537
538 m_fields[0]->GetCoords(x, y, z);
539
540 outfield = Array<OneD, NekDouble>(nq);
541 for (int k = 0; k < nq; k++)
542 {
543 outfield[k] = exp(-1.0 * m_pi * m_pi * time) * sin(m_pi * x[k]) *
544 sin(m_pi * y[k]) * sin(m_pi * z[k]);
545 }
546}

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_pi, and Nektar::UnitTests::z().

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

◆ TestPlaneProblem()

void Nektar::MMFDiffusion::TestPlaneProblem ( const NekDouble  time,
Array< OneD, NekDouble > &  outfield 
)
protected

Definition at line 508 of file MMFDiffusion.cpp.

511{
512 int nq = GetTotPoints();
513
517
518 m_fields[0]->GetCoords(x, y, z);
519
520 outfield = Array<OneD, NekDouble>(nq);
521 for (int k = 0; k < nq; k++)
522 {
523 outfield[k] = exp(-1.0 * m_pi * m_pi * time) * sin(m_pi * x[k]) *
524 cos(m_pi * y[k]);
525 }
526}

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_pi, and Nektar::UnitTests::z().

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

◆ v_EvaluateExactSolution()

void Nektar::MMFDiffusion::v_EvaluateExactSolution ( unsigned int  field,
Array< OneD, NekDouble > &  outfield,
const NekDouble  time 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 810 of file MMFDiffusion.cpp.

813{
814 switch (m_TestType)
815 {
816 case eTestPlane:
817 {
818 TestPlaneProblem(time, outfield);
819 }
820 break;
821
822 case eTestCube:
823 {
824 TestCubeProblem(time, outfield);
825 }
826 break;
827
830 {
831 Morphogenesis(time, field, outfield);
832 }
833 break;
834
835 case eFHNStandard:
836 case eFHNRogers:
837 case eFHNAlievPanf:
838 {
839 int nq = GetTotPoints();
840 outfield = Array<OneD, NekDouble>(nq, 0.0);
841 }
842 /* Falls through. */
843 default:
844 {
846 }
847 break;
848 }
849}
void Morphogenesis(const NekDouble time, unsigned int field, Array< OneD, NekDouble > &outfield)
void TestCubeProblem(const NekDouble time, Array< OneD, NekDouble > &outfield)
void TestPlaneProblem(const NekDouble time, Array< OneD, NekDouble > &outfield)
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)

References Nektar::eFHNAlievPanf, Nektar::eFHNRogers, Nektar::eFHNStandard, Nektar::eTestCube, Nektar::eTestLinearSphere, Nektar::eTestNonlinearSphere, Nektar::eTestPlane, FilterPython_Function::field, Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_TestType, Morphogenesis(), TestCubeProblem(), TestPlaneProblem(), and Nektar::SolverUtils::EquationSystem::v_EvaluateExactSolution().

◆ v_GenerateSummary()

void Nektar::MMFDiffusion::v_GenerateSummary ( SolverUtils::SummaryList s)
overrideprotectedvirtual

Prints a summary of the model parameters.

Reimplemented from Nektar::SolverUtils::MMFSystem.

Definition at line 851 of file MMFDiffusion.cpp.

852{
855 SolverUtils::AddSummaryItem(s, "epsilon0", m_epsilon[0]);
856 SolverUtils::AddSummaryItem(s, "epsilon1", m_epsilon[1]);
857 SolverUtils::AddSummaryItem(s, "epsilon2", m_epsilon[2]);
859 {
860 SolverUtils::AddSummaryItem(s, "epsilon for u", m_epsu[0]);
861 SolverUtils::AddSummaryItem(s, "epsilon for v", m_epsu[1]);
862 }
863}
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Virtual function for generating summary information.
Definition: MMFSystem.cpp:2463
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
const char *const TestTypeMap[]
Definition: MMFDiffusion.h:58

References Nektar::SolverUtils::AddSummaryItem(), Nektar::eTestLinearSphere, m_epsilon, m_epsu, m_TestType, Nektar::TestTypeMap, and Nektar::SolverUtils::MMFSystem::v_GenerateSummary().

◆ v_InitObject()

void Nektar::MMFDiffusion::v_InitObject ( bool  DeclareField = true)
overrideprotectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 63 of file MMFDiffusion.cpp.

64{
65 UnsteadySystem::v_InitObject(DeclareFields);
66
67 int nq = m_fields[0]->GetNpoints();
68 int nvar = m_fields.size();
69 int MFdim = 3;
70
71 // Diffusivity coefficient for e^j
73 m_session->LoadParameter("epsilon0", m_epsilon[0], 1.0);
74 m_session->LoadParameter("epsilon1", m_epsilon[1], 1.0);
75 m_session->LoadParameter("epsilon2", m_epsilon[2], 1.0);
76
77 // Diffusivity coefficient for u^j
79 m_session->LoadParameter("epsu0", m_epsu[0], 1.0);
80 m_session->LoadParameter("epsu1", m_epsu[1], 1.0);
81
82 m_session->LoadParameter("InitPtx", m_InitPtx, 0.0);
83 m_session->LoadParameter("InitPty", m_InitPty, 0.0);
84 m_session->LoadParameter("InitPtz", m_InitPtz, 0.0);
85
86 int shapedim = m_fields[0]->GetShapeDimension();
87 Array<OneD, Array<OneD, NekDouble>> Anisotropy(shapedim);
88 for (int j = 0; j < shapedim; ++j)
89 {
90 Anisotropy[j] = Array<OneD, NekDouble>(nq, 1.0);
91 Vmath::Fill(nq, sqrt(m_epsilon[j]), &Anisotropy[j][0], 1);
92 }
93
94 MMFSystem::MMFInitObject(Anisotropy);
95
96 // Define ProblemType
97 if (m_session->DefinesSolverInfo("TESTTYPE"))
98 {
99 std::string TestTypeStr = m_session->GetSolverInfo("TESTTYPE");
100 int i;
101 for (i = 0; i < (int)SIZE_TestType; ++i)
102 {
103 if (boost::iequals(TestTypeMap[i], TestTypeStr))
104 {
105 m_TestType = (TestType)i;
106 break;
107 }
108 }
109 }
110 else
111 {
112 m_TestType = (TestType)0;
113 }
114
115 if (m_session->DefinesSolverInfo("INITWAVETYPE"))
116 {
117 std::string InitWaveTypeStr = m_session->GetSolverInfo("INITWAVETYPE");
118 for (int i = 0; i < (int)SIZE_TestType; ++i)
119 {
120 if (boost::iequals(InitWaveTypeMap[i], InitWaveTypeStr))
121 {
123 break;
124 }
125 }
126 }
127 else
128 {
130 }
131
132 StdRegions::VarCoeffType MMFCoeffs[15] = {
141
142 int indx;
144 for (int k = 0; k < MFdim; ++k)
145 {
146 // For Moving Frames
147 indx = 5 * k;
148
149 for (int j = 0; j < m_spacedim; ++j)
150 {
151 Vmath::Vcopy(nq, &m_movingframes[k][j * nq], 1, &tmp[0], 1);
152 m_varcoeff[MMFCoeffs[indx + j]] = tmp;
153 }
154
155 // m_DivMF
156 Vmath::Vcopy(nq, &m_DivMF[k][0], 1, &tmp[0], 1);
157 m_varcoeff[MMFCoeffs[indx + 3]] = tmp;
158
159 // \| e^k \|
160 tmp = Array<OneD, NekDouble>(nq, 0.0);
161 for (int i = 0; i < m_spacedim; ++i)
162 {
163 Vmath::Vvtvp(nq, &m_movingframes[k][i * nq], 1,
164 &m_movingframes[k][i * nq], 1, &tmp[0], 1, &tmp[0], 1);
165 }
166
167 m_varcoeff[MMFCoeffs[indx + 4]] = tmp;
168 }
169
171 {
173 }
174
176}
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Computes the reaction terms and .
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Solve for the diffusion term.
int m_spacedim
Spatial dimension (>= expansion dim).
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_DivMF
Definition: MMFSystem.h:192
Array< OneD, Array< OneD, NekDouble > > m_movingframes
Definition: MMFSystem.h:183
SOLVER_UTILS_EXPORT void MMFInitObject(const Array< OneD, const Array< OneD, NekDouble > > &Anisotropy, const int TangentXelem=-1)
Definition: MMFSystem.cpp:51
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
const char *const InitWaveTypeMap[]
Definition: MMFDiffusion.h:74
@ SIZE_TestType
Length of enum list.
Definition: MMFDiffusion.h:55
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), DoImplicitSolve(), DoOdeRhs(), Nektar::StdRegions::eVarCoeffMF1Div, Nektar::StdRegions::eVarCoeffMF1Mag, Nektar::StdRegions::eVarCoeffMF1x, Nektar::StdRegions::eVarCoeffMF1y, Nektar::StdRegions::eVarCoeffMF1z, Nektar::StdRegions::eVarCoeffMF2Div, Nektar::StdRegions::eVarCoeffMF2Mag, Nektar::StdRegions::eVarCoeffMF2x, Nektar::StdRegions::eVarCoeffMF2y, Nektar::StdRegions::eVarCoeffMF2z, Nektar::StdRegions::eVarCoeffMF3Div, Nektar::StdRegions::eVarCoeffMF3Mag, Nektar::StdRegions::eVarCoeffMF3x, Nektar::StdRegions::eVarCoeffMF3y, Nektar::StdRegions::eVarCoeffMF3z, Vmath::Fill(), Nektar::InitWaveTypeMap, Nektar::SolverUtils::MMFSystem::m_DivMF, m_epsilon, m_epsu, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, m_InitPtx, m_InitPty, m_InitPtz, m_InitWaveType, Nektar::SolverUtils::MMFSystem::m_movingframes, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_TestType, m_varcoeff, Nektar::SolverUtils::MMFSystem::MMFInitObject(), Nektar::SIZE_TestType, tinysimd::sqrt(), Nektar::TestTypeMap, Nektar::SolverUtils::UnsteadySystem::v_InitObject(), Vmath::Vcopy(), and Vmath::Vvtvp().

◆ v_SetInitialConditions()

void Nektar::MMFDiffusion::v_SetInitialConditions ( NekDouble  initialtime = 0.0,
bool  dumpInitialConditions = true,
const int  domain = 0 
)
overrideprotectedvirtual

Sets a custom initial condition.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 436 of file MMFDiffusion.cpp.

439{
440 int nq = GetTotPoints();
441
442 switch (m_TestType)
443 {
444 case eTestPlane:
445 {
447
448 TestPlaneProblem(initialtime, u);
449 m_fields[0]->SetPhys(u);
450 }
451 break;
452
453 case eTestCube:
454 {
456
457 TestCubeProblem(initialtime, u);
458 m_fields[0]->SetPhys(u);
459 }
460 break;
461
464 {
467
468 Morphogenesis(initialtime, 0, u);
469 Morphogenesis(initialtime, 1, v);
470
471 m_fields[0]->SetPhys(u);
472 m_fields[1]->SetPhys(v);
473 }
474 break;
475
476 case eFHNStandard:
477 case eFHNRogers:
478 case eFHNAlievPanf:
479 {
481 m_fields[0]->SetPhys(PlanePhiWave());
482 m_fields[1]->SetPhys(Zero);
483 }
484 break;
485
486 default:
487 {
488 EquationSystem::v_SetInitialConditions(initialtime, false);
489 }
490 break;
491 }
492
493 // forward transform to fill the modal coeffs
494 for (int i = 0; i < m_fields.size(); ++i)
495 {
496 m_fields[i]->SetPhysState(true);
497 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
498 m_fields[i]->UpdateCoeffs());
499 }
500
501 if (dumpInitialConditions)
502 {
503 std::string outname = m_sessionName + "_initial.chk";
504 WriteFld(outname);
505 }
506}
Array< OneD, NekDouble > PlanePhiWave()
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
std::string m_sessionName
Name of the session.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::eFHNAlievPanf, Nektar::eFHNRogers, Nektar::eFHNStandard, Nektar::eTestCube, Nektar::eTestLinearSphere, Nektar::eTestNonlinearSphere, Nektar::eTestPlane, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_sessionName, m_TestType, Morphogenesis(), PlanePhiWave(), TestCubeProblem(), TestPlaneProblem(), Nektar::SolverUtils::EquationSystem::v_SetInitialConditions(), Nektar::SolverUtils::EquationSystem::WriteFld(), and Vmath::Zero().

Friends And Related Function Documentation

◆ MemoryManager< MMFDiffusion >

friend class MemoryManager< MMFDiffusion >
friend

Definition at line 74 of file MMFDiffusion.h.

Member Data Documentation

◆ className

string Nektar::MMFDiffusion::className
static
Initial value:
=
"MMFDiffusion", MMFDiffusion::create, "MMFDiffusion equation.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: MMFDiffusion.h:85
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 96 of file MMFDiffusion.h.

◆ m_epsilon

Array<OneD, NekDouble> Nektar::MMFDiffusion::m_epsilon
private

Definition at line 152 of file MMFDiffusion.h.

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

◆ m_epsu

Array<OneD, NekDouble> Nektar::MMFDiffusion::m_epsu
private

Definition at line 153 of file MMFDiffusion.h.

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

◆ m_InitPtx

NekDouble Nektar::MMFDiffusion::m_InitPtx
protected

Definition at line 146 of file MMFDiffusion.h.

Referenced by PlanePhiWave(), and v_InitObject().

◆ m_InitPty

NekDouble Nektar::MMFDiffusion::m_InitPty
protected

Definition at line 146 of file MMFDiffusion.h.

Referenced by PlanePhiWave(), and v_InitObject().

◆ m_InitPtz

NekDouble Nektar::MMFDiffusion::m_InitPtz
protected

Definition at line 146 of file MMFDiffusion.h.

Referenced by PlanePhiWave(), and v_InitObject().

◆ m_InitWaveType

InitWaveType Nektar::MMFDiffusion::m_InitWaveType
protected

Definition at line 108 of file MMFDiffusion.h.

Referenced by PlanePhiWave(), and v_InitObject().

◆ m_TestType

TestType Nektar::MMFDiffusion::m_TestType

◆ m_varcoeff

StdRegions::VarCoeffMap Nektar::MMFDiffusion::m_varcoeff
private

Variable diffusivity.

Definition at line 150 of file MMFDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().