|
OOFEM 3.0
|
#include <truss1d.h>
Public Member Functions | |
| Truss1d (int n, Domain *d) | |
| virtual | ~Truss1d () |
| FEInterpolation * | giveInterpolation () const override |
| void | computeLumpedMassMatrix (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeMassMatrix (FloatMatrix &answer, TimeStep *tStep) override |
| int | computeNumberOfDofs () override |
| void | giveDofManDofIDMask (int inode, IntArray &) const override |
| double | giveCharacteristicLength (const FloatArray &normalToCrackPlane) override |
| double | computeVolumeAround (GaussPoint *gp) override |
| void | computeStressVector (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override |
| void | computeConstitutiveMatrixAt (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override |
| void | computeConstitutiveMatrix_dPdF_At (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override |
| int | testElementExtension (ElementExtension ext) override |
| Interface * | giveInterface (InterfaceType it) override |
| void | drawRawGeometry (oofegGraphicContext &gc, TimeStep *tStep) override |
| void | drawDeformedGeometry (oofegGraphicContext &gc, TimeStep *tStep, UnknownType) override |
| void | drawScalar (oofegGraphicContext &gc, TimeStep *tStep) override |
| const char * | giveInputRecordName () const override |
| const char * | giveClassName () const override |
| MaterialMode | giveMaterialMode () override |
| Element_Geometry_Type | giveGeometryType () const override |
| void | NodalAveragingRecoveryMI_computeNodalValue (FloatArray &answer, int node, InternalStateType type, TimeStep *tStep) override |
| void | HuertaErrorEstimatorI_setupRefinedElementProblem (RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface::SetupMode sMode, TimeStep *tStep, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator::AnalysisMode aMode) override |
| void | HuertaErrorEstimatorI_computeNmatrixAt (GaussPoint *gp, FloatMatrix &answer) override |
| Public Member Functions inherited from oofem::NLStructuralElement | |
| NLStructuralElement (int n, Domain *d) | |
| virtual | ~NLStructuralElement () |
| Destructor. | |
| int | giveGeometryMode () |
| void | computeFirstPKStressVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
| void | computeCauchyStressVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
| void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override |
| void | computeInitialStressMatrix (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeStiffnessMatrix_withIRulesAsSubcells (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
| void | giveInternalForcesVector (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) override |
| void | giveInternalForcesVector_withIRulesAsSubcells (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) override |
| virtual void | computeDeformationGradientVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
| double | computeCurrentVolume (TimeStep *tStep) |
| void | initializeFrom (InputRecord &ir, int priority) override |
| void | giveInputRecord (DynamicInputRecord &input) override |
| 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 | computeConsistentMassMatrix (FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL) |
| virtual void | giveMassMtrxIntegrationgMask (IntArray &answer) |
| void | computeStiffnessMatrix_withIRulesAsSubcells (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
| virtual void | computeLumpedInitialStressMatrix (FloatMatrix &answer, TimeStep *tStep) |
| void | computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) override |
| 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 |
| 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 void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
| virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
| virtual int | giveNumberOfDofs () |
| virtual int | giveNumberOfInternalDofManagers () const |
| virtual DofManager * | giveInternalDofManager (int i) const |
| 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) |
| virtual int | computeNumberOfGlobalDofs () |
| int | computeNumberOfPrimaryMasterDofs () |
| virtual bool | computeGtoLRotationMatrix (FloatMatrix &answer) |
| virtual bool | giveRotationMatrix (FloatMatrix &answer) |
| virtual bool | computeDofTransformationMatrix (FloatMatrix &answer, const IntArray &nodes, bool includeInternal) |
| virtual void | giveInternalDofManDofIDMask (int inode, IntArray &answer) const |
| 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 double | computeLength () |
| 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 integrationDomain | giveIntegrationDomain () const |
| 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) |
| double | giveCharacteristicLengthForPlaneElements (const FloatArray &normalToCrackPlane) |
| double | giveCharacteristicLengthForAxisymmElements (const FloatArray &normalToCrackPlane) |
| virtual double | giveCharacteristicSize (GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method) |
| virtual double | giveParentElSize () const |
| virtual int | computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords) |
| virtual bool | computeLocalCoordinates (FloatArray &answer, const FloatArray &gcoords) |
| virtual int | giveLocalCoordinateSystem (FloatMatrix &answer) |
| 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 | 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 | 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::ZZNodalRecoveryModelInterface | |
| ZZNodalRecoveryModelInterface (Element *element) | |
| Constructor. | |
| virtual bool | ZZNodalRecoveryMI_computeNValProduct (FloatMatrix &answer, InternalStateType type, TimeStep *tStep) |
| virtual void | ZZNodalRecoveryMI_computeNNMatrix (FloatArray &answer, InternalStateType type) |
| Public Member Functions inherited from oofem::Interface | |
| Interface () | |
| Constructor. | |
| virtual | ~Interface () |
| Public Member Functions inherited from oofem::NodalAveragingRecoveryModelInterface | |
| NodalAveragingRecoveryModelInterface () | |
| Constructor. | |
| Public Member Functions inherited from oofem::SpatialLocalizerInterface | |
| SpatialLocalizerInterface (Element *element) | |
| virtual int | SpatialLocalizerI_containsPoint (const FloatArray &coords) |
| int | SpatialLocalizerI_BBoxContainsPoint (const FloatArray &coords) |
| virtual void | SpatialLocalizerI_giveBBox (FloatArray &bb0, FloatArray &bb1) |
| virtual double | SpatialLocalizerI_giveClosestPoint (FloatArray &lcoords, FloatArray &closest, const FloatArray &gcoords) |
| Public Member Functions inherited from oofem::ZZErrorEstimatorInterface | |
| ZZErrorEstimatorInterface (Element *element) | |
| Constructor. | |
| virtual void | ZZErrorEstimatorI_computeElementContributions (double &eNorm, double &sNorm, ZZErrorEstimator ::NormType norm, InternalStateType type, TimeStep *tStep) |
| virtual IntegrationRule * | ZZErrorEstimatorI_giveIntegrationRule () |
| virtual void | ZZErrorEstimatorI_computeLocalStress (FloatArray &answer, FloatArray &sig) |
| Public Member Functions inherited from oofem::HuertaErrorEstimatorInterface | |
| HuertaErrorEstimatorInterface () | |
| Constructor. | |
| virtual void | HuertaErrorEstimatorI_setupRefinedElementProblem (RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface ::SetupMode mode, TimeStep *tStep, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator ::AnalysisMode aMode)=0 |
Protected Member Functions | |
| void | computeBmatrixAt (GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override |
| void | computeBHmatrixAt (GaussPoint *gp, FloatMatrix &answer) override |
| void | computeNmatrixAt (const FloatArray &iLocCoord, FloatMatrix &answer) override |
| void | computeGaussPoints () override |
| Protected Member Functions inherited from oofem::NLStructuralElement | |
| int | checkConsistency () override |
| 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 | computeLoadGToLRotationMtrx (FloatMatrix &answer) |
| virtual int | computeLoadLEToLRotationMatrix (FloatMatrix &answer, int iEdge, GaussPoint *gp) |
| virtual int | computeLoadLSToLRotationMatrix (FloatMatrix &answer, int iSurf, GaussPoint *gp) |
| virtual int | giveNumberOfIPForMassMtrxIntegration () |
| void | condense (FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what) |
| virtual void | setupIRForMassMtrxIntegration (IntegrationRule &iRule) |
| virtual void | computeBodyLoadVectorAt (FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode) |
| Protected Member Functions inherited from oofem::HuertaErrorEstimatorInterface | |
| void | setupRefinedElementProblem1D (Element *element, RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface ::SetupMode mode, TimeStep *tStep, int nodes, FloatArray *corner, FloatArray &midNode, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator ::AnalysisMode aMode, const char *edgetype) |
| void | setupRefinedElementProblem2D (Element *element, RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface ::SetupMode mode, TimeStep *tStep, int nodes, FloatArray *corner, FloatArray *midSide, FloatArray &midNode, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator ::AnalysisMode aMode, const char *quadtype) |
| void | setupRefinedElementProblem3D (Element *element, RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface ::SetupMode mode, TimeStep *tStep, int nodes, FloatArray *corner, FloatArray *midSide, FloatArray *midFace, FloatArray &midNode, int &localNodeId, int &localElemId, int &localBcId, int hexaSideNode[1][3], int hexaFaceNode[1][3], IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator ::AnalysisMode aMode, const char *hexatype) |
Static Protected Attributes | |
| static FEI1dLin | interp |
| Static Protected Attributes inherited from oofem::NLStructuralElement | |
| static ParamKey | IPK_NLStructuralElement_nlgeoflag |
Additional Inherited Members | |
| Public Types inherited from oofem::HuertaErrorEstimatorInterface | |
| enum | SetupMode { CountMode = 0 , NodeMode = 1 , ElemMode = 2 , BCMode = 3 } |
| Mode for problem setup. More... | |
| 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 |
| Protected Attributes inherited from oofem::NLStructuralElement | |
| int | nlGeometry =0 |
| Flag indicating if geometrical nonlinearities apply. | |
| 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. | |
This class implements a two-node truss bar element for one-dimensional analysis.
| oofem::Truss1d::Truss1d | ( | int | n, |
| Domain * | d ) |
Definition at line 58 of file truss1d.C.
References oofem::HuertaErrorEstimatorInterface::HuertaErrorEstimatorInterface(), oofem::NLStructuralElement::NLStructuralElement(), oofem::NodalAveragingRecoveryModelInterface::NodalAveragingRecoveryModelInterface(), oofem::Element::numberOfDofMans, oofem::SpatialLocalizerInterface::SpatialLocalizerInterface(), oofem::ZZErrorEstimatorInterface::ZZErrorEstimatorInterface(), and oofem::ZZNodalRecoveryModelInterface::ZZNodalRecoveryModelInterface().
Referenced by oofem::Truss1dGradDamage::Truss1dGradDamage().
|
overrideprotectedvirtual |
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displacement gradient stored by columns. The components of this matrix are derivatives of the shape functions, but they are arranged in a somewhat different way from the usual B matrix.
| gp | Integration point. |
| answer | BF matrix at this point. |
Reimplemented from oofem::NLStructuralElement.
Definition at line 111 of file truss1d.C.
References computeBmatrixAt().
|
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 truss1d.C.
References oofem::FloatMatrix::beTranspositionOf(), oofem::GaussPoint::giveNaturalCoordinates(), and interp.
Referenced by computeBHmatrixAt().
|
overridevirtual |
Computes large strain 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::NLStructuralElement.
Definition at line 159 of file truss1d.C.
References oofem::StructuralCrossSection::giveStiffnessMatrix_dPdF_1d(), and oofem::StructuralElement::giveStructuralCrossSection().
|
overridevirtual |
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 153 of file truss1d.C.
References oofem::StructuralCrossSection::giveStiffnessMatrix_1d(), and oofem::StructuralElement::giveStructuralCrossSection().
|
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 69 of file truss1d.C.
References oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, and oofem::CrossSection::setupIntegrationPoints().
|
overridevirtual |
Computes lumped mass matrix of receiver. Default implementation returns lumped consistent mass matrix. Then returns lumped mass transformed into nodal coordinate system. The lumping procedure zeroes all off-diagonal members and zeroes also all diagonal members corresponding to non-displacement DOFs. Such diagonal matrix is then rescaled, to preserve the element mass. Requires the computeNmatrixAt and giveMassMtrxIntegrationgMask services to be implemented.
| answer | Lumped mass matrix. |
| tStep | Time step. |
Reimplemented from oofem::StructuralElement.
Definition at line 83 of file truss1d.C.
References oofem::FloatMatrix::at(), oofem::Element::computeLength(), oofem::CS_Area, oofem::CrossSection::give(), oofem::Element::giveCrossSection(), oofem::StructuralElement::giveStructuralCrossSection(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by computeMassMatrix().
|
inlineoverridevirtual |
Computes mass matrix of receiver. Default implementation returns consistent mass matrix and uses numerical integration. Returns result of this->computeConsistentMassMatrix service, transformed into nodal coordinate system. Requires the computeNmatrixAt and giveMassMtrxIntegrationgMask services to be implemented.
| answer | Mass matrix. |
| tStep | Time step. |
Reimplemented from oofem::StructuralElement.
Definition at line 71 of file truss1d.h.
References computeLumpedMassMatrix().
|
overrideprotectedvirtual |
Reimplemented from oofem::StructuralElement.
Definition at line 123 of file truss1d.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), interp, and oofem::FloatMatrix::resize().
Referenced by HuertaErrorEstimatorI_computeNmatrixAt().
|
inlineoverridevirtual |
Computes or simply returns total number of element's local DOFs. Must be defined by particular element.
Reimplemented from oofem::Element.
Reimplemented in oofem::Truss1dGradDamage.
|
overridevirtual |
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 147 of file truss1d.C.
References oofem::StructuralCrossSection::giveRealStress_1d(), and oofem::StructuralElement::giveStructuralCrossSection().
|
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 137 of file truss1d.C.
References oofem::CS_Area, oofem::CrossSection::give(), oofem::Element::giveCrossSection(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::GaussPoint::giveWeight(), and interp.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 247 of file truss1d.C.
References gc, oofem::Element::giveNode(), oofem::Node::giveUpdatedCoordinate(), OOFEG_DEFORMED_GEOMETRY_LAYER, and OOFEG_DEFORMED_GEOMETRY_WIDTH.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 222 of file truss1d.C.
References gc, oofem::DofManager::giveCoordinate(), oofem::Element::giveNode(), OOFEG_RAW_GEOMETRY_LAYER, and OOFEG_RAW_GEOMETRY_WIDTH.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 273 of file truss1d.C.
References oofem::FloatArray::at(), gc, oofem::DofManager::giveCoordinate(), oofem::StructuralElement::giveInternalStateAtNode(), oofem::StructuralElement::giveIPValue(), oofem::Element::giveNode(), oofem::Node::giveUpdatedCoordinate(), oofem::Element::integrationRulesArray, oofem::ISM_local, oofem::ISM_recovered, OOFEG_DEFORMED_GEOMETRY_WIDTH, OOFEG_VARPLOT_PATTERN_LAYER, oofem::SA_COLORZPROFILE, oofem::SA_ISO_LINE, oofem::SA_ISO_SURF, and oofem::SA_ZPROFILE.
|
inlineoverridevirtual |
Returns the size of element in the given direction, in some cases adjusted (e.g. if the direction is perpendicular to a planar element). Required by material models relying on the crack-band approach to achieve objectivity with respect to the mesh size.
| normalToCrackPlane | Normal to the expected crack band. |
Reimplemented from oofem::Element.
Definition at line 78 of file truss1d.h.
References oofem::Element::computeLength().
|
inlineoverridevirtual |
Reimplemented from oofem::NLStructuralElement.
Reimplemented in oofem::Truss1dGradDamage.
|
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.
Reimplemented in oofem::Truss1dGradDamage.
|
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.
Reimplemented in oofem::Truss1dGradDamage.
Definition at line 98 of file truss1d.h.
References _IFT_Truss1d_Name.
|
overridevirtual |
Interface requesting service.
Reimplemented from oofem::FEMComponent.
Definition at line 383 of file truss1d.C.
References oofem::HuertaErrorEstimatorInterface::HuertaErrorEstimatorInterface(), oofem::HuertaErrorEstimatorInterfaceType, oofem::NodalAveragingRecoveryModelInterface::NodalAveragingRecoveryModelInterface(), oofem::NodalAveragingRecoveryModelInterfaceType, oofem::SpatialLocalizerInterface::SpatialLocalizerInterface(), oofem::SpatialLocalizerInterfaceType, oofem::ZZErrorEstimatorInterface::ZZErrorEstimatorInterface(), oofem::ZZErrorEstimatorInterfaceType, oofem::ZZNodalRecoveryModelInterface::ZZNodalRecoveryModelInterface(), and oofem::ZZNodalRecoveryModelInterfaceType.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 80 of file truss1d.C.
References interp.
|
inlineoverridevirtual |
Returns material mode for receiver integration points. Should be specialized.
Reimplemented from oofem::Element.
Reimplemented in oofem::Truss1dGradDamage.
|
overridevirtual |
Implements oofem::HuertaErrorEstimatorInterface.
Definition at line 215 of file truss1d.C.
References computeNmatrixAt(), and oofem::GaussPoint::giveSubPatchCoordinates().
|
override |
Definition at line 174 of file truss1d.C.
References oofem::FloatArray::at(), oofem::HuertaErrorEstimatorInterface::BCMode, oofem::DofManager::giveCoordinates(), oofem::Element::giveNode(), oofem::HuertaErrorEstimator::HEE_linear, oofem::HuertaErrorEstimatorInterface::NodeMode, oofem::FloatArray::resize(), and oofem::HuertaErrorEstimatorInterface::setupRefinedElementProblem1D().
|
overridevirtual |
Computes the element value in given node.
| answer | Contains the result. |
| node | Element node number. |
| type | Determines the type of internal variable to be recovered. |
| tStep | Time step. |
Implements oofem::NodalAveragingRecoveryModelInterface.
Definition at line 402 of file truss1d.C.
References oofem::StructuralElement::giveIPValue(), 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.
|
staticprotected |
Definition at line 62 of file truss1d.h.
Referenced by oofem::Truss1dGradDamage::computeBdMatrixAt(), computeBmatrixAt(), oofem::Truss1dGradDamage::computeField(), oofem::Truss1dGradDamage::computeNdMatrixAt(), computeNmatrixAt(), computeVolumeAround(), and giveInterpolation().