|
OOFEM 3.0
|
#include <beam3d.h>
Public Member Functions | |
| Beam3d (int n, Domain *d) | |
| virtual | ~Beam3d () |
| FEInterpolation * | giveInterpolation () const override |
| void | computeConsistentMassMatrix (FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL) override |
| void | computeInitialStressMatrix (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeLumpedInitialStressMatrix (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override |
| int | giveLocalCoordinateSystem (FloatMatrix &answer) override |
| void | giveInternalForcesVector (FloatArray &answer, TimeStep *, int useUpdatedGpRecord=0) override |
| void | giveEndForcesVector (FloatArray &answer, TimeStep *tStep) |
| int | computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords) override |
| int | testElementExtension (ElementExtension ext) override |
| int | computeNumberOfDofs () override |
| int | computeNumberOfGlobalDofs () override |
| void | giveDofManDofIDMask (int inode, IntArray &) const override |
| int | giveNumberOfInternalDofManagers () const override |
| DofManager * | giveInternalDofManager (int i) const override |
| void | giveInternalDofManDofIDMask (int i, IntArray &answer) const override |
| void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) override |
| void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) override |
| double | computeVolumeAround (GaussPoint *gp) override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| void | printOutputAt (FILE *file, TimeStep *tStep) override |
| void | FiberedCrossSectionInterface_computeStrainVectorInFiber (FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *slaveGp, TimeStep *tStep) override |
| Interface * | giveInterface (InterfaceType it) override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| void | initializeFrom (InputRecord &ir, int priority) override |
| void | postInitialize () override |
| Performs post initialization steps. Called after all components are created and initialized. | |
| integrationDomain | giveIntegrationDomain () const override |
| Element_Geometry_Type | giveGeometryType () const override |
| void | updateLocalNumbering (EntityRenumberingFunctor &f) override |
| FloatMatrixF< 6, 6 > | B3SSMI_getUnknownsGtoLRotationMatrix () const override |
| Evaluate transformation matrix for reciver unknowns. | |
| void | giveCompositeExportData (std::vector< ExportRegion > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep) override |
| void | drawRawGeometry (oofegGraphicContext &gc, TimeStep *tStep) override |
| void | drawDeformedGeometry (oofegGraphicContext &gc, TimeStep *tStep, UnknownType) override |
| Public Member Functions inherited from oofem::BeamBaseElement | |
| BeamBaseElement (int n, Domain *d) | |
| virtual | ~BeamBaseElement () |
| Public Member Functions inherited from oofem::StructuralElement | |
| StructuralElement (int n, Domain *d) | |
| virtual | ~StructuralElement () |
| Destructor. | |
| void | giveCharacteristicMatrix (FloatMatrix &answer, CharType, TimeStep *tStep) override |
| void | giveCharacteristicVector (FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) override |
| virtual void | computeMassMatrix (FloatMatrix &answer, TimeStep *tStep) |
| virtual void | computeLumpedMassMatrix (FloatMatrix &answer, TimeStep *tStep) |
| virtual void | giveMassMtrxIntegrationgMask (IntArray &answer) |
| void | computeStiffnessMatrix_withIRulesAsSubcells (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
| void | computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) override |
| virtual void | giveInternalForcesVector_withIRulesAsSubcells (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) |
| virtual void | computeStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| virtual void | computeResultingIPTemperatureAt (FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode) |
| virtual void | computeResultingIPEigenstrainAt (FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode) |
| void | updateBeforeNonlocalAverage (TimeStep *tStep) override |
| virtual void | giveNonlocalLocationArray (IntArray &locationArray, const UnknownNumberingScheme &us) |
| virtual void | addNonlocalStiffnessContributions (SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep) |
| int | adaptiveUpdate (TimeStep *tStep) override |
| void | updateInternalState (TimeStep *tStep) override |
| void | updateYourself (TimeStep *tStep) override |
| int | checkConsistency () override |
| void | giveInputRecord (DynamicInputRecord &input) override |
| const char * | giveClassName () const override |
| int | giveInternalStateAtNode (FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep) override |
| void | showSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) override |
| Shows sparse structure. | |
| void | showExtendedSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) override |
| Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness). | |
| void | computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) override |
| void | computeBoundarySurfaceLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) override |
| void | computeBoundaryEdgeLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) override |
| virtual void | computeEdgeNMatrix (FloatMatrix &answer, int boundaryID, const FloatArray &lcoords) |
| computes edge interpolation matrix | |
| virtual void | computeSurfaceNMatrix (FloatMatrix &answer, int boundaryID, const FloatArray &lcoords) |
| StructuralCrossSection * | giveStructuralCrossSection () |
| Helper function which returns the structural cross-section for the element. | |
| virtual void | createMaterialStatus () |
| Public Member Functions inherited from oofem::Element | |
| Element (int n, Domain *aDomain) | |
| Element (const Element &src)=delete | |
| Element & | operator= (const Element &src)=delete |
| virtual | ~Element () |
| Virtual destructor. | |
| void | giveLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
| void | giveLocationArray (IntArray &locationArray, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
| virtual int | giveNumberOfDofs () |
| virtual void | setInternalDofManager (int num, std::unique_ptr< DofManager > dm) |
| virtual double | giveCharacteristicValue (CharType type, TimeStep *tStep) |
| virtual void | computeTangentFromSurfaceLoad (FloatMatrix &answer, BoundaryLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
| virtual void | computeTangentFromEdgeLoad (FloatMatrix &answer, BoundaryLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
| const IntArray & | giveBodyLoadList () const |
| const IntArray & | giveBoundaryLoadList () const |
| void | computeVectorOf (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
| void | computeVectorOf (const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
| void | computeBoundaryVectorOf (const IntArray &bNodes, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
| void | computeVectorOf (PrimaryField &field, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
| void | computeVectorOfPrescribed (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
| void | computeVectorOfPrescribed (const IntArray &dofIDMask, ValueModeType type, TimeStep *tStep, FloatArray &answer) |
| int | computeNumberOfPrimaryMasterDofs () |
| virtual bool | giveRotationMatrix (FloatMatrix &answer) |
| virtual bool | computeDofTransformationMatrix (FloatMatrix &answer, const IntArray &nodes, bool includeInternal) |
| virtual void | giveElementDofIDMask (IntArray &answer) const |
| virtual double | computeVolumeAreaOrLength () |
| Computes the volume, area or length of the element depending on its spatial dimension. | |
| double | computeMeanSize () |
| virtual double | computeVolume () |
| virtual double | computeArea () |
| virtual IntArray | giveBoundaryEdgeNodes (int boundary, bool includeHierarchical=false) const |
| virtual IntArray | giveBoundarySurfaceNodes (int boundary, bool includeHierarchical=false) const |
| virtual IntArray | giveBoundaryNodes (int boundary) const |
| virtual std::unique_ptr< IntegrationRule > | giveBoundaryEdgeIntegrationRule (int order, int boundary) |
| virtual std::unique_ptr< IntegrationRule > | giveBoundarySurfaceIntegrationRule (int order, int boundary) |
| int | giveDofManagerNumber (int i) const |
| const IntArray & | giveDofManArray () const |
| void | addDofManager (DofManager *dMan) |
| DofManager * | giveDofManager (int i) const |
| Node * | giveNode (int i) const |
| virtual ElementSide * | giveSide (int i) const |
| virtual FEInterpolation * | giveInterpolation (DofIDItem id) const |
| virtual const FEInterpolation * | getGeometryInterpolation () const |
| virtual Material * | giveMaterial () |
| int | giveMaterialNumber () const |
| CrossSection * | giveCrossSection () |
| int | getActivityTimeFunctionNumber () |
| void | setActivityTimeFunctionNumber (int funcIndx) |
| void | setMaterial (int matIndx) |
| virtual void | setCrossSection (int csIndx) |
| virtual int | giveNumberOfDofManagers () const |
| void | setNumberOfDofManagers (int i) |
| Sets number of element dof managers. | |
| virtual int | giveNumberOfNodes () const |
| void | setDofManagers (const IntArray &dmans) |
| void | setDofManager (int id, int dm) |
| void | setBodyLoads (const IntArray &bodyLoads) |
| void | setIntegrationRules (std ::vector< std ::unique_ptr< IntegrationRule > > irlist) |
| virtual int | giveIntegrationRuleLocalCodeNumbers (IntArray &answer, IntegrationRule &ie) |
| int | giveRegionNumber () |
| virtual void | initializeYourself (TimeStep *timeStepWhenICApply) |
| virtual bool | isActivated (TimeStep *tStep) |
| virtual bool | isCast (TimeStep *tStep) |
| virtual void | initForNewStep () |
| virtual Element_Geometry_Type | giveEdgeGeometryType (int id) const |
| Returns the receiver edge geometry type. | |
| virtual Element_Geometry_Type | giveSurfaceGeometryType (int id) const |
| Returns the receiver surface geometry type. | |
| virtual int | giveSpatialDimension () |
| virtual int | giveNumberOfBoundarySides () |
| Returns number of boundaries (entities of element_dimension-1: points, edges, surfaces). | |
| virtual int | giveNumberOfEdges () const |
| virtual int | giveNumberOfSurfaces () const |
| virtual int | giveDefaultIntegrationRule () const |
| virtual IntegrationRule * | giveDefaultIntegrationRulePtr () |
| int | giveNumberOfIntegrationRules () |
| virtual IntegrationRule * | giveIntegrationRule (int i) |
| std::vector< std ::unique_ptr< IntegrationRule > > & | giveIntegrationRulesArray () |
| int | giveGlobalIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
| virtual double | giveLengthInDir (const FloatArray &normalToCrackPlane) |
| virtual double | giveCharacteristicLength (const FloatArray &normalToCrackPlane) |
| double | giveCharacteristicLengthForPlaneElements (const FloatArray &normalToCrackPlane) |
| double | giveCharacteristicLengthForAxisymmElements (const FloatArray &normalToCrackPlane) |
| virtual double | giveCharacteristicSize (GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method) |
| virtual double | giveParentElSize () const |
| virtual bool | computeLocalCoordinates (FloatArray &answer, const FloatArray &gcoords) |
| virtual void | giveLocalCoordinateSystemVector (InternalStateType isttype, FloatArray &answer) |
| virtual void | computeMidPlaneNormal (FloatArray &answer, const GaussPoint *gp) |
| virtual int | adaptiveMap (Domain *oldd, TimeStep *tStep) |
| virtual int | mapStateVariables (Domain &iOldDom, const TimeStep &iTStep) |
| virtual int | adaptiveFinish (TimeStep *tStep) |
| void | updateLocalNumbering (EntityRenumberingFunctor &f) override |
| template<class T> | |
| void | ipEvaluator (T *src, void(T ::*f)(GaussPoint *gp)) |
| Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f. | |
| template<class T, class S> | |
| void | ipEvaluator (T *src, void(T ::*f)(GaussPoint *, S &), S &_val) |
| Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f as well as additional array. | |
| virtual void | drawYourself (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawAnnotation (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawScalar (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawSpecial (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | giveLocalIntVarMaxMin (oofegGraphicContext &gc, TimeStep *tStep, double &emin, double &emax) |
| virtual int | giveInternalStateAtSide (FloatArray &answer, InternalStateType type, InternalStateMode mode, int side, TimeStep *tStep) |
| int | giveLabel () const |
| int | giveGlobalNumber () const |
| void | setGlobalNumber (int num) |
| elementParallelMode | giveParallelMode () const |
| void | setParallelMode (elementParallelMode _mode) |
| Sets parallel mode of element. | |
| virtual elementParallelMode | giveKnotSpanParallelMode (int) const |
| int | packUnknowns (DataStream &buff, TimeStep *tStep) |
| int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep) |
| int | estimatePackSize (DataStream &buff) |
| const IntArray * | givePartitionList () const |
| void | setPartitionList (IntArray &pl) |
| virtual double | predictRelativeComputationalCost () |
| virtual double | giveRelativeSelfComputationalCost () |
| virtual double | predictRelativeRedistributionCost () |
| IntArray * | giveBodyLoadArray () |
| Returns array containing load numbers of loads acting on element. | |
| IntArray * | giveBoundaryLoadArray () |
| Returns array containing load numbers of boundary loads acting on element. | |
| void | initializeFrom (InputRecord &ir, int priority) override |
| void | initializeFinish () override |
| void | postInitialize () override |
| Performs post initialization steps. | |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| void | printOutputAt (FILE *file, TimeStep *tStep) override |
| virtual const IntArray | giveLocation () |
| virtual void | recalculateCoordinates (int nodeNumber, FloatArray &coords) |
| void | setSharedEdgeID (int iedge, int globalID) |
| void | setSharedSurfaceID (int isurf, int globalID) |
| const IntArray * | giveSharedEdgeIDs () const |
| const IntArray * | giveSharedSurfaceIDs () const |
| Public Member Functions inherited from oofem::FEMComponent | |
| FEMComponent (int n, Domain *d) | |
| virtual | ~FEMComponent ()=default |
| Virtual destructor. | |
| Domain * | giveDomain () const |
| virtual void | setDomain (Domain *d) |
| int | giveNumber () const |
| void | setNumber (int num) |
| virtual void | initializeFrom (InputRecord &ir) |
| virtual void | printYourself () |
| Prints receiver state on stdout. Useful for debugging. | |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). | |
| Public Member Functions inherited from oofem::FiberedCrossSectionInterface | |
| FiberedCrossSectionInterface () | |
| Public Member Functions inherited from oofem::Interface | |
| Interface () | |
| Constructor. | |
| virtual | ~Interface () |
| Public Member Functions inherited from oofem::Beam3dSubsoilMaterialInterface | |
| Beam3dSubsoilMaterialInterface () | |
| Constructor. | |
| virtual | ~Beam3dSubsoilMaterialInterface () |
| Public Member Functions inherited from oofem::VTKXMLExportModuleElementInterface | |
| VTKXMLExportModuleElementInterface () | |
| virtual void | giveCompositeExportData (ExportRegion &vtkPiece, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep) |
Protected Member Functions | |
| void | computeBoundaryEdgeLoadVector (FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) override |
| int | computeLoadGToLRotationMtrx (FloatMatrix &answer) override |
| void | computeBmatrixAt (GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override |
| void | computeNmatrixAt (const FloatArray &iLocCoord, FloatMatrix &) override |
| bool | computeGtoLRotationMatrix (FloatMatrix &answer) override |
| void | computeBodyLoadVectorAt (FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode) override |
| double | giveKappayCoeff (TimeStep *tStep) |
| double | giveKappazCoeff (TimeStep *tStep) |
| void | computeKappaCoeffs (TimeStep *tStep) |
| double | computeLength () override |
| void | computeGaussPoints () override |
| void | computeConstitutiveMatrixAt (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override |
| void | computeStressVector (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override |
| MaterialMode | giveMaterialMode () override |
| int | giveNumberOfIPForMassMtrxIntegration () override |
| bool | hasDofs2Condense () |
| void | computeSubSoilNMatrixAt (GaussPoint *gp, FloatMatrix &answer) |
| void | computeSubSoilStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
| void | giveInternalForcesVectorAtPoint (FloatArray &answer, TimeStep *tStep, FloatArray &coords) |
| void | computeInternalForcesFromBoundaryEdgeLoadVectorAtPoint (FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, FloatArray &pointCoords, double ds, bool global) |
| void | computeInternalForcesFromBodyLoadVectorAtPoint (FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode, FloatArray &pointCoords, double ds) |
| virtual int | computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords, const FloatArray &point) |
| Protected Member Functions inherited from oofem::BeamBaseElement | |
| virtual void | computeLocalForceLoadVector (FloatArray &answer, TimeStep *tStep, ValueModeType mode) |
| virtual void | computePointLoadVectorAt (FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global=true) |
| virtual void | giveEdgeDofMapping (IntArray &answer, int iEdge) const |
| virtual void | giveSurfaceDofMapping (IntArray &answer, int iSurf) const |
| virtual double | computeEdgeVolumeAround (GaussPoint *gp, int iEdge) |
| virtual double | computeSurfaceVolumeAround (GaussPoint *gp, int iSurf) |
| virtual int | computeLoadLEToLRotationMatrix (FloatMatrix &answer, int iEdge, GaussPoint *gp) |
| virtual int | computeLoadLSToLRotationMatrix (FloatMatrix &answer, int iSurf, GaussPoint *gp) |
| void | condense (FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what) |
| virtual void | setupIRForMassMtrxIntegration (IntegrationRule &iRule) |
Protected Attributes | |
| double | kappay |
| double | kappaz |
| double | length |
| int | referenceNode |
| FloatArray | yaxis |
| FloatArray | zaxis |
| double | referenceAngle = 0 |
| DofManager * | ghostNodes [2] |
| int | numberOfCondensedDofs |
| number of condensed DOFs | |
| int | subsoilMat |
| Subsoil material. | |
| Protected Attributes inherited from oofem::StructuralElement | |
| std::unique_ptr< FloatArray > | initialDisplacements |
| Initial displacement vector, describes the initial nodal displacements when element has been casted. | |
| Protected Attributes inherited from oofem::Element | |
| int | numberOfDofMans |
| Number of dofmanagers. | |
| IntArray | dofManArray |
| Array containing dofmanager numbers. | |
| int | material |
| Number of associated material. | |
| int | crossSection |
| Number of associated cross section. | |
| IntArray | bodyLoadArray |
| IntArray | boundaryLoadArray |
| std::vector< std ::unique_ptr< IntegrationRule > > | integrationRulesArray |
| FloatMatrix | elemLocalCS |
| Transformation material matrix, used in orthotropic and anisotropic materials, global->local transformation. | |
| int | activityTimeFunction |
| Element activity time function. If defined, nonzero value indicates active receiver, zero value inactive element. | |
| int | globalNumber |
| int | numberOfGaussPoints |
| elementParallelMode | parallel_mode |
| Determines the parallel mode of the element. | |
| IntArray | partitions |
| IntArray | globalEdgeIDs |
| IntArray | globalSurfaceIDs |
| Protected Attributes inherited from oofem::FEMComponent | |
| int | number |
| Component number. | |
| Domain * | domain |
| Link to domain object, useful for communicating with other FEM components. | |
Static Protected Attributes | |
| static FEI3dLineLin | interp |
| Geometry interpolator only. | |
| static ParamKey | IPK_Beam3d_dofsToCondense |
| [optional] DOFs to condense | |
| static ParamKey | IPK_Beam3d_refnode |
| [optional] Reference node for the beam | |
| static ParamKey | IPK_Beam3d_refangle |
| [optional] Reference angle for the beam | |
| static ParamKey | IPK_Beam3d_yaxis |
| [optional] Y axis for the beam | |
| static ParamKey | IPK_Beam3d_zaxis |
| [optional] Z axis for the beam | |
| static ParamKey | IPK_Beam3d_subsoilmat |
| [optional] Subsoil material for the beam | |
Additional Inherited Members | |
| Static Public Attributes inherited from oofem::Element | |
| static ParamKey | IPK_Element_mat |
| static ParamKey | IPK_Element_crosssect |
| static ParamKey | IPK_Element_nodes |
| static ParamKey | IPK_Element_bodyload |
| static ParamKey | IPK_Element_boundaryload |
| static ParamKey | IPK_Element_lcs |
| static ParamKey | IPK_Element_partitions |
| static ParamKey | IPK_Element_remote |
| static ParamKey | IPK_Element_activityTimeFunction |
| static ParamKey | IPK_Element_nip |
This class implements a 2-dimensional beam element with cubic lateral displacement interpolation (rotations are quadratic) and longitudial displacements are linear. This is an exact displacement approximation for beam with no nonnodal loading.
This class is not derived from liBeam3d or truss element, because it does not support any material nonlinearities (if should, stiffness must be integrated)
| oofem::Beam3d::Beam3d | ( | int | n, |
| Domain * | d ) |
Definition at line 74 of file beam3d.C.
References oofem::BeamBaseElement::BeamBaseElement(), ghostNodes, kappay, kappaz, length, numberOfCondensedDofs, oofem::Element::numberOfDofMans, oofem::Element::numberOfGaussPoints, referenceAngle, referenceNode, and subsoilMat.
Referenced by B3SSMI_getUnknownsGtoLRotationMatrix().
|
virtual |
Definition at line 89 of file beam3d.C.
References ghostNodes.
|
overridevirtual |
Evaluate transformation matrix for reciver unknowns.
Implements oofem::Beam3dSubsoilMaterialInterface.
Definition at line 345 of file beam3d.C.
References oofem::FloatMatrix::at(), oofem::FloatMatrixF< N, M >::at(), Beam3d(), and giveLocalCoordinateSystem().
Referenced by giveCompositeExportData().
|
overrideprotectedvirtual |
Computes the geometrical matrix of receiver in given integration point. The product of this matrix (assembled at given integration point) and element displacement vector is element strain vector. If lowerIndx and upperIndx parameters are specified, answer is formed only for strains within this interval. This will affects the size of answer.
| gp | Integration point for which answer is computed. |
| answer | Geometric matrix of receiver. |
| lowerIndx | If specified, answer is formed only for strain with index equal and greater than lowerIndx. This parameter has default value 1 (answer is formed from first strain). |
| upperIndx | If specified, answer is formed only for strain with index less and equal than upperIndx. This parameter has default value ALL_STRAINS (answer is formed for all strains). |
Implements oofem::StructuralElement.
Definition at line 98 of file beam3d.C.
References oofem::FloatMatrix::at(), computeLength(), oofem::FEMComponent::domain, giveKappayCoeff(), giveKappazCoeff(), oofem::GaussPoint::giveNaturalCoordinate(), kappay, kappaz, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by computeStiffnessMatrix().
|
overrideprotectedvirtual |
Computes the load vector due to body load acting on receiver, at given time step. Default implementation computes body load vector numerically as \( l=\int_V N^{\mathrm{T}} f \rho\;\mathrm{d}V \) using default integration rule. Result is transformed to global c.s.
| answer | Computed load vector due to body load |
| load | Body load which contribution is computed. |
| tStep | Time step. |
| mode | determines the response mode |
Reimplemented from oofem::StructuralElement.
Definition at line 725 of file beam3d.C.
References oofem::CS_Area, oofem::Element::giveCrossSection(), and oofem::FloatArray::times().
|
overrideprotectedvirtual |
Computes the contribution of the given load at the given boundary edge. In general, the answer should include only relevant DOFs at the edge. The related is giveBoundaryLocationArray method, which should return corresponding code numbers..
| answer | Requested contribution of load (in Global c.s.). |
| load | Load to compute contribution from. |
| edge | Edge number. |
| type | Type of the contribution. |
| mode | Determines mode of answer. |
| tStep | Time step when answer is computed. |
| global | if true (default) then contribution is in global c.s., when false then contribution is in element local c.s. |
Reimplemented from oofem::Element.
Definition at line 227 of file beam3d.C.
References oofem::FloatArray::clear(), computeGlobalCoordinates(), computeGtoLRotationMatrix(), computeLength(), computeLoadGToLRotationMtrx(), computeNmatrixAt(), oofem::Load::computeValues(), oofem::BoundaryLoad::giveCoordSystMode(), oofem::Element::giveDefaultIntegrationRulePtr(), oofem::Load::giveFormulationType(), N, OOFEM_ERROR, oofem::FloatArray::plusProduct(), and oofem::FloatArray::rotatedWith().
|
overridevirtual |
Computes consistent mass matrix of receiver using numerical integration over element volume. Mass matrix is computed as \( M = \int_V N^{\mathrm{T}} \rho N dV \), where \( N \) is displacement approximation matrix. The number of necessary integration points is determined using this->giveNumberOfIPForMassMtrxIntegration service. Only selected degrees of freedom participate in integration of mass matrix. This is described using dof mass integration mask. This mask is obtained from this->giveMassMtrxIntegrationgMask service. The nonzero mask value at i-th position indicates that i-th element DOF participates in mass matrix computation. The result is in element local coordinate system.
| answer | Mass matrix. |
| tStep | Time step. |
| mass | Total mass of receiver. |
Reimplemented from oofem::StructuralElement.
Definition at line 734 of file beam3d.C.
References oofem::FloatMatrix::at(), oofem::StructuralElement::computeConsistentMassMatrix(), computeLength(), oofem::CS_Area, oofem::CS_InertiaMomentY, oofem::CS_InertiaMomentZ, oofem::Element::giveCrossSection(), giveKappayCoeff(), giveKappazCoeff(), oofem::StructuralElement::giveStructuralCrossSection(), oofem::Element::integrationRulesArray, kappay, kappaz, oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
Computes constitutive matrix of receiver. Default implementation uses element cross section giveCharMaterialStiffnessMatrix service.
| answer | Constitutive matrix. |
| rMode | Material response mode of answer. |
| gp | Integration point for which constitutive matrix is computed. |
| tStep | Time step. |
Implements oofem::StructuralElement.
Definition at line 654 of file beam3d.C.
References oofem::StructuralElement::giveStructuralCrossSection().
Referenced by computeKappaCoeffs(), and computeStiffnessMatrix().
|
overrideprotectedvirtual |
Initializes the array of integration rules member variable. Element can have multiple integration rules for different tasks. For example structural element family class uses this feature to implement transparent support for reduced and selective integration of some strain components. Must be defined by terminator classes.
Reimplemented from oofem::Element.
Definition at line 143 of file beam3d.C.
References oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, and oofem::Element::numberOfGaussPoints.
|
overridevirtual |
Computes the global coordinates from given element's local coordinates.
| answer | Requested global coordinates. |
| lcoords | Local coordinates. |
Reimplemented from oofem::Element.
Definition at line 812 of file beam3d.C.
References oofem::FloatArray::at(), oofem::Element::giveNode(), and oofem::FloatArray::resize().
Referenced by computeBoundaryEdgeLoadVector(), and computeInternalForcesFromBoundaryEdgeLoadVectorAtPoint().
|
protectedvirtual |
Definition at line 1057 of file beam3d.C.
References oofem::FloatArray::at(), oofem::Element::giveNode(), and oofem::FloatArray::resize().
|
overrideprotectedvirtual |
Returns transformation matrix from global c.s. to local element c.s., i.e. \( r_l =T r_g \). If no transformation is necessary then answer is empty matrix and zero value is returned.
| answer | Computed rotation matrix. |
Reimplemented from oofem::Element.
Definition at line 291 of file beam3d.C.
References oofem::FloatMatrix::at(), oofem::FloatMatrix::beProductOf(), computeNumberOfGlobalDofs(), ghostNodes, giveLocalCoordinateSystem(), hasDofs2Condense(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by computeBoundaryEdgeLoadVector(), and computeInternalForcesFromBoundaryEdgeLoadVectorAtPoint().
|
overridevirtual |
Computes initial stress matrix for linear stability problem. Default implementation is not provided. Please note, that initial stress matrix depends on normal forces of element, corresponding engineering model must take this into account.
| answer | Computed initial stress matrix. |
| tStep | Time step. |
Reimplemented from oofem::StructuralElement.
Definition at line 830 of file beam3d.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computeLength(), giveEndForcesVector(), giveKappayCoeff(), giveKappazCoeff(), kappay, kappaz, oofem::min(), N, oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
|
protected |
Definition at line 1228 of file beam3d.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::BodyLoadBGT, oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), computeLoadGToLRotationMtrx(), computeNmatrixAt(), oofem::CS_Area, oofem::ForceLoadBVT, oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveCrossSection(), oofem::Element::giveDefaultIntegrationRulePtr(), oofem::FloatArray::giveSize(), OOFEM_ERROR, oofem::FloatArray::rotatedWith(), and oofem::FloatArray::times().
Referenced by giveInternalForcesVectorAtPoint().
|
protected |
Definition at line 1180 of file beam3d.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::FloatArray::clear(), computeGlobalCoordinates(), computeGtoLRotationMatrix(), computeLoadGToLRotationMtrx(), oofem::Load::computeValues(), oofem::BoundaryLoad::giveCoordSystMode(), oofem::Element::giveDefaultIntegrationRulePtr(), oofem::Load::giveFormulationType(), OOFEM_ERROR, and oofem::FloatArray::rotatedWith().
Referenced by giveInternalForcesVectorAtPoint().
|
protected |
Definition at line 440 of file beam3d.C.
References oofem::FloatMatrix::at(), computeConstitutiveMatrixAt(), computeLength(), oofem::Element::integrationRulesArray, kappay, and kappaz.
Referenced by giveKappayCoeff(), and giveKappazCoeff().
|
overrideprotectedvirtual |
Computes the length (zero for all but 1D geometries)
Reimplemented from oofem::Element.
Definition at line 421 of file beam3d.C.
References oofem::DofManager::giveCoordinate(), oofem::Element::giveNode(), and length.
Referenced by computeBmatrixAt(), computeBoundaryEdgeLoadVector(), computeConsistentMassMatrix(), computeInitialStressMatrix(), computeKappaCoeffs(), computeLumpedInitialStressMatrix(), computeNmatrixAt(), computeStiffnessMatrix(), computeSubSoilStiffnessMatrix(), and computeVolumeAround().
|
overrideprotectedvirtual |
Returns transformation matrix from global coordinate system to local element coordinate system for element load vector components. If no transformation is necessary, answer is empty matrix (default);
| answer | Transformation matrix. |
Reimplemented from oofem::StructuralElement.
Definition at line 272 of file beam3d.C.
References oofem::FloatMatrix::at(), giveLocalCoordinateSystem(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by computeBoundaryEdgeLoadVector(), computeInternalForcesFromBodyLoadVectorAtPoint(), and computeInternalForcesFromBoundaryEdgeLoadVectorAtPoint().
|
overridevirtual |
Computes lumped initial stress matrix of receiver.
| answer | Lumped initial stress matrix. |
| tStep | Time step. |
Reimplemented from oofem::StructuralElement.
Definition at line 902 of file beam3d.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computeLength(), giveEndForcesVector(), N, oofem::FloatMatrix::resize(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
Computes interpolation matrix for element unknowns. The order and meaning of unknowns is element dependent.
| iLocCoord | Local coordinates. |
| answer | Interpolation matrix evaluated at gp. |
Reimplemented from oofem::StructuralElement.
Definition at line 156 of file beam3d.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computeLength(), oofem::FEMComponent::domain, giveKappayCoeff(), giveKappazCoeff(), kappay, kappaz, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by computeBoundaryEdgeLoadVector(), computeInternalForcesFromBodyLoadVectorAtPoint(), computeSubSoilNMatrixAt(), and giveCompositeExportData().
|
inlineoverridevirtual |
Computes or simply returns total number of element's local DOFs. Must be defined by particular element.
Reimplemented from oofem::Element.
|
inlineoverridevirtual |
Computes the total number of element's global dofs. The transitions from global c.s. to nodal c.s. should NOT be included.
Reimplemented from oofem::Element.
Definition at line 130 of file beam3d.h.
Referenced by computeGtoLRotationMatrix().
|
overridevirtual |
Computes numerically stiffness matrix of receiver. Default implementation computes element stiffness using \( K=\int_v B^{\mathrm{T}} D B \mathrm{d}V \) formulae, where \( B \) is element geometric matrix and \( D \) is material stiffness matrix. No geometrical nonlinearity is taken into account. NUmerical integration procedure uses integrationRulesArray for numerical integration. Support for reduced or selected integration is implemented. The individual integration rules are assumed to correspond to different terms from which the overall matrix is assembled.
If there is one integration rule, the full integration of all coefficients is performed. Otherwise, integration is performed using following rules. Each integration rule can specify start and end strain index of strain vector components for which is valid. It is necessary to ensure that these start and end indexes, dividing geometrical matrix into blocks, are not overlapping and that each strain component is included.
Then stiffness matrix is obtained as summation of integrals \( I_{ij}=\int_v B^{\mathrm{T}}_i D_{ij} B_j \mathrm{d}V \) where \( B_i \) is i-th block of geometrical matrix and \( D_{ij} \) is corresponding constitutive sub-matrix. The geometrical matrix is obtained using computeBmatrixAt service and the constitutive matrix is obtained using computeConstitutiveMatrixAt service. The \( I_{ij} \) integral is evaluated using such integration rule, which is valid for i-th or j-th block and has smaller number of integration points.
For higher numerical performance, only one half of stiffness matrix is computed and answer is then symmetrized. Therefore, if element matrix will be generally nonsymmetric, one must specialize this method. Finally, the result is transformed into global coordinate system (or nodal coordinate system, if it is defined).
| answer | Computed stiffness matrix (symmetric). |
| rMode | Response mode. |
| tStep | Time step. |
Reimplemented from oofem::StructuralElement.
Definition at line 204 of file beam3d.C.
References oofem::FloatMatrix::add(), oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::clear(), computeBmatrixAt(), computeConstitutiveMatrixAt(), computeLength(), computeSubSoilStiffnessMatrix(), oofem::Element::giveDefaultIntegrationRulePtr(), oofem::FloatMatrix::plusProductSymmUpper(), subsoilMat, and oofem::FloatMatrix::symmetrized().
Referenced by giveInternalForcesVector().
|
overrideprotectedvirtual |
Computes the stress vector of receiver at given integration point, at time step tStep. The nature of these stresses depends on the element's type.
| answer | Stress vector. |
| strain | Strain vector. |
| gp | Integration point. |
| tStep | Time step. |
Implements oofem::StructuralElement.
Definition at line 661 of file beam3d.C.
References oofem::StructuralElement::giveStructuralCrossSection().
|
protected |
Definition at line 1032 of file beam3d.C.
References computeNmatrixAt(), and oofem::GaussPoint::giveNaturalCoordinates().
Referenced by computeSubSoilStiffnessMatrix().
|
protected |
Definition at line 1040 of file beam3d.C.
References oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::clear(), computeLength(), computeSubSoilNMatrixAt(), oofem::FEMComponent::domain, oofem::StructuralMaterial::give3dBeamSubSoilStiffMtrx(), oofem::Element::giveDefaultIntegrationRulePtr(), N, oofem::FloatMatrix::plusProductSymmUpper(), subsoilMat, and oofem::FloatMatrix::symmetrized().
Referenced by computeStiffnessMatrix(), giveInternalForcesVector(), and giveInternalForcesVectorAtPoint().
|
overridevirtual |
Returns volume related to given integration point. Used typically in subroutines, that perform integration over element volume. Should be implemented by particular elements.
| gp | Integration point for which volume is computed. |
Reimplemented from oofem::Element.
Definition at line 363 of file beam3d.C.
References computeLength(), and oofem::GaussPoint::giveWeight().
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 1004 of file beam3d.C.
References gc, oofem::Element::giveNode(), OOFEG_DEFORMED_GEOMETRY_LAYER, and OOFEG_DEFORMED_GEOMETRY_WIDTH.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 978 of file beam3d.C.
References gc, oofem::Element::giveNode(), OOFEG_RAW_GEOMETRY_LAYER, and OOFEG_RAW_GEOMETRY_WIDTH.
|
overridevirtual |
Computes full 3d strain vector in element fiber. This function is necessary if layered cross section is specified.
| answer | Full fiber strain vector. |
| masterGpStrain | Strain vector at master gauss point. |
| slaveGp | Slave integration point representing particular fiber. |
| tStep | Time step. |
Implements oofem::FiberedCrossSectionInterface.
Definition at line 936 of file beam3d.C.
References oofem::FloatArray::at(), oofem::GaussPoint::giveNaturalCoordinate(), and oofem::FloatArray::resize().
|
inlineoverridevirtual |
Reimplemented from oofem::Element.
Definition at line 159 of file beam3d.h.
References oofem::Element::giveLocationArray().
|
inlineoverridevirtual |
Returns the location array for the boundary of the element. Only takes into account nodes in the bNodes vector.
Reimplemented from oofem::Element.
Definition at line 155 of file beam3d.h.
References oofem::Element::giveLocationArray().
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
|
overridevirtual |
Reimplemented from oofem::VTKXMLExportModuleElementInterface.
Definition at line 1274 of file beam3d.C.
References oofem::__InternalStateTypeToString(), oofem::__UnknownTypeToString(), oofem::FloatArray::at(), oofem::IntArray::at(), B3SSMI_getUnknownsGtoLRotationMatrix(), Beam3d_nSubBeams, oofem::FloatArray::beProductOf(), oofem::FloatArray::beTProductOf(), computeNmatrixAt(), oofem::Element::computeVectorOf(), oofem::DofManager::giveCoordinate(), giveInternalForcesVectorAtPoint(), oofem::Element::giveLocalCoordinateSystemVector(), oofem::Element::giveNode(), and oofem::IntArray::giveSize().
|
overridevirtual |
Returns dofmanager dof mask for node. This mask defines the dofs which are used by element in node. Mask influences the code number ordering for particular node. Code numbers are ordered according to node order and dofs belonging to particular node are ordered according to this mask. If element requests dofs using node mask which are not in node then error is generated. This masking allows node to be shared by different elements with different dofs in same node. Elements local code numbers are extracted from node using this mask. Must be defined by particular element.
| inode | Mask is computed for local dofmanager with inode number. |
| answer | Mask for node. |
Reimplemented from oofem::Element.
| void oofem::Beam3d::giveEndForcesVector | ( | FloatArray & | answer, |
| TimeStep * | tStep ) |
Definition at line 668 of file beam3d.C.
References oofem::BeamBaseElement::computeLocalForceLoadVector(), giveInternalForcesVector(), oofem::FloatArray::giveSize(), and oofem::FloatArray::subtract().
Referenced by computeInitialStressMatrix(), computeLumpedInitialStressMatrix(), giveInternalForcesVectorAtPoint(), and printOutputAt().
|
inlineoverridevirtual |
Returns the element geometry type. This information is assumed to be of general interest, but it is required only for some specialized tasks.
Implements oofem::Element.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 179 of file beam3d.h.
References _IFT_Beam3d_Name.
|
inlineoverridevirtual |
Reimplemented from oofem::Element.
Definition at line 183 of file beam3d.h.
References oofem::_Line.
|
overridevirtual |
Interface requesting service.
Reimplemented from oofem::FEMComponent.
Definition at line 953 of file beam3d.C.
References oofem::Beam3dSubsoilMaterialInterface::Beam3dSubsoilMaterialInterface(), oofem::Beam3dSubsoilMaterialInterfaceType, oofem::FiberedCrossSectionInterface::FiberedCrossSectionInterface(), oofem::FiberedCrossSectionInterfaceType, oofem::VTKXMLExportModuleElementInterface::VTKXMLExportModuleElementInterface(), and oofem::VTKXMLExportModuleElementInterfaceType.
|
inlineoverridevirtual |
Returns i-th internal element dof manager of the receiver
| i | Internal number of DOF. |
Reimplemented from oofem::Element.
Definition at line 133 of file beam3d.h.
References ghostNodes, oofem::FEMComponent::number, and OOFEM_ERROR.
|
inlineoverridevirtual |
Returns internal dofmanager dof mask for node. This mask defines the dofs which are used by element in node. Mask influences the code number ordering for particular node. Code numbers are ordered according to node order and dofs belonging to particular node are ordered according to this mask. If element requests dofs using node mask which are not in node then error is generated. This masking allows node to be shared by different elements with different dofs in same node. Elements local code numbers are extracted from node using this mask. Must be defined by particular element.
| inode | Mask is computed for local dofmanager with inode number. |
| answer | mask for node. |
Reimplemented from oofem::Element.
Definition at line 142 of file beam3d.h.
References ghostNodes, oofem::FEMComponent::number, and OOFEM_ERROR.
|
overridevirtual |
Returns equivalent nodal forces vectors. Useful for nonlinear analysis. Default implementation computes result as \( F=\int_v B^{\mathrm{T}} \sigma \mathrm{d}V \), where \( \sigma \) is the real element stress vector obtained using computeStressVector service (if useUpdatedGpRecord=0) or (if useUpdatedGpRecord=1) from integration point status. The geometric matrix is obtained using computeBmatrixAt service. Integration is performed using default integration rule, which should produce always valid results, assuming that strains used for computation of stresses are valid.
| answer | Internal nodal forces vector. |
| tStep | Time step. |
| useUpdatedGpRecord | If equal to zero, the stresses in integration points are computed (slow but safe), else if nonzero the stresses are taken directly from integration point status (should be derived from StructuralMaterialStatus) (fast, but engineering model must ensure valid status data in each integration point). |
Reimplemented from oofem::StructuralElement.
Definition at line 628 of file beam3d.C.
References oofem::FloatArray::add(), oofem::FloatArray::beProductOf(), computeStiffnessMatrix(), computeSubSoilStiffnessMatrix(), oofem::Element::computeVectorOf(), and subsoilMat.
Referenced by giveEndForcesVector().
|
protected |
Definition at line 1074 of file beam3d.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::FloatArray::beProductOf(), oofem::FloatArray::beSubArrayOf(), oofem::Element::bodyLoadArray, oofem::BodyLoadBGT, oofem::Element::boundaryLoadArray, computeInternalForcesFromBodyLoadVectorAtPoint(), computeInternalForcesFromBoundaryEdgeLoadVectorAtPoint(), computeSubSoilStiffnessMatrix(), oofem::Element::computeVectorOf(), oofem::FEMComponent::domain, oofem::EdgeLoadBGT, oofem::EigenstrainBVT, oofem::ForceLoadBVT, oofem::BCTracker::getElementRecords(), oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveBodyLoadArray(), oofem::Element::giveBoundaryLoadArray(), oofem::DofManager::giveCoordinate(), giveEndForcesVector(), oofem::Element::giveNode(), oofem::FloatArray::giveSize(), oofem::GeneralBoundaryCondition::isImposed(), oofem::FEMComponent::number, OOFEM_ERROR, subsoilMat, oofem::FloatArray::subtract(), oofem::TemperatureBVT, and oofem::FloatArray::times().
Referenced by giveCompositeExportData().
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 95 of file beam3d.C.
References interp.
|
overridevirtual |
Returns the integration point corresponding value in full form.
| answer | Contain corresponding integration point value, zero sized if not available. |
| gp | Integration point to check. |
| type | Determines the type of internal variable. |
| tStep | Time step. |
Reimplemented from oofem::Element.
Definition at line 370 of file beam3d.C.
References oofem::FloatArray::at(), oofem::GaussPoint::giveMaterialStatus(), oofem::StructuralMaterialStatus::giveStrainVector(), oofem::StructuralMaterialStatus::giveStressVector(), and oofem::FloatArray::resize().
|
protected |
Definition at line 465 of file beam3d.C.
References computeKappaCoeffs(), and kappay.
Referenced by computeBmatrixAt(), computeConsistentMassMatrix(), computeInitialStressMatrix(), and computeNmatrixAt().
|
protected |
Definition at line 476 of file beam3d.C.
References computeKappaCoeffs(), and kappaz.
Referenced by computeBmatrixAt(), computeConsistentMassMatrix(), computeInitialStressMatrix(), and computeNmatrixAt().
|
overridevirtual |
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy. Returns a unit vectors of local coordinate system at element stored row-wise. If local system is equal to global one, set answer to empty matrix and return zero value.
Reimplemented from oofem::Element.
Definition at line 487 of file beam3d.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::beDifferenceOf(), oofem::FloatArray::beProductOf(), oofem::FloatArray::beVectorProductOf(), oofem::FloatArray::dotProduct(), oofem::DofManager::giveCoordinates(), oofem::FEMComponent::giveDomain(), oofem::Element::giveNode(), M_PI, oofem::FloatArray::normalize(), referenceAngle, referenceNode, oofem::FloatMatrix::resize(), yaxis, zaxis, and oofem::FloatMatrix::zero().
Referenced by B3SSMI_getUnknownsGtoLRotationMatrix(), computeGtoLRotationMatrix(), and computeLoadGToLRotationMtrx().
|
inlineoverrideprotectedvirtual |
Returns material mode for receiver integration points. Should be specialized.
Reimplemented from oofem::Element.
|
inlineoverridevirtual |
Reimplemented from oofem::Element.
Definition at line 132 of file beam3d.h.
References ghostNodes.
|
inlineoverrideprotectedvirtual |
Return desired number of integration points for consistent mass matrix computation, if required.
TODO this is without the jacobian and density
Reimplemented from oofem::StructuralElement.
|
inlineprotected |
Definition at line 229 of file beam3d.h.
References ghostNodes.
Referenced by computeGtoLRotationMatrix().
|
overridevirtual |
Reimplemented from oofem::FEMComponent.
Definition at line 564 of file beam3d.C.
References oofem::FEMComponent::domain, IPK_Beam3d_dofsToCondense, IPK_Beam3d_refangle, IPK_Beam3d_refnode, IPK_Beam3d_subsoilmat, IPK_Beam3d_yaxis, IPK_Beam3d_zaxis, oofem::FEMComponent::number, PM_UPDATE_PARAMETER, PM_UPDATE_TEMP_PARAMETER, referenceAngle, referenceNode, subsoilMat, yaxis, and zaxis.
|
overridevirtual |
Performs post initialization steps. Called after all components are created and initialized.
Reimplemented from oofem::FEMComponent.
Definition at line 578 of file beam3d.C.
References oofem::IntArray::at(), oofem::ParameterManager::checkIfSet(), oofem::ComponentInputException::ctElement, oofem::FEMComponent::domain, oofem::ParameterManager::getTempParam(), ghostNodes, oofem::FEMComponent::giveDomain(), oofem::IntArray::giveSize(), IPK_Beam3d_dofsToCondense, IPK_Beam3d_refangle, IPK_Beam3d_refnode, IPK_Beam3d_yaxis, IPK_Beam3d_zaxis, oofem::FEMComponent::number, numberOfCondensedDofs, OOFEM_WARNING, referenceNode, yaxis, and zaxis.
|
overridevirtual |
Prints output of receiver to stream, for given time step. This is used for output into the standard output file.
| file | File pointer to print to. |
| tStep | Time step to write for. |
Reimplemented from oofem::FEMComponent.
Definition at line 695 of file beam3d.C.
References oofem::Element::computeVectorOf(), giveEndForcesVector(), oofem::Element::giveLabel(), oofem::FEMComponent::giveNumber(), and oofem::Element::integrationRulesArray.
|
inlineoverridevirtual |
Tests if the element implements required extension. ElementExtension type defines the list of all available element extensions.
| ext | Extension to be tested. |
Reimplemented from oofem::Element.
Definition at line 124 of file beam3d.h.
References oofem::Element_EdgeLoadSupport.
|
overridevirtual |
Local renumbering support. For some tasks (parallel load balancing, for example) it is necessary to renumber the entities. The various FEM components (such as nodes or elements) typically contain links to other entities in terms of their local numbers, etc. This service allows to update these relations to reflect updated numbering. The renumbering function is passed, which is supposed to return an updated number of specified entity type based on old number.
Reimplemented from oofem::FEMComponent.
Definition at line 968 of file beam3d.C.
References oofem::ERS_DofManager, and referenceNode.
|
protected |
Definition at line 93 of file beam3d.h.
Referenced by Beam3d(), computeGtoLRotationMatrix(), giveInternalDofManager(), giveInternalDofManDofIDMask(), giveNumberOfInternalDofManagers(), hasDofs2Condense(), postInitialize(), and ~Beam3d().
|
staticprotected |
Geometry interpolator only.
Definition at line 80 of file beam3d.h.
Referenced by giveInterpolation().
|
staticprotected |
[optional] DOFs to condense
Definition at line 100 of file beam3d.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
[optional] Reference angle for the beam
Definition at line 102 of file beam3d.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
[optional] Reference node for the beam
Definition at line 101 of file beam3d.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
[optional] Subsoil material for the beam
Definition at line 105 of file beam3d.h.
Referenced by initializeFrom().
|
staticprotected |
[optional] Y axis for the beam
Definition at line 103 of file beam3d.h.
Referenced by initializeFrom(), and postInitialize().
|
staticprotected |
[optional] Z axis for the beam
Definition at line 104 of file beam3d.h.
Referenced by initializeFrom(), and postInitialize().
|
protected |
Definition at line 82 of file beam3d.h.
Referenced by Beam3d(), computeBmatrixAt(), computeConsistentMassMatrix(), computeInitialStressMatrix(), computeKappaCoeffs(), computeNmatrixAt(), and giveKappayCoeff().
|
protected |
Definition at line 82 of file beam3d.h.
Referenced by Beam3d(), computeBmatrixAt(), computeConsistentMassMatrix(), computeInitialStressMatrix(), computeKappaCoeffs(), computeNmatrixAt(), and giveKappazCoeff().
|
protected |
Definition at line 82 of file beam3d.h.
Referenced by Beam3d(), and computeLength().
|
protected |
number of condensed DOFs
Definition at line 95 of file beam3d.h.
Referenced by Beam3d(), and postInitialize().
|
protected |
Definition at line 85 of file beam3d.h.
Referenced by Beam3d(), giveLocalCoordinateSystem(), and initializeFrom().
|
protected |
Definition at line 83 of file beam3d.h.
Referenced by Beam3d(), giveLocalCoordinateSystem(), initializeFrom(), postInitialize(), and updateLocalNumbering().
|
protected |
Subsoil material.
Definition at line 98 of file beam3d.h.
Referenced by Beam3d(), computeStiffnessMatrix(), computeSubSoilStiffnessMatrix(), giveInternalForcesVector(), giveInternalForcesVectorAtPoint(), and initializeFrom().
|
protected |
Definition at line 84 of file beam3d.h.
Referenced by giveLocalCoordinateSystem(), initializeFrom(), and postInitialize().
|
protected |
Definition at line 84 of file beam3d.h.
Referenced by giveLocalCoordinateSystem(), initializeFrom(), and postInitialize().