21 #ifndef _FLEXIBLE_BODY_H_
22 #define _FLEXIBLE_BODY_H_
24 #include "mbsim/objects/body.h"
25 #include "namespace.h"
28 class FixedRelativeFrame;
32 namespace MBSimFlexibleBody {
34 class DiscretizationInterface;
66 virtual void updatedq() { dq = u*getStepSize(); }
67 virtual void updateqd() {
qd = u; }
68 virtual void updateh(
int k=0);
69 virtual void updateM();
70 virtual void updatedhdz();
72 virtual void updateVelocities(
NodeFrame* frame);
73 virtual void updateAccelerations(
NodeFrame* frame);
74 virtual void updateJacobians(
NodeFrame* frame,
int j=0);
75 virtual void updateGyroscopicAccelerations(
NodeFrame* frame);
79 virtual std::string getType()
const {
return "FlexibleBody"; }
80 virtual void initializeUsingXML(xercesc::DOMElement *element);
84 virtual void init(InitStage stage);
85 virtual double computeKineticEnergy();
86 virtual double computePotentialEnergy();
88 virtual void setq0(
fmatvec::Vec q0_) {
if(q0_.size()) MBSim::Body::setGeneralizedInitialPosition(q0_);
q<<
q0; }
89 virtual void setu0(
fmatvec::Vec u0_) {
if(u0_.size()) MBSim::Body::setGeneralizedInitialVelocity(u0_); u<<u0; }
102 virtual fmatvec::Vec3 getAngles(
int i) {
return fmatvec::Vec3(); }
103 virtual fmatvec::Vec3 getDerivativeOfAngles(
int i) {
return fmatvec::Vec3(); }
184 virtual void exportPositionVelocity(
const std::string & filenamePos,
const std::string & filenameVel = std::string(),
const int & deg = 3,
const bool &writePsFile =
false){
throw MBSim::MBSimError(
"exportPositionVelocity(const std::string& filenamePos, const std::string& filenameVel, const int& deg, const bool& writePsFile) is not implemented for " + this->getType()) ;}
191 virtual void importPositionVelocity(
const std::string& filenamePos,
const std::string& filenameVel = std::string()){
throw MBSim::MBSimError(
"importPositionVelocity(const std::string& filenamePos, const std::string& filenameVel) is not implemented for " + this->getType()) ;}
193 void resetUpToDate();
249 virtual std::string getType()
const {
return "FlexibleBodyContinuum"; }
252 void setContourNodes(
const std::vector<AT> nodes) {
userContourNodes = nodes; }
254 void setNodeOffset(
const AT nodeOffset_){
nodeOffset = nodeOffset_;}
255 AT getNodeOffset()
const {
return nodeOffset;}
virtual void BuildElements()=0
references finite element coordinates to assembled coordinates
upmost class for flexible body implementation
Definition: flexible_body.h:52
std::vector< DiscretizationInterface * > discretization
stl-vector of discretizations/finite elements
Definition: flexible_body.h:199
double d_massproportional
damping factor for mass proportion, see BodyFlexible::setMassProportionalDamping() ...
Definition: flexible_body.h:214
void addFrame(NodeFrame *frame)
Definition: flexible_body.cc:148
FlexibleBody(const std::string &name)
constructor
Definition: flexible_body.cc:44
AT nodeOffset
offset of the ROTNODE from the TRANSNODE
Definition: flexible_body.h:266
virtual void GlobalMatrixContribution(int CurrentElement, const fmatvec::Mat &locMat, fmatvec::Mat &gloMat)=0
insert 'local' information in global matrices
virtual void importPositionVelocity(const std::string &filenamePos, const std::string &filenameVel=std::string())
imports the interpolated position and optional the velocity files (created with exportPositionVelocit...
Definition: flexible_body.h:191
virtual ~FlexibleBody()
destructor
Definition: flexible_body.cc:46
flexible body entirely described within MBSim holding all informations about continuum approximations...
Definition: flexible_body.h:240
void setMassProportionalDamping(const double d_)
cartesian kinematic for contour or external frame (normal, tangent, binormal) is set by implementatio...
Definition: flexible_body.h:158
std::vector< fmatvec::Vec > uElement
stl-vector of finite element wise velocities
Definition: flexible_body.h:209
std::vector< fmatvec::Vec > qElement
stl-vector of finite element wise positions
Definition: flexible_body.h:204
std::vector< AT > userContourNodes
grid for contact point detection
Definition: flexible_body.h:261
bool updEle
vector of contour parameters each describing a frame
Definition: flexible_body.h:230
std::vector< Frame * > frame
virtual void GlobalVectorContribution(int CurrentElement, const fmatvec::Vec &locVec, fmatvec::Vec &gloVec)=0
insert 'local' information in global vectors
cartesian frame on nodes of flexible bodies
Definition: node_frame.h:31
virtual void exportPositionVelocity(const std::string &filenamePos, const std::string &filenameVel=std::string(), const int °=3, const bool &writePsFile=false)
interpolates the position and optional the velocity coordinates of the flexible body with Nurbs-packa...
Definition: flexible_body.h:184