OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::StructuralMaterial Class Reference

Abstract base class for all "structural" constitutive models. More...

#include <structuralmaterial.h>

+ Inheritance diagram for oofem::StructuralMaterial:
+ Collaboration diagram for oofem::StructuralMaterial:

Public Member Functions

 StructuralMaterial (int n, Domain *d)
 Constructor. More...
 
virtual ~StructuralMaterial ()
 Destructor. More...
 
virtual int hasMaterialModeCapability (MaterialMode mode)
 Tests if material supports material mode. More...
 
virtual const char * giveClassName () const
 
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 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_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. 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 (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 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 int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in Reduced form. 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...
 
Methods associated with the First PK stress tensor.

Computes the first Piola-Kirchhoff stress vector for given total deformation gradient and integration point.

The total deformation gradient is computed directly from displacement field at the given time step. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process. The First Piola-Kirchhoff stress vector is computed in Total Lagrangian mode

Parameters
answerContains result.
gpIntegration point.
reducedFDeformation gradient in in reduced form.
tStepCurrent time step (most models are able to respond only when tStep is current time step).
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 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...
 
Methods associated with the Cauchy stress tensor.

Computes the Cauchy stress vector for given increment of deformation gradient and given integration point.

The increment of deformation gradient is computed directly from displacement field at the given time step and it is computed wrt configuration which was reached in the last step. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process. The Cauchy stress vector is computed in Updated Lagrangian mode

Parameters
answerContains result.
gpIntegration point.
reducedFDeformation gradient in in reduced form.
tStepCurrent time step (most models are able to respond only when tStep is current time step).
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_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
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 hasNonLinearBehaviour ()
 Returns nonzero if receiver is non linear. More...
 
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 MaterialStatusgiveStatus (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 MaterialStatusCreateStatus (GaussPoint *gp) const
 Creates new copy of associated status and inserts it into given integration point. 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...
 
virtual const char * giveInputRecordName () const =0
 
DomaingiveDomain () 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 InterfacegiveInterface (InterfaceType t)
 Interface requesting service. More...
 
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros). More...
 

Static Public Member Functions

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

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

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...
 
Domaindomain
 Link to domain object, useful for communicating with other FEM components. More...
 

Friends

class CrossSection
 
class StructuralCrossSection
 
class SimpleCrossSection
 
class LayeredCrossSection
 

Detailed Description

Abstract base class for all "structural" constitutive models.

Todo:
Update documentation

It declares common services provided by all structural material models. The implementation of these services is partly left on derived classes, which will implement constitutive model dependent part. Some general purpose services are implemented on this level. For details, how to store material model related history variables in integration points, see base class Material documentation.

The constitutive model can in general support several material modes (plane stress, plane strain ,... modes). Its capabilities can be examined using hasMaterialModeCapability service. It is generally assumed, that results obtained from constitutive model services are according to valid material mode. This mode is determined from integration point, which is compulsory parameter of all material services. Structural material introduces several stress/strain modes. Full and reduced formats of stress/strain vectors are introduced for convenience. The full format includes all components, even if they are zero due to stress/strain mode nature, but in the reduced format, only generally nonzero components are stored. (full format must used only if absolutely necessary, to avoid wasting of space. It is used by output routines to print results in general form). Methods for converting vectors between full and reduced format are provided.

If in particular mode particular stress component is zero, the corresponding strain is not computed and not stored in reduced vector, and in full vector there is zero value on corresponding position. On the other hand, if some zero strain is imposed, On the other hand, if zero strain component is imposed this condition must be taken into account in geometrical relations (at element level), and corresponding component are included stress/strain reduced vectors.

Structural material introduces following basic stress/strain modes

  • 3d state - all components of general stress/strain vector are generally nonzero. General 3d strain vector has following components {sig_xx, sig_yy, sig_zz, tau_yz, tau_xz, tau_xy}
  • plane stress - sig_zz = tau_yz = tau_xz = 0.
  • plane strain - eps_z = gamma_xz = gamma_yz = 0. Note: as already described, if zero strain component is imposed (Plane strain, ..) this condition must be taken into account in geometrical relations, and corresponding component has to be included in reduced vector.
  • 1d uniaxial state - sigma_y = sigma_z = tau_yz = tau_zx = tau_xy = 0.
  • 2d beam layer - sigma_y=sigma_z=tau_zy=tau_xy = 0.
  • 3d shell layer, 2d plate layer - sigma_z = 0.

Derived classes can of course extend those modes. Generally speaking, there are following major tasks, covered by declared services.

  • Computing real/second PK stress vector (tensor) at integration point for given strain increment and updating its state (still temporary state, after overall equilibrium is reached).
  • Updating its state (final state), when equilibrium has been reached.
  • Returning its material stiffness (and/or flexibility) matrices for given material mode.
  • Storing/restoring its context to stream.
  • Returning its material properties.

Structural material services should not be called directly by elements. Instead, they always should pass their requests to corresponding cross section model. Cross section performs all necessary integration over its volume and invokes material model services.

Author
almost everyone
Jim Brouzoulis
Mikael Öhman

Definition at line 113 of file structuralmaterial.h.

Constructor & Destructor Documentation

oofem::StructuralMaterial::StructuralMaterial ( int  n,
Domain d 
)

Constructor.

Creates material with given number, belonging to given domain.

Parameters
nMaterial number.
dDomain to which new material will belong.

Definition at line 64 of file structuralmaterial.C.

Referenced by giveVI().

virtual oofem::StructuralMaterial::~StructuralMaterial ( )
inlinevirtual

Destructor.

Definition at line 136 of file structuralmaterial.h.

References hasMaterialModeCapability().

Member Function Documentation

void oofem::StructuralMaterial::applyDeviatoricElasticCompliance ( FloatArray strain,
const FloatArray stress,
double  EModulus,
double  nu 
)
static
void oofem::StructuralMaterial::applyDeviatoricElasticCompliance ( FloatArray strain,
const FloatArray stress,
double  GModulus 
)
static

Definition at line 1603 of file structuralmaterial.C.

References oofem::FloatArray::resize().

void oofem::StructuralMaterial::applyDeviatoricElasticStiffness ( FloatArray stress,
const FloatArray strain,
double  EModulus,
double  nu 
)
static
void oofem::StructuralMaterial::applyDeviatoricElasticStiffness ( FloatArray stress,
const FloatArray strain,
double  GModulus 
)
static

Definition at line 1622 of file structuralmaterial.C.

References oofem::FloatArray::resize().

void oofem::StructuralMaterial::applyElasticCompliance ( FloatArray strain,
const FloatArray stress,
double  EModulus,
double  nu 
)
static
void oofem::StructuralMaterial::applyElasticStiffness ( FloatArray stress,
const FloatArray strain,
double  EModulus,
double  nu 
)
static
double oofem::StructuralMaterial::computeDeviatoricVolumetricSplit ( FloatArray dev,
const FloatArray s 
)
static
void oofem::StructuralMaterial::computeDeviatoricVolumetricSum ( FloatArray s,
const FloatArray dev,
double  mean 
)
static
double oofem::StructuralMaterial::computeFirstCoordinate ( const FloatArray s)
static

Definition at line 1699 of file structuralmaterial.C.

References computeFirstInvariant().

Referenced by giveReferenceTemperature().

double oofem::StructuralMaterial::computeFirstInvariant ( const FloatArray s)
static
void oofem::StructuralMaterial::computePrincipalValDir ( FloatArray answer,
FloatMatrix dir,
const FloatArray s,
stressStrainPrincMode  mode 
)
static

Computes principal values and directions of stress or strain vector.

Parameters
answerComputed principal values.
dirPrincipal directions (stored column wise).
sStress/strain vector.
modeStress strain principal mode.

Definition at line 1423 of file structuralmaterial.C.

References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::giveSize(), oofem::FloatMatrix::jaco_(), OOFEM_SERROR, oofem::principal_deviatoricstress, oofem::principal_strain, oofem::principal_stress, oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), oofem::FloatArray::zero(), and oofem::FloatMatrix::zero().

Referenced by oofem::FCMMaterial::checkStrengthCriterion(), oofem::ConcreteDPM2::computeAlpha(), oofem::ConcreteDPM::computeDCosThetaDStress(), oofem::IsotropicDamageMaterial1::computeEta(), oofem::RankinePlasticMaterial::computeStressGradientVector(), oofem::DruckerPragerCutMat::computeStressGradientVector(), giveIPValue(), oofem::MDM::giveRawMDMParameters(), oofem::RCSDNLMaterial::giveRealStressVector(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), oofem::RCM2Material::giveRealStressVector(), oofem::FCMMaterial::giveRealStressVector(), oofem::Concrete2::giveRealStressVector_PlateLayer(), giveReferenceTemperature(), oofem::IDNLMaterial::giveRemoteNonlocalStiffnessContribution(), oofem::MazarsMaterial::initDamaged(), oofem::IsotropicDamageMaterial1::initDamaged(), oofem::ConcreteDPM::initDamaged(), oofem::ConcreteDPM2::initDamaged(), oofem::NCPrincipalStrain::nucleateEnrichmentItems(), oofem::NCPrincipalStress::nucleateEnrichmentItems(), oofem::ConcreteDPM::performRegularReturn(), oofem::ConcreteDPM2::performRegularReturn(), and oofem::PLPrincipalStrain::propagateInterface().

void oofem::StructuralMaterial::computePrincipalValues ( FloatArray answer,
const FloatArray s,
stressStrainPrincMode  mode 
)
static
void oofem::StructuralMaterial::computeStressIndependentStrainVector ( FloatArray answer,
GaussPoint gp,
TimeStep tStep,
ValueModeType  mode 
)
virtual

Computes reduced strain vector in given integration point, generated by internal processes in material, which are independent on loading in particular integration point.

Default implementation takes into account temperature induced strains and eigenstrains.

Parameters
answerReturned strain vector.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).
modeDetermines response mode (Total or incremental).
Todo:
Hack for loose gausspoints. We shouldn't ask for "gp->giveElement()". FIXME

Reimplemented in oofem::LatticeDamage2d, and oofem::RheoChainMaterial.

Definition at line 2269 of file structuralmaterial.C.

References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::Material::castingTime, oofem::FloatArray::clear(), oofem::Element::computeGlobalCoordinates(), oofem::StructuralElement::computeResultingIPEigenstrainAt(), oofem::StructuralElement::computeResultingIPTemperatureAt(), oofem::FEMComponent::domain, oofem::EngngModel::giveContext(), oofem::GaussPoint::giveElement(), oofem::Domain::giveEngngModel(), oofem::EngngModelContext::giveFieldManager(), oofem::GaussPoint::giveIntegrationRule(), oofem::TimeStep::giveIntrinsicTime(), oofem::GaussPoint::giveMaterialMode(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FEMComponent::giveNumber(), giveReducedSymVectorForm(), oofem::FloatArray::giveSize(), giveSizeOfVoigtSymVector(), giveThermalDilatationVector(), oofem::FloatArray::isEmpty(), oofem::FloatArray::isNotEmpty(), OOFEM_ERROR, referenceTemperature, oofem::FloatArray::times(), and oofem::FloatArray::zero().

Referenced by oofem::SteelRelaxMat::computeStressRelaxationStrainVector(), oofem::RheoChainMaterial::computeTrueStressIndependentStrainVector(), giveReferenceTemperature(), and giveStressDependentPartOfStrainVector().

double oofem::StructuralMaterial::computeThirdCoordinate ( const FloatArray s)
static
double oofem::StructuralMaterial::computeThirdStressInvariant ( const FloatArray s)
static

Definition at line 1688 of file structuralmaterial.C.

Referenced by computeThirdCoordinate(), and giveReferenceTemperature().

double oofem::StructuralMaterial::computeVonMisesStress ( const FloatArray currentStress)
static

Computes equivalent of von Mises stress.

Returns 0 if six stress components do not exist on the material.

Parameters
currentStressStress vector given by 6 components.

Definition at line 1737 of file structuralmaterial.C.

References oofem::FloatArray::at(), and oofem::FloatArray::giveSize().

Referenced by give3dMaterialStiffnessMatrix(), and giveIPValue().

void oofem::StructuralMaterial::convert_dSdE_2_dPdF ( FloatMatrix answer,
const FloatMatrix dSdE,
const FloatArray S,
const FloatArray F,
MaterialMode  matMode 
)
void oofem::StructuralMaterial::give1dStressStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing 1d stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 1d stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::ConcreteDPM2, oofem::MPlasticMaterial2, oofem::MPSDamMaterial, oofem::AnisotropicDamageMaterial, oofem::IsotropicDamageMaterial, oofem::RCM2Material, oofem::MPlasticMaterial, oofem::RheoChainMaterial, oofem::PlasticMaterial, oofem::PerfectlyPlasticMaterial, oofem::RankineMat, oofem::IsotropicLinearElasticMaterial, oofem::SteelRelaxMat, oofem::TrabBoneMaterial, oofem::MisesMatNl, oofem::MisesMat, oofem::MisesMatGrad, and oofem::IDGMaterial.

Definition at line 1109 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), give3dMaterialStiffnessMatrix(), and oofem::FloatMatrix::resize().

Referenced by oofem::MisesMat::give1dStressStiffMtrx(), oofem::RheoChainMaterial::give1dStressStiffMtrx(), oofem::MPlasticMaterial2::give1dStressStiffMtrx(), give1dStressStiffMtrx_dPdF(), oofem::SimpleCrossSection::give2dBeamStiffMtrx(), oofem::SimpleCrossSection::give3dBeamStiffMtrx(), give3dMaterialStiffnessMatrix(), oofem::SimpleCrossSection::giveCharMaterialStiffnessMatrix(), oofem::LinearElasticMaterial::giveRealStressVector_1d(), oofem::ConcreteDPM2::giveRealStressVector_1d(), giveStiffnessMatrix(), and oofem::SimpleCrossSection::giveStiffnessMatrix_1d().

void oofem::StructuralMaterial::give1dStressStiffMtrx_dCde ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::give1dStressStiffMtrx_dPdF ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::give2dBeamLayerStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing 2d beam layer stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d beam layer stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::MPlasticMaterial2, oofem::RCM2Material, oofem::MPlasticMaterial, oofem::PlasticMaterial, and oofem::PerfectlyPlasticMaterial.

Definition at line 1130 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), and give3dMaterialStiffnessMatrix().

Referenced by oofem::MPlasticMaterial2::give2dBeamLayerStiffMtrx(), oofem::LayeredCrossSection::give2dBeamStiffMtrx(), give3dMaterialStiffnessMatrix(), oofem::LinearElasticMaterial::giveRealStressVector_2dBeamLayer(), and giveStiffnessMatrix().

void oofem::StructuralMaterial::give2dLatticeStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing 2d lattice stiffness matrix of receiver.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::RheoChainMaterial.

Definition at line 1220 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by oofem::RheoChainMaterial::give2dLatticeStiffMtrx(), give3dMaterialStiffnessMatrix(), and giveStiffnessMatrix().

void oofem::StructuralMaterial::give2dPlateSubSoilStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing stiffness matrix of plate subsoil model.

Default method is emty; the implementation should be provided by the particular model.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::WinklerMaterial, and oofem::WinklerPasternakMaterial.

Definition at line 1245 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by oofem::SimpleCrossSection::give2dPlateSubSoilStiffMtrx(), and give3dMaterialStiffnessMatrix().

void oofem::StructuralMaterial::give2DStrainVectorTranformationMtrx ( FloatMatrix answer,
const FloatMatrix base,
bool  transpose = false 
)
static

Computes 2d strain vector transformation matrix from standard vector transformation matrix.

Parameters
answerTransformation matrix for strain vector.
baseA (2,2) matrix, where on each column are stored unit direction vectors of local coordinate axes to which we do transformation.
transposeDetermines if we transpose matrix before transforming.

Definition at line 1844 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beTranspositionOf(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().

Referenced by give3dMaterialStiffnessMatrix(), and oofem::FCMMaterial::initializeCrack().

void oofem::StructuralMaterial::give3dBeamSubSoilStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing stiffness matrix of beam3d subsoil model.

Default method is emty; the implementation should be provided by the particular model.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::WinklerMaterial.

Definition at line 1253 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by oofem::StructuralCrossSection::give3dBeamSubSoilStiffMtrx(), and give3dMaterialStiffnessMatrix().

void oofem::StructuralMaterial::give3dLatticeStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing 3d lattice stiffness matrix of receiver.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::RheoChainMaterial.

Definition at line 1232 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by oofem::RheoChainMaterial::give3dLatticeStiffMtrx(), give3dMaterialStiffnessMatrix(), and giveStiffnessMatrix().

virtual void oofem::StructuralMaterial::give3dMaterialStiffnessMatrix ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
inlinevirtual

Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point.

Parameters
answerComputed results.
modeMaterial response mode.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::ConcreteDPM2, oofem::ConcreteDPM, oofem::DustMaterial, oofem::MPSDamMaterial, oofem::MDM, oofem::AnisotropicDamageMaterial, oofem::DruckerPragerPlasticitySM, oofem::RheoChainMaterial, oofem::TrabBone3D, oofem::FCMMaterial, oofem::MPlasticMaterial2, oofem::RCM2Material, oofem::IsotropicDamageMaterial, oofem::MicroplaneMaterial, oofem::MPlasticMaterial, oofem::PerfectlyPlasticMaterial, oofem::CompoDamageMat, oofem::StructuralFE2Material, oofem::AbaqusUserMaterial, oofem::PlasticMaterial, oofem::OrthotropicLinearElasticMaterial, oofem::IsotropicLinearElasticMaterial, oofem::TrabBoneEmbed, oofem::TrabBoneNL3D, oofem::MisesMat, oofem::StructuralPythonMaterial, oofem::TrabBoneGrad3D, oofem::M1Material, oofem::MisesMatGrad, oofem::AnisotropicLinearElasticMaterial, oofem::MooneyRivlinMaterial, oofem::SimpleVitrificationMaterial, oofem::TutorialMaterial, oofem::StructuralMaterialSettable, and oofem::HyperElasticMaterial.

Definition at line 344 of file structuralmaterial.h.

References computeVonMisesStress(), give1dStressStiffMtrx(), give1dStressStiffMtrx_dCde(), give1dStressStiffMtrx_dPdF(), give2dBeamLayerStiffMtrx(), give2dLatticeStiffMtrx(), give2dPlateSubSoilStiffMtrx(), give2DStrainVectorTranformationMtrx(), give3dBeamSubSoilStiffMtrx(), give3dLatticeStiffMtrx(), give3dMaterialStiffnessMatrix_dCde(), give3dMaterialStiffnessMatrix_dPdF(), giveFiberStiffMtrx(), giveFullSymMatrixForm(), giveFullSymVectorForm(), giveFullVectorForm(), giveFullVectorFormF(), giveInvertedVoigtVectorMask(), giveIPValue(), givePlaneStrainStiffMtrx(), givePlaneStrainStiffMtrx_dCde(), givePlaneStrainStiffMtrx_dPdF(), givePlaneStressStiffMtrx(), givePlaneStressStiffMtrx_dCde(), givePlaneStressStiffMtrx_dPdF(), givePlaneStressVectorTranformationMtrx(), givePlateLayerStiffMtrx(), giveReducedMatrixForm(), giveReducedSymMatrixForm(), giveReducedSymVectorForm(), giveReducedVectorForm(), giveSizeOfVoigtSymVector(), giveSizeOfVoigtVector(), giveStrainVectorTranformationMtrx(), giveStressDependentPartOfStrainVector(), giveStressDependentPartOfStrainVector_3d(), giveStressVectorTranformationMtrx(), giveVoigtSymVectorMask(), giveVoigtVectorMask(), OOFEM_ERROR, setIPValue(), sortPrincDirAndValCloseTo(), transformStrainVectorTo(), and transformStressVectorTo().

Referenced by oofem::J2MPlasticMaterial::compute3dElasticModuli(), oofem::J2plasticMaterial::compute3dElasticModuli(), oofem::ConcreteDPM2::compute3dSecantStiffness(), oofem::XfemStructuralElementInterface::computeCohesiveTangent(), oofem::MDM::computeEffectiveStress(), give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), oofem::SimpleCrossSection::give3dDegeneratedShellStiffMtrx(), oofem::MisesMatGrad::give3dMaterialStiffnessMatrix(), oofem::MisesMat::give3dMaterialStiffnessMatrix(), oofem::PlasticMaterial::give3dMaterialStiffnessMatrix(), oofem::PerfectlyPlasticMaterial::give3dMaterialStiffnessMatrix(), oofem::MPlasticMaterial::give3dMaterialStiffnessMatrix(), oofem::IsotropicDamageMaterial::give3dMaterialStiffnessMatrix(), oofem::MPlasticMaterial2::give3dMaterialStiffnessMatrix(), oofem::RheoChainMaterial::give3dMaterialStiffnessMatrix(), oofem::ConcreteDPM::give3dMaterialStiffnessMatrix(), oofem::ConcreteDPM2::give3dMaterialStiffnessMatrix(), oofem::LargeStrainMasterMaterial::give3dMaterialStiffnessMatrix_dPdF(), give3dMaterialStiffnessMatrix_dPdF(), oofem::SimpleCrossSection::giveCharMaterialStiffnessMatrix(), giveFiberStiffMtrx(), oofem::IDNLMaterial::giveNormalElasticStiffnessMatrix(), oofem::LargeStrainMasterMaterialGrad::givePDGradMatrix_uu(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), givePlateLayerStiffMtrx(), oofem::LinearElasticMaterial::giveRealStressVector_3d(), oofem::ConcreteDPM2::giveRealStressVector_3d(), oofem::LinearElasticMaterial::giveRealStressVector_3dDegeneratedShell(), giveRealStressVector_ShellStressControl(), giveRealStressVector_StressControl(), giveStiffnessMatrix(), oofem::SimpleCrossSection::giveStiffnessMatrix_3d(), and oofem::StructuralMaterialEvaluator::solveYourself().

void oofem::StructuralMaterial::give3dMaterialStiffnessMatrix_dCde ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
virtual
Todo:
what should be default implementaiton?

Definition at line 731 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by give3dMaterialStiffnessMatrix(), and oofem::SimpleCrossSection::giveStiffnessMatrix_dCde().

virtual void oofem::StructuralMaterial::giveCauchyStressVector_1d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedF,
TimeStep tStep 
)
inlinevirtual
virtual void oofem::StructuralMaterial::giveCauchyStressVector_3d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedF,
TimeStep tStep 
)
inlinevirtual

Definition at line 245 of file structuralmaterial.h.

References OOFEM_ERROR.

Referenced by oofem::SimpleCrossSection::giveCauchyStresses().

virtual void oofem::StructuralMaterial::giveCauchyStressVector_PlaneStrain ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedF,
TimeStep tStep 
)
inlinevirtual

Definition at line 247 of file structuralmaterial.h.

References OOFEM_ERROR.

Referenced by oofem::SimpleCrossSection::giveCauchyStresses().

virtual void oofem::StructuralMaterial::giveCauchyStressVector_PlaneStress ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedF,
TimeStep tStep 
)
inlinevirtual

Definition at line 249 of file structuralmaterial.h.

References OOFEM_ERROR.

Referenced by oofem::SimpleCrossSection::giveCauchyStresses().

virtual const char* oofem::StructuralMaterial::giveClassName ( ) const
inlinevirtual
Returns
Class name of the receiver.

Implements oofem::FEMComponent.

Reimplemented in oofem::ConcreteDPM2, oofem::ConcreteDPM, oofem::DustMaterial, oofem::LatticeDamage2d, oofem::MDM, oofem::MPSMaterial, oofem::IsotropicDamageMaterial1, oofem::MPSDamMaterial, oofem::DruckerPragerPlasticitySM, oofem::TrabBone3D, oofem::RheoChainMaterial, oofem::AnisotropicDamageMaterial, oofem::Concrete2, oofem::FCMMaterial, oofem::MPlasticMaterial2, oofem::IsotropicDamageMaterial, oofem::RCM2Material, oofem::Eurocode2CreepMaterial, oofem::SteelRelaxMat, oofem::AbaqusUserMaterial, oofem::MPlasticMaterial, oofem::RankineMat, oofem::PerfectlyPlasticMaterial, oofem::TrabBoneMaterial, oofem::CompoDamageMat, oofem::StructuralFE2Material, oofem::PlasticMaterial, oofem::OrthotropicLinearElasticMaterial, oofem::B3SolidMaterial, oofem::RCSDMaterial, oofem::RCSDEMaterial, oofem::StructuralPythonMaterial, oofem::TrabBoneEmbed, oofem::RCSDNLMaterial, oofem::M4Material, oofem::IDNLMaterial, oofem::MazarsMaterial, oofem::DruckerPragerCutMat, oofem::M1Material, oofem::MazarsNLMaterial, oofem::FRCFCMNL, oofem::Masonry02, oofem::IsotropicLinearElasticMaterial, oofem::ConcreteFCM, oofem::MisesMat, oofem::B3Material, oofem::MisesMatNl, oofem::MicroMaterial, oofem::MooneyRivlinMaterial, oofem::RankineMatNl, oofem::TrabBoneNLEmbed, oofem::LinearElasticMaterial, oofem::RankineMatGrad, oofem::TrabBoneNL3D, oofem::SimpleVitrificationMaterial, oofem::TrabBoneNL, oofem::TrabBoneGrad3D, oofem::MisesMatGrad, oofem::WinklerMaterial, oofem::LargeStrainMasterMaterial, oofem::LargeStrainMasterMaterialGrad, oofem::WinklerPasternakMaterial, oofem::AnisotropicLinearElasticMaterial, oofem::MaxwellChainMaterial, oofem::KelvinChainSolidMaterial, oofem::HyperElasticMaterial, oofem::KelvinChainMaterial, oofem::TutorialMaterial, oofem::Concrete3, oofem::MicroplaneMaterial_Bazant, oofem::J2plasticMaterial, oofem::J2Mat, oofem::J2MPlasticMaterial, oofem::CebFip78Material, oofem::StructuralMaterialSettable, oofem::IDGMaterial, oofem::RankinePlasticMaterial, oofem::DoublePowerLawMaterial, oofem::Steel1, and oofem::DeformationTheoryMaterial.

Definition at line 139 of file structuralmaterial.h.

References giveFirstPKStressVector_1d(), giveFirstPKStressVector_3d(), giveFirstPKStressVector_PlaneStrain(), giveFirstPKStressVector_PlaneStress(), giveInputRecord(), giveRealStressVector(), giveRealStressVector_1d(), giveRealStressVector_2dBeamLayer(), giveRealStressVector_2dPlateSubSoil(), giveRealStressVector_3d(), giveRealStressVector_3dBeamSubSoil(), giveRealStressVector_Fiber(), giveRealStressVector_Lattice2d(), giveRealStressVector_Lattice3d(), giveRealStressVector_PlaneStrain(), giveRealStressVector_PlaneStress(), giveRealStressVector_PlateLayer(), giveRealStressVector_ShellStressControl(), giveRealStressVector_StressControl(), giveRealStressVector_Warping(), giveStiffnessMatrix(), and initializeFrom().

void oofem::StructuralMaterial::giveEshelbyStressVector_PlaneStrain ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedF,
TimeStep tStep 
)
virtual

Prototype for computation of Eshelby stress.

No default implementation is provided.

Reimplemented in oofem::LinearElasticMaterial.

Definition at line 611 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by giveCauchyStressVector_1d(), and oofem::SimpleCrossSection::giveEshelbyStresses().

void oofem::StructuralMaterial::giveFiberStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing 1d fiber stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 1d fiber stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::MPlasticMaterial2, oofem::MPlasticMaterial, and oofem::PlasticMaterial.

Definition at line 1192 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), and give3dMaterialStiffnessMatrix().

Referenced by give3dMaterialStiffnessMatrix(), oofem::MPlasticMaterial2::giveFiberStiffMtrx(), oofem::LinearElasticMaterial::giveRealStressVector_Fiber(), and giveStiffnessMatrix().

void oofem::StructuralMaterial::giveFirstPKStressVector_PlaneStrain ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedF,
TimeStep tStep 
)
virtual

Default implementation relies on giveFirstPKStressVector_3d.

Reimplemented in oofem::MooneyRivlinMaterial.

Definition at line 343 of file structuralmaterial.C.

References giveFirstPKStressVector_3d(), giveFullVectorFormF(), and giveReducedVectorForm().

Referenced by giveClassName(), and oofem::SimpleCrossSection::giveFirstPKStresses().

void oofem::StructuralMaterial::giveFullSymMatrixForm ( FloatMatrix answer,
const FloatMatrix red,
MaterialMode  matMode 
)
static
void oofem::StructuralMaterial::giveFullSymVectorForm ( FloatArray answer,
const FloatArray vec,
MaterialMode  matMode 
)
static

Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form.

Definition at line 2449 of file structuralmaterial.C.

References oofem::FloatArray::assemble(), oofem::FloatArray::giveSize(), giveVoigtSymVectorMask(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().

Referenced by oofem::RCSDEMaterial::computeCurrEquivStrain(), oofem::RCSDMaterial::computeCurrEquivStrain(), oofem::MDM::computeDamageOnPlane(), oofem::AnisotropicDamageMaterial::computeDamageTensor(), oofem::MazarsMaterial::computeEquivalentStrain(), oofem::IsotropicDamageMaterial1::computeEquivalentStrain(), oofem::AnisotropicDamageMaterial::computeEquivalentStrain(), oofem::B3Material::computeShrinkageStrainVector(), oofem::J2plasticMaterial::computeTrialStressIncrement(), oofem::MPlasticMaterial::computeTrialStressIncrement(), oofem::MPlasticMaterial2::computeTrialStressIncrement(), oofem::StructuralFE2Material::give3dMaterialStiffnessMatrix(), give3dMaterialStiffnessMatrix(), oofem::PlasticMaterial::giveConsistentStiffnessMatrix(), oofem::MPlasticMaterial::giveConsistentStiffnessMatrix(), oofem::MPlasticMaterial2::giveConsistentStiffnessMatrix(), oofem::MPlasticMaterial::giveElastoPlasticStiffnessMatrix(), oofem::MPlasticMaterial2::giveElastoPlasticStiffnessMatrix(), oofem::LinearElasticMaterial::giveEshelbyStressVector_PlaneStrain(), oofem::PlasticMaterial::giveIPValue(), oofem::PerfectlyPlasticMaterial::giveIPValue(), oofem::MPlasticMaterial::giveIPValue(), oofem::MPlasticMaterial2::giveIPValue(), oofem::FCMMaterial::giveIPValue(), oofem::MPSMaterial::giveIPValue(), giveIPValue(), oofem::PerfectlyPlasticMaterial::giveMaterialStiffnessMatrix(), oofem::MDM::giveRawMDMParameters(), oofem::RCSDNLMaterial::giveRealStressVector(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::PerfectlyPlasticMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), oofem::RCM2Material::giveRealStressVector(), oofem::AnisotropicDamageMaterial::giveRealStressVector(), oofem::StructuralFE2Material::giveRealStressVector_3d(), giveRealStressVector_PlaneStrain(), oofem::Concrete2::giveRealStressVector_PlateLayer(), oofem::IDNLMaterial::giveRemoteNonlocalStiffnessContribution(), oofem::IsotropicDamageMaterial1::initDamaged(), oofem::RCSDNLMaterialStatus::printOutputAt(), oofem::StructuralMaterialStatus::printOutputAt(), oofem::AnisotropicDamageMaterialStatus::printOutputAt(), oofem::M1MaterialStatus::restoreContext(), and oofem::MDM::transformStrainToPDC().

void oofem::StructuralMaterial::giveFullVectorForm ( FloatArray answer,
const FloatArray strainVector,
MaterialMode  matMode 
)
static

Converts the reduced symmetric Voigt vector (2nd order tensor) to full form.

Definition at line 2464 of file structuralmaterial.C.

References oofem::FloatArray::assemble(), giveVoigtVectorMask(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().

Referenced by oofem::MaterialForceEvaluator::computeMaterialForce(), and give3dMaterialStiffnessMatrix().

void oofem::StructuralMaterial::giveFullVectorFormF ( FloatArray answer,
const FloatArray strainVector,
MaterialMode  matMode 
)
static
void oofem::StructuralMaterial::giveInputRecord ( DynamicInputRecord input)
virtual

Setups the input record string of receiver.

Parameters
inputDynamic input record to be filled by receiver.

Reimplemented from oofem::Material.

Reimplemented in oofem::AnisotropicDamageMaterial, oofem::MDM, oofem::IsotropicDamageMaterial1, oofem::IsotropicDamageMaterial, oofem::MicroplaneMaterial, oofem::CompoDamageMat, oofem::AbaqusUserMaterial, oofem::StructuralFE2Material, oofem::OrthotropicLinearElasticMaterial, oofem::IDNLMaterial, oofem::IsotropicLinearElasticMaterial, oofem::MisesMatNl, oofem::StructuralPythonMaterial, oofem::RankineMatNl, oofem::TrabBoneNLEmbed, oofem::TrabBoneNL3D, oofem::TrabBoneNL, oofem::WinklerMaterial, oofem::WinklerPasternakMaterial, oofem::AnisotropicLinearElasticMaterial, oofem::SimpleVitrificationMaterial, oofem::LinearElasticMaterial, oofem::TutorialMaterial, and oofem::J2plasticMaterial.

Definition at line 2575 of file structuralmaterial.C.

References _IFT_StructuralMaterial_referencetemperature, oofem::Material::giveInputRecord(), referenceTemperature, and oofem::DynamicInputRecord::setField().

Referenced by giveClassName(), oofem::J2plasticMaterial::giveInputRecord(), oofem::TutorialMaterial::giveInputRecord(), oofem::LinearElasticMaterial::giveInputRecord(), oofem::SimpleVitrificationMaterial::giveInputRecord(), oofem::WinklerPasternakMaterial::giveInputRecord(), oofem::WinklerMaterial::giveInputRecord(), oofem::TrabBoneNL::giveInputRecord(), oofem::TrabBoneNL3D::giveInputRecord(), oofem::TrabBoneNLEmbed::giveInputRecord(), oofem::RankineMatNl::giveInputRecord(), oofem::StructuralPythonMaterial::giveInputRecord(), oofem::MisesMatNl::giveInputRecord(), oofem::IsotropicLinearElasticMaterial::giveInputRecord(), oofem::StructuralFE2Material::giveInputRecord(), oofem::AbaqusUserMaterial::giveInputRecord(), oofem::CompoDamageMat::giveInputRecord(), oofem::MicroplaneMaterial::giveInputRecord(), oofem::IsotropicDamageMaterial::giveInputRecord(), oofem::MDM::giveInputRecord(), and oofem::AnisotropicDamageMaterial::giveInputRecord().

void oofem::StructuralMaterial::giveInvertedVoigtVectorMask ( IntArray answer,
MaterialMode  mmode 
)
static
int oofem::StructuralMaterial::giveIPValue ( FloatArray answer,
GaussPoint gp,
InternalStateType  type,
TimeStep tStep 
)
virtual

Returns the integration point corresponding value in Reduced form.

Parameters
answerContain corresponding ip value, zero sized if not available.
gpIntegration point to which the value refers.
typeDetermines the type of internal variable.
tStepDetermines the time step.
Returns
Nonzero if the assignment can be done, zero if this type of variable is not supported.
Todo:
What about the stress meassure in large deformations here? The internal state type should specify "Cauchy" or something.
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!

Reimplemented from oofem::Material.

Reimplemented in oofem::ConcreteDPM2, oofem::ConcreteDPM, oofem::DustMaterial, oofem::LatticeDamage2d, oofem::DruckerPragerPlasticitySM, oofem::MDM, oofem::MPSMaterial, oofem::AnisotropicDamageMaterial, oofem::FCMMaterial, oofem::MPSDamMaterial, oofem::TrabBone3D, oofem::MPlasticMaterial2, oofem::RCM2Material, oofem::RheoChainMaterial, oofem::IsotropicDamageMaterial, oofem::MPlasticMaterial, oofem::RankineMat, oofem::RankineMatNl, oofem::SteelRelaxMat, oofem::PerfectlyPlasticMaterial, oofem::CompoDamageMat, oofem::PlasticMaterial, oofem::AbaqusUserMaterial, oofem::IDNLMaterial, oofem::DruckerPragerCutMat, oofem::MisesMat, oofem::TrabBoneEmbed, oofem::FRCFCMNL, oofem::RankineMatGrad, oofem::StructuralPythonMaterial, oofem::FRCFCM, oofem::M1Material, oofem::ConcreteFCM, oofem::LargeStrainMasterMaterial, oofem::LinearElasticMaterial, and oofem::TutorialMaterial.

Definition at line 2115 of file structuralmaterial.C.

References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::beColumnOf(), computePrincipalValDir(), computePrincipalValues(), oofem::StructuralElement::computeResultingIPEigenstrainAt(), computeVonMisesStress(), oofem::FEMComponent::domain, gc, oofem::EngngModel::giveContext(), oofem::GaussPoint::giveElement(), oofem::Domain::giveEngngModel(), oofem::FieldManager::giveField(), oofem::EngngModelContext::giveFieldManager(), giveFullSymVectorForm(), oofem::StructuralMaterialStatus::giveFVector(), oofem::Material::giveIPValue(), oofem::GaussPoint::giveMaterialMode(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FEMComponent::giveNumber(), oofem::FloatMatrix::giveNumberOfColumns(), oofem::StructuralMaterialStatus::givePVector(), oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::giveStrainVector(), oofem::StructuralMaterialStatus::giveStressVector(), oofem::StructuralMaterialStatus::giveTempStrainVector(), oofem::StructuralMaterialStatus::giveTempStressVector(), OOFEM_ERROR, oofem::principal_strain, oofem::principal_stress, oofem::FloatArray::resize(), transformStrainVectorTo(), transformStressVectorTo(), and oofem::FloatArray::zero().

Referenced by give3dMaterialStiffnessMatrix(), oofem::TutorialMaterial::giveIPValue(), oofem::LinearElasticMaterial::giveIPValue(), oofem::LargeStrainMasterMaterial::giveIPValue(), oofem::M1Material::giveIPValue(), oofem::StructuralPythonMaterial::giveIPValue(), oofem::TrabBoneEmbed::giveIPValue(), oofem::MisesMat::giveIPValue(), oofem::AbaqusUserMaterial::giveIPValue(), oofem::PlasticMaterial::giveIPValue(), oofem::CompoDamageMat::giveIPValue(), oofem::PerfectlyPlasticMaterial::giveIPValue(), oofem::SteelRelaxMat::giveIPValue(), oofem::RankineMat::giveIPValue(), oofem::MPlasticMaterial::giveIPValue(), oofem::IsotropicDamageMaterial::giveIPValue(), oofem::RheoChainMaterial::giveIPValue(), oofem::RCM2Material::giveIPValue(), oofem::MPlasticMaterial2::giveIPValue(), oofem::TrabBone3D::giveIPValue(), oofem::MPSDamMaterial::giveIPValue(), oofem::FCMMaterial::giveIPValue(), oofem::AnisotropicDamageMaterial::giveIPValue(), oofem::MDM::giveIPValue(), oofem::DruckerPragerPlasticitySM::giveIPValue(), oofem::LatticeDamage2d::giveIPValue(), oofem::DustMaterial::giveIPValue(), oofem::ConcreteDPM::giveIPValue(), oofem::ConcreteDPM2::giveIPValue(), and oofem::M1MaterialStatus::restoreContext().

void oofem::StructuralMaterial::givePlaneStrainStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing plane strain stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane strain stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix. Note: as already described, if zero strain component is imposed (Plane strain, ..) this condition must be taken into account in geometrical relations, and corresponding component has to be included in reduced vector. (So plane strain conditions are $ \epsilon_z = \gamma_{xz} = \gamma_{yz} = 0 $, but relations for $ \epsilon_z$ and $\sigma_z$ are included).

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::MPlasticMaterial2, oofem::MDM, oofem::MPSDamMaterial, oofem::AnisotropicDamageMaterial, oofem::IsotropicDamageMaterial, oofem::RCM2Material, oofem::MPlasticMaterial, oofem::RheoChainMaterial, oofem::PlasticMaterial, oofem::FCMMaterial, oofem::PerfectlyPlasticMaterial, oofem::IsotropicLinearElasticMaterial, oofem::MisesMatGrad, and oofem::IDGMaterial.

Definition at line 1075 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), give3dMaterialStiffnessMatrix(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().

Referenced by give3dMaterialStiffnessMatrix(), oofem::SimpleCrossSection::giveCharMaterialStiffnessMatrix(), oofem::RheoChainMaterial::givePlaneStrainStiffMtrx(), oofem::MPlasticMaterial2::givePlaneStrainStiffMtrx(), givePlaneStrainStiffMtrx_dPdF(), oofem::LinearElasticMaterial::giveRealStressVector_PlaneStrain(), giveStiffnessMatrix(), and oofem::SimpleCrossSection::giveStiffnessMatrix_PlaneStrain().

void oofem::StructuralMaterial::givePlaneStrainStiffMtrx_dCde ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::givePlaneStrainStiffMtrx_dPdF ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::givePlaneStressStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing plane stress stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane stress stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::MPlasticMaterial2, oofem::MDM, oofem::MPSDamMaterial, oofem::AnisotropicDamageMaterial, oofem::IsotropicDamageMaterial, oofem::RCM2Material, oofem::MPlasticMaterial, oofem::RheoChainMaterial, oofem::PlasticMaterial, oofem::FCMMaterial, oofem::PerfectlyPlasticMaterial, oofem::RankineMat, oofem::IsotropicLinearElasticMaterial, oofem::RankineMatNl, oofem::RankineMatGrad, and oofem::IDGMaterial.

Definition at line 1042 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), give3dMaterialStiffnessMatrix(), and oofem::FloatMatrix::resize().

Referenced by oofem::SimpleCrossSection::give2dPlateStiffMtrx(), give3dMaterialStiffnessMatrix(), oofem::SimpleCrossSection::give3dShellStiffMtrx(), oofem::SimpleCrossSection::giveCharMaterialStiffnessMatrix(), oofem::SimpleCrossSection::giveMembraneRotStiffMtrx(), oofem::RheoChainMaterial::givePlaneStressStiffMtrx(), oofem::MPlasticMaterial2::givePlaneStressStiffMtrx(), givePlaneStressStiffMtrx_dPdF(), oofem::LinearElasticMaterial::giveRealStressVector_PlaneStress(), giveStiffnessMatrix(), oofem::SimpleCrossSection::giveStiffnessMatrix_PlaneStress(), and oofem::M1MaterialStatus::restoreContext().

void oofem::StructuralMaterial::givePlaneStressStiffMtrx_dCde ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::givePlaneStressStiffMtrx_dPdF ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::givePlaneStressVectorTranformationMtrx ( FloatMatrix answer,
const FloatMatrix base,
bool  transpose = false 
)
static

Computes 2d stress vector transformation matrix from standard vector transformation matrix.

Parameters
answerTransformation matrix for stress vector.
baseA (2,2) matrix, where on each column are stored unit direction vectors of local coordinate axes to which we do transformation.
transposeDetermines if we transpose matrix before transforming.

Definition at line 1948 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beTranspositionOf(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().

Referenced by oofem::RankineMat::computeEta(), oofem::CrackExportModule::doOutput(), oofem::RankineMat::evaluatePlaneStressStiffMtrx(), give3dMaterialStiffnessMatrix(), oofem::RankineMatGrad::givePlaneStressKappaMatrix(), oofem::MPSDamMaterial::giveRealStressVector(), and oofem::FCMMaterial::initializeCrack().

void oofem::StructuralMaterial::givePlateLayerStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing 2d plate layer stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d plate layer stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::MPlasticMaterial2, oofem::RCM2Material, oofem::MPlasticMaterial, oofem::PlasticMaterial, oofem::PerfectlyPlasticMaterial, and oofem::Concrete2.

Definition at line 1154 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), and give3dMaterialStiffnessMatrix().

Referenced by oofem::LayeredCrossSection::give2dPlateStiffMtrx(), give3dMaterialStiffnessMatrix(), oofem::LayeredCrossSection::give3dShellStiffMtrx(), oofem::Concrete2::givePlateLayerStiffMtrx(), oofem::MPlasticMaterial2::givePlateLayerStiffMtrx(), oofem::LinearElasticMaterial::giveRealStressVector_PlateLayer(), and giveStiffnessMatrix().

void oofem::StructuralMaterial::giveRealStressVector ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
virtual

Computes the real stress vector for given total strain and integration point.

The total strain is defined as strain computed directly from displacement field at given time. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process.

Parameters
answerStress vector in reduced form. For large deformations it is treated as the second Piola-Kirchoff stress.
gpIntegration point.
reducedStrainStrain vector in reduced form. For large deformations it is treated as the Green-Lagrange strain.
tStepCurrent time step (most models are able to respond only when tStep is current time step).
Todo:
Move this to StructuralCrossSection ?

Reimplemented in oofem::LatticeDamage2d, oofem::MPSMaterial, oofem::MDM, oofem::AnisotropicDamageMaterial, oofem::MPSDamMaterial, oofem::FCMMaterial, oofem::MPlasticMaterial2, oofem::RCM2Material, oofem::IsotropicDamageMaterial, oofem::MPlasticMaterial, oofem::Eurocode2CreepMaterial, oofem::RheoChainMaterial, oofem::CompoDamageMat, oofem::PlasticMaterial, oofem::RCSDMaterial, oofem::SteelRelaxMat, oofem::PerfectlyPlasticMaterial, oofem::RCSDEMaterial, oofem::RCSDNLMaterial, oofem::B3SolidMaterial, oofem::FRCFCMNL, oofem::KelvinChainMaterial, oofem::MaxwellChainMaterial, and oofem::KelvinChainSolidMaterial.

Definition at line 79 of file structuralmaterial.C.

References oofem::GaussPoint::giveMaterialMode(), giveRealStressVector_1d(), giveRealStressVector_2dBeamLayer(), giveRealStressVector_2dPlateSubSoil(), giveRealStressVector_3d(), giveRealStressVector_3dBeamSubSoil(), giveRealStressVector_Fiber(), giveRealStressVector_Lattice2d(), giveRealStressVector_Lattice3d(), giveRealStressVector_PlaneStrain(), giveRealStressVector_PlaneStress(), and giveRealStressVector_PlateLayer().

Referenced by giveClassName(), and oofem::StructuralCrossSection::giveRealStresses().

void oofem::StructuralMaterial::giveRealStressVector_2dBeamLayer ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::giveRealStressVector_2dPlateSubSoil ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual

Default implementation is not provided.

Reimplemented in oofem::WinklerMaterial, and oofem::WinklerPasternakMaterial.

Definition at line 297 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by giveClassName(), oofem::SimpleCrossSection::giveGeneralizedStress_PlateSubSoil(), and giveRealStressVector().

void oofem::StructuralMaterial::giveRealStressVector_3d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::giveRealStressVector_3dBeamSubSoil ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::giveRealStressVector_Fiber ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual

Default implementation relies on giveRealStressVector_StressControl.

Reimplemented in oofem::PlasticMaterial, and oofem::LinearElasticMaterial.

Definition at line 275 of file structuralmaterial.C.

References giveRealStressVector_StressControl(), and giveVoigtSymVectorMask().

Referenced by giveClassName(), oofem::FiberedCrossSection::giveGeneralizedStress_Beam3d(), and giveRealStressVector().

void oofem::StructuralMaterial::giveRealStressVector_Lattice2d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual

Definition at line 284 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by giveClassName(), and giveRealStressVector().

void oofem::StructuralMaterial::giveRealStressVector_Lattice3d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual

Definition at line 290 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by giveClassName(), and giveRealStressVector().

void oofem::StructuralMaterial::giveRealStressVector_PlateLayer ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::giveRealStressVector_StressControl ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
const IntArray strainControl,
TimeStep tStep 
)
virtual
void oofem::StructuralMaterial::giveRealStressVector_Warping ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
virtual

Default implementation relies on giveRealStressVector_StressControl.

Reimplemented in oofem::LinearElasticMaterial.

Definition at line 119 of file structuralmaterial.C.

References OOFEM_ERROR.

Referenced by giveClassName(), and oofem::SimpleCrossSection::giveRealStress_Warping().

void oofem::StructuralMaterial::giveReducedMatrixForm ( FloatMatrix answer,
const FloatMatrix full,
MaterialMode  matMode 
)
static

Converts the full symmetric Voigt matrix (4th order tensor) to reduced form.

Definition at line 2522 of file structuralmaterial.C.

References oofem::FloatMatrix::beSubMatrixOf(), and giveVoigtVectorMask().

Referenced by give3dMaterialStiffnessMatrix(), and oofem::AbaqusUserMaterial::givePlaneStrainStiffMtrx_dPdF().

void oofem::StructuralMaterial::giveReducedSymVectorForm ( FloatArray answer,
const FloatArray vec,
MaterialMode  matMode 
)
static

Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form.

Definition at line 2497 of file structuralmaterial.C.

References oofem::FloatArray::beSubArrayOf(), oofem::IntArray::giveSize(), oofem::FloatArray::giveSize(), and giveVoigtSymVectorMask().

Referenced by oofem::MPSMaterial::computeB4AutogenousShrinkageStrainVector(), oofem::MPSMaterial::computeFibAutogenousShrinkageStrainVector(), oofem::PlasticMaterial::ComputeGradientVector(), oofem::MPlasticMaterial::computeGradientVector(), oofem::Eurocode2CreepMaterial::computeIncrementOfAutogenousShrinkageVector(), oofem::Eurocode2CreepMaterial::computeIncrementOfDryingShrinkageVector(), oofem::B3SolidMaterial::computePointShrinkageStrainVector(), oofem::MPSMaterial::computePointShrinkageStrainVector(), oofem::MPlasticMaterial2::computeReducedStressGradientVector(), oofem::B3Material::computeShrinkageStrainVector(), computeStressIndependentStrainVector(), oofem::J2MPlasticMaterial::computeStressSpaceHardeningVarsReducedGradient(), oofem::J2plasticMaterial::ComputeStressSpaceHardeningVarsReducedGradient(), oofem::B3Material::computeTotalAverageShrinkageStrainVector(), oofem::B3SolidMaterial::computeTotalAverageShrinkageStrainVector(), oofem::MPlasticMaterial::cuttingPlaneReturn(), give3dMaterialStiffnessMatrix(), give_dPdF_from(), oofem::NLStructuralElement::giveInternalForcesVector(), oofem::StructuralElement::giveInternalForcesVector(), oofem::PerfectlyPlasticMaterial::giveRealStressVector(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), oofem::PlasticMaterial::giveRealStressVector(), oofem::MPlasticMaterial::giveRealStressVector(), oofem::RCM2Material::giveRealStressVector(), oofem::MPlasticMaterial2::giveRealStressVector(), oofem::AnisotropicDamageMaterial::giveRealStressVector(), giveRealStressVector_PlaneStrain(), oofem::Concrete2::giveRealStressVector_PlateLayer(), oofem::IDNLMaterial::giveRemoteNonlocalStiffnessContribution(), oofem::MPSMaterial::giveShrinkageStrainVector(), oofem::IsotropicDamageMaterial1::MMI_map(), oofem::MDM::transformStressFromPDC(), and oofem::ZZErrorEstimatorInterface::ZZErrorEstimatorI_computeElementContributions().

int oofem::StructuralMaterial::giveSizeOfVoigtSymVector ( MaterialMode  mmode)
static

Returns the size of symmetric part of a reduced stress/strain vector according to given mode.

Parameters
mmodeMaterial response mode.

Definition at line 801 of file structuralmaterial.C.

References oofem::IntArray::giveSize(), and giveVoigtSymVectorMask().

Referenced by oofem::J2MPlasticMaterial::computeHardeningReducedModuli(), oofem::J2plasticMaterial::computeHardeningReducedModuli(), oofem::KelvinChainSolidMaterial::computeHiddenVars(), oofem::KelvinChainMaterial::computeHiddenVars(), oofem::J2MPlasticMaterial::computeReducedGradientMatrix(), oofem::J2plasticMaterial::computeReducedGradientMatrix(), oofem::J2Mat::computeReducedHardeningVarsSigmaGradient(), oofem::DruckerPragerCutMat::computeReducedHardeningVarsSigmaGradient(), oofem::J2Mat::computeReducedSKGradientMatrix(), oofem::DruckerPragerCutMat::computeReducedSKGradientMatrix(), oofem::J2Mat::computeReducedSSGradientMatrix(), oofem::B3Material::computeShrinkageStrainVector(), oofem::StructuralElementEvaluator::computeStrainVector(), oofem::StructuralElement::computeStrainVector(), computeStressIndependentStrainVector(), oofem::StressStrainBaseVector::convertFromFullForm(), give3dMaterialStiffnessMatrix(), oofem::NLStructuralElement::giveInternalForcesVector(), oofem::StructuralElement::giveInternalForcesVector(), oofem::StructuralElement::giveInternalForcesVector_withIRulesAsSubcells(), oofem::RheoChainMaterial::giveRealStressVector(), oofem::MPSDamMaterial::giveRealStressVector(), oofem::B3Material::giveShrinkageStrainVector(), oofem::B3SolidMaterial::giveShrinkageStrainVector(), oofem::Eurocode2CreepMaterial::giveShrinkageStrainVector(), oofem::MPSMaterial::giveShrinkageStrainVector(), oofem::J2Mat::giveSizeOfReducedHardeningVarsVector(), oofem::J2MPlasticMaterial::giveSizeOfReducedHardeningVarsVector(), oofem::J2plasticMaterial::giveSizeOfReducedHardeningVarsVector(), oofem::MisesMatGradStatus::initTempStatus(), oofem::RankineMatGradStatus::initTempStatus(), oofem::PlasticMaterialStatus::initTempStatus(), oofem::RCSDNLMaterialStatus::initTempStatus(), oofem::PerfectlyPlasticMaterialStatus::initTempStatus(), oofem::MPlasticMaterialStatus::initTempStatus(), oofem::StructuralMaterialStatus::initTempStatus(), oofem::MPlasticMaterial2Status::initTempStatus(), oofem::Concrete2MaterialStatus::initTempStatus(), oofem::RankineMatStatus::initTempStatus(), oofem::StressStrainBaseVector::letStressStrainModeBe(), oofem::MPSDamMaterialStatus::MPSDamMaterialStatus(), oofem::MPSMaterialStatus::MPSMaterialStatus(), oofem::RCSDNLMaterialStatus::RCSDNLMaterialStatus(), oofem::StressStrainBaseVector::StressStrainBaseVector(), and oofem::StructuralMaterialStatus::StructuralMaterialStatus().

int oofem::StructuralMaterial::giveSizeOfVoigtVector ( MaterialMode  mmode)
static

Returns the size of reduced stress/strain vector according to given mode.

Parameters
mmodeMaterial response mode.

Definition at line 810 of file structuralmaterial.C.

References oofem::IntArray::giveSize(), and giveVoigtVectorMask().

Referenced by give3dMaterialStiffnessMatrix().

void oofem::StructuralMaterial::giveStiffnessMatrix ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Computes the stiffness matrix for giveRealStressVector of receiver in given integration point, respecting its history.

The algorithm should use temporary or equilibrium history variables stored in integration point status to compute and return required result.

Parameters
answerContains result.
modeMaterial response mode.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::LatticeDamage2d, oofem::Masonry02, oofem::TrabBoneGrad3D, oofem::RankineMatGrad, oofem::MisesMatGrad, oofem::IDGMaterial, and oofem::LargeStrainMasterMaterialGrad.

Definition at line 638 of file structuralmaterial.C.

References oofem::__MaterialModeToString(), oofem::FloatMatrix::beUnitMatrix(), give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), give2dLatticeStiffMtrx(), give3dLatticeStiffMtrx(), give3dMaterialStiffnessMatrix(), giveFiberStiffMtrx(), oofem::GaussPoint::giveMaterialMode(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), givePlateLayerStiffMtrx(), OOFEM_ERROR, and oofem::FloatMatrix::resize().

Referenced by oofem::RCSDEMaterial::computeCurrEquivStrain(), oofem::RCSDMaterial::computeCurrEquivStrain(), oofem::MDM::computeEffectiveStress(), oofem::IsotropicDamageMaterial1::computeEquivalentStrain(), oofem::IsotropicDamageMaterial1::computeEta(), oofem::DruckerPragerCutMat::computeReducedElasticModuli(), oofem::Masonry02::computeReducedElasticModuli(), oofem::PlasticMaterial::computeReducedElasticModuli(), oofem::MPlasticMaterial::computeReducedElasticModuli(), oofem::MPlasticMaterial2::computeReducedElasticModuli(), oofem::J2plasticMaterial::computeTrialStressIncrement(), oofem::RankineMat::evaluatePlaneStressStiffMtrx(), oofem::PerfectlyPlasticMaterial::give1dStressStiffMtrx(), oofem::PlasticMaterial::give1dStressStiffMtrx(), oofem::MPlasticMaterial::give1dStressStiffMtrx(), oofem::IsotropicDamageMaterial::give1dStressStiffMtrx(), oofem::PerfectlyPlasticMaterial::give2dBeamLayerStiffMtrx(), oofem::PlasticMaterial::give2dBeamLayerStiffMtrx(), oofem::MPlasticMaterial::give2dBeamLayerStiffMtrx(), oofem::AnisotropicDamageMaterial::give3dMaterialStiffnessMatrix(), oofem::SimpleCrossSection::giveCharMaterialStiffnessMatrix(), oofem::LayeredCrossSection::giveCharMaterialStiffnessMatrix(), giveClassName(), oofem::MPlasticMaterial::giveConsistentStiffnessMatrix(), oofem::MPlasticMaterial2::giveConsistentStiffnessMatrix(), oofem::RCSDEMaterial::giveEffectiveMaterialStiffnessMatrix(), oofem::RCSDMaterial::giveEffectiveMaterialStiffnessMatrix(), oofem::PerfectlyPlasticMaterial::giveEffectiveMaterialStiffnessMatrix(), oofem::RCM2Material::giveEffectiveMaterialStiffnessMatrix(), oofem::MPlasticMaterial::giveElastoPlasticStiffnessMatrix(), oofem::MPlasticMaterial2::giveElastoPlasticStiffnessMatrix(), oofem::FiberedCrossSection::giveFiberMaterialStiffnessMatrix(), oofem::PlasticMaterial::giveFiberStiffMtrx(), oofem::MPlasticMaterial::giveFiberStiffMtrx(), oofem::IDNLMaterial::giveLocalNonlocalStiffnessContribution(), oofem::MDM::giveMaterialStiffnessMatrix(), oofem::FCMMaterial::giveMaterialStiffnessMatrix(), oofem::RCM2Material::giveNormalElasticStiffnessMatrix(), oofem::IDGMaterial::givePlaneStrainStiffMtrx(), oofem::MisesMatGrad::givePlaneStrainStiffMtrx(), oofem::PerfectlyPlasticMaterial::givePlaneStrainStiffMtrx(), oofem::PlasticMaterial::givePlaneStrainStiffMtrx(), oofem::MPlasticMaterial::givePlaneStrainStiffMtrx(), oofem::IsotropicDamageMaterial::givePlaneStrainStiffMtrx(), oofem::IDGMaterial::givePlaneStressStiffMtrx(), oofem::RankineMatNl::givePlaneStressStiffMtrx(), oofem::PerfectlyPlasticMaterial::givePlaneStressStiffMtrx(), oofem::PlasticMaterial::givePlaneStressStiffMtrx(), oofem::MPlasticMaterial::givePlaneStressStiffMtrx(), oofem::IsotropicDamageMaterial::givePlaneStressStiffMtrx(), oofem::AnisotropicDamageMaterial::givePlaneStressStiffMtrx(), oofem::PerfectlyPlasticMaterial::givePlateLayerStiffMtrx(), oofem::PlasticMaterial::givePlateLayerStiffMtrx(), oofem::MPlasticMaterial::givePlateLayerStiffMtrx(), oofem::RheoChainMaterial::giveRealStressVector(), oofem::IsotropicDamageMaterial::giveRealStressVector(), oofem::FCMMaterial::giveRealStressVector(), oofem::AnisotropicDamageMaterial::giveRealStressVector(), oofem::IDGMaterial::giveRealStressVectorGrad(), oofem::IDNLMaterial::giveRemoteNonlocalStiffnessContribution(), oofem::Masonry02::giveStiffnessMatrix(), oofem::RheoChainMaterial::giveUnitComplianceMatrix(), oofem::RheoChainMaterial::giveUnitStiffnessMatrix(), and oofem::IsotropicDamageMaterial1::initDamaged().

void oofem::StructuralMaterial::giveStrainVectorTranformationMtrx ( FloatMatrix answer,
const FloatMatrix base,
bool  transpose = false 
)
static

Computes 3d strain vector transformation matrix from standard vector transformation matrix.

Parameters
answerTransformation matrix for strain vector.
baseA (3,3) matrix, where on each column are stored unit direction vectors of local coordinate axes to which we do transformation.
transposeDetermines if we transpose matrix before transforming.

Definition at line 1777 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beTranspositionOf(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().

Referenced by give3dMaterialStiffnessMatrix(), oofem::CompoDamageMat::giveMatStiffRotationMatrix(), oofem::OrthotropicLinearElasticMaterial::giveRotationMatrix(), oofem::FCMMaterial::initializeCrack(), and transformStrainVectorTo().

void oofem::StructuralMaterial::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).

Calls computeStressIndependentStrainVector to obtain stress independent part of strain.

Parameters
answerComputed strain vector.
gpIntegration point.
reducedStrainVectorReduced strain vector.
tStepTime step.
modeDetermines value mode.

Definition at line 768 of file structuralmaterial.C.

References computeStressIndependentStrainVector(), oofem::FloatArray::giveSize(), and oofem::FloatArray::subtract().

Referenced by oofem::SteelRelaxMat::computeStressRelaxationStrainVector(), oofem::AnisotropicDamageMaterial::give3dMaterialStiffnessMatrix(), give3dMaterialStiffnessMatrix(), oofem::LinearElasticMaterial::giveIPValue(), oofem::IsotropicDamageMaterial::giveIPValue(), oofem::AnisotropicDamageMaterial::givePlaneStressStiffMtrx(), oofem::RCSDNLMaterial::giveRealStressVector(), oofem::PerfectlyPlasticMaterial::giveRealStressVector(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::SteelRelaxMat::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), oofem::PlasticMaterial::giveRealStressVector(), oofem::CompoDamageMat::giveRealStressVector(), oofem::RheoChainMaterial::giveRealStressVector(), oofem::MPlasticMaterial::giveRealStressVector(), oofem::IsotropicDamageMaterial::giveRealStressVector(), oofem::RCM2Material::giveRealStressVector(), oofem::MPlasticMaterial2::giveRealStressVector(), oofem::FCMMaterial::giveRealStressVector(), oofem::AnisotropicDamageMaterial::giveRealStressVector(), oofem::MDM::giveRealStressVector(), oofem::LatticeDamage2d::giveRealStressVector(), oofem::LinearElasticMaterial::giveRealStressVector_1d(), oofem::MisesMat::giveRealStressVector_1d(), oofem::LinearElasticMaterial::giveRealStressVector_2dBeamLayer(), oofem::SimpleVitrificationMaterial::giveRealStressVector_3d(), oofem::MisesMat::giveRealStressVector_3d(), oofem::LinearElasticMaterial::giveRealStressVector_3dDegeneratedShell(), oofem::LinearElasticMaterial::giveRealStressVector_Fiber(), oofem::LinearElasticMaterial::giveRealStressVector_PlaneStrain(), oofem::LinearElasticMaterial::giveRealStressVector_PlaneStress(), oofem::AnisotropicDamageMaterial::giveRealStressVector_PlaneStress(), oofem::LinearElasticMaterial::giveRealStressVector_PlateLayer(), oofem::Concrete2::giveRealStressVector_PlateLayer(), oofem::IDGMaterial::giveRealStressVectorGrad(), oofem::SteelRelaxMat::giveStressDependentPartOfStrainVector(), oofem::TrabBoneNL::updateBeforeNonlocAverage(), oofem::TrabBoneNLEmbed::updateBeforeNonlocAverage(), oofem::MazarsNLMaterial::updateBeforeNonlocAverage(), oofem::IDNLMaterial::updateBeforeNonlocAverage(), and oofem::TrabBoneNL3D::updateBeforeNonlocAverage().

void oofem::StructuralMaterial::giveStressVectorTranformationMtrx ( FloatMatrix answer,
const FloatMatrix base,
bool  transpose = false 
)
static

Computes 3d stress vector transformation matrix from standard vector transformation matrix.

Parameters
answerTransformation matrix for stress vector.
baseA (3,3) matrix, where on each column are stored unit direction vectors of local coordinate axes to which we do transformation.
transposeDetermines if we transpose matrix before transforming.

Definition at line 1881 of file structuralmaterial.C.

References oofem::FloatMatrix::at(), oofem::FloatMatrix::beTranspositionOf(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().

Referenced by give3dMaterialStiffnessMatrix(), oofem::RCM2Material::giveEffectiveMaterialStiffnessMatrix(), oofem::MPSDamMaterial::giveRealStressVector(), oofem::IDNLMaterial::giveRemoteNonlocalStiffnessContribution(), oofem::FCMMaterial::initializeCrack(), and transformStressVectorTo().

static int oofem::StructuralMaterial::giveSymVI ( int  ind1,
int  ind2 
)
inlinestatic

Definition at line 126 of file structuralmaterial.h.

Referenced by convert_dSdE_2_dPdF().

static int oofem::StructuralMaterial::giveVI ( int  ind1,
int  ind2 
)
inlinestatic

Definition at line 127 of file structuralmaterial.h.

References StructuralMaterial().

Referenced by convert_dSdE_2_dPdF().

int oofem::StructuralMaterial::giveVoigtSymVectorMask ( IntArray answer,
MaterialMode  mmode 
)
static

The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor.

Returns a mask of the vector indicies corresponding to components in a symmetric second order tensor of some stress/strain/deformation measure that performes work. Thus, components corresponding to imposed zero stress (e.g. plane stress etc.) are not included. On the other hand, if zero strain components are imposed( e.g. plane strain etc.) this condition must be taken into account in geometrical relations. Therefore, these corresponding components are included in the reduced vector. Which compnents to include are given by the particular MaterialMode. Ex: PlaneStress -> [1 2 6]

Parameters
answerReturned mask.
mmodeMaterial response mode.
Returns
The number of components in the corresponding full vector.
Todo:
This isn't actually fixed yet. Should be made compatible with 3dShell and 2dBeam
Todo:
This isn't actually fixed yet. Should be made compatible with 3dShell and 2dBeam

Definition at line 831 of file structuralmaterial.C.

References oofem::IntArray::clear(), and oofem::IntArray::enumerate().

Referenced by oofem::RCM2Material::checkForNewActiveCracks(), oofem::J2MPlasticMaterial::computeStressSpaceHardeningVars(), oofem::J2plasticMaterial::ComputeStressSpaceHardeningVars(), oofem::StressStrainBaseVector::convertFromFullForm(), oofem::StressStrainBaseVector::convertToFullForm(), give3dMaterialStiffnessMatrix(), oofem::RCM2Material::giveEffectiveMaterialStiffnessMatrix(), giveFullSymMatrixForm(), giveFullSymVectorForm(), giveInvertedVoigtVectorMask(), oofem::FCMMaterial::giveLocalCrackedStiffnessMatrix(), oofem::FCMMaterialStatus::giveMaxNumberOfCracks(), oofem::RCM2Material::giveNormalElasticStiffnessMatrix(), oofem::RCSDNLMaterial::giveRealStressVector(), giveRealStressVector_1d(), giveRealStressVector_2dBeamLayer(), giveRealStressVector_Fiber(), giveRealStressVector_PlaneStress(), giveRealStressVector_PlateLayer(), giveReducedSymMatrixForm(), giveReducedSymVectorForm(), giveSizeOfVoigtSymVector(), oofem::J2Mat::giveStressBackVector(), and oofem::RCM2Material::updateCrackStatus().

int oofem::StructuralMaterial::giveVoigtVectorMask ( IntArray answer,
MaterialMode  mmode 
)
static

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.

Thus, components corresponding to imposed zero stress (e.g. plane stress etc.) are not included. On the other hand, if zero strain components are imposed( e.g. plane strain etc.) this condition must be taken into account in geometrical relations. Therefore, these corresponding components are included in the reduced vector. Which compnents to include are given by the particular MaterialMode. Ex: PlaneStress -> [1 2 6 9]

Parameters
answerReturned mask.
mmodeMaterial response mode.
Returns
The number of components in the corresponding full vector.
Todo:
add additional modes if they relevant.

Definition at line 992 of file structuralmaterial.C.

References oofem::IntArray::at(), and oofem::IntArray::resize().

Referenced by oofem::PrescribedGradientBCWeak::computeTangent(), give3dMaterialStiffnessMatrix(), giveFirstPKStressVector_PlaneStress(), giveFullVectorForm(), giveFullVectorFormF(), giveReducedMatrixForm(), giveReducedVectorForm(), and giveSizeOfVoigtVector().

IRResultType oofem::StructuralMaterial::initializeFrom ( InputRecord ir)
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.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
Returns
IRResultType

Reimplemented from oofem::Material.

Reimplemented in oofem::ConcreteDPM2, oofem::ConcreteDPM, oofem::DustMaterial, oofem::LatticeDamage2d, oofem::AnisotropicDamageMaterial, oofem::MDM, oofem::MPSMaterial, oofem::IsotropicDamageMaterial1, oofem::IsotropicDamageMaterial, oofem::MPSDamMaterial, oofem::DruckerPragerPlasticitySM, oofem::TrabBone3D, oofem::RheoChainMaterial, oofem::Concrete2, oofem::FCMMaterial, oofem::RCM2Material, oofem::MicroplaneMaterial, oofem::Eurocode2CreepMaterial, oofem::SteelRelaxMat, oofem::PerfectlyPlasticMaterial, oofem::RankineMat, oofem::CompoDamageMat, oofem::TrabBoneMaterial, oofem::AbaqusUserMaterial, oofem::B3SolidMaterial, oofem::StructuralFE2Material, oofem::OrthotropicLinearElasticMaterial, oofem::RCSDMaterial, oofem::RCSDEMaterial, oofem::TrabBoneEmbed, oofem::RCSDNLMaterial, oofem::IDNLMaterial, oofem::FRCFCM, oofem::MazarsMaterial, oofem::IsotropicLinearElasticMaterial, oofem::M4Material, oofem::M1Material, oofem::FRCFCMNL, oofem::MazarsNLMaterial, oofem::DruckerPragerCutMat, oofem::Masonry02, oofem::ConcreteFCM, oofem::B3Material, oofem::MisesMatNl, oofem::StructuralPythonMaterial, oofem::MisesMat, oofem::RankineMatNl, oofem::TrabBoneNLEmbed, oofem::MicroMaterial, oofem::RankineMatGrad, oofem::TrabBoneNL3D, oofem::TrabBoneNL, oofem::TrabBoneGrad3D, oofem::MisesMatGrad, oofem::WinklerMaterial, oofem::WinklerPasternakMaterial, oofem::MooneyRivlinMaterial, oofem::AnisotropicLinearElasticMaterial, oofem::MaxwellChainMaterial, oofem::SimpleVitrificationMaterial, oofem::LargeStrainMasterMaterial, oofem::KelvinChainSolidMaterial, oofem::KelvinChainMaterial, oofem::LinearElasticMaterial, oofem::CebFip78Material, oofem::Concrete3, oofem::TutorialMaterial, oofem::IDGMaterial, oofem::J2Mat, oofem::J2MPlasticMaterial, oofem::J2plasticMaterial, oofem::DoublePowerLawMaterial, oofem::StructuralMaterialSettable, oofem::RankinePlasticMaterial, oofem::HyperElasticMaterial, oofem::LargeStrainMasterMaterialGrad, and oofem::Steel1.

Definition at line 2554 of file structuralmaterial.C.

References _IFT_StructuralMaterial_referencetemperature, _IFT_StructuralMaterial_talpha, oofem::Dictionary::add(), oofem::Dictionary::includes(), oofem::Material::initializeFrom(), IR_GIVE_OPTIONAL_FIELD, oofem::Material::propertyDictionary, referenceTemperature, and tAlpha.

Referenced by giveClassName(), oofem::HyperElasticMaterial::initializeFrom(), oofem::RankinePlasticMaterial::initializeFrom(), oofem::J2Mat::initializeFrom(), oofem::J2plasticMaterial::initializeFrom(), oofem::StructuralMaterialSettable::initializeFrom(), oofem::J2MPlasticMaterial::initializeFrom(), oofem::TutorialMaterial::initializeFrom(), oofem::LinearElasticMaterial::initializeFrom(), oofem::SimpleVitrificationMaterial::initializeFrom(), oofem::MooneyRivlinMaterial::initializeFrom(), oofem::WinklerPasternakMaterial::initializeFrom(), oofem::WinklerMaterial::initializeFrom(), oofem::MisesMat::initializeFrom(), oofem::StructuralPythonMaterial::initializeFrom(), oofem::Masonry02::initializeFrom(), oofem::DruckerPragerCutMat::initializeFrom(), oofem::TrabBoneEmbed::initializeFrom(), oofem::StructuralFE2Material::initializeFrom(), oofem::AbaqusUserMaterial::initializeFrom(), oofem::TrabBoneMaterial::initializeFrom(), oofem::RankineMat::initializeFrom(), oofem::SteelRelaxMat::initializeFrom(), oofem::MicroplaneMaterial::initializeFrom(), oofem::FCMMaterial::initializeFrom(), oofem::RheoChainMaterial::initializeFrom(), oofem::TrabBone3D::initializeFrom(), oofem::DruckerPragerPlasticitySM::initializeFrom(), oofem::IsotropicDamageMaterial::initializeFrom(), oofem::MDM::initializeFrom(), oofem::AnisotropicDamageMaterial::initializeFrom(), oofem::LatticeDamage2d::initializeFrom(), oofem::DustMaterial::initializeFrom(), oofem::ConcreteDPM::initializeFrom(), oofem::ConcreteDPM2::initializeFrom(), and oofem::M1MaterialStatus::restoreContext().

int oofem::StructuralMaterial::setIPValue ( const FloatArray value,
GaussPoint gp,
InternalStateType  type 
)
virtual

Sets the value of a certain variable at a given integration point to the given value.

Parameters
valueContains the value(s) to be set (in reduced form).
gpIntegration point.
typeDetermines the type of internal variable.
typeDetermines the type of internal variable.
Returns
Nonzero if ok, zero if var not supported.

Reimplemented from oofem::Material.

Reimplemented in oofem::ConcreteDPM, and oofem::DustMaterial.

Definition at line 2093 of file structuralmaterial.C.

References oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::letStrainVectorBe(), oofem::StructuralMaterialStatus::letStressVectorBe(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), and oofem::StructuralMaterialStatus::letTempStressVectorBe().

Referenced by give3dMaterialStiffnessMatrix(), oofem::DustMaterial::setIPValue(), and oofem::ConcreteDPM::setIPValue().

void oofem::StructuralMaterial::sortPrincDirAndValCloseTo ( FloatArray pVal,
FloatMatrix pDir,
FloatMatrix toPDir 
)
static

Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDir) to be closed to some (often previous) principal directions (toPDir).

pDir and toPDir should have eigenvectors stored in columns and normalized.

Parameters
pValNew eigenvalues.
pDirNew eigenvectors.
toPDirOld eigenvector.

Definition at line 2025 of file structuralmaterial.C.

References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::giveNumberOfRows(), oofem::FloatArray::giveSize(), oofem::FloatMatrix::isSquare(), and OOFEM_SERROR.

Referenced by give3dMaterialStiffnessMatrix(), and oofem::RCM2Material::giveRealPrincipalStressVector3d().

void oofem::StructuralMaterial::transformStrainVectorTo ( FloatArray answer,
const FloatMatrix base,
const FloatArray strainVector,
bool  transpose = false 
)
static

Transforms 3d strain vector into another coordinate system.

Parameters
answerTransformed strain vector
baseTransformation matrix. There are on each column stored unit vectors of coordinate system (so called base vectors) to which we do transformation. These vectors must be expressed in the same coordinate system as source strainVector.
strainVector3d strain.
transposeDetermines if we transpose matrix before transforming.

Definition at line 1985 of file structuralmaterial.C.

References oofem::FloatArray::beProductOf(), and giveStrainVectorTranformationMtrx().

Referenced by give3dMaterialStiffnessMatrix(), oofem::MITC4Shell::giveCharacteristicTensor(), giveIPValue(), oofem::MITC4Shell::giveMidplaneIPValue(), and oofem::CompoDamageMat::giveRealStressVector().

void oofem::StructuralMaterial::transformStressVectorTo ( FloatArray answer,
const FloatMatrix base,
const FloatArray stressVector,
bool  transpose = false 
)
static

Transforms 3d stress vector into another coordinate system.

Parameters
answerTransformed stress vector.
baseTransformation matrix. There are on each column stored unit vectors of coordinate system (so called base vectors) to which we do transformation. These vectors must be expressed in the same coordinate system as source stressVector.
stressVectorTransformed 3d strain.
transposeDetermines if we transpose matrix before transforming.

Definition at line 2004 of file structuralmaterial.C.

References oofem::FloatArray::beProductOf(), and giveStressVectorTranformationMtrx().

Referenced by oofem::ConcreteDPM2::computeAlpha(), give3dMaterialStiffnessMatrix(), oofem::MITC4Shell::giveCharacteristicTensor(), giveIPValue(), oofem::MITC4Shell::giveMidplaneIPValue(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), oofem::CompoDamageMat::giveRealStressVector(), oofem::RCM2Material::giveRealStressVector(), oofem::ConcreteDPM::performRegularReturn(), and oofem::ConcreteDPM2::performRegularReturn().

Friends And Related Function Documentation

friend class CrossSection
friend

Definition at line 700 of file structuralmaterial.h.

friend class LayeredCrossSection
friend

Definition at line 703 of file structuralmaterial.h.

friend class SimpleCrossSection
friend

Definition at line 702 of file structuralmaterial.h.

Definition at line 701 of file structuralmaterial.h.

Member Data Documentation

double oofem::StructuralMaterial::referenceTemperature
protected
std::vector< std::vector< int > > oofem::StructuralMaterial::svIndex
static
Initial value:
= {
{ 1, 6, 5 },
{ 6, 2, 4 },
{ 5, 4, 3 }
}

Symmetric Voigt index map.

Definition at line 124 of file structuralmaterial.h.

std::vector< std::vector< int > > oofem::StructuralMaterial::vIindex
static
Initial value:
= {
{ 1, 6, 5 },
{ 9, 2, 4 },
{ 8, 7, 3 }
}

Voigt index map.

Definition at line 121 of file structuralmaterial.h.


The documentation for this class was generated from the following files:

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:42 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011