OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
UEL interface from Abaqus user elements. More...
#include <AbaqusUserElement.h>
Public Member Functions | |
AbaqusUserElement (int n, Domain *d) | |
Constructor. More... | |
virtual | ~AbaqusUserElement () |
Destructor. More... | |
virtual void | computeConsistentMassMatrix (FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL) |
Computes consistent mass matrix of receiver using numerical integration over element volume. More... | |
virtual void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
Computes the stiffness matrix of receiver. More... | |
virtual void | giveInternalForcesVector (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) |
Evaluates nodal representation of real internal forces. More... | |
virtual void | giveInternalForcesVector (FloatArray &answer, TimeStep *tStep, FloatArray &U, FloatMatrix &DU, int useUpdatedGpRecord) |
virtual int | computeNumberOfDofs () |
Computes or simply returns total number of element's local DOFs. More... | |
virtual void | giveDofManDofIDMask (int inode, IntArray &answer) const |
Returns dofmanager dof mask for node. More... | |
virtual void | computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) |
Computes the unknown vector interpolated at the specified local coordinates. More... | |
virtual void | updateYourself (TimeStep *tStep) |
Updates element state after equilibrium in time step has been reached. More... | |
virtual void | updateInternalState (TimeStep *tStep) |
Updates element state after equilibrium in time step has been reached. More... | |
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 () |
virtual Interface * | giveInterface (InterfaceType it) |
Interface requesting service. More... | |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual void | giveInputRecord (DynamicInputRecord &input) |
Setups the input record string of receiver. More... | |
virtual void | postInitialize () |
Performs post initialization steps. More... | |
virtual const char * | giveClassName () const |
virtual const char * | giveInputRecordName () const |
virtual integrationDomain | giveIntegrationDomain () const |
Returns integration domain for receiver, used to initialize integration point over receiver volume. More... | |
virtual Element_Geometry_Type | giveGeometryType () const |
Returns the element geometry type. More... | |
Public Member Functions inherited from oofem::NLStructuralElement | |
NLStructuralElement (int n, Domain *d) | |
Constructor. More... | |
virtual | ~NLStructuralElement () |
Destructor. More... | |
int | giveGeometryMode () |
Returns the geometry mode describing the formulation used in the internal work 0 - Engineering (small deformation) stress-strain mode 1 - First Piola-Kirchhoff - Deformation gradient mode, P is defined as FS 2 - Second Piola-Kirchhoff - Green-Lagrange strain mode with deformation gradient as input (deprecated and not supported) More... | |
void | computeFirstPKStressVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
Computes the first Piola-Kirchhoff stress tensor on Voigt format. More... | |
void | computeCauchyStressVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
Computes the Cauchy stress tensor on Voigt format. More... | |
virtual void | computeInitialStressMatrix (FloatMatrix &answer, TimeStep *tStep) |
Computes the initial stiffness matrix of receiver. More... | |
void | computeStiffnessMatrix_withIRulesAsSubcells (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
Computes the stiffness matrix of receiver. More... | |
void | giveInternalForcesVector_withIRulesAsSubcells (FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0) |
Evaluates nodal representation of real internal forces. More... | |
virtual void | computeDeformationGradientVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
Computes the deformation gradient in Voigt form at integration point ip and at time step tStep. More... | |
double | computeCurrentVolume (TimeStep *tStep) |
Computes the current volume of element. More... | |
Public Member Functions inherited from oofem::StructuralElement | |
StructuralElement (int n, Domain *d) | |
Constructor. More... | |
virtual | ~StructuralElement () |
Destructor. More... | |
virtual void | giveCharacteristicMatrix (FloatMatrix &answer, CharType, TimeStep *tStep) |
Computes characteristic matrix of receiver of requested type in given time step. More... | |
virtual void | giveCharacteristicVector (FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) |
Computes characteristic vector of receiver of requested type in given time step. More... | |
virtual void | computeMassMatrix (FloatMatrix &answer, TimeStep *tStep) |
Computes mass matrix of receiver. More... | |
virtual void | computeLumpedMassMatrix (FloatMatrix &answer, TimeStep *tStep) |
Computes lumped mass matrix of receiver. More... | |
virtual void | giveMassMtrxIntegrationgMask (IntArray &answer) |
Returns mask indicating, which unknowns (their type and ordering is the same as element unknown vector) participate in mass matrix integration. More... | |
void | computeStiffnessMatrix_withIRulesAsSubcells (FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) |
virtual void | computeStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
Compute strain vector of receiver evaluated at given integration point at given time step from element displacement vector. More... | |
virtual void | computeResultingIPTemperatureAt (FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode) |
Computes at given time (tStep) the the resulting temperature component array. More... | |
virtual void | computeResultingIPEigenstrainAt (FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode) |
Computes at given time the resulting eigenstrain component array. More... | |
virtual int | adaptiveUpdate (TimeStep *tStep) |
Updates the internal state variables stored in all IPs according to already mapped state. More... | |
virtual int | giveInternalStateAtNode (FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep) |
Returns internal state variable (like stress,strain) at node of element in Reduced form, the way how is obtained is dependent on InternalValueType. More... | |
virtual void | showSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
Shows sparse structure. More... | |
virtual void | showExtendedSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness) More... | |
virtual void | computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) |
Computes the contribution of the given body load (volumetric). More... | |
virtual void | computeBoundarySurfaceLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
Computes the contribution of the given load at the given boundary surface in global coordinate system. More... | |
virtual void | computeBoundaryEdgeLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
Computes the contribution of the given load at the given boundary edge. More... | |
virtual void | computeEdgeNMatrix (FloatMatrix &answer, int boundaryID, const FloatArray &lcoords) |
computes edge interpolation matrix More... | |
virtual void | computeSurfaceNMatrix (FloatMatrix &answer, int boundaryID, const FloatArray &lcoords) |
Computes surface interpolation matrix. More... | |
StructuralCrossSection * | giveStructuralCrossSection () |
Helper function which returns the structural cross-section for the element. More... | |
virtual void | createMaterialStatus () |
virtual void | computeNmatrixAt (const FloatArray &iLocCoord, FloatMatrix &answer) |
Computes interpolation matrix for element unknowns. More... | |
virtual void | updateBeforeNonlocalAverage (TimeStep *tStep) |
Updates internal element state (in all integration points of receiver) before nonlocal averaging takes place. More... | |
virtual void | giveNonlocalLocationArray (IntArray &locationArray, const UnknownNumberingScheme &us) |
Returns the "nonlocal" location array of receiver. More... | |
virtual void | addNonlocalStiffnessContributions (SparseMtrx &dest, const UnknownNumberingScheme &s, TimeStep *tStep) |
Adds the "nonlocal" contribution to stiffness matrix, to account for nonlocality of material model. More... | |
Public Member Functions inherited from oofem::Element | |
Element (int n, Domain *aDomain) | |
Constructor. More... | |
Element (const Element &src)=delete | |
Element & | operator= (const Element &src)=delete |
virtual | ~Element () |
Virtual destructor. More... | |
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) |
Returns internal state variable (like stress,strain) at side of element in Reduced form If side is possessing DOFs, otherwise recover techniques will not work due to absence of side-shape functions. More... | |
int | giveLabel () const |
int | giveGlobalNumber () const |
void | setGlobalNumber (int num) |
Sets receiver globally unique number. More... | |
elementParallelMode | giveParallelMode () const |
Return elementParallelMode of receiver. More... | |
void | setParallelMode (elementParallelMode _mode) |
Sets parallel mode of element. More... | |
virtual elementParallelMode | giveKnotSpanParallelMode (int) const |
Returns the parallel mode for particular knot span of the receiver. More... | |
int | packUnknowns (DataStream &buff, TimeStep *tStep) |
Pack all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More... | |
int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep) |
Unpack and updates all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More... | |
int | estimatePackSize (DataStream &buff) |
Estimates the necessary pack size to hold all packed data of receiver. More... | |
const IntArray * | givePartitionList () const |
Returns partition list of receiver. More... | |
void | setPartitionList (IntArray &pl) |
Sets partition list of receiver. More... | |
virtual double | predictRelativeComputationalCost () |
Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element with linear approximation, single integration point and linear isotropic material. More... | |
virtual double | giveRelativeSelfComputationalCost () |
Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element. More... | |
virtual double | predictRelativeRedistributionCost () |
Returns the relative redistribution cost of the receiver. More... | |
IntArray * | giveBodyLoadArray () |
Returns array containing load numbers of loads acting on element. More... | |
IntArray * | giveBoundaryLoadArray () |
Returns array containing load numbers of boundary loads acting on element. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Stores receiver state to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Restores the receiver state previously written in stream. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to stream, for given time step. More... | |
void | giveLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
Returns the location array (array of code numbers) of receiver for given numbering scheme. More... | |
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) |
Returns the location array for the boundary of the element. More... | |
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 |
Returns i-th internal element dof manager of the receiver. More... | |
virtual double | giveCharacteristicValue (CharType type, TimeStep *tStep) |
Computes characteristic value of receiver of requested type in given time step. More... | |
virtual void | computeTangentFromSurfaceLoad (FloatMatrix &answer, SurfaceLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
Computes the tangent contribution of the given load at the given boundary. More... | |
virtual void | computeTangentFromEdgeLoad (FloatMatrix &answer, EdgeLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
Computes the tangent contribution of the given load at the given boundary. More... | |
const IntArray & | giveBodyLoadList () const |
Returns receiver list of bodyloads. More... | |
const IntArray & | giveBoundaryLoadList () const |
Returns receiver list of boundary loads. More... | |
void | computeVectorOf (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
Returns local vector of unknowns. More... | |
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) |
Boundary version of computeVectorOf. More... | |
void | computeVectorOf (PrimaryField &field, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
Returns local vector of unknowns. More... | |
void | computeVectorOfPrescribed (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
Returns local vector of prescribed unknowns. More... | |
void | computeVectorOfPrescribed (const IntArray &dofIDMask, ValueModeType type, TimeStep *tStep, FloatArray &answer) |
Returns local vector of prescribed unknowns. More... | |
virtual int | computeNumberOfGlobalDofs () |
Computes the total number of element's global dofs. More... | |
int | computeNumberOfPrimaryMasterDofs () |
Computes the total number of element's primary master DOFs. More... | |
virtual bool | computeGtoLRotationMatrix (FloatMatrix &answer) |
Returns transformation matrix from global c.s. More... | |
virtual bool | giveRotationMatrix (FloatMatrix &answer) |
Transformation matrices updates rotation matrix between element-local and primary DOFs, taking into account nodal c.s. More... | |
virtual bool | computeDofTransformationMatrix (FloatMatrix &answer, const IntArray &nodes, bool includeInternal) |
Returns transformation matrix for DOFs from global coordinate system to local coordinate system in nodes. More... | |
virtual void | giveInternalDofManDofIDMask (int inode, IntArray &answer) const |
Returns internal dofmanager dof mask for node. More... | |
virtual void | giveElementDofIDMask (IntArray &answer) const |
Returns element dof mask for node. More... | |
virtual double | computeVolumeAround (GaussPoint *gp) |
Returns volume related to given integration point. More... | |
virtual double | computeVolumeAreaOrLength () |
Computes the volume, area or length of the element depending on its spatial dimension. More... | |
double | computeMeanSize () |
Computes the size of the element defined as its length. More... | |
virtual double | computeVolume () |
Computes the volume. More... | |
virtual double | computeArea () |
Computes the area (zero for all but 2d geometries). More... | |
virtual double | computeLength () |
Computes the length (zero for all but 1D geometries) More... | |
virtual void | giveBoundaryEdgeNodes (IntArray &bNodes, int boundary) |
Returns list of receiver boundary nodes for given edge. More... | |
virtual void | giveBoundarySurfaceNodes (IntArray &bNodes, int boundary) |
Returns list of receiver boundary nodes for given surface. More... | |
virtual IntegrationRule * | giveBoundaryEdgeIntegrationRule (int order, int boundary) |
Returns boundary edge integration rule. More... | |
virtual IntegrationRule * | giveBoundarySurfaceIntegrationRule (int order, int boundary) |
Returns boundary surface integration rule. More... | |
int | giveDofManagerNumber (int i) const |
Translates local to global indices for dof managers. More... | |
const IntArray & | giveDofManArray () const |
void | addDofManager (DofManager *dMan) |
DofManager * | giveDofManager (int i) const |
Node * | giveNode (int i) const |
Returns reference to the i-th node of element. More... | |
virtual ElementSide * | giveSide (int i) const |
Returns reference to the i-th side of element. More... | |
virtual FEInterpolation * | giveInterpolation () const |
virtual FEInterpolation * | giveInterpolation (DofIDItem id) const |
Returns the interpolation for the specific dof id. More... | |
virtual Material * | giveMaterial () |
int | giveMaterialNumber () const |
CrossSection * | giveCrossSection () |
void | setMaterial (int matIndx) |
Sets the material of receiver. More... | |
virtual void | setCrossSection (int csIndx) |
Sets the cross section model of receiver. More... | |
virtual int | giveNumberOfDofManagers () const |
virtual int | giveNumberOfNodes () const |
Returns number of nodes of receiver. More... | |
void | setDofManagers (const IntArray &dmans) |
Sets receiver dofManagers. More... | |
void | setBodyLoads (const IntArray &bodyLoads) |
Sets receiver bodyLoadArray. More... | |
void | setIntegrationRules (std::vector< std::unique_ptr< IntegrationRule > > irlist) |
Sets integration rules. More... | |
virtual MaterialMode | giveMaterialMode () |
Returns material mode for receiver integration points. More... | |
virtual int | giveIntegrationRuleLocalCodeNumbers (IntArray &answer, IntegrationRule &ie) |
Assembles the code numbers of given integration element (sub-patch) This is done by obtaining list of nonzero shape functions and by collecting the code numbers of nodes corresponding to these shape functions. More... | |
int | giveRegionNumber () |
virtual void | initializeYourself (TimeStep *timeStepWhenICApply) |
Initialization according to state given by initial conditions. More... | |
virtual bool | isActivated (TimeStep *tStep) |
virtual bool | isCast (TimeStep *tStep) |
virtual void | initForNewStep () |
Initializes receivers state to new time step. More... | |
virtual int | giveSpatialDimension () |
Returns the element spatial dimension (1, 2, or 3). More... | |
virtual int | giveNumberOfBoundarySides () |
virtual int | giveDefaultIntegrationRule () const |
Returns id of default integration rule. More... | |
virtual IntegrationRule * | giveDefaultIntegrationRulePtr () |
Access method for default integration rule. More... | |
int | giveNumberOfIntegrationRules () |
virtual IntegrationRule * | giveIntegrationRule (int i) |
std::vector< std::unique_ptr< IntegrationRule > > & | giveIntegrationRulesArray () |
virtual int | testElementExtension (ElementExtension ext) |
Tests if the element implements required extension. More... | |
int | giveGlobalIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
virtual double | giveLengthInDir (const FloatArray &normalToCrackPlane) |
Default implementation returns length of element projection into specified direction. More... | |
virtual double | giveCharacteristicLength (const FloatArray &normalToCrackPlane) |
Returns the size of element in the given direction, in some cases adjusted (e.g. More... | |
double | giveCharacteristicLengthForPlaneElements (const FloatArray &normalToCrackPlane) |
Returns the size of element in the given direction if the direction is in the XY plane, otherwise gives the mean size defined as the square root of the element area. More... | |
double | giveCharacteristicLengthForAxisymmElements (const FloatArray &normalToCrackPlane) |
Returns the size of an axisymmetric element in the given direction if the direction is in the XY plane, otherwise gives the mean distance vrom the symmetry axis multiplied by pi. More... | |
virtual double | giveCharacteristicSize (GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method) |
Returns characteristic element size for a given integration point and given direction. More... | |
virtual double | giveParentElSize () const |
Returns the size (length, area or volume depending on element type) of the parent element. More... | |
virtual int | computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords) |
Computes the global coordinates from given element's local coordinates. More... | |
virtual bool | computeLocalCoordinates (FloatArray &answer, const FloatArray &gcoords) |
Computes the element local coordinates from given global coordinates. More... | |
virtual int | giveLocalCoordinateSystem (FloatMatrix &answer) |
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy. More... | |
virtual void | computeMidPlaneNormal (FloatArray &answer, const GaussPoint *gp) |
Computes mid-plane normal of receiver at integration point. More... | |
virtual int | adaptiveMap (Domain *oldd, TimeStep *tStep) |
Initializes the internal state variables stored in all IPs according to state in given domain. More... | |
virtual int | mapStateVariables (Domain &iOldDom, const TimeStep &iTStep) |
Maps the internal state variables stored in all IPs from the old domain to the new domain. More... | |
virtual int | adaptiveFinish (TimeStep *tStep) |
Finishes the mapping for given time step. More... | |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
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. More... | |
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. More... | |
Public Member Functions inherited from oofem::FEMComponent | |
FEMComponent (int n, Domain *d) | |
Regular constructor, creates component with given number and belonging to given domain. More... | |
virtual | ~FEMComponent () |
Virtual destructor. More... | |
Domain * | giveDomain () const |
virtual void | setDomain (Domain *d) |
Sets associated Domain. More... | |
int | giveNumber () const |
void | setNumber (int num) |
Sets number of receiver. More... | |
virtual void | printYourself () |
Prints receiver state on stdout. Useful for debugging. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Protected Member Functions | |
virtual void | computeStressVector (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) |
Computes the stress vector of receiver at given integration point, at time step tStep. More... | |
virtual void | computeConstitutiveMatrixAt (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) |
Computes constitutive matrix of receiver. More... | |
virtual void | computeBmatrixAt (GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS) |
Computes the geometrical matrix of receiver in given integration point. More... | |
virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
Returns the integration point corresponding value in full form. More... | |
Protected Member Functions inherited from oofem::NLStructuralElement | |
virtual int | checkConsistency () |
Performs consistency check. More... | |
virtual void | computeBHmatrixAt (GaussPoint *gp, FloatMatrix &answer) |
Computes a matrix which, multiplied by the column matrix of nodal displacements, gives the displacement gradient stored by columns. More... | |
Protected Member Functions inherited from oofem::StructuralElement | |
virtual void | computeBodyLoadVectorAt (FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode) |
Computes the load vector due to body load acting on receiver, at given time step. More... | |
virtual int | giveNumberOfIPForMassMtrxIntegration () |
Return desired number of integration points for consistent mass matrix computation, if required. More... | |
void | condense (FloatMatrix *stiff, FloatMatrix *mass, FloatArray *load, IntArray *what) |
General service for condensation of stiffness and optionally load vector and mass or initial stress matrices of receiver. More... | |
virtual void | setupIRForMassMtrxIntegration (IntegrationRule &iRule) |
Setup Integration Rule Gauss Points for Mass Matrix integration. More... | |
virtual void | computePointLoadVectorAt (FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode, bool global=true) |
Computes point load vector contribution of receiver for given load (should has BoundaryLoad Base). More... | |
virtual void | giveEdgeDofMapping (IntArray &answer, int iEdge) const |
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element DOFs. More... | |
virtual void | giveSurfaceDofMapping (IntArray &answer, int iSurf) const |
Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" element DOFs. More... | |
virtual IntegrationRule * | GetSurfaceIntegrationRule (int order) |
virtual double | computeEdgeVolumeAround (GaussPoint *gp, int iEdge) |
Computes volume related to integration point on local edge. More... | |
virtual double | computeSurfaceVolumeAround (GaussPoint *gp, int iSurf) |
Computes volume related to integration point on local surface. More... | |
virtual int | computeLoadGToLRotationMtrx (FloatMatrix &answer) |
Returns transformation matrix from global coordinate system to local element coordinate system for element load vector components. More... | |
virtual int | computeLoadLEToLRotationMatrix (FloatMatrix &answer, int iEdge, GaussPoint *gp) |
Returns transformation matrix from local edge c.s to element local coordinate system of load vector components. More... | |
virtual int | computeLoadLSToLRotationMatrix (FloatMatrix &answer, int iSurf, GaussPoint *gp) |
Returns transformation matrix from local surface c.s to element local coordinate system of load vector components. More... | |
Protected Member Functions inherited from oofem::Element | |
virtual void | computeGaussPoints () |
Initializes the array of integration rules member variable. More... | |
Private Attributes | |
void * | uelobj |
Dynamically loaded uel. More... | |
int | nCoords |
Coord transf. More... | |
IntArray | dofs |
FloatMatrix | coords |
int | numSvars |
Number of status variables. More... | |
FloatArray | props |
Element properties. More... | |
FloatArray | svars |
Status variables. More... | |
FloatArray | tempSvars |
int | jtype |
Element type. More... | |
int | nrhs = 2 |
Element rhs. More... | |
FloatMatrix | rhs |
FloatMatrix | tempRHS |
FloatMatrix | amatrx |
Element amatrx. More... | |
FloatMatrix | tempAmatrx |
int | ndofel = 0 |
ndofel More... | |
int | mcrd = 2 |
mcrd More... | |
int | mlvarx = 1 |
mlvarx More... | |
FloatArray | U |
Inputs to element routines. Velocity and Acceleration currently ignored. More... | |
FloatArray | V |
FloatArray | A |
FloatMatrix | DU |
IntArray | lFlags |
LFlags. More... | |
int | kinc = 1 |
kinc, kstep More... | |
int | kstep = 1 |
int | npredef = 1 |
predef More... | |
FloatArray | predef |
FloatArray | energy |
energy More... | |
int | ndLoad = 0 |
loads More... | |
int | mdLoad = 0 |
FloatMatrix | adlmag |
FloatMatrix | ddlmag |
int * | jdltype = NULL |
double | params [3] |
params More... | |
IntArray | jprops |
jprops More... | |
bool | hasTangentFlag |
Keeps track of whether the tangent has been obtained already. More... | |
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) More... | |
std::string | filename |
File containing the uel function. More... | |
Additional Inherited Members | |
Protected Attributes inherited from oofem::NLStructuralElement | |
int | nlGeometry |
Flag indicating if geometrical nonlinearities apply. More... | |
Protected Attributes inherited from oofem::StructuralElement | |
std::unique_ptr< FloatArray > | initialDisplacements |
Initial displacement vector, describes the initial nodal displacements when element has been casted. More... | |
Protected Attributes inherited from oofem::Element | |
int | numberOfDofMans |
Number of dofmanagers. More... | |
IntArray | dofManArray |
Array containing dofmanager numbers. More... | |
int | material |
Number of associated material. More... | |
int | crossSection |
Number of associated cross section. More... | |
IntArray | bodyLoadArray |
Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver. More... | |
IntArray | boundaryLoadArray |
std::vector< std::unique_ptr< IntegrationRule > > | integrationRulesArray |
List of integration rules of receiver (each integration rule contains associated integration points also). More... | |
FloatMatrix | elemLocalCS |
Transformation material matrix, used in orthotropic and anisotropic materials, global->local transformation. More... | |
int | activityTimeFunction |
Element activity time function. If defined, nonzero value indicates active receiver, zero value inactive element. More... | |
int | globalNumber |
In parallel mode, globalNumber contains globally unique DoFManager number. More... | |
int | numberOfGaussPoints |
Number of integration points as specified by nip. More... | |
elementParallelMode | parallel_mode |
Determines the parallel mode of the element. More... | |
IntArray | partitions |
List of partition sharing the shared element or remote partition containing remote element counterpart. More... | |
Protected Attributes inherited from oofem::FEMComponent | |
int | number |
Component number. More... | |
Domain * | domain |
Link to domain object, useful for communicating with other FEM components. More... | |
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 77 of file AbaqusUserElement.h.
oofem::AbaqusUserElement::AbaqusUserElement | ( | int | n, |
Domain * | d | ||
) |
Constructor.
Definition at line 54 of file AbaqusUserElement.C.
|
virtual |
|
inlineprotectedvirtual |
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 227 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
virtual |
Computes consistent mass matrix of receiver using numerical integration over element volume.
Mass matrix is computed as , where 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 303 of file AbaqusUserElement.C.
References ndofel, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
inlineprotectedvirtual |
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 223 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
inlinevirtual |
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::StructuralElement.
Definition at line 175 of file AbaqusUserElement.h.
References OOFEM_ERROR, updateInternalState(), and updateYourself().
|
inlinevirtual |
Computes or simply returns total number of element's local DOFs.
Must be defined by particular element.
Reimplemented from oofem::Element.
Definition at line 173 of file AbaqusUserElement.h.
References giveDofManDofIDMask(), and ndofel.
|
virtual |
Computes the stiffness matrix of receiver.
The response is evaluated using , where is the B-matrix which produces the displacement gradient vector when multiplied with the solution vector a. Reduced integration are taken into account.
answer | Computed stiffness matrix. |
rMode | Response mode. |
tStep | Time step. |
Reimplemented from oofem::NLStructuralElement.
Definition at line 184 of file AbaqusUserElement.C.
References DU, giveInternalForcesVector(), giveTempTangent(), hasTangent(), and U.
|
inlineprotectedvirtual |
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 220 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
inlinevirtual |
Reimplemented from oofem::NLStructuralElement.
Definition at line 210 of file AbaqusUserElement.h.
|
virtual |
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 179 of file AbaqusUserElement.C.
References dofs.
Referenced by computeNumberOfDofs().
|
inlinevirtual |
Returns the element geometry type.
This information is assumed to be of general interest, but it is required only for some specialized tasks.
Reimplemented from oofem::Element.
Definition at line 215 of file AbaqusUserElement.h.
|
virtual |
Setups the input record string of receiver.
input | Dynamic input record to be filled by receiver. |
Reimplemented from oofem::NLStructuralElement.
Definition at line 162 of file AbaqusUserElement.C.
References _IFT_AbaqusUserElement_dofs, _IFT_AbaqusUserElement_numcoords, _IFT_AbaqusUserElement_numsvars, _IFT_AbaqusUserElement_properties, _IFT_AbaqusUserElement_type, _IFT_AbaqusUserElement_userElement, coords, dofs, filename, oofem::StructuralElement::giveInputRecord(), jtype, numSvars, props, and oofem::DynamicInputRecord::setField().
Referenced by giveTempTangent().
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 211 of file AbaqusUserElement.h.
References _IFT_AbaqusUserElement_Name.
|
inlinevirtual |
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 212 of file AbaqusUserElement.h.
References oofem::_UnknownIntegrationDomain.
|
virtual |
Interface requesting service.
Reimplemented from oofem::FEMComponent.
Definition at line 174 of file AbaqusUserElement.C.
Referenced by giveTempTangent().
|
virtual |
Evaluates nodal representation of real internal forces.
Necessary transformations are taken into account.
answer | Equivalent nodal forces vector. |
tStep | Time step |
useUpdatedGpRecord | If equal to zero, the stresses in integration points are computed (slow but safe). |
Reimplemented from oofem::NLStructuralElement.
Definition at line 211 of file AbaqusUserElement.C.
References oofem::Element::computeVectorOf(), DU, oofem::Domain::giveClassName(), oofem::FEMComponent::giveDomain(), oofem::FloatMatrix::setColumn(), U, and oofem::FloatMatrix::zero().
Referenced by computeStiffnessMatrix(), and updateInternalState().
|
virtual |
Definition at line 225 of file AbaqusUserElement.C.
References A, adlmag, oofem::IntArray::at(), coords, oofem::FloatMatrix::copyColumn(), ddlmag, energy, oofem::IntArray::givePointer(), oofem::FloatArray::givePointer(), oofem::FloatMatrix::givePointer(), oofem::IntArray::giveSize(), oofem::FloatArray::giveSize(), 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, oofem::FEMComponent::number, oofem::Element::numberOfDofMans, numSvars, params, predef, props, rhs, uel, and V.
|
inlineprotectedvirtual |
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::StructuralElement.
Definition at line 231 of file AbaqusUserElement.h.
References OOFEM_ERROR.
|
inlinevirtual |
Definition at line 193 of file AbaqusUserElement.h.
References svars.
Referenced by giveInternalForcesVector().
|
inlinevirtual |
Definition at line 196 of file AbaqusUserElement.h.
References tempSvars.
|
inlinevirtual |
Definition at line 199 of file AbaqusUserElement.h.
References giveInputRecord(), giveInterface(), initializeFrom(), postInitialize(), and tempAmatrx.
Referenced by computeStiffnessMatrix().
|
inline |
Definition at line 180 of file AbaqusUserElement.h.
References hasTangentFlag.
Referenced by computeStiffnessMatrix().
|
virtual |
Initializes receiver according to object description stored in input record.
This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.
ir | Input record to initialize from. |
Reimplemented from oofem::NLStructuralElement.
Definition at line 72 of file AbaqusUserElement.C.
References _IFT_AbaqusUserElement_dofs, _IFT_AbaqusUserElement_name, _IFT_AbaqusUserElement_numcoords, _IFT_AbaqusUserElement_numsvars, _IFT_AbaqusUserElement_properties, _IFT_AbaqusUserElement_type, _IFT_AbaqusUserElement_userElement, oofem::Element::dofManArray, dofs, filename, oofem::IntArray::giveSize(), oofem::StructuralElement::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_OK, jtype, nCoords, oofem::Element::numberOfDofMans, numSvars, OOFEM_ERROR, props, uel, and uelobj.
Referenced by giveTempTangent().
|
inlinevirtual |
Definition at line 187 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
inlinevirtual |
Definition at line 190 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
inlinevirtual |
Definition at line 183 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
virtual |
Performs post initialization steps.
Reimplemented from oofem::Element.
Definition at line 132 of file AbaqusUserElement.C.
References A, amatrx, oofem::FloatMatrix::at(), coords, DU, energy, oofem::Node::giveCoordinate(), oofem::Element::giveNode(), oofem::FloatMatrix::isNotEmpty(), lFlags, mcrd, mlvarx, nCoords, ndofel, npredef, nrhs, oofem::Element::numberOfDofMans, numSvars, oofem::Element::postInitialize(), predef, oofem::IntArray::resize(), oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), rhs, svars, U, and V.
Referenced by giveTempTangent().
|
virtual |
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::StructuralElement.
Definition at line 205 of file AbaqusUserElement.C.
References giveInternalForcesVector().
Referenced by computeField().
|
virtual |
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::StructuralElement.
Definition at line 196 of file AbaqusUserElement.C.
References amatrx, hasTangentFlag, rhs, svars, tempAmatrx, tempRHS, tempSvars, and oofem::StructuralElement::updateYourself().
Referenced by computeField().
|
private |
Definition at line 117 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Definition at line 137 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Element amatrx.
Definition at line 105 of file AbaqusUserElement.h.
Referenced by postInitialize(), and updateYourself().
|
private |
Definition at line 86 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), and postInitialize().
|
private |
Definition at line 137 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Definition at line 85 of file AbaqusUserElement.h.
Referenced by giveDofManDofIDMask(), giveInputRecord(), and initializeFrom().
|
private |
Definition at line 118 of file AbaqusUserElement.h.
Referenced by computeStiffnessMatrix(), giveInternalForcesVector(), and postInitialize().
|
private |
energy
Definition at line 131 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
File containing the uel function.
Definition at line 159 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), and initializeFrom().
|
private |
Keeps track of whether the tangent has been obtained already.
Definition at line 147 of file AbaqusUserElement.h.
Referenced by hasTangent(), and updateYourself().
|
private |
Definition at line 138 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
jprops
Definition at line 144 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Element type.
Definition at line 98 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), and initializeFrom().
|
private |
kinc, kstep
Definition at line 124 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Definition at line 124 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
LFlags.
Definition at line 121 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
mcrd
Definition at line 111 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Definition at line 135 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
mlvarx
Definition at line 114 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Coord transf.
Definition at line 84 of file AbaqusUserElement.h.
Referenced by initializeFrom(), and postInitialize().
|
private |
|
private |
ndofel
Definition at line 108 of file AbaqusUserElement.h.
Referenced by computeConsistentMassMatrix(), computeNumberOfDofs(), giveInternalForcesVector(), and postInitialize().
|
private |
predef
Definition at line 127 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Element rhs.
Definition at line 101 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Number of status variables.
Definition at line 89 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), initializeFrom(), and postInitialize().
|
private |
params
Definition at line 141 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector().
|
private |
Definition at line 128 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().
|
private |
Element properties.
Definition at line 92 of file AbaqusUserElement.h.
Referenced by giveInputRecord(), giveInternalForcesVector(), and initializeFrom().
|
private |
Definition at line 102 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), postInitialize(), and updateYourself().
|
private |
Status variables.
Definition at line 95 of file AbaqusUserElement.h.
Referenced by giveStateVector(), postInitialize(), and updateYourself().
|
private |
Definition at line 105 of file AbaqusUserElement.h.
Referenced by giveTempTangent(), and updateYourself().
|
private |
Definition at line 102 of file AbaqusUserElement.h.
Referenced by updateYourself().
|
private |
Definition at line 95 of file AbaqusUserElement.h.
Referenced by giveTempStateVector(), and updateYourself().
|
private |
Inputs to element routines. Velocity and Acceleration currently ignored.
Definition at line 117 of file AbaqusUserElement.h.
Referenced by computeStiffnessMatrix(), giveInternalForcesVector(), and postInitialize().
|
private |
Pointer to the dynamically loaded uel-function (translated to C)
Definition at line 150 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and initializeFrom().
|
private |
Dynamically loaded uel.
Definition at line 81 of file AbaqusUserElement.h.
Referenced by initializeFrom(), and ~AbaqusUserElement().
|
private |
Definition at line 117 of file AbaqusUserElement.h.
Referenced by giveInternalForcesVector(), and postInitialize().