OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
SolveFactory.cpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2015, Chris Rogers
3 * All rights reserved.
4 * Redistribution and use in source and binary forms, with or without
5 * modification, are permitted provided that the following conditions are met:
6 * 1. Redistributions of source code must retain the above copyright notice,
7 * this list of conditions and the following disclaimer.
8 * 2. Redistributions in binary form must reproduce the above copyright notice,
9 * this list of conditions and the following disclaimer in the documentation
10 * and/or other materials provided with the distribution.
11 * 3. Neither the name of STFC nor the names of its contributors may be used to
12 * endorse or promote products derived from this software without specific
13 * prior written permission.
14 *
15 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
16 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
17 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
18 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
19 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
20 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
21 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
22 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
23 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
24 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
25 * POSSIBILITY OF SUCH DAMAGE.
26 */
27
28#include <sstream>
29
30#include <gsl/gsl_sf_pow_int.h>
31
33
36
37namespace interpolation {
38
39SolveFactory::SolveFactory(int smoothing_order,
40 int point_dim,
41 int value_dim,
42 std::vector< std::vector<double> > positions,
43 std::vector< std::vector<double> > deriv_positions,
44 std::vector< std::vector<int> >& deriv_indices) {
46 square_points_ = PPSolveFactory::getNearbyPointsSquares(point_dim, -1, smoothing_order);
47 if (positions.size() + deriv_positions.size() - n_poly_coeffs_ != 0) {
48 std::stringstream ss;
49 ss << "Total size of positions and deriv_positions ("
50 << positions.size() << "+" << deriv_positions.size() << ") should be "
53 "SolveFactory::SolveFactory",
54 ss.str());
55 }
56 BuildHInvMatrix(positions, deriv_positions, deriv_indices);
57 MMatrix<double> A_temp(value_dim, n_poly_coeffs_, 0.);
58 square_temp_ = SquarePolynomialVector(point_dim, A_temp);
59}
60
62 std::vector< std::vector<double> > positions,
63 std::vector< std::vector<double> > deriv_positions,
64 std::vector< std::vector<int> >& deriv_indices) {
65 int nCoeffs = positions.size();
67 for (int i = 0; i < nCoeffs; ++i) {
68 std::vector<double> poly_vec = MakeSquareVector(positions[i]);
69 for (int j = 0; j < int(poly_vec.size()); ++j) {
70 h_inv_(i+1, j+1) = poly_vec[j];
71 }
72 }
73
74 for (size_t i = 0; i < deriv_positions.size(); ++i) {
75 std::vector<double> deriv_vec = MakeSquareDerivVector(deriv_positions[i],
76 deriv_indices[i]);
77 for (int j = 0; j < n_poly_coeffs_; ++j) {
78 h_inv_(i+1+nCoeffs, j+1) = deriv_vec[j];
79 }
80 }
81 h_inv_.invert();
82}
83
84std::vector<double> SolveFactory::MakeSquareVector(std::vector<double> x) {
85 std::vector<double> square_vector(square_points_.size(), 1.);
86 for (size_t i = 0; i < square_points_.size(); ++i) {
87 std::vector<int>& point = square_points_[i];
88 for (size_t j = 0; j < point.size(); ++j)
89 square_vector[i] *= gsl_sf_pow_int(x[j], point[j]);
90 }
91 return square_vector;
92}
93
95 std::vector<double> x,
96 std::vector<int> deriv_indices) {
97 // vector like Product_i [x_i^{a_i - m_i}*m_i!/(a_i-m_i)]
98 // where:
99 // m_i is given by deriv_indices
100 // a_i are the polynomial vector indices
101
102 std::vector<double> deriv_vec(square_points_.size(), 1.);
103 int square_points_size = square_points_.size();
104 int dim = square_points_[0].size();
105 for (int i = 0; i < square_points_size; ++i) {
106 std::vector<int>& point = square_points_[i];
107 for (int j = 0; j < dim; ++j) {
108 int power = point[j] - deriv_indices[j]; // p_j - q_j
109 if (power < 0) {
110 deriv_vec[i] = 0.;
111 break;
112 } else {
113 // x^(p_j-q_j)
114 deriv_vec[i] *= gsl_sf_pow_int(x[j], power); // x_j^{power}
115 }
116 // p_j*(p_j-1)*(p_j-2)*...*(p_j-q_j)
117 for (int k = point[j]; k > power && k > 0; --k) {
118 deriv_vec[i] *= k;
119 }
120 }
121 }
122 return deriv_vec;
123}
124
126 const std::vector< std::vector<double> >& values,
127 const std::vector< std::vector<double> >& deriv_values) {
128 // Algorithm:
129 // define G_i = vector of F_i and d^pF/dF^p with values taken from coeffs
130 // and derivs respectively
131 // define H_ij = vector of f_i and d^pf/df^p)
132 // Then use G = H A => A = H^-1 G
133 // where A is vector of polynomial coefficients
134 // First set up G_i from coeffs and derivs; then calculate H_ij;
135 // then do the matrix inversion
136 // It is an error to include the same polynomial index more than once in
137 // coeffs or derivs or both; any polynomial indices that are missing will be
138 // assigned 0 value; polynomial order will be taken as the maximum
139 // polynomial order from coeffs and derivs.
140 // PointDimension and ValueDimension will be taken from coeffs and derivs;
141 // it is an error if these do not all have the same dimensions.
142
143 // OPTIMISATION - if we are doing this many times and only changing values,
144 // can reuse H
145 int nCoeffs = values.size();
146 int nDerivs = deriv_values.size();
147 if (values.size()+deriv_values.size() != size_t(n_poly_coeffs_)) {
149 "SolveFactory::PolynomialSolve",
150 "Values and derivatives over or under constrained"
151 );
152 }
153 for (int i = 1; i < nCoeffs && i < n_poly_coeffs_; ++i) {
154 if (values[i].size() < values[0].size()) {
155 throw GeneralClassicException("SolveFactory::PolynomialSolve",
156 "The vector of values is too short");
157 }
158 }
159 for (int i = 0; i < nDerivs; ++ i) {
160 if (deriv_values[i].size() < values[0].size()) {
161 throw GeneralClassicException("SolveFactory::PolynomialSolve",
162 "The vector of derivative values is too short");
163 }
164 }
165
166 int valueDim = 0;
167 if (!values.empty()) {
168 valueDim = values[0].size();
169 } else if (deriv_values.size() != 0) {
170 valueDim = deriv_values[0].size();
171 }
172
173 MMatrix<double> A(valueDim, n_poly_coeffs_, 0.);
174 for (size_t y_index = 0; y_index < values[0].size(); ++y_index) {
176 // First fill the values
177 for (int i = 0; i < nCoeffs && i < n_poly_coeffs_; ++i) {
178 G(i+1) = values[i][y_index];
179 }
180 // Now fill the derivatives
181 for (int i = 0; i < nDerivs; ++i) {
182 G(i+nCoeffs+1) = deriv_values[i][y_index];
183 }
184 MVector<double> A_vec = h_inv_*G;
185 for (int j = 0; j < n_poly_coeffs_; ++j) {
186 A(y_index+1, j+1) = A_vec(j+1);
187 }
188 }
190 temp->SetCoefficients(A);
191 return temp;
192}
193}
void invert()
turns this matrix into its inverse
static std::vector< std::vector< int > > getNearbyPointsSquares(int pointDim, int polyOrderLower, int polyOrderUpper)
SolveFactory(int smoothing_order, int point_dim, int value_dim, std::vector< std::vector< double > > positions, std::vector< std::vector< double > > deriv_positions, std::vector< std::vector< int > > &deriv_indices)
std::vector< double > MakeSquareDerivVector(std::vector< double > position, std::vector< int > derivIndices)
void BuildHInvMatrix(std::vector< std::vector< double > > positions, std::vector< std::vector< double > > deriv_positions, std::vector< std::vector< int > > &deriv_indices)
SquarePolynomialVector * PolynomialSolve(const std::vector< std::vector< double > > &values, const std::vector< std::vector< double > > &deriv_values)
SquarePolynomialVector square_temp_
Definition: SolveFactory.h:122
MMatrix< double > h_inv_
Definition: SolveFactory.h:123
std::vector< std::vector< int > > square_points_
Definition: SolveFactory.h:119
std::vector< double > MakeSquareVector(std::vector< double > position)
SquarePolynomialVector describes a vector of multivariate polynomials.
static unsigned int NumberOfPolynomialCoefficients(int pointDimension, int order)
void SetCoefficients(int pointDim, MMatrix< double > coeff)