OOFEM 3.0
Loading...
Searching...
No Matches
concretedpm.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 ConcreteDPM_h
36#define ConcreteDPM_h
37
39#include "floatarray.h"
40#include "floatmatrix.h"
41#include "cltypes.h"
44#include "gausspoint.h"
45#include "mathfem.h"
46
48
49#define _IFT_ConcreteDPM_Name "concretedpm"
50#define _IFT_ConcreteDPM_fc "fc"
51#define _IFT_ConcreteDPM_ft "ft"
52#define _IFT_ConcreteDPM_ecc "ecc"
53#define _IFT_ConcreteDPM_kinit "kinit"
54#define _IFT_ConcreteDPM_ahard "ahard"
55#define _IFT_ConcreteDPM_bhard "bhard"
56#define _IFT_ConcreteDPM_chard "chard"
57#define _IFT_ConcreteDPM_dhard "dhard"
58#define _IFT_ConcreteDPM_asoft "asoft"
59#define _IFT_ConcreteDPM_dilation "dilation"
60#define _IFT_ConcreteDPM_yieldtol "yieldtol"
61#define _IFT_ConcreteDPM_newtoniter "newtoniter"
62#define _IFT_ConcreteDPM_wf "wf"
63#define _IFT_ConcreteDPM_gf "gf"
64#define _IFT_ConcreteDPM_href "href"
65#define _IFT_ConcreteDPM_helem "helem"
67
68namespace oofem {
73// new approach to size-dependent adjustment of damage evolution
74// (can be deactivated by commenting out the following line)
75#define SOPHISTICATED_SIZEDEPENDENT_ADJUSTMENT
76
78{
79public:
92
98
99protected:
101
104
106
107 double dFDKappa = 0.;
108 double deltaLambda = 0.;
110
112
113 double kappaP = 0.;
114 double tempKappaP = 0.;
116
117 double le = 0.;
118
120
121 double equivStrain = 0.;
122 double tempEquivStrain = 0.;
123
124 double kappaD = 0.;
125 double tempKappaD = 0.;
126
127 double damage = 0.;
128 double tempDamage = 0.;
129
130 double deltaEquivStrain = 0.;
132
134
138
141
142
143#ifdef SOPHISTICATED_SIZEDEPENDENT_ADJUSTMENT
146 double epsloc = -1.;
147 double tempEpsloc = -1.;
148#endif
149
150public:
153
154 void initTempStatus() override;
155 void updateYourself(TimeStep *tStep) override;
156 void printOutputAt(FILE *file, TimeStep *tStep) const override;
157
158 void saveContext(DataStream &stream, ContextMode mode) override;
159 void restoreContext(DataStream &stream, ContextMode mode) override;
160
161 int setIPValue(const FloatArray &value, InternalStateType type);
162
163 const char *giveClassName() const override { return "ConcreteDPMStatus"; }
164
170
176 {
178 return sqrt(.5 * ( 2. * dev [ 0 ] * dev [ 0 ] + 2. * dev [ 1 ] * dev [ 1 ] + 2. * dev [ 2 ] * dev [ 2 ] +
179 dev [ 3 ] * dev [ 3 ] + dev [ 4 ] * dev [ 4 ] + dev [ 5 ] * dev [ 5 ] ) );
180 }
181
183 {
184 return 1. / 3. * ( plasticStrain [ 0 ] + plasticStrain [ 1 ] + plasticStrain [ 2 ] );
185 }
186
191 double giveKappaP() const { return kappaP; }
192
198 double giveKappaD() const { return kappaD; }
199
205 double giveEquivStrain() const { return equivStrain; }
206
207#ifdef SOPHISTICATED_SIZEDEPENDENT_ADJUSTMENT
213 double giveEpsLoc() const { return epsloc; }
214
220 void letTempEpslocBe(double v) { tempEpsloc = v; }
221
222#endif
223
229 double giveDamage() const { return damage; }
230
235 int giveStateFlag() const { return state_flag; }
236
242
248
249
255 double giveTempKappaP() const { return tempKappaP; }
256
262 double giveTempKappaD() const { return tempKappaD; }
263
269 double giveTempDamage() const { return tempDamage; }
270
276 double giveDeltaEquivStrain() const { return deltaEquivStrain; }
277
278
284 int giveTempStateFlag() const { return temp_state_flag; }
285
286 int giveTempVertexType() const { return tempVertexType; }
287
297
302 void letDeltaLambdaBe(double v) { deltaLambda = v; }
303
310
315 void letTempKappaPBe(double v)
316 { tempKappaP = v; }
317
318
323 void letTempKappaDBe(double v) { tempKappaD = v; }
327 void letKappaDBe(double v) { kappaD = v; }
332 void letTempEquivStrainBe(double v) { tempEquivStrain = v; }
336 void letEquivStrainBe(double v) { equivStrain = v; }
337
342 void letTempDamageBe(double v) { tempDamage = v; }
343
349
353 double giveLe() { return le; }
354
359 void setLe(double ls) { le = ls; }
360
367
368 void letTempVertexTypeBe(const int type) { tempVertexType = type; }
369};
370
371
383{
384protected:
385 // enum Concrete_VertexType { VT_Regular, VT_Tension, VT_Compression };
386 // mutable Concrete_VertexType vertexType; /// FIXME: This must be removed, material models must never have mutable state.
387
390
394 double fc = 0., ft = 0., ecc = 0.;
395
397 double AHard = 0.;
399 double BHard = 0.;
401 double CHard = 0.;
403 double DHard = 0.;
404
406 double ASoft = 0.;
407
409 double yieldHardInitial = 0.;
410
412 double dilationConst = 0.;
413
415 //double deltaLambda = 0.;
416
418 //double sig = 0.;
420 //double rho = 0.;
421
423 // double thetaTrial = 0.;
424
426 double m = 0.;
427
429 double mQ = 0.;
430
432 double helem = 0.;
433
435 double eM = 0.;
437 double gM = 0.;
439 double kM = 0.;
441 double nu = 0.;
442
444 //double kappaP = 0.;
445 //double tempKappaP = 0.;
446
448 //double kappaD = 0.;
449 //double tempKappaD = 0.;
450
452 //double damage = 0.;
453 //double tempDamage = 0.;
454
456 double ef = 0.;
457
459 double yieldTol = 0.;
460
462 int newtonIter = 0;
463
465 //FloatArray effectiveStress;
466
467#ifdef SOPHISTICATED_SIZEDEPENDENT_ADJUSTMENT
470 double href = 0.;
471#endif
472
473public:
475 ConcreteDPM(int n, Domain *d);
476
477 void initializeFrom(InputRecord &ir) override;
478
479 const char *giveClassName() const override { return "ConcreteDPM"; }
480 const char *giveInputRecordName() const override { return _IFT_ConcreteDPM_Name; }
481
483 { return static_cast< ConcreteDPMStatus * >( this->Material::giveStatus(gp) ); }
484
486 TimeStep *tStep) const override;
487
492 void performPlasticityReturn(GaussPoint *gp, FloatArrayF< 6 > &strain) const;
493
501
502 void checkForVertexCase(double &answer,
503 double sig,
504 double tempKappa,
505 GaussPoint *gp) const;
506
514 GaussPoint *gp,
515 double theta) const;
516
517
525 double apexStress,
526 GaussPoint *gp) const;
527
536 double computeYieldValue(double sig,
537 double rho,
538 double theta,
539 double tempKappa) const;
540
547 double computeHardeningOne(double tempKappa) const;
548
555 double computeHardeningOnePrime(double tempKappa) const;
556
557
567 double computeDFDKappa(double sig,
568 double rho,
569 double theta,
570 double tempKappa) const;
571
581 double computeDKappaDDeltaLambda(double sig,
582 double rho,
583 double theta,
584 double tempKappa) const;
585
593 virtual double computeDuctilityMeasure(double sig,
594 double rho,
595 double theta) const;
596
597
603 double rho,
604 double theta,
605 double tempKappa) const;
606
612 double rho,
613 double theta,
614 double tempKappa,
615 double deltaLambda) const;
616
622 double rho,
623 double tempKappa) const;
624
630 double computeRatioPotential(double sig,
631 double tempKappa) const;
632
638 double rho,
639 double tempKappa) const;
640
646 double rho,
647 double tempKappa) const;
648
655 double rho,
656 double theta,
657 double tempKappa) const;
658
663 double computeDDKappaDDeltaLambdaDKappa(double sig,
664 double rho,
665 double theta,
666 double tempKappa) const;
667
668
674 double rho,
675 double theta,
676 double tempKappa) const;
677
681 double computeTempKappa(double kappaInitial,
682 double sigTrial,
683 double rhoTrial,
684 double sig) const;
685
693 std::pair< double, double >computeDamage(const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const;
694
696 virtual double computeDamageParam(double kappa, GaussPoint *gp) const;
697
699 double computeInverseDamage(double dam, GaussPoint *gp) const;
700
702 virtual double computeEquivalentStrain(const FloatArrayF< 6 > &elasticStrain, GaussPoint *gp, TimeStep *tStep) const;
703
705 double computeDuctilityMeasureDamage(const FloatArrayF< 6 > &strain, GaussPoint *gp) const;
706
710 void initDamaged(double kappa,
711 const FloatArrayF< 6 > &elasticStrain,
712 GaussPoint *gp) const;
713
714
716 std::tuple< double, double, double >computeTrialCoordinates(const FloatArrayF< 6 > &stress, GaussPoint *gp) const;
717
719 void assignStateFlag(GaussPoint *gp) const;
720
723
725 void computeDSigDStress(FloatArrayF< 6 > &answer) const;
726
729
732
734 double computeDRDCosTheta(double theta, double ecc) const;
735
736 FloatMatrixF< 6, 6 >give3dMaterialStiffnessMatrix(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override;
737
738 bool isCharacteristicMtrxSymmetric(MatResponseMode rMode) const override { return false; }
739
740 int setIPValue(const FloatArray &value, GaussPoint *gp, InternalStateType type) override;
741
742 int giveIPValue(FloatArray &answer,
743 GaussPoint *gp,
745 TimeStep *tStep) override;
746
747 void restoreConsistency(GaussPoint *gp) override;
748 void saveContext(DataStream &stream, ContextMode mode) override;
749 void restoreContext(DataStream &stream, ContextMode mode) override;
750
751
752protected:
753 std::unique_ptr<MaterialStatus> CreateStatus(GaussPoint *gp) const override;
754};
755} // end namespace oofem
756#endif
void saveContext(DataStream &stream, ContextMode mode) override
double giveDeviatoricPlasticStrainNorm()
double giveTempKappaP() const
void letTempPlasticStrainBe(const FloatArrayF< 6 > &v)
void letTempDamageBe(double v)
const char * giveClassName() const override
ConcreteDPMStatus(GaussPoint *gp)
Constructor.
Definition concretedpm.C:54
void letDeltaEquivStrainBe(double v)
double giveDamage() const
state_flag_values
Values of history variable state_flag.
Definition concretedpm.h:81
FloatArrayF< 6 > tempPlasticStrain
void letTempKappaDBe(double v)
void letTempVertexTypeBe(const int type)
double giveTempVolumetricPlasticStrain() const
void printOutputAt(FILE *file, TimeStep *tStep) const override
Print receiver's output to given stream.
double giveKappaD() const
double giveVolumetricPlasticStrain() const
double giveDeltaEquivStrain() const
int setIPValue(const FloatArray &value, InternalStateType type)
double giveKappaP() const
void letEquivStrainBe(double v)
double giveTempKappaD() const
int giveTempVertexType() const
void initTempStatus() override
Definition concretedpm.C:65
void letTempEpslocBe(double v)
void letPlasticStrainBe(const FloatArrayF< 6 > &v)
const FloatArrayF< 6 > & giveTempPlasticStrain() const
void letTempStateFlagBe(int v)
void letTempKappaPBe(double v)
const FloatArrayF< 6 > & givePlasticStrain() const
void setLe(double ls)
void restoreContext(DataStream &stream, ContextMode mode) override
void letTempVolumetricPlasticStrainBe(double v)
double giveEpsLoc() const
double giveEquivStrain() const
void updateYourself(TimeStep *tStep) override
Definition concretedpm.C:80
FloatArrayF< 6 > plasticStrain
void letDeltaLambdaBe(double v)
int giveTempStateFlag() const
void letTempEquivStrainBe(double v)
void letKappaDBe(double v)
double giveTempDamage() const
double kM
Elastic bulk modulus.
void initDamaged(double kappa, const FloatArrayF< 6 > &elasticStrain, GaussPoint *gp) const
ConcreteDPMStatus * giveConcreteDPMStatus(GaussPoint *gp) const
IsotropicLinearElasticMaterial linearElasticMaterial
Linear elastic material.
double mQ
The dilation parameter of the plastic potential.
double ASoft
Parameter of the ductilityMeasure of the damage model.
double computeDDKappaDDeltaLambdaDKappa(double sig, double rho, double theta, double tempKappa) const
const char * giveClassName() const override
std::tuple< double, double, double > computeTrialCoordinates(const FloatArrayF< 6 > &stress, GaussPoint *gp) const
Compute the trial coordinates.
void restoreConsistency(GaussPoint *gp) override
double computeDFDKappa(double sig, double rho, double theta, double tempKappa) const
double yieldTol
Yield tolerance for the plasticity model.
FloatArrayF< 6 > computeDCosThetaDStress(const FloatArrayF< 6 > &stress) const
Computes the derivative of costheta with respect to the stress.
double helem
Element size (to be used in fracture energy approach (crack band).
double computeHardeningOne(double tempKappa) const
FloatArrayF< 6 > giveRealStressVector_3d(const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const override
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
void checkForVertexCase(double &answer, double sig, double tempKappa, GaussPoint *gp) const
FloatMatrixF< 6, 6 > computeDDRhoDDStress(const FloatArrayF< 6 > &stress) const
Computes the second derivative of rho with the respect to the stress.
std::unique_ptr< MaterialStatus > CreateStatus(GaussPoint *gp) const override
FloatMatrixF< 3, 3 > computeAMatrix(double sig, double rho, double theta, double tempKappa, double deltaLambda) const
double dilationConst
Control parameter for te volumetric plastic flow of the plastic potential.
FloatArrayF< 6 > computeDRhoDStress(const FloatArrayF< 6 > &stress) const
Computes the derivative of rho with respect to the stress.
void performVertexReturn(FloatArrayF< 6 > &stress, double apexStress, GaussPoint *gp) const
double gM
Elastic shear modulus.
double computeTempKappa(double kappaInitial, double sigTrial, double rhoTrial, double sig) const
FloatArrayF< 2 > computeDDKappaDDeltaLambdaDInv(double sig, double rho, double theta, double tempKappa) const
virtual double computeDuctilityMeasure(double sig, double rho, double theta) const
double computeYieldValue(double sig, double rho, double theta, double tempKappa) const
FloatMatrixF< 6, 6 > give3dMaterialStiffnessMatrix(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
double computeDKappaDDeltaLambda(double sig, double rho, double theta, double tempKappa) const
double computeInverseDamage(double dam, GaussPoint *gp) const
Compute the damage-driving variable from given damage.
void restoreContext(DataStream &stream, ContextMode mode) override
void computeDSigDStress(FloatArrayF< 6 > &answer) const
Computes the derivative of sig with respect to the stress.
ConcreteDPM(int n, Domain *d)
Constructor.
void initializeFrom(InputRecord &ir) override
double AHard
Parameter of the ductilityMeasure of the plasticity model.
double computeRatioPotential(double sig, double tempKappa) const
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
double nu
Elastic Poisson's ration.
double m
Plastic multiplier of the plasticity model.
double computeDuctilityMeasureDamage(const FloatArrayF< 6 > &strain, GaussPoint *gp) const
Compute the ductility measure for the damage model.
double ef
Hardening variable of plasticity model.
const char * giveInputRecordName() const override
FloatArrayF< 2 > computeDDGDInvDKappa(double sig, double rho, double tempKappa) const
std::pair< double, double > computeDamage(const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const
double CHard
Parameter of the ductilityMeasure of the plasticity model.
int newtonIter
Maximum number of iterations for stress return.
double computeDRDCosTheta(double theta, double ecc) const
Compute the derivative of R with respect to costheta.
FloatArrayF< 2 > computeDGDInv(double sig, double rho, double tempKappa) const
FloatArrayF< 2 > computeDDuctilityMeasureDInv(double sig, double rho, double theta, double tempKappa) const
FloatArrayF< 2 > computeDFDInv(double sig, double rho, double theta, double tempKappa) const
int setIPValue(const FloatArray &value, GaussPoint *gp, InternalStateType type) override
double href
Stress and its deviatoric part.
virtual double computeDamageParam(double kappa, GaussPoint *gp) const
Compute damage parameter.
void assignStateFlag(GaussPoint *gp) const
Assign state flag.
double computeHardeningOnePrime(double tempKappa) const
double DHard
Parameter of the ductilityMeasure of the plasticity model.
void performPlasticityReturn(GaussPoint *gp, FloatArrayF< 6 > &strain) const
virtual double computeEquivalentStrain(const FloatArrayF< 6 > &elasticStrain, GaussPoint *gp, TimeStep *tStep) const
Compute equivalent strain value.
double eM
Elastic Young's modulus.
double BHard
Parameter of the ductilityMeasure of the plasticity model.
double yieldHardInitial
Parameter of the hardening law of the plasticity model.
bool isCharacteristicMtrxSymmetric(MatResponseMode rMode) const override
void saveContext(DataStream &stream, ContextMode mode) override
FloatMatrixF< 2, 2 > computeDDGDDInv(double sig, double rho, double tempKappa) const
void performRegularReturn(FloatArrayF< 6 > &stress, GaussPoint *gp, double theta) const
GaussPoint * gp
Associated integration point.
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Definition material.C:206
StructuralMaterialStatus(GaussPoint *g)
Constructor. Creates new StructuralMaterialStatus with IntegrationPoint g.
StructuralMaterial(int n, Domain *d)
static FloatArrayF< 6 > computeDeviator(const FloatArrayF< 6 > &s)
#define _IFT_ConcreteDPM_Name
Definition concretedpm.h:49
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