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

virtual ~MMFDiffusion ()
 Desctructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::MMFSystem
SOLVER_UTILS_EXPORT MMFSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual SOLVER_UTILS_EXPORT ~MMFSystem ()
 
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
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...
 
- 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 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
 
- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 

Protected Member Functions

 MMFDiffusion (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
virtual void v_InitObject ()
 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 ()
 
virtual void v_SetInitialConditions (NekDouble initialtime, bool dumpInitialConditions, const int domain)
 Sets a custom initial condition. More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Prints a summary of the model parameters. More...
 
virtual void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
- 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 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 ()
 
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...
 
- 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_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

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
 
Array< OneD, NekDoublem_MFlength
 
- 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::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
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...
 
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...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
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...
 
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

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 []
 

Detailed Description

A model for cardiac conduction.

Definition at line 91 of file MMFDiffusion.h.

Constructor & Destructor Documentation

◆ ~MMFDiffusion()

Nektar::MMFDiffusion::~MMFDiffusion ( )
virtual

Desctructor.

Definition at line 195 of file MMFDiffusion.cpp.

196  {
197  }

◆ MMFDiffusion()

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

Constructor.

Definition at line 58 of file MMFDiffusion.cpp.

61  : UnsteadySystem(pSession, pGraph),
62  MMFSystem(pSession, pGraph)
63  {
64  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT MMFSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: MMFSystem.cpp:43

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 97 of file MMFDiffusion.h.

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

100  {
102  p->InitObject();
103  return p;
104  }
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.

◆ 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 206 of file MMFDiffusion.cpp.

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

Referenced by v_InitObject().

211  {
212  int nvariables = inarray.num_elements();
213  int nq = m_fields[0]->GetNpoints();
214 
215 
217  factors[StdRegions::eFactorTau] = 1.0;
218 
219  Array<OneD, Array< OneD, NekDouble> > F(nvariables);
220  factors[StdRegions::eFactorLambda] = 1.0/lambda;
221  F[0] = Array<OneD, NekDouble> (nq*nvariables);
222 
223  for (int n = 1; n < nvariables; ++n)
224  {
225  F[n] = F[n-1] + nq;
226  //cout << "F["<< n<<"=" << F[n][1] <<endl;
227  }
228 
229  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
230  // inarray = input: \hat{rhs} -> output: \hat{Y}
231  // outarray = output: nabla^2 \hat{Y}
232  // where \hat = modal coeffs
233  SetBoundaryConditions(time);
234 
235  for (int i = 0; i < nvariables; ++i)
236  {
237  factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsu[i];
238 
239  // Multiply 1.0/timestep
240  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda],inarray[i], 1, F[i], 1);
241 
242  /* for (int k = 0; k < 15; ++k)
243  cout << "inarray["<<i << "]"<< k<<"=" << inarray[i][k]<<endl;*/
244  // Solve a system of equations with Helmholtz solver and transform
245  // back into physical space.
246  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),NullFlagList, factors, m_varcoeff);
247 
248  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
249  /* Array<OneD, NekDouble> coefarray = m_fields[i]->GetCoeffs();
250  for (int k = 0; k < 15; ++k)
251  cout << "inarray["<< k<<"=" << coefarray[k]<<endl;*/
252  }
253  /* for (int kk = 0; kk < 15; ++kk)
254  cout << "inarray["<< kk<<"=" << m_varcoeff[StdRegions::eVarCoeffMF3Mag][kk]<<endl;*/
255 
256 
257  }
Array< OneD, NekDouble > m_epsu
Definition: MMFDiffusion.h:169
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
StdRegions::VarCoeffMap m_varcoeff
Variable diffusivity.
Definition: MMFDiffusion.h:166
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
static FlagList NullFlagList
An empty flag list.

◆ 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 263 of file MMFDiffusion.cpp.

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

Referenced by v_InitObject().

267  {
268  int nq = GetTotPoints();
269 
270  switch(m_TestType)
271  {
272  case eTestPlane:
273  {
274 
278 
279  m_fields[0]->GetCoords(x,y,z);
280 
281  for(int k=0; k<nq; k++)
282  {
283  outarray[0][k] = (m_epsilon[0]+m_epsilon[1]-1.0)*m_pi*m_pi*exp(-1.0*m_pi*m_pi*time)*sin(m_pi*x[k])*cos(m_pi*y[k]);
284  }
285  }
286  break;
287 
288  case eTestCube:
289  {
290 
294 
295  m_fields[0]->GetCoords(x,y,z);
296 
297  for(int k=0; k<nq; k++)
298  {
299  outarray[0][k] = (m_epsilon[0]+m_epsilon[1]+m_epsilon[2]-1.0)*m_pi*m_pi*exp(-1.0*m_pi*m_pi*time)*sin(m_pi*x[k])*sin(m_pi*y[k])*sin(m_pi*z[k]);
300 
301  }
302 
303  }
304  break;
305 
306  case eTestLinearSphere:
307  {
308  Array<OneD, NekDouble> temp(nq);
309 
310  NekDouble A = 2.0;
311  NekDouble B = 5.0;
312 
313  NekDouble m_a, m_b, m_c, m_d;
314  m_a = B-1.0;
315  m_b = A*A;
316  m_c = -1.0*B;
317  m_d = -1.0*A*A;
318 
319  temp = Array<OneD, NekDouble>(nq,0.0);
320  Vmath::Svtvp(nq,m_a,&inarray[0][0],1,&temp[0],1,&temp[0],1);
321  Vmath::Svtvp(nq,m_b,&inarray[1][0],1,&temp[0],1,&outarray[0][0],1);
322 
323  temp = Array<OneD, NekDouble>(nq,0.0);
324  Vmath::Svtvp(nq,m_c,&inarray[0][0],1,&temp[0],1,&temp[0],1);
325  Vmath::Svtvp(nq,m_d,&inarray[1][0],1,&temp[0],1,&outarray[1][0],1);
326  }
327  break;
328 
330  {
331  NekDouble A = 2.0;
332  NekDouble B = 5.0;
333 
334  Array<OneD, NekDouble> Aonevec(nq,A);
335 
336  // cube = phys0*phys0*phy1
337  Array<OneD, NekDouble> cube(nq);
338  Vmath::Vmul(nq,&inarray[0][0],1,&inarray[0][0],1,&cube[0],1);
339  Vmath::Vmul(nq,&inarray[1][0],1,&cube[0],1,&cube[0],1);
340 
341  // outarray[0] = A - B*phy0 + phy0*phy0*phy1 - phy0
342  NekDouble coeff = -1.0*B - 1.0;
343  Array<OneD, NekDouble> tmp(nq);
344  Vmath::Svtvp(nq,coeff,&inarray[0][0],1,&cube[0],1,&tmp[0],1);
345  Vmath::Vadd(nq,&Aonevec[0],1,&tmp[0],1,&outarray[0][0],1);
346 
347  // outarray[1] = B*phys0 - phy0*phy0*phy1
348  Vmath::Svtvm(nq,B,&inarray[0][0],1,&cube[0],1,&outarray[1][0],1);
349 
350  }
351  break;
352 
353  case eFHNStandard:
354  {
355  // \phi - \phi^3/3 - \psi
356  NekDouble a = 0.12;
357  NekDouble b = 0.011;
358  NekDouble c1 = 0.175;
359  NekDouble c2 = 0.03;
360  NekDouble d = 0.55;
361 
362  Array<OneD, NekDouble> tmp(nq);
363 
364  // Reaction for \phi = c1 \phi ( \phi - a)*(1 - \phi) - c2 v
365  Vmath::Smul(nq, -1.0*c1, inarray[0], 1, outarray[0], 1);
366  Vmath::Sadd(nq, -1.0*a, inarray[0], 1, tmp, 1);
367  Vmath::Vmul(nq, tmp, 1, inarray[0], 1, outarray[0], 1);
368  Vmath::Sadd(nq, -1.0, inarray[0], 1, tmp, 1);
369  Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
370 
371  Vmath::Smul(nq, -1.0*c2, inarray[1], 1, tmp, 1);
372  Vmath::Vadd(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
373 
374 
375  // Reaction for \psi = b (\phi - d \psi )
376  Vmath::Svtvp(nq, -1.0*d, inarray[1], 1, inarray[0], 1, outarray[1], 1);
377  Vmath::Smul(nq, b, outarray[1], 1, outarray[1], 1);
378  }
379  break;
380 
381  case eFHNRogers:
382  {
383  NekDouble a = 0.13;
384  NekDouble b = 0.013;
385  NekDouble c1 = 0.26;
386  NekDouble c2 = 0.1;
387  NekDouble d = 1.0;
388 
389  Array<OneD, NekDouble> tmp(nq);
390 
391  // Reaction for \phi = c1 \phi ( \phi - a)*(1 - \phi) - c2 u v
392  Vmath::Smul(nq, -1.0*c1, inarray[0], 1, outarray[0], 1);
393  Vmath::Sadd(nq, -1.0*a, inarray[0], 1, tmp, 1);
394  Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
395  Vmath::Sadd(nq, -1.0, inarray[0], 1, tmp, 1);
396  Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
397 
398  Vmath::Vmul(nq, inarray[0], 1, inarray[1], 1, tmp, 1);
399  Vmath::Smul(nq, -1.0*c2, tmp, 1, tmp, 1);
400  Vmath::Vadd(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
401 
402  // Reaction for \psi = b (\phi - d \psi )
403  Vmath::Svtvp(nq, -1.0*d, inarray[1], 1, inarray[0], 1, outarray[1], 1);
404  Vmath::Smul(nq, b, outarray[1], 1, outarray[1], 1);
405  }
406  break;
407 
408  case eFHNAlievPanf:
409  {
410 
411  NekDouble a = 0.15;
412  NekDouble c1 = 8.0;
413  NekDouble c2 = 1.0;
414  NekDouble c0 = 0.002;
415  NekDouble mu1 = 0.2;
416  NekDouble mu2 = 0.3;
417 
418  Array<OneD, NekDouble> tmp(nq);
419 
420  // Reaction for \phi = c1 \phi ( \phi - a)*(1 - \phi) - c2 u v
421  Vmath::Smul(nq, -1.0*c1, inarray[0], 1, outarray[0], 1);
422  Vmath::Sadd(nq, -1.0*a, inarray[0], 1, tmp, 1);
423  Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
424  Vmath::Sadd(nq, -1.0, inarray[0], 1, tmp, 1);
425  Vmath::Vmul(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
426 
427  Vmath::Vmul(nq, inarray[0], 1, inarray[1], 1, tmp, 1);
428  Vmath::Smul(nq, -1.0*c2, tmp, 1, tmp, 1);
429  Vmath::Vadd(nq, tmp, 1, outarray[0], 1, outarray[0], 1);
430 
431  // Reaction for \psi = (c0 + (\mu1 \psi/(\mu2+\phi) ) )*(-\psi - c1 * \phi*(\phi - a - 1) )
432 
433  Vmath::Smul(nq, mu1, inarray[1], 1, outarray[1], 1);
434  Vmath::Sadd(nq, mu2, inarray[0], 1, tmp, 1);
435  Vmath::Vdiv(nq, outarray[1], 1, tmp, 1, outarray[1], 1);
436  Vmath::Sadd(nq, c0, outarray[1], 1, outarray[1], 1);
437 
438  Vmath::Sadd(nq, (-a-1.0), inarray[0], 1, tmp, 1);
439  Vmath::Vmul(nq, inarray[0], 1, tmp, 1, tmp, 1);
440  Vmath::Smul(nq, c1, tmp, 1, tmp, 1);
441  Vmath::Vadd(nq, inarray[1], 1, tmp, 1, tmp, 1);
442  Vmath::Neg(nq, tmp, 1);
443 
444  Vmath::Vmul(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
445  }
446  break;
447 
448  default:
449  break;
450  }
451  }
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.cpp:488
Array< OneD, NekDouble > m_epsilon
Definition: MMFDiffusion.h:168
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.cpp:244
SOLVER_UTILS_EXPORT int GetTotPoints()
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
void Svtvm(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.cpp:521
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
double NekDouble
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.cpp:318
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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.cpp:302
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.cpp:186

◆ Morphogenesis()

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

Definition at line 580 of file MMFDiffusion.cpp.

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

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

583  {
584  int nq = GetTotPoints();
585 
586  int i, m, n, ind;
587  NekDouble a_n, d_n, gamma_n;
588  NekDouble A_mn, C_mn, theta, phi,radius;
589 
590  std::complex<double> Spericharmonic, delta_n, temp;
591  std::complex<double> varphi0, varphi1;
592  std::complex<double> B_mn, D_mn;
593 
594  // Set some parameter values
595  int Maxn = 6;
596  int Maxm = 2*Maxn-1;
597 
598  NekDouble A = 2.0;
599  NekDouble B = 5.0;
600 
601  NekDouble m_mu = 0.001;
602  NekDouble m_nu = 0.002;
603 
604  NekDouble m_a, m_b, m_c, m_d;
605 
606  m_a = B-1.0;
607  m_b = A*A;
608  m_c = -1.0*B;
609  m_d = -1.0*A*A;
610 
613 
614  for (i = 0; i < Maxn; ++i)
615  {
616  Ainit[i] = Array<OneD, NekDouble>(Maxm, 0.0);
617  Binit[i] = Array<OneD, NekDouble>(Maxm, 0.0);
618  }
619 
620  Ainit[5][0] = -0.5839;
621  Ainit[5][1] = -0.8436;
622  Ainit[5][2] = -0.4764;
623  Ainit[5][3] = 0.6475;
624  Ainit[5][4] = 0.1886;
625  Ainit[5][5] = 0.8709;
626  Ainit[5][6] = -0.8338;
627  Ainit[5][7] = 0.1795;
628  Ainit[5][8] = -0.7873;
629  Ainit[5][9] = 0.8842;
630  Ainit[5][10] = 0.2943;
631 
632  Binit[5][0] = -0.6263;
633  Binit[5][1] = 0.9803;
634  Binit[5][2] = 0.7222;
635  Binit[5][3] = 0.5945;
636  Binit[5][4] = 0.6026;
637  Binit[5][5] = -0.2076;
638  Binit[5][6] = 0.4556;
639  Binit[5][7] = 0.6024;
640  Binit[5][8] = 0.9695;
641  Binit[5][9] = -0.4936;
642  Binit[5][10] = 0.1098;
643 
649 
650  m_fields[0]->GetCoords(x,y,z);
651  for (int i = 0; i < nq; ++i)
652  {
653  radius = sqrt(x[i]*x[i] + y[i]*y[i] + z[i]*z[i]) ;
654 
655  // theta is in [0, pi]
656  theta = asin( z[i]/radius ) + 0.5*m_pi;
657 
658  // phi is in [0, 2*pi]
659  phi = atan2( y[i], x[i] ) + m_pi;
660 
661  varphi0 = 0.0*varphi0;
662  varphi1 = 0.0*varphi1;
663  for (n = 0; n < Maxn; ++n)
664  {
665  // Set up parameters
666  a_n = m_a - m_mu*( n*(n+1)/radius/radius );
667  d_n = m_d - m_nu*( n*(n+1)/radius/radius );
668 
669  gamma_n = 0.5*( a_n + d_n );
670 
671  temp = ( a_n + d_n )*( a_n + d_n ) - 4.0*( a_n*d_n - m_b*m_c );
672  delta_n = 0.5*sqrt( temp );
673 
674  for (m = -n; m <=n; ++m)
675  {
676  ind = m + n;
677  A_mn = Ainit[n][ind];
678  C_mn = Binit[n][ind];
679 
680  B_mn = ( (a_n - gamma_n)*Ainit[n][ind] + m_b*Binit[n][ind])/delta_n;
681  D_mn = ( m_c*Ainit[n][ind] + (d_n - gamma_n)*Binit[n][ind])/delta_n;
682 
683  Spericharmonic = boost::math::spherical_harmonic(n, m, theta, phi);
684  varphi0 += exp(gamma_n*time)*(A_mn*cosh(delta_n*time) + B_mn*sinh(delta_n*time))*Spericharmonic;
685  varphi1 += exp(gamma_n*time)*(C_mn*cosh(delta_n*time) + D_mn*sinh(delta_n*time))*Spericharmonic;
686  }
687  }
688 
689  u[i] = varphi0.real();
690  v[i] = varphi1.real();
691  }
692 
693  switch (field)
694  {
695  case 0:
696  {
697  outfield = u;
698  }
699  break;
700 
701  case 1:
702  {
703  outfield = v;
704  }
705  break;
706  }
707  }
NekDouble m_mu
SOLVER_UTILS_EXPORT int GetTotPoints()
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

◆ PlanePhiWave()

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

Definition at line 710 of file MMFDiffusion.cpp.

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

Referenced by v_SetInitialConditions().

711  {
712  int nq = GetTotPoints();
713  Array<OneD, NekDouble> outarray(nq,0.0);
714 
718 
719  m_fields[0]->GetCoords(x,y,z);
720 
721  NekDouble xmin, ymin, xmax;
722 
723  xmin = Vmath::Vmin(nq, x, 1);
724  xmax = Vmath::Vmax(nq, x, 1);
725  ymin = Vmath::Vmin(nq, y, 1);
726 
727  NekDouble xp, yp, xp2;
728  for (int i=0; i<nq; i++)
729  {
730  switch(m_InitWaveType)
731  {
732  case eLeft:
733  {
734  NekDouble radiusofinit = 4.0;
735  NekDouble frontstiff = 0.1;
736 
737  xp = x[i] - xmin;
738  outarray[i] = 1.0/( 1.0 + exp( ( xp - radiusofinit)/frontstiff ) );
739  }
740  break;
741 
742  case eBothEnds:
743  {
744  NekDouble radiusofinit = 3.0;
745  NekDouble frontstiff = 0.1;
746 
747  xp = x[i] - xmin;
748  xp2 = x[i] - xmax;
749 
750  outarray[i] = 1.0/( 1.0 + exp( ( sqrt(xp*xp) - radiusofinit)/frontstiff ) ) + 1.0/( 1.0 + exp( ( sqrt(xp2*xp2) - radiusofinit)/frontstiff ) );
751  }
752  break;
753 
754  case eCenter:
755  {
756  NekDouble radiusofinit = 6.0;
757  NekDouble frontstiff = 0.1;
758 
759  // NekDouble xc = 0.5*(Vmath::Vmax(nq, x, 1) + Vmath::Vmin(nq, x, 1));
760 
761  xp = x[i] - xmin;
762  outarray[i] =1.0/( 1.0 + exp( ( xp - radiusofinit)/frontstiff ) );
763  }
764  break;
765 
766  case eLeftBottomCorner:
767  {
768  NekDouble radiusofinit = 6.0;
769  NekDouble frontstiff = 0.1;
770  NekDouble bs = 2.0;
771 
772  xp = x[i] - xmin;
773  yp = y[i] - ymin;
774  outarray[i] = 1.0/( 1.0 + exp( ( sqrt(xp*xp+yp*yp)/bs - radiusofinit)/frontstiff ) );
775  }
776  break;
777 
778  case ePoint:
779  {
780  NekDouble xloc, yloc, zloc, rad;
781  NekDouble radiusofinit = 10.0;
782 
783  xloc = x[i]-m_InitPtx;
784  yloc = y[i]-m_InitPty;
785  zloc = z[i]-m_InitPtz;
786 
787  rad = sqrt(xloc*xloc + yloc*yloc + zloc*zloc);
788 
789  xloc = xloc/radiusofinit;
790  yloc = yloc/radiusofinit;
791  zloc = zloc/radiusofinit;
792 
793  if(rad<radiusofinit)
794  {
795  outarray[i] = exp( -(1.0/2.0)*( xloc*xloc + yloc*yloc + zloc*zloc) ) ;
796  }
797 
798  else
799  {
800  outarray[i] = 0.0;
801  }
802  }
803  break;
804 
805  case eSpiralDock:
806  {
807  NekDouble radiusofinit = 3.0;
808  NekDouble frontstiff = 0.1;
809  xp = x[i] - 4.0;
810  yp = y[i];
811  outarray[i] = (1.0/(1.0+exp(2.0*yp)))*(1.0/(1.0+exp(-2.0*xp)))*( 1.0/( 1.0 + exp( ( xp - radiusofinit)/frontstiff ) ) );
812  }
813  break;
814 
815  default:
816  break;
817  }
818 
819  }
820 
821  return outarray;
822  }
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.cpp:782
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.cpp:874
SOLVER_UTILS_EXPORT int GetTotPoints()
InitWaveType m_InitWaveType
Definition: MMFDiffusion.h:119
double NekDouble
static NekDouble rad(NekDouble x, NekDouble y)
Definition: Interpreter.cpp:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

◆ TestCubeProblem()

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

Definition at line 560 of file MMFDiffusion.cpp.

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

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

563  {
564  int nq = GetTotPoints();
565 
569 
570  m_fields[0]->GetCoords(x,y,z);
571 
572  outfield = Array<OneD, NekDouble> (nq);
573  for (int k=0; k<nq; k++)
574  {
575  outfield[k] = exp(-1.0*m_pi*m_pi*time)*sin(m_pi*x[k])*sin(m_pi*y[k])*sin(m_pi*z[k]);
576  }
577  }
SOLVER_UTILS_EXPORT int GetTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

◆ TestPlaneProblem()

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

Definition at line 540 of file MMFDiffusion.cpp.

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

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

543  {
544  int nq = GetTotPoints();
545 
549 
550  m_fields[0]->GetCoords(x,y,z);
551 
552  outfield = Array<OneD, NekDouble> (nq);
553  for (int k=0; k<nq; k++)
554  {
555  outfield[k] = exp(-1.0*m_pi*m_pi*time)*sin(m_pi*x[k])*cos(m_pi*y[k]);
556 
557  }
558  }
SOLVER_UTILS_EXPORT int GetTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

◆ v_EvaluateExactSolution()

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

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 824 of file MMFDiffusion.cpp.

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

827  {
828  switch(m_TestType)
829  {
830  case eTestPlane:
831  {
832  TestPlaneProblem(time,outfield);
833  }
834  break;
835 
836  case eTestCube:
837  {
838  TestCubeProblem(time,outfield);
839  }
840  break;
841 
842  case eTestLinearSphere:
844  {
845  Morphogenesis(time, field, outfield);
846  }
847  break;
848 
849  case eFHNStandard:
850  case eFHNRogers:
851  case eFHNAlievPanf:
852  {
853  int nq = GetTotPoints();
854  outfield = Array<OneD, NekDouble>(nq, 0.0);
855  }
856  /* Falls through. */
857  default:
858  {
859  EquationSystem::v_EvaluateExactSolution(field,outfield,time);
860  }
861  break;
862  }
863  }
void Morphogenesis(const NekDouble time, unsigned int field, Array< OneD, NekDouble > &outfield)
void TestCubeProblem(const NekDouble time, Array< OneD, NekDouble > &outfield)
SOLVER_UTILS_EXPORT int GetTotPoints()
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)

◆ v_GenerateSummary()

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

Prints a summary of the model parameters.

Reimplemented from Nektar::SolverUtils::MMFSystem.

Definition at line 865 of file MMFDiffusion.cpp.

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

866  {
869  SolverUtils::AddSummaryItem(s, "epsilon0", m_epsilon[0]);
870  SolverUtils::AddSummaryItem(s, "epsilon1", m_epsilon[1]);
871  SolverUtils::AddSummaryItem(s, "epsilon2", m_epsilon[2]);
873  {
874  SolverUtils::AddSummaryItem(s, "epsilon for u", m_epsu[0]);
875  SolverUtils::AddSummaryItem(s, "epsilon for v", m_epsu[1]);
876  }
877  }
const char *const TestTypeMap[]
Definition: MMFDiffusion.h:58
Array< OneD, NekDouble > m_epsu
Definition: MMFDiffusion.h:169
Array< OneD, NekDouble > m_epsilon
Definition: MMFDiffusion.h:168
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:49
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
Definition: MMFSystem.cpp:2492

◆ v_InitObject()

void Nektar::MMFDiffusion::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 66 of file MMFDiffusion.cpp.

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, Nektar::TestTypeMap, Nektar::SolverUtils::UnsteadySystem::v_InitObject(), Vmath::Vcopy(), and Vmath::Vvtvp().

67  {
69 
70  int nq = m_fields[0]->GetNpoints();
71  int nvar = m_fields.num_elements();
72  int MFdim = 3;
73 
74  // Diffusivity coefficient for e^j
76  m_session->LoadParameter("epsilon0", m_epsilon[0], 1.0);
77  m_session->LoadParameter("epsilon1", m_epsilon[1], 1.0);
78  m_session->LoadParameter("epsilon2", m_epsilon[2], 1.0);
79 
80  // Diffusivity coefficient for u^j
82  m_session->LoadParameter("epsu0", m_epsu[0], 1.0);
83  m_session->LoadParameter("epsu1", m_epsu[1], 1.0);
84 
85  m_session->LoadParameter("InitPtx", m_InitPtx, 0.0);
86  m_session->LoadParameter("InitPty", m_InitPty, 0.0);
87  m_session->LoadParameter("InitPtz", m_InitPtz, 0.0);
88 
89  int shapedim = m_fields[0]->GetShapeDimension();
90  Array<OneD, Array<OneD, NekDouble> > Anisotropy(shapedim);
91  for(int j=0; j<shapedim; ++j)
92  {
93  Anisotropy[j] = Array<OneD, NekDouble>(nq,1.0);
94  Vmath::Fill(nq, sqrt(m_epsilon[j]), &Anisotropy[j][0], 1);
95 
96  }
97 
98  MMFSystem::MMFInitObject(Anisotropy);
99 
100  // Define ProblemType
101  if(m_session->DefinesSolverInfo("TESTTYPE"))
102  {
103  std::string TestTypeStr = m_session->GetSolverInfo("TESTTYPE");
104  int i;
105  for(i = 0; i < (int) SIZE_TestType; ++i)
106  {
107  if(boost::iequals(TestTypeMap[i],TestTypeStr))
108  {
109  m_TestType = (TestType)i;
110  break;
111  }
112  }
113  }
114  else
115  {
116  m_TestType = (TestType)0;
117  }
118 
119  if(m_session->DefinesSolverInfo("INITWAVETYPE"))
120  {
121  std::string InitWaveTypeStr = m_session->GetSolverInfo("INITWAVETYPE");
122  for(int i = 0; i < (int) SIZE_TestType; ++i)
123  {
124  if(boost::iequals(InitWaveTypeMap[i],InitWaveTypeStr))
125  {
127  break;
128  }
129  }
130  }
131  else
132  {
134  }
135 
136 
152 
153  int indx;
154  Array<OneD, NekDouble> tmp(nq);
155  for (int k=0; k<MFdim; ++k)
156  {
157  // For Moving Frames
158  indx = 5*k;
159 
160  for (int j=0; j<m_spacedim; ++j)
161  {
162  m_varcoeff[MMFCoeffs[indx+j]] = Array<OneD, NekDouble>(nq, 0.0);
163  Vmath::Vcopy(nq, &m_movingframes[k][j*nq], 1, &m_varcoeff[MMFCoeffs[indx+j]][0], 1);
164  }
165 
166  // m_DivMF
167  m_varcoeff[MMFCoeffs[indx+3]] = Array<OneD, NekDouble>(nq, 0.0);
168  Vmath::Vcopy(nq, &m_DivMF[k][0], 1, &m_varcoeff[MMFCoeffs[indx+3]][0], 1);
169 
170  // \| e^k \|
171  m_varcoeff[MMFCoeffs[indx+4]] = Array<OneD, NekDouble>(nq,0.0);
172  tmp = Array<OneD, NekDouble>(nq,0.0);
173  for (int i=0; i<m_spacedim; ++i)
174  {
175  Vmath::Vvtvp(nq, &m_movingframes[k][i*nq], 1, &m_movingframes[k][i*nq], 1, &tmp[0], 1, &tmp[0], 1);
176  }
177 
178  Vmath::Vcopy(nq, &tmp[0], 1, &m_varcoeff[MMFCoeffs[indx+4]][0], 1);
179 
180 
181  }
182 
183 
184  if (!m_explicitDiffusion)
185  {
187  }
188 
190  }
const char *const TestTypeMap[]
Definition: MMFDiffusion.h:58
Array< OneD, NekDouble > m_epsu
Definition: MMFDiffusion.h:169
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
const char *const InitWaveTypeMap[]
Definition: MMFDiffusion.h:80
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45
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.cpp:445
Array< OneD, NekDouble > m_epsilon
Definition: MMFDiffusion.h:168
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.
Array< OneD, Array< OneD, NekDouble > > m_movingframes
Definition: MMFSystem.h:180
Array< OneD, Array< OneD, NekDouble > > m_DivMF
Definition: MMFSystem.h:189
InitWaveType m_InitWaveType
Definition: MMFDiffusion.h:119
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
int m_spacedim
Spatial dimension (>= expansion dim).
StdRegions::VarCoeffMap m_varcoeff
Variable diffusivity.
Definition: MMFDiffusion.h:166
Length of enum list.
Definition: MMFDiffusion.h:55
SOLVER_UTILS_EXPORT void MMFInitObject(const Array< OneD, const Array< OneD, NekDouble >> &Anisotropy, const int TangentXelem=-1)
Definition: MMFSystem.cpp:53
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Computes the reaction terms and .
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_SetInitialConditions()

void Nektar::MMFDiffusion::v_SetInitialConditions ( NekDouble  initialtime,
bool  dumpInitialConditions,
const int  domain 
)
protectedvirtual

Sets a custom initial condition.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 457 of file MMFDiffusion.cpp.

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

460  {
461  boost::ignore_unused(domain);
462 
463  int nq = GetTotPoints();
464 
465  switch(m_TestType)
466  {
467  case eTestPlane:
468  {
470 
471  TestPlaneProblem(initialtime,u);
472  m_fields[0]->SetPhys(u);
473 
474  }
475  break;
476 
477  case eTestCube:
478  {
480 
481  TestCubeProblem(initialtime,u);
482  m_fields[0]->SetPhys(u);
483  /*for (int k=0; k<nq; ++k)
484  {
485  //for (int j=0; j<m_spacedim; ++j)
486  //{
487  cout << "_varcoeff" << u[k] <<endl;
488  // }
489  }*/
490 
491  }
492  break;
493 
494  case eTestLinearSphere:
496  {
499 
500  Morphogenesis(initialtime,0,u);
501  Morphogenesis(initialtime,1,v);
502 
503  m_fields[0]->SetPhys(u);
504  m_fields[1]->SetPhys(v);
505  }
506  break;
507 
508  case eFHNStandard:
509  case eFHNRogers:
510  case eFHNAlievPanf:
511  {
513  m_fields[0]->SetPhys(PlanePhiWave());
514  m_fields[1]->SetPhys(Zero);
515  }
516  break;
517 
518  default:
519  {
520  EquationSystem::v_SetInitialConditions(initialtime,false);
521  }
522  break;
523  }
524 
525  // forward transform to fill the modal coeffs
526  for(int i = 0; i < m_fields.num_elements(); ++i)
527  {
528  m_fields[i]->SetPhysState(true);
529  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),m_fields[i]->UpdateCoeffs());
530  }
531 
532  if(dumpInitialConditions)
533  {
534  std::string outname = m_sessionName + "_initial.chk";
535  WriteFld(outname);
536  }
537  }
void Morphogenesis(const NekDouble time, unsigned int field, Array< OneD, NekDouble > &outfield)
void TestCubeProblem(const NekDouble time, Array< OneD, NekDouble > &outfield)
std::string m_sessionName
Name of the session.
SOLVER_UTILS_EXPORT int GetTotPoints()
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
void TestPlaneProblem(const NekDouble time, Array< OneD, NekDouble > &outfield)
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, NekDouble > PlanePhiWave()
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376

Friends And Related Function Documentation

◆ MemoryManager< MMFDiffusion >

friend class MemoryManager< MMFDiffusion >
friend

Definition at line 94 of file MMFDiffusion.h.

Member Data Documentation

◆ className

string Nektar::MMFDiffusion::className
static
Initial value:
RegisterCreatorFunction("MMFDiffusion",
"MMFDiffusion equation.")

Name of class.

Definition at line 107 of file MMFDiffusion.h.

◆ m_epsilon

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

Definition at line 168 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 169 of file MMFDiffusion.h.

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

◆ m_InitPtx

NekDouble Nektar::MMFDiffusion::m_InitPtx
protected

Definition at line 162 of file MMFDiffusion.h.

Referenced by PlanePhiWave(), and v_InitObject().

◆ m_InitPty

NekDouble Nektar::MMFDiffusion::m_InitPty
protected

Definition at line 162 of file MMFDiffusion.h.

Referenced by PlanePhiWave(), and v_InitObject().

◆ m_InitPtz

NekDouble Nektar::MMFDiffusion::m_InitPtz
protected

Definition at line 162 of file MMFDiffusion.h.

Referenced by PlanePhiWave(), and v_InitObject().

◆ m_InitWaveType

InitWaveType Nektar::MMFDiffusion::m_InitWaveType
protected

Definition at line 119 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 166 of file MMFDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().