|
OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This class allows for custom user materials from Abaqus (UMAT). More...
#include <abaqususermaterial.h>
Inheritance diagram for oofem::AbaqusUserMaterial:
Collaboration diagram for oofem::AbaqusUserMaterial:Public Member Functions | |
| AbaqusUserMaterial (int n, Domain *d) | |
| Constructor. More... | |
| virtual | ~AbaqusUserMaterial () |
| Destructor. More... | |
| virtual IRResultType | initializeFrom (InputRecord *ir) |
| Reads the following values;. More... | |
| virtual void | giveInputRecord (DynamicInputRecord &input) |
| Setups the input record string of receiver. More... | |
| virtual MaterialStatus * | CreateStatus (GaussPoint *gp) const |
| Creates new copy of associated status and inserts it into given integration point. More... | |
| virtual void | give3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point. More... | |
| virtual void | give3dMaterialStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| virtual void | givePlaneStrainStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual void | giveRealStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
| virtual void | giveFirstPKStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
| virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
| Returns the integration point corresponding value in Reduced form. More... | |
| virtual int | hasNonLinearBehaviour () |
| Returns nonzero if receiver is non linear. More... | |
| virtual const char * | giveClassName () const |
| virtual const char * | giveInputRecordName () const |
Public Member Functions inherited from oofem::StructuralMaterial | |
| StructuralMaterial (int n, Domain *d) | |
| Constructor. More... | |
| virtual | ~StructuralMaterial () |
| Destructor. More... | |
| virtual int | hasMaterialModeCapability (MaterialMode mode) |
| Tests if material supports material mode. More... | |
| virtual void | giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| Computes the stiffness matrix for giveRealStressVector of receiver in given integration point, respecting its history. More... | |
| virtual void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) |
| Computes the real stress vector for given total strain and integration point. More... | |
| virtual void | giveRealStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector_3d. More... | |
| virtual void | giveRealStressVector_StressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep) |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· More... | |
| virtual void | giveRealStressVector_ShellStressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep) |
| virtual void | giveRealStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector_StressControl. More... | |
| virtual void | giveRealStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector_StressControl. More... | |
| virtual void | giveRealStressVector_Warping (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector_StressControl. More... | |
| virtual void | giveRealStressVector_2dBeamLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector_StressControl. More... | |
| virtual void | giveRealStressVector_PlateLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector_StressControl. More... | |
| virtual void | giveRealStressVector_Fiber (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation relies on giveRealStressVector_StressControl. More... | |
| virtual void | giveRealStressVector_Lattice2d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| virtual void | giveRealStressVector_Lattice3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| virtual void | giveRealStressVector_2dPlateSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| Default implementation is not provided. More... | |
| virtual void | giveRealStressVector_3dBeamSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
| virtual void | giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| Prototype for computation of Eshelby stress. More... | |
| void | give_dPdF_from (const FloatMatrix &dSdE, FloatMatrix &answer, GaussPoint *gp) |
| void | convert_dSdE_2_dPdF (FloatMatrix &answer, const FloatMatrix &dSdE, const FloatArray &S, const FloatArray &F, MaterialMode matMode) |
| virtual void | giveThermalDilatationVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
| Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis. More... | |
| double | giveReferenceTemperature () |
| Returns the reference temperature of receiver. More... | |
| virtual void | computeStressIndependentStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) |
| Computes reduced strain vector in given integration point, generated by internal processes in material, which are independent on loading in particular integration point. More... | |
| virtual void | computeStressIndependentStrainVector_3d (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode) |
| virtual void | give3dMaterialStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| void | giveStressDependentPartOfStrainVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) |
| Method for subtracting from reduced space strain vector its stress-independent parts (caused by temperature, shrinkage, creep and possibly by other phenomena). More... | |
| void | giveStressDependentPartOfStrainVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) |
| virtual int | setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) |
| Sets the value of a certain variable at a given integration point to the given value. More... | |
| virtual void | give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing 2d beam layer stiffness matrix of receiver. More... | |
| virtual void | givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing 2d plate layer stiffness matrix of receiver. More... | |
| virtual void | giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing 1d fiber stiffness matrix of receiver. More... | |
| virtual void | give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing 2d lattice stiffness matrix of receiver. More... | |
| virtual void | give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing 3d lattice stiffness matrix of receiver. More... | |
| virtual void | give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing stiffness matrix of plate subsoil model. More... | |
| virtual void | give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing stiffness matrix of beam3d subsoil model. More... | |
| virtual void | giveFirstPKStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| Default implementation relies on giveFirstPKStressVector_3d. More... | |
| virtual void | giveFirstPKStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| Default implementation relies on giveFirstPKStressVector_3d. More... | |
| virtual void | giveFirstPKStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| Default implementation relies on giveFirstPKStressVector_3d. More... | |
| virtual void | giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing plane stress stiffness matrix of receiver. More... | |
| virtual void | givePlaneStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual void | givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual void | givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing plane strain stiffness matrix of receiver. More... | |
| virtual void | givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual void | give1dStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| Method for computing 1d stiffness matrix of receiver. More... | |
| virtual void | give1dStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Public Member Functions inherited from oofem::Material | |
| Material (int n, Domain *d) | |
| Constructor. More... | |
| virtual | ~Material () |
| Destructor. More... | |
| virtual bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) |
| Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More... | |
| virtual double | give (int aProperty, GaussPoint *gp) |
| Returns the value of material property 'aProperty'. More... | |
| virtual bool | hasProperty (int aProperty, GaussPoint *gp) |
| Returns true if 'aProperty' exists on material. More... | |
| virtual void | modifyProperty (int aProperty, double value, GaussPoint *gp) |
| Modify 'aProperty', which already exists on material. More... | |
| double | giveCastingTime () |
| virtual bool | isActivated (TimeStep *tStep) |
| virtual int | hasCastingTimeSupport () |
| Tests if material supports casting time. More... | |
| virtual void | printYourself () |
| Prints receiver state on stdout. Useful for debugging. More... | |
| virtual contextIOResultType | saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| Stores integration point state to output stream. More... | |
| virtual contextIOResultType | restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| Reads integration point state to output stream. More... | |
| virtual int | checkConsistency () |
| Allows programmer to test some internal data, before computation begins. More... | |
| virtual int | initMaterial (Element *element) |
| Optional function to call specific procedures when initializing a material. More... | |
| virtual MaterialStatus * | giveStatus (GaussPoint *gp) const |
| Returns material status of receiver in given integration point. More... | |
| virtual int | packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
| Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More... | |
| virtual int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
| Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More... | |
| virtual int | estimatePackSize (DataStream &buff, GaussPoint *ip) |
| Estimates the necessary pack size to hold all packed data of receiver. More... | |
| virtual double | predictRelativeComputationalCost (GaussPoint *gp) |
| Returns the weight representing relative computational cost of receiver The reference material model is linear isotropic material - its weight is set to 1.0 The other material models should compare to this reference model. More... | |
| virtual double | predictRelativeRedistributionCost (GaussPoint *gp) |
| Returns the relative redistribution cost of the receiver. More... | |
| virtual void | initTempStatus (GaussPoint *gp) |
| Initializes temporary variables stored in integration point status at the beginning of new time step. 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 | updateLocalNumbering (EntityRenumberingFunctor &f) |
| Local renumbering support. 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... | |
| virtual Interface * | giveInterface (InterfaceType t) |
| Interface requesting service. More... | |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). More... | |
Private Attributes | |
| void * | umatobj |
| Dynamically loaded umat. More... | |
| void(* | umat )(double *stress, double *statev, double *ddsdde, double *sse, double *spd, double *scd, double *rpl, double *ddsddt, double *drplde, double *drpldt, double *stran, double *dstran, double time[2], double *dtime, double *temp, double *dtemp, double predef[1], double dpred[1], char cmname[80], int *ndi, int *nshr, int *ntens, int *nstatv, double *props, int *nprops, double coords[3], double *drot, double *pnewdt, double *celent, double *dfgrd0, double *dfgrd1, int *noel, int *npt, int *layer, int *kspt, int *kstep, int *kinc) |
| Pointer to the dynamically loaded umat-function (translated to C) More... | |
| char | cmname [80] |
| Name for material routine. More... | |
| int | numState |
| Size of the state vector. More... | |
| FloatArray | properties |
| Material properties. More... | |
| FloatArray | initialStress |
| Initial stress. More... | |
| int | mStressInterpretation |
| Flag to determine how the stress and Jacobian are interpreted. More... | |
| bool | mUseNumericalTangent |
| Flag to determine if numerical tangent should be used. More... | |
| double | mPerturbation |
| Size of perturbation if numerical tangent is used. More... | |
| std::string | filename |
| Name of the file that contains the umat function. More... | |
Static Private Attributes | |
| static int | n = 1 |
| Gausspoint counter. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from oofem::StructuralMaterial | |
| static int | giveSymVI (int ind1, int ind2) |
| static int | giveVI (int ind1, int ind2) |
| static int | giveVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| Returns a mask of the vector indicies corresponding to components in a general (non-symmetric) second order tensor of some stress/strain/deformation measure that performes work. More... | |
| static int | giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode) |
| The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor. More... | |
| static void | giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| Gives the inverted version of giveVoigtVectorMask. More... | |
| static int | giveSizeOfVoigtVector (MaterialMode mmode) |
| Returns the size of reduced stress/strain vector according to given mode. More... | |
| static int | giveSizeOfVoigtSymVector (MaterialMode mmode) |
| Returns the size of symmetric part of a reduced stress/strain vector according to given mode. More... | |
| static void | giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. More... | |
| static void | giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced deformation gradient Voigt vector (2nd order tensor). More... | |
| static void | giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. More... | |
| static void | giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. More... | |
| static void | giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. More... | |
| static void | giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More... | |
| static void | giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. More... | |
| static void | giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More... | |
| static void | transformStrainVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &strainVector, bool transpose=false) |
| Transforms 3d strain vector into another coordinate system. More... | |
| static void | transformStressVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &stressVector, bool transpose=false) |
| Transforms 3d stress vector into another coordinate system. More... | |
| static double | computeVonMisesStress (const FloatArray *currentStress) |
| Computes equivalent of von Mises stress. More... | |
| static void | giveStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
| Computes 3d strain vector transformation matrix from standard vector transformation matrix. More... | |
| static void | give2DStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
| Computes 2d strain vector transformation matrix from standard vector transformation matrix. More... | |
| static void | giveStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
| Computes 3d stress vector transformation matrix from standard vector transformation matrix. More... | |
| static void | givePlaneStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false) |
| Computes 2d stress vector transformation matrix from standard vector transformation matrix. More... | |
| static void | sortPrincDirAndValCloseTo (FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir) |
| Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDir) to be closed to some (often previous) principal directions (toPDir). More... | |
| static void | computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode) |
| Common functions for convenience. More... | |
| static void | computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode) |
| Computes principal values and directions of stress or strain vector. More... | |
| static double | computeDeviatoricVolumetricSplit (FloatArray &dev, const FloatArray &s) |
| Computes split of receiver into deviatoric and volumetric part. More... | |
| static void | computeDeviatoricVolumetricSum (FloatArray &s, const FloatArray &dev, double mean) |
| static void | applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu) |
| static void | applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double GModulus) |
| static void | applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu) |
| static void | applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double GModulus) |
| static void | applyElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu) |
| static void | applyElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu) |
| static double | computeStressNorm (const FloatArray &stress) |
| static double | computeFirstInvariant (const FloatArray &s) |
| static double | computeSecondStressInvariant (const FloatArray &s) |
| static double | computeThirdStressInvariant (const FloatArray &s) |
| static double | computeFirstCoordinate (const FloatArray &s) |
| static double | computeSecondCoordinate (const FloatArray &s) |
| static double | computeThirdCoordinate (const FloatArray &s) |
Static Public Attributes inherited from oofem::StructuralMaterial | |
| static std::vector< std::vector< int > > | vIindex |
| Voigt index map. More... | |
| static std::vector< std::vector< int > > | svIndex |
| Symmetric Voigt index map. More... | |
Protected Attributes inherited from oofem::StructuralMaterial | |
| double | referenceTemperature |
| Reference temperature (temperature, when material has been built into structure). More... | |
Protected Attributes inherited from oofem::Material | |
| Dictionary | propertyDictionary |
| Property dictionary. More... | |
| double | castingTime |
| Casting time. 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... | |
This class allows for custom user materials from Abaqus (UMAT).
The umat material should be compiled as a shared library. For example, compile like this;
* gfortran -fPIC -fdefault-real-8 -shared -Wl,-soname,umat.so -o umat.so my_umat_code.f *
The interface for a user material file for Abaqus is as follows: SUBROUTINE UMAT(STRESS,STATEV,DDSDDE,SSE,SPD, SCD,RPL,DDSDDT,DRPLDE,DRPLDT, STRAN,DSTRAN,TIME,DTIME,TEMP,DTEMP,PREDEF,DPRED,CMNAME, NDI,NSHR,NTENS,NSTATV,PROPS,NPROPS,COORDS,DROT,PNEWDT, CELENT,DFGRD0,DFGRD1,NOEL,NPT,LAYER,KSPT,KSTEP,KINC)
CHARACTER*80 CMNAME DOUBLE PRECISION STRESS(NTENS),STATEV(NSTATV), DDSDDE(NTENS,NTENS),DDSDDT(NTENS),DRPLDE(NTENS), STRAN(NTENS),DSTRAN(NTENS),TIME(2),PREDEF(1),DPRED(1), PROPS(NPROPS),COORDS(3),DROT(3,3),DFGRD0(3,3),DFGRD1(3,3)
Definition at line 80 of file abaqususermaterial.h.
| oofem::AbaqusUserMaterial::AbaqusUserMaterial | ( | int | n, |
| Domain * | d | ||
| ) |
Constructor.
Definition at line 53 of file abaqususermaterial.C.
|
virtual |
|
virtual |
Creates new copy of associated status and inserts it into given integration point.
| gp | Integration point where newly created status will be stored. |
Reimplemented from oofem::Material.
Definition at line 144 of file abaqususermaterial.C.
References oofem::FEMComponent::giveDomain(), n, and numState.
|
virtual |
Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point.
| answer | Computed results. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 150 of file abaqususermaterial.C.
References oofem::FloatArray::at(), oofem::GaussPoint::giveMaterialStatus(), giveRealStressVector_3d(), oofem::FloatArray::giveSize(), oofem::Material::giveStatus(), oofem::AbaqusUserMaterialStatus::giveTempTangent(), oofem::AbaqusUserMaterialStatus::hasTangent(), oofem::FloatMatrix::printYourself(), oofem::FloatArray::subtract(), oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
virtual |
Reimplemented from oofem::StructuralMaterial.
Definition at line 187 of file abaqususermaterial.C.
References oofem::FloatArray::at(), giveFirstPKStressVector_3d(), oofem::GaussPoint::giveMaterialStatus(), oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::giveTempFVector(), oofem::AbaqusUserMaterialStatus::giveTempTangent(), oofem::AbaqusUserMaterialStatus::hasTangent(), mPerturbation, mUseNumericalTangent, oofem::FloatMatrix::setColumn(), oofem::FloatArray::subtract(), and oofem::FloatArray::times().
Referenced by givePlaneStrainStiffMtrx_dPdF().
|
inlinevirtual |
Reimplemented from oofem::StructuralMaterial.
Definition at line 162 of file abaqususermaterial.h.
|
virtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 412 of file abaqususermaterial.C.
References oofem::FloatMatrix::at(), oofem::FloatArray::beDifferenceOf(), oofem::FloatMatrix::beInverseOf(), oofem::FloatMatrix::beMatrixForm(), oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::beProductTOf(), oofem::FloatArray::beSymVectorForm(), oofem::FloatArray::beSymVectorFormOfStrain(), oofem::FloatMatrix::beTProductOf(), oofem::FloatMatrix::beUnitMatrix(), oofem::FloatArray::beVectorForm(), oofem::FloatArray::changeComponentOrder(), oofem::FloatMatrix::changeComponentOrder(), cmname, oofem::Element::computeGlobalCoordinates(), E, oofem::FloatMatrix::giveDeterminant(), oofem::GaussPoint::giveElement(), oofem::StructuralMaterialStatus::giveFVector(), oofem::TimeStep::giveMetaStepNumber(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FEMComponent::giveNumber(), oofem::FloatArray::givePointer(), oofem::FloatMatrix::givePointer(), oofem::FloatArray::giveSize(), oofem::AbaqusUserMaterialStatus::giveStateVector(), oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::giveStrainVector(), oofem::TimeStep::giveTargetTime(), oofem::TimeStep::giveTimeIncrement(), oofem::StructuralMaterialStatus::letTempFVectorBe(), oofem::StructuralMaterialStatus::letTempPVectorBe(), oofem::AbaqusUserMaterialStatus::letTempStateVectorBe(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), oofem::AbaqusUserMaterialStatus::letTempTangentBe(), mStressInterpretation, numState, OOFEM_LOG_DEBUG, properties, S, oofem::FloatMatrix::times(), and umat.
Referenced by give3dMaterialStiffnessMatrix_dPdF().
|
virtual |
Setups the input record string of receiver.
| input | Dynamic input record to be filled by receiver. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 134 of file abaqususermaterial.C.
References _IFT_AbaqusUserMaterial_name, _IFT_AbaqusUserMaterial_numState, _IFT_AbaqusUserMaterial_properties, _IFT_AbaqusUserMaterial_userMaterial, cmname, filename, oofem::StructuralMaterial::giveInputRecord(), numState, properties, and oofem::DynamicInputRecord::setField().
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 163 of file abaqususermaterial.h.
References _IFT_AbaqusUserMaterial_Name.
|
virtual |
Returns the integration point corresponding value in Reduced form.
| answer | Contain corresponding ip value, zero sized if not available. |
| gp | Integration point to which the value refers. |
| type | Determines the type of internal variable. |
| tStep | Determines the time step. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 620 of file abaqususermaterial.C.
References oofem::FloatArray::add(), oofem::StructuralMaterial::giveIPValue(), oofem::AbaqusUserMaterialStatus::giveStateVector(), oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::giveStressVector(), and initialStress.
|
virtual |
Reimplemented from oofem::StructuralMaterial.
Definition at line 246 of file abaqususermaterial.C.
References give3dMaterialStiffnessMatrix_dPdF(), and oofem::StructuralMaterial::giveReducedMatrixForm().
|
virtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 255 of file abaqususermaterial.C.
References oofem::FloatArray::add(), oofem::FloatArray::beDifferenceOf(), oofem::FloatMatrix::beUnitMatrix(), oofem::FloatArray::changeComponentOrder(), oofem::FloatMatrix::changeComponentOrder(), cmname, oofem::Element::computeGlobalCoordinates(), oofem::GaussPoint::giveElement(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FEMComponent::giveNumber(), oofem::FloatArray::givePointer(), oofem::FloatMatrix::givePointer(), oofem::FloatArray::giveSize(), oofem::AbaqusUserMaterialStatus::giveStateVector(), oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::giveStrainVector(), oofem::StructuralMaterialStatus::giveStressVector(), oofem::TimeStep::giveTargetTime(), oofem::TimeStep::giveTimeIncrement(), initialStress, oofem::AbaqusUserMaterialStatus::letTempStateVectorBe(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), oofem::AbaqusUserMaterialStatus::letTempTangentBe(), numState, OOFEM_LOG_DEBUG, properties, oofem::FloatArray::subtract(), and umat.
Referenced by give3dMaterialStiffnessMatrix().
|
inlinevirtual |
Returns nonzero if receiver is non linear.
Reimplemented from oofem::Material.
Definition at line 161 of file abaqususermaterial.h.
|
virtual |
Reads the following values;.
Reimplemented from oofem::StructuralMaterial.
Definition at line 75 of file abaqususermaterial.C.
References _IFT_AbaqusUserMaterial_initialStress, _IFT_AbaqusUserMaterial_name, _IFT_AbaqusUserMaterial_numericalTangent, _IFT_AbaqusUserMaterial_numericalTangentPerturbation, _IFT_AbaqusUserMaterial_numState, _IFT_AbaqusUserMaterial_properties, _IFT_AbaqusUserMaterial_userMaterial, cmname, filename, oofem::InputRecord::hasField(), oofem::StructuralMaterial::initializeFrom(), initialStress, IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_OK, mPerturbation, mUseNumericalTangent, numState, OOFEM_ERROR, properties, umat, and umatobj.
|
private |
Name for material routine.
Definition at line 97 of file abaqususermaterial.h.
Referenced by giveFirstPKStressVector_3d(), giveInputRecord(), giveRealStressVector_3d(), and initializeFrom().
|
private |
Name of the file that contains the umat function.
Definition at line 121 of file abaqususermaterial.h.
Referenced by giveInputRecord(), and initializeFrom().
|
private |
Initial stress.
Definition at line 103 of file abaqususermaterial.h.
Referenced by giveIPValue(), giveRealStressVector_3d(), and initializeFrom().
|
private |
Size of perturbation if numerical tangent is used.
Definition at line 118 of file abaqususermaterial.h.
Referenced by give3dMaterialStiffnessMatrix_dPdF(), and initializeFrom().
|
private |
Flag to determine how the stress and Jacobian are interpreted.
0 implies that the P and dPdF are returned from the umat routine.
Definition at line 108 of file abaqususermaterial.h.
Referenced by giveFirstPKStressVector_3d().
|
private |
Flag to determine if numerical tangent should be used.
Definition at line 113 of file abaqususermaterial.h.
Referenced by give3dMaterialStiffnessMatrix_dPdF(), and initializeFrom().
|
staticprivate |
Gausspoint counter.
Definition at line 84 of file abaqususermaterial.h.
Referenced by CreateStatus().
|
private |
Size of the state vector.
Definition at line 99 of file abaqususermaterial.h.
Referenced by CreateStatus(), giveFirstPKStressVector_3d(), giveInputRecord(), giveRealStressVector_3d(), and initializeFrom().
|
private |
Material properties.
Definition at line 101 of file abaqususermaterial.h.
Referenced by giveFirstPKStressVector_3d(), giveInputRecord(), giveRealStressVector_3d(), and initializeFrom().
|
private |
Pointer to the dynamically loaded umat-function (translated to C)
Definition at line 89 of file abaqususermaterial.h.
Referenced by giveFirstPKStressVector_3d(), giveRealStressVector_3d(), and initializeFrom().
|
private |
Dynamically loaded umat.
Definition at line 86 of file abaqususermaterial.h.
Referenced by initializeFrom(), and ~AbaqusUserMaterial().