Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
prandtl_reuss_impl.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2014-2015: G-CSC, Goethe University Frankfurt
3 * Author: Raphael Prohl
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifndef PRANDTL_REUSS_IMPL_H_
34#define PRANDTL_REUSS_IMPL_H_
35
36#include "prandtl_reuss.h"
37
38#define PROFILE_PRANDTL_REUSS
39#ifdef PROFILE_PRANDTL_REUSS
40 #define PRANDTL_REUSS_PROFILE_FUNC() PROFILE_FUNC_GROUP("Prandtl Reuss")
41 #define PRANDTL_REUSS_PROFILE_BEGIN(name) PROFILE_BEGIN_GROUP(name, "Prandtl Reuss")
42 #define PRANDTL_REUSS_PROFILE_END() PROFILE_END()
43#else
44 #define PRANDTL_REUSS_PROFILE_FUNC()
45 #define PRANDTL_REUSS_PROFILE_BEGIN(name)
46 #define PRANDTL_REUSS_PROFILE_END()
47#endif
48
49namespace ug{
50namespace SmallStrainMechanics{
51
52template <typename TDomain>
54 IMaterialLaw<TDomain>(),
55 m_MaxHardIter(100), m_HardAccuracy(0.0), m_tangentAccur(1e-08),
56 m_bHardModulus(false), m_bHardExp(false)
57{
58 // set default material constants
59 matConsts.mu = 0.0; matConsts.kappa = 0.0;
60
61 matConsts.K_0 = 0.0; m_hardening = 0;
62 matConsts.K_inf = 0.0; matConsts.Hard = 0.0; matConsts.omega = 0.0;
63
64 m_max_k_tan = 0.0;
65 m_min_k_tan = 100.0;
66 m_plasticIPs = 0;
67
68 std::stringstream ss;
69 ss << "Prandtl Reuss Plasticity: \n";
70 m_materialConfiguration = ss.str();
71
73}
74
75template <typename TDomain>
76void
79{
80 std::stringstream ss;
81
82 switch (hard){
83 // perfect hardening
84 case 0: m_hardening = 0;
85 ss << m_materialConfiguration << "perfect hardening \n";
86 m_materialConfiguration = ss.str();
87 break;
88
89 // linear hardening
90 case 1: m_hardening = 1;
91 ss << m_materialConfiguration << "linear hardening \n";
92 m_materialConfiguration = ss.str();
93 break;
94
95 // exponential hardening
96 case 2: m_hardening = 2;
97 m_MaxHardIter = 100;
98 m_HardAccuracy = 1e-10;
99 ss << m_materialConfiguration << "exponential hardening \n"
100 << " max. hardening iterations = " << m_MaxHardIter << "\n"
101 << " hardening accuracy = " << m_HardAccuracy << "\n";
102 m_materialConfiguration = ss.str();
103 break;
104
105 default: UG_THROW(hard << " is not a valid hardening behavior! "
106 "Choose 0 (perfect), 1 (linear) or 2 (exponential) ! \n");
107 }
108}
109
110template <typename TDomain>
111void
113init()
114{
115 // check, if all parameter for a hardening behavior are set
116 if (m_hardening == 1){ // linear hardening law
117 if (!m_bHardModulus)
118 UG_THROW("No hardening modulus set! This is necessary for a "
119 "linear hardening law in PrandtlReuss::init() \n");
120 }
121
122 if (m_hardening == 2){ // exponential hardening law
123 if (!m_bHardModulus)
124 UG_THROW("No hardening modulus set! This is necessary for an "
125 "exponential hardening law in PrandtlReuss::init() \n");
126 if (!m_bHardExp)
127 UG_THROW("No hardening exponent set! This is necessary for an "
128 "exponential hardening law in PrandtlReuss::init() \n");
129 }
130}
131
132template <typename TDomain>
133void
136 const MathMatrix<dim, dim>& strain_p_old_t, const number alpha)
137{
138 MathMatrix<dim, dim> strial, normal, strain, strain_p_new;
139 number gamma;
140
141 Flowrule(strain_p_new, strain, gamma, strial, normal, GradU, strain_p_old_t, alpha);
142 ConstLaw(stressTens, strain, strial, gamma, normal);
143}
144
145template <typename TDomain>
146void
148stressTensor(MathMatrix<dim,dim>& stressTens, const size_t ip,
149 const MathMatrix<dim, dim>& GradU)
150{
151 PRANDTL_REUSS_PROFILE_BEGIN(PrandtlReuss_stressTensor);
152
153 // get internal variables
154 MathMatrix<dim, dim>& strain_p_old_t = m_pElemData->internalVars[ip].strain_p_old_t;
155 number alpha = m_pElemData->internalVars[ip].alpha;
156
157 StressTensor(stressTens, GradU, strain_p_old_t, alpha);
158}
159
160
161template <typename TDomain>
164elasticityTensor(const size_t ip, const MathMatrix<dim, dim>& GradU_const)
165{
166 PRANDTL_REUSS_PROFILE_BEGIN(PrandtlReuss_elasticityTensor);
167
168 // get internal variables
169 MathMatrix<dim, dim>& strain_p_old_t = m_pElemData->internalVars[ip].strain_p_old_t;
170 number alpha = m_pElemData->internalVars[ip].alpha;
171
172 MathMatrix<dim,dim>& GradU = const_cast<MathMatrix<dim,dim>&>(GradU_const);
173
174 // computes the elasticity Tensor by means of numerical differentiation
176 MathMatrix<dim, dim> stressT, stressTT;
177
178 // for formulation in reference config
179 StressTensor(stressT, GradU, strain_p_old_t, alpha);
180
181 for (size_t i = 0; i < (size_t) dim; ++i)
182 for (size_t j = 0; j < (size_t) dim; ++j)
183 {
184 GradU[i][j] += m_tangentAccur;
185
186 StressTensor(stressTT, GradU, strain_p_old_t, alpha);
187
188 for (size_t k = 0; k < (size_t) dim; ++k)
189 for (size_t l = 0; l < (size_t) dim; ++l)
190 elastTens[i][j][k][l] = 1.0/ m_tangentAccur * (stressTT[k][l] - stressT[k][l]);
191
192 GradU[i][j] -= m_tangentAccur;
193 }
194
195 // compute maximal relative error of numerical differentiation
196 // (ref: Deuflhard Numerische Mathematik 1)
197 MathMatrix<dim, dim> diffStress;
198 MatSubtract(diffStress, stressTT, stressT);
199 if (MatFrobeniusNorm(stressT) > 0.0)
200 {
201 const number k_tan = MatFrobeniusNorm(diffStress)/MatFrobeniusNorm(stressT);
202 if (k_tan > m_max_k_tan)
203 m_max_k_tan = k_tan;
204 if (k_tan < m_min_k_tan)
205 m_min_k_tan = k_tan;
206 }
207
208 // TODO: change this smart pointer to a member variable
209 // do the same with m_GradU!
211 return spElastTens;
212}
213
214template <typename TDomain>
215void
217init_internal_vars(TBaseElem* elem, const size_t numIP)
218{
219 m_aaElemData[elem].internalVars.resize(numIP);
220
221 // set plastic strain (eps_p) and hardening variable (alpha)
222 // to zero at every ip (in ElemData-struct)
223 for (size_t ip = 0; ip < numIP; ++ip)
224 {
225 m_aaElemData[elem].internalVars[ip].strain_p_old_t = 0.0;
226 m_aaElemData[elem].internalVars[ip].alpha = 0.0;
227 }
228
229 if (!base_type::m_bInit)
230 base_type::m_bInit = true;
231}
232
233template <typename TDomain>
234inline
235void
238{
239 m_pElemData = &m_aaElemData[elem];
240}
241
242template <typename TDomain>
243void
245update_internal_vars(const size_t ip, const MathMatrix<dim, dim>& GradU)
246{
247 MathMatrix<dim, dim>& strain_p_old_t = m_pElemData->internalVars[ip].strain_p_old_t;
248 // TODO use a reference for alpha here!
249 number alpha = m_pElemData->internalVars[ip].alpha;
250
251 // update the internal variables: strain_p_old_t (plastic strain tensor)
252 // and alpha (hardening parameter) at ip
253 Update_internal_vars(strain_p_old_t, alpha, GradU, strain_p_old_t);
254 m_pElemData->internalVars[ip].alpha = alpha;
255}
256
257template <typename TDomain>
258void
261 number& alpha,
262 const MathMatrix<dim, dim>& GradU,
263 const MathMatrix<dim, dim>& strain_p_old_t)
264{
265 // compute trial strain tensor (eps)
266 MathMatrix<dim, dim> strainTensTrial;
267 for(size_t i = 0; i < (size_t) dim; ++i)
268 for(size_t j = 0; j < (size_t) dim; ++j)
269 strainTensTrial[i][j] = 0.5 * (GradU[i][j] + GradU[j][i]) - strain_p_old_t[i][j];
270
271 MathMatrix<dim, dim> dev_strainTrial;
272 MatDeviatorTrace(strainTensTrial, dev_strainTrial);
273
275 for(size_t i = 0; i < (size_t) dim; ++i)
276 for(size_t j = 0; j < (size_t) dim; ++j)
277 strial[i][j] = 2.0 * matConsts.mu * dev_strainTrial[i][j];
278
279 number strialnorm = MatFrobeniusNorm(strial);
280 number flowcondtrial = strialnorm - sqrt(2.0 / 3.0) * Hardening(alpha);
281
282 if (flowcondtrial <= 0){
283 strain_p_new = strain_p_old_t; return;
284 }
285
286 m_plasticIPs += 1;
287
288 number gamma = 0.0;
289
290 UG_ASSERT(strialnorm > 0.0, "norm of strial needs to be > 0.0");
291
292 // computation of gamma (plastic corrector/multiplicator)
293 switch (m_hardening)
294 {
295 case 0: gamma = PerfectPlasticity(flowcondtrial); break;
296 case 1: gamma = LinearHardening(flowcondtrial); break;
297 case 2: gamma = ExponentialHardening(strialnorm, alpha); break;
298 default:
299 UG_THROW(m_hardening << " in 'Update' is not a valid hardening behavior! \n");
300 }
301
302 MathMatrix<dim, dim> normaltrial;
303 MatScale(normaltrial, 1.0 / strialnorm, strial);
304
305 UG_ASSERT(gamma > 0.0, "gamma needs to be > 0.0");
306
307 MathMatrix<dim, dim> d_strain_p;
308 MatScale(d_strain_p, gamma, normaltrial);
309 MatAdd(strain_p_new, strain_p_old_t, d_strain_p);
310
311 alpha += sqrt(2.0 / 3.0) * gamma;
312}
313
314template <typename TDomain>
315inline
316void
319{
320 for(size_t i = 0; i < (size_t) dim; ++i)
321 for(size_t j = 0; j < (size_t) dim; ++j)
322 strainTens[i][j] = 0.5 * (GradU[i][j] + GradU[j][i]);
323}
324
325template <typename TDomain>
326void
328Flowrule(MathMatrix<dim, dim>& strain_p_new, MathMatrix<dim, dim>& strain, number& gamma,
330 const MathMatrix<dim, dim>& strain_p_old_t, const number alpha)
331{
333 // TRIAL ELASTIC STEP
335
336 MathMatrix<dim, dim> strain_trial;
337
338 // compute linearized strain tensor (eps)
339 for(size_t i = 0; i < (size_t) dim; ++i)
340 for(size_t j = 0; j < (size_t) dim; ++j)
341 {
342 strain[i][j] = 0.5 * (GradU[i][j] + GradU[j][i]);
343 // get eps_trial := eps_{n+1} - eps^p_{n}
344 strain_trial[i][j] = strain[i][j] - strain_p_old_t[i][j];
345 }
346
347 MathMatrix<dim, dim> dev_strain_trial;
348 MatDeviatorTrace(strain_trial, dev_strain_trial);
349
350 // compute trial strain deviator
351 MatScale(strial, 2.0 * matConsts.mu, dev_strain_trial);
352
353 number strialnorm = MatFrobeniusNorm(strial);
354 number flowcondtrial = strialnorm - sqrt(2.0 / 3.0) * Hardening(alpha);
355
357 // CHECKING YIELD-CONDITION:
359
360 if (flowcondtrial <= 0){
361 strain_p_new = strain_p_old_t; gamma = 0.0;
362 return;
363 }
364
366 // RETURN-MAPPING (corrector-step)
368
369 UG_ASSERT(strialnorm > 0.0, "norm of strial needs to be > 0.0");
370
371 // computation of gamma (plastic corrector/multiplicator)
372 // accordingly to Simo/Hughes 98 p.121/122
373 switch (m_hardening)
374 {
375 case 0: gamma = PerfectPlasticity(flowcondtrial); break;
376 case 1: gamma = LinearHardening(flowcondtrial); break;
377 case 2: gamma = ExponentialHardening(strialnorm, alpha); break;
378 default:
379 UG_THROW(m_hardening << " in 'Flowrule' is not a valid hardening behavior! \n");
380 }
381
382 UG_ASSERT(gamma > 0.0, "gamma needs to be > 0.0");
383
384 MatScale(normal, 1.0 / strialnorm, strial);
385
386 MathMatrix<dim, dim> d_strain_p;
387 MatScale(d_strain_p, gamma, normal);
388 MatAdd(strain_p_new, strain_p_old_t, d_strain_p);
389}
390
391template <typename TDomain>
392void
394ConstLaw(MathMatrix<dim, dim>& stressTens, const MathMatrix<dim, dim>& strain, const MathMatrix<dim, dim>& strial,
395 const number& gamma, const MathMatrix<dim, dim>& normal)
396{
397 // constitutive law taken from Simo/Hughes 98 'Computational Inelasticity' p. 124
398 number trStrain = Trace(strain);
399
400 // compute sigma = kappa * tr[eps] * id + strial - 2 * mu * gamma * normal
401 for(size_t i = 0; i < (size_t) dim; ++i)
402 {
403 for(size_t j = 0; j < (size_t) dim; ++j)
404 stressTens[i][j] = strial[i][j] - 2.0 * matConsts.mu * gamma * normal[i][j];
405
406 stressTens[i][i] += matConsts.kappa * trStrain;
407 }
408
409}
410
411template <typename TDomain>
412number
414plastic_multiplier(const size_t ip, const MathMatrix<dim, dim>& GradU)
415{
416 // get internal variables
417 MathMatrix<dim, dim>& strain_p_old_t = m_pElemData->internalVars[ip].strain_p_old_t;
418 number alpha = m_pElemData->internalVars[ip].alpha;
419
421 // TRIAL ELASTIC STEP
423
424 MathMatrix<dim, dim> strain_trial;
425
426 // compute linearized strain tensor (eps)
427 for(size_t i = 0; i < (size_t) dim; ++i)
428 for(size_t j = 0; j < (size_t) dim; ++j)
429 {
430 // get eps_trial := eps_{n+1} - eps^p_{n}
431 strain_trial[i][j] = 0.5 * (GradU[i][j] + GradU[j][i]) - strain_p_old_t[i][j];
432 }
433
434 MathMatrix<dim, dim> dev_strain_trial;
435 MatDeviatorTrace(strain_trial, dev_strain_trial);
436
437 // compute trial strain deviator
439 MatScale(strial, 2.0 * matConsts.mu, dev_strain_trial);
440
441 number strialnorm = MatFrobeniusNorm(strial);
442 number flowcondtrial = strialnorm - sqrt(2.0 / 3.0) * Hardening(alpha);
443
445 // CHECKING YIELD-CONDITION:
447
448 if (flowcondtrial <= 0){
449 return 0.0;
450 }
451
453 // RETURN-MAPPING (corrector-step)
455
456 UG_ASSERT(strialnorm > 0.0, "norm of strial needs to be > 0.0");
457
458 number gamma = 0.0;
459
460 // computation of gamma (plastic corrector/multiplicator)
461 // accordingly to Simo/Hughes 98 p.121/122
462 switch (m_hardening)
463 {
464 case 0: gamma = PerfectPlasticity(flowcondtrial); break;
465 case 1: gamma = LinearHardening(flowcondtrial); break;
466 case 2: gamma = ExponentialHardening(strialnorm, alpha); break;
467 default:
468 UG_THROW(m_hardening << " in 'Flowrule' is not a valid hardening behavior! \n");
469 }
470
471 UG_ASSERT(gamma > 0.0, "gamma needs to be > 0.0");
472
473 return gamma;
474}
475
476template <typename TDomain>
477inline
478number
480Hardening(const number alpha)
481{
482 return (matConsts.K_0 + matConsts.Hard * alpha +
483 (matConsts.K_inf - matConsts.K_0) * (1.0 - exp(-matConsts.omega * alpha)));
484}
485
486// "Hardening_d": derivative of nonlinear hardening law "Hardening"
487template <typename TDomain>
488inline
489number
491Hardening_d(const number alpha)
492{
493 return (matConsts.Hard + (matConsts.K_inf - matConsts.K_0) * matConsts.omega
494 * exp(-matConsts.omega * alpha));
495}
496
497template <typename TDomain>
498inline
499number
501PerfectPlasticity(const number flowcondtrial)
502{
503 return flowcondtrial / (2.0 * matConsts.mu);
504}
505
506template <typename TDomain>
507inline
508number
510LinearHardening(const number flowcondtrial)
511{
512 return flowcondtrial / (2.0 * (matConsts.mu + matConsts.Hard/3.0));
513}
514
515template <typename TDomain>
516number
518ExponentialHardening(const number strialnorm, const number alpha)
519{
520 number gamma = 0.0;
521
522 for (size_t i = 0; i < m_MaxHardIter; ++i) {
523 // f_cap has to be equal to 0, due to the consistency condition
524 number f_cap = strialnorm - 2.0 * gamma * matConsts.mu - sqrt(2.0 / 3.0) * Hardening(
525 alpha + sqrt(2.0 / 3.0) * gamma);
526
527 //abs(f_cap) < HardAccuracy * strialnorm vs. - f_cap > - HardAccuracy * strialnorm
528 if (f_cap < m_HardAccuracy * strialnorm){ break;}
529
530 gamma += f_cap / (2.0 * matConsts.mu * (1.0 + Hardening_d(
531 alpha + sqrt(2.0 / 3.0) * gamma) / (3.0 * matConsts.mu)));
532 }
533
534 return gamma;
535}
536
537template <typename TDomain>
538void
541{
542 UG_LOG("maximal k_tan:" << m_max_k_tan << "\n");
543 UG_LOG("minimal k_tan:" << m_min_k_tan << "\n");
544
545 // print: at how many gauss points we are in the plastic zone,...
546 UG_LOG("# of plastic IPs in this timestep:" << m_plasticIPs << "\n");
547 // reset data
548 m_plasticIPs = 0;
549}
550
551}// end of namespace SmallStrainMechanics
552}// end of namespace ug
553
554#endif /* PRANDTL_REUSS_IMPL_H_ */
Definition mat_law_interface.h:44
virtual SmartPtr< MathTensor4< dim, dim, dim, dim > > elasticityTensor()
Definition mat_law_interface.h:84
bool m_bConstElastTens
flag indicating, if elasticity tensor is constant
Definition mat_law_interface.h:131
std::string m_materialConfiguration
Definition mat_law_interface.h:124
matrix_t::value_type MatDeviatorTrace(const matrix_t &m, matrix_t &dev)
void MatSubtract(matrix_t &mOut, const matrix_t &m, typename matrix_t::value_type s)
MathMatrix< 1, 1, T >::value_type Trace(const MathMatrix< 1, 1, T > &m)
void MatAdd(matrix_t &mOut, const matrix_t &m, typename matrix_t::value_type s)
void MatScale(matrix_t &mOut, typename matrix_t::value_type s, const matrix_t &m)
matrix_t::value_type MatFrobeniusNorm(matrix_t &m)
PrandtlReuss()
constructor
Definition prandtl_reuss_impl.h:53
number m_min_k_tan
Definition prandtl_reuss.h:266
number Hardening(const number alpha)
Definition prandtl_reuss_impl.h:480
void strainTensor(MathMatrix< dim, dim > &strainTens, const MathMatrix< dim, dim > &GradU)
Definition prandtl_reuss_impl.h:318
void ConstLaw(MathMatrix< dim, dim > &stressTens, const MathMatrix< dim, dim > &strainTens, const MathMatrix< dim, dim > &strial, const number &gamma, const MathMatrix< dim, dim > &normal)
Definition prandtl_reuss_impl.h:394
void set_hardening_behavior(int hard)
Definition prandtl_reuss_impl.h:78
void StressTensor(MathMatrix< dim, dim > &stressTens, const MathMatrix< dim, dim > &GradU, const MathMatrix< dim, dim > &strain_p_old_t, const number alpha)
Definition prandtl_reuss_impl.h:135
struct ug::SmallStrainMechanics::PrandtlReuss::MaterialConstants matConsts
virtual void update_internal_vars(const size_t ip, const MathMatrix< dim, dim > &GradU)
Definition prandtl_reuss_impl.h:245
void Flowrule(MathMatrix< dim, dim > &strain_p_new, MathMatrix< dim, dim > &strain, number &gamma, MathMatrix< dim, dim > &strial, MathMatrix< dim, dim > &normal, const MathMatrix< dim, dim > &GradU, const MathMatrix< dim, dim > &strain_p_old_t, const number alpha)
Definition prandtl_reuss_impl.h:328
number ExponentialHardening(const number strialnorm, const number alpha)
Definition prandtl_reuss_impl.h:518
number PerfectPlasticity(const number flowcondtrial)
Definition prandtl_reuss_impl.h:501
virtual void write_data_to_console(const number t)
Definition prandtl_reuss_impl.h:540
base_type::TBaseElem TBaseElem
base element type
Definition prandtl_reuss.h:83
void stressTensor(MathMatrix< dim, dim > &stressTens, const size_t ip, const MathMatrix< dim, dim > &GradU)
computes the cauchy stress tensor sigma at an integration point 'ip'
Definition prandtl_reuss_impl.h:148
number Hardening_d(const number alpha)
Definition prandtl_reuss_impl.h:491
number m_max_k_tan
max condition number for numerical approximated matrix
Definition prandtl_reuss.h:265
virtual void internal_vars(TBaseElem *elem)
Definition prandtl_reuss_impl.h:237
virtual number plastic_multiplier(const size_t ip, const MathMatrix< dim, dim > &GradU)
Definition prandtl_reuss_impl.h:414
number LinearHardening(const number flowcondtrial)
Definition prandtl_reuss_impl.h:510
int m_hardening
hardening behavior
Definition prandtl_reuss.h:244
virtual void init_internal_vars(TBaseElem *elem, const size_t numIP)
Definition prandtl_reuss_impl.h:217
void init()
Definition prandtl_reuss_impl.h:113
size_t m_plasticIPs
Definition prandtl_reuss.h:267
void Update_internal_vars(MathMatrix< dim, dim > &strain_p_new, number &alpha, const MathMatrix< dim, dim > &GradU, const MathMatrix< dim, dim > &strain_p_old_t)
Definition prandtl_reuss_impl.h:260
#define UG_ASSERT(expr, msg)
#define UG_THROW(msg)
#define UG_LOG(msg)
double number
#define PRANDTL_REUSS_PROFILE_BEGIN(name)
Definition prandtl_reuss_impl.h:41