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

This class implements a ConcreteFCM material in a finite element problem. More...

#include <concretefcm.h>

+ Inheritance diagram for oofem::ConcreteFCM:
+ Collaboration diagram for oofem::ConcreteFCM:

Public Member Functions

 ConcreteFCM (int n, Domain *d)
 
virtual ~ConcreteFCM ()
 
virtual IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
virtual int hasNonLinearBehaviour ()
 Returns nonzero if receiver is non linear. More...
 
virtual const char * giveClassName () const
 
virtual const char * giveInputRecordName () const
 
virtual MaterialStatusCreateStatus (GaussPoint *gp) const
 Creates new copy of associated status and inserts it into given integration point. More...
 
virtual double give (int aProperty, GaussPoint *gp)
 Returns the value of material property 'aProperty'. More...
 
virtual int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in Reduced form. More...
 
virtual MaterialStatusgiveStatus (GaussPoint *gp) const
 Returns material status of receiver in given integration point. More...
 
- Public Member Functions inherited from oofem::FCMMaterial
 FCMMaterial (int n, Domain *d)
 
virtual ~FCMMaterial ()
 
virtual int hasMaterialModeCapability (MaterialMode mode)
 Tests if material supports material mode. More...
 
IsotropicLinearElasticMaterialgiveLinearElasticMaterial ()
 
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 givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane stress stiffness matrix of receiver. More...
 
virtual void givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane strain stiffness matrix of receiver. 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 initializeCrack (GaussPoint *gp, FloatMatrix &base, int nCrack)
 
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_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_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 double computeNormalCrackOpening (GaussPoint *gp, int i)
 uses temporary cracking strain and characteristic length to obtain the crack opening More...
 
virtual double computeMaxNormalCrackOpening (GaussPoint *gp, int i)
 uses maximum equilibrated cracking strain and characteristic length to obtain the maximum reached crack opening More...
 
virtual double computeShearSlipOnCrack (GaussPoint *gp, int i)
 computes total shear slip on a given crack plane (i = 1, 2, 3); the slip is computed from the temporary cracking strain More...
 
- Public Member Functions inherited from oofem::StructuralMaterial
 StructuralMaterial (int n, Domain *d)
 Constructor. More...
 
virtual ~StructuralMaterial ()
 Destructor. 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_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_Warping (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_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 void give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d beam layer stiffness matrix of receiver. More...
 
virtual void givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d plate layer stiffness matrix of receiver. More...
 
virtual void giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 1d fiber stiffness matrix of receiver. More...
 
virtual void give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d lattice stiffness matrix of receiver. More...
 
virtual void give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 3d lattice stiffness matrix of receiver. More...
 
virtual void give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of plate subsoil model. More...
 
virtual void give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of beam3d subsoil model. More...
 
virtual void giveFirstPKStressVector_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...
 
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_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
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 bool hasProperty (int aProperty, GaussPoint *gp)
 Returns true if 'aProperty' exists on material. More...
 
virtual void modifyProperty (int aProperty, double value, GaussPoint *gp)
 Modify 'aProperty', which already exists on material. More...
 
double giveCastingTime ()
 
virtual bool isActivated (TimeStep *tStep)
 
virtual int hasCastingTimeSupport ()
 Tests if material supports casting time. More...
 
virtual void printYourself ()
 Prints receiver state on stdout. Useful for debugging. More...
 
virtual contextIOResultType saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Stores integration point state to output stream. More...
 
virtual contextIOResultType restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Reads integration point state to output stream. More...
 
virtual int checkConsistency ()
 Allows programmer to test some internal data, before computation begins. More...
 
virtual int initMaterial (Element *element)
 Optional function to call specific procedures when initializing a material. More...
 
virtual int packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
 Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
 Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int estimatePackSize (DataStream &buff, GaussPoint *ip)
 Estimates the necessary pack size to hold all packed data of receiver. More...
 
virtual double predictRelativeComputationalCost (GaussPoint *gp)
 Returns the weight representing relative computational cost of receiver The reference material model is linear isotropic material - its weight is set to 1.0 The other material models should compare to this reference model. More...
 
virtual double predictRelativeRedistributionCost (GaussPoint *gp)
 Returns the relative redistribution cost of the receiver. More...
 
virtual void initTempStatus (GaussPoint *gp)
 Initializes temporary variables stored in integration point status at the beginning of new time step. More...
 
- Public Member Functions inherited from oofem::FEMComponent
 FEMComponent (int n, Domain *d)
 Regular constructor, creates component with given number and belonging to given domain. More...
 
virtual ~FEMComponent ()
 Virtual destructor. More...
 
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...
 
- Public Member Functions inherited from oofem::RandomMaterialExtensionInterface
 RandomMaterialExtensionInterface ()
 Constructor. More...
 
virtual ~RandomMaterialExtensionInterface ()
 Destructor. More...
 
IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
void giveInputRecord (DynamicInputRecord &ir)
 
bool give (int key, GaussPoint *gp, double &value)
 Returns the property in associated status of given integration point if defined. More...
 
- Public Member Functions inherited from oofem::Interface
 Interface ()
 Constructor. More...
 
virtual ~Interface ()
 

Protected Types

enum  SofteningType {
  ST_NONE, ST_Exponential, ST_Linear, ST_Hordijk,
  ST_UserDefinedCrack, ST_LinearHardeningStrain, ST_UserDefinedStrain, ST_Unknown
}
 type of post-peak behavior in the normal direction to the crack plane More...
 
enum  ShearRetentionType {
  SHR_NONE, SHR_Const_ShearRetFactor, SHR_Const_ShearFactorCoeff, SHR_UserDefined_ShearRetFactor,
  SHR_Unknown
}
 type of reduction of the shear stiffness caused by cracking More...
 
enum  ShearStrengthType { SHS_NONE, SHS_Const_Ft, SHS_Collins_Interlock, SHS_Unknown }
 defines the maximum value of shear stress More...
 

Protected Member Functions

virtual double giveTensileStrength (GaussPoint *gp)
 returns tensile strength (can be random) More...
 
virtual double giveFractureEnergy (GaussPoint *gp)
 returns fracture energy (can be random) More...
 
virtual double giveCrackingModulus (MatResponseMode rMode, GaussPoint *gp, int i)
 returns stiffness in the normal direction of the i-th crack More...
 
virtual double computeEffectiveShearModulus (GaussPoint *gp, int i)
 returns Geff which is necessary in the global stiffness matrix More...
 
virtual double computeD2ModulusForCrack (GaussPoint *gp, int icrack)
 shear modulus for a given crack plane (1, 2, 3) More...
 
virtual double giveNormalCrackingStress (GaussPoint *gp, double eps_cr, int i)
 computes normal stress associated with i-th crack direction More...
 
virtual double maxShearStress (GaussPoint *gp, int i)
 computes the maximum value of the shear stress; if the shear stress exceeds this value, it is cropped More...
 
virtual void checkSnapBack (GaussPoint *gp, int crack)
 checks possible snap-back More...
 
- Protected Member Functions inherited from oofem::FCMMaterial
virtual void updateCrackStatus (GaussPoint *gp)
 updates crack statuses More...
 
virtual double giveCharacteristicElementLength (GaussPoint *gp, const FloatArray &crackPlaneNormal)
 returns characteristic element length in given direction More...
 
virtual double computeTotalD2Modulus (GaussPoint *gp, int i)
 shear modulus for a given shear direction (4, 5, 6) More...
 
virtual bool isIntactForShear (GaussPoint *gp, int i)
 returns true for closed or no cracks associated to given shear direction (i = 4, 5, 6) More...
 
virtual bool isIntact (GaussPoint *gp, int icrack)
 returns true for closed or no crack (i = 1, 2, 3) More...
 
virtual bool checkStrengthCriterion (FloatMatrix &newBase, const FloatArray &globalStress, GaussPoint *gp, TimeStep *tStep, int nCrack)
 checks if the globalStress does not exceed strength in the direction of newBase for n-th crack More...
 
virtual bool isStrengthExceeded (const FloatMatrix &base, GaussPoint *gp, TimeStep *tStep, int iCrack, double trialStress)
 compares trial stress with strength. Returns true if the strength is exceeded. Function oveloaded in the nonlocal approach for the fiber reinforced composites More...
 
virtual double computeShearStiffnessRedistributionFactor (GaussPoint *gp, int ithCrackPlane, int jthCrackDirection)
 function calculating ratio used to split shear slips on two crack planes More...
 
virtual double giveCrackSpacing (void)
 returns either user-provided value of crack spacing or a value computed from composition More...
 
virtual double giveNumberOfCracksInDirection (GaussPoint *gp, int iCrack)
 returns number of fictiotious parallel cracks in the direction of i-th crack More...
 
virtual double giveNumberOfCracksForShearDirection (GaussPoint *gp, int i)
 returns number of cracks for given shear direction (i = 4, 5, 6) which is treated as the maximum of the two associated normal directions More...
 
virtual void giveMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep)
 
virtual void giveLocalCrackedStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
 returns local stiffness matrix of the crack More...
 
virtual double computeOverallElasticStiffness (void)
 returns overall Young's modulus More...
 
virtual double computeOverallElasticShearModulus (void)
 returns overall shear modulus More...
 
- Protected Member Functions inherited from oofem::RandomMaterialExtensionInterface
void _generateStatusVariables (GaussPoint *) const
 Sets up (generates) the variables identified in randVariables array using generators given in randomVariableGenerators and stores them in given status. More...
 

Protected Attributes

double Gf
 Fracture energy. More...
 
double Ft
 Tensile strenght. More...
 
double beta
 shear retention factor More...
 
double sf
 shear factor More...
 
double fc
 Collins' aggregate interlock: compressive strength in MPa. More...
 
double ag
 Collins' aggregate interlock: aggregate diameter in appropriate units (same as FE mesh) More...
 
double lengthScale
 Collins' aggregate interlock: 1 for meter, 1000 for analysis in mm. More...
 
FloatArray soft_w
 user-defined softening (traction-COD) More...
 
FloatArray soft_function_w
 
FloatArray soft_eps
 user-defined softening (traction-strain) More...
 
FloatArray soft_function_eps
 
FloatArray beta_w
 user-defined shear retention factor (with respect to crack opening) More...
 
FloatArray beta_function
 
double H
 hardening modulus More...
 
double eps_f
 strain at failure More...
 
SofteningType softType
 
ShearRetentionType shearType
 
ShearStrengthType shearStrengthType
 
- Protected Attributes inherited from oofem::FCMMaterial
IsotropicLinearElasticMateriallinearElasticMaterial
 
int nAllowedCracks
 allowed number of cracks (user-defined) More...
 
ElementCharSizeMethod ecsMethod
 Method used for evaluation of characteristic element size. More...
 
bool multipleCrackShear
 if true = takes shear compliance of all cracks, false = only dominant crack contribution, default value is false More...
 
double crackSpacing
 value of crack spacing (allows to "have" more parallel cracks in one direction if the element size exceeds user-defined or computed crack spacing). More...
 
- Protected Attributes inherited from oofem::StructuralMaterial
double referenceTemperature
 Reference temperature (temperature, when material has been built into structure). More...
 
- Protected Attributes inherited from oofem::Material
Dictionary propertyDictionary
 Property dictionary. More...
 
double castingTime
 Casting time. More...
 
- Protected Attributes inherited from oofem::FEMComponent
int number
 Component number. More...
 
Domaindomain
 Link to domain object, useful for communicating with other FEM components. More...
 
- Protected Attributes inherited from oofem::RandomMaterialExtensionInterface
IntArray randVariables
 Array of randomized variables (identified by a key). More...
 
IntArray randomVariableGenerators
 Array of generators id's for corresponding randomized variables. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from oofem::StructuralMaterial
static int giveSymVI (int ind1, int ind2)
 
static int giveVI (int ind1, int ind2)
 
static int giveVoigtVectorMask (IntArray &answer, MaterialMode mmode)
 Returns a mask of the vector indicies corresponding to components in a general (non-symmetric) second order tensor of some stress/strain/deformation measure that performes work. More...
 
static int giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode)
 The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor. More...
 
static void giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode)
 Gives the inverted version of giveVoigtVectorMask. More...
 
static int giveSizeOfVoigtVector (MaterialMode mmode)
 Returns the size of reduced stress/strain vector according to given mode. More...
 
static int giveSizeOfVoigtSymVector (MaterialMode mmode)
 Returns the size of symmetric part of a reduced stress/strain vector according to given mode. More...
 
static void giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
 Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. More...
 
static void giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
 Converts the reduced deformation gradient Voigt vector (2nd order tensor). More...
 
static void giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. More...
 
static void giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. More...
 
static void giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. More...
 
static void giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode)
 Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More...
 
static void giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
 Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. More...
 
static void giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
 Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More...
 
static void transformStrainVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &strainVector, bool transpose=false)
 Transforms 3d strain vector into another coordinate system. More...
 
static void transformStressVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &stressVector, bool transpose=false)
 Transforms 3d stress vector into another coordinate system. More...
 
static double computeVonMisesStress (const FloatArray *currentStress)
 Computes equivalent of von Mises stress. More...
 
static void giveStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 3d strain vector transformation matrix from standard vector transformation matrix. More...
 
static void give2DStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 2d strain vector transformation matrix from standard vector transformation matrix. More...
 
static void giveStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 3d stress vector transformation matrix from standard vector transformation matrix. More...
 
static void givePlaneStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 2d stress vector transformation matrix from standard vector transformation matrix. More...
 
static void sortPrincDirAndValCloseTo (FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir)
 Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDir) to be closed to some (often previous) principal directions (toPDir). More...
 
static void computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode)
 Common functions for convenience. More...
 
static void computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
 Computes principal values and directions of stress or strain vector. More...
 
static double computeDeviatoricVolumetricSplit (FloatArray &dev, const FloatArray &s)
 Computes split of receiver into deviatoric and volumetric part. More...
 
static void computeDeviatoricVolumetricSum (FloatArray &s, const FloatArray &dev, double mean)
 
static void applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu)
 
static void applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double GModulus)
 
static void applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu)
 
static void applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double GModulus)
 
static void applyElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu)
 
static void applyElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu)
 
static double computeStressNorm (const FloatArray &stress)
 
static double computeFirstInvariant (const FloatArray &s)
 
static double computeSecondStressInvariant (const FloatArray &s)
 
static double computeThirdStressInvariant (const FloatArray &s)
 
static double computeFirstCoordinate (const FloatArray &s)
 
static double computeSecondCoordinate (const FloatArray &s)
 
static double computeThirdCoordinate (const FloatArray &s)
 
- Static Public Attributes inherited from oofem::StructuralMaterial
static std::vector< std::vector< int > > vIindex
 Voigt index map. More...
 
static std::vector< std::vector< int > > svIndex
 Symmetric Voigt index map. More...
 

Detailed Description

This class implements a ConcreteFCM material in a finite element problem.

Cracking is described using Fixed Crack Model based on crack band approach. The model offers several predefined functions for softening and for reduction of shear stiffness caused by cracking.

Definition at line 98 of file concretefcm.h.

Member Enumeration Documentation

type of reduction of the shear stiffness caused by cracking

Enumerator
SHR_NONE 
SHR_Const_ShearRetFactor 
SHR_Const_ShearFactorCoeff 
SHR_UserDefined_ShearRetFactor 
SHR_Unknown 

Definition at line 181 of file concretefcm.h.

defines the maximum value of shear stress

Enumerator
SHS_NONE 
SHS_Const_Ft 
SHS_Collins_Interlock 
SHS_Unknown 

Definition at line 185 of file concretefcm.h.

type of post-peak behavior in the normal direction to the crack plane

Enumerator
ST_NONE 
ST_Exponential 
ST_Linear 
ST_Hordijk 
ST_UserDefinedCrack 
ST_LinearHardeningStrain 
ST_UserDefinedStrain 
ST_Unknown 

Definition at line 177 of file concretefcm.h.

Constructor & Destructor Documentation

oofem::ConcreteFCM::ConcreteFCM ( int  n,
Domain d 
)
virtual oofem::ConcreteFCM::~ConcreteFCM ( )
inlinevirtual

Definition at line 102 of file concretefcm.h.

References oofem::MaterialStatus::initializeFrom().

Member Function Documentation

virtual MaterialStatus* oofem::ConcreteFCM::CreateStatus ( GaussPoint gp) const
inlinevirtual

Creates new copy of associated status and inserts it into given integration point.

Parameters
gpIntegration point where newly created status will be stored.
Returns
Reference to new status.

Reimplemented from oofem::FCMMaterial.

Reimplemented in oofem::FRCFCMNL, and oofem::FRCFCM.

Definition at line 112 of file concretefcm.h.

References oofem::ConcreteFCMStatus::ConcreteFCMStatus(), oofem::FEMComponent::domain, and oofem::IntegrationPointStatus::gp.

Referenced by giveStatus().

double oofem::ConcreteFCM::give ( int  aProperty,
GaussPoint gp 
)
virtual

Returns the value of material property 'aProperty'.

Property must be identified by unique int id. Integration point also passed to allow for materials with spatially varying properties

Parameters
aPropertyID of property requested.
gpIntegration point,
Returns
Property value.

Reimplemented from oofem::FCMMaterial.

Definition at line 767 of file concretefcm.C.

References Ft, ft_strength, Gf, gf_ID, oofem::RandomMaterialExtensionInterface::give(), and oofem::FCMMaterial::give().

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

Reimplemented from oofem::FCMMaterial.

Reimplemented in oofem::FRCFCMNL.

Definition at line 109 of file concretefcm.h.

virtual double oofem::ConcreteFCM::giveFractureEnergy ( GaussPoint gp)
inlineprotectedvirtual

returns fracture energy (can be random)

Definition at line 156 of file concretefcm.h.

References gf_ID.

Referenced by checkSnapBack(), giveCrackingModulus(), and giveNormalCrackingStress().

virtual const char* oofem::ConcreteFCM::giveInputRecordName ( ) const
inlinevirtual
Returns
Input record name of the receiver.

Implements oofem::FEMComponent.

Reimplemented in oofem::FRCFCMNL.

Definition at line 110 of file concretefcm.h.

References _IFT_ConcreteFCM_Name.

int oofem::ConcreteFCM::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.

Reimplemented from oofem::FCMMaterial.

Reimplemented in oofem::FRCFCMNL, and oofem::FRCFCM.

Definition at line 782 of file concretefcm.C.

References oofem::FCMMaterial::giveIPValue().

Referenced by oofem::FRCFCM::giveIPValue().

virtual double oofem::ConcreteFCM::giveTensileStrength ( GaussPoint gp)
inlineprotectedvirtual

returns tensile strength (can be random)

Implements oofem::FCMMaterial.

Definition at line 153 of file concretefcm.h.

References ft_strength.

Referenced by checkSnapBack(), giveCrackingModulus(), giveNormalCrackingStress(), and oofem::FRCFCM::isStrengthExceeded().

virtual int oofem::ConcreteFCM::hasNonLinearBehaviour ( )
inlinevirtual

Returns nonzero if receiver is non linear.

Reimplemented from oofem::FCMMaterial.

Definition at line 108 of file concretefcm.h.

IRResultType oofem::ConcreteFCM::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::FCMMaterial.

Reimplemented in oofem::FRCFCM, and oofem::FRCFCMNL.

Definition at line 53 of file concretefcm.C.

References _IFT_ConcreteFCM_ag, _IFT_ConcreteFCM_beta, _IFT_ConcreteFCM_beta_function, _IFT_ConcreteFCM_beta_w, _IFT_ConcreteFCM_eps_f, _IFT_ConcreteFCM_fc, _IFT_ConcreteFCM_ft, _IFT_ConcreteFCM_gf, _IFT_ConcreteFCM_H, _IFT_ConcreteFCM_lengthScale, _IFT_ConcreteFCM_sf, _IFT_ConcreteFCM_shearStrengthType, _IFT_ConcreteFCM_shearType, _IFT_ConcreteFCM_soft_eps, _IFT_ConcreteFCM_soft_function_eps, _IFT_ConcreteFCM_soft_function_w, _IFT_ConcreteFCM_soft_w, _IFT_ConcreteFCM_softType, ag, oofem::FloatArray::at(), beta, beta_function, beta_w, eps_f, fc, Ft, Gf, oofem::FCMMaterial::giveCrackSpacing(), oofem::FloatArray::giveSize(), H, oofem::RandomMaterialExtensionInterface::initializeFrom(), oofem::IsotropicLinearElasticMaterial::initializeFrom(), oofem::FCMMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_BAD_FORMAT, oofem::IRRT_OK, lengthScale, oofem::FCMMaterial::linearElasticMaterial, OOFEM_ERROR, OOFEM_WARNING, sf, shearStrengthType, shearType, SHR_Const_ShearFactorCoeff, SHR_Const_ShearRetFactor, SHR_NONE, SHR_UserDefined_ShearRetFactor, SHS_Collins_Interlock, SHS_Const_Ft, SHS_NONE, soft_eps, soft_function_eps, soft_function_w, soft_w, softType, ST_Exponential, ST_Hordijk, ST_Linear, ST_LinearHardeningStrain, ST_NONE, ST_UserDefinedCrack, and ST_UserDefinedStrain.

Referenced by oofem::FRCFCM::initializeFrom().

double oofem::ConcreteFCM::maxShearStress ( GaussPoint gp,
int  i 
)
protectedvirtual

computes the maximum value of the shear stress; if the shear stress exceeds this value, it is cropped

Implements oofem::FCMMaterial.

Reimplemented in oofem::FRCFCM.

Definition at line 720 of file concretefcm.C.

References ag, oofem::FCMMaterial::computeMaxNormalCrackOpening(), oofem::FCMMaterial::computeNormalCrackOpening(), fc, fcm_BIGNUMBER, Ft, oofem::FCMMaterial::isIntactForShear(), lengthScale, oofem::max(), OOFEM_ERROR, shearStrengthType, SHS_Collins_Interlock, SHS_Const_Ft, and SHS_NONE.

Referenced by oofem::FRCFCM::maxShearStress().

Member Data Documentation

double oofem::ConcreteFCM::ag
protected

Collins' aggregate interlock: aggregate diameter in appropriate units (same as FE mesh)

Definition at line 136 of file concretefcm.h.

Referenced by initializeFrom(), and maxShearStress().

double oofem::ConcreteFCM::beta
protected

shear retention factor

Definition at line 128 of file concretefcm.h.

Referenced by computeD2ModulusForCrack(), computeEffectiveShearModulus(), and initializeFrom().

FloatArray oofem::ConcreteFCM::beta_function
protected

Definition at line 145 of file concretefcm.h.

Referenced by computeD2ModulusForCrack(), and initializeFrom().

FloatArray oofem::ConcreteFCM::beta_w
protected

user-defined shear retention factor (with respect to crack opening)

Definition at line 145 of file concretefcm.h.

Referenced by computeD2ModulusForCrack(), and initializeFrom().

double oofem::ConcreteFCM::eps_f
protected

strain at failure

Definition at line 150 of file concretefcm.h.

Referenced by giveCrackingModulus(), giveNormalCrackingStress(), and initializeFrom().

double oofem::ConcreteFCM::fc
protected

Collins' aggregate interlock: compressive strength in MPa.

Definition at line 134 of file concretefcm.h.

Referenced by initializeFrom(), and maxShearStress().

double oofem::ConcreteFCM::Ft
protected
double oofem::ConcreteFCM::Gf
protected

Fracture energy.

Definition at line 124 of file concretefcm.h.

Referenced by checkSnapBack(), ConcreteFCM(), give(), giveCrackingModulus(), giveNormalCrackingStress(), and initializeFrom().

double oofem::ConcreteFCM::H
protected

hardening modulus

Definition at line 148 of file concretefcm.h.

Referenced by giveCrackingModulus(), giveNormalCrackingStress(), and initializeFrom().

double oofem::ConcreteFCM::lengthScale
protected

Collins' aggregate interlock: 1 for meter, 1000 for analysis in mm.

Definition at line 138 of file concretefcm.h.

Referenced by initializeFrom(), and maxShearStress().

double oofem::ConcreteFCM::sf
protected

shear factor

Definition at line 130 of file concretefcm.h.

Referenced by computeD2ModulusForCrack(), and initializeFrom().

ShearStrengthType oofem::ConcreteFCM::shearStrengthType
protected

Definition at line 186 of file concretefcm.h.

Referenced by initializeFrom(), and maxShearStress().

FloatArray oofem::ConcreteFCM::soft_eps
protected

user-defined softening (traction-strain)

Definition at line 143 of file concretefcm.h.

Referenced by checkSnapBack(), giveCrackingModulus(), giveNormalCrackingStress(), and initializeFrom().

FloatArray oofem::ConcreteFCM::soft_function_eps
protected
FloatArray oofem::ConcreteFCM::soft_function_w
protected
FloatArray oofem::ConcreteFCM::soft_w
protected

user-defined softening (traction-COD)

Definition at line 141 of file concretefcm.h.

Referenced by checkSnapBack(), giveCrackingModulus(), giveNormalCrackingStress(), and initializeFrom().

SofteningType oofem::ConcreteFCM::softType
protected

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:34 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011