OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
PPSolveFactory.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 <vector>
29 #include <iostream>
30 #include <string>
31 #include <sstream>
32 
33 #include <gsl/gsl_sf_pow_int.h>
34 #include <gsl/gsl_sf_gamma.h>
35 
36 #include "Utility/Inform.h" // ippl
37 
39 
43 
44 extern Inform* gmsg;
45 
46 namespace interpolation {
47 
49 
50 namespace {
51  std::vector<double> getDeltaPos(Mesh* mesh) {
52  Mesh::Iterator it = mesh->begin();
53  std::vector<double> pos1 = it.getPosition();
54  for (size_t i = 0; i < it.getState().size(); ++i)
55  it[i]++;
56  std::vector<double> pos2 = it.getPosition();
57  for (size_t i = 0; i < pos2.size(); ++i)
58  pos2[i] -= pos1[i];
59  return pos2;
60  }
61 }
62 
64  std::vector<std::vector<double> > values,
65  int polyPatchOrder,
66  int smoothingOrder)
67  : polyPatchOrder_m(polyPatchOrder),
68  smoothingOrder_m(smoothingOrder),
69  polyDim_m(0),
70  polyMesh_m(),
71  points_m(points),
72  values_m(values),
73  polynomials_m(),
74  thisPoints_m(),
75  thisValues_m(),
76  derivPoints_m(),
77  derivValues_m(),
78  derivIndices_m(),
79  edgePoints_m(),
80  smoothingPoints_m() {
81  if (points == NULL) {
83  "PPSolveFactory::Solve",
84  "NULL Mesh input to Solve");
85  }
86  if (points->end().toInteger() != int(values.size())) {
87  std::stringstream ss;
88  ss << "Mismatch between Mesh and values in Solve. Mesh size was "
89  << points->end().toInteger() << " but field values size was "
90  << values.size();
92  "PPSolveFactory::Solve",
93  ss.str());
94  }
95  if (smoothingOrder < polyPatchOrder) {
97  "PPSolveFactory::Solve",
98  "Polynomial order must be <= smoothing order in Solve");
99  }
100  polyMesh_m = points->dual();
101  if (polyMesh_m == NULL)
103  "PPSolveFactory::Solve",
104  "Dual of Mesh was NULL");
105  polyDim_m = values[0].size();
106  int posDim = points->getPositionDimension();
107 
108  edgePoints_m = std::vector< std::vector< std::vector<int> > >(posDim+1);
109  for (int i = 1; i <= posDim; ++i)
111  (i, 0, smoothingOrder-polyPatchOrder);
113  (posDim, polyPatchOrder_m-1, polyPatchOrder_m);
114 }
115 
117  Mesh::Iterator outOfBoundsIt) {
118  const Mesh* mesh = outOfBoundsIt.getMesh();
119  std::vector<int> state = outOfBoundsIt.getState();
120  std::vector<int> delta = std::vector<int>(state.size(), 2);
121  Mesh::Iterator end = mesh->end()-1;
122  std::vector<double> pos1 = mesh->begin().getPosition();
123  std::vector<double> pos2 = Mesh::Iterator(delta, mesh).getPosition();
124  for (size_t i = 0; i < state.size(); ++i)
125  if (state[i] < 1) {
126  state[i] = 1;
127  } else if (state[i] > end[i]) {
128  state[i] = end[i];
129  }
130  std::vector<double> position = Mesh::Iterator(state, mesh).getPosition();
131  state = outOfBoundsIt.getState();
132  for (size_t i = 0; i < state.size(); ++i)
133  if (state[i] < 1) {
134  position[i] = (pos2[i] - pos1[i])*(state[i]-1) + position[i];
135  } else if (state[i] > end[i]) {
136  position[i] = (pos2[i] - pos1[i])*(state[i]-end[i]) + position[i];
137  }
138  return position;
139 }
140 
142  int posDim = points_m->getPositionDimension();
143  std::vector< std::vector<int> > dataPoints =
145  thisPoints_m = std::vector< std::vector<double> >(dataPoints.size());
146  std::vector<double> deltaPos = getDeltaPos(points_m);
147  // now iterate over the index_by_power_list elements - offset; each of
148  // these makes a point in the polynomial fit
149  for (size_t i = 0; i < dataPoints.size(); ++i) {
150  thisPoints_m[i] = std::vector<double>(posDim);
151  for (int j = 0; j < posDim; ++j)
152  thisPoints_m[i][j] = (0.5-dataPoints[i][j])*deltaPos[j];
153  }
154 }
155 
157  thisValues_m = std::vector< std::vector<double> >();
158  int posDim = it.getState().size();
159  std::vector< std::vector<int> > dataPoints =
161  size_t dataPointsSize = dataPoints.size();
162  Mesh::Iterator end = points_m->end()-1;
163  // now iterate over the indexByPowerList elements - offset; each of
164  // these makes a point in the polynomial fit
165  for (size_t i = 0; i < dataPointsSize; ++i) {
166  std::vector<int> itState = it.getState();
167  for (int j = 0; j < posDim; ++j)
168  itState[j]++;
169  Mesh::Iterator itCurrent = Mesh::Iterator(itState, points_m);
170  bool outOfBounds = false; // element is off the edge of the mesh
171  for (int j = 0; j < posDim; ++j) {
172  itCurrent[j] -= dataPoints[i][j];
173  outOfBounds = outOfBounds ||
174  itCurrent[j] < 1 ||
175  itCurrent[j] > end[j];
176  }
177  if (outOfBounds) { // if off the edge, then just constrain to zero
178  thisValues_m.push_back(values_m[it.toInteger()]);
179  } else { // else fit using values
180  thisValues_m.push_back(values_m[itCurrent.toInteger()]);
181  }
182  }
183 }
184 
186  derivPoints_m = std::vector< std::vector<double> >();
187  derivIndices_m = std::vector< std::vector<int> >();
188  derivOrigins_m = std::vector< std::vector<int> >();
189  derivPolyVec_m = std::vector<MVector<double> >();
190  int posDim = points_m->getPositionDimension();
191  // get the outer layer of points
192  int deltaOrder = smoothingOrder_m - polyPatchOrder_m;
193  if (deltaOrder <= 0)
194  return;
195  std::vector<double> deltaPos = getDeltaPos(points_m);
196  // smoothingPoints_m is the last layer of points used for polynomial fit
197  // walk around smoothingPoints_m finding the points used for smoothing
198  for (size_t i = 0; i < smoothingPoints_m.size(); ++i) {
199  // make a list of the axes that are on the edge of the space
200  std::vector<int> equalAxes;
201  for (size_t j = 0; j < smoothingPoints_m[i].size(); ++j)
202  if (smoothingPoints_m[i][j] == polyPatchOrder_m)
203  equalAxes.push_back(j);
204  std::vector<std::vector<int> > edgePoints =
205  edgePoints_m[equalAxes.size()];
206  // note the first point, 0,0, is ignored
207  for (size_t j = 0; j < edgePoints.size(); ++j) {
208  derivIndices_m.push_back(std::vector<int>(posDim));
209  derivOrigins_m.push_back(smoothingPoints_m[i]);
210  derivPoints_m.push_back(std::vector<double>(posDim));
211 
212  for (size_t k = 0; k < smoothingPoints_m[i].size(); ++k) {
213  derivPoints_m.back()[k] = (smoothingPoints_m[i][k]-0.5)*deltaPos[k];
214  }
215  for (size_t k = 0; k < edgePoints[j].size(); ++k) {
216  derivIndices_m.back()[equalAxes[k]] = edgePoints[j][k];
217  }
218  int index = 0;
219  // potential bug ... derivIndexByPower_m is never used, is it
220  // needed?
221  while (true) {
222  bool isEqual = true;
223  std::vector<int> polyIndex =
225  for (int k = 0; k < posDim && isEqual; ++k) {
226  isEqual = polyIndex[k] == derivIndices_m.back()[k];
227  }
228  if (isEqual) {
229  derivIndexByPower_m.push_back(index);
230  break;
231  }
232  ++index;
233  }
234  }
235  }
236 
238  (posDim, smoothingOrder_m);
239  for (size_t i = 0; i < derivPoints_m.size(); ++i) {
240  derivPolyVec_m.push_back(MVector<double>(nPolyCoeffs, 1.));
241  // Product[x^{p_j-q_j}*p_j!/(p_j-q_j)!]
242  for (int j = 0; j < nPolyCoeffs; ++j) {
243  std::vector<int> pVec =
245  std::vector<int> qVec = derivIndices_m[i];
246  for (int l = 0; l < posDim; ++l) {
247  int pow = pVec[l]-qVec[l];
248  if (pow < 0) {
249  derivPolyVec_m.back()(j+1) = 0.;
250  break;
251  }
252  derivPolyVec_m.back()(j+1) *= gsl_sf_fact(pVec[l])/gsl_sf_fact(pow)*gsl_sf_pow_int(-deltaPos[l]/2., pow);
253  }
254  }
255  }
256 }
257 
259  derivValues_m = std::vector< std::vector<double> >(derivPoints_m.size());
260 
261  std::vector<double> itPos = it.getPosition();
262  int posDim = it.getState().size();
263  // get the outer layer of points
264  Mesh::Iterator end = it.getMesh()->end()-1;
265  int deltaOrder = smoothingOrder_m - polyPatchOrder_m;
266  if (deltaOrder <= 0)
267  return;
268  for (size_t i = 0; i < derivPoints_m.size(); ++i) {
269  std::vector<double> point = std::vector<double>(posDim, 0.);
270  Mesh::Iterator nearest = it;
271  for (int j = 0; j < posDim; ++j) {
272  point[j] = itPos[j]-derivPoints_m[i][j];
273  }
274  for (int j = 0; j < posDim; ++j) {
275  nearest[j] += derivOrigins_m[i][j];
276  if (nearest[j] > end[j]) {
277  nearest = it;
278  break;
279  }
280  }
281  if (polynomials_m[nearest.toInteger()] == NULL) {
282  derivValues_m[i] = std::vector<double>(polyDim_m, 0.);
283  } else {
284  MMatrix<double> coeffs =
285  polynomials_m[nearest.toInteger()]->GetCoefficientsAsMatrix();
286  MVector<double> values = coeffs*derivPolyVec_m[i];
287  derivValues_m[i] = std::vector<double>(polyDim_m);
288  for(int j = 0; j < posDim; ++j) {
289  derivValues_m[i][j] = values(j+1);
290  }
291  }
292  }
293 }
294 
296  Mesh::Iterator begin = points_m->begin();
297  Mesh::Iterator end = points_m->end()-1;
298  int meshSize = polyMesh_m->end().toInteger();
299  polynomials_m = std::vector<SquarePolynomialVector*>(meshSize, NULL);
300  // get the list of points that are needed to make a given poly vector
301  getPoints();
302  getDerivPoints();
306  polyDim_m,
307  thisPoints_m,
310  int total = (polyMesh_m->end()-1).toInteger();
311  double oldPercentage = 0.;
312  for (Mesh::Iterator it = polyMesh_m->end()-1;
313  it >= polyMesh_m->begin();
314  --it) {
315  double newPercentage = (total-it.toInteger())/double(total)*100.;
316  if (newPercentage - oldPercentage > 10.) {
317  *gmsg << " Done " << newPercentage << " %" << endl;
318  oldPercentage = newPercentage;
319  }
320  // find the set of values and derivatives for this point in the mesh
321  getValues(it);
322  getDerivs(it);
323  // The polynomial is found using simultaneous equation solve
324  polynomials_m[it.toInteger()] = solver.PolynomialSolve
326  }
328 }
329 
331  std::vector<int> check,
332  size_t checkIndex,
333  size_t polyPower,
334  std::vector<std::vector<int> >& nearbyPoints) {
335  check[checkIndex] = polyPower;
336  nearbyPoints.push_back(check);
337  if (checkIndex+1 == check.size())
338  return;
339  for (int i = 1; i < int(polyPower); ++i)
340  nearbyPointsRecursive(check, checkIndex+1, i, nearbyPoints);
341  for (int i = 0; i <= int(polyPower); ++i) {
342  check[checkIndex] = i;
343  nearbyPointsRecursive(check, checkIndex+1, polyPower, nearbyPoints);
344  }
345 }
346 
347 std::vector<std::vector<int> > PPSolveFactory::getNearbyPointsSquares(
348  int pointDim,
349  int polyOrderLower,
350  int polyOrderUpper) {
351  if (pointDim < 1)
352  throw(GeneralClassicException("PPSolveFactory::getNearbyPointsSquares",
353  "Point dimension must be > 0"));
354  if (polyOrderLower > polyOrderUpper)
355  throw(GeneralClassicException("PPSolveFactory::getNearbyPointsSquares",
356  "Polynomial lower bound must be <= polynomial upper bound"));
357  // polyOrder -1 (or less) means no terms
358  if (polyOrderUpper == polyOrderLower || polyOrderUpper < 0)
359  return std::vector<std::vector<int> >();
360  // polyOrder 0 means constant term
361  std::vector<std::vector<int> > nearbyPoints(1, std::vector<int>(pointDim, 0));
362  // polyOrder > 0 means corresponding poly terms (1 is linear, 2 is quadratic...)
363  for (size_t thisPolyOrder = 0;
364  thisPolyOrder < size_t(polyOrderUpper);
365  ++thisPolyOrder)
366  nearbyPointsRecursive(nearbyPoints[0], 0, thisPolyOrder+1, nearbyPoints);
367  if (polyOrderLower < 0)
368  return nearbyPoints; // no terms in lower bound
369  int nPointsLower = 1;
370  for (int dim = 0; dim < pointDim; ++dim)
371  nPointsLower *= polyOrderLower+1;
372  nearbyPoints.erase(nearbyPoints.begin(), nearbyPoints.begin()+nPointsLower);
373  return nearbyPoints;
374 }
375 }
static std::vector< int > IndexByPower(int index, int nInputVariables)
std::vector< double > outOfBoundsPosition(Mesh::Iterator outOfBoundsIt)
static void nearbyPointsRecursive(std::vector< int > check, size_t checkIndex, size_t polyPower, std::vector< std::vector< int > > &nearbyPoints)
std::vector< std::vector< double > > thisPoints_m
PolynomialCoefficient represents a coefficient in a multi-dimensional polynomial. ...
PolynomialCoefficient Coeff
Inform * gmsg
Definition: Main.cpp:21
virtual Mesh::Iterator end() const =0
virtual Mesh * dual() const =0
void getDerivs(Mesh::Iterator it)
void getValues(Mesh::Iterator it)
PPSolveFactory(Mesh *points, std::vector< std::vector< double > > values, int polyPatchOrder, int smoothingOrder)
std::vector< std::vector< int > > derivIndices_m
std::vector< int > derivIndexByPower_m
SolveFactory is a factory class for solving a set of linear equations to generate a polynomial based ...
Definition: SolveFactory.h:44
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Definition: TpsMath.h:76
virtual Mesh::Iterator begin() const =0
Base class for meshing routines.
Definition: Mesh.h:49
std::vector< std::vector< double > > values_m
std::vector< std::vector< double > > derivValues_m
static std::vector< std::vector< int > > getNearbyPointsSquares(int pointDim, int polyOrderLower, int polyOrderUpper)
std::vector< SquarePolynomialVector * > polynomials_m
virtual void getPosition(double *point) const
std::vector< int > getState() const
Patches together many SquarePolynomialVectors to make a multidimensional polynomial spline...
std::vector< std::vector< double > > derivPoints_m
std::vector< std::vector< double > > thisValues_m
std::vector< std::vector< int > > smoothingPoints_m
std::vector< MVector< double > > derivPolyVec_m
const Mesh * getMesh() const
Definition: Inform.h:41
virtual int getPositionDimension() const =0
std::vector< std::vector< std::vector< int > > > edgePoints_m
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
static unsigned int NumberOfPolynomialCoefficients(int pointDimension, int order)
std::vector< std::vector< int > > derivOrigins_m