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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::LinearElasticSystem:
Collaboration graph
[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 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 >
boost::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 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 void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. 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 WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $. More...
 
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form. 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 ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to 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 GetNumElmVelocity ()
 
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 SetStepsToOne ()
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 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)
 Default constructor. More...
 
virtual void v_InitObject ()
 Set up the linear elasticity system. More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Generate summary at runtime. More...
 
virtual void v_DoSolve ()
 Solve elliptic linear elastic system. More...
 
virtual void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const string &s1, const string &s2)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Virtual function for initialisation 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_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)
 
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
 
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 

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...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
map< std::string, Array< OneD, Array< OneD, float > > > m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD, Array< OneD, unsigned int > > > m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_base
 Base fields. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_derivedfields
 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...
 
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_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, Array< OneD, Array< OneD, NekDouble > > > m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_tanbasis
 2 x m_spacedim x nq 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...
 
int m_NumMode
 Mode to use in case of single mode analysis. 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...
 

Detailed Description

Base class for linear elastic system.

Definition at line 50 of file LinearElasticSystem.h.

Constructor & Destructor Documentation

Nektar::LinearElasticSystem::LinearElasticSystem ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Default constructor.

Definition at line 119 of file LinearElasticSystem.cpp.

121  : EquationSystem(pSession)
122 {
123 }
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises EquationSystem class members.

Member Function Documentation

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 842 of file LinearElasticSystem.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::eFULL, and Vmath::Smul().

Referenced by BuildMatrixSystem().

847 {
848  const int nCoeffs = exp->GetNcoeffs();
849  const int nPhys = exp->GetTotPoints();
850  int i;
851 
853  nCoeffs, nCoeffs, 0.0, eFULL);
854 
855  Array<OneD, NekDouble> tmp2(nPhys);
856  Array<OneD, NekDouble> tmp3(nPhys);
857 
858  for (i = 0; i < nCoeffs; ++i)
859  {
860  Array<OneD, NekDouble> tmp1(nCoeffs, 0.0);
861  tmp1[i] = 1.0;
862 
863  exp->BwdTrans ( tmp1, tmp2);
864  exp->PhysDeriv (k1, tmp2, tmp3);
865  exp->IProductWRTDerivBase(k2, tmp3, tmp1);
866 
867  Vmath::Smul(
868  nCoeffs, scale, &tmp1[0], 1, &(ret->GetPtr())[0]+i*nCoeffs, 1);
869  }
870 
871  return ret;
872 }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
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:199
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.

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

Referenced by v_DoSolve().

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(m_fields[0]->GetOffset_Elmt_Id(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 
265  exp->DetShapeType(),
266  *exp, factors, varcoeffA);
268  exp->DetShapeType(),
269  *exp, factors, varcoeffD);
270 
271  /*
272  * mat holds the linear operator [ A B ] acting on [ u ].
273  * [ C D ] [ v ]
274  */
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: "
287  << (int)(100.0 * n / nEl) << "%" << flush;
288  }
289  }
290  }
291  else if (nVel == 3)
292  {
293  for (n = 0; n < nEl; ++n)
294  {
295  exp = m_fields[0]->GetExp(m_fields[0]->GetOffset_Elmt_Id(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 
312  exp->DetShapeType(),
313  *exp, factors, varcoeffA);
315  exp->DetShapeType(),
316  *exp, factors, varcoeffE);
318  exp->DetShapeType(),
319  *exp, factors, 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  */
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: "
345  << (int)(100.0 * n / nEl) << "%" << flush;
346  }
347  }
348  }
349 
350  if (verbose && root)
351  {
352  cout << "\rBuilding matrix system: done." << endl;
353  }
354 }
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_nu
Poisson ratio.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:248
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:225
double NekDouble
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.
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
static EquationSystemSharedPtr Nektar::LinearElasticSystem::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 57 of file LinearElasticSystem.h.

59  {
61  LinearElasticSystem>::AllocateSharedPtr(pSession);
62  p->InitObject();
63  return p;
64  }
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
Base class for linear elastic system.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
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 754 of file LinearElasticSystem.cpp.

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

Referenced by BuildMatrixSystem().

758 {
759  int i, j, k, l;
760  const int nVel = mat.GetRows();
761  const int nB = exp->NumBndryCoeffs();
762  const int nI = exp->GetNcoeffs() - nB;
763  const int nBnd = exp->NumBndryCoeffs() * nVel;
764  const int nInt = exp->GetNcoeffs() * nVel - nBnd;
765  const MatrixStorage s = eFULL; // Maybe look into doing symmetric
766  // version of this?
767 
768  DNekMatSharedPtr A =
769  MemoryManager<DNekMat>::AllocateSharedPtr(nBnd, nBnd, 0.0, s);
770  DNekMatSharedPtr B =
771  MemoryManager<DNekMat>::AllocateSharedPtr(nBnd, nInt, 0.0, s);
772  DNekMatSharedPtr C =
773  MemoryManager<DNekMat>::AllocateSharedPtr(nInt, nBnd, 0.0, s);
774  DNekMatSharedPtr D =
775  MemoryManager<DNekMat>::AllocateSharedPtr(nInt, nInt, 0.0, s);
776 
777  for (i = 0; i < nVel; ++i)
778  {
779  for (j = 0; j < nVel; ++j)
780  {
781  // Boundary-boundary and boundary-interior
782  for (k = 0; k < nB; ++k)
783  {
784  for (l = 0; l < nB; ++l)
785  {
786  (*A)(k + i*nB, l + j*nB) =
787  (*mat[i][j])(m_bmap[n][k], m_bmap[n][l]);
788  }
789 
790  for (l = 0; l < nI; ++l)
791  {
792  (*B)(k + i*nB, l + j*nI) =
793  (*mat[i][j])(m_bmap[n][k], m_imap[n][l]);
794  }
795  }
796 
797  // Interior-boundary / interior-interior
798  for (k = 0; k < nI; ++k)
799  {
800  for (l = 0; l < nB; ++l)
801  {
802  (*C)(k + i*nI, l + j*nB) =
803  (*mat[i][j])(m_imap[n][k], m_bmap[n][l]);
804  }
805 
806  for (l = 0; l < nI; ++l)
807  {
808  (*D)(k + i*nI, l + j*nI) =
809  (*mat[i][j])(m_imap[n][k], m_imap[n][l]);
810  }
811  }
812  }
813  }
814 
815  // Construct static condensation matrices.
816  D->Invert();
817  (*B) = (*B)*(*D);
818  (*A) = (*A) - (*B)*(*C);
819 
820  DNekScalMatSharedPtr tmp_mat;
821  m_schurCompl->SetBlock(
823  1.0, A));
824  m_BinvD ->SetBlock(
826  1.0, B));
827  m_C ->SetBlock(
829  1.0, C));
830  m_Dinv ->SetBlock(
832  1.0, D));
833 }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
DNekScalBlkMatSharedPtr m_BinvD
Matrix of elemental components.
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
boost::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
DNekScalBlkMatSharedPtr m_Dinv
Matrix of elemental components.
Array< OneD, Array< OneD, unsigned int > > m_imap
Interior maps for each of the fields.
DNekScalBlkMatSharedPtr m_schurCompl
Schur complement boundary-boundary matrix.
DNekScalBlkMatSharedPtr m_C
Matrix of elemental components.
Array< OneD, Array< OneD, unsigned int > > m_bmap
Boundary maps for each of the fields.
void Nektar::LinearElasticSystem::v_DoSolve ( void  )
protectedvirtual

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.

Reimplemented in Nektar::IterativeElasticSystem.

Definition at line 381 of file LinearElasticSystem.cpp.

References ASSERTL0, BuildMatrixSystem(), Nektar::SpatialDomains::eDeformed, Nektar::eDIAGONAL, Nektar::MultiRegions::eDirectStaticCond, Nektar::eFULL, Nektar::MultiRegions::eIterativeStaticCond, Nektar::StdRegions::eLinearAdvectionReaction, Nektar::SolverUtils::EquationSystem::EvaluateFunction(), Vmath::Fill(), Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::iterator, 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(), Vmath::Neg(), Nektar::MultiRegions::NullPreconditionerSharedPtr, sign, Vmath::Smul(), Nektar::Transpose(), Vmath::Vadd(), and Vmath::Vcopy().

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

382 {
383 
384  int i, j, k, l, nv;
385  const int nVel = m_fields[0]->GetCoordim(0);
386 
387  // Build initial matrix system.
389 
390  // Now we've got the matrix system set up, create a GlobalLinSys
391  // object. We mask ourselves as LinearAdvectionReaction to create a full
392  // matrix instead of symmetric storage.
396 
397  // Currently either direct or iterative static condensation is
398  // supported.
399  if (m_assemblyMap->GetGlobalSysSolnType() == MultiRegions::eDirectStaticCond)
400  {
401  linSys = MemoryManager<
402  MultiRegions::GlobalLinSysDirectStaticCond>::AllocateSharedPtr(
403  key, m_fields[0], m_schurCompl, m_BinvD, m_C, m_Dinv,
404  m_assemblyMap);
405  }
406  else if (m_assemblyMap->GetGlobalSysSolnType() ==
408  {
409  linSys = MemoryManager<
411  key, m_fields[0], m_schurCompl, m_BinvD, m_C, m_Dinv,
413  }
414 
415  linSys->Initialise(m_assemblyMap);
416 
417  const int nCoeffs = m_fields[0]->GetNcoeffs();
418  const int nGlobDofs = m_assemblyMap->GetNumGlobalCoeffs();
419 
420  //
421  // -- Evaluate forcing functions
422  //
423 
424  // Evaluate the forcing function from the XML file.
425  Array<OneD, Array<OneD, NekDouble> > forcing(nVel);
426  EvaluateFunction(forcing, "Forcing");
427 
428  // Add temperature term
429  string tempEval;
430  m_session->LoadSolverInfo("Temperature", tempEval, "None");
431 
432  if (tempEval == "Jacobian")
433  {
434  // Allocate storage
436 
437  for (nv = 0; nv < nVel; ++nv)
438  {
441  m_fields[nv]->GetNpoints());
442 
443  for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
444  {
445  // Calculate element area
447  m_fields[0]->GetExp(i);
449  exp->GetPointsKeys();
451  exp->GetMetricInfo()->GetJac(pkey);
452 
453  int offset = m_fields[0]->GetPhys_Offset(i);
454 
455  if (exp->GetMetricInfo()->GetGtype() ==
457  {
458  Vmath::Smul(exp->GetTotPoints(), m_beta, jac, 1,
459  tmp = m_temperature[nv] + offset, 1);
460  }
461  else
462  {
463  Vmath::Fill(exp->GetTotPoints(), m_beta*jac[0],
464  tmp = m_temperature[nv] + offset, 1);
465  }
466  }
467  m_fields[nv]->PhysDeriv(nv, m_temperature[nv], forcing[nv]);
468  }
469  }
470  else if (tempEval == "Metric")
471  {
472  ASSERTL0((m_fields[0]->GetCoordim(0) == 2 &&
473  m_graph->GetAllQuadGeoms().size() == 0) ||
474  (m_fields[0]->GetCoordim(0) == 3 &&
475  m_graph->GetAllPrismGeoms().size() == 0 &&
476  m_graph->GetAllPyrGeoms ().size() == 0 &&
477  m_graph->GetAllHexGeoms ().size() == 0),
478  "LinearIdealMetric temperature only implemented for "
479  "two-dimensional triangular meshes or three-dimensional "
480  "tetrahedral meshes.");
481 
484 
485  for (nv = 0; nv < nVel; ++nv)
486  {
488  m_fields[nv]->GetNpoints(), 0.0);
489  m_stress[nv] = Array<OneD, Array<OneD, NekDouble> >(nVel);
490 
491  for (i = 0; i < nVel; ++i)
492  {
493  m_stress[nv][i] = Array<OneD, NekDouble>(
494  m_fields[nv]->GetNpoints(), 0.0);
495  }
496  }
497 
498  for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
499  {
500  // Calculate element area
502  m_fields[0]->GetExp(i);
503  LibUtilities::PointsKeyVector pkey = exp->GetPointsKeys();
505  exp->GetMetricInfo()->GetDeriv(pkey);
506  int offset = m_fields[0]->GetPhys_Offset(i);
507 
508  DNekMat i2rm = MappingIdealToRef(exp->GetGeom());
509 
510  // Compute metric tensor
511  DNekMat jac (nVel, nVel, 0.0, eFULL);
512  DNekMat jacIdeal(nVel, nVel, 0.0, eFULL);
513  DNekMat metric (nVel, nVel, 0.0, eFULL);
514 
515  for (j = 0; j < deriv[0][0].num_elements(); ++j)
516  {
517  for (k = 0; k < nVel; ++k)
518  {
519  for (l = 0; l < nVel; ++l)
520  {
521  jac(l,k) = deriv[k][l][j];
522  }
523  }
524 
525  jacIdeal = jac * i2rm;
526  metric = Transpose(jacIdeal) * jacIdeal;
527 
528  // Compute eigenvalues/eigenvectors of metric tensor using
529  // ideal mapping.
530  char jobvl = 'N', jobvr = 'V';
531  int worklen = 8*nVel, info;
532 
533  DNekMat eval (nVel, nVel, 0.0, eDIAGONAL);
534  DNekMat evec (nVel, nVel, 0.0, eFULL);
535  DNekMat evecinv(nVel, nVel, 0.0, eFULL);
536  Array<OneD, NekDouble> vl (nVel*nVel);
537  Array<OneD, NekDouble> work(worklen);
538  Array<OneD, NekDouble> wi (nVel);
539 
540  Lapack::Dgeev(jobvl, jobvr, nVel, metric.GetRawPtr(), nVel,
541  &(eval.GetPtr())[0], &wi[0], &vl[0], nVel,
542  &(evec.GetPtr())[0], nVel,
543  &work[0], worklen, info);
544 
545  evecinv = evec;
546  evecinv.Invert();
547 
548  // rescaling of the eigenvalues
549  for (nv = 0; nv < nVel; ++nv)
550  {
551  eval(nv,nv) = m_beta * (sqrt(eval(nv,nv)) - 1.0);
552  }
553 
554  DNekMat beta = evec * eval * evecinv;
555 
556  for (k = 0; k < nVel; ++k)
557  {
558  for (l = 0; l < nVel; ++l)
559  {
560  m_stress[k][l][offset+j] = beta(k,l);
561  }
562  }
563  }
564 
565  if (deriv[0][0].num_elements() != exp->GetTotPoints())
566  {
568  for (k = 0; k < nVel; ++k)
569  {
570  for (l = 0; l < nVel; ++l)
571  {
572  Vmath::Fill(
573  exp->GetTotPoints(), m_stress[k][l][offset],
574  tmp = m_stress[k][l] + offset, 1);
575  }
576  }
577  }
578  }
579 
580  // Calculate divergence of stress tensor.
582  for (i = 0; i < nVel; ++i)
583  {
584  for (j = 0; j < nVel; ++j)
585  {
586  m_fields[i]->PhysDeriv(j, m_stress[i][j], tmpderiv);
587  Vmath::Vadd (m_fields[i]->GetNpoints(), tmpderiv, 1,
588  m_temperature[i], 1, m_temperature[i], 1);
589  }
590 
592  m_temperature[i], 1, forcing[i], 1);
593  }
594  }
595  else if (tempEval != "None")
596  {
597  ASSERTL0(false, "Unknown temperature form: " + tempEval);
598  }
599 
600  // Set up some temporary storage.
601  //
602  // - forCoeffs holds the forcing coefficients in a local ordering;
603  // however note that the ordering is different and dictated by the
604  // assembly map. We loop over each element, then the boundary degrees of
605  // freedom for u, boundary for v, followed by the interior for u and then
606  // interior for v.
607  // - rhs is the global assembly of forCoeffs.
608  // - inout holds the Dirichlet degrees of freedom in the global ordering,
609  // which have been assembled from the boundary expansion.
610  Array<OneD, NekDouble> forCoeffs(nVel * nCoeffs, 0.0);
611  Array<OneD, NekDouble> inout (nGlobDofs, 0.0);
612  Array<OneD, NekDouble> rhs (nGlobDofs, 0.0);
613 
614  for (nv = 0; nv < nVel; ++nv)
615  {
616  // Take the inner product of the forcing function.
617  Array<OneD, NekDouble> tmp(nCoeffs);
618  m_fields[nv]->IProductWRTBase_IterPerExp(forcing[nv], tmp);
619 
620  // Scatter forcing into RHS vector according to the ordering dictated in
621  // the comment above.
622  for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
623  {
624  int nBnd = m_bmap[i].num_elements();
625  int nInt = m_imap[i].num_elements();
626  int offset = m_fields[nv]->GetCoeff_Offset(i);
627 
628  for (j = 0; j < nBnd; ++j)
629  {
630  forCoeffs[nVel*offset + nv*nBnd + j] =
631  tmp[offset+m_bmap[i][j]];
632  }
633  for (j = 0; j < nInt; ++j)
634  {
635  forCoeffs[nVel*(offset + nBnd) + nv*nInt + j] =
636  tmp[offset+m_imap[i][j]];
637  }
638  }
639  }
640 
641  // -- Impose Dirichlet boundary conditions.
642 
643  // First try to do parallel assembly: the intention here is that Dirichlet
644  // values at some edges/vertices need to be communicated to processes which
645  // don't contain the entire boundary region. See
646  // ContField2D::v_ImposeDirichletConditions for more detail.
647  map<int, vector<MultiRegions::ExtraDirDof> > &extraDirDofs =
648  m_assemblyMap->GetExtraDirDofs();
649  map<int, vector<MultiRegions::ExtraDirDof> >::iterator it;
650 
651  for (nv = 0; nv < nVel; ++nv)
652  {
654  = m_fields[nv]->GetBndCondExpansions();
655 
656  // First try to do parallel stuff
657  for (it = extraDirDofs.begin(); it != extraDirDofs.end(); ++it)
658  {
659  for (i = 0; i < it->second.size(); ++i)
660  {
661  inout[it->second.at(i).get<1>()*nVel + nv] =
662  bndCondExp[it->first]->GetCoeffs()[
663  it->second.at(i).get<0>()]*it->second.at(i).get<2>();
664  }
665  }
666  }
667 
668  m_assemblyMap->UniversalAssemble(inout);
669 
670  // Counter for the local Dirichlet boundary to global ordering.
671  int bndcnt = 0;
672 
673  // Now assemble local boundary contributions.
674  for (nv = 0; nv < nVel; ++nv)
675  {
677  = m_fields[nv]->GetBndCondExpansions();
678  const Array<OneD, const int> &bndMap
679  = m_assemblyMap->GetBndCondCoeffsToGlobalCoeffsMap();
680 
681  for (i = 0; i < bndCondExp.num_elements(); ++i)
682  {
683  const Array<OneD,const NekDouble> &bndCoeffs =
684  bndCondExp[i]->GetCoeffs();
685 
686  for (j = 0; j < bndCondExp[i]->GetNcoeffs(); ++j)
687  {
688  NekDouble sign =
689  m_assemblyMap->GetBndCondCoeffsToGlobalCoeffsSign(
690  bndcnt);
691  inout[bndMap[bndcnt++]] = sign * bndCoeffs[j];
692  }
693  }
694  }
695 
696  //
697  // -- Perform solve
698  //
699 
700  // Assemble forcing into the RHS.
701  m_assemblyMap->Assemble(forCoeffs, rhs);
702 
703  // Negate RHS to be consistent with matrix definition.
704  Vmath::Neg(rhs.num_elements(), rhs, 1);
705 
706  // Solve.
707  linSys->Solve(rhs, inout, m_assemblyMap);
708 
709  // Scatter the global ordering back to the alternate local ordering.
710  Array<OneD, NekDouble> tmp(nVel * nCoeffs);
711  m_assemblyMap->GlobalToLocal(inout, tmp);
712 
713  //
714  // -- Postprocess
715  //
716 
717  // Scatter back to field degrees of freedom
718  for (nv = 0; nv < nVel; ++nv)
719  {
720  for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
721  {
722  int nBnd = m_bmap[i].num_elements();
723  int nInt = m_imap[i].num_elements();
724  int offset = m_fields[nv]->GetCoeff_Offset(i);
725 
726  for (j = 0; j < nBnd; ++j)
727  {
728  m_fields[nv]->UpdateCoeffs()[offset+m_bmap[i][j]] =
729  tmp[nVel*offset + nv*nBnd + j];
730  }
731  for (j = 0; j < nInt; ++j)
732  {
733  m_fields[nv]->UpdateCoeffs()[offset+m_imap[i][j]] =
734  tmp[nVel*(offset + nBnd) + nv*nInt + j];
735  }
736  }
737  m_fields[nv]->BwdTrans(m_fields[nv]->GetCoeffs(),
738  m_fields[nv]->UpdatePhys());
739  }
740 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:220
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
NekDouble m_beta
Parameter dictating amount of thermal stress to add.
DNekScalBlkMatSharedPtr m_BinvD
Matrix of elemental components.
Array< OneD, Array< OneD, NekDouble > > m_temperature
Storage for the temperature terms.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_stress
Storage for the thermal stress terms.
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:199
static PreconditionerSharedPtr NullPreconditionerSharedPtr
DNekScalBlkMatSharedPtr m_Dinv
Matrix of elemental components.
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
double NekDouble
Array< OneD, Array< OneD, unsigned int > > m_imap
Interior maps for each of the fields.
SOLVER_UTILS_EXPORT void EvaluateFunction(Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
Evaluates a function as specified in the session file.
Describe a linear system.
DNekScalBlkMatSharedPtr m_schurCompl
Schur complement boundary-boundary matrix.
DNekScalBlkMatSharedPtr m_C
Matrix of elemental components.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
CoupledAssemblyMapSharedPtr m_assemblyMap
Assembly map for the coupled (u,v,w) system.
Array< OneD, Array< OneD, unsigned int > > m_bmap
Boundary maps for each of the fields.
boost::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:51
void BuildMatrixSystem()
Build matrix system for linear elasticity equations.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
DNekMat MappingIdealToRef(SpatialDomains::GeometrySharedPtr geom)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
Geometry is curved or has non-constant factors.
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:285
void Nektar::LinearElasticSystem::v_ExtraFldOutput ( std::vector< Array< OneD, NekDouble > > &  fieldcoeffs,
std::vector< std::string > &  variables 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 874 of file LinearElasticSystem.cpp.

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

877 {
878  const int nVel = m_fields[0]->GetCoordim(0);
879  const int nCoeffs = m_fields[0]->GetNcoeffs();
880 
881  if (m_temperature.num_elements() == 0)
882  {
883  return;
884  }
885 
886  for (int i = 0; i < nVel; ++i)
887  {
888  Array<OneD, NekDouble> tFwd(nCoeffs);
889  m_fields[i]->FwdTrans(m_temperature[i], tFwd);
890  fieldcoeffs.push_back(tFwd);
891  variables.push_back(
892  "ThermStressDiv" + boost::lexical_cast<std::string>(i));
893  }
894 
895  if (m_stress.num_elements() == 0)
896  {
897  return;
898  }
899 
900  for (int i = 0; i < nVel; ++i)
901  {
902  for (int j = 0; j < nVel; ++j)
903  {
904  Array<OneD, NekDouble> tFwd(nCoeffs);
905  m_fields[i]->FwdTrans(m_stress[i][j], tFwd);
906  fieldcoeffs.push_back(tFwd);
907  variables.push_back(
908  "ThermStress"
909  + boost::lexical_cast<std::string>(i)
910  + boost::lexical_cast<std::string>(j));
911  }
912  }
913 }
Array< OneD, Array< OneD, NekDouble > > m_temperature
Storage for the temperature terms.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_stress
Storage for the thermal stress terms.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::LinearElasticSystem::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Generate summary at runtime.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::IterativeElasticSystem.

Definition at line 359 of file LinearElasticSystem.cpp.

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

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

360 {
361  EquationSystem::SessionSummary(s);
362 
363  AddSummaryItem(s, "Young's modulus", m_E);
364  AddSummaryItem(s, "Poisson ratio", m_nu);
365 }
NekDouble m_nu
Poisson ratio.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
NekDouble m_E
Young's modulus.
void Nektar::LinearElasticSystem::v_InitObject ( )
protectedvirtual

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.

Reimplemented in Nektar::IterativeElasticSystem.

Definition at line 132 of file LinearElasticSystem.cpp.

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

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

133 {
134  EquationSystem::v_InitObject();
135 
136  const int nVel = m_fields[0]->GetCoordim(0);
137  int n;
138 
139  ASSERTL0(nVel > 1, "Linear elastic solver not set up for"
140  " this dimension (only 2D/3D supported).");
141 
142  // Make sure that we have Young's modulus and Poisson ratio set.
143  m_session->LoadParameter("E", m_E, 1.00);
144  m_session->LoadParameter("nu", m_nu, 0.25);
145  m_session->LoadParameter("Beta", m_beta, 1.00);
146 
147  // Create a coupled assembly map which allows us to tie u, v and w
148  // fields together.
149  if (nVel == 2)
150  {
151  MultiRegions::ContField2DSharedPtr u = boost::dynamic_pointer_cast<
155  m_graph,
156  u->GetLocalToGlobalMap(),
158  m_fields);
159  }
160 
161  if (nVel == 3)
162  {
163  MultiRegions::ContField3DSharedPtr u = boost::dynamic_pointer_cast<
165  m_assemblyMap = MemoryManager<CoupledAssemblyMap>
167  m_graph,
168  u->GetLocalToGlobalMap(),
170  m_fields);
171  }
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.
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(m_fields[0]->GetOffset_Elmt_Id(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  sizeBnd, sizeInt, blkmatStorage);
204  sizeInt, sizeBnd, blkmatStorage);
206  sizeInt, sizeInt, blkmatStorage);
207 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
NekDouble m_nu
Poisson ratio.
NekDouble m_beta
Parameter dictating amount of thermal stress to add.
boost::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:291
DNekScalBlkMatSharedPtr m_BinvD
Matrix of elemental components.
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField2D.h:56
DNekScalBlkMatSharedPtr m_Dinv
Matrix of elemental components.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
Array< OneD, Array< OneD, unsigned int > > m_imap
Interior maps for each of the fields.
DNekScalBlkMatSharedPtr m_schurCompl
Schur complement boundary-boundary matrix.
DNekScalBlkMatSharedPtr m_C
Matrix of elemental components.
NekDouble m_E
Young's modulus.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
CoupledAssemblyMapSharedPtr m_assemblyMap
Assembly map for the coupled (u,v,w) system.
Array< OneD, Array< OneD, unsigned int > > m_bmap
Boundary maps for each of the fields.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
boost::shared_ptr< ContField3D > ContField3DSharedPtr
Definition: ContField3D.h:190

Friends And Related Function Documentation

friend class MemoryManager< LinearElasticSystem >
friend

Class may only be instantiated through the MemoryManager.

Definition at line 54 of file LinearElasticSystem.h.

Member Data Documentation

string Nektar::LinearElasticSystem::className
static
Initial value:
RegisterCreatorFunction("LinearElasticSystem",

Name of class.

Definition at line 67 of file LinearElasticSystem.h.

CoupledAssemblyMapSharedPtr Nektar::LinearElasticSystem::m_assemblyMap
protected

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

Definition at line 93 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().

NekDouble Nektar::LinearElasticSystem::m_beta
protected

Parameter dictating amount of thermal stress to add.

Definition at line 91 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_InitObject().

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_BinvD
protected

Matrix of elemental $ B^{-1}D $ components.

Definition at line 97 of file LinearElasticSystem.h.

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

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

Boundary maps for each of the fields.

Definition at line 103 of file LinearElasticSystem.h.

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

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_C
protected

Matrix of elemental $ C $ components.

Definition at line 99 of file LinearElasticSystem.h.

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

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_Dinv
protected

Matrix of elemental $ D^{-1} $ components.

Definition at line 101 of file LinearElasticSystem.h.

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

NekDouble Nektar::LinearElasticSystem::m_E
protected

Young's modulus.

Definition at line 89 of file LinearElasticSystem.h.

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

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

Interior maps for each of the fields.

Definition at line 105 of file LinearElasticSystem.h.

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

NekDouble Nektar::LinearElasticSystem::m_nu
protected

Poisson ratio.

Definition at line 87 of file LinearElasticSystem.h.

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

DNekScalBlkMatSharedPtr Nektar::LinearElasticSystem::m_schurCompl
protected

Schur complement boundary-boundary matrix.

Definition at line 95 of file LinearElasticSystem.h.

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

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

Storage for the thermal stress terms.

Definition at line 109 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_ExtraFldOutput().

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

Storage for the temperature terms.

Definition at line 107 of file LinearElasticSystem.h.

Referenced by v_DoSolve(), and v_ExtraFldOutput().