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 "SingularTripletReporter.h" 11 : #include "JsonIO.h" 12 : 13 : registerMooseObject("StochasticToolsApp", SingularTripletReporter); 14 : 15 : InputParameters 16 100 : SingularTripletReporter::validParams() 17 : { 18 100 : InputParameters params = GeneralReporter::validParams(); 19 100 : params += MappingInterface::validParams(); 20 : 21 100 : params.addClassDescription("Tool for accessing and outputting the singular triplets of a " 22 : "singular value decomposition in PODMapping."); 23 : 24 200 : params.addRequiredParam<UserObjectName>( 25 : "pod_mapping", "The PODMapping object whose singular triplets should be printed."); 26 : 27 200 : params.addRequiredParam<std::vector<VariableName>>( 28 : "variables", "The names of the variables whose SVD should be printed."); 29 : 30 100 : return params; 31 0 : } 32 : 33 50 : SingularTripletReporter::SingularTripletReporter(const InputParameters & parameters) 34 : : GeneralReporter(parameters), 35 : MappingInterface(this), 36 50 : _variable_names(getParam<std::vector<VariableName>>("variables")), 37 50 : _left_singular_vectors( 38 200 : declareValueByName<std::map<VariableName, std::vector<DenseVector<Real>>>>( 39 : "left_singular_vectors", REPORTER_MODE_ROOT)), 40 50 : _right_singular_vectors( 41 150 : declareValueByName<std::map<VariableName, std::vector<DenseVector<Real>>>>( 42 : "right_singular_vectors", REPORTER_MODE_ROOT)), 43 50 : _singular_values(declareValueByName<std::map<VariableName, std::vector<Real>>>( 44 50 : "singular_values", REPORTER_MODE_ROOT)) 45 : 46 : { 47 100 : for (const auto & vname : _variable_names) 48 : { 49 50 : _left_singular_vectors.emplace(vname, std::vector<DenseVector<Real>>()); 50 50 : _right_singular_vectors.emplace(vname, std::vector<DenseVector<Real>>()); 51 100 : _singular_values.emplace(vname, std::vector<Real>()); 52 : } 53 50 : } 54 : 55 : void 56 50 : SingularTripletReporter::initialSetup() 57 : { 58 50 : _pod_mapping = dynamic_cast<PODMapping *>(&getMapping("pod_mapping")); 59 : 60 50 : if (!_pod_mapping) 61 0 : paramError("pod_mapping", "The given mapping is not a PODMapping!"); 62 : 63 50 : const auto & pod_mapping_variables = _pod_mapping->getVariableNames(); 64 100 : for (const auto & vname : _variable_names) 65 50 : if (std::find(pod_mapping_variables.begin(), pod_mapping_variables.end(), vname) == 66 : pod_mapping_variables.end()) 67 0 : paramError("variables", 68 0 : "The SVD of the requested variable " + vname + 69 : " is not in the PODMapping object!"); 70 : 71 100 : for (const auto & vname : _variable_names) 72 : { 73 50 : _left_singular_vectors[vname].clear(); 74 50 : _right_singular_vectors[vname].clear(); 75 50 : _singular_values[vname].clear(); 76 : } 77 50 : } 78 : 79 : void 80 50 : SingularTripletReporter::execute() 81 : { 82 100 : for (const auto & vname : _variable_names) 83 : { 84 50 : _pod_mapping->buildMapping(vname); 85 50 : _singular_values[vname] = _pod_mapping->singularValues(vname); 86 : 87 50 : const auto & left_modes = _pod_mapping->leftBasis(vname); 88 50 : const auto & right_modes = _pod_mapping->rightBasis(vname); 89 : 90 150 : for (auto mode_i : index_range(left_modes)) 91 : { 92 100 : _left_singular_vectors[vname].push_back(left_modes[mode_i]); 93 100 : _right_singular_vectors[vname].push_back(right_modes[mode_i]); 94 : } 95 : } 96 50 : }