Plugins
prandtl_reuss.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_H_
34 #define PRANDTL_REUSS_H_
35 
36 namespace ug{
37 namespace SmallStrainMechanics{
38 
43 
67 template <typename TDomain>
69  : public IMaterialLaw<TDomain>
70 {
71  private:
74 
77 
78  public:
80  static const int dim = base_type::dim;
81 
83  typedef typename base_type::TBaseElem TBaseElem;
84 
85  public:
87  PrandtlReuss();
88 
91 
93  void set_bulk_modulus(const number bulkModulus)
94  {matConsts.kappa = bulkModulus;
95  std::stringstream ss; ss << m_materialConfiguration << " bulk modulus kappa = " << bulkModulus << "\n";
96  m_materialConfiguration = ss.str();}
97 
98  void set_shear_modulus(const number shearModulus)
99  {matConsts.mu = shearModulus;
100  std::stringstream ss; ss << m_materialConfiguration << " shear modulus mu = " << shearModulus << "\n";
101  m_materialConfiguration = ss.str();}
102 
103  void set_initial_flow_stress(const number initialFlowStress)
104  {matConsts.K_0 = initialFlowStress;
105  std::stringstream ss; ss << m_materialConfiguration << " initial flow-stress K_0 = " << initialFlowStress << "\n";
106  m_materialConfiguration = ss.str();}
107 
108  void set_residual_flow_stress(const number resFlowStress)
109  {matConsts.K_inf = resFlowStress;
110  std::stringstream ss; ss << m_materialConfiguration << " residual flow-stress (saturation stress) K_inf = " << resFlowStress << "\n";
111  m_materialConfiguration = ss.str();}
112 
113  void set_hardening_modulus(const number hardModulus)
114  {matConsts.Hard = hardModulus; m_bHardModulus = true;
115  std::stringstream ss; ss << m_materialConfiguration << " hardening modulus = " << hardModulus << "\n";
116  m_materialConfiguration = ss.str();}
117 
118  void set_hardening_exponent(const number hardExponent)
119  {matConsts.omega = hardExponent; m_bHardExp = true;
120  std::stringstream ss; ss << m_materialConfiguration << " hardening exponent = " << hardExponent << "\n";
121  m_materialConfiguration = ss.str();}
122 
123  void set_hardening_behavior(int hard);
124 
126  void set_tangent_precision(const number tanAccur)
127  {m_tangentAccur = tanAccur;
128  std::stringstream ss; ss << m_materialConfiguration << " accuracy of the tangent approximation = " << tanAccur << "\n";
129  m_materialConfiguration = ss.str();}
130 
131 
132  public:
134  // INTERFACE-METHODS
136  void init();
137 
139  void stressTensor(MathMatrix<dim,dim>& stressTens, const size_t ip,
140  const MathMatrix<dim, dim>& GradU);
141 
144  elasticityTensor(const size_t ip, const MathMatrix<dim, dim>& GradU);
145 
146  virtual bool needs_to_add_jac_m(){return false;}
147 
148 
149  virtual void init_internal_vars(TBaseElem* elem, const size_t numIP);
150  virtual void internal_vars(TBaseElem* elem);
151  virtual void update_internal_vars(const size_t ip, const MathMatrix<dim, dim>& GradU);
152 
153  virtual void write_data_to_console(const number t);
154 
155 
157  MathMatrix<dim, dim>& inelastStrain = m_pElemData->internalVars[ip].strain_p_old_t;
158  SmartPtr<MathMatrix<dim, dim> > spInelasticStrain (new MathMatrix<dim, dim>(inelastStrain));
159  return spInelasticStrain;
160  }
161  virtual number hardening_parameter(const size_t ip){
162  return m_pElemData->internalVars[ip].alpha;
163  }
164  virtual number plastic_multiplier(const size_t ip, const MathMatrix<dim, dim>& GradU);
165 
167 
169  {
170  //todo: Move this to an initialization routine. Best would probably to make
171  // set_domain(...) of the base class virtual and to attach m_aElemData there.
172 
173  grid.template attach_to<TBaseElem>(m_aElemData);
174  m_aaElemData.access(grid, m_aElemData);
175  }
176 
177  virtual void clear_attachments(typename TDomain::grid_type& grid)
178  {
179  if(grid.template has_attachment<TBaseElem>(m_aElemData)){
180  grid.template detach_from<TBaseElem>(m_aElemData);
181  m_aaElemData.invalidate();
182  }
183  }
184 
185  private:
186  void strainTensor(MathMatrix<dim,dim>& strainTens, const MathMatrix<dim, dim>& GradU);
187 
188  void Flowrule(MathMatrix<dim, dim>& strain_p_new,
189  MathMatrix<dim, dim>& strain,
190  number& gamma,
191  MathMatrix<dim, dim>& strial,
193  const MathMatrix<dim, dim>& GradU,
194  const MathMatrix<dim, dim>& strain_p_old_t,
195  const number alpha);
196 
197  void ConstLaw(MathMatrix<dim, dim>& stressTens, const MathMatrix<dim, dim>& strainTens,
198  const MathMatrix<dim, dim>& strial, const number& gamma,
200 
201  void StressTensor(MathMatrix<dim,dim>& stressTens, const MathMatrix<dim, dim>& GradU,
202  const MathMatrix<dim, dim>& strain_p_old_t, const number alpha);
203 
204  number Hardening(const number alpha);
205 
206  number Hardening_d(const number alpha);
207 
208  number PerfectPlasticity(const number flowcondtrial);
209 
210  number LinearHardening(const number flowcondtrial);
211 
212  number ExponentialHardening(const number strialnorm, const number alpha);
213 
214  void Update_internal_vars(MathMatrix<dim, dim>& strain_p_new,
215  number& alpha,
216  const MathMatrix<dim, dim>& GradU,
217  const MathMatrix<dim, dim>& strain_p_old_t);
218 
219  private:
221  struct InternalVars {
222  // plastic part of strain tensor wrt to last time step
224 
225  // hardening variable
227  };
228 
229  // std-vector of InternalVars
230  struct ElemData{
231  std::vector<InternalVars> internalVars;
232  };
233 
235 
236  // attachment type: attachment of ElemData
239  // the attachment accessor
240  typedef Grid::AttachmentAccessor<TBaseElem, AElemData> ElemDataAccessor;
242 
247 
249  number mu; // shear modulus
250  number kappa; // bulk modulus
251  number K_0; // initial flow-stress:
252  number K_inf; // residual flow-stress (saturation stress)
253  number Hard; // linear hardening modulus
254  number omega; // hardening exponent
256 
259 
263 
267  size_t m_plasticIPs;
268 
269  public:
271 };
272 
273 }// end of namespace SmallStrainMechanics
274 }// end of namespace ug
275 
276 #include "prandtl_reuss_impl.h"
277 
278 #endif /* PRANDTL_REUSS_H_ */
Definition: mat_law_interface.h:44
static const int dim
World dimension.
Definition: mat_law_interface.h:47
virtual SmartPtr< MathTensor4< dim, dim, dim, dim > > elasticityTensor()
Definition: mat_law_interface.h:84
std::string m_materialConfiguration
Definition: mat_law_interface.h:124
domain_traits< TDomain::dim >::grid_base_object TBaseElem
base element type of associated domain
Definition: mat_law_interface.h:50
Definition: prandtl_reuss.h:70
std::string m_materialConfiguration
Definition: mat_law_interface.h:124
SmartPtr< TGrid > grid()
TGrid grid_type
PrandtlReuss()
constructor
Definition: prandtl_reuss_impl.h:53
void set_tangent_precision(const number tanAccur)
set precision of numerical approximation of the tangent
Definition: prandtl_reuss.h:126
MathMatrix< dim, dim > strain_p_old_t
Definition: prandtl_reuss.h:223
virtual number hardening_parameter(const size_t ip)
Definition: prandtl_reuss.h:161
void set_initial_flow_stress(const number initialFlowStress)
Definition: prandtl_reuss.h:103
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
ElemDataAccessor m_aaElemData
Definition: prandtl_reuss.h:241
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
~PrandtlReuss()
Destructor.
Definition: prandtl_reuss.h:90
number Hard
Definition: prandtl_reuss.h:253
Attachment< ElemData > AElemData
Definition: prandtl_reuss.h:237
void set_residual_flow_stress(const number resFlowStress)
Definition: prandtl_reuss.h:108
number m_HardAccuracy
Definition: prandtl_reuss.h:246
number m_tangentAccur
tangent accuracy
Definition: prandtl_reuss.h:258
std::vector< InternalVars > internalVars
Definition: prandtl_reuss.h:231
number omega
Definition: prandtl_reuss.h:254
size_t m_MaxHardIter
Definition: prandtl_reuss.h:245
bool m_bHardModulus
flags indicating if hardening variables are set
Definition: prandtl_reuss.h:261
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
bool m_bHardExp
Definition: prandtl_reuss.h:262
virtual void attach_internal_vars(typename TDomain::grid_type &grid)
use this method to make sure that all required attachments are attached
Definition: prandtl_reuss.h:168
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
void set_hardening_exponent(const number hardExponent)
Definition: prandtl_reuss.h:118
void set_shear_modulus(const number shearModulus)
Definition: prandtl_reuss.h:98
ElemData * m_pElemData
Definition: prandtl_reuss.h:234
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
AElemData m_aElemData
Definition: prandtl_reuss.h:238
virtual void write_data_to_console(const number t)
Definition: prandtl_reuss_impl.h:540
number alpha
Definition: prandtl_reuss.h:226
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
virtual void clear_attachments(typename TDomain::grid_type &grid)
Definition: prandtl_reuss.h:177
number Hardening_d(const number alpha)
Definition: prandtl_reuss_impl.h:491
void set_hardening_modulus(const number hardModulus)
Definition: prandtl_reuss.h:113
PrandtlReuss< TDomain > this_type
own type
Definition: prandtl_reuss.h:76
number kappa
Definition: prandtl_reuss.h:250
number K_inf
Definition: prandtl_reuss.h:252
virtual bool needs_to_add_jac_m()
Definition: prandtl_reuss.h:146
number m_max_k_tan
max condition number for numerical approximated matrix
Definition: prandtl_reuss.h:265
virtual SmartPtr< MathMatrix< dim, dim > > inelastic_strain_tensor(const size_t ip)
Definition: prandtl_reuss.h:156
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
static const int dim
World dimension.
Definition: prandtl_reuss.h:80
IMaterialLaw< TDomain > base_type
Base class type.
Definition: prandtl_reuss.h:73
number LinearHardening(const number flowcondtrial)
Definition: prandtl_reuss_impl.h:510
void set_bulk_modulus(const number bulkModulus)
set-methods for material constants
Definition: prandtl_reuss.h:93
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
Grid::AttachmentAccessor< TBaseElem, AElemData > ElemDataAccessor
Definition: prandtl_reuss.h:240
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
normal
double number
attached ElemData
Definition: prandtl_reuss.h:221