mbsimflexiblebody  4.0.0
MBSim Flexible Body Module
2s_13_disk.h
1/* Copyright (C) 2004-2015 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 _FINITE_ELEMENT_2S_13_DISK_H_
21#define _FINITE_ELEMENT_2S_13_DISK_H_
22
23#include "mbsimFlexibleBody/discretization_interface.h"
24#include <cmath>
25
26namespace MBSimFlexibleBody {
27
39 public:
46 FiniteElement2s13Disk(double E_, double nu_, double rho_);
47
51 ~FiniteElement2s13Disk() override = default;
52
53 /* INTERFACE OF DISCRETIZATIONINTERFACE */
54 const fmatvec::SymMat& getM() const override { return M; }
55 const fmatvec::Vec& geth() const override;
56 const fmatvec::SqrMat& getdhdq() const override;
57 const fmatvec::SqrMat& getdhdu() const override;
58 int getqSize() const override { return RefDofs + 4*NodeDofs; }
59 int getuSize() const override { return RefDofs + 4*NodeDofs; }
60 void computeM(const fmatvec::Vec& q) override;
61 void computeh(const fmatvec::Vec& q,const fmatvec::Vec& u) override;
62 void computedhdz(const fmatvec::Vec& q,const fmatvec::Vec& u) override;
63 double computeKineticEnergy(const fmatvec::Vec& q,const fmatvec::Vec& u) override;
64 double computeGravitationalEnergy(const fmatvec::Vec& q) override;
65 double computeElasticEnergy(const fmatvec::Vec& q) override;
66 virtual fmatvec::Vec3 getPosition(const fmatvec::Vec& qElement, const fmatvec::Vec2 &s);
67 virtual fmatvec::SqrMat3 getOrientation(const fmatvec::Vec& qElement, const fmatvec::Vec2 &s);
68 virtual fmatvec::Vec3 getVelocity (const fmatvec::Vec& qElement, const fmatvec::Vec& qpElement, const fmatvec::Vec2 &s);
69 virtual fmatvec::Vec3 getAngularVelocity(const fmatvec::Vec& qElement, const fmatvec::Vec& qpElement, const fmatvec::Vec2 &s);
70 virtual fmatvec::Mat getJacobianOfMotion(const fmatvec::Vec& qElement, const fmatvec::Vec2 &s);
71 /***************************************************/
72
73 /* GETTER / SETTER */
74 const fmatvec::SymMat& getK() const { return K; }
75 void setEModul(double E_) { E = E_; }
76 void setPoissonRatio(double nu_) { nu = nu_; }
77 void setDensity(double rho_) { rho = rho_; }
78 void setShearCorrectionFactor(double alphaS_) { alphaS = alphaS_; }
79 /***************************************************/
80
87 void computeConstantSystemMatrices(const fmatvec::Vec &NodeCoordinates, double d1, double d2);
88
98 fmatvec::Vector<fmatvec::Fixed<6>, double> getPositions(const fmatvec::Vec &NodeCoordinates, const fmatvec::Vec &qElement, const fmatvec::Vec2 &s, double d1, double d2);
99
100 fmatvec::Vector<fmatvec::Fixed<6>, double> getVelocities(const fmatvec::Vec &NodeCoordinates, const fmatvec::Vec &qElement, const fmatvec::Vec &qpElement, const fmatvec::Vec2 &s, double d1, double d2);
101
107 fmatvec::Mat JGeneralized(const fmatvec::Vec &NodeCoordinates, const fmatvec::Vec2 &s);
108
109 private:
113 double E;
114
118 double nu;
119
123 double G;
124
128 double rho;
129
133 double alphaS;
134
139
144
148 int Nodes;
149
153 fmatvec::SymMat M, K;
154 };
155
156}
157
158#endif /* _FINITE_ELEMENT_2S_13_RCM_H_ */
discretization interface for flexible systems
Definition: discretization_interface.h:36
FE for Reissner-Mindlin Plate using MFR.
Definition: 2s_13_disk.h:38
int getqSize() const override
Definition: 2s_13_disk.h:58
double computeKineticEnergy(const fmatvec::Vec &q, const fmatvec::Vec &u) override
compute kinetic energy
Definition: 2s_13_disk.cc:457
int RefDofs
reference dof
Definition: 2s_13_disk.h:138
fmatvec::Vector< fmatvec::Fixed< 6 >, double > getPositions(const fmatvec::Vec &NodeCoordinates, const fmatvec::Vec &qElement, const fmatvec::Vec2 &s, double d1, double d2)
Definition: 2s_13_disk.cc:275
const fmatvec::SymMat & getM() const override
Definition: 2s_13_disk.h:54
FiniteElement2s13Disk(double E_, double nu_, double rho_)
constructor
Definition: 2s_13_disk.cc:35
double G
shear modulus
Definition: 2s_13_disk.h:123
double nu
Poisson ratio.
Definition: 2s_13_disk.h:118
int NodeDofs
elastic dof per node
Definition: 2s_13_disk.h:143
void computedhdz(const fmatvec::Vec &q, const fmatvec::Vec &u) override
compute Jacobian for implicit integration
Definition: 2s_13_disk.cc:453
const fmatvec::SqrMat & getdhdu() const override
Definition: 2s_13_disk.cc:441
double rho
density
Definition: 2s_13_disk.h:128
double computeElasticEnergy(const fmatvec::Vec &q) override
compute elastic energy
Definition: 2s_13_disk.cc:465
double alphaS
shear correction factor
Definition: 2s_13_disk.h:133
void computeh(const fmatvec::Vec &q, const fmatvec::Vec &u) override
compute smooth right hand side
Definition: 2s_13_disk.cc:449
void computeM(const fmatvec::Vec &q) override
compute mass matrix
Definition: 2s_13_disk.cc:445
double computeGravitationalEnergy(const fmatvec::Vec &q) override
compute gravitational energy
Definition: 2s_13_disk.cc:461
int Nodes
number of nodes
Definition: 2s_13_disk.h:148
~FiniteElement2s13Disk() override=default
destructor
const fmatvec::Vec & geth() const override
Definition: 2s_13_disk.cc:433
double E
Young's modulus.
Definition: 2s_13_disk.h:113
int getuSize() const override
Definition: 2s_13_disk.h:59
fmatvec::SymMat M
mass and stiffness matrix
Definition: 2s_13_disk.h:153
void computeConstantSystemMatrices(const fmatvec::Vec &NodeCoordinates, double d1, double d2)
computes mass and stiffness matrix
Definition: 2s_13_disk.cc:37
const fmatvec::SqrMat & getdhdq() const override
Definition: 2s_13_disk.cc:437
fmatvec::Mat JGeneralized(const fmatvec::Vec &NodeCoordinates, const fmatvec::Vec2 &s)
compute Jacobian of contact description at contour point
Definition: 2s_13_disk.cc:355