00001 #include "SundanceFeketeTriangleQuadrature.hpp"
00002 #include "SundanceExceptions.hpp"
00003 #include "SundanceOut.hpp"
00004 #include "SundanceTabs.hpp"
00005
00006 using namespace Sundance;
00007 using namespace Teuchos;
00008
00009
00010
00011
00012
00013
00014 void FeketeTriangleQuadrature::getPoints(int order, Array<double>& wgt, Array<
00015 double>& x, Array<double>& y)
00016 {
00017 Array<double> w;
00018 Array<int> multiplicity;
00019 Array<Array<double> > q;
00020
00021 if (order == 1)
00022 {
00023
00024 multiplicity = tuple(3);
00025 q.resize(1);
00026 w = tuple(1.0 / 3.0);
00027 q[0] = tuple(1.0, 0.0, 0.0);
00028 }
00029 else if (order == 2)
00030 {
00031
00032 multiplicity = tuple(3, 3);
00033 q.resize(2);
00034 w = tuple(0.0, 1.0 / 3.0);
00035 q[0] = tuple(1.0, 0.0, 0.0);
00036 q[1] = tuple(0.0, 0.5, 0.5);
00037 }
00038 else if (order == 3)
00039 {
00040 multiplicity = tuple(1, 3, 6);
00041 q.resize(3);
00042 w = tuple(0.45, 1.0 / 60.0, 1.0 / 12.0);
00043 q[0] = tuple(1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0);
00044 q[1] = tuple(1.0, 0.0, 0.0);
00045 q[2] = tuple(0.000000000000000, 0.276393202250021, 0.723606797749979);
00046 }
00047 else if (order == 4)
00048 {
00049 multiplicity = tuple(3, 3, 3, 6);
00050 q.resize(4);
00051 w = tuple(-0.002443433685378, 0.040684938686721, 0.200360939516545,
00052 0.047365444407723);
00053 q[0] = tuple(1.0, 0.0, 0.0);
00054 q[1] = tuple(0.0, 0.5, 0.5);
00055 q[2] = tuple(0.551583507555305, 0.224208246222347, 0.224208246222347);
00056 q[3] = tuple(0.000000000000000, 0.172673164646012, 0.827326835353988);
00057 }
00058 else if (order == 5)
00059 {
00060 multiplicity = tuple(3, 3, 3, 6, 6);
00061 q.resize(5);
00062 w = tuple(0.005695992854379, 0.125094845896272, 0.116460019007601,
00063 0.012270695896559, 0.030770541890982);
00064 q[0] = tuple(1.0, 0.0, 0.0);
00065 q[1] = tuple(0.684472514501908, 0.157763742749046, 0.157763742749046);
00066 q[2] = tuple(0.171245477332074, 0.414377261333963, 0.414377261333963);
00067 q[3] = tuple(0.000000000000000, 0.117472338035268, 0.882527661964732);
00068 q[4] = tuple(0.000000000000000, 0.357384241759678, 0.642615758240322);
00069 }
00070 else if (order == 6)
00071 {
00072 multiplicity = tuple(1, 3, 3, 3, 6, 6, 6);
00073 q.resize(7);
00074 w = tuple(0.078877036160935, -0.001814501354491, 0.021979435329947,
00075 0.055816419204784, 0.013312944349338, 0.013197985920130,
00076 0.089018887113590);
00077 q[0] = tuple(1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0);
00078 q[1] = tuple(1.0, 0.0, 0.0);
00079 q[2] = tuple(0.0, 0.5, 0.5);
00080 q[3] = tuple(0.769520861253251, 0.115239569373374, 0.115239569373375);
00081 q[4] = tuple(0.000000000000000, 0.084888051860717, 0.915111948139283);
00082 q[5] = tuple(0.000000000000000, 0.265575603264643, 0.734424396735357);
00083 q[6] = tuple(0.127944523240232, 0.317617605129902, 0.554437871629866);
00084 }
00085 else if (order == 9)
00086 {
00087 multiplicity = tuple(1, 3, 3, 3, 3, 6, 6, 6, 6, 6, 6, 6);
00088 q.resize(12);
00089 w = tuple(0.054830037550851, 0.001827879484114, 0.019001036234018,
00090 0.036420799722133, 0.030554508665561, -0.000151089387317,
00091 0.005590374483628, 0.004964725734335, 0.006760863782205,
00092 0.020422905832759, 0.035099032364350, 0.040939402211986);
00093 q[0] = tuple(1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0);
00094 q[1] = tuple(1.0, 0.0, 0.0);
00095 q[2] = tuple(0.888511356254118, 0.055744321872941, 0.055744321872941);
00096 q[3] = tuple(0.643273196352075, 0.178363401823962, 0.178363401823962);
00097 q[4] = tuple(0.067157296821421, 0.466421351589290, 0.466421351589290);
00098 q[5] = tuple(0.000000000000000, 0.040233045916771, 0.959766954083229);
00099 q[6] = tuple(0.000000000000000, 0.130613067447248, 0.869386932552752);
00100 q[7] = tuple(0.000000000000000, 0.261037525094778, 0.738962474905222);
00101 q[8] = tuple(0.000000000000000, 0.417360521166807, 0.582639478833193);
00102 q[9] = tuple(0.063288094133481, 0.162084689378845, 0.774627216487674);
00103 q[10] = tuple(0.066372381175690, 0.304692583561054, 0.628935035263256);
00104 q[11] = tuple(0.185067879630320, 0.326702931348863, 0.488229189020817);
00105 }
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144 else
00145 {
00146 #ifndef TRILINOS_7
00147 SUNDANCE_ERROR("symmetric Fekete quadrature rule order "
00148 << order <<
00149 " for triangles not available");
00150 #else
00151 SUNDANCE_ERROR7("symmetric Fekete quadrature rule order "
00152 << order <<
00153 " for triangles not available");
00154 #endif
00155 }
00156
00157 for (int i = 0; i < q.length(); i++)
00158 {
00159 Array<Array<double> > qPerm;
00160 permute(multiplicity[i], q[i], qPerm);
00161 for (int j = 0; j < multiplicity[i]; j++)
00162 {
00163 x.append(qPerm[j][0]);
00164 y.append(qPerm[j][1]);
00165 wgt.append(w[i]);
00166 }
00167 }
00168
00169 }
00170
00171 bool FeketeTriangleQuadrature::supportsOrder(int order)
00172 {
00173 if ((order >= 1 && order <= 6) || order == 9)
00174 return true;
00175 return false;
00176 }
00177
00178
00179
00180
00181
00182
00183 void FeketeTriangleQuadrature::evalPKDpolynomials(int order, double x,
00184 double y, double* resultPtr)
00185 {
00186 int i = 0;
00187 double Legendre = 1.0;
00188 double Legendre1 = 0.0;
00189 for (int k = 0; k <= order; k++)
00190 {
00191 double Jacobi = 1.0;
00192 double Jacobi1 = 0.0;
00193
00194
00195 for (int l = 0; l <= order - k; l++)
00196 {
00197
00198 resultPtr[i++] = Legendre * pow(x + y, k) * Jacobi;
00199
00200
00201 double Jacobi2 = Jacobi1;
00202 Jacobi1 = Jacobi;
00203 int c1 = 2 * (l + 1) * (l + 2 * k + 2) * (2 * l + 2 * k + 1);
00204 int c2 = (2 * l + 2 * k + 2) * (2 * k + 1) * (2 * k + 1);
00205 int c3 = (2 * l + 2 * k + 1) * (2 * l + 2 * k + 2) * (2 * l + 2 * k
00206 + 3);
00207 int c4 = 2 * (l + 2 * k + 1) * l * (2 * l + 2 * k + 3);
00208 Jacobi = ((c2 + c3 * (1.0 - 2.0 * x - 2.0 * y)) * Jacobi1 - c4
00209 * Jacobi2) / c1;
00210 }
00211
00212
00213 double Legendre2 = Legendre1;
00214 Legendre1 = Legendre;
00215 double xy = x + y;
00216 if (::fabs(xy) > 0.0)
00217 {
00218 xy = (2.0 * k + 1.0) * ((x - y) / xy) * Legendre1;
00219 }
00220 Legendre = (xy - k * Legendre2) / (k + 1);
00221 }
00222 }
00223
00224 void FeketeTriangleQuadrature::permute(int m, const Array<double>& q, Array<
00225 Array<double> >& qPerm)
00226 {
00227 qPerm.resize(m);
00228 if (m == 1)
00229 {
00230 qPerm[0] = q;
00231 }
00232 else if (m == 3)
00233 {
00234 qPerm[0] = tuple(q[0], q[1], q[2]);
00235 qPerm[1] = tuple(q[1], q[0], q[2]);
00236 qPerm[2] = tuple(q[2], q[1], q[0]);
00237 }
00238 else if (m == 6)
00239 {
00240 qPerm[0] = tuple(q[0], q[1], q[2]);
00241 qPerm[1] = tuple(q[0], q[2], q[1]);
00242 qPerm[2] = tuple(q[1], q[0], q[2]);
00243 qPerm[3] = tuple(q[1], q[2], q[0]);
00244 qPerm[4] = tuple(q[2], q[1], q[0]);
00245 qPerm[5] = tuple(q[2], q[0], q[1]);
00246 }
00247 else
00248 {
00249 #ifndef TRILINOS_7
00250 SUNDANCE_ERROR("invalid multiplicity "
00251 << m <<
00252 " in FeketeTriangleQuadrature::permute()");
00253 #else
00254 SUNDANCE_ERROR7("invalid multiplicity "
00255 << m <<
00256 " in FeketeTriangleQuadrature::permute()");
00257 #endif
00258 }
00259 }
00260
00261 bool FeketeTriangleQuadrature::test(int p)
00262 {
00263 Array<double> w;
00264 Array<double> x;
00265 Array<double> y;
00266
00267 getPoints(p, w, x, y);
00268 bool pass = true;
00269
00270 for (int a = 0; a <= p; a++)
00271 {
00272 for (int b = 0; b < p - a; b++)
00273 {
00274 int cMax = p - a - b;
00275 for (int c = 0; c <= cMax; c++)
00276 {
00277 double sum = 0.0;
00278 for (int q = 0; q < w.length(); q++)
00279 {
00280 sum += 0.5 * w[q] * pow(x[q], (double) a) * pow(y[q],
00281 (double) b) * pow(1.0 - x[q] - y[q], (double) c);
00282 }
00283 double err = fabs(sum - exact(a, b, c));
00284 bool localPass = err < 1.0e-14;
00285 pass = pass && localPass;
00286 if (!localPass)
00287 {
00288 fprintf(
00289 stderr,
00290 "order=%d m (%d, %d, %d) q=%22.15g exact=%22.15g\n",
00291 p, a, b, c, sum, exact(a, b, c));
00292 std::cerr << "error = " << err << std::endl;
00293 }
00294 }
00295 }
00296 }
00297 return pass;
00298 }
00299
00300 double FeketeTriangleQuadrature::exact(int a, int b, int c)
00301 {
00302 return fact(a) * fact(b) * fact(c) / fact(a + b + c + 2);
00303 }
00304
00305 double FeketeTriangleQuadrature::fact(int x)
00306 {
00307 if (x == 0)
00308 return 1.0;
00309 return x * fact(x - 1);
00310 }
00311