Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 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 : 10 : #include "CrystalPlasticityStateVariable.h" 11 : 12 : #include <fstream> 13 : 14 : registerMooseObject("SolidMechanicsApp", CrystalPlasticityStateVariable); 15 : 16 : InputParameters 17 156 : CrystalPlasticityStateVariable::validParams() 18 : { 19 156 : InputParameters params = CrystalPlasticityUOBase::validParams(); 20 312 : params.addParam<FileName>( 21 : "state_variable_file_name", 22 : "", 23 : "Name of the file containing the initial values of slip system resistances"); 24 312 : MooseEnum intvar_read_options("file_input inline_input user_input", "inline_input"); 25 312 : params.addParam<MooseEnum>( 26 : "intvar_read_type", 27 : intvar_read_options, 28 : "Read from options for initial value of internal variables: Default from .i file"); 29 312 : params.addParam<std::vector<unsigned int>>("groups", 30 : {}, 31 : "To group the initial values on different " 32 : "slip systems 'format: [start end)', i.e.'0 " 33 : "4 8 11' groups 0-3, 4-7 and 8-11 "); 34 312 : params.addParam<std::vector<Real>>("group_values", 35 : {}, 36 : "The initial values corresponding to each " 37 : "group, i.e. '0.0 1.0 2.0' means 0-4 = 0.0, " 38 : "4-8 = 1.0 and 8-12 = 2.0 "); 39 312 : params.addParam<std::vector<std::string>>("uo_state_var_evol_rate_comp_name", 40 : "Name of state variable evolution rate component " 41 : "property: Same as state variable evolution rate " 42 : "component user object specified in input file."); 43 312 : params.addParam<Real>("zero", 0.0, "Numerical zero for interval variable"); 44 312 : params.addParam<std::vector<Real>>("scale_factor", "Scale factor of individual component."); 45 156 : params.addClassDescription( 46 : "Crystal plasticity state variable class. Override the virtual functions in your class"); 47 156 : return params; 48 156 : } 49 : 50 78 : CrystalPlasticityStateVariable::CrystalPlasticityStateVariable(const InputParameters & parameters) 51 : : CrystalPlasticityUOBase(parameters), 52 78 : _num_mat_state_var_evol_rate_comps( 53 78 : parameters.get<std::vector<std::string>>("uo_state_var_evol_rate_comp_name").size()), 54 78 : _mat_prop_state_var(getMaterialProperty<std::vector<Real>>(_name)), 55 78 : _state_variable_file_name(getParam<FileName>("state_variable_file_name")), 56 156 : _intvar_read_type(getParam<MooseEnum>("intvar_read_type")), 57 156 : _groups(getParam<std::vector<unsigned int>>("groups")), 58 156 : _group_values(getParam<std::vector<Real>>("group_values")), 59 156 : _zero(getParam<Real>("zero")), 60 312 : _scale_factor(getParam<std::vector<Real>>("scale_factor")) 61 : { 62 78 : if (_scale_factor.size() != _num_mat_state_var_evol_rate_comps) 63 0 : mooseError("CrystalPlasticityStateVariable: Scale factor should be have the same size of " 64 : "evolution rate components."); 65 : 66 78 : _mat_prop_state_var_evol_rate_comps.resize(_num_mat_state_var_evol_rate_comps); 67 : 68 156 : for (unsigned int i = 0; i < _num_mat_state_var_evol_rate_comps; ++i) 69 78 : _mat_prop_state_var_evol_rate_comps[i] = &getMaterialProperty<std::vector<Real>>( 70 78 : parameters.get<std::vector<std::string>>("uo_state_var_evol_rate_comp_name")[i]); 71 78 : } 72 : 73 : void 74 12416 : CrystalPlasticityStateVariable::initSlipSysProps(std::vector<Real> & val, 75 : const Point & q_point) const 76 : { 77 6208 : switch (_intvar_read_type) 78 : { 79 64 : case 0: 80 64 : readInitialValueFromFile(val); 81 64 : break; 82 6144 : case 1: 83 6144 : readInitialValueFromInline(val); 84 6144 : break; 85 0 : case 2: 86 0 : provideInitialValueByUser(val, q_point); 87 0 : break; 88 0 : default: 89 0 : mooseError("CrystalPlasticityStateVariable: Read option for initial value of internal " 90 : "variables is not supported."); 91 : } 92 : 93 85312 : for (unsigned int i = 0; i < _variable_size; ++i) 94 79104 : if (val[i] <= 0.0) 95 0 : mooseError("CrystalPlasticityStateVariable: Value of state variables ", i, " non positive"); 96 6208 : } 97 : 98 : void 99 64 : CrystalPlasticityStateVariable::readInitialValueFromFile(std::vector<Real> & val) const 100 : { 101 64 : MooseUtils::checkFileReadable(_state_variable_file_name); 102 : 103 64 : std::ifstream file; 104 64 : file.open(_state_variable_file_name.c_str()); 105 : 106 832 : for (unsigned int i = 0; i < _variable_size; ++i) 107 768 : if (!(file >> val[i])) 108 0 : mooseError("Error CrystalPlasticityStateVariable: Premature end of state_variable file"); 109 : 110 64 : file.close(); 111 64 : } 112 : 113 : void 114 6144 : CrystalPlasticityStateVariable::readInitialValueFromInline(std::vector<Real> & val) const 115 : { 116 6144 : if (_groups.size() <= 0) 117 0 : mooseError("CrystalPlasticityStateVariable: Error in reading initial state variable values: " 118 : "Specify input in .i file or in state_variable file"); 119 6144 : else if (_groups.size() != (_group_values.size() + 1)) 120 0 : mooseError( 121 : "CrystalPlasticityStateVariable: The size of the groups and group_values does not match."); 122 : 123 24576 : for (unsigned int i = 0; i < _groups.size() - 1; ++i) 124 : { 125 : unsigned int is, ie; 126 : 127 18432 : is = _groups[i]; 128 18432 : ie = _groups[i + 1] - 1; 129 : 130 18432 : if (is > ie) 131 0 : mooseError("CrystalPlasticityStateVariable: Start index is = ", 132 : is, 133 : " should be greater than end index ie = ", 134 : ie, 135 : " in state variable read"); 136 : 137 96768 : for (unsigned int j = is; j <= ie; ++j) 138 78336 : val[j] = _group_values[i]; 139 : } 140 6144 : } 141 : 142 : void 143 0 : CrystalPlasticityStateVariable::provideInitialValueByUser(std::vector<Real> & /*val*/, 144 : const Point & /*q_point*/) const 145 : { 146 0 : mooseError("Error CrystalPlasticityStateVariable: User has to overwrite " 147 : "'provideInitialValueByUser' function" 148 : "in order to provide specific initial values based on quadrature point location."); 149 : } 150 : 151 : bool 152 319872 : CrystalPlasticityStateVariable::updateStateVariable(unsigned int qp, 153 : Real dt, 154 : std::vector<Real> & val, 155 : std::vector<Real> & val_old) const 156 : { 157 4725120 : for (unsigned int i = 0; i < _variable_size; ++i) 158 : { 159 4405248 : val[i] = 0.0; 160 8810496 : for (unsigned int j = 0; j < _num_mat_state_var_evol_rate_comps; j++) 161 4405248 : val[i] += (*_mat_prop_state_var_evol_rate_comps[j])[qp][i] * dt * _scale_factor[j]; 162 : } 163 : 164 4725120 : for (unsigned int i = 0; i < _variable_size; ++i) 165 : { 166 4405248 : if (val_old[i] < _zero && val[i] < 0.0) 167 0 : val[i] = val_old[i]; 168 : else 169 4405248 : val[i] = val_old[i] + val[i]; 170 : 171 4405248 : if (val[i] < 0.0) 172 : return false; 173 : } 174 : return true; 175 : }