Line data Source code
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 : 10 : #include "CrystalPlasticitySlipRateGSS.h" 11 : 12 : #include <fstream> 13 : 14 : registerMooseObject("TensorMechanicsApp", CrystalPlasticitySlipRateGSS); 15 : 16 : InputParameters 17 78 : CrystalPlasticitySlipRateGSS::validParams() 18 : { 19 78 : InputParameters params = CrystalPlasticitySlipRate::validParams(); 20 156 : params.addParam<std::string>("uo_state_var_name", 21 : "Name of state variable property: Same as " 22 : "state variable user object specified in input " 23 : "file."); 24 78 : params.addClassDescription("Phenomenological constitutive model slip rate class. Override the " 25 : "virtual functions in your class"); 26 78 : return params; 27 0 : } 28 : 29 39 : CrystalPlasticitySlipRateGSS::CrystalPlasticitySlipRateGSS(const InputParameters & parameters) 30 : : CrystalPlasticitySlipRate(parameters), 31 39 : _mat_prop_state_var( 32 39 : getMaterialProperty<std::vector<Real>>(parameters.get<std::string>("uo_state_var_name"))), 33 39 : _pk2(getMaterialPropertyByName<RankTwoTensor>("pk2")), 34 39 : _a0(_variable_size), 35 39 : _xm(_variable_size), 36 117 : _flow_direction(getMaterialProperty<std::vector<RankTwoTensor>>(_name + "_flow_direction")) 37 : { 38 39 : if (_slip_sys_flow_prop_file_name.length() != 0) 39 0 : readFileFlowRateParams(); 40 : else 41 39 : getFlowRateParams(); 42 39 : } 43 : 44 : void 45 0 : CrystalPlasticitySlipRateGSS::readFileFlowRateParams() 46 : { 47 0 : MooseUtils::checkFileReadable(_slip_sys_flow_prop_file_name); 48 : 49 0 : std::ifstream file; 50 0 : file.open(_slip_sys_flow_prop_file_name.c_str()); 51 : 52 : std::vector<Real> vec; 53 0 : vec.resize(_num_slip_sys_flowrate_props); 54 : 55 0 : for (unsigned int i = 0; i < _variable_size; ++i) 56 : { 57 0 : for (unsigned int j = 0; j < _num_slip_sys_flowrate_props; ++j) 58 0 : if (!(file >> vec[j])) 59 0 : mooseError( 60 : "Error CrystalPlasticitySlipRateGSS: Premature end of slip_sys_flow_rate_param file"); 61 : 62 0 : _a0(i) = vec[0]; 63 0 : _xm(i) = vec[1]; 64 : } 65 : 66 0 : file.close(); 67 0 : } 68 : 69 : void 70 39 : CrystalPlasticitySlipRateGSS::getFlowRateParams() 71 : { 72 39 : if (_flowprops.size() <= 0) 73 0 : mooseError("CrystalPlasticitySlipRateGSS: Error in reading flow rate properties: Specify " 74 : "input in .i file or a slip_sys_flow_prop_file_name"); 75 : 76 39 : _a0.resize(_variable_size); 77 39 : _xm.resize(_variable_size); 78 : 79 39 : unsigned int num_data_grp = 2 + _num_slip_sys_flowrate_props; // Number of data per group e.g. 80 : // start_slip_sys, end_slip_sys, 81 : // value1, value2, .. 82 : 83 156 : for (unsigned int i = 0; i < _flowprops.size() / num_data_grp; ++i) 84 : { 85 : Real vs, ve; 86 : unsigned int is, ie; 87 : 88 117 : vs = _flowprops[i * num_data_grp]; 89 117 : ve = _flowprops[i * num_data_grp + 1]; 90 : 91 117 : if (vs <= 0 || ve <= 0) 92 0 : mooseError("CrystalPlasticitySlipRateGSS: Indices in flow rate parameter read must be " 93 : "positive integers: is = ", 94 : vs, 95 : " ie = ", 96 : ve); 97 : 98 117 : if (vs != std::floor(vs) || ve != std::floor(ve)) 99 0 : mooseError("CrystalPlasticitySlipRateGSS: Error in reading flow props: Values specifying " 100 : "start and end number of slip system groups should be integer"); 101 : 102 117 : is = static_cast<unsigned int>(vs); 103 117 : ie = static_cast<unsigned int>(ve); 104 : 105 117 : if (is > ie) 106 0 : mooseError("CrystalPlasticitySlipRateGSS: Start index is = ", 107 : is, 108 : " should be greater than end index ie = ", 109 : ie, 110 : " in flow rate parameter read"); 111 : 112 693 : for (unsigned int j = is; j <= ie; ++j) 113 : { 114 576 : _a0(j - 1) = _flowprops[i * num_data_grp + 2]; 115 576 : _xm(j - 1) = _flowprops[i * num_data_grp + 3]; 116 : } 117 : } 118 : 119 615 : for (unsigned int i = 0; i < _variable_size; ++i) 120 : { 121 576 : if (!(_a0(i) > 0.0 && _xm(i) > 0.0)) 122 : { 123 0 : mooseWarning( 124 : "CrystalPlasticitySlipRateGSS: Non-positive flow rate parameters ", _a0(i), ",", _xm(i)); 125 0 : break; 126 : } 127 : } 128 39 : } 129 : 130 : void 131 171015 : CrystalPlasticitySlipRateGSS::calcFlowDirection(unsigned int qp, 132 : std::vector<RankTwoTensor> & flow_direction) const 133 : { 134 171015 : DenseVector<Real> mo(LIBMESH_DIM * _variable_size), no(LIBMESH_DIM * _variable_size); 135 : 136 : // Update slip direction and normal with crystal orientation 137 2518107 : for (unsigned int i = 0; i < _variable_size; ++i) 138 : { 139 9388368 : for (const auto j : make_range(Moose::dim)) 140 : { 141 7041276 : mo(i * LIBMESH_DIM + j) = 0.0; 142 28165104 : for (const auto k : make_range(Moose::dim)) 143 21123828 : mo(i * LIBMESH_DIM + j) = 144 21123828 : mo(i * LIBMESH_DIM + j) + _crysrot[qp](j, k) * _mo(i * LIBMESH_DIM + k); 145 : } 146 : 147 9388368 : for (const auto j : make_range(Moose::dim)) 148 : { 149 7041276 : no(i * LIBMESH_DIM + j) = 0.0; 150 28165104 : for (const auto k : make_range(Moose::dim)) 151 21123828 : no(i * LIBMESH_DIM + j) = 152 21123828 : no(i * LIBMESH_DIM + j) + _crysrot[qp](j, k) * _no(i * LIBMESH_DIM + k); 153 : } 154 : } 155 : 156 : // Calculate Schmid tensor and resolved shear stresses 157 2518107 : for (unsigned int i = 0; i < _variable_size; ++i) 158 9388368 : for (const auto j : make_range(Moose::dim)) 159 28165104 : for (const auto k : make_range(Moose::dim)) 160 21123828 : flow_direction[i](j, k) = mo(i * LIBMESH_DIM + j) * no(i * LIBMESH_DIM + k); 161 171015 : } 162 : 163 : bool 164 619900 : CrystalPlasticitySlipRateGSS::calcSlipRate(unsigned int qp, Real dt, std::vector<Real> & val) const 165 : { 166 619900 : DenseVector<Real> tau(_variable_size); 167 : 168 9012376 : for (unsigned int i = 0; i < _variable_size; ++i) 169 8392476 : tau(i) = _pk2[qp].doubleContraction(_flow_direction[qp][i]); 170 : 171 9003460 : for (unsigned int i = 0; i < _variable_size; ++i) 172 : { 173 8384303 : val[i] = _a0(i) * std::pow(std::abs(tau(i) / _mat_prop_state_var[qp][i]), 1.0 / _xm(i)) * 174 8384303 : std::copysign(1.0, tau(i)); 175 8384303 : if (std::abs(val[i] * dt) > _slip_incr_tol) 176 : { 177 : #ifdef DEBUG 178 : mooseWarning("Maximum allowable slip increment exceeded ", std::abs(val[i]) * dt); 179 : #endif 180 : return false; 181 : } 182 : } 183 : 184 : return true; 185 : } 186 : 187 : bool 188 619157 : CrystalPlasticitySlipRateGSS::calcSlipRateDerivative(unsigned int qp, 189 : Real /*dt*/, 190 : std::vector<Real> & val) const 191 : { 192 619157 : DenseVector<Real> tau(_variable_size); 193 : 194 9002717 : for (unsigned int i = 0; i < _variable_size; ++i) 195 8383560 : tau(i) = _pk2[qp].doubleContraction(_flow_direction[qp][i]); 196 : 197 9002717 : for (unsigned int i = 0; i < _variable_size; ++i) 198 8383560 : val[i] = _a0(i) / _xm(i) * 199 8383560 : std::pow(std::abs(tau(i) / _mat_prop_state_var[qp][i]), 1.0 / _xm(i) - 1.0) / 200 : _mat_prop_state_var[qp][i]; 201 : 202 619157 : return true; 203 : }