www.mooseframework.org
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
MultiPlasticityRawComponentAssembler Class Reference

MultiPlasticityRawComponentAssembler holds and computes yield functions, flow directions, etc, for use in FiniteStrainMultiPlasticity. More...

#include <MultiPlasticityRawComponentAssembler.h>

Inheritance diagram for MultiPlasticityRawComponentAssembler:
[legend]

Public Member Functions

 MultiPlasticityRawComponentAssembler (const MooseObject *moose_object)
 
virtual ~MultiPlasticityRawComponentAssembler ()
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual void yieldFunction (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< Real > &f)
 The active yield function(s) More...
 
virtual void dyieldFunction_dstress (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankTwoTensor > &df_dstress)
 The derivative of the active yield function(s) with respect to stress. More...
 
virtual void dyieldFunction_dintnl (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< Real > &df_dintnl)
 The derivative of active yield function(s) with respect to their internal parameters (the user objects assume there is exactly one internal param per yield function) More...
 
virtual void flowPotential (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankTwoTensor > &r)
 The active flow potential(s) - one for each yield function. More...
 
virtual void dflowPotential_dstress (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankFourTensor > &dr_dstress)
 The derivative of the active flow potential(s) with respect to stress. More...
 
virtual void dflowPotential_dintnl (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankTwoTensor > &dr_dintnl)
 The derivative of the active flow potentials with respect to the active internal parameters The UserObjects explicitly assume that r[alpha] is only dependent on intnl[alpha]. More...
 
virtual void hardPotential (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< Real > &h)
 The active hardening potentials (one for each internal parameter and for each yield function) by assumption in the Userobjects, the h[a][alpha] is nonzero only if the surface alpha is part of model a, so we only calculate those here. More...
 
virtual void dhardPotential_dstress (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankTwoTensor > &dh_dstress)
 The derivative of the active hardening potentials with respect to stress By assumption in the Userobjects, the h[a][alpha] is nonzero only for a = alpha, so we only calculate those here. More...
 
virtual void dhardPotential_dintnl (const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< Real > &dh_dintnl)
 The derivative of the active hardening potentials with respect to the active internal parameters. More...
 
virtual void buildActiveConstraints (const std::vector< Real > &f, const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &Eijkl, std::vector< bool > &act)
 Constructs a set of active constraints, given the yield functions, f. More...
 
unsigned int modelNumber (unsigned int surface)
 returns the model number, given the surface number More...
 
bool anyActiveSurfaces (int model, const std::vector< bool > &active)
 returns true if any internal surfaces of the given model are active according to 'active' More...
 
void activeModelSurfaces (int model, const std::vector< bool > &active, std::vector< unsigned int > &active_surfaces_of_model)
 Returns the internal surface number(s) of the active surfaces of the given model This may be of size=0 if there are no active surfaces of the given model. More...
 
void activeSurfaces (int model, const std::vector< bool > &active, std::vector< unsigned int > &active_surfaces)
 Returns the external surface number(s) of the active surfaces of the given model This may be of size=0 if there are no active surfaces of the given model. More...
 
bool returnMapAll (const RankTwoTensor &trial_stress, const std::vector< Real > &intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &stress, std::vector< Real > &intnl, std::vector< Real > &pm, std::vector< Real > &cumulative_pm, RankTwoTensor &delta_dp, std::vector< Real > &yf, unsigned &num_successful_plastic_returns, unsigned &custom_model)
 Performs a returnMap for each plastic model using their inbuilt returnMap functions. More...
 

Protected Attributes

const InputParameters & _params
 
unsigned int _num_models
 Number of plastic models for this material. More...
 
unsigned int _num_surfaces
 Number of surfaces within the plastic models. More...
 
std::vector< std::vector< unsigned int > > _surfaces_given_model
 _surfaces_given_model[model_number] = vector of surface numbers for this model More...
 
MooseEnum _specialIC
 Allows initial set of active constraints to be chosen optimally. More...
 
std::vector< const TensorMechanicsPlasticModel * > _f
 User objects that define the yield functions, flow potentials, etc. More...
 

Private Member Functions

void buildActiveConstraintsRock (const std::vector< Real > &f, const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &Eijkl, std::vector< bool > &act)
 "Rock" version Constructs a set of active constraints, given the yield functions, f. More...
 
void buildActiveConstraintsJoint (const std::vector< Real > &f, const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &Eijkl, std::vector< bool > &act)
 "Joint" version Constructs a set of active constraints, given the yield functions, f. More...
 

Private Attributes

std::vector< unsigned int > _model_given_surface
 given a surface number, this returns the model number More...
 
std::vector< unsigned int > _model_surface_given_surface
 given a surface number, this returns the corresponding-model's internal surface number More...
 

Detailed Description

MultiPlasticityRawComponentAssembler holds and computes yield functions, flow directions, etc, for use in FiniteStrainMultiPlasticity.

Here there are a number of numbering systems.

There are _num_models of plastic models, read from the input file. Each of these models can, in principal, contain multiple 'internal surfaces'. Models are numbered from 0 to _num_models - 1.

There are num_surfaces surfaces. This = sum{plastic_models} (numberSurfaces in model) Evidently _num_surface >= _num_models Surfaces are numbered from 0 to _num_surfaces - 1.

The std::vectors _model_given_surface, _model_surface_given_surface and _surfaces_given_model allow translation between these

Definition at line 44 of file MultiPlasticityRawComponentAssembler.h.

Constructor & Destructor Documentation

◆ MultiPlasticityRawComponentAssembler()

MultiPlasticityRawComponentAssembler::MultiPlasticityRawComponentAssembler ( const MooseObject *  moose_object)

Definition at line 38 of file MultiPlasticityRawComponentAssembler.C.

40  : UserObjectInterface(moose_object),
41  _params(moose_object->parameters()),
42  _num_models(_params.get<std::vector<UserObjectName>>("plastic_models").size()),
43  _num_surfaces(0),
44  _specialIC(_params.get<MooseEnum>("specialIC"))
45 {
46  _f.resize(_num_models);
47  for (unsigned model = 0; model < _num_models; ++model)
48  _f[model] = &getUserObjectByName<TensorMechanicsPlasticModel>(
49  _params.get<std::vector<UserObjectName>>("plastic_models")[model]);
50 
51  for (unsigned model = 0; model < _num_models; ++model)
52  _num_surfaces += _f[model]->numberSurfaces();
53 
56  unsigned int surface = 0;
57  for (unsigned model = 0; model < _num_models; ++model)
58  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
59  {
60  _model_given_surface[surface] = model;
61  _model_surface_given_surface[surface] = model_surface;
62  surface++;
63  }
64 
66  surface = 0;
67  for (unsigned model = 0; model < _num_models; ++model)
68  {
69  _surfaces_given_model[model].resize(_f[model]->numberSurfaces());
70  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
71  _surfaces_given_model[model][model_surface] = surface++;
72  }
73 
74  // check the plastic_models for specialIC
75  if (_specialIC == "rock")
76  {
77  if (_num_models != 2)
78  mooseError("Choosing specialIC=rock, you must have plasticity models of type 'TensileMulti "
79  "MohrCoulombMulti'\n");
80  if (!(_f[0]->modelName().compare("TensileMulti") == 0 &&
81  _f[1]->modelName().compare("MohrCoulombMulti") == 0))
82  mooseError("Choosing specialIC=rock, you must have plasticity models of type 'TensileMulti "
83  "MohrCoulombMulti'\n");
84  }
85  if (_specialIC == "joint")
86  {
87  if (_num_models != 2)
88  mooseError("Choosing specialIC=joint, you must have plasticity models of type "
89  "'WeakPlaneTensile WeakPlaneShear'\n");
90  if (!(_f[0]->modelName().compare("WeakPlaneTensile") == 0 &&
91  _f[1]->modelName().compare("WeakPlaneShear") == 0))
92  mooseError("Choosing specialIC=joint, you must have plasticity models of type "
93  "'WeakPlaneTensile WeakPlaneShear'\n");
94  }
95 }

◆ ~MultiPlasticityRawComponentAssembler()

virtual MultiPlasticityRawComponentAssembler::~MultiPlasticityRawComponentAssembler ( )
inlinevirtual

Definition at line 51 of file MultiPlasticityRawComponentAssembler.h.

51 {}

Member Function Documentation

◆ activeModelSurfaces()

void MultiPlasticityRawComponentAssembler::activeModelSurfaces ( int  model,
const std::vector< bool > &  active,
std::vector< unsigned int > &  active_surfaces_of_model 
)
protected

Returns the internal surface number(s) of the active surfaces of the given model This may be of size=0 if there are no active surfaces of the given model.

Parameters
modelthe model number
activearray with entries being 'true' if the surface is active
[out]active_surfaces_of_modelthe output

Definition at line 811 of file MultiPlasticityRawComponentAssembler.C.

815 {
816  active_surfaces_of_model.resize(0);
817  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
818  if (active[_surfaces_given_model[model][model_surface]])
819  active_surfaces_of_model.push_back(model_surface);
820 }

Referenced by dflowPotential_dintnl(), dflowPotential_dstress(), dhardPotential_dintnl(), dhardPotential_dstress(), dyieldFunction_dintnl(), dyieldFunction_dstress(), flowPotential(), hardPotential(), and yieldFunction().

◆ activeSurfaces()

void MultiPlasticityRawComponentAssembler::activeSurfaces ( int  model,
const std::vector< bool > &  active,
std::vector< unsigned int > &  active_surfaces 
)
protected

Returns the external surface number(s) of the active surfaces of the given model This may be of size=0 if there are no active surfaces of the given model.

Parameters
modelthe model number
activearray with entries being 'true' if the surface is active
[out]active_surfacesthe output

Definition at line 800 of file MultiPlasticityRawComponentAssembler.C.

803 {
804  active_surfaces.resize(0);
805  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
806  if (active[_surfaces_given_model[model][model_surface]])
807  active_surfaces.push_back(_surfaces_given_model[model][model_surface]);
808 }

Referenced by MultiPlasticityLinearSystem::calculateConstraints().

◆ anyActiveSurfaces()

bool MultiPlasticityRawComponentAssembler::anyActiveSurfaces ( int  model,
const std::vector< bool > &  active 
)
protected

returns true if any internal surfaces of the given model are active according to 'active'

Definition at line 791 of file MultiPlasticityRawComponentAssembler.C.

792 {
793  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
794  if (active[_surfaces_given_model[model][model_surface]])
795  return true;
796  return false;
797 }

Referenced by MultiPlasticityLinearSystem::calculateJacobian(), MultiPlasticityLinearSystem::calculateRHS(), MultiPlasticityDebugger::checkSolution(), MultiPlasticityDebugger::dof_included(), MultiPlasticityLinearSystem::nrStep(), and ComputeMultiPlasticityStress::residual2().

◆ buildActiveConstraints()

void MultiPlasticityRawComponentAssembler::buildActiveConstraints ( const std::vector< Real > &  f,
const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const RankFourTensor Eijkl,
std::vector< bool > &  act 
)
protectedvirtual

Constructs a set of active constraints, given the yield functions, f.

This uses TensorMechanicsPlasticModel::activeConstraints to identify the active constraints for each model.

Parameters
fyield functions (should be _num_surfaces of these)
stressstress tensor
intnlinternal parameters
Eijklelasticity tensor (stress = Eijkl*strain)
[out]actthe set of active constraints (will be resized to _num_surfaces)

Definition at line 344 of file MultiPlasticityRawComponentAssembler.C.

349 {
350  mooseAssert(f.size() == _num_surfaces,
351  "buildActiveConstraints called with f.size = " << f.size() << " while there are "
352  << _num_surfaces << " surfaces");
353  mooseAssert(intnl.size() == _num_models,
354  "buildActiveConstraints called with intnl.size = "
355  << intnl.size() << " while there are " << _num_models << " models");
356 
357  if (_specialIC == "rock")
358  buildActiveConstraintsRock(f, stress, intnl, Eijkl, act);
359  else if (_specialIC == "joint")
360  buildActiveConstraintsJoint(f, stress, intnl, Eijkl, act);
361  else // no specialIC
362  {
363  act.resize(0);
364  unsigned ind = 0;
365  for (unsigned model = 0; model < _num_models; ++model)
366  {
367  std::vector<Real> model_f(0);
368  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
369  model_f.push_back(f[ind++]);
370  std::vector<bool> model_act;
371  RankTwoTensor returned_stress;
372  _f[model]->activeConstraints(
373  model_f, stress, intnl[model], Eijkl, model_act, returned_stress);
374  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
375  act.push_back(model_act[model_surface]);
376  }
377  }
378 }

Referenced by ComputeMultiPlasticityStress::returnMap().

◆ buildActiveConstraintsJoint()

void MultiPlasticityRawComponentAssembler::buildActiveConstraintsJoint ( const std::vector< Real > &  f,
const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const RankFourTensor Eijkl,
std::vector< bool > &  act 
)
private

"Joint" version Constructs a set of active constraints, given the yield functions, f.

This uses TensorMechanicsPlasticModel::activeConstraints to identify the active constraints for each model.

Parameters
fyield functions (should be _num_surfaces of these)
stressstress tensor
intnlinternal parameters
Eijklelasticity tensor (stress = Eijkl*strain)
[out]actthe set of active constraints (will be resized to _num_surfaces)

Definition at line 381 of file MultiPlasticityRawComponentAssembler.C.

386 {
387  act.assign(2, false);
388 
389  RankTwoTensor returned_stress;
390  std::vector<bool> active_tensile;
391  std::vector<bool> active_shear;
392  std::vector<Real> f_single;
393 
394  // first try tensile alone
395  f_single.assign(1, 0);
396  f_single[0] = f[0];
397  _f[0]->activeConstraints(f_single, stress, intnl[0], Eijkl, active_tensile, returned_stress);
398  _f[1]->yieldFunctionV(returned_stress, intnl[1], f_single);
399  if (f_single[0] <= _f[1]->_f_tol)
400  {
401  act[0] = active_tensile[0];
402  return;
403  }
404 
405  // next try shear alone
406  f_single.assign(1, 0);
407  f_single[0] = f[1];
408  _f[1]->activeConstraints(f_single, stress, intnl[1], Eijkl, active_shear, returned_stress);
409  _f[0]->yieldFunctionV(returned_stress, intnl[0], f_single);
410  if (f_single[0] <= _f[0]->_f_tol)
411  {
412  act[1] = active_shear[0];
413  return;
414  }
415 
416  // must be mixed
417  act[0] = act[1] = true;
418  return;
419 }

Referenced by buildActiveConstraints().

◆ buildActiveConstraintsRock()

void MultiPlasticityRawComponentAssembler::buildActiveConstraintsRock ( const std::vector< Real > &  f,
const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const RankFourTensor Eijkl,
std::vector< bool > &  act 
)
private

"Rock" version Constructs a set of active constraints, given the yield functions, f.

This uses TensorMechanicsPlasticModel::activeConstraints to identify the active constraints for each model.

Parameters
fyield functions (should be _num_surfaces of these)
stressstress tensor
intnlinternal parameters
Eijklelasticity tensor (stress = Eijkl*strain)
[out]actthe set of active constraints (will be resized to _num_surfaces)

Definition at line 422 of file MultiPlasticityRawComponentAssembler.C.

427 {
428  act.assign(9, false);
429 
430  RankTwoTensor returned_stress;
431  std::vector<bool> active_tensile;
432  std::vector<bool> active_MC;
433  std::vector<Real> f_single;
434 
435  // first try tensile alone
436  f_single.assign(3, 0);
437  f_single[0] = f[0];
438  f_single[1] = f[1];
439  f_single[2] = f[2];
440  _f[0]->activeConstraints(f_single, stress, intnl[0], Eijkl, active_tensile, returned_stress);
441  _f[1]->yieldFunctionV(returned_stress, intnl[1], f_single);
442  if (f_single[0] <= _f[1]->_f_tol && f_single[1] <= _f[1]->_f_tol &&
443  f_single[2] <= _f[1]->_f_tol && f_single[3] <= _f[1]->_f_tol &&
444  f_single[4] <= _f[1]->_f_tol && f_single[5] <= _f[1]->_f_tol)
445  {
446  act[0] = active_tensile[0];
447  act[1] = active_tensile[1];
448  act[2] = active_tensile[2];
449  return;
450  }
451 
452  // next try MC alone
453  f_single.assign(6, 0);
454  f_single[0] = f[3];
455  f_single[1] = f[4];
456  f_single[2] = f[5];
457  f_single[3] = f[6];
458  f_single[4] = f[7];
459  f_single[5] = f[8];
460  _f[1]->activeConstraints(f_single, stress, intnl[1], Eijkl, active_MC, returned_stress);
461  _f[0]->yieldFunctionV(returned_stress, intnl[0], f_single);
462  if (f_single[0] <= _f[0]->_f_tol && f_single[1] <= _f[0]->_f_tol && f_single[2] <= _f[0]->_f_tol)
463  {
464  act[3] = active_MC[0];
465  act[4] = active_MC[1];
466  act[5] = active_MC[2];
467  act[6] = active_MC[3];
468  act[7] = active_MC[4];
469  act[8] = active_MC[5];
470  return;
471  }
472 
473  // must be a mix.
474  // The possibilities are enumerated below.
475 
476  // tensile=edge, MC=tip (two possibilities)
477  if (active_tensile[0] == false && active_tensile[1] == true && active_tensile[2] == true &&
478  active_MC[0] == true && active_MC[1] == true && active_MC[2] == false &&
479  active_MC[3] == true && active_MC[4] == false && active_MC[5] == false)
480  {
481  act[1] = act[2] = act[6] = true;
482  act[4] = true;
483  return;
484  }
485  if (active_tensile[0] == false && active_tensile[1] == true && active_tensile[2] == true &&
486  active_MC[0] == false && active_MC[1] == true && active_MC[2] == false &&
487  active_MC[3] == true && active_MC[4] == false && active_MC[5] == true)
488  {
489  act[1] = act[2] = act[6] = true; // i don't think act[4] is necessary, is it?!
490  return;
491  }
492 
493  // tensile = edge, MC=edge (two possibilities)
494  if (active_tensile[0] == false && active_tensile[1] == true && active_tensile[2] == true &&
495  active_MC[0] == false && active_MC[1] == true && active_MC[2] == false &&
496  active_MC[3] == true && active_MC[4] == false && active_MC[5] == false)
497  {
498  act[1] = act[2] = act[4] = act[6] = true;
499  return;
500  }
501  if (active_tensile[0] == false && active_tensile[1] == true && active_tensile[2] == true &&
502  active_MC[0] == false && active_MC[1] == false && active_MC[2] == false &&
503  active_MC[3] == true && active_MC[4] == false && active_MC[5] == true)
504  {
505  act[1] = act[2] = act[4] = act[6] = true;
506  return;
507  }
508 
509  // tensile = edge, MC=face
510  if (active_tensile[0] == false && active_tensile[1] == true && active_tensile[2] == true &&
511  active_MC[0] == false && active_MC[1] == false && active_MC[2] == false &&
512  active_MC[3] == true && active_MC[4] == false && active_MC[5] == false)
513  {
514  act[1] = act[2] = act[6] = true;
515  return;
516  }
517 
518  // tensile = face, MC=tip (two possibilities)
519  if (active_tensile[0] == false && active_tensile[1] == false && active_tensile[2] == true &&
520  active_MC[0] == true && active_MC[1] == true && active_MC[2] == false &&
521  active_MC[3] == true && active_MC[4] == false && active_MC[5] == false)
522  {
523  act[2] = act[6] = true;
524  act[4] = true;
525  act[8] = true;
526  return;
527  }
528  if (active_tensile[0] == false && active_tensile[1] == false && active_tensile[2] == true &&
529  active_MC[0] == false && active_MC[1] == true && active_MC[2] == false &&
530  active_MC[3] == true && active_MC[4] == false && active_MC[5] == true)
531  {
532  act[2] = act[6] = true;
533  act[8] = true;
534  return;
535  }
536 
537  // tensile = face, MC=face
538  if (active_tensile[0] == false && active_tensile[1] == false && active_tensile[2] == true &&
539  active_MC[0] == false && active_MC[1] == false && active_MC[2] == false &&
540  active_MC[3] == true && active_MC[4] == false && active_MC[5] == false)
541  {
542  act[1] = act[2] = act[6] = true;
543  return;
544  }
545 
546  // tensile = face, MC=edge (two possibilites).
547  act[2] = true; // tensile face
548  act[3] = active_MC[0];
549  act[4] = active_MC[1];
550  act[5] = active_MC[2];
551  act[6] = active_MC[3];
552  act[7] = active_MC[4];
553  act[8] = active_MC[5];
554  return;
555 }

Referenced by buildActiveConstraints().

◆ dflowPotential_dintnl()

void MultiPlasticityRawComponentAssembler::dflowPotential_dintnl ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< RankTwoTensor > &  dr_dintnl 
)
protectedvirtual

The derivative of the active flow potentials with respect to the active internal parameters The UserObjects explicitly assume that r[alpha] is only dependent on intnl[alpha].

Parameters
stressthe stress at which to calculate the flow potential
intnlvector of internal parameters
activeset of active constraints - only the active derivatives are put into "dr_dintnl"
[out]dr_dintnlthe derivatives. dr_dintnl[alpha](i, j) = dr[alpha](i, j)/dintnl[alpha]

Definition at line 235 of file MultiPlasticityRawComponentAssembler.C.

239 {
240  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
241  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
242 
243  dr_dintnl.resize(0);
244  std::vector<unsigned int> active_surfaces_of_model;
245  std::vector<unsigned int>::iterator active_surface;
246  std::vector<RankTwoTensor> model_dr_dintnl;
247  for (unsigned model = 0; model < _num_models; ++model)
248  {
249  activeModelSurfaces(model, active, active_surfaces_of_model);
250  if (active_surfaces_of_model.size() > 0)
251  {
252  _f[model]->dflowPotential_dintnlV(stress, intnl[model], model_dr_dintnl);
253  for (active_surface = active_surfaces_of_model.begin();
254  active_surface != active_surfaces_of_model.end();
255  ++active_surface)
256  dr_dintnl.push_back(model_dr_dintnl[*active_surface]);
257  }
258  }
259 }

Referenced by MultiPlasticityLinearSystem::calculateJacobian(), MultiPlasticityDebugger::checkDerivatives(), and ComputeMultiPlasticityStress::consistentTangentOperator().

◆ dflowPotential_dstress()

void MultiPlasticityRawComponentAssembler::dflowPotential_dstress ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< RankFourTensor > &  dr_dstress 
)
protectedvirtual

The derivative of the active flow potential(s) with respect to stress.

Parameters
stressthe stress at which to calculate the flow potential
intnlvector of internal parameters
activeset of active constraints - only the active derivatives are put into "dr_dstress"
[out]dr_dstressthe derivative. dr_dstress[alpha](i, j, k, l) = dr[alpha](i, j)/dstress(k, l)

Definition at line 207 of file MultiPlasticityRawComponentAssembler.C.

212 {
213  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
214  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
215 
216  dr_dstress.resize(0);
217  std::vector<unsigned int> active_surfaces_of_model;
218  std::vector<unsigned int>::iterator active_surface;
219  std::vector<RankFourTensor> model_dr_dstress;
220  for (unsigned model = 0; model < _num_models; ++model)
221  {
222  activeModelSurfaces(model, active, active_surfaces_of_model);
223  if (active_surfaces_of_model.size() > 0)
224  {
225  _f[model]->dflowPotential_dstressV(stress, intnl[model], model_dr_dstress);
226  for (active_surface = active_surfaces_of_model.begin();
227  active_surface != active_surfaces_of_model.end();
228  ++active_surface)
229  dr_dstress.push_back(model_dr_dstress[*active_surface]);
230  }
231  }
232 }

Referenced by MultiPlasticityLinearSystem::calculateJacobian(), MultiPlasticityDebugger::checkDerivatives(), and ComputeMultiPlasticityStress::consistentTangentOperator().

◆ dhardPotential_dintnl()

void MultiPlasticityRawComponentAssembler::dhardPotential_dintnl ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< Real > &  dh_dintnl 
)
protectedvirtual

The derivative of the active hardening potentials with respect to the active internal parameters.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlvector of internal parameters
activeset of active constraints - only the active derivatives are put into "dh_dintnl"
[out]dh_dintnlthe derivatives. dh_dintnl[a][alpha][b] = dh[a][alpha]/dintnl[b]. Note that the userobjects assume that there is exactly one internal parameter per yield function, so the derivative is only nonzero for a=alpha=b, so that is all we calculate

Definition at line 317 of file MultiPlasticityRawComponentAssembler.C.

321 {
322  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
323  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
324 
325  dh_dintnl.resize(0);
326  std::vector<unsigned int> active_surfaces_of_model;
327  std::vector<unsigned int>::iterator active_surface;
328  std::vector<Real> model_dh_dintnl;
329  for (unsigned model = 0; model < _num_models; ++model)
330  {
331  activeModelSurfaces(model, active, active_surfaces_of_model);
332  if (active_surfaces_of_model.size() > 0)
333  {
334  _f[model]->dhardPotential_dintnlV(stress, intnl[model], model_dh_dintnl);
335  for (active_surface = active_surfaces_of_model.begin();
336  active_surface != active_surfaces_of_model.end();
337  ++active_surface)
338  dh_dintnl.push_back(model_dh_dintnl[*active_surface]);
339  }
340  }
341 }

Referenced by MultiPlasticityLinearSystem::calculateJacobian().

◆ dhardPotential_dstress()

void MultiPlasticityRawComponentAssembler::dhardPotential_dstress ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< RankTwoTensor > &  dh_dstress 
)
protectedvirtual

The derivative of the active hardening potentials with respect to stress By assumption in the Userobjects, the h[a][alpha] is nonzero only for a = alpha, so we only calculate those here.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlvector of internal parameters
activeset of active constraints - only the active derivatives are put into "dh_dstress"
[out]dh_dstressthe derivative. dh_dstress[a](i, j) = dh[a]/dstress(k, l)

Definition at line 289 of file MultiPlasticityRawComponentAssembler.C.

294 {
295  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
296  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
297 
298  dh_dstress.resize(0);
299  std::vector<unsigned int> active_surfaces_of_model;
300  std::vector<unsigned int>::iterator active_surface;
301  std::vector<RankTwoTensor> model_dh_dstress;
302  for (unsigned model = 0; model < _num_models; ++model)
303  {
304  activeModelSurfaces(model, active, active_surfaces_of_model);
305  if (active_surfaces_of_model.size() > 0)
306  {
307  _f[model]->dhardPotential_dstressV(stress, intnl[model], model_dh_dstress);
308  for (active_surface = active_surfaces_of_model.begin();
309  active_surface != active_surfaces_of_model.end();
310  ++active_surface)
311  dh_dstress.push_back(model_dh_dstress[*active_surface]);
312  }
313  }
314 }

Referenced by MultiPlasticityLinearSystem::calculateJacobian().

◆ dyieldFunction_dintnl()

void MultiPlasticityRawComponentAssembler::dyieldFunction_dintnl ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< Real > &  df_dintnl 
)
protectedvirtual

The derivative of active yield function(s) with respect to their internal parameters (the user objects assume there is exactly one internal param per yield function)

Parameters
stressthe stress at which to calculate the yield function
intnlvector of internal parameters
activeset of active constraints - only the active derivatives are put into "df_dintnl"
[out]df_dintnlthe derivatives. df_dstress[alpha] = dyieldFunction[alpha]/dintnl[alpha]

Definition at line 153 of file MultiPlasticityRawComponentAssembler.C.

157 {
158  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
159  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
160 
161  df_dintnl.resize(0);
162  std::vector<unsigned int> active_surfaces_of_model;
163  std::vector<unsigned int>::iterator active_surface;
164  std::vector<Real> model_df_dintnl;
165  for (unsigned model = 0; model < _num_models; ++model)
166  {
167  activeModelSurfaces(model, active, active_surfaces_of_model);
168  if (active_surfaces_of_model.size() > 0)
169  {
170  _f[model]->dyieldFunction_dintnlV(stress, intnl[model], model_df_dintnl);
171  for (active_surface = active_surfaces_of_model.begin();
172  active_surface != active_surfaces_of_model.end();
173  ++active_surface)
174  df_dintnl.push_back(model_df_dintnl[*active_surface]);
175  }
176  }
177 }

Referenced by MultiPlasticityLinearSystem::calculateJacobian(), MultiPlasticityDebugger::checkDerivatives(), and ComputeMultiPlasticityStress::consistentTangentOperator().

◆ dyieldFunction_dstress()

void MultiPlasticityRawComponentAssembler::dyieldFunction_dstress ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< RankTwoTensor > &  df_dstress 
)
protectedvirtual

The derivative of the active yield function(s) with respect to stress.

Parameters
stressthe stress at which to calculate the yield function
intnlvector of internal parameters
activeset of active constraints - only the active derivatives are put into "df_dstress"
[out]df_dstressthe derivative (or derivatives in the case of multisurface plasticity). df_dstress[alpha](i, j) = dyieldFunction[alpha]/dstress(i, j)

Definition at line 125 of file MultiPlasticityRawComponentAssembler.C.

130 {
131  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
132  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
133 
134  df_dstress.resize(0);
135  std::vector<unsigned int> active_surfaces_of_model;
136  std::vector<unsigned int>::iterator active_surface;
137  std::vector<RankTwoTensor> model_df_dstress;
138  for (unsigned model = 0; model < _num_models; ++model)
139  {
140  activeModelSurfaces(model, active, active_surfaces_of_model);
141  if (active_surfaces_of_model.size() > 0)
142  {
143  _f[model]->dyieldFunction_dstressV(stress, intnl[model], model_df_dstress);
144  for (active_surface = active_surfaces_of_model.begin();
145  active_surface != active_surfaces_of_model.end();
146  ++active_surface)
147  df_dstress.push_back(model_df_dstress[*active_surface]);
148  }
149  }
150 }

Referenced by ComputeMultiPlasticityStress::buildDumbOrder(), MultiPlasticityLinearSystem::calculateJacobian(), MultiPlasticityDebugger::checkDerivatives(), ComputeMultiPlasticityStress::consistentTangentOperator(), and MultiPlasticityLinearSystem::eliminateLinearDependence().

◆ flowPotential()

void MultiPlasticityRawComponentAssembler::flowPotential ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< RankTwoTensor > &  r 
)
protectedvirtual

The active flow potential(s) - one for each yield function.

Parameters
stressthe stress at which to calculate the flow potential
intnlvector of internal parameters
activeset of active constraints - only the active flow potentials are put into "r"
[out]rthe flow potential (flow potentials in the multi-surface case)

Definition at line 180 of file MultiPlasticityRawComponentAssembler.C.

184 {
185  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
186  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
187 
188  r.resize(0);
189  std::vector<unsigned int> active_surfaces_of_model;
190  std::vector<unsigned int>::iterator active_surface;
191  std::vector<RankTwoTensor> model_r;
192  for (unsigned model = 0; model < _num_models; ++model)
193  {
194  activeModelSurfaces(model, active, active_surfaces_of_model);
195  if (active_surfaces_of_model.size() > 0)
196  {
197  _f[model]->flowPotentialV(stress, intnl[model], model_r);
198  for (active_surface = active_surfaces_of_model.begin();
199  active_surface != active_surfaces_of_model.end();
200  ++active_surface)
201  r.push_back(model_r[*active_surface]);
202  }
203  }
204 }

Referenced by MultiPlasticityLinearSystem::calculateConstraints(), MultiPlasticityLinearSystem::calculateJacobian(), ComputeMultiPlasticityStress::consistentTangentOperator(), MultiPlasticityDebugger::fddflowPotential_dintnl(), and MultiPlasticityDebugger::fddflowPotential_dstress().

◆ hardPotential()

void MultiPlasticityRawComponentAssembler::hardPotential ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< Real > &  h 
)
protectedvirtual

The active hardening potentials (one for each internal parameter and for each yield function) by assumption in the Userobjects, the h[a][alpha] is nonzero only if the surface alpha is part of model a, so we only calculate those here.

Parameters
stressthe stress at which to calculate the hardening potential
intnlvector of internal parameters
activeset of active constraints - only the active hardening potentials are put into "h"
[out]hthe hardening potentials. h[alpha] = hardening potential for yield fcn alpha (and, by the above assumption we know which hardening parameter, a, this belongs to)

Definition at line 262 of file MultiPlasticityRawComponentAssembler.C.

266 {
267  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
268  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
269 
270  h.resize(0);
271  std::vector<unsigned int> active_surfaces_of_model;
272  std::vector<unsigned int>::iterator active_surface;
273  std::vector<Real> model_h;
274  for (unsigned model = 0; model < _num_models; ++model)
275  {
276  activeModelSurfaces(model, active, active_surfaces_of_model);
277  if (active_surfaces_of_model.size() > 0)
278  {
279  _f[model]->hardPotentialV(stress, intnl[model], model_h);
280  for (active_surface = active_surfaces_of_model.begin();
281  active_surface != active_surfaces_of_model.end();
282  ++active_surface)
283  h.push_back(model_h[*active_surface]);
284  }
285  }
286 }

Referenced by MultiPlasticityLinearSystem::calculateConstraints(), MultiPlasticityLinearSystem::calculateJacobian(), and ComputeMultiPlasticityStress::consistentTangentOperator().

◆ modelNumber()

unsigned int MultiPlasticityRawComponentAssembler::modelNumber ( unsigned int  surface)
protected

◆ returnMapAll()

bool MultiPlasticityRawComponentAssembler::returnMapAll ( const RankTwoTensor trial_stress,
const std::vector< Real > &  intnl_old,
const RankFourTensor E_ijkl,
Real  ep_plastic_tolerance,
RankTwoTensor stress,
std::vector< Real > &  intnl,
std::vector< Real > &  pm,
std::vector< Real > &  cumulative_pm,
RankTwoTensor delta_dp,
std::vector< Real > &  yf,
unsigned &  num_successful_plastic_returns,
unsigned &  custom_model 
)
protected

Performs a returnMap for each plastic model using their inbuilt returnMap functions.

Performs a returnMap for each plastic model.

This may be used to quickly ascertain whether a (trial_stress, intnl_old) configuration is admissible, or whether a single model's customized returnMap function can provide a solution to the return-map problem, or whether a full Newton-Raphson approach such as implemented in ComputeMultiPlasticityStress is needed.

There are three cases mentioned below: (A) The (trial_stress, intnl_old) configuration is admissible according to all plastic models (B) The (trial_stress, intnl_old) configuration is inadmissible to exactly one plastic model, and that model can successfully use its customized returnMap function to provide a returned (stress, intnl) configuration, and that configuration is admissible according to all plastic models (C) All other cases. This includes customized returnMap functions failing, or more than one plastic_model being inadmissible, etc

Parameters
trial_stressthe trial stress
intnl_oldthe old values of the internal parameters
E_ijklthe elasticity tensor
ep_plastic_tolerancethe tolerance on the plastic strain
[out]stressis set to trial_stress in case (A) or (C), and the returned value of stress in case (B).
[out]intnlis set to intnl_old in case (A) or (C), and the returned value of intnl in case (B)
[out]pmZero in case (A) or (C), otherwise the plastic multipliers needed to bring about the returnMap in case (B)
[in/out]cumulative_pm cumulative plastic multipliers, updated in case (B), otherwise left untouched
[out]delta_dpis unchanged in case (A) or (C), and is set to the change in plastic strain in case(B)
[out]yfwill contain the yield function values at (stress, intnl)
[out]num_successful_plastic_returnswill be 0 for (A) and (C), and 1 for (B)
Returns
true in case (A) and (B), and false in case (C)

If all models actually signal "elastic" by returning true from their returnMap, and by returning model_plastically_active=0, then yf will contain the yield function values num_successful_plastic_returns will be zero intnl = intnl_old delta_dp will be unchanged from its input value stress will be set to trial_stress pm will be zero cumulative_pm will be unchanged return value will be true num_successful_plastic_returns = 0

If only one model signals "plastically active" by returning true from its returnMap, and by returning model_plastically_active=1, then yf will contain the yield function values num_successful_plastic_returns will be one intnl will be set by the returnMap algorithm delta_dp will be set by the returnMap algorithm stress will be set by the returnMap algorithm pm will be nonzero for the single model, and zero for other models cumulative_pm will be updated return value will be true num_successful_plastic_returns = 1

If >1 model signals "plastically active" or if >=1 model's returnMap fails, then yf will contain the yield function values num_successful_plastic_returns will be set appropriately intnl = intnl_old delta_dp will be unchanged from its input value stress will be set to trial_stress pm will be zero cumulative_pm will be unchanged return value will be true if all returnMap functions returned true, otherwise it will be false num_successful_plastic_returns is set appropriately

Definition at line 599 of file MultiPlasticityRawComponentAssembler.C.

611 {
612  mooseAssert(intnl_old.size() == _num_models,
613  "returnMapAll: Incorrect size of internal parameters");
614  mooseAssert(intnl.size() == _num_models, "returnMapAll: Incorrect size of internal parameters");
615  mooseAssert(pm.size() == _num_surfaces, "returnMapAll: Incorrect size of pm");
616 
617  num_successful_plastic_returns = 0;
618  yf.resize(0);
619  pm.assign(_num_surfaces, 0.0);
620 
621  RankTwoTensor returned_stress; // each model will give a returned_stress. if only one model is
622  // plastically active, i set stress=returned_stress, so as to
623  // record this returned value
624  std::vector<Real> model_f;
625  RankTwoTensor model_delta_dp;
626  std::vector<Real> model_pm;
627  bool trial_stress_inadmissible;
628  bool successful_return = true;
629  unsigned the_single_plastic_model = 0;
630  bool using_custom_return_map = true;
631 
632  // run through all the plastic models, performing their
633  // returnMap algorithms.
634  // If one finds (trial_stress, intnl) inadmissible and
635  // successfully returns, break from the loop to evaluate
636  // all the others at that returned stress
637  for (unsigned model = 0; model < _num_models; ++model)
638  {
639  if (using_custom_return_map)
640  {
641  model_pm.assign(_f[model]->numberSurfaces(), 0.0);
642  bool model_returned = _f[model]->returnMap(trial_stress,
643  intnl_old[model],
644  E_ijkl,
645  ep_plastic_tolerance,
646  returned_stress,
647  intnl[model],
648  model_pm,
649  model_delta_dp,
650  model_f,
651  trial_stress_inadmissible);
652  if (!trial_stress_inadmissible)
653  {
654  // in the elastic zone: record the yield-function values (returned_stress, intnl, model_pm
655  // and model_delta_dp are undefined)
656  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces();
657  ++model_surface)
658  yf.push_back(model_f[model_surface]);
659  }
660  else if (trial_stress_inadmissible && !model_returned)
661  {
662  // in the plastic zone, and the customized returnMap failed
663  // for some reason (or wasn't implemented). The coder
664  // should have correctly returned model_f(trial_stress, intnl_old)
665  // so record them
666  // (returned_stress, intnl, model_pm and model_delta_dp are undefined)
667  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces();
668  ++model_surface)
669  yf.push_back(model_f[model_surface]);
670  // now there's almost zero point in using the custom
671  // returnMap functions
672  using_custom_return_map = false;
673  successful_return = false;
674  }
675  else
676  {
677  // in the plastic zone, and the customized returnMap
678  // succeeded.
679  // record the first returned_stress and delta_dp if everything is going OK
680  // as they could be the actual answer
681  if (trial_stress_inadmissible)
682  num_successful_plastic_returns++;
683  the_single_plastic_model = model;
684  stress = returned_stress;
685  // note that i break here, and don't push_back
686  // model_f to yf. So now yf contains only the values of
687  // model_f from previous models to the_single_plastic_model
688  // also i don't set delta_dp = model_delta_dp yet, because
689  // i might find problems later on
690  // also, don't increment cumulative_pm for the same reason
691 
692  break;
693  }
694  }
695  else
696  {
697  // not using custom returnMap functions because one
698  // has already failed and that one said trial_stress
699  // was inadmissible. So now calculate the yield functions
700  // at the trial stress
701  _f[model]->yieldFunctionV(trial_stress, intnl_old[model], model_f);
702  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
703  yf.push_back(model_f[model_surface]);
704  }
705  }
706 
707  if (num_successful_plastic_returns == 0)
708  {
709  // here either all the models were elastic (successful_return=true),
710  // or some were plastic and either the customized returnMap failed
711  // or wasn't implemented (successful_return=false).
712  // In either case, have to set the following:
713  stress = trial_stress;
714  for (unsigned model = 0; model < _num_models; ++model)
715  intnl[model] = intnl_old[model];
716  return successful_return;
717  }
718 
719  // Now we know that num_successful_plastic_returns == 1 and all the other
720  // models (with model number < the_single_plastic_model) must have been
721  // admissible at (trial_stress, intnl). However, all models might
722  // not be admissible at (trial_stress, intnl), so must check that
723  std::vector<Real> yf_at_returned_stress(0);
724  bool all_admissible = true;
725  for (unsigned model = 0; model < _num_models; ++model)
726  {
727  if (model == the_single_plastic_model)
728  {
729  // no need to spend time calculating the yield function: we know its admissible
730  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
731  yf_at_returned_stress.push_back(model_f[model_surface]);
732  continue;
733  }
734  _f[model]->yieldFunctionV(stress, intnl_old[model], model_f);
735  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
736  {
737  if (model_f[model_surface] > _f[model]->_f_tol)
738  // bummer, this model is not admissible at the returned_stress
739  all_admissible = false;
740  yf_at_returned_stress.push_back(model_f[model_surface]);
741  }
742  if (!all_admissible)
743  // no point in continuing computing yield functions
744  break;
745  }
746 
747  if (!all_admissible)
748  {
749  // we tried using the returned value of stress predicted by
750  // the_single_plastic_model, but it wasn't admissible according
751  // to other plastic models. We need to set:
752  stress = trial_stress;
753  for (unsigned model = 0; model < _num_models; ++model)
754  intnl[model] = intnl_old[model];
755  // and calculate the remainder of the yield functions at trial_stress
756  for (unsigned model = the_single_plastic_model; model < _num_models; ++model)
757  {
758  _f[model]->yieldFunctionV(trial_stress, intnl[model], model_f);
759  for (unsigned model_surface = 0; model_surface < _f[model]->numberSurfaces(); ++model_surface)
760  yf.push_back(model_f[model_surface]);
761  }
762  num_successful_plastic_returns = 0;
763  return false;
764  }
765 
766  // So the customized returnMap algorithm can provide a returned
767  // (stress, intnl) configuration, and that is admissible according
768  // to all plastic models
769  yf.resize(0);
770  for (unsigned surface = 0; surface < yf_at_returned_stress.size(); ++surface)
771  yf.push_back(yf_at_returned_stress[surface]);
772  delta_dp = model_delta_dp;
773  for (unsigned model_surface = 0; model_surface < _f[the_single_plastic_model]->numberSurfaces();
774  ++model_surface)
775  {
776  cumulative_pm[_surfaces_given_model[the_single_plastic_model][model_surface]] +=
777  model_pm[model_surface];
778  pm[_surfaces_given_model[the_single_plastic_model][model_surface]] = model_pm[model_surface];
779  }
780  custom_model = the_single_plastic_model;
781  return true;
782 }

Referenced by ComputeMultiPlasticityStress::quickStep().

◆ validParams()

InputParameters MultiPlasticityRawComponentAssembler::validParams ( )
static

Definition at line 16 of file MultiPlasticityRawComponentAssembler.C.

17 {
18  InputParameters params = emptyInputParameters();
19  MooseEnum specialIC("none rock joint", "none");
20  params.addParam<MooseEnum>("specialIC",
21  specialIC,
22  "For certain combinations of plastic models, the set of active "
23  "constraints can be initialized optimally. 'none': no special "
24  "initialization is performed. For all other choices, the "
25  "plastic_models must be chosen to have the following types. 'rock': "
26  "'TensileMulti MohrCoulombMulti'. 'joint': 'WeakPlaneTensile "
27  "WeakPlaneShear'.");
28  params.addParam<std::vector<UserObjectName>>(
29  "plastic_models",
30  "List of names of user objects that define the plastic models that could "
31  "be active for this material. If no plastic_models are provided, only "
32  "elasticity will be used.");
33  params.addClassDescription("RawComponentAssembler class to calculate yield functions, etc, used "
34  "in multi-surface finite-strain plasticity");
35  return params;
36 }

Referenced by MultiPlasticityLinearSystem::validParams().

◆ yieldFunction()

void MultiPlasticityRawComponentAssembler::yieldFunction ( const RankTwoTensor stress,
const std::vector< Real > &  intnl,
const std::vector< bool > &  active,
std::vector< Real > &  f 
)
protectedvirtual

The active yield function(s)

Parameters
stressthe stress at which to calculate the yield function
intnlvector of internal parameters
activeset of active constraints - only the active yield functions are put into "f"
[out]fthe yield function (or functions in the case of multisurface plasticity)

Definition at line 98 of file MultiPlasticityRawComponentAssembler.C.

102 {
103  mooseAssert(intnl.size() == _num_models, "Incorrect size of internal parameters");
104  mooseAssert(active.size() == _num_surfaces, "Incorrect size of active");
105 
106  f.resize(0);
107  std::vector<unsigned int> active_surfaces_of_model;
108  std::vector<unsigned int>::iterator active_surface;
109  std::vector<Real> model_f;
110  for (unsigned model = 0; model < _num_models; ++model)
111  {
112  activeModelSurfaces(model, active, active_surfaces_of_model);
113  if (active_surfaces_of_model.size() > 0)
114  {
115  _f[model]->yieldFunctionV(stress, intnl[model], model_f);
116  for (active_surface = active_surfaces_of_model.begin();
117  active_surface != active_surfaces_of_model.end();
118  ++active_surface)
119  f.push_back(model_f[*active_surface]);
120  }
121  }
122 }

Referenced by ComputeMultiPlasticityStress::buildDumbOrder(), MultiPlasticityLinearSystem::calculateConstraints(), ComputeMultiPlasticityStress::checkAdmissible(), MultiPlasticityDebugger::fddyieldFunction_dintnl(), MultiPlasticityDebugger::fddyieldFunction_dstress(), and ComputeMultiPlasticityStress::returnMap().

Member Data Documentation

◆ _f

std::vector<const TensorMechanicsPlasticModel *> MultiPlasticityRawComponentAssembler::_f
protected

◆ _model_given_surface

std::vector<unsigned int> MultiPlasticityRawComponentAssembler::_model_given_surface
private

given a surface number, this returns the model number

Definition at line 295 of file MultiPlasticityRawComponentAssembler.h.

Referenced by modelNumber(), and MultiPlasticityRawComponentAssembler().

◆ _model_surface_given_surface

std::vector<unsigned int> MultiPlasticityRawComponentAssembler::_model_surface_given_surface
private

given a surface number, this returns the corresponding-model's internal surface number

Definition at line 298 of file MultiPlasticityRawComponentAssembler.h.

Referenced by MultiPlasticityRawComponentAssembler().

◆ _num_models

unsigned int MultiPlasticityRawComponentAssembler::_num_models
protected

◆ _num_surfaces

unsigned int MultiPlasticityRawComponentAssembler::_num_surfaces
protected

Number of surfaces within the plastic models.

For many situations this will be = _num_models since each model will contain just one surface. More generally it is >= _num_models. For instance, Mohr-Coulomb is a single model with 6 surfaces

Definition at line 66 of file MultiPlasticityRawComponentAssembler.h.

Referenced by ComputeMultiPlasticityStress::activeCombinationNumber(), ComputeMultiPlasticityStress::applyKuhnTucker(), buildActiveConstraints(), ComputeMultiPlasticityStress::buildDumbOrder(), MultiPlasticityLinearSystem::calculateConstraints(), MultiPlasticityLinearSystem::calculateJacobian(), MultiPlasticityLinearSystem::calculateRHS(), ComputeMultiPlasticityStress::canAddConstraints(), ComputeMultiPlasticityStress::canIncrementDumb(), ComputeMultiPlasticityStress::changeScheme(), ComputeMultiPlasticityStress::checkAdmissible(), MultiPlasticityDebugger::checkDerivatives(), MultiPlasticityDebugger::checkJacobian(), ComputeMultiPlasticityStress::checkKuhnTucker(), MultiPlasticityDebugger::checkSolution(), ComputeMultiPlasticityStress::ComputeMultiPlasticityStress(), ComputeMultiPlasticityStress::computeQpStress(), ComputeMultiPlasticityStress::consistentTangentOperator(), dflowPotential_dintnl(), dflowPotential_dstress(), dhardPotential_dintnl(), dhardPotential_dstress(), MultiPlasticityDebugger::dof_included(), dyieldFunction_dintnl(), dyieldFunction_dstress(), MultiPlasticityLinearSystem::eliminateLinearDependence(), MultiPlasticityDebugger::fddflowPotential_dintnl(), MultiPlasticityDebugger::fddflowPotential_dstress(), MultiPlasticityDebugger::fddyieldFunction_dintnl(), MultiPlasticityDebugger::fddyieldFunction_dstress(), MultiPlasticityDebugger::fdJacobian(), flowPotential(), hardPotential(), ComputeMultiPlasticityStress::incrementDumb(), ComputeMultiPlasticityStress::initQpStatefulProperties(), MultiPlasticityRawComponentAssembler(), MultiPlasticityLinearSystem::nrStep(), ComputeMultiPlasticityStress::numberActive(), MultiPlasticityDebugger::outputAndCheckDebugParameters(), ComputeMultiPlasticityStress::plasticStep(), ComputeMultiPlasticityStress::reinstateLinearDependentConstraints(), ComputeMultiPlasticityStress::residual2(), ComputeMultiPlasticityStress::returnMap(), returnMapAll(), ComputeMultiPlasticityStress::singleStep(), and yieldFunction().

◆ _params

const InputParameters& MultiPlasticityRawComponentAssembler::_params
protected

◆ _specialIC

MooseEnum MultiPlasticityRawComponentAssembler::_specialIC
protected

Allows initial set of active constraints to be chosen optimally.

Definition at line 72 of file MultiPlasticityRawComponentAssembler.h.

Referenced by buildActiveConstraints(), and MultiPlasticityRawComponentAssembler().

◆ _surfaces_given_model

std::vector<std::vector<unsigned int> > MultiPlasticityRawComponentAssembler::_surfaces_given_model
protected

_surfaces_given_model[model_number] = vector of surface numbers for this model

Definition at line 69 of file MultiPlasticityRawComponentAssembler.h.

Referenced by activeModelSurfaces(), activeSurfaces(), anyActiveSurfaces(), MultiPlasticityRawComponentAssembler(), ComputeMultiPlasticityStress::quickStep(), and returnMapAll().


The documentation for this class was generated from the following files:
MultiPlasticityRawComponentAssembler::_model_given_surface
std::vector< unsigned int > _model_given_surface
given a surface number, this returns the model number
Definition: MultiPlasticityRawComponentAssembler.h:295
MultiPlasticityRawComponentAssembler::_specialIC
MooseEnum _specialIC
Allows initial set of active constraints to be chosen optimally.
Definition: MultiPlasticityRawComponentAssembler.h:72
MultiPlasticityRawComponentAssembler::_num_surfaces
unsigned int _num_surfaces
Number of surfaces within the plastic models.
Definition: MultiPlasticityRawComponentAssembler.h:66
MultiPlasticityRawComponentAssembler::_f
std::vector< const TensorMechanicsPlasticModel * > _f
User objects that define the yield functions, flow potentials, etc.
Definition: MultiPlasticityRawComponentAssembler.h:75
MultiPlasticityRawComponentAssembler::activeModelSurfaces
void activeModelSurfaces(int model, const std::vector< bool > &active, std::vector< unsigned int > &active_surfaces_of_model)
Returns the internal surface number(s) of the active surfaces of the given model This may be of size=...
Definition: MultiPlasticityRawComponentAssembler.C:811
MultiPlasticityRawComponentAssembler::_num_models
unsigned int _num_models
Number of plastic models for this material.
Definition: MultiPlasticityRawComponentAssembler.h:57
MultiPlasticityRawComponentAssembler::_params
const InputParameters & _params
Definition: MultiPlasticityRawComponentAssembler.h:54
MultiPlasticityRawComponentAssembler::buildActiveConstraintsJoint
void buildActiveConstraintsJoint(const std::vector< Real > &f, const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &Eijkl, std::vector< bool > &act)
"Joint" version Constructs a set of active constraints, given the yield functions,...
Definition: MultiPlasticityRawComponentAssembler.C:381
RankTwoTensorTempl< Real >
MultiPlasticityRawComponentAssembler::_model_surface_given_surface
std::vector< unsigned int > _model_surface_given_surface
given a surface number, this returns the corresponding-model's internal surface number
Definition: MultiPlasticityRawComponentAssembler.h:298
MultiPlasticityRawComponentAssembler::_surfaces_given_model
std::vector< std::vector< unsigned int > > _surfaces_given_model
_surfaces_given_model[model_number] = vector of surface numbers for this model
Definition: MultiPlasticityRawComponentAssembler.h:69
MultiPlasticityRawComponentAssembler::buildActiveConstraintsRock
void buildActiveConstraintsRock(const std::vector< Real > &f, const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &Eijkl, std::vector< bool > &act)
"Rock" version Constructs a set of active constraints, given the yield functions, f.
Definition: MultiPlasticityRawComponentAssembler.C:422