www.mooseframework.org
FiniteStrainCrystalPlasticity.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "petscblaslapack.h"
12 #include "libmesh/utility.h"
13 
14 #include <fstream>
15 #include <cmath>
16 
18 
20 
21 InputParameters
23 {
24  InputParameters params = ComputeStressBase::validParams();
25  params.addClassDescription(
26  "Crystal Plasticity base class: FCC system with power law flow rule implemented");
27  params.addRequiredParam<int>("nss", "Number of slip systems");
28  params.addParam<std::vector<Real>>("gprops", "Initial values of slip system resistances");
29  params.addParam<std::vector<Real>>("hprops", "Hardening properties");
30  params.addParam<std::vector<Real>>("flowprops", "Parameters used in slip rate equations");
31  params.addRequiredParam<FileName>("slip_sys_file_name",
32  "Name of the file containing the slip system");
33  params.addParam<FileName>(
34  "slip_sys_res_prop_file_name",
35  "",
36  "Name of the file containing the initial values of slip system resistances");
37  params.addParam<FileName>(
38  "slip_sys_flow_prop_file_name",
39  "",
40  "Name of the file containing the values of slip rate equation parameters");
41  params.addParam<FileName>(
42  "slip_sys_hard_prop_file_name",
43  "",
44  "Name of the file containing the values of hardness evolution parameters");
45  params.addParam<Real>("rtol", 1e-6, "Constitutive stress residue relative tolerance");
46  params.addParam<Real>("abs_tol", 1e-6, "Constitutive stress residue absolute tolerance");
47  params.addParam<Real>("gtol", 1e2, "Constitutive slip system resistance residual tolerance");
48  params.addParam<Real>("slip_incr_tol", 2e-2, "Maximum allowable slip in an increment");
49  params.addParam<unsigned int>("maxiter", 100, "Maximum number of iterations for stress update");
50  params.addParam<unsigned int>(
51  "maxitergss", 100, "Maximum number of iterations for slip system resistance update");
52  params.addParam<unsigned int>(
53  "num_slip_sys_flowrate_props",
54  2,
55  "Number of flow rate properties for a slip system"); // Used for reading flow rate parameters
56  params.addParam<UserObjectName>("read_prop_user_object",
57  "The ElementReadPropertyFile "
58  "GeneralUserObject to read element "
59  "specific property values from file");
60  MooseEnum tan_mod_options("exact none", "none"); // Type of read
61  params.addParam<MooseEnum>("tan_mod_type",
62  tan_mod_options,
63  "Type of tangent moduli for preconditioner: default elastic");
64  MooseEnum intvar_read_options("slip_sys_file slip_sys_res_file none", "none");
65  params.addParam<MooseEnum>(
66  "intvar_read_type",
67  intvar_read_options,
68  "Read from options for initial value of internal variables: Default from .i file");
69  params.addParam<unsigned int>("num_slip_sys_props",
70  0,
71  "Number of slip system specific properties provided in the file "
72  "containing slip system normals and directions");
73  params.addParam<bool>(
74  "gen_random_stress_flag",
75  false,
76  "Flag to generate random stress to perform time cutback on constitutive failure");
77  params.addParam<bool>("input_random_scaling_var",
78  false,
79  "Flag to input scaling variable: _Cijkl(0,0,0,0) when false");
80  params.addParam<Real>("random_scaling_var",
81  1e9,
82  "Random scaling variable: Large value can cause non-positive definiteness");
83  params.addParam<unsigned int>(
84  "random_seed",
85  2000,
86  "Random integer used to generate random stress when constitutive failure occurs");
87  params.addParam<unsigned int>(
88  "maximum_substep_iteration", 1, "Maximum number of substep iteration");
89  params.addParam<bool>("use_line_search", false, "Use line search in constitutive update");
90  params.addParam<Real>("min_line_search_step_size", 0.01, "Minimum line search step size");
91  params.addParam<Real>("line_search_tol", 0.5, "Line search bisection method tolerance");
92  params.addParam<unsigned int>(
93  "line_search_maxiter", 20, "Line search bisection method maximum number of iteration");
94  MooseEnum line_search_method("CUT_HALF BISECTION", "CUT_HALF");
95  params.addParam<MooseEnum>(
96  "line_search_method", line_search_method, "The method used in line search");
97 
98  return params;
99 }
100 
102  : ComputeStressBase(parameters),
103  _nss(getParam<int>("nss")),
104  _gprops(getParam<std::vector<Real>>("gprops")),
105  _hprops(getParam<std::vector<Real>>("hprops")),
106  _flowprops(getParam<std::vector<Real>>("flowprops")),
107  _slip_sys_file_name(getParam<FileName>("slip_sys_file_name")),
108  _slip_sys_res_prop_file_name(getParam<FileName>("slip_sys_res_prop_file_name")),
109  _slip_sys_flow_prop_file_name(getParam<FileName>("slip_sys_flow_prop_file_name")),
110  _slip_sys_hard_prop_file_name(getParam<FileName>("slip_sys_hard_prop_file_name")),
111  _rtol(getParam<Real>("rtol")),
112  _abs_tol(getParam<Real>("abs_tol")),
113  _gtol(getParam<Real>("gtol")),
114  _slip_incr_tol(getParam<Real>("slip_incr_tol")),
115  _maxiter(getParam<unsigned int>("maxiter")),
116  _maxiterg(getParam<unsigned int>("maxitergss")),
117  _num_slip_sys_flowrate_props(getParam<unsigned int>("num_slip_sys_flowrate_props")),
118  _tan_mod_type(getParam<MooseEnum>("tan_mod_type")),
119  _intvar_read_type(getParam<MooseEnum>("intvar_read_type")),
120  _num_slip_sys_props(getParam<unsigned int>("num_slip_sys_props")),
121  _gen_rndm_stress_flag(getParam<bool>("gen_random_stress_flag")),
122  _input_rndm_scale_var(getParam<bool>("input_random_scaling_var")),
123  _rndm_scale_var(getParam<Real>("random_scaling_var")),
124  _rndm_seed(getParam<unsigned int>("random_seed")),
125  _max_substep_iter(getParam<unsigned int>("maximum_substep_iteration")),
126  _use_line_search(getParam<bool>("use_line_search")),
127  _min_lsrch_step(getParam<Real>("min_line_search_step_size")),
128  _lsrch_tol(getParam<Real>("line_search_tol")),
129  _lsrch_max_iter(getParam<unsigned int>("line_search_maxiter")),
130  _lsrch_method(getParam<MooseEnum>("line_search_method")),
131  _fp(declareProperty<RankTwoTensor>("fp")), // Plastic deformation gradient
132  _fp_old(getMaterialPropertyOld<RankTwoTensor>(
133  "fp")), // Plastic deformation gradient of previous increment
134  _pk2(declareProperty<RankTwoTensor>("pk2")), // 2nd Piola-Kirchoff Stress
135  _pk2_old(getMaterialPropertyOld<RankTwoTensor>(
136  "pk2")), // 2nd Piola Kirchoff Stress of previous increment
137  _lag_e(declareProperty<RankTwoTensor>("lage")), // Lagrangian strain
138  _lag_e_old(
139  getMaterialPropertyOld<RankTwoTensor>("lage")), // Lagrangian strain of previous increment
140  _gss(declareProperty<std::vector<Real>>("gss")), // Slip system resistances
141  _gss_old(getMaterialPropertyOld<std::vector<Real>>(
142  "gss")), // Slip system resistances of previous increment
143  _acc_slip(declareProperty<Real>("acc_slip")), // Accumulated slip
144  _acc_slip_old(
145  getMaterialPropertyOld<Real>("acc_slip")), // Accumulated slip of previous increment
146  _update_rot(declareProperty<RankTwoTensor>(
147  "update_rot")), // Rotation tensor considering material rotation and crystal orientation
148  _deformation_gradient(getMaterialProperty<RankTwoTensor>("deformation_gradient")),
149  _deformation_gradient_old(getMaterialPropertyOld<RankTwoTensor>("deformation_gradient")),
150  _elasticity_tensor_name(_base_name + "elasticity_tensor"),
151  _elasticity_tensor(getMaterialPropertyByName<RankFourTensor>(_elasticity_tensor_name)),
152  _crysrot(getMaterialProperty<RankTwoTensor>("crysrot")),
153  _mo(_nss * LIBMESH_DIM),
154  _no(_nss * LIBMESH_DIM),
155  _slip_incr(_nss),
156  _tau(_nss),
157  _dslipdtau(_nss),
158  _s0(_nss),
159  _gss_tmp(_nss),
160  _gss_tmp_old(_nss),
161  _dgss_dsliprate(_nss, _nss)
162 {
163  _err_tol = false;
164 
165  if (_num_slip_sys_props > 0)
167 
168  _pk2_tmp.zero();
169  _delta_dfgrd.zero();
170 
171  _first_step_iter = false;
172  _last_step_iter = false;
173  // Initialize variables in the first iteration of substepping
174  _first_substep = true;
175 
176  _read_from_slip_sys_file = false;
177  if (_intvar_read_type == "slip_sys_file")
179 
181  mooseError("Crystal Plasticity Error: Specify number of internal variable's initial values to "
182  "be read from slip system file");
183 
184  getSlipSystems();
185 
186  RankTwoTensor::initRandom(_rndm_seed);
187 }
188 
189 void
191 {
192  _stress[_qp].zero();
193 
194  _fp[_qp].setToIdentity();
195 
196  _pk2[_qp].zero();
197  _acc_slip[_qp] = 0.0;
198  _lag_e[_qp].zero();
199 
200  _update_rot[_qp].setToIdentity();
201 
202  initSlipSysProps(); // Initializes slip system related properties
204 }
205 
206 void
208 {
209  switch (_intvar_read_type)
210  {
211  case 0:
213  break;
214  case 1:
216  break;
217  default:
219  }
220 
221  if (_slip_sys_flow_prop_file_name.length() != 0)
223  else
225 
226  if (_slip_sys_hard_prop_file_name.length() != 0)
228  else
230 }
231 
232 void
234 {
235  _gss[_qp].resize(_nss);
236 
237  for (unsigned int i = 0; i < _nss; ++i)
238  _gss[_qp][i] = _slip_sys_props(i);
239 }
240 
241 // Read initial slip system resistances from .txt file. See test.
242 void
244 {
245  _gss[_qp].resize(_nss);
246 
247  MooseUtils::checkFileReadable(_slip_sys_res_prop_file_name);
248 
249  std::ifstream file;
250  file.open(_slip_sys_res_prop_file_name.c_str());
251 
252  for (unsigned int i = 0; i < _nss; ++i)
253  if (!(file >> _gss[_qp][i]))
254  mooseError("Error FiniteStrainCrystalPlasticity: Premature end of slip_sys_res_prop file");
255 
256  file.close();
257 }
258 
259 // Read initial slip system resistances from .i file
260 void
262 {
263  if (_gprops.size() <= 0)
264  mooseError("FiniteStrainCrystalPLasticity: Error in reading slip system resistance properties: "
265  "Specify input in .i file or in slip_sys_res_prop_file or in slip_sys_file");
266 
267  _gss[_qp].resize(_nss, 0.0);
268 
269  unsigned int num_data_grp = 3; // Number of data per group e.g. start_slip_sys, end_slip_sys,
270  // value
271 
272  for (unsigned int i = 0; i < _gprops.size() / num_data_grp; ++i)
273  {
274  Real vs, ve;
275  unsigned int is, ie;
276 
277  vs = _gprops[i * num_data_grp];
278  ve = _gprops[i * num_data_grp + 1];
279 
280  if (vs <= 0 || ve <= 0)
281  mooseError("FiniteStrainCrystalPLasticity: Indices in gss property read must be positive "
282  "integers: is = ",
283  vs,
284  " ie = ",
285  ve);
286 
287  if (vs != floor(vs) || ve != floor(ve))
288  mooseError("FiniteStrainCrystalPLasticity: Error in reading slip system resistances: Values "
289  "specifying start and end number of slip system groups should be integer");
290 
291  is = static_cast<unsigned int>(vs);
292  ie = static_cast<unsigned int>(ve);
293 
294  if (is > ie)
295  mooseError("FiniteStrainCrystalPLasticity: Start index is = ",
296  is,
297  " should be greater than end index ie = ",
298  ie,
299  " in slip system resistance property read");
300 
301  for (unsigned int j = is; j <= ie; ++j)
302  _gss[_qp][j - 1] = _gprops[i * num_data_grp + 2];
303  }
304 
305  for (unsigned int i = 0; i < _nss; ++i)
306  if (_gss[_qp][i] <= 0.0)
307  mooseError("FiniteStrainCrystalPLasticity: Value of resistance for slip system ",
308  i + 1,
309  " non positive");
310 }
311 
312 // Read flow rate parameters from .txt file. See test.
313 void
315 {
316  _a0.resize(_nss);
317  _xm.resize(_nss);
318 
319  MooseUtils::checkFileReadable(_slip_sys_flow_prop_file_name);
320 
321  std::ifstream file;
322  file.open(_slip_sys_flow_prop_file_name.c_str());
323 
324  std::vector<Real> vec;
325  vec.resize(_num_slip_sys_flowrate_props);
326 
327  for (unsigned int i = 0; i < _nss; ++i)
328  {
329  for (unsigned int j = 0; j < _num_slip_sys_flowrate_props; ++j)
330  if (!(file >> vec[j]))
331  mooseError(
332  "Error FiniteStrainCrystalPlasticity: Premature end of slip_sys_flow_rate_param file");
333 
334  _a0(i) = vec[0];
335  _xm(i) = vec[1];
336  }
337 
338  file.close();
339 }
340 
341 // Read flow rate parameters from .i file
342 void
344 {
345  if (_flowprops.size() <= 0)
346  mooseError("FiniteStrainCrystalPLasticity: Error in reading flow rate properties: Specify "
347  "input in .i file or a slip_sys_flow_prop_file_name");
348 
349  _a0.resize(_nss);
350  _xm.resize(_nss);
351 
352  unsigned int num_data_grp = 2 + _num_slip_sys_flowrate_props; // Number of data per group e.g.
353  // start_slip_sys, end_slip_sys,
354  // value1, value2, ..
355 
356  for (unsigned int i = 0; i < _flowprops.size() / num_data_grp; ++i)
357  {
358  Real vs, ve;
359  unsigned int is, ie;
360 
361  vs = _flowprops[i * num_data_grp];
362  ve = _flowprops[i * num_data_grp + 1];
363 
364  if (vs <= 0 || ve <= 0)
365  mooseError("FiniteStrainCrystalPLasticity: Indices in flow rate parameter read must be "
366  "positive integers: is = ",
367  vs,
368  " ie = ",
369  ve);
370 
371  if (vs != floor(vs) || ve != floor(ve))
372  mooseError("FiniteStrainCrystalPLasticity: Error in reading flow props: Values specifying "
373  "start and end number of slip system groups should be integer");
374 
375  is = static_cast<unsigned int>(vs);
376  ie = static_cast<unsigned int>(ve);
377 
378  if (is > ie)
379  mooseError("FiniteStrainCrystalPLasticity: Start index is = ",
380  is,
381  " should be greater than end index ie = ",
382  ie,
383  " in flow rate parameter read");
384 
385  for (unsigned int j = is; j <= ie; ++j)
386  {
387  _a0(j - 1) = _flowprops[i * num_data_grp + 2];
388  _xm(j - 1) = _flowprops[i * num_data_grp + 3];
389  }
390  }
391 
392  for (unsigned int i = 0; i < _nss; ++i)
393  {
394  if (!(_a0(i) > 0.0 && _xm(i) > 0.0))
395  {
396  mooseWarning(
397  "FiniteStrainCrystalPlasticity: Non-positive flow rate parameters ", _a0(i), ",", _xm(i));
398  break;
399  }
400  }
401 }
402 
403 // Read hardness parameters from .txt file
404 void
406 {
407 }
408 
409 // Read hardness parameters from .i file
410 void
412 {
413  if (_hprops.size() <= 0)
414  mooseError("FiniteStrainCrystalPLasticity: Error in reading hardness properties: Specify input "
415  "in .i file or a slip_sys_hard_prop_file_name");
416 
417  _r = _hprops[0];
418  _h0 = _hprops[1];
419  _tau_init = _hprops[2];
420  _tau_sat = _hprops[3];
421 }
422 
423 // Read slip systems from file
424 void
426 {
427  Real vec[LIBMESH_DIM];
428  std::ifstream fileslipsys;
429 
430  MooseUtils::checkFileReadable(_slip_sys_file_name);
431 
432  fileslipsys.open(_slip_sys_file_name.c_str());
433 
434  for (unsigned int i = 0; i < _nss; ++i)
435  {
436  // Read the slip normal
437  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
438  if (!(fileslipsys >> vec[j]))
439  mooseError("Crystal Plasticity Error: Premature end of file reading slip system file \n");
440 
441  // Normalize the vectors
442  Real mag;
443  mag = Utility::pow<2>(vec[0]) + Utility::pow<2>(vec[1]) + Utility::pow<2>(vec[2]);
444  mag = std::sqrt(mag);
445 
446  for (unsigned j = 0; j < LIBMESH_DIM; ++j)
447  _no(i * LIBMESH_DIM + j) = vec[j] / mag;
448 
449  // Read the slip direction
450  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
451  if (!(fileslipsys >> vec[j]))
452  mooseError("Crystal Plasticity Error: Premature end of file reading slip system file \n");
453 
454  // Normalize the vectors
455  mag = Utility::pow<2>(vec[0]) + Utility::pow<2>(vec[1]) + Utility::pow<2>(vec[2]);
456  mag = std::sqrt(mag);
457 
458  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
459  _mo(i * LIBMESH_DIM + j) = vec[j] / mag;
460 
461  mag = 0.0;
462  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
463  mag += _mo(i * LIBMESH_DIM + j) * _no(i * LIBMESH_DIM + j);
464 
465  if (std::abs(mag) > 1e-8)
466  mooseError(
467  "Crystal Plasicity Error: Slip direction and normal not orthonormal, System number = ",
468  i,
469  "\n");
470 
472  for (unsigned int j = 0; j < _num_slip_sys_props; ++j)
473  if (!(fileslipsys >> _slip_sys_props(i * _num_slip_sys_props + j)))
474  mooseError("Crystal Plasticity Error: Premature end of file reading slip system file - "
475  "check in slip system file read input options/values\n");
476  }
477 
478  fileslipsys.close();
479 }
480 
481 // Initialize addtional stateful material properties
482 void
484 {
485 }
486 
491 void
493 {
494  unsigned int substep_iter = 1; // Depth of substepping; Limited to maximum substep iteration
495  unsigned int num_substep = 1; // Calculated from substep_iter as 2^substep_iter
496  Real dt_original = _dt; // Stores original _dt; Reset at the end of solve
497  _first_substep = true; // Initialize variables at substep_iter = 1
498 
499  if (_max_substep_iter > 1)
500  {
502  if (_dfgrd_tmp_old.det() == 0)
503  _dfgrd_tmp_old.addIa(1.0);
504 
506  _err_tol = true; // Indicator to continue substepping
507  }
508 
509  // Substepping loop
510  while (_err_tol && _max_substep_iter > 1)
511  {
512  _dt = dt_original / num_substep;
513 
514  for (unsigned int istep = 0; istep < num_substep; ++istep)
515  {
516  _first_step_iter = false;
517  if (istep == 0)
518  _first_step_iter = true;
519 
520  _last_step_iter = false;
521  if (istep == num_substep - 1)
522  _last_step_iter = true;
523 
524  _dfgrd_scale_factor = (static_cast<Real>(istep) + 1) / num_substep;
525 
526  preSolveQp();
527  solveQp();
528 
529  if (_err_tol)
530  {
531  substep_iter++;
532  num_substep *= 2;
533  break;
534  }
535  }
536 
537  _first_substep = false; // Prevents reinitialization
538  _dt = dt_original; // Resets dt
539 
540 #ifdef DEBUG
541  if (substep_iter > _max_substep_iter)
542  mooseWarning("FiniteStrainCrystalPlasticity: Failure with substepping");
543 #endif
544 
545  if (!_err_tol || substep_iter > _max_substep_iter)
546  postSolveQp(); // Evaluate variables after successful solve or indicate failure
547  }
548 
549  // No substepping
550  if (_max_substep_iter == 1)
551  {
552  preSolveQp();
553  solveQp();
554  postSolveQp();
555  }
556 }
557 
558 void
560 {
561  // Initialize variable
562  if (_first_substep)
563  {
564  _Jacobian_mult[_qp].zero(); // Initializes jacobian for preconditioner
566  }
567 
568  if (_max_substep_iter == 1)
569  _dfgrd_tmp = _deformation_gradient[_qp]; // Without substepping
570  else
572 
573  _err_tol = false;
574 }
575 
576 void
578 {
580  solveStatevar();
581  if (_err_tol)
582  return;
584 }
585 
586 void
588 {
589  if (_err_tol)
590  {
591  _err_tol = false;
593  {
595  _rndm_scale_var = _elasticity_tensor[_qp](0, 0, 0, 0);
596 
597  _stress[_qp] = RankTwoTensor::genRandomSymmTensor(_rndm_scale_var, 1.0);
598  }
599  else
600  mooseError("FiniteStrainCrystalPlasticity: Constitutive failure");
601  }
602  else
603  {
604  _stress[_qp] = _fe * _pk2[_qp] * _fe.transpose() / _fe.det();
605 
606  _Jacobian_mult[_qp] += calcTangentModuli(); // Calculate jacobian for preconditioner
607 
608  RankTwoTensor iden(RankTwoTensor::initIdentity);
609 
610  _lag_e[_qp] = _deformation_gradient[_qp].transpose() * _deformation_gradient[_qp] - iden;
611  _lag_e[_qp] = _lag_e[_qp] * 0.5;
612 
613  RankTwoTensor rot;
614  rot = get_current_rotation(_deformation_gradient[_qp]); // Calculate material rotation
615  _update_rot[_qp] = rot * _crysrot[_qp];
616  }
617 }
618 
619 void
621 {
622  if (_max_substep_iter == 1) // No substepping
623  {
624  _gss_tmp = _gss_old[_qp];
626  }
627  else
628  {
629  if (_first_step_iter)
630  {
631  _gss_tmp = _gss_tmp_old = _gss_old[_qp];
633  }
634  else
636  }
637 }
638 
639 void
641 {
642  Real gmax, gdiff;
643  unsigned int iterg;
644  std::vector<Real> gss_prev(_nss);
645 
646  gmax = 1.1 * _gtol;
647  iterg = 0;
648 
649  while (gmax > _gtol && iterg < _maxiterg) // Check for slip system resistance update tolerance
650  {
651  preSolveStress();
652  solveStress();
653  if (_err_tol)
654  return;
655  postSolveStress();
656 
657  gss_prev = _gss_tmp;
658 
659  update_slip_system_resistance(); // Update slip system resistance
660 
661  gmax = 0.0;
662  for (unsigned i = 0; i < _nss; ++i)
663  {
664  gdiff = std::abs(gss_prev[i] - _gss_tmp[i]); // Calculate increment size
665 
666  if (gdiff > gmax)
667  gmax = gdiff;
668  }
669  iterg++;
670  }
671 
672  if (iterg == _maxiterg)
673  {
674 #ifdef DEBUG
675  mooseWarning("FiniteStrainCrystalPLasticity: Hardness Integration error gmax", gmax, "\n");
676 #endif
677  _err_tol = true;
678  }
679 }
680 
681 void
683 {
684  if (_max_substep_iter == 1) // No substepping
685  {
686  _gss[_qp] = _gss_tmp;
687  _acc_slip[_qp] = _accslip_tmp;
688  }
689  else
690  {
691  if (_last_step_iter)
692  {
693  _gss[_qp] = _gss_tmp;
694  _acc_slip[_qp] = _accslip_tmp;
695  }
696  else
697  {
700  }
701  }
702 }
703 
704 void
706 {
707  if (_max_substep_iter == 1) // No substepping
708  {
709  _pk2_tmp = _pk2_old[_qp];
710  _fp_old_inv = _fp_old[_qp].inverse();
713  }
714  else
715  {
716  if (_first_step_iter)
717  {
718  _pk2_tmp = _pk2_tmp_old = _pk2_old[_qp];
719  _fp_old_inv = _fp_old[_qp].inverse();
720  }
721  else
723 
726  }
727 }
728 
729 void
731 {
732  unsigned int iter = 0;
733  RankTwoTensor resid, dpk2;
734  RankFourTensor jac;
735  Real rnorm, rnorm0, rnorm_prev;
736 
737  calc_resid_jacob(resid, jac); // Calculate stress residual
738  if (_err_tol)
739  {
740 #ifdef DEBUG
741  mooseWarning(
742  "FiniteStrainCrystalPLasticity: Slip increment exceeds tolerance - Element number ",
743  _current_elem->id(),
744  " Gauss point = ",
745  _qp);
746 #endif
747  return;
748  }
749 
750  rnorm = resid.L2norm();
751  rnorm0 = rnorm;
752 
753  while (rnorm > _rtol * rnorm0 && rnorm0 > _abs_tol &&
754  iter < _maxiter) // Check for stress residual tolerance
755  {
756  dpk2 = -jac.invSymm() * resid; // Calculate stress increment
757  _pk2_tmp = _pk2_tmp + dpk2; // Update stress
758  calc_resid_jacob(resid, jac);
759  internalVariableUpdateNRiteration(); // update _fp_prev_inv
760 
761  if (_err_tol)
762  {
763 #ifdef DEBUG
764  mooseWarning(
765  "FiniteStrainCrystalPLasticity: Slip increment exceeds tolerance - Element number ",
766  _current_elem->id(),
767  " Gauss point = ",
768  _qp);
769 #endif
770  return;
771  }
772 
773  rnorm_prev = rnorm;
774  rnorm = resid.L2norm();
775 
776  if (_use_line_search && rnorm > rnorm_prev && !line_search_update(rnorm_prev, dpk2))
777  {
778 #ifdef DEBUG
779  mooseWarning("FiniteStrainCrystalPLasticity: Failed with line search");
780 #endif
781  _err_tol = true;
782  return;
783  }
784 
785  if (_use_line_search)
786  rnorm = resid.L2norm();
787 
788  iter++;
789  }
790 
791  if (iter >= _maxiter)
792  {
793 #ifdef DEBUG
794  mooseWarning("FiniteStrainCrystalPLasticity: Stress Integration error rmax = ", rnorm);
795 #endif
796  _err_tol = true;
797  }
798 }
799 
800 void
802 {
803  if (_max_substep_iter == 1) // No substepping
804  {
805  _fp[_qp] = _fp_inv.inverse();
806  _pk2[_qp] = _pk2_tmp;
807  }
808  else
809  {
810  if (_last_step_iter)
811  {
812  _fp[_qp] = _fp_inv.inverse();
813  _pk2[_qp] = _pk2_tmp;
814  }
815  else
816  {
819  }
820  }
821 }
822 
823 // Update slip system resistance. Overide to incorporate new slip system resistance laws
824 void
826 {
827  updateGss();
828 }
829 
834 void
836 {
837  DenseVector<Real> hb(_nss);
838  Real qab;
839 
840  Real a = _hprops[4]; // Kalidindi
841 
843  for (unsigned int i = 0; i < _nss; ++i)
844  _accslip_tmp += std::abs(_slip_incr(i));
845 
846  // Real val = std::cosh(_h0 * _accslip_tmp / (_tau_sat - _tau_init)); // Karthik
847  // val = _h0 * std::pow(1.0/val,2.0); // Kalidindi
848 
849  for (unsigned int i = 0; i < _nss; ++i)
850  // hb(i)=val;
851  hb(i) = _h0 * std::pow(std::abs(1.0 - _gss_tmp[i] / _tau_sat), a) *
852  std::copysign(1.0, 1.0 - _gss_tmp[i] / _tau_sat);
853 
854  for (unsigned int i = 0; i < _nss; ++i)
855  {
856  if (_max_substep_iter == 1) // No substepping
857  _gss_tmp[i] = _gss_old[_qp][i];
858  else
859  _gss_tmp[i] = _gss_tmp_old[i];
860 
861  for (unsigned int j = 0; j < _nss; ++j)
862  {
863  unsigned int iplane, jplane;
864  iplane = i / 3;
865  jplane = j / 3;
866 
867  if (iplane == jplane) // Kalidindi
868  qab = 1.0;
869  else
870  qab = _r;
871 
872  _gss_tmp[i] += qab * hb(j) * std::abs(_slip_incr(j));
873  _dgss_dsliprate(i, j) = qab * hb(j) * std::copysign(1.0, _slip_incr(j)) * _dt;
874  }
875  }
876 }
877 
878 // Calculates stress residual equation and jacobian
879 void
881 {
882  calcResidual(resid);
883  if (_err_tol)
884  return;
885  calcJacobian(jac);
886 }
887 
888 void
890 {
891  RankTwoTensor iden(RankTwoTensor::initIdentity), ce, ee, ce_pk2, eqv_slip_incr, pk2_new;
892 
893  _fe = _dfgrd_tmp * _fp_prev_inv; // _fp_inv ==> _fp_prev_inv
894 
895  ce = _fe.transpose() * _fe;
896  ce_pk2 = ce * _pk2_tmp;
897  ce_pk2 = ce_pk2 / _fe.det();
898 
899  // Calculate Schmid tensor and resolved shear stresses
900  for (unsigned int i = 0; i < _nss; ++i)
901  _tau(i) = ce_pk2.doubleContraction(_s0[i]);
902 
903  getSlipIncrements(); // Calculate dslip,dslipdtau
904 
905  if (_err_tol)
906  return;
907 
908  eqv_slip_incr.zero();
909  for (unsigned int i = 0; i < _nss; ++i)
910  eqv_slip_incr += _s0[i] * _slip_incr(i);
911 
912  eqv_slip_incr = iden - eqv_slip_incr;
913  _fp_inv = _fp_old_inv * eqv_slip_incr;
914  _fe = _dfgrd_tmp * _fp_inv;
915 
916  ce = _fe.transpose() * _fe;
917  ee = ce - iden;
918  ee *= 0.5;
919 
920  pk2_new = _elasticity_tensor[_qp] * ee;
921 
922  resid = _pk2_tmp - pk2_new;
923 }
924 
925 void
927 {
928  RankFourTensor dfedfpinv, deedfe, dfpinvdpk2;
929 
930  std::vector<RankTwoTensor> dtaudpk2(_nss), dfpinvdslip(_nss);
931 
932  for (unsigned int i = 0; i < _nss; ++i)
933  {
934  dtaudpk2[i] = _s0[i];
935  dfpinvdslip[i] = -_fp_old_inv * _s0[i];
936  }
937 
938  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
939  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
940  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
941  dfedfpinv(i, j, k, j) = _dfgrd_tmp(i, k);
942 
943  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
944  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
945  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
946  {
947  deedfe(i, j, k, i) = deedfe(i, j, k, i) + _fe(k, j) * 0.5;
948  deedfe(i, j, k, j) = deedfe(i, j, k, j) + _fe(k, i) * 0.5;
949  }
950 
951  for (unsigned int i = 0; i < _nss; ++i)
952  dfpinvdpk2 += (dfpinvdslip[i] * _dslipdtau(i)).outerProduct(dtaudpk2[i]);
953 
954  jac =
955  RankFourTensor::IdentityFour() - (_elasticity_tensor[_qp] * deedfe * dfedfpinv * dfpinvdpk2);
956 }
957 
958 // Calculate slip increment,dslipdtau. Override to modify.
959 void
961 {
962  for (unsigned int i = 0; i < _nss; ++i)
963  {
964  _slip_incr(i) = _a0(i) * std::pow(std::abs(_tau(i) / _gss_tmp[i]), 1.0 / _xm(i)) *
965  std::copysign(1.0, _tau(i)) * _dt;
966  if (std::abs(_slip_incr(i)) > _slip_incr_tol)
967  {
968  _err_tol = true;
969 #ifdef DEBUG
970  mooseWarning("Maximum allowable slip increment exceeded ", std::abs(_slip_incr(i)));
971 #endif
972  return;
973  }
974  }
975 
976  for (unsigned int i = 0; i < _nss; ++i)
977  _dslipdtau(i) = _a0(i) / _xm(i) *
978  std::pow(std::abs(_tau(i) / _gss_tmp[i]), 1.0 / _xm(i) - 1.0) / _gss_tmp[i] *
979  _dt;
980 }
981 
982 // Calls getMatRot to perform RU factorization of a tensor.
985 {
986  return getMatRot(a);
987 }
988 
989 // Performs RU factorization of a tensor
992 {
993  RankTwoTensor rot;
994  RankTwoTensor c, diag, evec;
995  PetscScalar cmat[LIBMESH_DIM][LIBMESH_DIM], work[10];
996  PetscReal w[LIBMESH_DIM];
997  PetscBLASInt nd = LIBMESH_DIM, lwork = 10, info;
998 
999  c = a.transpose() * a;
1000 
1001  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1002  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1003  cmat[i][j] = c(i, j);
1004 
1005  LAPACKsyev_("V", "U", &nd, &cmat[0][0], &nd, w, work, &lwork, &info);
1006 
1007  if (info != 0)
1008  mooseError("FiniteStrainCrystalPLasticity: DSYEV function call in getMatRot function failed");
1009 
1010  diag.zero();
1011 
1012  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1013  diag(i, i) = std::sqrt(w[i]);
1014 
1015  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1016  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1017  evec(i, j) = cmat[i][j];
1018 
1019  rot = a * ((evec.transpose() * diag * evec).inverse());
1020 
1021  return rot;
1022 }
1023 
1024 // Calculates tangent moduli which is used for global solve
1025 void
1027 {
1028 }
1029 
1032 {
1033  RankFourTensor tan_mod;
1034 
1035  switch (_tan_mod_type)
1036  {
1037  case 0:
1038  tan_mod = elastoPlasticTangentModuli();
1039  break;
1040  default:
1041  tan_mod = elasticTangentModuli();
1042  }
1043 
1044  return tan_mod;
1045 }
1046 
1047 void
1049 {
1050  DenseVector<Real> mo(LIBMESH_DIM * _nss), no(LIBMESH_DIM * _nss);
1051 
1052  // Update slip direction and normal with crystal orientation
1053  for (unsigned int i = 0; i < _nss; ++i)
1054  {
1055  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1056  {
1057  mo(i * LIBMESH_DIM + j) = 0.0;
1058  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1059  mo(i * LIBMESH_DIM + j) =
1060  mo(i * LIBMESH_DIM + j) + _crysrot[_qp](j, k) * _mo(i * LIBMESH_DIM + k);
1061  }
1062 
1063  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1064  {
1065  no(i * LIBMESH_DIM + j) = 0.0;
1066  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1067  no(i * LIBMESH_DIM + j) =
1068  no(i * LIBMESH_DIM + j) + _crysrot[_qp](j, k) * _no(i * LIBMESH_DIM + k);
1069  }
1070  }
1071 
1072  // Calculate Schmid tensor and resolved shear stresses
1073  for (unsigned int i = 0; i < _nss; ++i)
1074  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1075  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1076  _s0[i](j, k) = mo(i * LIBMESH_DIM + j) * no(i * LIBMESH_DIM + k);
1077 }
1078 
1081 {
1082  RankFourTensor tan_mod;
1083  RankTwoTensor pk2fet, fepk2;
1084  RankFourTensor deedfe, dsigdpk2dfe;
1085 
1086  // Fill in the matrix stiffness material property
1087 
1088  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1089  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1090  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
1091  {
1092  deedfe(i, j, k, i) = deedfe(i, j, k, i) + _fe(k, j) * 0.5;
1093  deedfe(i, j, k, j) = deedfe(i, j, k, j) + _fe(k, i) * 0.5;
1094  }
1095 
1096  dsigdpk2dfe = _fe.mixedProductIkJl(_fe) * _elasticity_tensor[_qp] * deedfe;
1097 
1098  pk2fet = _pk2_tmp * _fe.transpose();
1099  fepk2 = _fe * _pk2_tmp;
1100 
1101  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
1102  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
1103  for (unsigned int l = 0; l < LIBMESH_DIM; ++l)
1104  {
1105  tan_mod(i, j, i, l) = tan_mod(i, j, i, l) + pk2fet(l, j);
1106  tan_mod(i, j, j, l) = tan_mod(i, j, j, l) + fepk2(i, l);
1107  }
1108 
1109  tan_mod += dsigdpk2dfe;
1110 
1111  Real je = _fe.det();
1112  if (je > 0.0)
1113  tan_mod /= je;
1114 
1115  return tan_mod;
1116 }
1117 
1120 {
1121  return _elasticity_tensor[_qp]; // update jacobian_mult
1122 }
1123 
1124 bool
1126 {
1127  if (_lsrch_method == "CUT_HALF")
1128  {
1129  Real rnorm;
1130  RankTwoTensor resid;
1131  Real step = 1.0;
1132 
1133  do
1134  {
1135  _pk2_tmp = _pk2_tmp - step * dpk2;
1136  step /= 2.0;
1137  _pk2_tmp = _pk2_tmp + step * dpk2;
1138 
1139  calcResidual(resid);
1140  rnorm = resid.L2norm();
1141  } while (rnorm > rnorm_prev && step > _min_lsrch_step);
1142 
1143  if (rnorm > rnorm_prev && step <= _min_lsrch_step)
1144  return false;
1145 
1146  return true;
1147  }
1148  else if (_lsrch_method == "BISECTION")
1149  {
1150  unsigned int count = 0;
1151  Real step_a = 0.0;
1152  Real step_b = 1.0;
1153  Real step = 1.0;
1154  Real s_m = 1000.0;
1155  Real rnorm = 1000.0;
1156 
1157  RankTwoTensor resid;
1158  calcResidual(resid);
1159  Real s_b = resid.doubleContraction(dpk2);
1160  Real rnorm1 = resid.L2norm();
1161  _pk2_tmp = _pk2_tmp - dpk2;
1162  calcResidual(resid);
1163  Real s_a = resid.doubleContraction(dpk2);
1164  Real rnorm0 = resid.L2norm();
1165  _pk2_tmp = _pk2_tmp + dpk2;
1166 
1167  if ((rnorm1 / rnorm0) < _lsrch_tol || s_a * s_b > 0)
1168  {
1169  calcResidual(resid);
1170  return true;
1171  }
1172 
1173  while ((rnorm / rnorm0) > _lsrch_tol && count < _lsrch_max_iter)
1174  {
1175  _pk2_tmp = _pk2_tmp - step * dpk2;
1176  step = 0.5 * (step_b + step_a);
1177  _pk2_tmp = _pk2_tmp + step * dpk2;
1178  calcResidual(resid);
1179  s_m = resid.doubleContraction(dpk2);
1180  rnorm = resid.L2norm();
1181 
1182  if (s_m * s_a < 0.0)
1183  {
1184  step_b = step;
1185  s_b = s_m;
1186  }
1187  if (s_m * s_b < 0.0)
1188  {
1189  step_a = step;
1190  s_a = s_m;
1191  }
1192  count++;
1193  }
1194 
1195  if ((rnorm / rnorm0) < _lsrch_tol && count < _lsrch_max_iter)
1196  return true;
1197 
1198  return false;
1199  }
1200  else
1201  {
1202  mooseError("Line search meothod is not provided.");
1203  return false;
1204  }
1205 }
1206 
1207 void
1209 {
1210  _fp_prev_inv = _fp_inv; // update _fp_prev_inv
1211 }
FiniteStrainCrystalPlasticity::preSolveQp
virtual void preSolveQp()
This function set variables for stress and internal variable solve.
Definition: FiniteStrainCrystalPlasticity.C:559
FiniteStrainCrystalPlasticity::_min_lsrch_step
Real _min_lsrch_step
Minimum line search step size.
Definition: FiniteStrainCrystalPlasticity.h:293
FiniteStrainCrystalPlasticity::_input_rndm_scale_var
bool _input_rndm_scale_var
Input option for scaling variable to generate random stress when convergence fails.
Definition: FiniteStrainCrystalPlasticity.h:278
FiniteStrainCrystalPlasticity::_h0
Real _h0
Definition: FiniteStrainCrystalPlasticity.h:330
FiniteStrainCrystalPlasticity::_err_tol
bool _err_tol
Definition: FiniteStrainCrystalPlasticity.h:351
FiniteStrainCrystalPlasticity::readFileInitSlipSysRes
virtual void readFileInitSlipSysRes()
This function read slip system resistances from file - see test.
Definition: FiniteStrainCrystalPlasticity.C:243
FiniteStrainCrystalPlasticity::_fp_prev_inv
RankTwoTensor _fp_prev_inv
Definition: FiniteStrainCrystalPlasticity.h:336
FiniteStrainCrystalPlasticity::_pk2_tmp
RankTwoTensor _pk2_tmp
Definition: FiniteStrainCrystalPlasticity.h:340
FiniteStrainCrystalPlasticity::_gss_old
const MaterialProperty< std::vector< Real > > & _gss_old
Definition: FiniteStrainCrystalPlasticity.h:311
FiniteStrainCrystalPlasticity::_slip_sys_flow_prop_file_name
std::string _slip_sys_flow_prop_file_name
File should contain values of the flow rate equation parameters.
Definition: FiniteStrainCrystalPlasticity.h:244
FiniteStrainCrystalPlasticity::getSlipIncrements
virtual void getSlipIncrements()
This function updates the slip increments.
Definition: FiniteStrainCrystalPlasticity.C:960
FiniteStrainCrystalPlasticity::calcJacobian
virtual void calcJacobian(RankFourTensor &)
This function calculate jacobian.
Definition: FiniteStrainCrystalPlasticity.C:926
FiniteStrainCrystalPlasticity::get_current_rotation
RankTwoTensor get_current_rotation(const RankTwoTensor &a)
This function perform RU decomposition to obtain the rotation tensor.
Definition: FiniteStrainCrystalPlasticity.C:984
FiniteStrainCrystalPlasticity::_deformation_gradient
const MaterialProperty< RankTwoTensor > & _deformation_gradient
Definition: FiniteStrainCrystalPlasticity.h:316
ComputeStressBase::_stress
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
Definition: ComputeStressBase.h:50
FiniteStrainCrystalPlasticity::update_slip_system_resistance
virtual void update_slip_system_resistance()
This function updates the slip system resistances.
Definition: FiniteStrainCrystalPlasticity.C:825
FiniteStrainCrystalPlasticity::getHardnessParams
virtual void getHardnessParams()
This function assign flow rate parameters from .i file - see test.
Definition: FiniteStrainCrystalPlasticity.C:411
FiniteStrainCrystalPlasticity::_slip_incr_tol
Real _slip_incr_tol
Slip increment tolerance.
Definition: FiniteStrainCrystalPlasticity.h:256
FiniteStrainCrystalPlasticity
FiniteStrainCrystalPlasticity uses the multiplicative decomposition of deformation gradient and solve...
Definition: FiniteStrainCrystalPlasticity.h:26
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
FiniteStrainCrystalPlasticity::_pk2_old
const MaterialProperty< RankTwoTensor > & _pk2_old
Definition: FiniteStrainCrystalPlasticity.h:307
ComputeStressBase::_Jacobian_mult
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition: ComputeStressBase.h:61
FiniteStrainCrystalPlasticity::_rtol
Real _rtol
Stress residual equation relative tolerance.
Definition: FiniteStrainCrystalPlasticity.h:250
FiniteStrainCrystalPlasticity::_rndm_seed
unsigned int _rndm_seed
Seed value.
Definition: FiniteStrainCrystalPlasticity.h:284
FiniteStrainCrystalPlasticity::FiniteStrainCrystalPlasticity
FiniteStrainCrystalPlasticity(const InputParameters &parameters)
Definition: FiniteStrainCrystalPlasticity.C:101
FiniteStrainCrystalPlasticity::_slip_sys_hard_prop_file_name
std::string _slip_sys_hard_prop_file_name
The hardening parameters in this class are read from .i file. The user can override to read from file...
Definition: FiniteStrainCrystalPlasticity.h:247
FiniteStrainCrystalPlasticity::_last_step_iter
bool _last_step_iter
Definition: FiniteStrainCrystalPlasticity.h:358
FiniteStrainCrystalPlasticity::_num_slip_sys_flowrate_props
unsigned int _num_slip_sys_flowrate_props
Number of slip system flow rate parameters.
Definition: FiniteStrainCrystalPlasticity.h:264
FiniteStrainCrystalPlasticity::assignSlipSysRes
virtual void assignSlipSysRes()
This function assign initial values of slip system resistances/internal variables read from getSlipSy...
Definition: FiniteStrainCrystalPlasticity.C:233
FiniteStrainCrystalPlasticity::_fp_inv
RankTwoTensor _fp_inv
Definition: FiniteStrainCrystalPlasticity.h:336
FiniteStrainCrystalPlasticity::computeQpStress
virtual void computeQpStress()
This function updates the stress at a quadrature point.
Definition: FiniteStrainCrystalPlasticity.C:492
FiniteStrainCrystalPlasticity::_use_line_search
bool _use_line_search
Flag to activate line serach.
Definition: FiniteStrainCrystalPlasticity.h:290
FiniteStrainCrystalPlasticity::_pk2_tmp_old
RankTwoTensor _pk2_tmp_old
Definition: FiniteStrainCrystalPlasticity.h:340
FiniteStrainCrystalPlasticity::_pk2
MaterialProperty< RankTwoTensor > & _pk2
Definition: FiniteStrainCrystalPlasticity.h:306
FiniteStrainCrystalPlasticity::initQpStatefulProperties
virtual void initQpStatefulProperties()
This function initializes the stateful properties such as stress, plastic deformation gradient,...
Definition: FiniteStrainCrystalPlasticity.C:190
FiniteStrainCrystalPlasticity::readFileHardnessParams
virtual void readFileHardnessParams()
This function read hardness parameters from file.
Definition: FiniteStrainCrystalPlasticity.C:405
FiniteStrainCrystalPlasticity::solveStress
virtual void solveStress()
This function solves for stress, updates plastic deformation gradient.
Definition: FiniteStrainCrystalPlasticity.C:730
FiniteStrainCrystalPlasticity::elastoPlasticTangentModuli
virtual RankFourTensor elastoPlasticTangentModuli()
This function calculate the exact tangent moduli for preconditioner.
Definition: FiniteStrainCrystalPlasticity.C:1080
FiniteStrainCrystalPlasticity::_elasticity_tensor
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
Definition: FiniteStrainCrystalPlasticity.h:321
FiniteStrainCrystalPlasticity::elasticTangentModuli
virtual RankFourTensor elasticTangentModuli()
This function calculate the elastic tangent moduli for preconditioner.
Definition: FiniteStrainCrystalPlasticity.C:1119
registerMooseObject
registerMooseObject("TensorMechanicsApp", FiniteStrainCrystalPlasticity)
defineLegacyParams
defineLegacyParams(FiniteStrainCrystalPlasticity)
FiniteStrainCrystalPlasticity::_fe
RankTwoTensor _fe
Definition: FiniteStrainCrystalPlasticity.h:336
FiniteStrainCrystalPlasticity::calcResidual
virtual void calcResidual(RankTwoTensor &)
This function calculate stress residual.
Definition: FiniteStrainCrystalPlasticity.C:889
ComputeStressBase
ComputeStressBase is the base class for stress tensors.
Definition: ComputeStressBase.h:26
FiniteStrainCrystalPlasticity::updateGss
virtual void updateGss()
This function updates the slip system resistances.
Definition: FiniteStrainCrystalPlasticity.C:835
FiniteStrainCrystalPlasticity::_delta_dfgrd
RankTwoTensor _delta_dfgrd
Flag to check whether convergence is achieved.
Definition: FiniteStrainCrystalPlasticity.h:354
FiniteStrainCrystalPlasticity::line_search_update
bool line_search_update(const Real rnorm_prev, const RankTwoTensor)
This function performs the line search update.
Definition: FiniteStrainCrystalPlasticity.C:1125
FiniteStrainCrystalPlasticity::solveStatevar
virtual void solveStatevar()
This function solves internal variables.
Definition: FiniteStrainCrystalPlasticity.C:640
ComputeStressBase::validParams
static InputParameters validParams()
Definition: ComputeStressBase.C:17
FiniteStrainCrystalPlasticity::_flowprops
std::vector< Real > _flowprops
Definition: FiniteStrainCrystalPlasticity.h:232
FiniteStrainCrystalPlasticity.h
FiniteStrainCrystalPlasticity::_crysrot
const MaterialProperty< RankTwoTensor > & _crysrot
Definition: FiniteStrainCrystalPlasticity.h:322
FiniteStrainCrystalPlasticity::_dfgrd_scale_factor
Real _dfgrd_scale_factor
Scales the substepping increment to obtain deformation gradient at a substep iteration.
Definition: FiniteStrainCrystalPlasticity.h:356
FiniteStrainCrystalPlasticity::postSolveStatevar
virtual void postSolveStatevar()
This function update internal variable after solve.
Definition: FiniteStrainCrystalPlasticity.C:682
FiniteStrainCrystalPlasticity::validParams
static InputParameters validParams()
Definition: FiniteStrainCrystalPlasticity.C:22
FiniteStrainCrystalPlasticity::_acc_slip
MaterialProperty< Real > & _acc_slip
Definition: FiniteStrainCrystalPlasticity.h:312
FiniteStrainCrystalPlasticity::_fp_old
const MaterialProperty< RankTwoTensor > & _fp_old
Definition: FiniteStrainCrystalPlasticity.h:305
FiniteStrainCrystalPlasticity::_xm
DenseVector< Real > _xm
Definition: FiniteStrainCrystalPlasticity.h:328
FiniteStrainCrystalPlasticity::getMatRot
RankTwoTensor getMatRot(const RankTwoTensor &a)
This function perform RU decomposition to obtain the rotation tensor.
Definition: FiniteStrainCrystalPlasticity.C:991
FiniteStrainCrystalPlasticity::_mo
DenseVector< Real > _mo
Definition: FiniteStrainCrystalPlasticity.h:324
FiniteStrainCrystalPlasticity::_tau_sat
Real _tau_sat
Definition: FiniteStrainCrystalPlasticity.h:331
FiniteStrainCrystalPlasticity::postSolveQp
virtual void postSolveQp()
This function update stress and internal variable after solve.
Definition: FiniteStrainCrystalPlasticity.C:587
FiniteStrainCrystalPlasticity::_first_substep
bool _first_substep
Definition: FiniteStrainCrystalPlasticity.h:358
FiniteStrainCrystalPlasticity::calc_schmid_tensor
void calc_schmid_tensor()
This function calculate the Schmid tensor.
Definition: FiniteStrainCrystalPlasticity.C:1048
FiniteStrainCrystalPlasticity::_slip_sys_props
DenseVector< Real > _slip_sys_props
Definition: FiniteStrainCrystalPlasticity.h:345
FiniteStrainCrystalPlasticity::_maxiter
unsigned int _maxiter
Maximum number of iterations for stress update.
Definition: FiniteStrainCrystalPlasticity.h:259
FiniteStrainCrystalPlasticity::_fp
MaterialProperty< RankTwoTensor > & _fp
Definition: FiniteStrainCrystalPlasticity.h:304
FiniteStrainCrystalPlasticity::_nss
const unsigned int _nss
Number of slip system resistance.
Definition: FiniteStrainCrystalPlasticity.h:228
FiniteStrainCrystalPlasticity::_a0
DenseVector< Real > _a0
Definition: FiniteStrainCrystalPlasticity.h:327
FiniteStrainCrystalPlasticity::_slip_sys_res_prop_file_name
std::string _slip_sys_res_prop_file_name
File should contain initial values of the slip system resistances.
Definition: FiniteStrainCrystalPlasticity.h:238
FiniteStrainCrystalPlasticity::preSolveStress
virtual void preSolveStress()
This function set variables for stress solve.
Definition: FiniteStrainCrystalPlasticity.C:705
FiniteStrainCrystalPlasticity::_gen_rndm_stress_flag
bool _gen_rndm_stress_flag
Definition: FiniteStrainCrystalPlasticity.h:275
FiniteStrainCrystalPlasticity::_lsrch_max_iter
unsigned int _lsrch_max_iter
Line search bisection method maximum iteration number.
Definition: FiniteStrainCrystalPlasticity.h:299
FiniteStrainCrystalPlasticity::_s0
std::vector< RankTwoTensor > _s0
Definition: FiniteStrainCrystalPlasticity.h:338
FiniteStrainCrystalPlasticity::_gss
MaterialProperty< std::vector< Real > > & _gss
Definition: FiniteStrainCrystalPlasticity.h:310
FiniteStrainCrystalPlasticity::_accslip_tmp
Real _accslip_tmp
Definition: FiniteStrainCrystalPlasticity.h:341
FiniteStrainCrystalPlasticity::_gtol
Real _gtol
Internal variable update equation tolerance.
Definition: FiniteStrainCrystalPlasticity.h:254
FiniteStrainCrystalPlasticity::_no
DenseVector< Real > _no
Definition: FiniteStrainCrystalPlasticity.h:325
FiniteStrainCrystalPlasticity::_dslipdtau
DenseVector< Real > _dslipdtau
Definition: FiniteStrainCrystalPlasticity.h:337
FiniteStrainCrystalPlasticity::postSolveStress
virtual void postSolveStress()
This function update stress and plastic deformation gradient after solve.
Definition: FiniteStrainCrystalPlasticity.C:801
FiniteStrainCrystalPlasticity::calcTangentModuli
virtual RankFourTensor calcTangentModuli()
This function calculate the tangent moduli for preconditioner.
Definition: FiniteStrainCrystalPlasticity.C:1031
FiniteStrainCrystalPlasticity::initAdditionalProps
virtual void initAdditionalProps()
This function initializes additional parameters.
Definition: FiniteStrainCrystalPlasticity.C:483
FiniteStrainCrystalPlasticity::getInitSlipSysRes
virtual void getInitSlipSysRes()
This function assign slip system resistances - see test.
Definition: FiniteStrainCrystalPlasticity.C:261
FiniteStrainCrystalPlasticity::_deformation_gradient_old
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
Definition: FiniteStrainCrystalPlasticity.h:317
FiniteStrainCrystalPlasticity::_lsrch_tol
Real _lsrch_tol
Line search bisection method tolerance.
Definition: FiniteStrainCrystalPlasticity.h:296
FiniteStrainCrystalPlasticity::preSolveStatevar
virtual void preSolveStatevar()
This function set variables for internal variable solve.
Definition: FiniteStrainCrystalPlasticity.C:620
FiniteStrainCrystalPlasticity::_update_rot
MaterialProperty< RankTwoTensor > & _update_rot
Definition: FiniteStrainCrystalPlasticity.h:314
FiniteStrainCrystalPlasticity::getFlowRateParams
virtual void getFlowRateParams()
This function assign flow rate parameters - see test.
Definition: FiniteStrainCrystalPlasticity.C:343
FiniteStrainCrystalPlasticity::_r
Real _r
Definition: FiniteStrainCrystalPlasticity.h:333
FiniteStrainCrystalPlasticity::_dgss_dsliprate
DenseMatrix< Real > _dgss_dsliprate
Definition: FiniteStrainCrystalPlasticity.h:347
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
FiniteStrainCrystalPlasticity::_maxiterg
unsigned int _maxiterg
Maximum number of iterations for internal variable update.
Definition: FiniteStrainCrystalPlasticity.h:261
FiniteStrainCrystalPlasticity::_gss_tmp
std::vector< Real > _gss_tmp
Definition: FiniteStrainCrystalPlasticity.h:342
FiniteStrainCrystalPlasticity::readFileFlowRateParams
virtual void readFileFlowRateParams()
This function read flow rate parameters from file - see test.
Definition: FiniteStrainCrystalPlasticity.C:314
FiniteStrainCrystalPlasticity::_dfgrd_tmp
RankTwoTensor _dfgrd_tmp
Definition: FiniteStrainCrystalPlasticity.h:335
FiniteStrainCrystalPlasticity::initSlipSysProps
virtual void initSlipSysProps()
This function initializes slip system resistances.
Definition: FiniteStrainCrystalPlasticity.C:207
FiniteStrainCrystalPlasticity::solveQp
virtual void solveQp()
This function solves stress and internal variables.
Definition: FiniteStrainCrystalPlasticity.C:577
FiniteStrainCrystalPlasticity::_lsrch_method
MooseEnum _lsrch_method
Definition: FiniteStrainCrystalPlasticity.h:302
FiniteStrainCrystalPlasticity::_slip_incr
DenseVector< Real > _slip_incr
Definition: FiniteStrainCrystalPlasticity.h:337
FiniteStrainCrystalPlasticity::_intvar_read_type
MooseEnum _intvar_read_type
Read from options for initial values of internal variables.
Definition: FiniteStrainCrystalPlasticity.h:270
FiniteStrainCrystalPlasticity::_tan_mod_type
MooseEnum _tan_mod_type
Type of tangent moduli calculation.
Definition: FiniteStrainCrystalPlasticity.h:267
FiniteStrainCrystalPlasticity::calc_resid_jacob
virtual void calc_resid_jacob(RankTwoTensor &, RankFourTensor &)
This function calls the residual and jacobian functions used in the stress update algorithm.
Definition: FiniteStrainCrystalPlasticity.C:880
FiniteStrainCrystalPlasticity::_tau_init
Real _tau_init
Definition: FiniteStrainCrystalPlasticity.h:332
FiniteStrainCrystalPlasticity::_dfgrd_tmp_old
RankTwoTensor _dfgrd_tmp_old
Definition: FiniteStrainCrystalPlasticity.h:354
FiniteStrainCrystalPlasticity::getSlipSystems
virtual void getSlipSystems()
This function reads slip system from file - see test.
Definition: FiniteStrainCrystalPlasticity.C:425
FiniteStrainCrystalPlasticity::_gprops
std::vector< Real > _gprops
Definition: FiniteStrainCrystalPlasticity.h:230
FiniteStrainCrystalPlasticity::internalVariableUpdateNRiteration
void internalVariableUpdateNRiteration()
This function updates internal variables after each NewTon Raphson iteration (_fp_inv)
Definition: FiniteStrainCrystalPlasticity.C:1208
FiniteStrainCrystalPlasticity::_num_slip_sys_props
unsigned int _num_slip_sys_props
Number of slip system specific properties provided in the file containing slip system normals and dir...
Definition: FiniteStrainCrystalPlasticity.h:273
RankTwoTensorTempl< Real >
FiniteStrainCrystalPlasticity::_accslip_tmp_old
Real _accslip_tmp_old
Definition: FiniteStrainCrystalPlasticity.h:341
FiniteStrainCrystalPlasticity::_lag_e
MaterialProperty< RankTwoTensor > & _lag_e
Definition: FiniteStrainCrystalPlasticity.h:308
FiniteStrainCrystalPlasticity::_rndm_scale_var
Real _rndm_scale_var
Scaling value.
Definition: FiniteStrainCrystalPlasticity.h:281
FiniteStrainCrystalPlasticity::_read_from_slip_sys_file
bool _read_from_slip_sys_file
Definition: FiniteStrainCrystalPlasticity.h:349
FiniteStrainCrystalPlasticity::_first_step_iter
bool _first_step_iter
Flags to reset variables and reinitialize variables.
Definition: FiniteStrainCrystalPlasticity.h:358
FiniteStrainCrystalPlasticity::_gss_tmp_old
std::vector< Real > _gss_tmp_old
Definition: FiniteStrainCrystalPlasticity.h:343
FiniteStrainCrystalPlasticity::_slip_sys_file_name
std::string _slip_sys_file_name
File should contain slip plane normal and direction. See test.
Definition: FiniteStrainCrystalPlasticity.h:235
FiniteStrainCrystalPlasticity::computeQpElasticityTensor
virtual void computeQpElasticityTensor()
This function updates the elasticity tensor at a quadrature point.
Definition: FiniteStrainCrystalPlasticity.C:1026
FiniteStrainCrystalPlasticity::_abs_tol
Real _abs_tol
Stress residual equation absolute tolerance.
Definition: FiniteStrainCrystalPlasticity.h:252
FiniteStrainCrystalPlasticity::_fp_old_inv
RankTwoTensor _fp_old_inv
Definition: FiniteStrainCrystalPlasticity.h:336
FiniteStrainCrystalPlasticity::_hprops
std::vector< Real > _hprops
Definition: FiniteStrainCrystalPlasticity.h:231
FiniteStrainCrystalPlasticity::_tau
DenseVector< Real > _tau
Definition: FiniteStrainCrystalPlasticity.h:337
FiniteStrainCrystalPlasticity::_max_substep_iter
unsigned int _max_substep_iter
Maximum number of substep iterations.
Definition: FiniteStrainCrystalPlasticity.h:287
FiniteStrainCrystalPlasticity::_acc_slip_old
const MaterialProperty< Real > & _acc_slip_old
Definition: FiniteStrainCrystalPlasticity.h:313