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 "Sampler1DReal.h" 11 : 12 : registerMooseObject("ThermalHydraulicsApp", Sampler1DReal); 13 : registerMooseObject("ThermalHydraulicsApp", ADSampler1DReal); 14 : 15 : template <bool is_ad> 16 : InputParameters 17 126 : Sampler1DRealTempl<is_ad>::validParams() 18 : { 19 126 : InputParameters params = GeneralVectorPostprocessor::validParams(); 20 126 : params += SamplerBase::validParams(); 21 126 : params += BlockRestrictable::validParams(); 22 : 23 252 : params.addRequiredParam<std::vector<std::string>>( 24 : "property", "Names of the material properties to be output along a line"); 25 : 26 : // This parameter exists in BlockRestrictable, but it is made required here 27 : // because it is undesirable to use the default, which is to add all blocks. 28 252 : params.addRequiredParam<std::vector<SubdomainName>>( 29 : "block", "The list of block ids (SubdomainID) for which this object will be applied"); 30 : 31 126 : params.addClassDescription( 32 : "Samples material properties at all quadrature points in mesh block(s)"); 33 : 34 126 : return params; 35 0 : } 36 : 37 : template <bool is_ad> 38 64 : Sampler1DRealTempl<is_ad>::Sampler1DRealTempl(const InputParameters & parameters) 39 : : GeneralVectorPostprocessor(parameters), 40 : SamplerBase(parameters, this, _communicator), 41 : BlockRestrictable(this), 42 64 : _mesh(_subproblem.mesh()), 43 64 : _qrule(_assembly.qRule()), 44 64 : _q_point(_assembly.qPoints()) 45 : { 46 128 : std::vector<std::string> material_property_names = getParam<std::vector<std::string>>("property"); 47 126 : for (unsigned int i = 0; i < material_property_names.size(); ++i) 48 : { 49 64 : if (!hasGenericMaterialProperty<Real, is_ad>(material_property_names[i])) 50 4 : mooseError("The material property '" + material_property_names[i] + "' does not exist."); 51 62 : _material_properties.push_back( 52 62 : &getGenericMaterialProperty<Real, is_ad>(material_property_names[i])); 53 : } 54 : 55 62 : SamplerBase::setupVariables(material_property_names); 56 62 : } 57 : 58 : template <bool is_ad> 59 : void 60 59 : Sampler1DRealTempl<is_ad>::initialize() 61 : { 62 59 : SamplerBase::initialize(); 63 59 : } 64 : 65 : template <bool is_ad> 66 : void 67 59 : Sampler1DRealTempl<is_ad>::execute() 68 : { 69 59 : std::vector<Real> values(_material_properties.size()); 70 : 71 : std::unordered_set<unsigned int> needed_mat_props; 72 59 : const auto & mp_deps = getMatPropDependencies(); 73 59 : needed_mat_props.insert(mp_deps.begin(), mp_deps.end()); 74 59 : _fe_problem.setActiveMaterialProperties(needed_mat_props, _tid); 75 : 76 59 : ConstElemRange & elem_range = *(_mesh.getActiveLocalElementRange()); 77 6569 : for (typename ConstElemRange::const_iterator el = elem_range.begin(); el != elem_range.end(); 78 : ++el) 79 : { 80 3255 : const Elem * elem = *el; 81 : 82 3255 : if (elem->processor_id() != processor_id()) 83 0 : continue; 84 : 85 3255 : if (!hasBlocks(elem->subdomain_id())) 86 0 : continue; 87 : 88 3255 : _subproblem.setCurrentSubdomainID(elem, _tid); 89 3255 : _subproblem.prepare(elem, _tid); 90 3255 : _subproblem.reinitElem(elem, _tid); 91 : 92 : // Set up Sentinel class so that, even if reinitMaterials() throws, we 93 : // still remember to swap back during stack unwinding. 94 3255 : SwapBackSentinel sentinel(_fe_problem, &FEProblem::swapBackMaterials, _tid); 95 3255 : _fe_problem.reinitMaterials(elem->subdomain_id(), _tid); 96 : 97 8540 : for (unsigned int qp = 0; qp < _qrule->n_points(); ++qp) 98 : { 99 10570 : for (unsigned int j = 0; j < _material_properties.size(); ++j) 100 5285 : values[j] = MetaPhysicL::raw_value((*_material_properties[j])[qp]); 101 : 102 : // use the "x" coordinate as the "id"; at this time, it is not used for anything 103 5285 : addSample(_q_point[qp], _q_point[qp](0), values); 104 : } 105 : } 106 59 : _fe_problem.clearActiveMaterialProperties(_tid); 107 59 : } 108 : 109 : template <bool is_ad> 110 : void 111 59 : Sampler1DRealTempl<is_ad>::finalize() 112 : { 113 59 : SamplerBase::finalize(); 114 59 : }