Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 NekDouble &pTime=0.0, 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 NekDouble &pTime=0.0, 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, NekDouble
ErrorExtraPoints (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::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 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 SetSteps (const int steps)
 
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...
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

static 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)
 
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::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_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...
 
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...
 
std::set< std::string > m_loadedFields
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, 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 127 of file LinearElasticSystem.cpp.

129  : EquationSystem(pSession)
130 {
131 }
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 878 of file LinearElasticSystem.cpp.

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

Referenced by BuildMatrixSystem().

883 {
884  const int nCoeffs = exp->GetNcoeffs();
885  const int nPhys = exp->GetTotPoints();
886  int i;
887 
889  nCoeffs, nCoeffs, 0.0, eFULL);
890 
891  Array<OneD, NekDouble> tmp2(nPhys);
892  Array<OneD, NekDouble> tmp3(nPhys);
893 
894  for (i = 0; i < nCoeffs; ++i)
895  {
896  Array<OneD, NekDouble> tmp1(nCoeffs, 0.0);
897  tmp1[i] = 1.0;
898 
899  exp->BwdTrans ( tmp1, tmp2);
900  exp->PhysDeriv (k1, tmp2, tmp3);
901  exp->IProductWRTDerivBase(k2, tmp3, tmp1);
902 
903  Vmath::Smul(
904  nCoeffs, scale, &tmp1[0], 1, &(ret->GetPtr())[0]+i*nCoeffs, 1);
905  }
906 
907  return ret;
908 }
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 238 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().

239 {
240  const int nEl = m_fields[0]->GetExpSize();
241  const int nVel = m_fields[0]->GetCoordim(0);
242 
244  int n;
245 
246  // Factors map for matrix keys.
248 
249  // Calculate various constants
250  NekDouble mu = m_E * 0.5 / (1.0 + m_nu);
251  NekDouble lambda = m_E * m_nu / (1.0 + m_nu) / (1.0 - 2.0*m_nu);
252 
253  bool verbose = m_session->DefinesCmdLineArgument("verbose");
254  bool root = m_session->GetComm()->GetRank() == 0;
255 
256  // Loop over each element and construct matrices.
257  if (nVel == 2)
258  {
259  for (n = 0; n < nEl; ++n)
260  {
261  exp = m_fields[0]->GetExp(m_fields[0]->GetOffset_Elmt_Id(n));
262  const int nPhys = exp->GetTotPoints();
263  Array<OneD, NekDouble> aArr(nPhys, lambda + 2*mu);
264  Array<OneD, NekDouble> bArr(nPhys, mu);
265 
266  StdRegions::VarCoeffMap varcoeffA, varcoeffD;
267  varcoeffA[StdRegions::eVarCoeffD00] = aArr;
268  varcoeffA[StdRegions::eVarCoeffD11] = bArr;
269  varcoeffD[StdRegions::eVarCoeffD00] = bArr;
270  varcoeffD[StdRegions::eVarCoeffD11] = aArr;
271 
273  exp->DetShapeType(),
274  *exp, factors, varcoeffA);
276  exp->DetShapeType(),
277  *exp, factors, varcoeffD);
278 
279  /*
280  * mat holds the linear operator [ A B ] acting on [ u ].
281  * [ C D ] [ v ]
282  */
284  mat[0][0] = exp->GenMatrix(matkeyA);
285  mat[0][1] = BuildLaplacianIJMatrix(1, 0, mu+lambda, exp);
286  mat[1][0] = mat[0][1];
287  mat[1][1] = exp->GenMatrix(matkeyD);
288 
289  // Set up the statically condensed block for this element.
290  SetStaticCondBlock(n, exp, mat);
291 
292  if (verbose && root)
293  {
294  cout << "\rBuilding matrix system: "
295  << (int)(100.0 * n / nEl) << "%" << flush;
296  }
297  }
298  }
299  else if (nVel == 3)
300  {
301  for (n = 0; n < nEl; ++n)
302  {
303  exp = m_fields[0]->GetExp(m_fields[0]->GetOffset_Elmt_Id(n));
304  const int nPhys = exp->GetTotPoints();
305  Array<OneD, NekDouble> aArr(nPhys, lambda + 2*mu);
306  Array<OneD, NekDouble> bArr(nPhys, mu);
307 
308  StdRegions::VarCoeffMap varcoeffA, varcoeffE, varcoeffI;
309  varcoeffA[StdRegions::eVarCoeffD00] = aArr;
310  varcoeffA[StdRegions::eVarCoeffD11] = bArr;
311  varcoeffA[StdRegions::eVarCoeffD22] = bArr;
312  varcoeffE[StdRegions::eVarCoeffD00] = bArr;
313  varcoeffE[StdRegions::eVarCoeffD11] = aArr;
314  varcoeffE[StdRegions::eVarCoeffD22] = bArr;
315  varcoeffI[StdRegions::eVarCoeffD00] = bArr;
316  varcoeffI[StdRegions::eVarCoeffD11] = bArr;
317  varcoeffI[StdRegions::eVarCoeffD22] = aArr;
318 
320  exp->DetShapeType(),
321  *exp, factors, varcoeffA);
323  exp->DetShapeType(),
324  *exp, factors, varcoeffE);
326  exp->DetShapeType(),
327  *exp, factors, varcoeffI);
328 
329  /*
330  * mat holds the linear operator [ A B C ] acting on [ u ].
331  * [ D E F ] [ v ]
332  * [ G H I ] [ w ]
333  */
335  mat[0][0] = exp->GenMatrix(matkeyA);
336  mat[0][1] = BuildLaplacianIJMatrix(1, 0, mu + lambda, exp);
337  mat[0][2] = BuildLaplacianIJMatrix(2, 0, mu + lambda, exp);
338 
339  mat[1][0] = mat[0][1];
340  mat[1][1] = exp->GenMatrix(matkeyE);
341  mat[1][2] = BuildLaplacianIJMatrix(2, 1, mu + lambda, exp);
342 
343  mat[2][0] = mat[0][2];
344  mat[2][1] = mat[1][2];
345  mat[2][2] = exp->GenMatrix(matkeyI);
346 
347  // Set up the statically condensed block for this element.
348  SetStaticCondBlock(n, exp, mat);
349 
350  if (verbose && root)
351  {
352  cout << "\rBuilding matrix system: "
353  << (int)(100.0 * n / nEl) << "%" << flush;
354  }
355  }
356  }
357 
358  if (verbose && root)
359  {
360  cout << "\rBuilding matrix system: done." << endl;
361  }
362 }
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:251
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
double NekDouble
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
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.
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 790 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().

794 {
795  int i, j, k, l;
796  const int nVel = mat.GetRows();
797  const int nB = exp->NumBndryCoeffs();
798  const int nI = exp->GetNcoeffs() - nB;
799  const int nBnd = exp->NumBndryCoeffs() * nVel;
800  const int nInt = exp->GetNcoeffs() * nVel - nBnd;
801  const MatrixStorage s = eFULL; // Maybe look into doing symmetric
802  // version of this?
803 
804  DNekMatSharedPtr A =
805  MemoryManager<DNekMat>::AllocateSharedPtr(nBnd, nBnd, 0.0, s);
806  DNekMatSharedPtr B =
807  MemoryManager<DNekMat>::AllocateSharedPtr(nBnd, nInt, 0.0, s);
808  DNekMatSharedPtr C =
809  MemoryManager<DNekMat>::AllocateSharedPtr(nInt, nBnd, 0.0, s);
810  DNekMatSharedPtr D =
811  MemoryManager<DNekMat>::AllocateSharedPtr(nInt, nInt, 0.0, s);
812 
813  for (i = 0; i < nVel; ++i)
814  {
815  for (j = 0; j < nVel; ++j)
816  {
817  // Boundary-boundary and boundary-interior
818  for (k = 0; k < nB; ++k)
819  {
820  for (l = 0; l < nB; ++l)
821  {
822  (*A)(k + i*nB, l + j*nB) =
823  (*mat[i][j])(m_bmap[n][k], m_bmap[n][l]);
824  }
825 
826  for (l = 0; l < nI; ++l)
827  {
828  (*B)(k + i*nB, l + j*nI) =
829  (*mat[i][j])(m_bmap[n][k], m_imap[n][l]);
830  }
831  }
832 
833  // Interior-boundary / interior-interior
834  for (k = 0; k < nI; ++k)
835  {
836  for (l = 0; l < nB; ++l)
837  {
838  (*C)(k + i*nI, l + j*nB) =
839  (*mat[i][j])(m_imap[n][k], m_bmap[n][l]);
840  }
841 
842  for (l = 0; l < nI; ++l)
843  {
844  (*D)(k + i*nI, l + j*nI) =
845  (*mat[i][j])(m_imap[n][k], m_imap[n][l]);
846  }
847  }
848  }
849  }
850 
851  // Construct static condensation matrices.
852  D->Invert();
853  (*B) = (*B)*(*D);
854  (*A) = (*A) - (*B)*(*C);
855 
856  DNekScalMatSharedPtr tmp_mat;
857  m_schurCompl->SetBlock(
859  1.0, A));
860  m_BinvD ->SetBlock(
862  1.0, B));
863  m_C ->SetBlock(
865  1.0, C));
866  m_Dinv ->SetBlock(
868  1.0, D));
869 }
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 389 of file LinearElasticSystem.cpp.

References ASSERTL0, BuildMatrixSystem(), Nektar::SpatialDomains::eDeformed, Nektar::eDIAGONAL, Nektar::MultiRegions::eDirectStaticCond, Nektar::SpatialDomains::eDirichlet, Nektar::eFULL, Nektar::MultiRegions::eIterativeStaticCond, Nektar::StdRegions::eLinearAdvectionReaction, Nektar::MultiRegions::ePETScStaticCond, Nektar::SolverUtils::EquationSystem::EvaluateFunction(), Nektar::MultiRegions::eXxtStaticCond, 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().

390 {
391 
392  int i, j, k, l, nv;
393  const int nVel = m_fields[0]->GetCoordim(0);
394 
395  // Build initial matrix system.
397 
398  // Now we've got the matrix system set up, create a GlobalLinSys
399  // object. We mask ourselves as LinearAdvectionReaction to create a full
400  // matrix instead of symmetric storage.
404 
405  // Currently either direct or iterative static condensation is
406  // supported.
407  if (m_assemblyMap->GetGlobalSysSolnType() == MultiRegions::eDirectStaticCond)
408  {
409  linSys = MemoryManager<
410  MultiRegions::GlobalLinSysDirectStaticCond>::AllocateSharedPtr(
411  key, m_fields[0], m_schurCompl, m_BinvD, m_C, m_Dinv,
412  m_assemblyMap);
413  }
414  else if (m_assemblyMap->GetGlobalSysSolnType() ==
416  {
417  linSys = MemoryManager<
419  key, m_fields[0], m_schurCompl, m_BinvD, m_C, m_Dinv,
421  }
422 #ifdef NEKTAR_USE_PETSC
423  else if (m_assemblyMap->GetGlobalSysSolnType() ==
425  {
426  linSys = MemoryManager<
427  MultiRegions::GlobalLinSysPETScStaticCond>::AllocateSharedPtr(
428  key, m_fields[0], m_schurCompl, m_BinvD, m_C, m_Dinv,
429  m_assemblyMap);
430  }
431 #endif
432 #ifdef NEKTAR_USE_MPI
433  else if (m_assemblyMap->GetGlobalSysSolnType() ==
435  {
436  linSys = MemoryManager<
437  MultiRegions::GlobalLinSysXxtStaticCond>::AllocateSharedPtr(
438  key, m_fields[0], m_schurCompl, m_BinvD, m_C, m_Dinv,
439  m_assemblyMap);
440  }
441 #endif
442 
443  linSys->Initialise(m_assemblyMap);
444 
445  const int nCoeffs = m_fields[0]->GetNcoeffs();
446  const int nGlobDofs = m_assemblyMap->GetNumGlobalCoeffs();
447 
448  //
449  // -- Evaluate forcing functions
450  //
451 
452  // Evaluate the forcing function from the XML file.
453  Array<OneD, Array<OneD, NekDouble> > forcing(nVel);
454  EvaluateFunction(forcing, "Forcing");
455 
456  // Add temperature term
457  string tempEval;
458  m_session->LoadSolverInfo("Temperature", tempEval, "None");
459 
460  if (tempEval == "Jacobian")
461  {
462  // Allocate storage
464 
465  for (nv = 0; nv < nVel; ++nv)
466  {
469  m_fields[nv]->GetNpoints());
470 
471  for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
472  {
473  // Calculate element area
475  m_fields[0]->GetExp(i);
477  exp->GetPointsKeys();
479  exp->GetMetricInfo()->GetJac(pkey);
480 
481  int offset = m_fields[0]->GetPhys_Offset(i);
482 
483  if (exp->GetMetricInfo()->GetGtype() ==
485  {
486  Vmath::Smul(exp->GetTotPoints(), m_beta, jac, 1,
487  tmp = m_temperature[nv] + offset, 1);
488  }
489  else
490  {
491  Vmath::Fill(exp->GetTotPoints(), m_beta*jac[0],
492  tmp = m_temperature[nv] + offset, 1);
493  }
494  }
495  m_fields[nv]->PhysDeriv(nv, m_temperature[nv], forcing[nv]);
496  }
497  }
498  else if (tempEval == "Metric")
499  {
500  ASSERTL0((m_fields[0]->GetCoordim(0) == 2 &&
501  m_graph->GetAllQuadGeoms().size() == 0) ||
502  (m_fields[0]->GetCoordim(0) == 3 &&
503  m_graph->GetAllPrismGeoms().size() == 0 &&
504  m_graph->GetAllPyrGeoms ().size() == 0 &&
505  m_graph->GetAllHexGeoms ().size() == 0),
506  "LinearIdealMetric temperature only implemented for "
507  "two-dimensional triangular meshes or three-dimensional "
508  "tetrahedral meshes.");
509 
512 
513  for (nv = 0; nv < nVel; ++nv)
514  {
516  m_fields[nv]->GetNpoints(), 0.0);
517  m_stress[nv] = Array<OneD, Array<OneD, NekDouble> >(nVel);
518 
519  for (i = 0; i < nVel; ++i)
520  {
521  m_stress[nv][i] = Array<OneD, NekDouble>(
522  m_fields[nv]->GetNpoints(), 0.0);
523  }
524  }
525 
526  for (i = 0; i < m_fields[0]->GetExpSize(); ++i)
527  {
528  // Calculate element area
530  m_fields[0]->GetExp(i);
531  LibUtilities::PointsKeyVector pkey = exp->GetPointsKeys();
533  exp->GetMetricInfo()->GetDeriv(pkey);
534  int offset = m_fields[0]->GetPhys_Offset(i);
535 
536  DNekMat i2rm = MappingIdealToRef(exp->GetGeom());
537 
538  // Compute metric tensor
539  DNekMat jac (nVel, nVel, 0.0, eFULL);
540  DNekMat jacIdeal(nVel, nVel, 0.0, eFULL);
541  DNekMat metric (nVel, nVel, 0.0, eFULL);
542 
543  for (j = 0; j < deriv[0][0].num_elements(); ++j)
544  {
545  for (k = 0; k < nVel; ++k)
546  {
547  for (l = 0; l < nVel; ++l)
548  {
549  jac(l,k) = deriv[k][l][j];
550  }
551  }
552 
553  jacIdeal = jac * i2rm;
554  metric = Transpose(jacIdeal) * jacIdeal;
555 
556  // Compute eigenvalues/eigenvectors of metric tensor using
557  // ideal mapping.
558  char jobvl = 'N', jobvr = 'V';
559  int worklen = 8*nVel, info;
560 
561  DNekMat eval (nVel, nVel, 0.0, eDIAGONAL);
562  DNekMat evec (nVel, nVel, 0.0, eFULL);
563  DNekMat evecinv(nVel, nVel, 0.0, eFULL);
564  Array<OneD, NekDouble> vl (nVel*nVel);
565  Array<OneD, NekDouble> work(worklen);
566  Array<OneD, NekDouble> wi (nVel);
567 
568  Lapack::Dgeev(jobvl, jobvr, nVel, metric.GetRawPtr(), nVel,
569  &(eval.GetPtr())[0], &wi[0], &vl[0], nVel,
570  &(evec.GetPtr())[0], nVel,
571  &work[0], worklen, info);
572 
573  evecinv = evec;
574  evecinv.Invert();
575 
576  // rescaling of the eigenvalues
577  for (nv = 0; nv < nVel; ++nv)
578  {
579  eval(nv,nv) = m_beta * (sqrt(eval(nv,nv)) - 1.0);
580  }
581 
582  DNekMat beta = evec * eval * evecinv;
583 
584  for (k = 0; k < nVel; ++k)
585  {
586  for (l = 0; l < nVel; ++l)
587  {
588  m_stress[k][l][offset+j] = beta(k,l);
589  }
590  }
591  }
592 
593  if (deriv[0][0].num_elements() != exp->GetTotPoints())
594  {
596  for (k = 0; k < nVel; ++k)
597  {
598  for (l = 0; l < nVel; ++l)
599  {
600  Vmath::Fill(
601  exp->GetTotPoints(), m_stress[k][l][offset],
602  tmp = m_stress[k][l] + offset, 1);
603  }
604  }
605  }
606  }
607 
608  // Calculate divergence of stress tensor.
610  for (i = 0; i < nVel; ++i)
611  {
612  for (j = 0; j < nVel; ++j)
613  {
614  m_fields[i]->PhysDeriv(j, m_stress[i][j], tmpderiv);
615  Vmath::Vadd (m_fields[i]->GetNpoints(), tmpderiv, 1,
616  m_temperature[i], 1, m_temperature[i], 1);
617  }
618 
620  m_temperature[i], 1, forcing[i], 1);
621  }
622  }
623  else if (tempEval != "None")
624  {
625  ASSERTL0(false, "Unknown temperature form: " + tempEval);
626  }
627 
628  // Set up some temporary storage.
629  //
630  // - forCoeffs holds the forcing coefficients in a local ordering;
631  // however note that the ordering is different and dictated by the
632  // assembly map. We loop over each element, then the boundary degrees of
633  // freedom for u, boundary for v, followed by the interior for u and then
634  // interior for v.
635  // - rhs is the global assembly of forCoeffs.
636  // - inout holds the Dirichlet degrees of freedom in the global ordering,
637  // which have been assembled from the boundary expansion.
638  Array<OneD, NekDouble> forCoeffs(nVel * nCoeffs, 0.0);
639  Array<OneD, NekDouble> inout (nGlobDofs, 0.0);
640  Array<OneD, NekDouble> rhs (nGlobDofs, 0.0);
641 
642  for (nv = 0; nv < nVel; ++nv)
643  {
644  // Take the inner product of the forcing function.
645  Array<OneD, NekDouble> tmp(nCoeffs);
646  m_fields[nv]->IProductWRTBase_IterPerExp(forcing[nv], tmp);
647 
648  // Scatter forcing into RHS vector according to the ordering dictated in
649  // the comment above.
650  for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
651  {
652  int nBnd = m_bmap[i].num_elements();
653  int nInt = m_imap[i].num_elements();
654  int offset = m_fields[nv]->GetCoeff_Offset(i);
655 
656  for (j = 0; j < nBnd; ++j)
657  {
658  forCoeffs[nVel*offset + nv*nBnd + j] =
659  tmp[offset+m_bmap[i][j]];
660  }
661  for (j = 0; j < nInt; ++j)
662  {
663  forCoeffs[nVel*(offset + nBnd) + nv*nInt + j] =
664  tmp[offset+m_imap[i][j]];
665  }
666  }
667  }
668 
669  // -- Impose Dirichlet boundary conditions.
670 
671  // First try to do parallel assembly: the intention here is that Dirichlet
672  // values at some edges/vertices need to be communicated to processes which
673  // don't contain the entire boundary region. See
674  // ContField2D::v_ImposeDirichletConditions for more detail.
675  map<int, vector<MultiRegions::ExtraDirDof> > &extraDirDofs =
676  m_assemblyMap->GetExtraDirDofs();
677  map<int, vector<MultiRegions::ExtraDirDof> >::iterator it;
678 
679  for (nv = 0; nv < nVel; ++nv)
680  {
682  = m_fields[nv]->GetBndCondExpansions();
683 
684  // First try to do parallel stuff
685  for (it = extraDirDofs.begin(); it != extraDirDofs.end(); ++it)
686  {
687  for (i = 0; i < it->second.size(); ++i)
688  {
689  inout[it->second.at(i).get<1>()*nVel + nv] =
690  bndCondExp[it->first]->GetCoeffs()[
691  it->second.at(i).get<0>()]*it->second.at(i).get<2>();
692  }
693  }
694  }
695 
696  m_assemblyMap->UniversalAssemble(inout);
697 
698  // Counter for the local Dirichlet boundary to global ordering.
699  int bndcnt = 0;
700 
701  // Now assemble local boundary contributions.
702  for (nv = 0; nv < nVel; ++nv)
703  {
705  = m_fields[nv]->GetBndCondExpansions();
706  const Array<OneD, const int> &bndMap
707  = m_assemblyMap->GetBndCondCoeffsToGlobalCoeffsMap();
708 
709  for (i = 0; i < bndCondExp.num_elements(); ++i)
710  {
711  if (m_fields[nv]->GetBndConditions()[i]->GetBoundaryConditionType()
713  {
714  const Array<OneD,const NekDouble> &bndCoeffs =
715  bndCondExp[i]->GetCoeffs();
716 
717  for (j = 0; j < bndCondExp[i]->GetNcoeffs(); ++j)
718  {
719  NekDouble sign =
720  m_assemblyMap->GetBndCondCoeffsToGlobalCoeffsSign(
721  bndcnt);
722  inout[bndMap[bndcnt++]] = sign * bndCoeffs[j];
723  }
724  }
725  else
726  {
727  bndcnt += bndCondExp[i]->GetNcoeffs();
728  }
729  }
730  }
731 
732  //
733  // -- Perform solve
734  //
735 
736  // Assemble forcing into the RHS.
737  m_assemblyMap->Assemble(forCoeffs, rhs);
738 
739  // Negate RHS to be consistent with matrix definition.
740  Vmath::Neg(rhs.num_elements(), rhs, 1);
741 
742  // Solve.
743  linSys->Solve(rhs, inout, m_assemblyMap);
744 
745  // Scatter the global ordering back to the alternate local ordering.
746  Array<OneD, NekDouble> tmp(nVel * nCoeffs);
747  m_assemblyMap->GlobalToLocal(inout, tmp);
748 
749  //
750  // -- Postprocess
751  //
752 
753  // Scatter back to field degrees of freedom
754  for (nv = 0; nv < nVel; ++nv)
755  {
756  for (i = 0; i < m_fields[nv]->GetExpSize(); ++i)
757  {
758  int nBnd = m_bmap[i].num_elements();
759  int nInt = m_imap[i].num_elements();
760  int offset = m_fields[nv]->GetCoeff_Offset(i);
761 
762  for (j = 0; j < nBnd; ++j)
763  {
764  m_fields[nv]->UpdateCoeffs()[offset+m_bmap[i][j]] =
765  tmp[nVel*offset + nv*nBnd + j];
766  }
767  for (j = 0; j < nInt; ++j)
768  {
769  m_fields[nv]->UpdateCoeffs()[offset+m_imap[i][j]] =
770  tmp[nVel*(offset + nBnd) + nv*nInt + j];
771  }
772  }
773  m_fields[nv]->BwdTrans(m_fields[nv]->GetCoeffs(),
774  m_fields[nv]->UpdatePhys());
775  }
776 }
#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
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
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.
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:52
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:1047
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 910 of file LinearElasticSystem.cpp.

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

913 {
914  const int nVel = m_fields[0]->GetCoordim(0);
915  const int nCoeffs = m_fields[0]->GetNcoeffs();
916 
917  if (m_temperature.num_elements() == 0)
918  {
919  return;
920  }
921 
922  for (int i = 0; i < nVel; ++i)
923  {
924  Array<OneD, NekDouble> tFwd(nCoeffs);
925  m_fields[i]->FwdTrans(m_temperature[i], tFwd);
926  fieldcoeffs.push_back(tFwd);
927  variables.push_back(
928  "ThermStressDiv" + boost::lexical_cast<std::string>(i));
929  }
930 
931  if (m_stress.num_elements() == 0)
932  {
933  return;
934  }
935 
936  for (int i = 0; i < nVel; ++i)
937  {
938  for (int j = 0; j < nVel; ++j)
939  {
940  Array<OneD, NekDouble> tFwd(nCoeffs);
941  m_fields[i]->FwdTrans(m_stress[i][j], tFwd);
942  fieldcoeffs.push_back(tFwd);
943  variables.push_back(
944  "ThermStress"
945  + boost::lexical_cast<std::string>(i)
946  + boost::lexical_cast<std::string>(j));
947  }
948  }
949 }
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 367 of file LinearElasticSystem.cpp.

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

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

368 {
369  EquationSystem::SessionSummary(s);
370 
371  AddSummaryItem(s, "Young's modulus", m_E);
372  AddSummaryItem(s, "Poisson ratio", m_nu);
373 }
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 140 of file LinearElasticSystem.cpp.

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, and Nektar::SolverUtils::EquationSystem::m_session.

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

141 {
142  EquationSystem::v_InitObject();
143 
144  const int nVel = m_fields[0]->GetCoordim(0);
145  int n;
146 
147  ASSERTL0(nVel > 1, "Linear elastic solver not set up for"
148  " this dimension (only 2D/3D supported).");
149 
150  // Make sure that we have Young's modulus and Poisson ratio set.
151  m_session->LoadParameter("E", m_E, 1.00);
152  m_session->LoadParameter("nu", m_nu, 0.25);
153  m_session->LoadParameter("Beta", m_beta, 1.00);
154 
155  // Create a coupled assembly map which allows us to tie u, v and w
156  // fields together.
157  if (nVel == 2)
158  {
159  MultiRegions::ContField2DSharedPtr u = boost::dynamic_pointer_cast<
163  m_graph,
164  u->GetLocalToGlobalMap(),
165  m_fields[0]->GetBndConditions(),
166  m_fields);
167  }
168 
169  if (nVel == 3)
170  {
171  MultiRegions::ContField3DSharedPtr u = boost::dynamic_pointer_cast<
173  m_assemblyMap = MemoryManager<CoupledAssemblyMap>
175  m_graph,
176  u->GetLocalToGlobalMap(),
177  m_fields[0]->GetBndConditions(),
178  m_fields);
179  }
180 
181  // Figure out size of our new matrix systems by looping over all expansions
182  // and multiply number of coefficients by velocity components.
183  const int nEl = m_fields[0]->GetExpSize();
185 
186  Array<OneD, unsigned int> sizeBnd(nEl);
187  Array<OneD, unsigned int> sizeInt(nEl);
188 
189  // Allocate storage for boundary and interior map caches.
192 
193  // Cache interior and boundary maps to avoid recalculating
194  // constantly. Should really be handled by manager in StdRegions but not
195  // really working yet.
196  for (n = 0; n < nEl; ++n)
197  {
198  exp = m_fields[0]->GetExp(m_fields[0]->GetOffset_Elmt_Id(n));
199  sizeBnd[n] = nVel * exp->NumBndryCoeffs();
200  sizeInt[n] = nVel * exp->GetNcoeffs() - sizeBnd[n];
201  exp->GetBoundaryMap(m_bmap[n]);
202  exp->GetInteriorMap(m_imap[n]);
203  }
204 
205  // Create block matrix storage for the statically condensed system.
206  MatrixStorage blkmatStorage = eDIAGONAL;
208  sizeBnd, sizeBnd, blkmatStorage);
210  sizeBnd, sizeInt, blkmatStorage);
212  sizeInt, sizeBnd, blkmatStorage);
214  sizeInt, sizeInt, blkmatStorage);
215 }
#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:293
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.
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
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.
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:191

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