00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043 #include "SundanceFeketeQuadQuadrature.hpp"
00044 #include "SundanceOut.hpp"
00045 #include "SundancePoint.hpp"
00046 #include "SundanceGaussLobatto1D.hpp"
00047 #include "PlayaTabs.hpp"
00048
00049 using namespace Sundance;
00050 using namespace Teuchos;
00051
00052 extern "C"
00053 {
00054
00055
00056 void dgetrf_(const int* M, const int* N, double* A, const int* lda,
00057 const int* iPiv, int* info);
00058
00059
00060 void dgetri_(const int* n, double* a, const int* lda, const int* iPiv,
00061 double* work, const int* lwork, int* info);
00062 }
00063
00064 void FeketeQuadQuadrature::getPoints(int order, Array<double>& wgt, Array<
00065 double>& x, Array<double>& y)
00066 {
00067 int p = order + 3;
00068 p = p + (p % 2);
00069 int nNodes = p / 2;
00070 GaussLobatto1D rule(nNodes, 0.0, 1.0);
00071 Array<double> s = rule.nodes();
00072 Array<double> t = s;
00073 Array<double> w = rule.weights();
00074 int n = rule.nPoints();
00075
00076 wgt.resize(n * n);
00077 x.resize(n * n);
00078 y.resize(n * n);
00079
00080 int k = 0;
00081 for (int i = 0; i < n; i++)
00082 {
00083 double p = s[i];
00084 for (int j = 0; j < n; j++, k++)
00085 {
00086 double q = t[j];
00087 x[k] = p;
00088 y[k] = q;
00089 wgt[k] = w[i] * w[j];
00090 }
00091 }
00092 }
00093
00094 void FeketeQuadQuadrature::computeBasisCoeffs(const int order, Array<double>& basisCoeffs)
00095 {
00096
00097 Array<Point> feketePts;
00098
00099 Array<double> x;
00100 Array<double> y;
00101 Array<double> w;
00102 getPoints(order, w, x, y);
00103 int nFeketePts = w.length();
00104 feketePts.resize(nFeketePts);
00105 for (int i = 0; i < nFeketePts; i++)
00106 feketePts[i] = Point(x[i], y[i]);
00107
00108
00109 basisCoeffs.resize(nFeketePts * nFeketePts);
00110
00111
00112
00113 for (int n = 0; n < nFeketePts; n++)
00114 {
00115
00116
00117 double* start = &(basisCoeffs[n * nFeketePts]);
00118 evalPolynomials(nFeketePts, feketePts[n][0], feketePts[n][1], start);
00119 }
00120
00121
00122
00123 int lapack_err = 0;
00124 Array<int> pivot;
00125 pivot.resize(nFeketePts);
00126 Array<double> work;
00127 work.resize(1);
00128 int lwork = -1;
00129
00130
00131 ::dgetrf_(&nFeketePts, &nFeketePts, &(basisCoeffs[0]), &nFeketePts,
00132 &(pivot[0]), &lapack_err);
00133
00134 TEUCHOS_TEST_FOR_EXCEPTION(
00135 lapack_err != 0,
00136 std::runtime_error,
00137 "FeketeQuadQuadrature::computeBasisCoeffs(): factorization of generalized Vandermonde matrix failed");
00138
00139
00140 ::dgetri_(&nFeketePts, &(basisCoeffs[0]), &nFeketePts, &(pivot[0]),
00141 &(work[0]), &lwork, &lapack_err);
00142 lwork = (int) work[0];
00143 work.resize(lwork);
00144 ::dgetri_(&nFeketePts, &(basisCoeffs[0]), &nFeketePts, &(pivot[0]),
00145 &(work[0]), &lwork, &lapack_err);
00146
00147 TEUCHOS_TEST_FOR_EXCEPTION(
00148 lapack_err != 0,
00149 std::runtime_error,
00150 "FeketeQuadQuadrature::computeBasisCoeffs(): inversion of generalized Vandermonde matrix failed");
00151 }
00152
00153 void FeketeQuadQuadrature::evalPolynomials(int nPts, double x,
00154 double y, double* resultPtr)
00155 {
00156
00157 int i = 0;
00158 int order = (int) sqrt((double) nPts);
00159
00160 double xLegendre = 1.;
00161 double xLegendre1 = 0.;
00162
00163 for (int k = 0; k < order; k++)
00164 {
00165 double yLegendre = 1.;
00166 double yLegendre1 = 0.;
00167 double normxLegendre = sqrt( 2. * k + 1 ) * xLegendre;
00168
00169
00170 for (int l = 0; l < order ; l++)
00171 {
00172 double normyLegendre = sqrt( 2. * l + 1 ) * yLegendre;
00173
00174
00175 resultPtr[i++] = normxLegendre * normyLegendre;
00176
00177 double yLegendre2 = yLegendre1;
00178 yLegendre1 = yLegendre;
00179 yLegendre = ( (2 * l + 1) * (2 * y - 1) * yLegendre1 - l * yLegendre2 ) / ( l + 1 );
00180 }
00181
00182 double xLegendre2 = xLegendre1;
00183 xLegendre1 = xLegendre;
00184 xLegendre = ( (2 * k + 1) * (2 * x - 1) * xLegendre1 - k * xLegendre2 ) / ( k + 1 );
00185 }
00186 }
00187
00188 bool FeketeQuadQuadrature::test(int p)
00189 {
00190 Array<double> w;
00191 Array<double> x;
00192 Array<double> y;
00193
00194 getPoints(p, w, x, y);
00195 bool pass = true;
00196 for (int a = 0; a <= p; a++)
00197 {
00198 int bMax = p - a;
00199 for (int b = 0; b <= bMax; b++)
00200 {
00201 double sum = 0.0;
00202 for (int q = 0; q < w.length(); q++)
00203 {
00204 sum += w[q] * pow(x[q], (double) a) * pow(y[q], (double) b);
00205 }
00206 double err = fabs(sum - exact(a, b));
00207 bool localPass = err < 1.0e-14;
00208 pass = pass && localPass;
00209 if (!localPass)
00210 {
00211 fprintf(stderr,
00212 "order=%d m (%d, %d) q=%22.15g exact=%22.15g\n", p, a,
00213 b, sum, exact(a, b));
00214 std::cerr << "error = " << err << std::endl;
00215 }
00216 }
00217 }
00218 return pass;
00219 }
00220
00221 double FeketeQuadQuadrature::exact(int a, int b)
00222 {
00223 return 1.0 / (a + 1) / (b + 1);
00224 }
00225