20#ifndef _FLEXIBLE_BODY_2S_13_DISK_H_
21#define _FLEXIBLE_BODY_2S_13_DISK_H_
23#include "mbsimFlexibleBody/flexible_body/2s_13.h"
24#include "mbsimFlexibleBody/flexible_body/fe/2s_13_disk.h"
26namespace MBSimFlexibleBody {
61 fmatvec::Vec3 evalPosition()
override;
62 fmatvec::SqrMat3 evalOrientation()
override;
68 void updateGyroscopicAccelerations(
Frame2s*
frame)
override;
70 void updatePositions(
int node)
override;
71 void updateVelocities(
int node)
override;
72 void updateAccelerations(
int node)
override;
73 void updateJacobians(
int node,
int j=0)
override;
74 void updateGyroscopicAccelerations(
int node)
override;
77 void init(
InitStage stage,
const MBSim::InitConfigSet &config)
override;
84 fmatvec::Vec
transformCW(
const fmatvec::Vec& WrPoint)
override;
plate according to Reissner-Mindlin with axial moving frame of reference
Definition: 2s_13_disk.h:42
void GlobalVectorContribution(int CurrentElement, const fmatvec::Vec &locVec, fmatvec::Vec &gloVec) override
insert 'local' information in global vectors
Definition: 2s_13_disk.cc:473
void BuildElements() override
references finite element coordinates to assembled coordinates
Definition: 2s_13_disk.cc:38
fmatvec::Vec transformCW(const fmatvec::Vec &WrPoint) override
transform Cartesian to cylinder system
Definition: 2s_13_disk.cc:337
~FlexibleBody2s13Disk() override=default
destructor
void initMatrices() override
calculate the matrices for the first time
Definition: 2s_13_disk.cc:352
void GlobalMatrixContribution(int CurrentElement, const fmatvec::SymMat &locMat, fmatvec::SymMat &gloMat) override
insert 'local' information in global matrices
void updateAG() override
update the transformation matrices A and G
Definition: 2s_13_disk.cc:463
void GlobalMatrixContribution(int CurrentElement, const fmatvec::Mat &locMat, fmatvec::Mat &gloMat) override
insert 'local' information in global matrices
FlexibleBody2s13Disk(const std::string &name)
constructor
Definition: 2s_13_disk.cc:34
plate according to Reissner-Mindlin with moving frame of reference
Definition: 2s_13.h:75
Definition: frame_2s.h:27
std::vector< Frame * > frame