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 "ADComputeShellStress.h" 11 : #include "RankTwoTensor.h" 12 : #include "RankFourTensor.h" 13 : #include "ADComputeIsotropicElasticityTensorShell.h" 14 : 15 : #include "libmesh/quadrature.h" 16 : #include "libmesh/utility.h" 17 : #include "libmesh/enum_quadrature_type.h" 18 : #include "libmesh/fe_type.h" 19 : #include "libmesh/string_to_enum.h" 20 : #include "libmesh/quadrature_gauss.h" 21 : 22 : registerMooseObject("TensorMechanicsApp", ADComputeShellStress); 23 : 24 : InputParameters 25 160 : ADComputeShellStress::validParams() 26 : { 27 160 : InputParameters params = Material::validParams(); 28 160 : params.addClassDescription("Compute in-plane stress using elasticity for shell"); 29 320 : params.addRequiredParam<std::string>("through_thickness_order", 30 : "Quadrature order in out of plane direction"); 31 160 : return params; 32 0 : } 33 : 34 120 : ADComputeShellStress::ADComputeShellStress(const InputParameters & parameters) 35 120 : : Material(parameters) 36 : { 37 : // get number of quadrature points along thickness based on order 38 : std::unique_ptr<QGauss> t_qrule = std::make_unique<QGauss>( 39 240 : 1, Utility::string_to_enum<Order>(getParam<std::string>("through_thickness_order"))); 40 120 : _t_points = t_qrule->get_points(); 41 120 : _elasticity_tensor.resize(_t_points.size()); 42 120 : _stress.resize(_t_points.size()); 43 120 : _stress_old.resize(_t_points.size()); 44 120 : _strain_increment.resize(_t_points.size()); 45 120 : _covariant_transformation_matrix.resize(_t_points.size()); 46 120 : _global_stress.resize(_t_points.size()); 47 360 : for (unsigned int t = 0; t < _t_points.size(); ++t) 48 : { 49 240 : _elasticity_tensor[t] = 50 240 : &getADMaterialProperty<RankFourTensor>("elasticity_tensor_t_points_" + std::to_string(t)); 51 240 : _stress[t] = &declareADProperty<RankTwoTensor>("stress_t_points_" + std::to_string(t)); 52 240 : _stress_old[t] = 53 240 : &getMaterialPropertyOldByName<RankTwoTensor>("stress_t_points_" + std::to_string(t)); 54 240 : _strain_increment[t] = 55 240 : &getADMaterialProperty<RankTwoTensor>("strain_increment_t_points_" + std::to_string(t)); 56 : // rotation matrix and stress for output purposes only 57 240 : _covariant_transformation_matrix[t] = &getMaterialProperty<RankTwoTensor>( 58 240 : "covariant_transformation_t_points_" + std::to_string(t)); 59 240 : _global_stress[t] = 60 480 : &declareProperty<RankTwoTensor>("global_stress_t_points_" + std::to_string(t)); 61 : } 62 120 : } 63 : 64 : void 65 2212 : ADComputeShellStress::initQpStatefulProperties() 66 : { 67 : // initialize stress tensor to zero 68 6636 : for (unsigned int i = 0; i < _t_points.size(); ++i) 69 4424 : (*_stress[i])[_qp].zero(); 70 2212 : } 71 : 72 : void 73 19880 : ADComputeShellStress::computeQpProperties() 74 : { 75 59640 : for (unsigned int i = 0; i < _t_points.size(); ++i) 76 : { 77 39760 : (*_stress[i])[_qp] = 78 39760 : (*_stress_old[i])[_qp] + (*_elasticity_tensor[i])[_qp] * (*_strain_increment[i])[_qp]; 79 : 80 159040 : for (unsigned int ii = 0; ii < 3; ++ii) 81 477120 : for (unsigned int jj = 0; jj < 3; ++jj) 82 357840 : _unrotated_stress(ii, jj) = MetaPhysicL::raw_value((*_stress[i])[_qp](ii, jj)); 83 39760 : (*_global_stress[i])[_qp] = (*_covariant_transformation_matrix[i])[_qp].transpose() * 84 39760 : _unrotated_stress * (*_covariant_transformation_matrix[i])[_qp]; 85 : } 86 19880 : }