OOFEM 3.0
Loading...
Searching...
No Matches
frcfcm.h
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35#ifndef frcfcm_h
36#define frcfcm_h
37
38
39#include "concretefcm.h"
40
42
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"
67
69
70namespace oofem {
75{
76protected:
78 double damage = 0.;
80 double tempDamage = 0.;
81
82public:
84
86 double giveDamage() const { return damage; }
88 double giveTempDamage() const { return tempDamage; }
90 void setTempDamage(double newDamage) { tempDamage = newDamage; }
91
92 void printOutputAt(FILE *file, TimeStep *tStep) const override;
93
94 const char *giveClassName() const override { return "FRCFCMStatus"; }
95
96 void initTempStatus() override;
97 void updateYourself(TimeStep *tStep) override;
98
99 void saveContext(DataStream &stream, ContextMode mode) override;
100 void restoreContext(DataStream &stream, ContextMode mode) override;
101};
102
103
110class FRCFCM : public ConcreteFCM
111{
112public:
113 FRCFCM(int n, Domain *d);
114
115 void initializeFrom(InputRecord &ir) override;
116
117 std::unique_ptr<MaterialStatus> CreateStatus(GaussPoint *gp) const override { return std::make_unique<FRCFCMStatus>(gp); }
118
119 int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override;
120
121protected:
123 double tau_0 = 0.;
124
126 double b0 = 0.;
128 double b1 = 0., b2 = 0., b3 = 0.;
129
131 double f = 0.;
132
134 double g = 0.;
135
137 double Vf = 0.;
138
140 double Lf = 0.;
141
143 double Df = 0.;
144
146 double Ef = 0.;
147
149 double Gfib = 0.;
150
152 double kfib = 0.;
153
155 double w_star = 0.;
156
158 double eta = 0.;
159
161 double gammaCrackFail = 0.;
162
164 double minDamageOpening = 0.;
165
172 int M = 0;
173
176
179
186 double dw0 = 0., dw1 = 0.;
187 bool smoothen = false;
188
195
199
209
210 double giveCrackingModulus(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep, int i) const override;
212 virtual double computeFiberBond(double w) const;
213 double giveNormalCrackingStress(GaussPoint *gp, TimeStep *tStep, double eps_cr, int i) const override;
215 virtual double computeStressInFibersInCracked(GaussPoint *gp, TimeStep *tStep, double eps_cr, int i) const;
216 double computeEffectiveShearModulus(GaussPoint *gp, TimeStep *tStep, int i) const override;
217 double computeD2ModulusForCrack(GaussPoint *gp, TimeStep *tStep, int icrack) const override;
219 virtual double estimateD2ModulusForCrack(GaussPoint *gp, TimeStep *tStep, int icrack) const;
220 double maxShearStress(GaussPoint *gp, TimeStep *tStep, int i) const override;
222 virtual double computeTempDamage(GaussPoint *gp, TimeStep *tStep) const;
224 virtual double computeCrackSpacing();
226 virtual double computeCrackFibreAngle(GaussPoint *gp, int i) const;
227 void checkSnapBack(GaussPoint *gp, TimeStep *tStep, int crack) const override { }
228 bool isStrengthExceeded(const FloatMatrix &base, GaussPoint *gp, TimeStep *tStep, int iCrack, double trialStress) const override;
229 double computeShearStiffnessRedistributionFactor(GaussPoint *gp, TimeStep *tStep, int ithCrackPlane, int jthCrackDirection) const override;
230 double computeOverallElasticStiffness(GaussPoint *gp, TimeStep *tStep) const override;
231 double computeOverallElasticShearModulus(GaussPoint *gp, TimeStep *tStep) const override { return this->computeOverallElasticStiffness(gp, tStep) / ( 2. * ( 1. + linearElasticMaterial.givePoissonsRatio() ) ); }
232};
233} // end namespace oofem
234#endif // frcfcm_h
ConcreteFCMStatus(GaussPoint *g)
ConcreteFCM(int n, Domain *d)
Definition concretefcm.C:44
IsotropicLinearElasticMaterial linearElasticMaterial
Definition fcm.h:199
void printOutputAt(FILE *file, TimeStep *tStep) const override
Print receiver's output to given stream.
Definition frcfcm.C:1215
void initTempStatus() override
Definition frcfcm.C:1232
void saveContext(DataStream &stream, ContextMode mode) override
Definition frcfcm.C:1260
void setTempDamage(double newDamage)
Sets the temp damage level to given value.
Definition frcfcm.h:90
double giveTempDamage() const
Returns the temporary damage level.
Definition frcfcm.h:88
FRCFCMStatus(GaussPoint *g)
Definition frcfcm.C:1209
void updateYourself(TimeStep *tStep) override
Definition frcfcm.C:1245
double damage
Damage level of material.
Definition frcfcm.h:78
const char * giveClassName() const override
Definition frcfcm.h:94
double giveDamage() const
Returns the last equilibrated damage level.
Definition frcfcm.h:86
void restoreContext(DataStream &stream, ContextMode mode) override
Definition frcfcm.C:1270
double tempDamage
Non-equilibrated damage level of material.
Definition frcfcm.h:80
FiberDamageType
Type of fibre damage which is triggered by crack shearing strain = w / u.
Definition frcfcm.h:197
@ FDAM_GammaCrackExp
Definition frcfcm.h:197
@ FDAM_GammaCrackLin
Definition frcfcm.h:197
FiberDamageType fiberDamageType
Definition frcfcm.h:198
virtual double computeFiberBond(double w) const
evaluates the fiber bond if w > w*
Definition frcfcm.C:563
double computeOverallElasticStiffness(GaussPoint *gp, TimeStep *tStep) const override
returns overall Young's modulus
Definition frcfcm.C:1183
double f
snubbing factor "f"
Definition frcfcm.h:131
FRCFCM(int n, Domain *d)
Definition frcfcm.C:45
virtual double computeStressInFibersInCracked(GaussPoint *gp, TimeStep *tStep, double eps_cr, int i) const
compute the nominal stress in fibers in the i-th crack
Definition frcfcm.C:620
virtual double estimateD2ModulusForCrack(GaussPoint *gp, TimeStep *tStep, int icrack) const
estimate shear modulus for a given crack plane (1, 2, 3). Uses equilibrated value of damage.
Definition frcfcm.C:915
double b2
Definition frcfcm.h:128
double kfib
fiber cross-sectional shear factor
Definition frcfcm.h:152
double tau_0
fiber shear strength at zero slip
Definition frcfcm.h:123
FiberShearStrengthType
Definition frcfcm.h:193
bool smoothen
Definition frcfcm.h:187
void initializeFrom(InputRecord &ir) override
Definition frcfcm.C:49
double gammaCrackFail
shear strain at full fibers rupture
Definition frcfcm.h:161
double computeOverallElasticShearModulus(GaussPoint *gp, TimeStep *tStep) const override
returns overall shear modulus
Definition frcfcm.h:231
std::unique_ptr< MaterialStatus > CreateStatus(GaussPoint *gp) const override
Definition frcfcm.h:117
double b0
micromechanical parameter for fiber shear according to Sajdlova
Definition frcfcm.h:126
double b1
micromechanical parameter for fiber shear according to Kabele
Definition frcfcm.h:128
double dw1
Definition frcfcm.h:186
virtual double computeTempDamage(GaussPoint *gp, TimeStep *tStep) const
evaluates temporary value of damage caused by fibre shearing
Definition frcfcm.C:955
double g
auxiliary parameter computed from snubbing factor "f"
Definition frcfcm.h:134
void checkSnapBack(GaussPoint *gp, TimeStep *tStep, int crack) const override
checks possible snap-back
Definition frcfcm.h:227
double Ef
fiber Young's modulus
Definition frcfcm.h:146
double giveNormalCrackingStress(GaussPoint *gp, TimeStep *tStep, double eps_cr, int i) const override
computes normal stress associated with i-th crack direction
Definition frcfcm.C:602
double Df
fiber diameter
Definition frcfcm.h:143
double Lf
fiber length
Definition frcfcm.h:140
double w_star
transitional opening
Definition frcfcm.h:155
FiberShearStrengthType fiberShearStrengthType
Definition frcfcm.h:194
double computeD2ModulusForCrack(GaussPoint *gp, TimeStep *tStep, int icrack) const override
shear modulus for a given crack plane (1, 2, 3)
Definition frcfcm.C:874
double computeEffectiveShearModulus(GaussPoint *gp, TimeStep *tStep, int i) const override
returns Geff which is necessary in the global stiffness matrix
Definition frcfcm.C:817
double maxShearStress(GaussPoint *gp, TimeStep *tStep, int i) const override
computes the maximum value of the shear stress; if the shear stress exceeds this value,...
Definition frcfcm.C:1004
double computeShearStiffnessRedistributionFactor(GaussPoint *gp, TimeStep *tStep, int ithCrackPlane, int jthCrackDirection) const override
function calculating ratio used to split shear slips on two crack planes
Definition frcfcm.C:1144
virtual double computeCrackSpacing()
computes crack spacing based on composition of the fibre composite
Definition frcfcm.C:1093
bool isStrengthExceeded(const FloatMatrix &base, GaussPoint *gp, TimeStep *tStep, int iCrack, double trialStress) const override
compares trial stress with strength. Returns true if the strength is exceeded. Function oveloaded in ...
Definition frcfcm.C:1117
virtual double computeCrackFibreAngle(GaussPoint *gp, int i) const
compute the angle between the fibre and i-th crack normal
Definition frcfcm.C:258
double minDamageOpening
minimum opening at which damage can start
Definition frcfcm.h:164
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
Definition frcfcm.C:1159
double giveCrackingModulus(MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep, int i) const override
returns stiffness in the normal direction of the i-th crack
Definition frcfcm.C:290
double dw0
Definition frcfcm.h:186
double fibreActivationOpening
crack opening at which the crossing fibers begin to be activated
Definition frcfcm.h:178
double b3
Definition frcfcm.h:128
double eta
aux. factor
Definition frcfcm.h:158
FiberType fiberType
Definition frcfcm.h:208
FloatArray orientationVector
orientation of fibres
Definition frcfcm.h:175
double Gfib
fiber shear modulus
Definition frcfcm.h:149
double Vf
volume fraction of fibers
Definition frcfcm.h:137
long ContextMode
Definition contextmode.h:43

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011