43 #define _IFT_FRCFCM_Name "frcfcm" 44 #define _IFT_FRCFCM_Vf "vf" 45 #define _IFT_FRCFCM_Lf "lf" 46 #define _IFT_FRCFCM_Df "df" 47 #define _IFT_FRCFCM_Ef "ef" 48 #define _IFT_FRCFCM_nuf "nuf" 49 #define _IFT_FRCFCM_Gfib "gfib" 50 #define _IFT_FRCFCM_kfib "kfib" 51 #define _IFT_FRCFCM_tau_0 "tau_0" 52 #define _IFT_FRCFCM_b0 "b0" 53 #define _IFT_FRCFCM_b1 "b1" 54 #define _IFT_FRCFCM_b2 "b2" 55 #define _IFT_FRCFCM_b3 "b3" 56 #define _IFT_FRCFCM_f "f" 57 #define _IFT_FRCFCM_M "m" 58 #define _IFT_FRCFCM_orientationVector "fibreorientationvector" 59 #define _IFT_FRCFCM_fssType "fsstype" 60 #define _IFT_FRCFCM_fDamType "fdamtype" 61 #define _IFT_FRCFCM_fiberType "fibertype" 62 #define _IFT_FRCFCM_gammaCrack "gammacrack" 63 #define _IFT_FRCFCM_computeCrackSpacing "computecrackspacing" 64 #define _IFT_FRCFCM_fibreActivationOpening "fibreactivationopening" 65 #define _IFT_FRCFCM_dw0 "dw0" 66 #define _IFT_FRCFCM_dw1 "dw1" 203 enum FiberDamageType { FDAM_NONE, FDAM_GammaCrackLin, FDAM_GammaCrackExp, FDAM_Unknown };
218 virtual double computeFiberBond(
double w);
220 virtual double giveNormalCrackingStress(
GaussPoint *gp,
double eps_cr,
int i);
222 virtual double computeStressInFibersInCracked(
GaussPoint *gp,
double eps_cr,
int i);
224 virtual double computeEffectiveShearModulus(
GaussPoint *gp,
int i);
226 virtual double computeD2ModulusForCrack(
GaussPoint *gp,
int icrack);
228 virtual double estimateD2ModulusForCrack(
GaussPoint *gp,
int icrack);
230 virtual double maxShearStress(
GaussPoint *gp,
int i);
232 virtual double computeTempDamage(
GaussPoint *gp);
234 virtual double computeCrackSpacing(
void);
236 virtual double computeCrackFibreAngle(
GaussPoint *gp,
int i);
242 virtual double computeShearStiffnessRedistributionFactor(
GaussPoint *gp,
int ithCrackPlane,
int jthCrackDirection);
244 virtual double computeOverallElasticStiffness(
void);
246 virtual double computeOverallElasticShearModulus(
void) {
return this->computeOverallElasticStiffness() / ( 2. * ( 1. + this->giveLinearElasticMaterial()->givePoissonsRatio() ) ); }
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual const char * giveClassName() const
GaussPoint * gp
Associated integration point.
double gammaCrackFail
shear strain at full fibers rupture
double minDamageOpening
minimum opening at which damage can start
Domain * domain
Link to domain object, useful for communicating with other FEM components.
double kfib
fiber cross-sectional shear factor
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
FRCFCMStatus(int n, Domain *d, GaussPoint *g)
This class manages the status of ConcreteFCM.
void setTempDamage(double newDamage)
Sets the temp damage level to given value.
double f
snubbing factor "f"
FiberShearStrengthType fiberShearStrengthType
virtual void initTempStatus()
initializes temporary status
MatResponseMode
Describes the character of characteristic material matrix.
double fibreActivationOpening
crack opening at which the crossing fibers begin to be activated
FiberDamageType fiberDamageType
int M
Exponent in the unloading-reloading constitutive law.
FiberShearStrengthType
Type strength of the shear bond.
double g
auxiliary parameter computed from snubbing factor "f"
double b0
micromechanical parameter for fiber shear according to Sajdlova
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
saves current context(state) into stream
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
double w_star
transitional opening
This class implements a ConcreteFCM material in a finite element problem.
FiberType
Type fo fibers in the composite.
This class implements a FRCFCM material (Fiber Reinforced Concrete base on Fixed Crack Model) in a fi...
FiberDamageType
Type of fibre damage which is triggered by crack shearing strain = w / u.
double giveDamage()
Returns the last equilibrated damage level.
virtual void checkSnapBack(GaussPoint *gp, int crack)
overrides real checking from concretefcm.C, here we assume that the fibers provide sufficient strengt...
virtual void updateYourself(TimeStep *tStep)
replaces equilibrated values with temporary values
double giveTempDamage()
Returns the temporary damage level.
double tempDamage
Non-equilibrated damage level of material.
FloatArray orientationVector
orientation of fibres
double Ef
fiber Young's modulus
double tau_0
fiber shear strength at zero slip
Abstract base class representing a material status information.
Class representing vector of real numbers.
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
restores context(state) from stream
Implementation of matrix containing floating point numbers.
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Writes information into the output file.
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
double Gfib
fiber shear modulus
This class manages the status of FRCFCM.
the oofem namespace is to define a context or scope in which all oofem names are defined.
Class representing integration point in finite element program.
virtual double computeOverallElasticShearModulus(void)
from the Poisson's ratio of matrix and the overall stiffness estimates G
Class representing solution step.
double Vf
volume fraction of fibers
double damage
Damage level of material.