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 "NearestPointSurrogate.h" 11 : 12 : registerMooseObject("StochasticToolsApp", NearestPointSurrogate); 13 : 14 : InputParameters 15 170 : NearestPointSurrogate::validParams() 16 : { 17 170 : InputParameters params = SurrogateModel::validParams(); 18 170 : params.addClassDescription("Surrogate that evaluates the value from the nearest point from data " 19 : "in [NearestPointTrainer.md]"); 20 170 : return params; 21 0 : } 22 : 23 85 : NearestPointSurrogate::NearestPointSurrogate(const InputParameters & parameters) 24 : : SurrogateModel(parameters), 25 85 : _sample_points(getModelData<std::vector<std::vector<Real>>>("_sample_points")), 26 255 : _sample_results(getModelData<std::vector<std::vector<Real>>>("_sample_results")) 27 : { 28 85 : } 29 : 30 : Real 31 10700 : NearestPointSurrogate::evaluate(const std::vector<Real> & x) const 32 : { 33 : // Check whether input point has same dimensionality as training data 34 : mooseAssert(_sample_points.size() == x.size(), 35 : "Input point does not match dimensionality of training data."); 36 : 37 10700 : return _sample_results[0][findNearestPoint(x)]; 38 : } 39 : 40 : void 41 50 : NearestPointSurrogate::evaluate(const std::vector<Real> & x, std::vector<Real> & y) const 42 : { 43 : mooseAssert(_sample_points.size() == x.size(), 44 : "Input point does not match dimensionality of training data."); 45 : 46 50 : y.assign(_sample_results.size(), 0.0); 47 : 48 50 : unsigned int idx = findNearestPoint(x); 49 : 50 550 : for (const auto & r : index_range(y)) 51 500 : y[r] = _sample_results[r][idx]; 52 50 : } 53 : 54 : unsigned int 55 10750 : NearestPointSurrogate::findNearestPoint(const std::vector<Real> & x) const 56 : { 57 : unsigned int idx = 0; 58 : 59 : // Container of current minimum distance during training sample loop 60 : Real dist_min = std::numeric_limits<Real>::max(); 61 : 62 10030750 : for (dof_id_type p = 0; p < _sample_points[0].size(); ++p) 63 : { 64 : // Sum over the distance of each point dimension 65 : Real dist = 0; 66 40095200 : for (unsigned int i = 0; i < x.size(); ++i) 67 : { 68 30075200 : Real diff = (x[i] - _sample_points[i][p]); 69 30075200 : dist += diff * diff; 70 : } 71 : 72 : // Check if this training point distance is smaller than the current minimum 73 10020000 : if (dist < dist_min) 74 : { 75 336750 : idx = p; 76 : dist_min = dist; 77 : } 78 : } 79 10750 : return idx; 80 : }