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

Base class for linear elastic system. More...

#include <LinearElasticSystem.h>

Inheritance diagram for Nektar::LinearElasticSystem:
[legend]

Public Member Functions

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

Static Public Member Functions

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

Static Public Attributes

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

Protected Member Functions

 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

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_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...
 

Friends

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

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

Base class for linear elastic system.

Definition at line 49 of file LinearElasticSystem.h.

Constructor & Destructor Documentation

◆ LinearElasticSystem()

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

Default constructor.

Definition at line 136 of file LinearElasticSystem.cpp.

139 : EquationSystem(pSession, pGraph)
140{
141}
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises EquationSystem class members.

Member Function Documentation

◆ BuildLaplacianIJMatrix()

DNekMatSharedPtr Nektar::LinearElasticSystem::BuildLaplacianIJMatrix ( const int  k1,
const int  k2,
const NekDouble  scale,
LocalRegions::ExpansionSharedPtr  exp 
)

Construct a LaplacianIJ matrix for a given expansion.

This routine constructs matrices whose entries contain the evaluation of

\[ \partial_{\tt k1} \phi_i \partial_{\tt k2} \phi_j \,dx \]

Definition at line 786 of file LinearElasticSystem.cpp.

789{
790 const int nCoeffs = exp->GetNcoeffs();
791 const int nPhys = exp->GetTotPoints();
792 int i;
793
794 DNekMatSharedPtr ret =
795 MemoryManager<DNekMat>::AllocateSharedPtr(nCoeffs, nCoeffs, 0.0, eFULL);
796
797 Array<OneD, NekDouble> tmp2(nPhys);
798 Array<OneD, NekDouble> tmp3(nPhys);
799
800 for (i = 0; i < nCoeffs; ++i)
801 {
802 Array<OneD, NekDouble> tmp1(nCoeffs, 0.0);
803 tmp1[i] = 1.0;
804
805 exp->BwdTrans(tmp1, tmp2);
806 exp->PhysDeriv(k1, tmp2, tmp3);
807 exp->IProductWRTDerivBase(k2, tmp3, tmp1);
808
809 Vmath::Smul(nCoeffs, scale, &tmp1[0], 1,
810 &(ret->GetPtr())[0] + i * nCoeffs, 1);
811 }
812
813 return ret;
814}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eFULL, and Vmath::Smul().

Referenced by BuildMatrixSystem().

◆ BuildMatrixSystem()

void Nektar::LinearElasticSystem::BuildMatrixSystem ( )

Build matrix system for linear elasticity equations.

This routine constructs the matrix discretisation arising from the weak formulation of the linear elasticity equations. The resulting matrix system is then passed to LinearElasticSystem::SetStaticCondBlock in order to construct the statically condensed system.

All of the matrices involved in the construction of the divergence of the stress tensor are Laplacian-like. We use the variable coefficient functionality within the library to multiply the Laplacian by the appropriate constants for all diagonal terms, and off-diagonal terms use LinearElasticSystem::BuildLaplacianIJMatrix to construct matrices containing

\[ \int \partial_{x_i} \phi_k \partial_{x_j} \phi_l \]

where mixed derivative terms are present. Symmetry (in terms of k,l) is exploited to avoid constructing this matrix repeatedly.

Todo:
Make static condensation optional and construct full system instead.

Definition at line 230 of file LinearElasticSystem.cpp.

231{
232 const int nEl = m_fields[0]->GetExpSize();
233 const int nVel = m_fields[0]->GetCoordim(0);
234
236 int n;
237
238 // Factors map for matrix keys.
240
241 // Calculate various constants
242 NekDouble mu = m_E * 0.5 / (1.0 + m_nu);
243 NekDouble lambda = m_E * m_nu / (1.0 + m_nu) / (1.0 - 2.0 * m_nu);
244
245 bool verbose = m_session->DefinesCmdLineArgument("verbose");
246 bool root = m_session->GetComm()->GetRank() == 0;
247
248 // Loop over each element and construct matrices.
249 if (nVel == 2)
250 {
251 for (n = 0; n < nEl; ++n)
252 {
253 exp = m_fields[0]->GetExp(n);
254 const int nPhys = exp->GetTotPoints();
255 Array<OneD, NekDouble> aArr(nPhys, lambda + 2 * mu);
256 Array<OneD, NekDouble> bArr(nPhys, mu);
257
258 StdRegions::VarCoeffMap varcoeffA, varcoeffD;
259 varcoeffA[StdRegions::eVarCoeffD00] = aArr;
260 varcoeffA[StdRegions::eVarCoeffD11] = bArr;
261 varcoeffD[StdRegions::eVarCoeffD00] = bArr;
262 varcoeffD[StdRegions::eVarCoeffD11] = aArr;
263
264 LocalRegions::MatrixKey matkeyA(StdRegions::eLaplacian,
265 exp->DetShapeType(), *exp, factors,
266 varcoeffA);
267 LocalRegions::MatrixKey matkeyD(StdRegions::eLaplacian,
268 exp->DetShapeType(), *exp, factors,
269 varcoeffD);
270
271 /*
272 * mat holds the linear operator [ A B ] acting on [ u ].
273 * [ C D ] [ v ]
274 */
275 Array<TwoD, DNekMatSharedPtr> mat(2, 2);
276 mat[0][0] = exp->GenMatrix(matkeyA);
277 mat[0][1] = BuildLaplacianIJMatrix(1, 0, mu + lambda, exp);
278 mat[1][0] = mat[0][1];
279 mat[1][1] = exp->GenMatrix(matkeyD);
280
281 // Set up the statically condensed block for this element.
282 SetStaticCondBlock(n, exp, mat);
283
284 if (verbose && root)
285 {
286 cout << "\rBuilding matrix system: " << (int)(100.0 * n / nEl)
287 << "%" << flush;
288 }
289 }
290 }
291 else if (nVel == 3)
292 {
293 for (n = 0; n < nEl; ++n)
294 {
295 exp = m_fields[0]->GetExp(n);
296 const int nPhys = exp->GetTotPoints();
297 Array<OneD, NekDouble> aArr(nPhys, lambda + 2 * mu);
298 Array<OneD, NekDouble> bArr(nPhys, mu);
299
300 StdRegions::VarCoeffMap varcoeffA, varcoeffE, varcoeffI;
301 varcoeffA[StdRegions::eVarCoeffD00] = aArr;
302 varcoeffA[StdRegions::eVarCoeffD11] = bArr;
303 varcoeffA[StdRegions::eVarCoeffD22] = bArr;
304 varcoeffE[StdRegions::eVarCoeffD00] = bArr;
305 varcoeffE[StdRegions::eVarCoeffD11] = aArr;
306 varcoeffE[StdRegions::eVarCoeffD22] = bArr;
307 varcoeffI[StdRegions::eVarCoeffD00] = bArr;
308 varcoeffI[StdRegions::eVarCoeffD11] = bArr;
309 varcoeffI[StdRegions::eVarCoeffD22] = aArr;
310
311 LocalRegions::MatrixKey matkeyA(StdRegions::eLaplacian,
312 exp->DetShapeType(), *exp, factors,
313 varcoeffA);
314 LocalRegions::MatrixKey matkeyE(StdRegions::eLaplacian,
315 exp->DetShapeType(), *exp, factors,
316 varcoeffE);
317 LocalRegions::MatrixKey matkeyI(StdRegions::eLaplacian,
318 exp->DetShapeType(), *exp, factors,
319 varcoeffI);
320
321 /*
322 * mat holds the linear operator [ A B C ] acting on [ u ].
323 * [ D E F ] [ v ]
324 * [ G H I ] [ w ]
325 */
326 Array<TwoD, DNekMatSharedPtr> mat(3, 3);
327 mat[0][0] = exp->GenMatrix(matkeyA);
328 mat[0][1] = BuildLaplacianIJMatrix(1, 0, mu + lambda, exp);
329 mat[0][2] = BuildLaplacianIJMatrix(2, 0, mu + lambda, exp);
330
331 mat[1][0] = mat[0][1];
332 mat[1][1] = exp->GenMatrix(matkeyE);
333 mat[1][2] = BuildLaplacianIJMatrix(2, 1, mu + lambda, exp);
334
335 mat[2][0] = mat[0][2];
336 mat[2][1] = mat[1][2];
337 mat[2][2] = exp->GenMatrix(matkeyI);
338
339 // Set up the statically condensed block for this element.
340 SetStaticCondBlock(n, exp, mat);
341
342 if (verbose && root)
343 {
344 cout << "\rBuilding matrix system: " << (int)(100.0 * n / nEl)
345 << "%" << flush;
346 }
347 }
348 }
349
350 if (verbose && root)
351 {
352 cout << "\rBuilding matrix system: done." << endl;
353 }
354}
NekDouble m_nu
Poisson ratio.
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.
NekDouble m_E
Young's modulus.
DNekMatSharedPtr BuildLaplacianIJMatrix(const int k1, const int k2, const NekDouble scale, LocalRegions::ExpansionSharedPtr exp)
Construct a LaplacianIJ matrix for a given expansion.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
double NekDouble

References BuildLaplacianIJMatrix(), Nektar::StdRegions::eLaplacian, Nektar::StdRegions::eVarCoeffD00, Nektar::StdRegions::eVarCoeffD11, Nektar::StdRegions::eVarCoeffD22, Nektar::VarcoeffHashingTest::factors, m_E, Nektar::SolverUtils::EquationSystem::m_fields, m_nu, Nektar::SolverUtils::EquationSystem::m_session, and SetStaticCondBlock().

Referenced by v_DoSolve().

◆ create()

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

Creates an instance of this class.

Definition at line 56 of file LinearElasticSystem.h.

59 {
62 pGraph);
63 p->InitObject();
64 return p;
65 }
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

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

◆ SetStaticCondBlock()

void Nektar::LinearElasticSystem::SetStaticCondBlock ( const int  n,
const LocalRegions::ExpansionSharedPtr  exp,
Array< TwoD, DNekMatSharedPtr > &  mat 
)

Given a block matrix for an element, construct its static condensation matrices.

This routine essentially duplicates the logic present in many of the LocalRegions matrix generation routines to construct the statically condensed equivalent of mat to pass to the GlobalLinSys solver.

Parameters
nElement/block number.
expPointer to expansion.
matBlock matrix containing matrix operator.

Definition at line 703 of file LinearElasticSystem.cpp.

706{
707 int i, j, k, l;
708 const int nVel = mat.GetRows();
709 const int nB = exp->NumBndryCoeffs();
710 const int nI = exp->GetNcoeffs() - nB;
711 const int nBnd = exp->NumBndryCoeffs() * nVel;
712 const int nInt = exp->GetNcoeffs() * nVel - nBnd;
713 const MatrixStorage s = eFULL; // Maybe look into doing symmetric
714 // version of this?
715
724
725 for (i = 0; i < nVel; ++i)
726 {
727 for (j = 0; j < nVel; ++j)
728 {
729 // Boundary-boundary and boundary-interior
730 for (k = 0; k < nB; ++k)
731 {
732 for (l = 0; l < nB; ++l)
733 {
734 (*A)(k + i * nB, l + j * nB) =
735 (*mat[i][j])(m_bmap[n][k], m_bmap[n][l]);
736 }
737
738 for (l = 0; l < nI; ++l)
739 {
740 (*B)(k + i * nB, l + j * nI) =
741 (*mat[i][j])(m_bmap[n][k], m_imap[n][l]);
742 }
743 }
744
745 // Interior-boundary / interior-interior
746 for (k = 0; k < nI; ++k)
747 {
748 for (l = 0; l < nB; ++l)
749 {
750 (*C)(k + i * nI, l + j * nB) =
751 (*mat[i][j])(m_imap[n][k], m_bmap[n][l]);
752 }
753
754 for (l = 0; l < nI; ++l)
755 {
756 (*D)(k + i * nI, l + j * nI) =
757 (*mat[i][j])(m_imap[n][k], m_imap[n][l]);
758 }
759 }
760 }
761 }
762
763 // Construct static condensation matrices.
764 D->Invert();
765 (*B) = (*B) * (*D);
766 (*A) = (*A) - (*B) * (*C);
767
768 DNekScalMatSharedPtr tmp_mat;
769 m_schurCompl->SetBlock(
771 m_BinvD->SetBlock(
772 n, n, tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, B));
773 m_C->SetBlock(
774 n, n, tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, C));
775 m_Dinv->SetBlock(
776 n, n, tmp_mat = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, D));
777}
Array< OneD, Array< OneD, unsigned int > > m_imap
Interior maps for each of the fields.
DNekScalBlkMatSharedPtr m_Dinv
Matrix of elemental components.
DNekScalBlkMatSharedPtr m_BinvD
Matrix of elemental components.
DNekScalBlkMatSharedPtr m_schurCompl
Schur complement boundary-boundary matrix.
Array< OneD, Array< OneD, unsigned int > > m_bmap
Boundary maps for each of the fields.
DNekScalBlkMatSharedPtr m_C
Matrix of elemental components.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eFULL, m_BinvD, m_bmap, m_C, m_Dinv, m_imap, and m_schurCompl.

Referenced by BuildMatrixSystem().

◆ v_DoSolve()

void Nektar::LinearElasticSystem::v_DoSolve ( void  )
overrideprotectedvirtual

Solve elliptic linear elastic system.

The solve proceeds as follows:

  • Create a MultiRegions::GlobalLinSys object.
  • Evaluate a forcing function from the session file.
  • If the temperature term is enabled, evaluate this and add to forcing.
  • Apply Dirichlet boundary conditions.
  • Scatter forcing into the correct ordering according to the coupled assembly map.
  • Do the solve.
  • Scatter solution back to fields and backwards transform to physical space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 381 of file LinearElasticSystem.cpp.

382{
383 int i, j, k, l, nv;
384 const int nVel = m_fields[0]->GetCoordim(0);
385
386 // Build initial matrix system.
388
389 // Now we've got the matrix system set up, create a GlobalLinSys
390 // object. We mask ourselves as LinearAdvectionReaction to create a full
391 // matrix instead of symmetric storage.
392 MultiRegions::GlobalLinSysKey key(StdRegions::eLinearAdvectionReaction,
395
396 // Currently either direct or iterative static condensation is
397 // supported.
398 if (m_assemblyMap->GetGlobalSysSolnType() ==
400 {
404 }
405 else if (m_assemblyMap->GetGlobalSysSolnType() ==
407 {
412 }
413#ifdef NEKTAR_USE_PETSC
414 else if (m_assemblyMap->GetGlobalSysSolnType() ==
416 {
420 }
421#endif
422#ifdef NEKTAR_USE_MPI
423 else if (m_assemblyMap->GetGlobalSysSolnType() ==
425 {
429 }
430#endif
431
432 linSys->Initialise(m_assemblyMap);
433
434 const int nCoeffs = m_fields[0]->GetNcoeffs();
435
436 //
437 // -- Evaluate forcing functions
438 //
439
440 // Evaluate the forcing function from the XML file.
441 Array<OneD, Array<OneD, NekDouble>> forcing(nVel);
442 GetFunction("Forcing")->Evaluate(forcing);
443
444 // Add temperature term
445 string tempEval;
446 m_session->LoadSolverInfo("Temperature", tempEval, "None");
447
448 if (tempEval == "Jacobian")
449 {
450 // Allocate storage
451 m_temperature = Array<OneD, Array<OneD, NekDouble>>(nVel);
452
453 for (nv = 0; nv < nVel; ++nv)
454 {
455 Array<OneD, NekDouble> tmp;
456 m_temperature[nv] =
457 Array<OneD, NekDouble>(m_fields[nv]->GetNpoints());
458
459 for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
460 {
461 // Calculate element area
462 LocalRegions::ExpansionSharedPtr exp = m_fields[0]->GetExp(i);
463 LibUtilities::PointsKeyVector pkey = exp->GetPointsKeys();
464 Array<OneD, NekDouble> jac = exp->GetMetricInfo()->GetJac(pkey);
465
466 int offset = m_fields[0]->GetPhys_Offset(i);
467
468 if (exp->GetMetricInfo()->GetGtype() ==
470 {
471 Vmath::Smul(exp->GetTotPoints(), m_beta, jac, 1,
472 tmp = m_temperature[nv] + offset, 1);
473 }
474 else
475 {
476 Vmath::Fill(exp->GetTotPoints(), m_beta * jac[0],
477 tmp = m_temperature[nv] + offset, 1);
478 }
479 }
480 m_fields[nv]->PhysDeriv(nv, m_temperature[nv], forcing[nv]);
481 }
482 }
483 else if (tempEval == "Metric")
484 {
485 ASSERTL0((m_fields[0]->GetCoordim(0) == 2 &&
486 m_graph->GetAllQuadGeoms().size() == 0) ||
487 (m_fields[0]->GetCoordim(0) == 3 &&
488 m_graph->GetAllPrismGeoms().size() == 0 &&
489 m_graph->GetAllPyrGeoms().size() == 0 &&
490 m_graph->GetAllHexGeoms().size() == 0),
491 "LinearIdealMetric temperature only implemented for "
492 "two-dimensional triangular meshes or three-dimensional "
493 "tetrahedral meshes.");
494
495 m_temperature = Array<OneD, Array<OneD, NekDouble>>(nVel);
496 m_stress = Array<OneD, Array<OneD, Array<OneD, NekDouble>>>(nVel);
497
498 for (nv = 0; nv < nVel; ++nv)
499 {
500 m_temperature[nv] =
501 Array<OneD, NekDouble>(m_fields[nv]->GetNpoints(), 0.0);
502 m_stress[nv] = Array<OneD, Array<OneD, NekDouble>>(nVel);
503
504 for (i = 0; i < nVel; ++i)
505 {
506 m_stress[nv][i] =
507 Array<OneD, NekDouble>(m_fields[nv]->GetNpoints(), 0.0);
508 }
509 }
510
511 for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
512 {
513 // Calculate element area
514 LocalRegions::ExpansionSharedPtr exp = m_fields[0]->GetExp(i);
515 LibUtilities::PointsKeyVector pkey = exp->GetPointsKeys();
516 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> deriv =
517 exp->GetMetricInfo()->GetDeriv(pkey);
518 int offset = m_fields[0]->GetPhys_Offset(i);
519
520 DNekMat i2rm = MappingIdealToRef(exp->GetGeom());
521
522 // Compute metric tensor
523 DNekMat jac(nVel, nVel, 0.0, eFULL);
524 DNekMat jacIdeal(nVel, nVel, 0.0, eFULL);
525 DNekMat metric(nVel, nVel, 0.0, eFULL);
526
527 for (j = 0; j < deriv[0][0].size(); ++j)
528 {
529 for (k = 0; k < nVel; ++k)
530 {
531 for (l = 0; l < nVel; ++l)
532 {
533 jac(l, k) = deriv[k][l][j];
534 }
535 }
536
537 jacIdeal = jac * i2rm;
538 metric = Transpose(jacIdeal) * jacIdeal;
539
540 // Compute eigenvalues/eigenvectors of metric tensor using
541 // ideal mapping.
542 char jobvl = 'N', jobvr = 'V';
543 int worklen = 8 * nVel, info;
544
545 DNekMat eval(nVel, nVel, 0.0, eDIAGONAL);
546 DNekMat evec(nVel, nVel, 0.0, eFULL);
547 DNekMat evecinv(nVel, nVel, 0.0, eFULL);
548 Array<OneD, NekDouble> vl(nVel * nVel);
549 Array<OneD, NekDouble> work(worklen);
550 Array<OneD, NekDouble> wi(nVel);
551
552 Lapack::Dgeev(jobvl, jobvr, nVel, metric.GetRawPtr(), nVel,
553 &(eval.GetPtr())[0], &wi[0], &vl[0], nVel,
554 &(evec.GetPtr())[0], nVel, &work[0], worklen,
555 info);
556
557 evecinv = evec;
558 evecinv.Invert();
559
560 // rescaling of the eigenvalues
561 for (nv = 0; nv < nVel; ++nv)
562 {
563 eval(nv, nv) = m_beta * (sqrt(eval(nv, nv)) - 1.0);
564 }
565
566 DNekMat beta = evec * eval * evecinv;
567
568 for (k = 0; k < nVel; ++k)
569 {
570 for (l = 0; l < nVel; ++l)
571 {
572 m_stress[k][l][offset + j] = beta(k, l);
573 }
574 }
575 }
576
577 if (deriv[0][0].size() != exp->GetTotPoints())
578 {
579 Array<OneD, NekDouble> tmp;
580 for (k = 0; k < nVel; ++k)
581 {
582 for (l = 0; l < nVel; ++l)
583 {
584 Vmath::Fill(exp->GetTotPoints(), m_stress[k][l][offset],
585 tmp = m_stress[k][l] + offset, 1);
586 }
587 }
588 }
589 }
590
591 // Calculate divergence of stress tensor.
592 Array<OneD, NekDouble> tmpderiv(m_fields[0]->GetNpoints());
593 for (i = 0; i < nVel; ++i)
594 {
595 for (j = 0; j < nVel; ++j)
596 {
597 m_fields[i]->PhysDeriv(j, m_stress[i][j], tmpderiv);
598 Vmath::Vadd(m_fields[i]->GetNpoints(), tmpderiv, 1,
599 m_temperature[i], 1, m_temperature[i], 1);
600 }
601
603 forcing[i], 1);
604 }
605 }
606 else if (tempEval != "None")
607 {
608 ASSERTL0(false, "Unknown temperature form: " + tempEval);
609 }
610
611 // Set up some temporary storage.
612 //
613 // - forCoeffs holds the forcing coefficients in a local ordering;
614 // however note that the ordering is different and dictated by the
615 // assembly map. We loop over each element, then the boundary degrees of
616 // freedom for u, boundary for v, followed by the interior for u and then
617 // interior for v.
618 // - inout holds the Dirichlet degrees of freedom in the local ordering,
619 // which have the boundary conditions imposed
620 Array<OneD, NekDouble> forCoeffs(nVel * nCoeffs, 0.0);
621 Array<OneD, NekDouble> inout(nVel * nCoeffs, 0.0);
622 Array<OneD, NekDouble> tmp(nCoeffs);
623
624 for (nv = 0; nv < nVel; ++nv)
625 {
626 // Take the inner product of the forcing function.
627 m_fields[nv]->IProductWRTBase(forcing[nv], tmp);
628
629 // Impose Dirichlet condition on field which should be initialised
630 Array<OneD, NekDouble> loc_inout = m_fields[nv]->UpdateCoeffs();
631 m_fields[nv]->ImposeDirichletConditions(loc_inout);
632
633 // Scatter forcing into RHS vector according to the ordering dictated in
634 // the comment above.
635 for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
636 {
637 int nBnd = m_bmap[i].size();
638 int nInt = m_imap[i].size();
639 int offset = m_fields[nv]->GetCoeff_Offset(i);
640
641 for (j = 0; j < nBnd; ++j)
642 {
643 forCoeffs[nVel * offset + nv * nBnd + j] =
644 -1.0 * tmp[offset + m_bmap[i][j]];
645 inout[nVel * offset + nv * nBnd + j] =
646 loc_inout[offset + m_bmap[i][j]];
647 }
648 for (j = 0; j < nInt; ++j)
649 {
650 forCoeffs[nVel * (offset + nBnd) + nv * nInt + j] =
651 -1.0 * tmp[offset + m_imap[i][j]];
652 }
653 }
654 }
655
656 //
657 // -- Perform solve
658 //
659 // Solve.
660 linSys->Solve(forCoeffs, inout, m_assemblyMap);
661
662 //
663 // -- Postprocess
664 //
665
666 // Scatter back to field degrees of freedom
667 for (nv = 0; nv < nVel; ++nv)
668 {
669 for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
670 {
671 int nBnd = m_bmap[i].size();
672 int nInt = m_imap[i].size();
673 int offset = m_fields[nv]->GetCoeff_Offset(i);
674
675 for (j = 0; j < nBnd; ++j)
676 {
677 m_fields[nv]->UpdateCoeffs()[offset + m_bmap[i][j]] =
678 inout[nVel * offset + nv * nBnd + j];
679 }
680 for (j = 0; j < nInt; ++j)
681 {
682 m_fields[nv]->UpdateCoeffs()[offset + m_imap[i][j]] =
683 inout[nVel * (offset + nBnd) + nv * nInt + j];
684 }
685 }
686 m_fields[nv]->BwdTrans(m_fields[nv]->GetCoeffs(),
687 m_fields[nv]->UpdatePhys());
688 }
689}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void BuildMatrixSystem()
Build matrix system for linear elasticity equations.
Array< OneD, Array< OneD, NekDouble > > m_temperature
Storage for the temperature terms.
CoupledAssemblyMapSharedPtr m_assemblyMap
Assembly map for the coupled (u,v,w) system.
NekDouble m_beta
Parameter dictating amount of thermal stress to add.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_stress
Storage for the thermal stress terms.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
static void Dgeev(const char &uplo, const char &lrev, const int &n, const double *a, const int &lda, double *wr, double *wi, double *rev, const int &ldr, double *lev, const int &ldv, double *work, const int &lwork, int &info)
Solve general real matrix eigenproblem.
Definition: Lapack.hpp:329
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
static PreconditionerSharedPtr NullPreconditionerSharedPtr
std::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:51
@ eDeformed
Geometry is curved or has non-constant factors.
NekMatrix< NekDouble, StandardMatrixTag > DNekMat
Definition: NekTypeDefs.hpp:50
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
DNekMat MappingIdealToRef(SpatialDomains::GeometrySharedPtr geom)
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 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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::beta, BuildMatrixSystem(), Lapack::Dgeev(), Nektar::SpatialDomains::eDeformed, Nektar::eDIAGONAL, Nektar::MultiRegions::eDirectStaticCond, Nektar::eFULL, Nektar::MultiRegions::eIterativeStaticCond, Nektar::StdRegions::eLinearAdvectionReaction, Nektar::MultiRegions::ePETScStaticCond, Nektar::MultiRegions::eXxtStaticCond, Vmath::Fill(), Nektar::SolverUtils::EquationSystem::GetFunction(), Nektar::SolverUtils::EquationSystem::GetNpoints(), m_assemblyMap, m_beta, m_BinvD, m_bmap, m_C, m_Dinv, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_graph, m_imap, m_schurCompl, Nektar::SolverUtils::EquationSystem::m_session, m_stress, m_temperature, Nektar::MappingIdealToRef(), Nektar::MultiRegions::NullPreconditionerSharedPtr, Vmath::Smul(), tinysimd::sqrt(), Nektar::Transpose(), Vmath::Vadd(), and Vmath::Vcopy().

Referenced by Nektar::IterativeElasticSystem::v_DoSolve().

◆ v_ExtraFldOutput()

void Nektar::LinearElasticSystem::v_ExtraFldOutput ( std::vector< Array< OneD, NekDouble > > &  fieldcoeffs,
std::vector< std::string > &  variables 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 816 of file LinearElasticSystem.cpp.

819{
820 const int nVel = m_fields[0]->GetCoordim(0);
821 const int nCoeffs = m_fields[0]->GetNcoeffs();
822
823 if (m_temperature.size() == 0)
824 {
825 return;
826 }
827
828 for (int i = 0; i < nVel; ++i)
829 {
830 Array<OneD, NekDouble> tFwd(nCoeffs);
831 m_fields[i]->FwdTrans(m_temperature[i], tFwd);
832 fieldcoeffs.push_back(tFwd);
833 variables.push_back("ThermStressDiv" + std::to_string(i));
834 }
835
836 if (m_stress.size() == 0)
837 {
838 return;
839 }
840
841 for (int i = 0; i < nVel; ++i)
842 {
843 for (int j = 0; j < nVel; ++j)
844 {
845 Array<OneD, NekDouble> tFwd(nCoeffs);
846 m_fields[i]->FwdTrans(m_stress[i][j], tFwd);
847 fieldcoeffs.push_back(tFwd);
848 variables.push_back("ThermStress" + std::to_string(i) +
849 std::to_string(j));
850 }
851 }
852}

References Nektar::SolverUtils::EquationSystem::m_fields, m_stress, and m_temperature.

◆ v_GenerateSummary()

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

Generate summary at runtime.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 359 of file LinearElasticSystem.cpp.

360{
362
363 AddSummaryItem(s, "Young's modulus", m_E);
364 AddSummaryItem(s, "Poisson ratio", m_nu);
365}
SOLVER_UTILS_EXPORT void SessionSummary(SummaryList &vSummary)
Write out a session summary.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47

References Nektar::SolverUtils::AddSummaryItem(), m_E, m_nu, and Nektar::SolverUtils::EquationSystem::SessionSummary().

Referenced by Nektar::IterativeElasticSystem::v_GenerateSummary().

◆ v_InitObject()

void Nektar::LinearElasticSystem::v_InitObject ( bool  DeclareFields = true)
overrideprotectedvirtual

Set up the linear elasticity system.

This routine loads the E and nu variables from the session file, creates a coupled assembly map and creates containers for the statically condensed block matrix system.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 150 of file LinearElasticSystem.cpp.

151{
152 EquationSystem::v_InitObject(DeclareFields);
153
154 const int nVel = m_fields[0]->GetCoordim(0);
155 int n;
156
157 ASSERTL0(nVel > 1, "Linear elastic solver not set up for"
158 " this dimension (only 2D/3D supported).");
159
160 // Make sure that we have Young's modulus and Poisson ratio set.
161 m_session->LoadParameter("E", m_E, 1.00);
162 m_session->LoadParameter("nu", m_nu, 0.25);
163 m_session->LoadParameter("Beta", m_beta, 1.00);
164
165 // Create a coupled assembly map which allows us to tie u, v and w
166 // fields together.
168 std::dynamic_pointer_cast<MultiRegions::ContField>(m_fields[0]);
170 m_session, m_graph, u->GetLocalToGlobalMap(),
171 m_fields[0]->GetBndConditions(), m_fields);
172
173 // Figure out size of our new matrix systems by looping over all expansions
174 // and multiply number of coefficients by velocity components.
175 const int nEl = m_fields[0]->GetExpSize();
177
178 Array<OneD, unsigned int> sizeBnd(nEl);
179 Array<OneD, unsigned int> sizeInt(nEl);
180
181 // Allocate storage for boundary and interior map caches.
182 m_bmap = Array<OneD, Array<OneD, unsigned int>>(nEl);
183 m_imap = Array<OneD, Array<OneD, unsigned int>>(nEl);
184
185 // Cache interior and boundary maps to avoid recalculating
186 // constantly. Should really be handled by manager in StdRegions but not
187 // really working yet.
188 for (n = 0; n < nEl; ++n)
189 {
190 exp = m_fields[0]->GetExp(n);
191 sizeBnd[n] = nVel * exp->NumBndryCoeffs();
192 sizeInt[n] = nVel * exp->GetNcoeffs() - sizeBnd[n];
193 exp->GetBoundaryMap(m_bmap[n]);
194 exp->GetInteriorMap(m_imap[n]);
195 }
196
197 // Create block matrix storage for the statically condensed system.
198 MatrixStorage blkmatStorage = eDIAGONAL;
200 sizeBnd, sizeBnd, blkmatStorage);
202 blkmatStorage);
204 blkmatStorage);
206 blkmatStorage);
207}
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:268

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::eDIAGONAL, m_assemblyMap, m_beta, m_BinvD, m_bmap, m_C, m_Dinv, m_E, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_graph, m_imap, m_nu, m_schurCompl, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::SolverUtils::EquationSystem::v_InitObject().

Referenced by Nektar::IterativeElasticSystem::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< LinearElasticSystem >

friend class MemoryManager< LinearElasticSystem >
friend

Class may only be instantiated through the MemoryManager.

Definition at line 1 of file LinearElasticSystem.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 68 of file LinearElasticSystem.h.

◆ m_assemblyMap

CoupledAssemblyMapSharedPtr Nektar::LinearElasticSystem::m_assemblyMap
protected

Assembly map for the coupled (u,v,w) system.

Definition at line 91 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_beta

NekDouble Nektar::LinearElasticSystem::m_beta
protected

Parameter dictating amount of thermal stress to add.

Definition at line 89 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_BinvD

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_BinvD
protected

Matrix of elemental \( B^{-1}D \) components.

Definition at line 95 of file LinearElasticSystem.h.

Referenced by SetStaticCondBlock(), v_DoSolve(), and v_InitObject().

◆ m_bmap

Array<OneD, Array<OneD, unsigned int> > Nektar::LinearElasticSystem::m_bmap
protected

Boundary maps for each of the fields.

Definition at line 101 of file LinearElasticSystem.h.

Referenced by SetStaticCondBlock(), v_DoSolve(), and v_InitObject().

◆ m_C

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_C
protected

Matrix of elemental \( C \) components.

Definition at line 97 of file LinearElasticSystem.h.

Referenced by SetStaticCondBlock(), v_DoSolve(), and v_InitObject().

◆ m_Dinv

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_Dinv
protected

Matrix of elemental \( D^{-1} \) components.

Definition at line 99 of file LinearElasticSystem.h.

Referenced by SetStaticCondBlock(), v_DoSolve(), and v_InitObject().

◆ m_E

NekDouble Nektar::LinearElasticSystem::m_E
protected

Young's modulus.

Definition at line 87 of file LinearElasticSystem.h.

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

◆ m_imap

Array<OneD, Array<OneD, unsigned int> > Nektar::LinearElasticSystem::m_imap
protected

Interior maps for each of the fields.

Definition at line 103 of file LinearElasticSystem.h.

Referenced by SetStaticCondBlock(), v_DoSolve(), and v_InitObject().

◆ m_nu

NekDouble Nektar::LinearElasticSystem::m_nu
protected

Poisson ratio.

Definition at line 85 of file LinearElasticSystem.h.

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

◆ m_schurCompl

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_schurCompl
protected

Schur complement boundary-boundary matrix.

Definition at line 93 of file LinearElasticSystem.h.

Referenced by SetStaticCondBlock(), v_DoSolve(), and v_InitObject().

◆ m_stress

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::LinearElasticSystem::m_stress
protected

Storage for the thermal stress terms.

Definition at line 107 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_ExtraFldOutput().

◆ m_temperature

Array<OneD, Array<OneD, NekDouble> > Nektar::LinearElasticSystem::m_temperature
protected

Storage for the temperature terms.

Definition at line 105 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_ExtraFldOutput().