OOFEM 3.0
Loading...
Searching...
No Matches
mplasticmaterial2.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 mplasticmaterial2_h
36#define mplasticmaterial2_h
37
40#include "dictionary.h"
41#include "intarray.h"
42#include "floatarray.h"
43#include "floatmatrix.h"
45
46#include <vector>
47#include <set>
48
49namespace oofem {
50class GaussPoint;
51
68{
69public:
71
72protected:
76
83
85 state_flag_values state_flag = MPlasticMaterial2Status :: PM_Elastic;
86 state_flag_values temp_state_flag = MPlasticMaterial2Status :: PM_Elastic;
87
89 double damage = 0., tempDamage = 0.;
90
95
96public:
97 MPlasticMaterial2Status(GaussPoint * g, int statusSize);
98
99 void printOutputAt(FILE *file, TimeStep *tStep) const override;
100
101 void initTempStatus() override;
102 void updateYourself(TimeStep *tStep) override;
103
104 void saveContext(DataStream &stream, ContextMode mode) override;
105 void restoreContext(DataStream &stream, ContextMode mode) override;
106
115
120
121 void letTempDamageBe(double v) { tempDamage = v; }
122 double giveDamage() { return damage; }
123 double giveTempDamage() { return tempDamage; }
124
128
131 const FloatArray &giveTempGamma() { return tempGamma; }
132 void setTempGamma(FloatArray v) { tempGamma = std :: move(v); }
133
134 const char *giveClassName() const override { return "MPlasticMaterial2Status"; }
135};
136
174{
175protected:
179 int nsurf;
188 mutable std :: set< long >populationSet; // FIXME: Race condition when threaded.
189
190public:
191 MPlasticMaterial2(int n, Domain * d);
192 virtual ~MPlasticMaterial2();
193
194 bool hasMaterialModeCapability(MaterialMode mode) const override;
195 const char *giveClassName() const override { return "MPlasticMaterial2"; }
196
199
200 bool isCharacteristicMtrxSymmetric(MatResponseMode rMode) const override { return true; }
201
203 {
204 double alpha = this->linearElasticMaterial->give(tAlpha, gp);
205 return {alpha, alpha, alpha, 0., 0., 0.};
206 }
207
208 FloatMatrixF<6,6> give3dMaterialStiffnessMatrix(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override;
209
211 const FloatArray &strain, TimeStep *tStep) const override;
212
214 {
215 FloatArray answer;
216 const_cast<MPlasticMaterial2*>(this)->giveRealStressVector(answer, gp, strain, tStep);
217 return answer;
218 }
220 {
221 FloatArray answer;
222 const_cast<MPlasticMaterial2*>(this)->giveRealStressVector(answer, gp, strain, tStep);
223 return answer;
224 }
226 {
227 FloatArray answer;
228 const_cast<MPlasticMaterial2*>(this)->giveRealStressVector(answer, gp, strain, tStep);
229 return answer;
230 }
232 {
233 FloatArray answer;
234 const_cast<MPlasticMaterial2*>(this)->giveRealStressVector(answer, gp, strain, tStep);
235 return answer;
236 }
237
238 virtual double computeDamage(GaussPoint *gp, const FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const;
239
240 int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override;
241
242 // auxiliary functions
243 virtual int giveSizeOfFullHardeningVarsVector() const { return 0; }
244 virtual int giveSizeOfReducedHardeningVarsVector(GaussPoint *) const { return 0; }
245
246 std::unique_ptr<MaterialStatus> CreateStatus(GaussPoint *gp) const override;
247
248protected:
249
250 virtual int giveMaxNumberOfActiveYieldConds(GaussPoint *gp) const = 0;
251 void closestPointReturn(FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma,
252 GaussPoint *gp,
253 const FloatArray &totalStrain, FloatArray &plasticStrainR,
254 FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const;
255
256 void cuttingPlaneReturn(FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma,
257 GaussPoint *gp,
258 const FloatArray &totalStrain, FloatArray &plasticStrainR,
259 FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const;
260
261 // add here some auxiliary functions if needed
262 /* void computeGradientVector (FloatArray& answer, functType ftype, int isurf, GaussPoint* gp, const FloatArray& fullStressVector,
263 * const FloatArray& fullStressSpaceHardeningVars);*/
264 void computeResidualVector(FloatArray &answer, GaussPoint *gp, const FloatArray &gamma,
265 const IntArray &activeConditionMap, const FloatArray &plasticStrainVectorR,
266 const FloatArray &strainSpaceHardeningVariables, std :: vector< FloatArray > &gradVec) const;
267 virtual FloatMatrix giveConsistentStiffnessMatrix(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const;
268
269 virtual FloatMatrix giveElastoPlasticStiffnessMatrix(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const;
270
272 GaussPoint *gp, const FloatMatrix &elasticModuliInverse,
273 const FloatArray &gamma, const IntArray &activeConditionMap,
274 const FloatArray &fullStressVector,
275 const FloatArray &strainSpaceHardeningVariables) const;
276
277 /*void computeDiagModuli(FloatMatrix& answer,
278 * GaussPoint *gp, FloatMatrix &elasticModuliInverse,
279 * FloatMatrix &hardeningModuliInverse);*/
280
282 virtual double computeYieldValueAt(GaussPoint *gp, int isurf, const FloatArray &stressVector,
283 const FloatArray &strainSpaceHardeningVariables) const = 0;
284
286 virtual void computeStressGradientVector(FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector,
287 const FloatArray &strainSpaceHardeningVariables) const = 0;
288 void computeReducedStressGradientVector(FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector,
289 const FloatArray &strainSpaceHardeningVariables) const;
290
301 const FloatArray &stress, const FloatArray &dlambda,
302 const FloatArray &dplasticStrain, const IntArray &activeConditionMap) const = 0;
306 virtual void computeKGradientVector(FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, FloatArray &fullStressVector,
307 const FloatArray &strainSpaceHardeningVariables) const = 0;
308
312 virtual void computeReducedHardeningVarsSigmaGradient(FloatMatrix &answer, GaussPoint *gp, const IntArray &activeConditionMap,
313 const FloatArray &fullStressVector,
314 const FloatArray &strainSpaceHardeningVars,
315 const FloatArray &gamma) const = 0;
317 virtual void computeReducedHardeningVarsLamGradient(FloatMatrix &answer, GaussPoint *gp, int actSurf,
318 const IntArray &activeConditionMap,
319 const FloatArray &fullStressVector,
320 const FloatArray &strainSpaceHardeningVars,
321 const FloatArray &gamma) const = 0;
325 virtual int hasHardening() const = 0;
326 /* virtual void computeReducedGradientMatrix (FloatMatrix& answer, int isurf,
327 * GaussPoint *gp,
328 * const FloatArray& stressVector,
329 * const FloatArray& stressSpaceHardeningVars) = 0;*/
331 virtual void computeReducedSSGradientMatrix(FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector,
332 const FloatArray &strainSpaceHardeningVariables) const = 0;
334 virtual void computeReducedSKGradientMatrix(FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector,
335 const FloatArray &strainSpaceHardeningVariables) const = 0;
336
344 virtual void computeTrialStressIncrement(FloatArray &answer, GaussPoint *gp,
345 const FloatArray &strainIncrement, TimeStep *tStep) const;
346 virtual void computeReducedElasticModuli(FloatMatrix &answer, GaussPoint *gp,
347 TimeStep *tStep) const;
348 //virtual void compute3dElasticModuli(FloatMatrix& answer, GaussPoint *gp,
349 // TimeStep *tStep) = 0;
350
351 // next functions overloaded rom structural material level
352 FloatMatrixF<3,3> givePlaneStressStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override;
353 FloatMatrixF<4,4> givePlaneStrainStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override;
354 FloatMatrixF<1,1> give1dStressStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override;
355 FloatMatrixF<2,2> give2dBeamLayerStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override;
356 FloatMatrixF<5,5> givePlateLayerStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override;
357 FloatMatrixF<3,3> giveFiberStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override;
358
359protected:
360 long getPopulationSignature(IntArray &mask) const;
361 int testPopulation(long pop) const;
362 void clearPopulationSet() const;
363 void addNewPopulation(IntArray &mask) const;
364 int getNewPopulation(IntArray &result, IntArray &candidateMask, int degree, int size) const;
365};
366} // end namespace oofem
367#endif // mplasticmaterial2_h
void letTempPlasticStrainVectorBe(FloatArray v)
void letTempStateFlagBe(state_flag_values v)
double damage
Isotropic damage variables.
const FloatArray & givePlasticStrainVector() const
Returns the equilibrated strain vector.
const FloatArray & giveTempPlasticStrainVector() const
Returns the actual (temp) strain vector.
const FloatArray & giveTempStrainSpaceHardeningVarsVector() const
Returns the actual (temp) hardening variable vector.
MPlasticMaterial2Status(GaussPoint *g, int statusSize)
void letStrainSpaceHardeningVarsVectorBe(FloatArray v)
void letTempStrainSpaceHardeningVarsVectorBe(FloatArray v)
state_flag_values giveTempStateFlag()
void letPlasticStrainVectorBe(FloatArray v)
IntArray activeConditionMap
Active set of yield functions (needed for algorithmic stiffness).
const char * giveClassName() const override
state_flag_values state_flag
Yield function status indicator.
void printOutputAt(FILE *file, TimeStep *tStep) const override
Print receiver's output to given stream.
const FloatArray & giveStrainSpaceHardeningVars() const
Returns the equilibrated hardening variable vector.
void saveContext(DataStream &stream, ContextMode mode) override
void updateYourself(TimeStep *tStep) override
void restoreContext(DataStream &stream, ContextMode mode) override
FloatArray plasticStrainVector
Plastic strain vector.
FloatArray gamma
Consistency parameter values (needed for algorithmic stiffness).
const IntArray & giveTempActiveConditionMap()
virtual int hasHardening() const =0
virtual void computeTrialStressIncrement(FloatArray &answer, GaussPoint *gp, const FloatArray &strainIncrement, TimeStep *tStep) const
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
bool isCharacteristicMtrxSymmetric(MatResponseMode rMode) const override
void computeReducedStressGradientVector(FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const
std::unique_ptr< MaterialStatus > CreateStatus(GaussPoint *gp) const override
enum oofem::MPlasticMaterial2::ReturnMappingAlgoType rmType
LinearElasticMaterial * linearElasticMaterial
Reference to bulk (undamaged) material.
int getNewPopulation(IntArray &result, IntArray &candidateMask, int degree, int size) const
FloatMatrixF< 3, 3 > givePlaneStressStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
const char * giveClassName() const override
void addNewPopulation(IntArray &mask) const
virtual double computeDamage(GaussPoint *gp, const FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const
FloatArrayF< 1 > giveRealStressVector_1d(const FloatArrayF< 1 > &strain, GaussPoint *gp, TimeStep *tStep) const override
Default implementation relies on giveRealStressVector_StressControl.
functType
Type that allows to distinguish between yield function and loading function.
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 computeResidualVector(FloatArray &answer, GaussPoint *gp, const FloatArray &gamma, const IntArray &activeConditionMap, const FloatArray &plasticStrainVectorR, const FloatArray &strainSpaceHardeningVariables, std ::vector< FloatArray > &gradVec) const
MPlasticMaterial2(int n, Domain *d)
FloatMatrixF< 3, 3 > giveFiberStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
enum oofem::MPlasticMaterial2::plastType plType
virtual void computeReducedSSGradientMatrix(FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const =0
Computes second derivative of loading function with respect to stress.
virtual double computeYieldValueAt(GaussPoint *gp, int isurf, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const =0
Computes the value of yield function.
LinearElasticMaterial * giveLinearElasticMaterial()
Returns reference to undamaged (bulk) material.
virtual int giveMaxNumberOfActiveYieldConds(GaussPoint *gp) const =0
void cuttingPlaneReturn(FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma, GaussPoint *gp, const FloatArray &totalStrain, FloatArray &plasticStrainR, FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const
ReturnMappingAlgoType
Protected type to determine the return mapping algorithm.
FloatMatrixF< 6, 6 > give3dMaterialStiffnessMatrix(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 4 > giveRealStressVector_PlaneStrain(const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const override
Default implementation relies on giveRealStressVector_3d.
bool iterativeUpdateOfActiveConds
Flag indicating whether iterative update of a set of active yield conditions takes place.
FloatArrayF< 3 > giveRealStressVector_PlaneStress(const FloatArrayF< 3 > &strain, GaussPoint *gp, TimeStep *tStep) const override
Default implementation relies on giveRealStressVector_StressControl.
bool hasMaterialModeCapability(MaterialMode mode) const override
int nsurf
Number of yield surfaces.
FloatMatrixF< 1, 1 > give1dStressStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
void giveRealStressVector(FloatArray &answer, GaussPoint *gp, const FloatArray &strain, TimeStep *tStep) const override
virtual void computeReducedHardeningVarsSigmaGradient(FloatMatrix &answer, GaussPoint *gp, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &gamma) const =0
FloatMatrixF< 2, 2 > give2dBeamLayerStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
virtual void computeReducedHardeningVarsLamGradient(FloatMatrix &answer, GaussPoint *gp, int actSurf, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &gamma) const =0
computes derivative of vector with respect to lambda vector
virtual FloatMatrix giveElastoPlasticStiffnessMatrix(MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const
virtual void computeStressGradientVector(FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const =0
Computes the stress gradient of yield/loading function (df/d_sigma).
virtual int giveSizeOfFullHardeningVarsVector() const
virtual void computeReducedElasticModuli(FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep) const
virtual void computeKGradientVector(FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const =0
long getPopulationSignature(IntArray &mask) const
FloatArrayF< 6 > giveThermalDilatationVector(GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 5, 5 > givePlateLayerStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 4, 4 > givePlaneStrainStiffMtrx(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
virtual void computeReducedSKGradientMatrix(FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const =0
Computes second derivative of loading function with respect to stress and hardening vars.
void computeAlgorithmicModuli(FloatMatrix &answer, GaussPoint *gp, const FloatMatrix &elasticModuliInverse, const FloatArray &gamma, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const
virtual FloatMatrix giveConsistentStiffnessMatrix(MatResponseMode, GaussPoint *gp, TimeStep *tStep) const
int testPopulation(long pop) const
virtual int giveSizeOfReducedHardeningVarsVector(GaussPoint *) const
virtual void computeStrainHardeningVarsIncrement(FloatArray &answer, GaussPoint *gp, const FloatArray &stress, const FloatArray &dlambda, const FloatArray &dplasticStrain, const IntArray &activeConditionMap) const =0
std ::set< long > populationSet
Set for keeping record of generated populations of active yield conditions during return.
void closestPointReturn(FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma, GaussPoint *gp, const FloatArray &totalStrain, FloatArray &plasticStrainR, FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const
virtual double give(int aProperty, GaussPoint *gp) const
Definition material.C:51
StructuralMaterialStatus(GaussPoint *g)
Constructor. Creates new StructuralMaterialStatus with IntegrationPoint g.
StructuralMaterial(int n, Domain *d)
#define tAlpha
Definition matconst.h:67
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