107 double l, ksi,
kappa, c1;
113 c1 = 1. + 2. *
kappa;
118 answer.
at(1, 1) = -1. / l;
119 answer.
at(1, 4) = 1. / l;
120 answer.
at(2, 2) = ( 6. - 12. * ksi ) / ( l * l * c1 );
121 answer.
at(2, 3) = ( -2. * ( 2. +
kappa ) + 6. * ksi ) / ( l * c1 );
122 answer.
at(2, 5) = ( -6. + 12. * ksi ) / ( l * l * c1 );
123 answer.
at(2, 6) = ( -2. * ( 1. -
kappa ) + 6. * ksi ) / ( l * c1 );
124 answer.
at(3, 2) = ( -2. *
kappa ) / ( l * c1 );
125 answer.
at(3, 3) =
kappa / ( c1 );
126 answer.
at(3, 5) = 2. *
kappa / ( l * c1 );
127 answer.
at(3, 6) =
kappa / ( c1 );
132Beam2d :: computeGaussPoints()
152 double l, ksi, ksi2, ksi3,
kappa, c1;
153 TimeStep *tStep = this->
domain->giveEngngModel()->giveCurrentStep();
156 ksi = 0.5 + 0.5 * iLocCoord.
at(1);
158 c1 = 1. + 2. *
kappa;
165 answer.
at(1, 1) = 1. - ksi;
166 answer.
at(1, 4) = ksi;
167 answer.
at(2, 2) = ( c1 - 2. *
kappa * ksi - 3. * ksi2 + 2. * ksi3 ) / c1;
168 answer.
at(2, 3) = l * ( -( 1. +
kappa ) * ksi + ( 2. +
kappa ) * ksi2 - ksi3 ) / c1;
169 answer.
at(2, 5) = ( 2. *
kappa * ksi + 3. * ksi2 - 2. * ksi3 ) / c1;
170 answer.
at(2, 6) = l * (
kappa * ksi + ( 1. -
kappa ) * ksi2 - ksi3 ) / c1;
171 answer.
at(3, 2) = ( 6. * ksi - 6. * ksi2 ) / ( l * c1 );
172 answer.
at(3, 3) = ( c1 - 2. * ( 2. +
kappa ) * ksi + 3. * ksi2 ) / c1;
173 answer.
at(3, 5) = ( -6. * ksi + 6. * ksi2 ) / ( l * c1 );
174 answer.
at(3, 6) = ( -2. * ( 1. -
kappa ) * ksi + 3. * ksi2 ) / c1;
187 double dV = gp->giveWeight() * 0.5 * l;
215 answer.
resize(ndofs, ndofs);
220 answer.
at(1, 1) = cosine;
221 answer.
at(1, 2) = sine;
222 answer.
at(2, 1) = -sine;
223 answer.
at(2, 2) = cosine;
224 answer.
at(3, 3) = 1.;
225 answer.
at(4, 4) = cosine;
226 answer.
at(4, 5) = sine;
227 answer.
at(5, 4) = -sine;
228 answer.
at(5, 5) = cosine;
229 answer.
at(6, 6) = 1.;
231 for (
int i = 7; i <= ndofs; i++ ) {
232 answer.
at(i, i) = 1.0;
236 int condensedDofCounter = 0;
243 for (
int inode = 0; inode < 2; inode++ ) {
245 for (
int idof = 0; idof < 3; idof++ ) {
246 int eq = inode * 3 + idof + 1;
248 if (
ghostNodes [ inode ]->hasDofID(dofids [ idof ]) ) {
249 condensedDofCounter++;
250 l2p.
at(eq, 6 + condensedDofCounter) = 1.0;
254 l2p.
at(eq, eq) = 1.0;
280 double layerZeta, layerZCoord, top, bottom;
285 layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
289 answer.
at(1) = masterGpStrain.
at(1) + masterGpStrain.
at(2) * layerZCoord;
290 answer.
at(2) = masterGpStrain.
at(3);
295Beam2d :: giveDofManDofIDMask(
int inode,
IntArray &answer)
const
304Beam2d :: computeLength()
314 length = sqrt(dx * dx + dy * dy);
324 double xA, xB, yA, yB;
327 if (
pitch == 10. ) {
334 pitch = atan2(yB - yA, xB - xA);
351 kappa = 6. * d.
at(2, 2) / ( d.
at(3, 3) * l * l );
373 answer.
at(1, 1) = cosine;
374 answer.
at(1, 2) = sine;
375 answer.
at(2, 1) = -sine;
376 answer.
at(2, 2) = cosine;
377 answer.
at(3, 3) = 1.0;
387 BeamBaseElement :: initializeFrom(ir, priority);
393Beam2d :: postInitialize()
395 BeamBaseElement :: postInitialize();
399 IntArray val (std::get<IntArray>(*_val));
409 for (
int i = 1; i <= val.
giveSize(); i++ ) {
410 if ( val.
at(i) <= 3 ) {
428 BeamBaseElement :: giveInternalForcesVector(answer, tStep, useUpdatedGpRecord);
454 OOFEM_ERROR(
"Beam2D only has 1 edge (the midline) that supports loads. Attempted to apply load to edge %d", edge);
457 if ( type != ExternalForcesVector ) {
467 const FloatArray &lcoords = gp->giveNaturalCoordinates();
470 load->
computeValues(t, tStep, lcoords, { D_u, D_w, R_v }, mode);
473 load->
computeValues(t, tStep, coords, { D_u, D_w, R_v }, mode);
482 double dl = gp->giveWeight() * 0.5 * l;
498 BeamBaseElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
506 if ( type == IST_BeamForceMomentTensor ) {
509 }
else if ( type == IST_BeamStrainCurvatureTensor ) {
512 }
else if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ) {
515 const FloatArray &help = type == IST_ShellForceTensor ?
520 answer.
at(1) = help.
at(1);
524 answer.
at(5) = help.
at(2);
527 }
else if ( type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
528 const FloatArray &help = type == IST_ShellMomentTensor ?
537 answer.
at(6) = help.
at(3);
540 return BeamBaseElement :: giveIPValue(answer, gp, type, tStep);
546Beam2d :: printOutputAt(FILE *File,
TimeStep *tStep)
558 fprintf(File,
" local displacements ");
559 for (
auto &val : rl ) {
560 fprintf(File,
" %.4e", val);
563 fprintf(File,
"\n local end forces ");
564 for (
auto &val : Fl ) {
565 fprintf(File,
" %.4e", val);
571 iRule->printOutputAt(File, tStep);
577Beam2d :: computeConsistentMassMatrix(
FloatMatrix &answer,
TimeStep *tStep,
double &mass,
const double *ipDensity)
593 if ( ipDensity != NULL ) {
595 density = * ipDensity;
599 double c2 = ( area * density ) / ( ( 1. + 2. *
kappa ) * ( 1. + 2. *
kappa ) );
600 double c1 = ( area * density );
605 answer.
at(1, 1) = c1 * l / 3.0;
606 answer.
at(1, 4) = c1 * l / 6.0;
607 answer.
at(2, 2) = c2 * l * ( 13. / 35. + 7. *
kappa / 5. + 4. * kappa2 / 3. );
608 answer.
at(2, 3) = -c2 * l * l * ( 11. / 210. +
kappa * 11. / 60. + kappa2 / 6. );
609 answer.
at(2, 5) = c2 * l * ( 9. / 70. +
kappa * 3. / 5. + kappa2 * 2. / 3. );
610 answer.
at(2, 6) = c2 * l * l * ( 13. / 420. +
kappa * 3. / 20. + kappa2 / 6. );
611 answer.
at(3, 3) = c2 * l * l * l * ( 1. / 105. +
kappa / 30. + kappa2 / 30. );
612 answer.
at(3, 5) = -c2 * l * l * ( 13. / 420. +
kappa * 3. / 20. + kappa2 / 6. );
613 answer.
at(3, 6) = -c2 * l * l * l * ( 1. / 140. +
kappa / 30. + kappa2 / 30. );
615 answer.
at(4, 4) = c1 * l / 3.0;
616 answer.
at(5, 5) = c2 * l * ( 13. / 35. +
kappa * 7. / 5. + kappa2 * 4. / 3. );
617 answer.
at(5, 6) = c2 * l * l * ( 11. / 210. +
kappa * 11. / 60. + kappa2 / 6. );
618 answer.
at(6, 6) = c2 * l * l * l * ( 1. / 105. +
kappa / 30. + kappa2 / 30. );
622 mass = l * area * density;
642 answer.
at(2, 2) = 4. * kappa2 + 4. *
kappa + 6. / 5.;
643 answer.
at(2, 3) = -l / 10.;
644 answer.
at(2, 5) = -4. * kappa2 - 4. *
kappa - 6. / 5.;
645 answer.
at(2, 6) = -l / 10.;
646 answer.
at(3, 3) = l * l * ( kappa2 / 3. +
kappa / 3. + 2. / 15. );
647 answer.
at(3, 5) = l / 10.;
648 answer.
at(3, 6) = -l * l * ( kappa2 / 3. +
kappa / 3. + 1. / 30. );
650 answer.
at(5, 5) = 4. * kappa2 + 4. *
kappa + 6. / 5.;
651 answer.
at(5, 6) = l / 10.;
652 answer.
at(6, 6) = l * l * ( kappa2 / 3. +
kappa / 3. + 2. / 15. );
654 answer.
at(1, 1) =
min( fabs( answer.
at(2, 2) ), fabs( answer.
at(3, 3) ) ) / 1000.;
655 answer.
at(1, 4) = -answer.
at(1, 1);
656 answer.
at(4, 4) = answer.
at(1, 1);
662 N = ( -endForces.
at(1) + endForces.
at(4) ) / 2.;
682 answer.
at(2, 2) = 1.;
683 answer.
at(2, 5) =-1.;
684 answer.
at(5, 2) =-1.;
685 answer.
at(5, 5) = 1.;
689 N = ( -endForces.
at(1) + endForces.
at(4) ) / 2.;
700 if ( !
gc.testElementGraphicActivity(
this) ) {
707 EASValsSetColor(
gc.getElementColor() );
709 p [ 0 ].x = ( FPNum ) this->
giveNode(1)->giveCoordinate(1);
711 p [ 0 ].z = ( FPNum ) this->
giveNode(1)->giveCoordinate(3);
712 p [ 1 ].x = ( FPNum ) this->
giveNode(2)->giveCoordinate(1);
714 p [ 1 ].z = ( FPNum ) this->
giveNode(2)->giveCoordinate(3);
715 go = CreateLine3D(p);
716 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
717 EGAttachObject(go, ( EObjectP )
this);
718 EMAddGraphicsToModel(ESIModel(), go);
726 if ( !
gc.testElementGraphicActivity(
this) ) {
730 double defScale =
gc.getDefScale();
734 EASValsSetColor(
gc.getDeformedElementColor() );
736 p [ 0 ].x = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
738 p [ 0 ].z = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
740 p [ 1 ].x = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
742 p [ 1 ].z = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
743 go = CreateLine3D(p);
744 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
745 EMAddGraphicsToModel(ESIModel(), go);
#define REGISTER_Element(class)
double giveKappaCoeff(TimeStep *tStep)
int computeNumberOfGlobalDofs() override
DofManager * ghostNodes[2]
double computeLength() override
void giveEndForcesVector(FloatArray &answer, TimeStep *tStep)
bool computeGtoLRotationMatrix(FloatMatrix &answer) override
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) override
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &) override
static FEI2dLineLin interp_geom
static ParamKey IPK_Beam2d_dofsToCondense
[optional] DOFs to condense
int numberOfCondensedDofs
number of condensed DOFs
void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override
void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override
virtual void computeLocalForceLoadVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
BeamBaseElement(int n, Domain *d)
CoordSystType giveCoordSystMode() override
double giveCoordinate(int i) const
Node * giveNode(int i) const
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
int numberOfDofMans
Number of dofmanagers.
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
CrossSection * giveCrossSection()
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Domain * giveDomain() const
Domain * domain
Link to domain object, useful for communicating with other FEM components.
int number
Component number.
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
void rotatedWith(FloatMatrix &r, char mode)
void subtract(const FloatArray &src)
bool isNotEmpty() const
Returns true if receiver is not empty.
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
void resize(Index rows, Index cols)
*Sets size of receiver to be an empty matrix It will have zero rows and zero columns size void clear()
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
IntegrationPointStatus * giveMaterialStatus(IntegrationPointStatusIDType key=IPSID_Default)
double giveWeight()
Returns integration weight of receiver.
LayeredCrossSectionInterface()
virtual void computeValues(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, const IntArray &dofids, ValueModeType mode)
virtual FormulationType giveFormulationType()
bool checkIfSet(size_t componentIndex, size_t paramIndex)
std::optional< paramValue > getTempParam(size_t componentIndex, size_t paramIndex) const
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver's strain vector.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
@ CS_BottomZCoord
Bottom z coordinate.
@ CS_TopZCoord
Top z coordinate.
FloatArrayF< N > min(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ LayeredCrossSectionInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_DEFORMED_GEOMETRY_LAYER
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_LAYER
#define PM_UPDATE_TEMP_PARAMETER(_type, _pm, _ir, _componentnum, _paramkey, _prio)