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

Class for iterative elastic system, in which linear elasticity is applied in substeps to attain a large deformation. More...

#include <IterativeElasticSystem.h>

Inheritance diagram for Nektar::IterativeElasticSystem:
[legend]

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::LinearElasticSystem
static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::LinearElasticSystem
static std::string className
 Name of class. More...
 

Protected Member Functions

 IterativeElasticSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void v_InitObject (bool DeclareFields=true) override
 Initialisation object for EquationSystem. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Virtual function for generating summary information. More...
 
void v_DoSolve () override
 Virtual function for solve implementation. More...
 
void WriteGeometry (const int i)
 Write out a file in serial or directory in parallel containing new mesh geometry. More...
 
- Protected Member Functions inherited from Nektar::LinearElasticSystem
 LinearElasticSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Default constructor. More...
 
void v_InitObject (bool DeclareFields=true) override
 Set up the linear elasticity system. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Generate summary at runtime. More...
 
void v_DoSolve () override
 Solve elliptic linear elastic system. More...
 
void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) override
 
- 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

int m_numSteps
 Number of steps to split deformation across. More...
 
bool m_repeatBCs
 Flag determining whether to repeat boundary conditions. More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_savedBCs
 Storage for boundary conditions. More...
 
std::vector< int > m_toDeform
 Vector of boundary regions to deform. More...
 
- Protected Attributes inherited from Nektar::LinearElasticSystem
NekDouble m_nu
 Poisson ratio. More...
 
NekDouble m_E
 Young's modulus. More...
 
NekDouble m_beta
 Parameter dictating amount of thermal stress to add. More...
 
CoupledAssemblyMapSharedPtr m_assemblyMap
 Assembly map for the coupled (u,v,w) system. More...
 
DNekScalBlkMatSharedPtr m_schurCompl
 Schur complement boundary-boundary matrix. More...
 
DNekScalBlkMatSharedPtr m_BinvD
 Matrix of elemental \( B^{-1}D \) components. More...
 
DNekScalBlkMatSharedPtr m_C
 Matrix of elemental \( C \) components. More...
 
DNekScalBlkMatSharedPtr m_Dinv
 Matrix of elemental \( D^{-1} \) components. More...
 
Array< OneD, Array< OneD, unsigned int > > m_bmap
 Boundary maps for each of the fields. More...
 
Array< OneD, Array< OneD, unsigned int > > m_imap
 Interior maps for each of the fields. More...
 
Array< OneD, Array< OneD, NekDouble > > m_temperature
 Storage for the temperature terms. More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_stress
 Storage for the thermal stress terms. 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_movingFrameVelsxyz
 Moving frame of reference velocities (u, v, w, omega_x, omega_y, omega_z, a_x, a_y, a_z, domega_x, domega_y, domega_z) More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. 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...
 

Friends

class MemoryManager< IterativeElasticSystem >
 Class may only be instantiated through the MemoryManager. More...
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::LinearElasticSystem
void BuildMatrixSystem ()
 Build matrix system for linear elasticity equations. More...
 
void SetStaticCondBlock (const int n, const LocalRegions::ExpansionSharedPtr exp, Array< TwoD, DNekMatSharedPtr > &mat)
 Given a block matrix for an element, construct its static condensation matrices. More...
 
DNekMatSharedPtr BuildLaplacianIJMatrix (const int k1, const int k2, const NekDouble scale, LocalRegions::ExpansionSharedPtr exp)
 Construct a LaplacianIJ matrix for a given expansion. 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...
 
- 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

Class for iterative elastic system, in which linear elasticity is applied in substeps to attain a large deformation.

Definition at line 50 of file IterativeElasticSystem.h.

Constructor & Destructor Documentation

◆ IterativeElasticSystem()

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

Definition at line 57 of file IterativeElasticSystem.cpp.

60 : LinearElasticSystem(pSession, pGraph)
61{
62}
LinearElasticSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Default constructor.

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 57 of file IterativeElasticSystem.h.

60 {
63 pGraph);
64 p->InitObject();
65 return p;
66 }
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.

◆ v_DoSolve()

void Nektar::IterativeElasticSystem::v_DoSolve ( void  )
overrideprotectedvirtual

Virtual function for solve implementation.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 127 of file IterativeElasticSystem.cpp.

128{
129 int i, j, k;
130
131 // Write initial geometry for consistency/script purposes
132 WriteGeometry(0);
133
134 Array<OneD, Array<OneD, NekDouble>> physvals(m_fields.size());
135
136 for (int i = 0; i < m_fields.size(); ++i)
137 {
138 physvals[i] = m_fields[i]->UpdatePhys();
139 }
140
141 // Now loop over desired number of steps
142 for (i = 1; i <= m_numSteps; ++i)
143 {
144 int invalidElmtId = -1;
145
146 // Perform solve for this iteration and update geometry accordingly.
149 WriteGeometry(i);
150
151 // Check for invalid elements.
152 for (j = 0; j < m_fields[0]->GetExpSize(); ++j)
153 {
155 m_fields[0]->GetExp(j)->GetGeom()->GetGeomFactors();
156
157 if (!geomFac->IsValid())
158 {
159 invalidElmtId =
160 m_fields[0]->GetExp(j)->GetGeom()->GetGlobalID();
161 break;
162 }
163 }
164
165 m_session->GetComm()->AllReduce(invalidElmtId, LibUtilities::ReduceMax);
166
167 // If we found an invalid element, exit loop without writing output.
168 if (invalidElmtId >= 0)
169 {
170 if (m_session->GetComm()->GetRank() == 0)
171 {
172 cout << "- Detected negative Jacobian in element "
173 << invalidElmtId
174 << "; terminating at"
175 " step: "
176 << i << endl;
177 }
178
179 break;
180 }
181
182 if (m_session->GetComm()->GetRank() == 0)
183 {
184 cout << "Step: " << i << endl;
185 }
186
187 // Update boundary conditions
188 if (m_repeatBCs)
189 {
190 for (j = 0; j < m_fields.size(); ++j)
191 {
192 string varName = m_session->GetVariable(j);
193 m_fields[j]->EvaluateBoundaryConditions(m_time, varName);
194 }
195 }
196 else
197 {
198 for (j = 0; j < m_fields.size(); ++j)
199 {
200 const Array<OneD, const MultiRegions::ExpListSharedPtr>
201 &bndCondExp = m_fields[j]->GetBndCondExpansions();
202
203 for (k = 0; k < m_toDeform.size(); ++k)
204 {
205 const int id = m_toDeform[k];
206 const int nCoeffs = bndCondExp[id]->GetNcoeffs();
207 Vmath::Vcopy(nCoeffs, m_savedBCs[k][j], 1,
208 bndCondExp[id]->UpdateCoeffs(), 1);
209 }
210 }
211 }
212 }
213}
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_savedBCs
Storage for boundary conditions.
void WriteGeometry(const int i)
Write out a file in serial or directory in parallel containing new mesh geometry.
bool m_repeatBCs
Flag determining whether to repeat boundary conditions.
std::vector< int > m_toDeform
Vector of boundary regions to deform.
int m_numSteps
Number of steps to split deformation across.
void v_DoSolve() override
Solve elliptic linear elastic system.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void UpdateGeometry(SpatialDomains::MeshGraphSharedPtr graph, Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, Array< OneD, NekDouble > > &PhysVals, bool modal)
Update geometry according to displacement that is in current fields.
Definition: Deform.cpp:59
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:60
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_graph, m_numSteps, m_repeatBCs, m_savedBCs, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_time, m_toDeform, Nektar::LibUtilities::ReduceMax, Nektar::GlobalMapping::UpdateGeometry(), Nektar::LinearElasticSystem::v_DoSolve(), Vmath::Vcopy(), and WriteGeometry().

◆ v_GenerateSummary()

void Nektar::IterativeElasticSystem::v_GenerateSummary ( SolverUtils::SummaryList l)
overrideprotectedvirtual

Virtual function for generating summary information.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 122 of file IterativeElasticSystem.cpp.

123{
125}
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Generate summary at runtime.

References Nektar::LinearElasticSystem::v_GenerateSummary().

◆ v_InitObject()

void Nektar::IterativeElasticSystem::v_InitObject ( bool  DeclareFeld = true)
overrideprotectedvirtual

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 64 of file IterativeElasticSystem.cpp.

65{
67
68 const int nVel = m_fields[0]->GetCoordim(0);
69
70 // Read in number of steps to take.
71 m_session->LoadParameter("NumSteps", m_numSteps, 0);
72 ASSERTL0(m_numSteps > 0, "You must specify at least one step");
73
74 // Read in whether to repeatedly apply boundary conditions (for e.g.
75 // rotation purposes).
76 string bcType;
77 m_session->LoadSolverInfo("BCType", bcType, "Normal");
78 m_repeatBCs = bcType != "Normal";
79
80 if (!m_repeatBCs)
81 {
82 // Loop over BCs, identify which ones we need to deform.
83 const Array<OneD, const SpatialDomains::BoundaryConditionShPtr>
84 &bndCond = m_fields[0]->GetBndConditions();
85
86 for (int i = 0; i < bndCond.size(); ++i)
87 {
88 if (boost::iequals(bndCond[i]->GetUserDefined(), "Wall"))
89 {
90 m_toDeform.push_back(i);
91 }
92 }
93
94 int numDeform = m_toDeform.size();
95 m_comm->AllReduce(numDeform, LibUtilities::ReduceMax);
96 ASSERTL0(numDeform > 0, "You must specify at least one WALL tag on"
97 "a boundary condition");
98
100 Array<OneD, Array<OneD, Array<OneD, NekDouble>>>(m_toDeform.size());
101
102 for (int i = 0; i < m_toDeform.size(); ++i)
103 {
104 m_savedBCs[i] = Array<OneD, Array<OneD, NekDouble>>(nVel);
105 for (int j = 0; j < nVel; ++j)
106 {
107 const int id = m_toDeform[i];
109 m_fields[j]->GetBndCondExpansions()[id];
110 int nCoeffs = bndCondExp->GetNcoeffs();
111
112 m_savedBCs[i][j] = Array<OneD, NekDouble>(nCoeffs);
113 Vmath::Smul(nCoeffs, 1.0 / m_numSteps, bndCondExp->GetCoeffs(),
114 1, bndCondExp->UpdateCoeffs(), 1);
115 Vmath::Vcopy(nCoeffs, bndCondExp->GetCoeffs(), 1,
116 m_savedBCs[i][j], 1);
117 }
118 }
119 }
120}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void v_InitObject(bool DeclareFields=true) override
Set up the linear elasticity system.
LibUtilities::CommSharedPtr m_comm
Communicator.
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
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 ASSERTL0, Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, m_numSteps, m_repeatBCs, m_savedBCs, Nektar::SolverUtils::EquationSystem::m_session, m_toDeform, Nektar::LibUtilities::ReduceMax, Vmath::Smul(), Nektar::LinearElasticSystem::v_InitObject(), and Vmath::Vcopy().

◆ WriteGeometry()

void Nektar::IterativeElasticSystem::WriteGeometry ( const int  i)
protected

Write out a file in serial or directory in parallel containing new mesh geometry.

Definition at line 219 of file IterativeElasticSystem.cpp.

220{
221 fs::path filename;
222 stringstream s;
223 s << m_session->GetSessionName() << "-" << i;
224
225 if (m_session->GetComm()->GetSize() > 1)
226 {
227 s << "_xml";
228
229 string ss = s.str();
230 if (!fs::is_directory(ss))
231 {
232 fs::create_directory(ss);
233 }
234
235 boost::format pad("P%1$07d.xml");
236 pad % m_session->GetComm()->GetRank();
237 filename = fs::path(ss) / fs::path(pad.str());
238 }
239 else
240 {
241 s << ".xml";
242 filename = fs::path(s.str());
243 }
244
245 string fname = LibUtilities::PortablePath(filename);
246 m_fields[0]->GetGraph()->WriteGeometry(fname);
247}
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56

References CellMLToNektar.pycml::format, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::LibUtilities::PortablePath().

Referenced by v_DoSolve().

Friends And Related Function Documentation

◆ MemoryManager< IterativeElasticSystem >

friend class MemoryManager< IterativeElasticSystem >
friend

Class may only be instantiated through the MemoryManager.

Definition at line 1 of file IterativeElasticSystem.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 69 of file IterativeElasticSystem.h.

◆ m_numSteps

int Nektar::IterativeElasticSystem::m_numSteps
protected

Number of steps to split deformation across.

Definition at line 73 of file IterativeElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_repeatBCs

bool Nektar::IterativeElasticSystem::m_repeatBCs
protected

Flag determining whether to repeat boundary conditions.

Definition at line 75 of file IterativeElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_savedBCs

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::IterativeElasticSystem::m_savedBCs
protected

Storage for boundary conditions.

Definition at line 77 of file IterativeElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_toDeform

std::vector<int> Nektar::IterativeElasticSystem::m_toDeform
protected

Vector of boundary regions to deform.

Definition at line 79 of file IterativeElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().