|
OOFEM 3.0
|
#include <AbaqusUserElement.h>
Public Member Functions | |
| AbaqusUserElement (int n, Domain *d) | |
| Constructor. | |
| virtual | ~AbaqusUserElement () |
| Destructor. | |
| void | computeConsistentMassMatrix (FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL) override |
| void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) override |
| void | giveInternalForcesVector (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) override |
| virtual void | giveInternalForcesVector (FloatArray &answer, TimeStep *tStep, FloatArray &U, FloatMatrix &DU, int useUpdatedGpRecord) |
| int | computeNumberOfDofs () override |
| void | giveDofManDofIDMask (int inode, IntArray &answer) const override |
| void | computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) override |
| void | updateYourself (TimeStep *tStep) override |
| void | updateInternalState (TimeStep *tStep) override |
| bool | hasTangent () const |
| virtual void | letTempTangentBe (FloatMatrix &src) |
| virtual FloatMatrix & | letTempRhsBe (FloatMatrix &src) |
| virtual FloatArray & | letTempSvarsBe (FloatArray &src) |
| virtual const FloatArray & | giveStateVector () const |
| virtual const FloatArray & | giveTempStateVector () const |
| virtual const FloatMatrix & | giveTempTangent () |
| void | initializeFrom (InputRecord &ir, int priority) override |
| void | giveInputRecord (DynamicInputRecord &input) override |
| void | postInitialize () override |
| Performs post initialization steps. Called after all components are created and initialized. | |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| integrationDomain | giveIntegrationDomain () const override |
| Element_Geometry_Type | giveGeometryType () const 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 | 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) |
| virtual void | computeInitialStressMatrix (FloatMatrix &answer, TimeStep *tStep) |
| virtual void | computeLumpedInitialStressMatrix (FloatMatrix &answer, 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) |
| virtual void | computeNmatrixAt (const FloatArray &iLocCoord, FloatMatrix &answer) |
| 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 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 | computeVolumeAround (GaussPoint *gp) |
| 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 () 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 MaterialMode | giveMaterialMode () |
| 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 () |
| virtual int | testElementExtension (ElementExtension ext) |
| 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 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 | drawRawGeometry (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawDeformedGeometry (oofegGraphicContext &gc, TimeStep *tStep, UnknownType) |
| 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. | |
| virtual Interface * | giveInterface (InterfaceType t) |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). | |
Protected Member Functions | |
| void | computeStressVector (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override |
| void | computeConstitutiveMatrixAt (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) override |
| void | computeBmatrixAt (GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) 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) |
| virtual void | computeGaussPoints () |
Private Attributes | |
| void * | uelobj |
| Dynamically loaded uel. | |
| int | nCoords |
| Coord transf. | |
| IntArray | dofs |
| FloatMatrix | coords |
| int | numSvars |
| Number of status variables. | |
| FloatArray | props |
| Element properties. | |
| FloatArray | svars |
| Status variables. | |
| FloatArray | tempSvars |
| int | jtype |
| Element type. | |
| int | nrhs = 2 |
| Element rhs. | |
| FloatMatrix | rhs |
| FloatMatrix | tempRHS |
| FloatMatrix | amatrx |
| Element amatrx. | |
| FloatMatrix | tempAmatrx |
| int | ndofel = 0 |
| ndofel | |
| int | mcrd = 2 |
| mcrd | |
| int | mlvarx = 1 |
| mlvarx | |
| FloatArray | U |
| Inputs to element routines. Velocity and Acceleration currently ignored. | |
| FloatArray | V |
| FloatArray | A |
| FloatMatrix | DU |
| IntArray | lFlags |
| LFlags. | |
| int | kinc = 1 |
| kinc, kstep | |
| int | kstep = 1 |
| int | npredef = 1 |
| predef | |
| FloatArray | predef |
| FloatArray | energy |
| energy | |
| int | ndLoad = 0 |
| loads | |
| int | mdLoad = 0 |
| FloatMatrix | adlmag |
| FloatMatrix | ddlmag |
| int * | jdltype = nullptr |
| double | params [3] |
| params | |
| IntArray | jprops |
| jprops | |
| bool | hasTangentFlag |
| Keeps track of whether the tangent has been obtained already. | |
| void(* | uel )(double *rhs, double *amatrx, double *svars, double energy[8], int *ndofel, int *nrhs, int *nsvars, double *props, int *nprops, double *coords, int *mcrd, int *nnode, double *u, double *du, double *v, double *a, int *jtype, double time[2], double *dtime, int *kstep, int *kinc, int *jelem, double params[3], int *ndload, int *jdltyp, double *adlmag, double *predef, int *npredef, int *lflags, int *mvarx, double *ddlmag, int *mdload, double *pnewdt, int *jprops, int *njprop, double *period) |
| Pointer to the dynamically loaded uel-function (translated to C). | |
| std::string | filename |
| File containing the uel function. | |
Static Private Attributes | |
| static ParamKey | IPK_AbaqusUserElement_userElement |
| static ParamKey | IPK_AbaqusUserElement_numcoords |
| static ParamKey | IPK_AbaqusUserElement_dofs |
| static ParamKey | IPK_AbaqusUserElement_numsvars |
| static ParamKey | IPK_AbaqusUserElement_properties |
| static ParamKey | IPK_AbaqusUserElement_type |
| static ParamKey | IPK_AbaqusUserElement_name |
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 |
| 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. | |
UEL interface from Abaqus user elements.
The function prototype for UEL is: SUBROUTINE UEL(RHS,AMATRX,SVARS,ENERGY,NDOFEL,NRHS,NSVARS, PROPS,NPROPS,COORDS,MCRD,NNODE,U,DU,V,A,JTYPE,TIME, DTIME,KSTEP,KINC,JELEM,PARAMS,NDLOAD,JDLTYP,ADLMAG, PREDEF,NPREDF,LFLAGS,MLVARX,DDLMAG,MDLOAD,PNEWDT, JPROPS,NJPROP,PERIOD)
DIMENSION RHS(MLVARX,*),AMATRX(NDOFEL,NDOFEL), SVARS(NSVARS),ENERGY(8),PROPS(*),COORDS(MCRD,NNODE), U(NDOFEL),DU(MLVARX,*),V(NDOFEL),A(NDOFEL),TIME(2), PARAMS(3),JDLTYP(MDLOAD,*),ADLMAG(MDLOAD,*), DDLMAG(MDLOAD,*),PREDEF(2,NPREDF,NNODE),LFLAGS(*), JPROPS(*)
Definition at line 76 of file AbaqusUserElement.h.
| oofem::AbaqusUserElement::AbaqusUserElement | ( | int | n, |
| Domain * | d ) |
Constructor.
Definition at line 63 of file AbaqusUserElement.C.
References hasTangentFlag, oofem::StructuralElement::StructuralElement(), uel, and uelobj.
|
virtual |
|
inlineoverrideprotectedvirtual |
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 228 of file AbaqusUserElement.h.
References ALL_STRAINS, and OOFEM_ERROR.
|
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 306 of file AbaqusUserElement.C.
References ndofel, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
inlineoverrideprotectedvirtual |
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 224 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
Computes the unknown vector interpolated at the specified local coordinates. Used for exporting data and mapping fields.
| mode | Mode (total, increment, etc) of the output |
| tStep | Time step to evaluate at |
| lcoords | Local coordinates to evaluate at |
| answer | Results |
Reimplemented from oofem::Element.
Definition at line 182 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
Computes or simply returns total number of element's local DOFs. Must be defined by particular element.
Reimplemented from oofem::Element.
Definition at line 180 of file AbaqusUserElement.h.
|
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 187 of file AbaqusUserElement.C.
References DU, giveInternalForcesVector(), giveTempTangent(), hasTangent(), and U.
|
inlineoverrideprotectedvirtual |
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 221 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 215 of file AbaqusUserElement.h.
|
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.
Definition at line 182 of file AbaqusUserElement.C.
References dofs.
|
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.
Definition at line 218 of file AbaqusUserElement.h.
|
overridevirtual |
Setups the input record string of receiver.
| input | Dynamic input record to be filled by receiver. |
Reimplemented from oofem::FEMComponent.
Definition at line 170 of file AbaqusUserElement.C.
References coords, dofs, filename, oofem::StructuralElement::giveInputRecord(), IPK_AbaqusUserElement_dofs, IPK_AbaqusUserElement_numcoords, IPK_AbaqusUserElement_numsvars, IPK_AbaqusUserElement_properties, IPK_AbaqusUserElement_type, IPK_AbaqusUserElement_userElement, jtype, numSvars, props, and oofem::DynamicInputRecord::setField().
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 216 of file AbaqusUserElement.h.
References _IFT_AbaqusUserElement_Name.
|
inlineoverridevirtual |
Returns integration domain for receiver, used to initialize integration point over receiver volume. Default behavior is taken from the default interpolation.
Reimplemented from oofem::Element.
Definition at line 217 of file AbaqusUserElement.h.
References oofem::_UnknownIntegrationDomain.
|
virtual |
Definition at line 228 of file AbaqusUserElement.C.
References A, adlmag, coords, oofem::FloatMatrix::copyColumn(), ddlmag, DU, energy, oofem::FloatArray::givePointer(), oofem::FloatMatrix::givePointer(), giveStateVector(), oofem::TimeStep::giveTargetTime(), oofem::TimeStep::giveTimeIncrement(), jdltype, jprops, jtype, kinc, kstep, letTempRhsBe(), letTempSvarsBe(), letTempTangentBe(), lFlags, mcrd, mdLoad, mlvarx, ndLoad, ndofel, oofem::FloatMatrix::negated(), npredef, nrhs, numSvars, params, predef, props, rhs, U, uel, and V.
|
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 214 of file AbaqusUserElement.C.
References oofem::Element::computeVectorOf(), DU, oofem::Domain::giveClassName(), oofem::FEMComponent::giveDomain(), giveInternalForcesVector(), and U.
Referenced by computeStiffnessMatrix(), giveInternalForcesVector(), and updateInternalState().
|
inlineoverrideprotectedvirtual |
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 232 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
inlinevirtual |
Definition at line 200 of file AbaqusUserElement.h.
References svars.
Referenced by giveInternalForcesVector().
|
inlinevirtual |
Definition at line 203 of file AbaqusUserElement.h.
References tempSvars.
|
inlinevirtual |
Definition at line 206 of file AbaqusUserElement.h.
References tempAmatrx.
Referenced by computeStiffnessMatrix().
|
inline |
Definition at line 187 of file AbaqusUserElement.h.
References hasTangentFlag.
Referenced by computeStiffnessMatrix().
|
overridevirtual |
Reimplemented from oofem::FEMComponent.
Definition at line 81 of file AbaqusUserElement.C.
References dofs, oofem::Domain::elementPPM, filename, oofem::FEMComponent::giveDomain(), oofem::Element::initializeFrom(), IPK_AbaqusUserElement_dofs, IPK_AbaqusUserElement_numcoords, IPK_AbaqusUserElement_numsvars, IPK_AbaqusUserElement_properties, IPK_AbaqusUserElement_type, IPK_AbaqusUserElement_userElement, jtype, nCoords, oofem::FEMComponent::number, numSvars, PM_UPDATE_PARAMETER, and props.
|
inlinevirtual |
Definition at line 194 of file AbaqusUserElement.h.
References tempRHS.
Referenced by giveInternalForcesVector().
|
inlinevirtual |
Definition at line 197 of file AbaqusUserElement.h.
References tempSvars.
Referenced by giveInternalForcesVector().
|
inlinevirtual |
Definition at line 190 of file AbaqusUserElement.h.
References hasTangentFlag, and tempAmatrx.
Referenced by giveInternalForcesVector().
|
overridevirtual |
Performs post initialization steps. Called after all components are created and initialized.
Reimplemented from oofem::FEMComponent.
Definition at line 95 of file AbaqusUserElement.C.
References _IFT_AbaqusUserElement_name, A, amatrx, coords, oofem::Element::dofManArray, DU, oofem::Domain::elementPPM, energy, filename, oofem::DofManager::giveCoordinate(), oofem::FEMComponent::giveDomain(), oofem::Element::giveNode(), IPK_AbaqusUserElement_numcoords, IPK_AbaqusUserElement_numsvars, IPK_AbaqusUserElement_properties, IPK_AbaqusUserElement_type, IPK_AbaqusUserElement_userElement, IR_GIVE_OPTIONAL_FIELD, jtype, lFlags, mcrd, mlvarx, nCoords, ndofel, npredef, nrhs, oofem::FEMComponent::number, oofem::Element::numberOfDofMans, numSvars, OOFEM_ERROR, PM_ELEMENT_ERROR_IFNOTSET, oofem::Element::postInitialize(), predef, rhs, svars, U, uel, uelobj, and V.
|
overridevirtual |
Updates element state after equilibrium in time step has been reached. Default implementation updates all integration rules defined by integrationRulesArray member variable. Doing this, all integration points and their material statuses are updated also. All temporary history variables, which now describe equilibrium state are copied into equilibrium ones. The existing internal state is used for update.
| tStep | Time step for newly reached state. |
Reimplemented from oofem::Element.
Definition at line 208 of file AbaqusUserElement.C.
References giveInternalForcesVector().
|
overridevirtual |
Updates element state after equilibrium in time step has been reached. Default implementation updates all integration rules defined by integrationRulesArray member variable. Doing this, all integration points and their material statuses are updated also. All temporary history variables, which now describe equilibrium state are copied into equilibrium ones. The existing internal state is used for update.
| tStep | Time step for newly reached state. |
Reimplemented from oofem::Element.
Definition at line 199 of file AbaqusUserElement.C.
References amatrx, hasTangentFlag, rhs, svars, tempAmatrx, tempRHS, tempSvars, and oofem::StructuralElement::updateYourself().
|
private |
Definition at line 116 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Definition at line 136 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Element amatrx.
Definition at line 104 of file AbaqusUserElement.h.
Referenced by postInitialize(), and updateYourself().
|
private |
Definition at line 85 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), and postInitialize().
|
private |
Definition at line 136 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Definition at line 84 of file AbaqusUserElement.h.
Referenced by giveDofManDofIDMask(), giveInputRecord(), and initializeFrom().
|
private |
Definition at line 117 of file AbaqusUserElement.h.
Referenced by computeStiffnessMatrix(), giveInternalForcesVector(), giveInternalForcesVector(), and postInitialize().
|
private |
energy
Definition at line 130 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
File containing the uel function.
Definition at line 158 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), initializeFrom(), and postInitialize().
|
private |
Keeps track of whether the tangent has been obtained already.
Definition at line 146 of file AbaqusUserElement.h.
Referenced by AbaqusUserElement(), hasTangent(), letTempTangentBe(), and updateYourself().
|
staticprivate |
Definition at line 162 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), and initializeFrom().
|
staticprivate |
Definition at line 166 of file AbaqusUserElement.h.
|
staticprivate |
Definition at line 161 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), initializeFrom(), and postInitialize().
|
staticprivate |
Definition at line 163 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), initializeFrom(), and postInitialize().
|
staticprivate |
Definition at line 164 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), initializeFrom(), and postInitialize().
|
staticprivate |
Definition at line 165 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), initializeFrom(), and postInitialize().
|
staticprivate |
Definition at line 160 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), initializeFrom(), and postInitialize().
|
private |
Definition at line 137 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
jprops
Definition at line 143 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Element type.
Definition at line 97 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), initializeFrom(), and postInitialize().
|
private |
kinc, kstep
Definition at line 123 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Definition at line 123 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
LFlags.
Definition at line 120 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
mcrd
Definition at line 110 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Definition at line 134 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
mlvarx
Definition at line 113 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Coord transf.
Definition at line 83 of file AbaqusUserElement.h.
Referenced by initializeFrom(), and postInitialize().
|
private |
|
private |
ndofel
Definition at line 107 of file AbaqusUserElement.h.
Referenced by computeConsistentMassMatrix(), giveInternalForcesVector(), and postInitialize().
|
private |
predef
Definition at line 126 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Element rhs.
Definition at line 100 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Number of status variables.
Definition at line 88 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), initializeFrom(), and postInitialize().
|
private |
params
Definition at line 140 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Definition at line 127 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Element properties.
Definition at line 91 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), and initializeFrom().
|
private |
Definition at line 101 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), postInitialize(), and updateYourself().
|
private |
Status variables.
Definition at line 94 of file AbaqusUserElement.h.
Referenced by giveStateVector(), postInitialize(), and updateYourself().
|
private |
Definition at line 104 of file AbaqusUserElement.h.
Referenced by giveTempTangent(), letTempTangentBe(), and updateYourself().
|
private |
Definition at line 101 of file AbaqusUserElement.h.
Referenced by letTempRhsBe(), and updateYourself().
|
private |
Definition at line 94 of file AbaqusUserElement.h.
Referenced by giveTempStateVector(), letTempSvarsBe(), and updateYourself().
|
private |
Inputs to element routines. Velocity and Acceleration currently ignored.
Definition at line 116 of file AbaqusUserElement.h.
Referenced by computeStiffnessMatrix(), giveInternalForcesVector(), giveInternalForcesVector(), and postInitialize().
|
private |
Pointer to the dynamically loaded uel-function (translated to C).
Definition at line 149 of file AbaqusUserElement.h.
Referenced by AbaqusUserElement(), giveInternalForcesVector(), and postInitialize().
|
private |
Dynamically loaded uel.
Definition at line 80 of file AbaqusUserElement.h.
Referenced by AbaqusUserElement(), postInitialize(), and ~AbaqusUserElement().
|
private |
Definition at line 116 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().