-
Notifications
You must be signed in to change notification settings - Fork 10
/
ElasticBase.h
150 lines (124 loc) · 6.5 KB
/
ElasticBase.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
// $Id$
//==============================================================================
//!
//! \file ElasticBase.h
//!
//! \date Jul 04 2014
//!
//! \author Knut Morten Okstad / SINTEF
//!
//! \brief Base class representing FEM integrands for elasticity problems.
//!
//==============================================================================
#ifndef _ELASTIC_BASE_H
#define _ELASTIC_BASE_H
#include "IntegrandBase.h"
#include "Vec3.h"
class Material;
namespace TimeIntegration { class BDFD2; }
/*!
\brief Base class representing the FEM integrand of elasticity problems.
*/
class ElasticBase : public IntegrandBase
{
protected:
//! \brief The default constructor is protected to allow sub-classes only.
ElasticBase();
public:
//! \brief The destructor deletes the BDFD2 object, if any.
virtual ~ElasticBase();
//! \brief Parses material properties from a character string.
virtual Material* parseMatProp(char*, bool) { return nullptr; }
//! \brief Parses material properties from an XML-element.
virtual Material* parseMatProp(const tinyxml2::XMLElement*, bool)
{ return nullptr; }
//! \brief Parses a data section from an XML-element.
virtual bool parse(const tinyxml2::XMLElement* elem);
//! \brief Defines the material properties.
virtual void setMaterial(Material*) {}
//! \brief Returns the current material object.
virtual Material* getMaterial() const { return nullptr; }
//! \brief Initializes time integration parameters.
virtual bool init(const TimeDomain&) { return true; }
//! \brief Defines the gravitation vector.
void setGravity(double gx, double gy = 0.0, double gz = 0.0)
{ gravity.x = gx; gravity.y = gy; gravity.z = gz; }
//! \brief Returns the gravitation vector.
const Vec3& getGravity() const { return gravity; }
//! \brief Defines the number solution vectors.
void setNoSolutions(size_t n, bool dual = false) { nSV = n; dS = dual; }
//! \brief Defines the solution mode before the element assembly is started.
//! \param[in] mode The solution mode to use
virtual void setMode(SIM::SolutionMode mode);
//! \brief Updates the external load vector index for gradient calculation.
void setLoadGradientMode() { eS = 2; }
//! \brief Initializes an integration parameter for the integrand.
//! \param[in] i Index of the integration parameter to define
//! \param[in] prm The parameter value to assign
virtual void setIntegrationPrm(unsigned short int i, double prm);
//! \brief Returns an integration parameter for the integrand.
//! \param[in] i Index of the integration parameter to return
virtual double getIntegrationPrm(unsigned short int i) const;
//! \brief Advances the %BDF time step scheme one step forward.
void advanceStep(double dt, double dtn);
//! \brief Returns whether this integrand has explicit boundary contributions.
virtual bool hasBoundaryTerms() const { return false; }
//! \brief Returns the number of primary/secondary solution field components.
//! \param[in] fld which field set to consider (1=primary, 2=secondary)
virtual size_t getNoFields(int fld) const { return fld == 1 ? npv : 0; }
//! \brief Returns the name of a primary solution field component.
//! \param[in] i Field component index
//! \param[in] prefix Name prefix for all components
virtual std::string getField1Name(size_t i, const char* prefix) const;
//! \brief Evaluates the point load integrand at a specified point.
//! \param elmInt The local integral object to receive the contributions
//! \param[in] fe Finite element data of current integration point
//! \param[in] pval Load value at the specified point
virtual bool evalPoint(LocalIntegral& elmInt, const FiniteElement& fe,
const Vec3& pval);
using IntegrandBase::finalizeElement;
//! \brief Finalizes the element matrices after the numerical integration.
//! \param elmInt The local integral object to receive the contributions
//! \param[in] time Parameters for nonlinear and time-dependent simulations
//!
//! \details This method is used to pass time step size parameters to the
//! integrand in case of a dynamics simulation, where it is needed to compute
//! the effective stiffness/mass matrix used in the Newton iterations.
virtual bool finalizeElement(LocalIntegral& elmInt,
const TimeDomain& time, size_t = 0);
protected:
Vec3 gravity; //!< Gravitation vector
// Finite element quantities, i.e., indices into element matrices and vectors.
// These indices will be identical for all elements in a model and can thus
// be stored here, even when doing multi-threading. Note that the indices are
// 1-based, since the value zero is used to signal non-existing matrix/vector.
short int eKm; //!< Index to element material stiffness matrix
short int eKg; //!< Index to element geometric stiffness matrix
short int eM; //!< Index to element mass matrix
short int eS; //!< Index to element load vector
short int gS; //!< Index to element load gradient vector
short int iS; //!< Index to element internal force vector
short int dS; //!< Index to element dual force vector
short int nSV; //!< Number of consequtive solution vectors in core
std::vector<const char*> matNames; //!< Element matrix names (for debug print)
std::vector<const char*> vecNames; //!< Element vector names (for debug print)
//! \brief Newmark time integration parameters.
//! \details The interpretation of each parameter
//! depends on the actual simulator drivers, as follows: <UL>
//! <LI>[0]: Mass-proportional damping coefficient (Rayleigh damping).
//! <LI>[1]: Stiffness-proportional damping coefficient (Rayleigh damping).
//! <LI>[2]: α<SUB>H</SUB> for nonlinear Newmark drivers.
//! β or α<SUB>m</SUB> for linear Newmark drivers.
//! For linear drivers, a negative value signals that displacement increments
//! are used as primary unknowns, otherwise accelerations are used.
//! <LI>[3]: A positive value indicates that the solution driver is linear,
//! and the actual value is then the γ or α<SUB>f</SUB> parameter.
//! A zero value indicates a nonlinear driver, with stiffness-proportional
//! damping (if any) depending on both material and geometric stiffness.
//! A negative value means a nonlinear driver, with stiffness-proportional
//! damping (if any), depending on material stiffness only.
//! <LI>[4]: 1.0 if HHTSIM is used, 2.0 if GenAlphaSIM is used, otherwise 0.0.
double intPrm[5]; //!< </UL>
TimeIntegration::BDFD2* bdf; //!< BDF time discretization parameters
};
#endif