MBDyn-1.7.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups
strnode.cc File Reference
#include "mbconfig.h"
#include <cerrno>
#include <cstring>
#include <sstream>
#include "mynewmem.h"
#include "strnode.h"
#include "body.h"
#include "autostr.h"
#include "dataman.h"
#include "matvecexp.h"
#include "Rot.hh"
Include dependency graph for strnode.cc:

Go to the source code of this file.

Functions

NodeReadStructNode (DataManager *pDM, MBDynParser &HP, DofOwner *pDO, unsigned int uLabel)
 

Variables

static const char xyz [] = "xyz"
 
static const char * sdn_dof []
 
static const char * sdn_eq []
 
static const char * sdn_initial_dof []
 
static const char * sdn_initial_eq []
 
static const char * sn_dof []
 
static const char * sn_eq []
 
static const char * sn_modal_eq []
 
static const char * sn_initial_dof []
 
static const char * sn_initial_eq []
 

Function Documentation

Node* ReadStructNode ( DataManager pDM,
MBDynParser HP,
DofOwner pDO,
unsigned int  uLabel 
)

Definition at line 3978 of file strnode.cc.

References AbsRefFrame, ASSERT, Elem::AUTOMATICSTRUCTURAL, DataManager::bDoesOmegaRotate(), DataManager::bIsStaticModel(), DataManager::bOutputAccelerations(), DEBUGCOUT, DataManager::dGetInitialPositionStiffness(), DataManager::dGetInitialVelocityStiffness(), dRaDegr, DataManager::dReadScale(), EULER_123, EULER_313, EULER_321, Eye3, DataManager::fReadOutput(), IncludeParser::GetLineData(), DataManager::GetLogFile(), MBDynParser::GetOmeAbs(), MBDynParser::GetPosAbs(), MBDynParser::GetPosRel(), StructNode::GetRCurr(), HighParser::GetReal(), MBDynParser::GetRotAbs(), MBDynParser::GetRotRel(), MBDynParser::GetVelAbs(), HighParser::GetWord(), StructDispNode::GetXCurr(), HighParser::GetYesNoOrBool(), DataManager::IncElemCount(), HighParser::IsArg(), HighParser::IsKeyWord(), LASTKEYWORD, MatR2EulerAngles123(), MatR2EulerAngles313(), MatR2EulerAngles321(), MBDYN_EXCEPT_ARGS, ORIENTATION_MATRIX, ORIENTATION_VECTOR, StructDispNode::OUTPUT_ACCELERATIONS, StructDispNode::OUTPUT_INERTIA, DataManager::pGetRBK(), R, DataManager::ReadNode(), ReadOptionalOrientationDescription(), MBDynParser::RF, SAFENEWWITHCONSTRUCTOR, DofOwner::SetScale(), Node::STRUCTURAL, DofOwner::STRUCTURALNODE, HighParser::UNKNOWN, UNKNOWN_ORIENTATION_DESCRIPTION, RotManip::VecRot(), Write(), Vec3::Write(), Mat3x3::Write(), and Zero3.

Referenced by DataManager::ReadNodes().

3982 {
3983  DEBUGCOUT("Entering ReadStructNode(" << uLabel << ")" << std::endl);
3984 
3985  const char* sKeyWords[] = {
3986  "static" "displacement",
3987  "dynamic" "displacement",
3988  "static",
3989  "dynamic",
3990  "modal",
3991  "dummy",
3992  "offset",
3993  "relative" "frame",
3994  0
3995  };
3996 
3997  /* enum delle parole chiave */
3998  enum KeyWords {
3999  UNKNOWN = -1,
4000 
4001  STATIC_DISP = 0,
4002  DYNAMIC_DISP,
4003  STATIC,
4004  DYNAMIC,
4005  MODAL,
4006  DUMMY,
4007 
4008  OFFSET,
4009  RELATIVEFRAME,
4010 
4011  LASTKEYWORD
4012  };
4013 
4014  /* tabella delle parole chiave */
4015  KeyTable K(HP, sKeyWords);
4016 
4017  /* lettura dati specifici */
4018  KeyWords CurrType((KeyWords)HP.IsKeyWord());
4019 
4020  /*
4021  * explicit node type required; default is no longer "DYNAMIC"
4022  */
4023  if (CurrType == UNKNOWN) {
4024  silent_cerr("StructNode(" << uLabel << "): "
4025  "missing node type at line " << HP.GetLineData()
4026  << std::endl);
4028  }
4029 
4030 #ifdef DEBUG
4031  switch (CurrType) {
4032  case STATIC_DISP:
4033  std::cout << "Static structural displacement node" << std::endl;
4034  break;
4035  case DYNAMIC_DISP:
4036  std::cout << "Dynamic structural displacement node" << std::endl;
4037  break;
4038  case STATIC:
4039  std::cout << "Static structural node" << std::endl;
4040  break;
4041  case DYNAMIC:
4042  std::cout << "Dynamic structural node" << std::endl;
4043  break;
4044  case DUMMY:
4045  std::cout << "Dummy structural node" << std::endl;
4046  break;
4047  case MODAL:
4048  std::cout << "Modal node" << std::endl;
4049  break;
4050  default:
4051  std::cout << "Unknown structural node" << std::endl;
4052  break;
4053  }
4054 #endif /* DEBUG */
4055 
4056  StructDispNode* pNd = NULL;
4058  KeyWords DummyType = UNKNOWN;
4059  if (CurrType == DUMMY) {
4060  const StructNode* pNode = pDM->ReadNode<const StructNode, Node::STRUCTURAL>(HP);
4061 
4062  DummyType = KeyWords(HP.GetWord());
4063  switch (DummyType) {
4064  case OFFSET: {
4065  ReferenceFrame RF(pNode);
4066  Vec3 f(HP.GetPosRel(RF));
4067  Mat3x3 R(HP.GetRotRel(RF));
4068 
4069  od = ReadOptionalOrientationDescription(pDM, HP);
4070 
4071  flag fOut = pDM->fReadOutput(HP, Node::STRUCTURAL);
4074  OffsetDummyStructNode(uLabel, pDO, pNode,
4075  f, R, od, fOut));
4076  } break;
4077 
4078  case RELATIVEFRAME: {
4079  const StructNode* pNodeRef = pDM->ReadNode<const StructNode, Node::STRUCTURAL>(HP);
4080 
4081  ReferenceFrame RF(pNodeRef);
4082 
4083  Vec3 fh(Zero3);
4084  if (HP.IsKeyWord("position")) {
4085  fh = HP.GetPosRel(RF);
4086  }
4087 
4088  Mat3x3 Rh(Eye3);
4089  if (HP.IsKeyWord("orientation")) {
4090  Rh = HP.GetRotRel(RF);
4091 
4092  od = ReadOptionalOrientationDescription(pDM, HP);
4093  }
4094 
4095  const StructNode *pNodeRef2 = 0;
4096  Vec3 fh2(Zero3);
4097  Mat3x3 Rh2(Eye3);
4098  if (HP.IsKeyWord("pivot" "node")) {
4099  pNodeRef2 = pDM->ReadNode<const StructNode, Node::STRUCTURAL>(HP);
4100 
4101  if (HP.IsKeyWord("position")) {
4102  fh2 = HP.GetPosRel(RF);
4103  }
4104 
4105  if (HP.IsKeyWord("orientation")) {
4106  Rh2 = HP.GetRotRel(RF);
4107  }
4108  }
4109 
4110  flag fOut = pDM->fReadOutput(HP, Node::STRUCTURAL);
4111 
4112  if (pNodeRef2) {
4115  PivotRelFrameDummyStructNode(uLabel, pDO,
4116  pNode, pNodeRef, fh, Rh,
4117  pNodeRef2, fh2, Rh2, od, fOut));
4118 
4119  } else {
4122  RelFrameDummyStructNode(uLabel, pDO,
4123  pNode, pNodeRef, fh, Rh, od, fOut));
4124  }
4125  } break;
4126 
4127  default:
4128  silent_cerr("StructNode(" << uLabel << "): "
4129  "unknown dummy node type "
4130  "at line " << HP.GetLineData() << std::endl);
4132  }
4133 
4134  } else {
4135  bool bDisp(false);
4136 
4137  if (CurrType == STATIC_DISP || CurrType == DYNAMIC_DISP) {
4138  bDisp = true;
4139  }
4140 
4141  /* posizione (vettore di 3 elementi) */
4142  if (!HP.IsKeyWord("position")) {
4143  pedantic_cerr("StructNode(" << uLabel << "): "
4144  "missing keyword \"position\" at line "
4145  << HP.GetLineData() << std::endl);
4146  }
4147  Vec3 X0(HP.GetPosAbs(::AbsRefFrame));
4148  DEBUGCOUT("X0 =" << std::endl << X0 << std::endl);
4149 
4150  /* sistema di riferimento (trucco dei due vettori) */
4151  Mat3x3 R0;
4152  if (!bDisp) {
4153  if (!HP.IsKeyWord("orientation")) {
4154  pedantic_cerr("StructNode(" << uLabel << "): "
4155  "missing keyword \"orientation\" at line "
4156  << HP.GetLineData() << std::endl);
4157  }
4158  R0 = HP.GetRotAbs(::AbsRefFrame);
4159 
4160  od = ReadOptionalOrientationDescription(pDM, HP);
4161 
4162  DEBUGCOUT("R0 =" << std::endl << R0 << std::endl);
4163  }
4164 
4165  /* Velocita' iniziali (due vettori di 3 elementi, con la possibilita'
4166  * di usare "null" per porli uguali a zero) */
4167  if (!HP.IsKeyWord("velocity")) {
4168  pedantic_cerr("StructNode(" << uLabel << "): "
4169  "missing keyword \"velocity\" at line "
4170  << HP.GetLineData() << std::endl);
4171  }
4172  Vec3 XPrime0(HP.GetVelAbs(::AbsRefFrame, X0));
4173  DEBUGCOUT("Xprime0 =" << std::endl << XPrime0 << std::endl);
4174 
4175  Vec3 Omega0;
4176  if (!bDisp) {
4177  if (!HP.IsKeyWord("angular" "velocity")) {
4178  pedantic_cerr("StructNode(" << uLabel << "): "
4179  "missing keyword \"angular velocity\" at line "
4180  << HP.GetLineData() << std::endl);
4181  }
4182  Omega0 = HP.GetOmeAbs(::AbsRefFrame);
4183  DEBUGCOUT("Omega0 =" << std::endl << Omega0 << std::endl);
4184  }
4185 
4186  const StructNode *pRefNode = 0;
4187  if (HP.IsKeyWord("prediction" "node")) {
4188  switch (CurrType) {
4189  case STATIC_DISP:
4190  case DYNAMIC_DISP:
4191  case STATIC:
4192  case DYNAMIC:
4193  break;
4194 
4195  default:
4196  silent_cerr("StructNode(" << uLabel << "): "
4197  "prediction node allowed "
4198  "for static and dynamic nodes only, "
4199  "at line " << HP.GetLineData()
4200  << std::endl);
4202  }
4203  pRefNode = pDM->ReadNode<const StructNode, Node::STRUCTURAL>(HP);
4204 
4205 #ifndef MBDYN_X_RELATIVE_PREDICTION
4206  silent_cerr("warning, relative prediction disabled; "
4207  "absolute prediction will be used" << std::endl);
4208 #endif /* ! MBDYN_X_RELATIVE_PREDICTION */
4209  }
4210 
4211  const RigidBodyKinematics *pRBK = pDM->pGetRBK();
4212 
4213  /* Rigidezza in assemblaggio diversa da quella di default
4214  * e flag di output */
4215  doublereal dPosStiff = pDM->dGetInitialPositionStiffness();
4216  doublereal dVelStiff = pDM->dGetInitialVelocityStiffness();
4217  bool bOmRot = pDM->bDoesOmegaRotate();
4218 
4219  if (HP.IsArg()) {
4220  if (HP.IsKeyWord("assembly")) {
4221  dPosStiff = HP.GetReal(dPosStiff);
4222  dVelStiff = HP.GetReal(dVelStiff);
4223 
4224  DEBUGCOUT("Initial position stiffness: " << dPosStiff << std::endl);
4225  DEBUGCOUT("Initial velocity stiffness: " << dVelStiff << std::endl);
4226 
4227  if (!bDisp) {
4228  bOmRot = HP.GetYesNoOrBool(bOmRot);
4229 
4230  DEBUGCOUT("Omega rotates? : " << (bOmRot ? "yes" : "no") << std::endl);
4231  }
4232  }
4233  }
4234 
4236 
4237  flag fOut = pDM->fReadOutput(HP, Node::STRUCTURAL);
4238  switch (CurrType) {
4239  case DYNAMIC_DISP:
4240  case DYNAMIC:
4241  case MODAL:
4242  {
4243  bool bGotAccels(false);
4244  bool bAccels(pDM->bOutputAccelerations());
4245 
4246  bool bGotInertia(false);
4247  bool bInertia(false);
4248 
4249  while (HP.IsArg()) {
4250  if (HP.IsKeyWord("accelerations")) {
4251  if (bGotAccels) {
4252  silent_cerr("StructNode(" << uLabel << "): "
4253  "\"accelerations\" already set, "
4254  "repeated at line " << HP.GetLineData()
4255  << std::endl);
4257  }
4258  bGotAccels = true;
4259 
4260  if (HP.IsArg()) {
4261  if (HP.GetYesNoOrBool(false)) {
4262  bAccels = true;
4263 
4264  } else {
4265  bAccels = false;
4266  }
4267 
4268  } else {
4269  // deprecated
4270  silent_cout("StructNode(" << uLabel << "): "
4271  "warning, \"accelerations\" needs \"yes\" or \"no\" "
4272  "at line " << HP.GetLineData() << std::endl);
4273  bAccels = true;
4274  }
4275 
4276  } else if (HP.IsKeyWord("output" "inertia")) {
4277  if (bGotInertia) {
4278  silent_cerr("StructNode(" << uLabel << "): "
4279  "\"inertia\" already set, "
4280  "repeated at line " << HP.GetLineData()
4281  << std::endl);
4283  }
4284  bGotInertia = true;
4285 
4286  if (!HP.IsArg()) {
4287  silent_cerr("StructNode(" << uLabel << "): "
4288  "\"inertia\" needs \"yes\" or \"no\" "
4289  "at line " << HP.GetLineData() << std::endl);
4291  }
4292 
4293  if (HP.GetYesNoOrBool(true)) {
4294  bInertia = true;
4295 
4296  } else {
4297  bInertia = false;
4298  }
4299 
4300  } else {
4301  break;
4302  }
4303  }
4304 
4305  if (fOut) {
4306  // restore legacy behavior
4307  if (!bGotInertia) {
4308  bInertia = true;
4309  }
4310 
4311  if (bInertia) {
4313  }
4314 
4315  if (bAccels) {
4317  }
4318  }
4319  } break;
4320 
4321  default:
4322  break;
4323  }
4324 
4325  if (CurrType == DYNAMIC && pDM->bIsStaticModel()) {
4326  pedantic_cout("DynamicStructNode(" << uLabel << ") turned into static" << std::endl);
4327  CurrType = STATIC;
4328 
4329  } else if (CurrType == DYNAMIC_DISP && pDM->bIsStaticModel()) {
4330  pedantic_cout("DynamicStructDispNode(" << uLabel << ") turned into static" << std::endl);
4331  CurrType = STATIC_DISP;
4332  }
4333 
4334  /* Se non c'e' il punto e virgola finale */
4335  if (HP.IsArg()) {
4336  silent_cerr("ReadStructNode(" << uLabel << "): semicolon expected "
4337  "at line " << HP.GetLineData() << std::endl);
4339  }
4340 
4341  /* costruzione del nodo */
4342  switch (CurrType) {
4343  case STATIC_DISP:
4345  StaticStructDispNode(uLabel, pDO,
4346  X0,
4347  XPrime0,
4348  pRefNode,
4349  pRBK,
4350  dPosStiff, dVelStiff,
4351  od, fOut));
4352  break;
4353 
4354  case DYNAMIC_DISP:
4356  DynamicStructDispNode(uLabel, pDO,
4357  X0,
4358  XPrime0,
4359  pRefNode,
4360  pRBK,
4361  dPosStiff, dVelStiff,
4362  od, fOut));
4363 
4364  /* Incrementa il numero di elementi automatici dei nodi dinamici */
4366  break;
4367 
4368  case STATIC:
4370  StaticStructNode(uLabel, pDO,
4371  X0, R0,
4372  XPrime0, Omega0,
4373  pRefNode,
4374  pRBK,
4375  dPosStiff, dVelStiff,
4376  bOmRot, od, fOut));
4377  break;
4378 
4379  case DYNAMIC:
4381  DynamicStructNode(uLabel, pDO,
4382  X0, R0,
4383  XPrime0, Omega0,
4384  pRefNode,
4385  pRBK,
4386  dPosStiff, dVelStiff,
4387  bOmRot, od, fOut));
4388 
4389  /* Incrementa il numero di elementi automatici dei nodi dinamici */
4391  break;
4392 
4393  case MODAL:
4395  ModalNode(uLabel, pDO,
4396  X0, R0,
4397  XPrime0, Omega0,
4398  pRBK,
4399  dPosStiff, dVelStiff,
4400  bOmRot, od, fOut));
4401  break;
4402 
4403  default:
4404  ASSERT(false);
4406  }
4407  }
4408 
4409  std::ostream& out = pDM->GetLogFile();
4410 
4411  const char *description = "structural node: ";
4412 
4413  switch (CurrType){
4414  case DUMMY:
4415  switch (DummyType){
4416  case RELATIVEFRAME:
4417  description = "relative frame structural node: ";
4418  break;
4419 
4420  default:
4421  break;
4422  }
4423  break;
4424 
4425  default:
4426  break;
4427  }
4428 
4429  out << description << uLabel
4430  << " ", pNd->GetXCurr().Write(out, " ")
4431  << " ";
4432  const StructNode *pSN = dynamic_cast<const StructNode *>(pNd);
4433  if (pSN) {
4434  switch (od) {
4435  case EULER_123:
4436  out << "euler123 ",
4437  (MatR2EulerAngles123(pSN->GetRCurr())*dRaDegr).Write(out, " ");
4438  break;
4439 
4440  case EULER_313:
4441  out << "euler313 ",
4442  (MatR2EulerAngles313(pSN->GetRCurr())*dRaDegr).Write(out, " ");
4443  break;
4444 
4445  case EULER_321:
4446  out << "euler321 ",
4447  (MatR2EulerAngles321(pSN->GetRCurr())*dRaDegr).Write(out, " ");
4448  break;
4449 
4450  case ORIENTATION_VECTOR:
4451  out << "phi ",
4452  RotManip::VecRot(pSN->GetRCurr()).Write(out, " ");
4453  break;
4454 
4455  case ORIENTATION_MATRIX:
4456  out << "mat ",
4457  pSN->GetRCurr().Write(out, " ");
4458  break;
4459 
4460  default:
4461  /* impossible */
4462  break;
4463  }
4464 
4465  } else {
4466  switch (od) {
4467  case EULER_123:
4468  out << "euler123 ",
4469  ::Zero3.Write(out, " ");
4470  break;
4471 
4472  case EULER_313:
4473  out << "euler313 ",
4474  ::Zero3.Write(out, " ");
4475  break;
4476 
4477  case EULER_321:
4478  out << "euler321 ",
4479  ::Zero3.Write(out, " ");
4480  break;
4481 
4482  case ORIENTATION_VECTOR:
4483  out << "phi ",
4484  ::Zero3.Write(out, " ");
4485  break;
4486 
4487  case ORIENTATION_MATRIX:
4488  out << "mat ",
4489  ::Eye3.Write(out, " ");
4490  break;
4491 
4492  default:
4493  /* impossible */
4494  break;
4495  }
4496  }
4497 
4498  out << std::endl;
4499 
4500  ASSERT(pNd != NULL);
4501 
4502  return pNd;
4503 } /* End of ReadStructNode() */
flag fReadOutput(MBDynParser &HP, const T &t) const
Definition: dataman.h:1064
Mat3x3 GetRotRel(const ReferenceFrame &rf)
Definition: mbpar.cc:1795
const Vec3 Zero3(0., 0., 0.)
std::ostream & Write(std::ostream &out, const char *sFill=" ") const
Definition: matvec3.cc:738
long int flag
Definition: mbdyn.h:43
#define MBDYN_EXCEPT_ARGS
Definition: except.h:63
Definition: matvec3.h:98
const doublereal & dGetInitialVelocityStiffness(void) const
Definition: dataman2.cc:117
virtual const Mat3x3 & GetRCurr(void) const
Definition: strnode.h:1012
Mat3x3 GetRotAbs(const ReferenceFrame &rf)
Definition: mbpar.cc:1857
std::ostream & Write(std::ostream &out, const FullMatrixHandler &m, const char *s, const char *s2)
Definition: fullmh.cc:376
Vec3 GetVelAbs(const ReferenceFrame &rf, const Vec3 &x)
Definition: mbpar.cc:1482
const Mat3x3 Eye3(1., 0., 0., 0., 1., 0., 0., 0., 1.)
OrientationDescription
Definition: matvec3.h:1597
OrientationDescription ReadOptionalOrientationDescription(DataManager *pDM, MBDynParser &HP)
Definition: dataman3.cc:2531
bool bIsStaticModel(void) const
Definition: dataman.h:482
Vec3 VecRot(const Mat3x3 &Phi)
Definition: Rot.cc:136
virtual bool GetYesNoOrBool(bool bDefval=false)
Definition: parser.cc:1038
doublereal dReadScale(MBDynParser &HP, enum DofOwner::Type t) const
Definition: dataman3.cc:1491
void SetScale(const doublereal &d)
Definition: dofown.cc:44
const ReferenceFrame AbsRefFrame(0, Vec3(0., 0., 0), Mat3x3(1., 0., 0., 0., 1., 0., 0., 0., 1.), Vec3(0., 0., 0), Vec3(0., 0., 0), EULER_123)
Vec3 MatR2EulerAngles313(const Mat3x3 &R)
Definition: matvec3.cc:927
void IncElemCount(Elem::Type type)
Definition: dataman2.cc:129
Vec3 GetPosRel(const ReferenceFrame &rf)
Definition: mbpar.cc:1331
virtual bool IsKeyWord(const char *sKeyWord)
Definition: parser.cc:910
Vec3 GetOmeAbs(const ReferenceFrame &rf)
Definition: mbpar.cc:1551
bool bDoesOmegaRotate(void) const
Definition: dataman2.cc:123
#define DEBUGCOUT(msg)
Definition: myassert.h:232
Vec3 GetPosAbs(const ReferenceFrame &rf)
Definition: mbpar.cc:1401
Vec3 MatR2EulerAngles123(const Mat3x3 &R)
Definition: matvec3.cc:893
const doublereal dRaDegr
Definition: matvec3.cc:884
#define ASSERT(expression)
Definition: colamd.c:977
KeyWords
Definition: dataman4.cc:94
#define SAFENEWWITHCONSTRUCTOR(pnt, item, constructor)
Definition: mynewmem.h:698
virtual const Vec3 & GetXCurr(void) const
Definition: strnode.h:310
const RigidBodyKinematics * pGetRBK(void) const
Definition: dataman.h:485
virtual bool IsArg(void)
Definition: parser.cc:807
const doublereal & dGetInitialPositionStiffness(void) const
Definition: dataman2.cc:111
virtual int GetWord(void)
Definition: parser.cc:1083
std::ostream & GetLogFile(void) const
Definition: dataman.h:326
Vec3 MatR2EulerAngles321(const Mat3x3 &R)
Definition: matvec3.cc:941
virtual HighParser::ErrOut GetLineData(void) const
Definition: parsinc.cc:697
std::ostream & Write(std::ostream &out, const char *sFill=" ", const char *sFill2=NULL) const
Definition: matvec3.cc:722
Node * ReadNode(MBDynParser &HP, Node::Type type) const
Definition: dataman3.cc:2309
Mat3x3 R
bool bOutputAccelerations(void) const
Definition: dataman.cc:878
virtual doublereal GetReal(const doublereal &dDefval=0.0)
Definition: parser.cc:1056

Here is the call graph for this function:

Variable Documentation

const char* sdn_dof[]
static
Initial value:
= {
"position P",
"momentum B"
}

Definition at line 49 of file strnode.cc.

Referenced by DynamicStructDispNode::DescribeDof().

const char* sdn_eq[]
static
Initial value:
= {
"momentum definition B",
"force equilibrium F"
}

Definition at line 53 of file strnode.cc.

Referenced by StructDispNode::DescribeEq(), and DynamicStructDispNode::DescribeEq().

const char* sdn_initial_dof[]
static
Initial value:
= {
"position P",
"velocity v"
}

Definition at line 57 of file strnode.cc.

Referenced by StructDispNode::DescribeDof().

const char* sdn_initial_eq[]
static
Initial value:
= {
"position constraint P",
"position constraint derivative v"
}

Definition at line 61 of file strnode.cc.

Referenced by StructDispNode::DescribeEq().

const char* sn_dof[]
static
Initial value:
= {
sdn_dof[0],
"incremental rotation parameter g",
sdn_dof[1],
"momenta moment G"
}
static const char * sdn_dof[]
Definition: strnode.cc:49

Definition at line 66 of file strnode.cc.

Referenced by DynamicStructNode::DescribeDof().

const char* sn_eq[]
static
Initial value:
= {
sdn_eq[0],
"momenta moment definition G",
sdn_eq[1],
"moment equilibrium M"
}
static const char * sdn_eq[]
Definition: strnode.cc:53

Definition at line 72 of file strnode.cc.

Referenced by StructNode::DescribeEq(), and DynamicStructNode::DescribeEq().

const char* sn_initial_dof[]
static
Initial value:
= {
"incremental rotation parameter g",
"angular velocity w"
}
static const char * sdn_initial_dof[]
Definition: strnode.cc:57

Definition at line 84 of file strnode.cc.

Referenced by StructNode::DescribeDof(), and ModalNode::DescribeDof().

const char* sn_initial_eq[]
static
Initial value:
= {
"orientation constraint g",
"orientation constraint derivative w"
}
static const char * sdn_initial_eq[]
Definition: strnode.cc:61

Definition at line 90 of file strnode.cc.

Referenced by StructNode::DescribeEq(), and ModalNode::DescribeEq().

const char* sn_modal_eq[]
static
Initial value:
= {
"linear velocity definition v",
"angular velocity definition w",
"force equilibrium F",
"moment equilibrium M"
}

Definition at line 78 of file strnode.cc.

Referenced by ModalNode::DescribeEq().