72 double l, ksi, n1, n2, n1x, n2x;
80 n1 = 0.5 * ( 1 - ksi );
81 n2 = 0.5 * ( 1. + ksi );
85 answer.
at(1, 1) = -1. / l;
86 answer.
at(1, 7) = 1. / l;
88 answer.
at(2, 3) = n1x;
90 answer.
at(2, 9) = n2x;
91 answer.
at(2, 11) = n2;
93 answer.
at(3, 2) = n1x;
94 answer.
at(3, 6) = -n1;
95 answer.
at(3, 8) = n2x;
96 answer.
at(3, 12) = -n2;
98 answer.
at(4, 4) = -1. / l;
99 answer.
at(4, 10) = 1. / l;
101 answer.
at(5, 5) = n1x;
102 answer.
at(5, 11) = n2x;
104 answer.
at(6, 6) = n1x;
105 answer.
at(6, 12) = n2x;
110LIBeam3d :: computeGaussPoints()
145 answer.
at(1, 1) = answer.
at(2, 2) = answer.
at(3, 3) = halfMass;
146 answer.
at(7, 7) = answer.
at(8, 8) = answer.
at(9, 9) = halfMass;
157 ksi = iLocCoord.
at(1);
158 n1 = ( 1. - ksi ) * 0.5;
159 n2 = ( 1. + ksi ) * 0.5;
164 answer.
at(1, 1) = n1;
165 answer.
at(1, 7) = n2;
166 answer.
at(2, 2) = n1;
167 answer.
at(2, 8) = n2;
168 answer.
at(3, 3) = n1;
169 answer.
at(3, 9) = n2;
171 answer.
at(4, 4) = n1;
172 answer.
at(4, 10) = n2;
173 answer.
at(5, 5) = n1;
174 answer.
at(5, 11) = n2;
175 answer.
at(6, 6) = n1;
176 answer.
at(6, 12) = n2;
184 StructuralElement :: computeStiffnessMatrix(answer, rMode, tStep);
197 for (
int i = 1; i <= 3; i++ ) {
198 for (
int j = 1; j <= 3; j++ ) {
199 answer.
at(i, j) = lcs.
at(i, j);
200 answer.
at(i + 3, j + 3) = lcs.
at(i, j);
201 answer.
at(i + 6, j + 6) = lcs.
at(i, j);
202 answer.
at(i + 9, j + 9) = lcs.
at(i, j);
224 if ( type == IST_BeamForceMomentTensor ) {
227 }
else if ( type == IST_BeamStrainCurvatureTensor ) {
231 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
237LIBeam3d :: giveDofManDofIDMask(
int inode,
IntArray &answer)
const
239 answer = {D_u, D_v, D_w, R_u, R_v, R_w};
249 n1 = ( 1. - ksi ) * 0.5;
250 n2 = ( 1. + ksi ) * 0.5;
253 answer.
at(1) = n1 * this->
giveNode(1)->giveCoordinate(1) + n2 *this->
giveNode(2)->giveCoordinate(1);
254 answer.
at(2) = n1 * this->
giveNode(1)->giveCoordinate(2) + n2 *this->
giveNode(2)->giveCoordinate(2);
255 answer.
at(3) = n1 * this->
giveNode(1)->giveCoordinate(3) + n2 *this->
giveNode(2)->giveCoordinate(3);
269LIBeam3d :: computeLength()
281 length = sqrt(dx * dx + dy * dy + dz * dz);
291 StructuralElement :: initializeFrom(ir, priority);
297LIBeam3d :: postInitialize()
301 StructuralElement :: postInitialize();
309LIBeam3d :: giveEdgeDofMapping(
IntArray &answer,
int iEdge)
const
321 for (
int i = 1; i <= 12; i++ ) {
327LIBeam3d :: computeEdgeVolumeAround(
GaussPoint *gp,
int iEdge)
354 for (
int i = 1; i <= 3; i++ ) {
355 for (
int j = 1; j <= 3; j++ ) {
356 answer.
at(i, j) = lcs.
at(i, j);
357 answer.
at(3 + i, 3 + j) = lcs.
at(i, j);
382 StructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
386Node* LIBeam3d :: giveReferenceNode(
int refNode)
393 int nnodes = this->
giveDomain()->giveNumberOfDofManagers();
395 for (
int i = 1; i <= nnodes; i++ ) {
403 OOFEM_ERROR(
"Could not find the reference node. Check numbering.");
415 Node *nodeA, *nodeB, *refNode;
423 for (
int i = 1; i <= 3; i++ ) {
428 lz.beVectorProductOf(lx, help);
430 ly.beVectorProductOf(lz, lx);
433 for (
int i = 1; i <= 3; i++ ) {
434 answer.
at(1, i) = lx.at(i);
435 answer.
at(2, i) = ly.at(i);
436 answer.
at(3, i) = lz.at(i);
443LIBeam3d :: FiberedCrossSectionInterface_computeStrainVectorInFiber(
FloatArray &answer,
const FloatArray &masterGpStrain,
446 double layerYCoord, layerZCoord;
453 answer.
at(1) = masterGpStrain.
at(1) + masterGpStrain.
at(5) * layerZCoord - masterGpStrain.
at(6) * layerYCoord;
454 answer.
at(2) = masterGpStrain.
at(2);
455 answer.
at(3) = masterGpStrain.
at(3);
459LIBeam3d :: NodalAveragingRecoveryMI_computeNodalValue(
FloatArray &answer,
int node,
485 if ( !
gc.testElementGraphicActivity(
this) ) {
492 EASValsSetColor(
gc.getElementColor() );
494 p [ 0 ].x = ( FPNum ) this->
giveNode(1)->giveCoordinate(1);
495 p [ 0 ].y = ( FPNum ) this->
giveNode(1)->giveCoordinate(2);
496 p [ 0 ].z = ( FPNum ) this->
giveNode(1)->giveCoordinate(3);
497 p [ 1 ].x = ( FPNum ) this->
giveNode(2)->giveCoordinate(1);
498 p [ 1 ].y = ( FPNum ) this->
giveNode(2)->giveCoordinate(2);
499 p [ 1 ].z = ( FPNum ) this->
giveNode(2)->giveCoordinate(3);
500 go = CreateLine3D(p);
501 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
502 EGAttachObject(go, ( EObjectP )
this);
503 EMAddGraphicsToModel(ESIModel(), go);
512 if ( !
gc.testElementGraphicActivity(
this) ) {
516 double defScale =
gc.getDefScale();
520 EASValsSetColor(
gc.getDeformedElementColor() );
522 p [ 0 ].x = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
523 p [ 0 ].y = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
524 p [ 0 ].z = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
526 p [ 1 ].x = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
527 p [ 1 ].y = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
528 p [ 1 ].z = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
529 go = CreateLine3D(p);
530 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
531 EMAddGraphicsToModel(ESIModel(), go);
#define REGISTER_Element(class)
double giveCoordinate(int i) const
Node * giveNode(int i) const
int numberOfDofMans
Number of dofmanagers.
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
CrossSection * giveCrossSection()
Domain * giveDomain() const
Domain * domain
Link to domain object, useful for communicating with other FEM components.
int number
Component number.
FiberedCrossSectionInterface()
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 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.
static ParamKey IPK_LIBeam3d_refnode
double computeLength() override
int giveLocalCoordinateSystem(FloatMatrix &answer) override
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Node * giveReferenceNode(int refNode)
NodalAveragingRecoveryModelInterface()
Constructor.
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
StructuralElement(int n, Domain *d)
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.
@ NodalAveragingRecoveryModelInterfaceType
@ FiberedCrossSectionInterfaceType
@ Element_EdgeLoadSupport
Element extension for edge loads.
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_ELEMENT_ERROR_IFNOTSET(_pm, _componentnum, _paramkey)
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)