OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
TwoPolynomial.cpp
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2018, Martin Duy Tat
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 <cmath>
29 #include <iostream>
30 #include <stdexcept>
31 #include "gsl/gsl_fft_complex.h"
32 #include "TwoPolynomial.h"
33 
34 namespace polynomial {
35 
36 TwoPolynomial::TwoPolynomial(): maxXorder_m(0), maxSorder_m(0) {
37  std::vector<int> temp;
38  temp.push_back(0);
39  coefficients_m.push_back(temp);
40 }
41 
42 TwoPolynomial::TwoPolynomial(const std::vector<std::vector<int>> &coefficients):
43  maxXorder_m(coefficients.size() - 1),
44  coefficients_m(coefficients) {
45  if (coefficients.size() == 0) {
46  maxXorder_m = 0;
47  maxSorder_m = 0;
48  std::vector<int> temp;
49  temp.push_back(0);
50  coefficients_m.push_back(temp);
51  return;
52  }
53  maxSorder_m = coefficients[0].size() - 1;
54  std::size_t len = maxSorder_m + 1;
55  for (std::size_t i = 0; i < coefficients.size(); i++) {
56  if (coefficients[i].size() != len) {
57  throw std::length_error("2D vector not rectangular");
58  }
59  if (coefficients[i].size() == 0) {
60  maxXorder_m = 0;
61  maxSorder_m = 0;
62  std::vector<int> temp;
63  temp.push_back(0);
64  coefficients_m.resize(0);
65  coefficients_m.push_back(temp);
66  return;
67  }
68  }
69 }
70 
71 TwoPolynomial::TwoPolynomial(const TwoPolynomial &poly):
72  maxXorder_m(poly.maxXorder_m), maxSorder_m(poly.maxSorder_m),
73  coefficients_m(poly.coefficients_m), dSfactors_m(poly.dSfactors_m) {
74 }
75 
76 TwoPolynomial::~TwoPolynomial() {
77 }
78 
79 TwoPolynomial& TwoPolynomial::operator= (const TwoPolynomial &poly) {
80  maxXorder_m = poly.maxXorder_m;
81  maxSorder_m = poly.maxSorder_m;
82  coefficients_m = poly.coefficients_m;
83  dSfactors_m = poly.dSfactors_m;
84  return *this;
85 }
86 
87 void TwoPolynomial::differentiateX() {
88  if (maxXorder_m == 0 && maxSorder_m == 0 && coefficients_m[0][0] == 0) {
89  return;
90  }
91  if (maxXorder_m == 0) {
92  std::vector<std::vector<int>> temp {{0}};
93  coefficients_m = temp;
94  maxSorder_m = 0;
95  return;
96  }
97  for (std::size_t i = 0; i < maxXorder_m; i++) {
98  for (std::size_t j = 0; j <= maxSorder_m; j++) {
99  coefficients_m[i][j] = coefficients_m[i + 1][j] * (i + 1);
100  }
101  }
102  coefficients_m.pop_back();
103  maxXorder_m--;
104 }
105 
106 void TwoPolynomial::differentiateS() {
107  if (maxXorder_m == 0 && maxSorder_m == 0 && coefficients_m[0][0] == 0) {
108  return;
109  }
110  if (maxSorder_m == 0) {
111  std::vector<std::vector<int>> temp {{0}};
112  coefficients_m = temp;
113  maxXorder_m = 0;
114  return;
115  }
116  for (std::size_t i = 0; i <= maxXorder_m; i++) {
117  for (std::size_t j = 0; j < maxSorder_m; j++) {
118  coefficients_m[i][j] = coefficients_m[i][j + 1] * (j + 1);
119  }
120  coefficients_m[i].pop_back();
121  }
122  if (dSfactors_m.size() != 0) {
123  dSfactors_m[0] = dSfactors_m[0] + 1;
124  } else {
125  dSfactors_m.push_back(1);
126  }
127  maxSorder_m--;
128 }
129 
130 void TwoPolynomial::multiplyConstant (const int &constant) {
131  if (maxXorder_m == 0 && maxSorder_m == 0 && coefficients_m[0][0] == 0) {
132  return;
133  }
134  if (constant == 0) {
135  setZero();
136  return;
137  }
138  for (std::size_t i = 0; i <= maxXorder_m; i++) {
139  for (std::size_t j = 0; j <= maxSorder_m; j++) {
140  coefficients_m[i][j] = constant * coefficients_m[i][j];
141  }
142  }
143 }
144 
145 void TwoPolynomial::multiplydSfactors(const TwoPolynomial &poly) {
146  if (poly.dSfactors_m.size() > dSfactors_m.size()) {
147  dSfactors_m.resize(poly.dSfactors_m.size(), 0);
148  }
149  for (std::size_t k = 0; k < poly.dSfactors_m.size(); k++) {
150  dSfactors_m[k] = dSfactors_m[k] + poly.dSfactors_m[k];
151  }
152 }
153 
154 void TwoPolynomial::convert2Dto1Darray(double *vec1, double *vec2,
155  const vectorLengths &nn,
156  const TwoPolynomial &poly) const {
157  for (std::size_t i = 0; i < nn.nx; i++) {
158  for (std::size_t j = 0; j < nn.ny; j++) {
159  if (i < nn.nx1 && j < nn.ny1) {
160  vec1[2 * (nn.ny * i + j)] = coefficients_m[i][j];
161  vec1[2 * (nn.ny * i + j) + 1] = 0;
162  } else {
163  vec1[2 * (nn.ny * i + j)] = 0;
164  vec1[2 * (nn.ny * i + j) + 1] = 0;
165  }
166  if (i < nn.nx2 && j < nn.ny2) {
167  vec2[2 * (nn.ny * i + j)] = poly.coefficients_m[i][j];
168  vec2[2 * (nn.ny * i + j) + 1] = 0;
169  } else {
170  vec2[2 * (nn.ny * i + j)] = 0;
171  vec2[2 * (nn.ny * i + j) + 1] = 0;
172  }
173  }
174  }
175 }
176 
177 void TwoPolynomial::convolution(double *vec1,
178  double *vec2,
179  const std::size_t &nx,
180  const std::size_t &ny) const {
182  for (std::size_t i = 0; i < nx; i++) {
183  double *p1 = &vec1[0] + 2 * i * ny;
184  gsl_fft_complex_radix2_forward(p1, 1, ny);
185  double *p2 = &vec2[0] + 2 * i * ny;
186  gsl_fft_complex_radix2_forward(p2, 1, ny);
187  }
188  for (std::size_t j = 0; j < ny; j++) {
189  double *p1 = &vec1[0] + 2 * j;
190  gsl_fft_complex_radix2_forward(p1, ny, nx);
191  double *p2 = &vec2[0] + 2 * j;
192  gsl_fft_complex_radix2_forward(p2, ny, nx);
193  }
195  for (std::size_t i = 0; i < nx; i++) {
196  for (std::size_t j = 0; j < ny; j++) {
197  double real1 = vec1[2 * (ny * i + j)];
198  double real2 = vec2[2 * (ny * i + j)];
199  double imag1 = vec1[2 * (ny * i + j) + 1];
200  double imag2 = vec2[2 * (ny * i + j) + 1];
201  vec1[2 * (ny * i + j)] = real1 * real2 - imag1 * imag2;
202  vec1[2 * (ny * i + j) + 1] = real1 * imag2 + imag1 * real2;
203  }
204  }
206  for (std::size_t i = 0; i < nx; i++) {
207  double *p1 = &vec1[0] + 2 * i * ny;
208  gsl_fft_complex_radix2_inverse(p1, 1, ny);
209  }
210  for (std::size_t j = 0; j < ny; j++) {
211  double *p1 = &vec1[0] + 2 * j;
212  gsl_fft_complex_radix2_inverse(p1, ny, nx);
213  }
214 }
215 
216 void TwoPolynomial::convert1Dto2Darray(double *vec1,
217  const vectorLengths &nn) {
218  coefficients_m.resize(nn.nx1 + nn.nx2 - 1);
219  for (std::size_t i = 0; i < (nn.nx1 + nn.nx2 - 1); i++) {
220  coefficients_m[i].resize(nn.ny1 + nn.ny2 - 1);
221  for (std::size_t j = 0; j < (nn.ny1 + nn.ny2 - 1); j++) {
222  coefficients_m[i][j] = std::round(vec1[2 * (nn.ny * i + j)]);
223  }
224  }
225 }
226 
227 void TwoPolynomial::multiplyPolynomial(const TwoPolynomial &poly) {
228  if (maxXorder_m == 0 && maxSorder_m == 0 && coefficients_m[0][0] == 0) {
229  return;
230  }
231  if (poly.maxXorder_m == 0 &&
232  poly.maxSorder_m == 0 &&
233  poly.coefficients_m[0][0] == 0) {
234  setZero();
235  return;
236  }
238  multiplydSfactors(poly);
240  vectorLengths nn;
241  nn.nx = 1;
242  nn.ny = 1;
243  nn.nx1 = coefficients_m.size();
244  nn.nx2 = poly.coefficients_m.size();
245  nn.ny1 = coefficients_m[0].size();
246  nn.ny2 = poly.coefficients_m[0].size();
247  while ((nn.nx1 + nn.nx2 - 1) > nn.nx) {
248  nn.nx *= 2;
249  }
250  while ((nn.ny1 + nn.ny2 - 1) > nn.ny) {
251  nn.ny *= 2;
252  }
254  double *vec1 = new double[2 * nn.nx * nn.ny];
255  double *vec2 = new double[2 * nn.nx * nn.ny];
257  convert2Dto1Darray(vec1, vec2, nn, poly);
259  convolution(vec1, vec2, nn.nx, nn.ny);
261  convert1Dto2Darray(vec1, nn);
262  delete[] vec1;
263  delete[] vec2;
264  maxXorder_m = coefficients_m.size() - 1;
265  maxSorder_m = coefficients_m[0].size() - 1;
266 }
267 
268 void TwoPolynomial::printPolynomial() const {
269  if (maxXorder_m == 0 && maxSorder_m == 0 && coefficients_m[0][0] == 0) {
270  return;
271  }
272  std::cout << "(";
273  for (std::size_t i = 0; i <= maxXorder_m; i++) {
274  for (std::size_t j = 0; j <= maxSorder_m; j++) {
275  if (coefficients_m[i][j] > 0) {
276  std::cout << " + " << coefficients_m[i][j];
277  std::cout << "x^" << i << "s^" << j;
278  } else if (coefficients_m[i][j] < 0) {
279  std::cout << " - " << -coefficients_m[i][j];
280  std::cout << "x^" << i << "s^" << j;
281  }
282  }
283  }
284  std::cout << ")";
285  for (std::size_t i = 0; i < dSfactors_m.size(); i++) {
286  std::cout << "(d^" << i + 1 << "S/ds)^" << dSfactors_m[i];
287  }
288 }
289 
290 int TwoPolynomial::getCoefficient(const std::size_t &Xorder,
291  const std::size_t &Sorder) const {
292  if (Xorder > maxXorder_m || Sorder > maxSorder_m) {
293  return 0;
294  }
295  return coefficients_m[Xorder][Sorder];
296 }
297 
298 void TwoPolynomial::setCoefficient(const int &coefficient,
299  const std::size_t &Xorder,
300  const std::size_t &Sorder) {
301  if (Xorder > maxXorder_m) {
302  coefficients_m.resize(Xorder + 1);
303  for (std::size_t i = (maxXorder_m + 1); i <= Xorder; i++) {
304  coefficients_m[i].resize(maxSorder_m + 1, 0);
305  }
306  maxXorder_m = coefficients_m.size() - 1;
307  }
308  if (Sorder > maxSorder_m) {
309  for (std::size_t i = 0; i < maxXorder_m; i++) {
310  coefficients_m[i].resize(Sorder + 1, 0);
311  }
312  maxSorder_m = coefficients_m[0].size() - 1;
313  }
314  coefficients_m[Xorder][Sorder] = coefficient;
315 }
316 
317 void TwoPolynomial::setMaxXorder(const std::size_t &maxXorder) {
318  coefficients_m.resize(maxXorder + 1);
319  for (std::size_t i = (maxXorder_m + 1); i <= maxXorder; i++) {
320  coefficients_m[i].resize(maxSorder_m + 1);
321  }
322  maxXorder_m = maxXorder;
323 }
324 
325 void TwoPolynomial::setMaxSorder(const std::size_t &maxSorder) {
326  for (std::size_t i = 0; i <= maxXorder_m; i++) {
327  coefficients_m[i].resize(maxSorder + 1);
328  }
329  maxSorder_m = maxSorder;
330 }
331 
332 void TwoPolynomial::setZero() {
333  maxXorder_m = 0;
334  maxSorder_m = 0;
335  coefficients_m.resize(1);
336  coefficients_m[0].resize(1);
337  coefficients_m[0][0] = 0;
338 }
339 
340 double TwoPolynomial::evaluatePolynomial(const double &x,
341  const double &s) const {
342  if (maxXorder_m == 0 && maxSorder_m == 0 &&
343  coefficients_m[0][0] == 0) {
344  return 0;
345  }
346  double result = 0.0;
347  std::size_t i = maxXorder_m + 1;
348  while (i != 0) {
349  i--;
350  std::size_t j = maxSorder_m + 1;
351  double temp = 0.0;
352  while (j != 0) {
353  j--;
354  temp = temp * s + coefficients_m[i][j];
355  }
356  result = result * x + temp;
357  }
358  return result;
359 }
360 
361 void TwoPolynomial::addPolynomial(const TwoPolynomial &poly) {
362  if (maxXorder_m < poly.maxXorder_m) {
363  setMaxXorder(poly.maxXorder_m);
364  }
365  if (maxSorder_m < poly.maxSorder_m) {
366  setMaxSorder(poly.maxSorder_m);
367  }
368  for (std::size_t i = 0; i <= poly.maxXorder_m; i++) {
369  for (std::size_t j = 0; j <= poly.maxSorder_m; j++) {
370  coefficients_m[i][j] += poly.coefficients_m[i][j];
371  }
372  }
373 }
374 
375 bool operator < (const TwoPolynomial &left, const TwoPolynomial &right) {
376  std::size_t nleft = left.dSfactors_m.size();
377  std::size_t nright = right.dSfactors_m.size();
378  if (nleft < nright) {
379  return true;
380  } else if (nright < nleft) {
381  return false;
382  } else {
383  for (std::size_t n = 0; n < nleft; n++) {
384  std::size_t elemleft = left.dSfactors_m[n];
385  std::size_t elemright = right.dSfactors_m[n];
386  if (elemleft < elemright) {
387  return true;
388  } else if (elemright < elemleft) {
389  return false;
390  }
391  }
392  return false;
393  }
394 }
395 
396 bool operator == (const TwoPolynomial &left, const TwoPolynomial &right) {
397  if (left < right || right < left) {
398  return false;
399  } else {
400  return true;
401  }
402 }
403 
404 }
bool operator<(const TwoPolynomial &left, const TwoPolynomial &right)
bool operator==(const TwoPolynomial &left, const TwoPolynomial &right)