All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Pages
flexible_body_1s_33_ancf.h
1 /* Copyright (C) 2004-2014 MBSim Development Team
2  *
3  * This library is free software; you can redistribute it and/or
4  * modify it under the terms of the GNU Lesser General Public
5  * License as published by the Free Software Foundation; either
6  * version 2.1 of the License, or (at your option) any later version.
7  *
8  * This library is distributed in the hope that it will be useful,
9  * but WITHOUT ANY WARRANTY; without even the implied warranty of
10  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11  * Lesser General Public License for more details.
12  *
13  * You should have received a copy of the GNU Lesser General Public
14  * License along with this library; if not, write to the Free Software
15  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
16  *
17  * Contact: thorsten.schindler@mytum.de
18  */
19 
20 #ifndef _FLEXIBLE_BODY_1S_33_ANCF_H_
21 #define _FLEXIBLE_BODY_1S_33_ANCF_H_
22 
23 #include "mbsimFlexibleBody/flexible_body/flexible_body_1s.h"
24 
25 namespace MBSimFlexibleBody {
26 
42 
43  public:
49  FlexibleBody1s33ANCF(const std::string &name, bool openStructure);
50 
51  /* INHERITED INTERFACE OF FLEXIBLE BODY */
52  virtual void updateM() { }
53  virtual void updateLLM() { }
54 
55  virtual void BuildElements();
56 
57  virtual void GlobalVectorContribution(int n, const fmatvec::Vec& locVec, fmatvec::Vec& gloVec);
58  virtual void GlobalMatrixContribution(int n, const fmatvec::Mat& locMat, fmatvec::Mat& gloMat);
59  virtual void GlobalMatrixContribution(int n, const fmatvec::SymMat& locMat, fmatvec::SymMat& gloMat);
60 
61  virtual void updatePositions(Frame1s* frame);
62  virtual void updateVelocities(Frame1s* frame);
63  virtual void updateAccelerations(Frame1s* frame);
64  virtual void updateJacobians(Frame1s* frame, int j=0);
65  virtual void updateGyroscopicAccelerations(Frame1s* frame);
66 
67  virtual void updatePositions(NodeFrame* frame);
68  virtual void updateVelocities(NodeFrame* frame);
69  virtual void updateAccelerations(NodeFrame* frame);
70  virtual void updateJacobians(NodeFrame* frame, int j=0);
71  virtual void updateGyroscopicAccelerations(NodeFrame* frame);
72  /****************************************/
73 
74  /* INHERITED INTERFACE OF OBJECT */
75  virtual void init(InitStage stage);
76  /***************************************************/
77 
78  /* INHERITED INTERFACE OF ELEMENT */
79  virtual void plot();
80  virtual std::string getType() const { return "FlexibleBody1s33ANCF"; }
81  /***************************************************/
82 
83  /* GETTER / SETTER */
84  void setNumberElements(int n);
85  void setEModul(double E_) { E = E_; }
86  void setShearModul(double G_) { G = G_; }
87  void setCrossSectionalArea(double A_) { A = A_; }
88  void setMomentInertia(double I0_,double I1_,double I2_) { I0 = I0_; I1 = I1_; I2 = I2_; }
89  void setDensity(double rho_) { rho = rho_; }
90  void setCurlRadius(double rc1_,double rc2_);
91  int getNumberElements(){ return Elements; }
92  /***************************************************/
93 
97  void initInfo();
98 
103  void initRelaxed(double alpha);
104 
105  private:
109  int Elements;
110 
114  double l0;
115 
119  double E;
120 
124  double G;
125 
129  double A;
130 
134  double I0;
135 
139  double I1;
140 
144  double I2;
145 
149  double rho;
150 
154  double rc1;
155 
159  double rc2;
160 
165 
169  void initM();
170 
177  void BuildElement(const double& sGlobal, double& sLocal, int& currentElement);
178 
183 
188  };
189 
190 }
191 #endif /* _FLEXIBLE_BODY_1S_33_ANCF_H_ */
void BuildElement(const double &sGlobal, double &sLocal, int &currentElement)
detect current finite element
Definition: flexible_body_1s_33_ancf.cc:234
double I0
polar moment of inertia of cross-section
Definition: flexible_body_1s_33_ancf.h:134
void initInfo()
initialise beam only for giving information with respect to state, number elements, length, (not for simulation)
Definition: flexible_body_1s_33_ancf.cc:263
FlexibleBody1s33ANCF & operator=(const FlexibleBody1s33ANCF &)
assignment operator is declared private
FlexibleBody1s33ANCF(const std::string &name, bool openStructure)
constructor:
double I2
moment of inertia of cross-section
Definition: flexible_body_1s_33_ancf.h:144
double rc1
radius of undeformed shape
Definition: flexible_body_1s_33_ancf.h:154
double I1
moment of inertia of cross-section
Definition: flexible_body_1s_33_ancf.h:139
tbd
Definition: flexible_body_1s.h:33
bool initialised
flag for testing if beam is initialised
Definition: flexible_body_1s_33_ancf.h:164
Definition: frame_1s.h:27
double rho
material density
Definition: flexible_body_1s_33_ancf.h:149
int Elements
number of finite elements used for discretisation
Definition: flexible_body_1s_33_ancf.h:109
Absolute Nodal Coordinate Formulation for flexible spatial beams.
Definition: flexible_body_1s_33_ancf.h:41
bool openStructure
flag for open (cantilever beam) or closed (rings) structures
Definition: flexible_body_1s.h:73
virtual void BuildElements()
references finite element coordinates to assembled coordinates
Definition: flexible_body_1s_33_ancf.cc:216
void initRelaxed(double alpha)
initialise beam state concerning a straight cantilever setting or a circle shaped ring ...
Definition: flexible_body_1s_33_ancf.cc:274
void initM()
initialize mass matrix and calculate Cholesky decomposition
Definition: flexible_body_1s_33_ancf.cc:250
double l0
length of one finite element
Definition: flexible_body_1s_33_ancf.h:114
double rc2
radius of undeformed shape
Definition: flexible_body_1s_33_ancf.h:159
std::vector< Frame * > frame
double A
cross-section area
Definition: flexible_body_1s_33_ancf.h:129
virtual void GlobalVectorContribution(int n, const fmatvec::Vec &locVec, fmatvec::Vec &gloVec)
insert &#39;local&#39; information in global vectors
Definition: flexible_body_1s_33_ancf.cc:40
double G
shear modulus
Definition: flexible_body_1s_33_ancf.h:124
cartesian frame on nodes of flexible bodies
Definition: node_frame.h:31
virtual void GlobalMatrixContribution(int n, const fmatvec::Mat &locMat, fmatvec::Mat &gloMat)
insert &#39;local&#39; information in global matrices
Definition: flexible_body_1s_33_ancf.cc:52
double E
modulus of linear elasticity
Definition: flexible_body_1s_33_ancf.h:119

Impressum / Disclaimer / Datenschutz Generated by doxygen 1.8.5 Valid HTML