71LatticeBeam3d :: ~LatticeBeam3d()
88LatticeBeam3d :: giveLength()
105LatticeBeam3d :: computeStiffnessMatrix(
FloatMatrix &answer, MatResponseMode rMode,
126 answer.
at(1, 1) = a / l;
127 answer.
at(1, 7) = -a / l;
129 answer.
at(2, 2) = 12. * iz / pow(l, 3.);
130 answer.
at(2, 6) = 6. * iz / pow(l, 2.);
131 answer.
at(2, 8) = -12. * iz / pow(l, 3.);
132 answer.
at(2, 12) = 6. * iz / pow(l, 2.);
134 answer.
at(3, 3) = 12 * iy / pow(l, 3.);
135 answer.
at(3, 5) = -6. * iy / pow(l, 2.);
136 answer.
at(3, 9) = -12. * iy / pow(l, 3.);
137 answer.
at(3, 11) = -6. * iy / pow(l, 2.);
139 answer.
at(4, 4) = j / ( 2. * ( 1. + nu ) * l );
140 answer.
at(4, 10) = -j / ( 2. * ( 1. + nu ) * l );
142 answer.
at(5, 3) = answer.
at(3, 5);
143 answer.
at(5, 5) = 4 * iy / l;
144 answer.
at(5, 9) = 6. * iy / pow(l, 2.);
145 answer.
at(5, 11) = 2. * iy / l;
147 answer.
at(6, 2) = answer.
at(2, 6);
148 answer.
at(6, 6) = 4. * iz / l;
149 answer.
at(6, 8) = -6. * iz / pow(l, 2.);
150 answer.
at(6, 12) = 2. * iz / l;
152 answer.
at(7, 1) = answer.
at(1, 7);
153 answer.
at(7, 7) = a / l;
155 answer.
at(8, 2) = answer.
at(2, 8);
156 answer.
at(8, 6) = answer.
at(6, 8);
157 answer.
at(8, 8) = 12. * iz / pow(l, 3.);
158 answer.
at(8, 12) = -6. * iz / pow(l, 2.);
160 answer.
at(9, 3) = answer.
at(3, 9);
161 answer.
at(9, 5) = answer.
at(5, 9);
162 answer.
at(9, 9) = 12 * iy / pow(l, 3.);
163 answer.
at(9, 11) = 6 * iy / pow(l, 2.);
165 answer.
at(10, 4) = answer.
at(4, 10.);
166 answer.
at(10, 10) = j / ( 2 * ( 1 + nu ) * l );
168 answer.
at(11, 3) = answer.
at(3, 11);
169 answer.
at(11, 5) = answer.
at(5, 11);
170 answer.
at(11, 9) = answer.
at(9, 11);
171 answer.
at(11, 11) = 4 * iy / l;
173 answer.
at(12, 2) = answer.
at(2, 12);
174 answer.
at(12, 6) = answer.
at(6, 12);
175 answer.
at(12, 8) = answer.
at(8, 12);
176 answer.
at(12, 12) = 4 * iz / l;
183void LatticeBeam3d :: computeGaussPoints()
193double LatticeBeam3d :: giveArea() {
212 for ( i = 1; i <= 3; i++ ) {
213 for ( j = 1; j <= 3; j++ ) {
214 answer.
at(i, j) = lcs.
at(i, j);
215 answer.
at(i + 3, j + 3) = lcs.
at(i, j);
216 answer.
at(i + 6, j + 6) = lcs.
at(i, j);
217 answer.
at(i + 9, j + 9) = lcs.
at(i, j);
249LatticeBeam3d :: giveDofManDofIDMask(
int inode,
IntArray &answer)
const
252 D_u, D_v, D_w, R_u, R_v, R_w
260 LatticeStructuralElement :: initializeFrom(ir, priority);
265LatticeBeam3d :: postInitialize()
268 this->LatticeStructuralElement :: postInitialize();
287LatticeBeam3d :: computeGeometryProperties()
296 for (
int i = 0; i < 3; i++ ) {
306 for (
int i = 0; i < 3; i++ ) {
307 this->
normal.at(i + 1) = coordsB.
at(i + 1) - coordsA.at(i + 1);
314 for (
int i = 0; i < 3; i++ ) {
315 this->
midPoint.at(i + 1) = 0.5 * ( coordsB.
at(i + 1) + coordsA.at(i + 1) );
318 for (
int i = 0; i < 3; i++ ) {
332LatticeBeam3d :: computeCrossSectionProperties() {
336 if ( this->
normal.at(1) == 0 ) {
338 s.at(2) = this->
normal.at(3);
339 s.at(3) = -this->
normal.at(2);
340 }
else if ( this->
normal.at(2) == 0 ) {
341 s.at(1) = this->
normal.at(3);
343 s.at(3) = -this->
normal.at(1);
345 s.at(1) = this->
normal.at(2);
346 s.at(2) = -this->
normal.at(1);
359 for (
int i = 1; i <= 3; i++ ) {
375 TimeStep *tStep,
int useUpdatedGpRecord)
406 answer.
at(1) = -stress.
at(1) *
area;
407 answer.
at(7) = -answer.
at(1);
447 if ( !
gc.testElementGraphicActivity(
this) ) {
452 EASValsSetColor(
gc.getElementColor() );
455 p [ 0 ].x = ( FPNum ) this->
giveNode(1)->giveCoordinate(1);
456 p [ 0 ].y = ( FPNum ) this->
giveNode(1)->giveCoordinate(2);
457 p [ 0 ].z = ( FPNum ) this->
giveNode(1)->giveCoordinate(3);
458 p [ 1 ].x = ( FPNum ) this->
giveNode(2)->giveCoordinate(1);
459 p [ 1 ].y = ( FPNum ) this->
giveNode(2)->giveCoordinate(2);
460 p [ 1 ].z = ( FPNum ) this->
giveNode(2)->giveCoordinate(3);
462 go = CreateLine3D(p);
463 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
464 EGAttachObject(go, ( EObjectP )
this);
465 EMAddGraphicsToModel(ESIModel(), go);
473 if ( !
gc.testElementGraphicActivity(
this) ) {
477 double defScale =
gc.getDefScale();
482 EASValsSetColor(
gc.getDeformedElementColor() );
485 p [ 0 ].x = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
486 p [ 0 ].y = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
487 p [ 0 ].z = ( FPNum ) this->
giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
489 p [ 1 ].x = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
490 p [ 1 ].y = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
491 p [ 1 ].z = ( FPNum ) this->
giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
493 go = CreateLine3D(p);
494 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
495 EMAddGraphicsToModel(ESIModel(), go);
#define REGISTER_Element(class)
double giveCoordinate(int i) const
Node * giveNode(int i) const
virtual bool isActivated(TimeStep *tStep)
virtual Material * giveMaterial()
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
int numberOfDofMans
Number of dofmanagers.
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
CrossSection * giveCrossSection()
virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
Domain * giveDomain() const
int number
Component number.
void zero()
Zeroes all coefficients of receiver.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
void subtract(const FloatArray &src)
void resize(Index rows, Index cols)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
virtual double giveLength() override
FloatArray globalCentroid
FloatMatrix localCoordinateSystem
virtual int giveLocalCoordinateSystem(FloatMatrix &answer) override
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
virtual void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
virtual void computeGeometryProperties()
static ParamKey IPK_LatticeBeam3d_diameter
virtual void drawDeformedGeometry(oofegGraphicContext &, TimeStep *tStep, UnknownType) override
virtual void computeCrossSectionProperties()
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override
LatticeStructuralElement(int n, Domain *d)
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
@ OGC_eigenVectorGeometry
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)