1 /* Discretization tools */ 2 3 #include <petscdt.h> /*I "petscdt.h" I*/ 4 #include <petscblaslapack.h> 5 #include <petsc/private/petscimpl.h> 6 #include <petsc/private/dtimpl.h> 7 #include <petscviewer.h> 8 #include <petscdmplex.h> 9 #include <petscdmshell.h> 10 11 #if defined(PETSC_HAVE_MPFR) 12 #include <mpfr.h> 13 #endif 14 15 const char *const PetscDTNodeTypes[] = {"gaussjacobi", "equispaced", "tanhsinh", "PETSCDTNODES_", NULL}; 16 17 static PetscBool GolubWelschCite = PETSC_FALSE; 18 const char GolubWelschCitation[] = "@article{GolubWelsch1969,\n" 19 " author = {Golub and Welsch},\n" 20 " title = {Calculation of Quadrature Rules},\n" 21 " journal = {Math. Comp.},\n" 22 " volume = {23},\n" 23 " number = {106},\n" 24 " pages = {221--230},\n" 25 " year = {1969}\n}\n"; 26 27 /* Numerical tests in src/dm/dt/tests/ex1.c show that when computing the nodes and weights of Gauss-Jacobi 28 quadrature rules: 29 30 - in double precision, Newton's method and Golub & Welsch both work for moderate degrees (< 100), 31 - in single precision, Newton's method starts producing incorrect roots around n = 15, but 32 the weights from Golub & Welsch become a problem before then: they produces errors 33 in computing the Jacobi-polynomial Gram matrix around n = 6. 34 35 So we default to Newton's method (required fewer dependencies) */ 36 PetscBool PetscDTGaussQuadratureNewton_Internal = PETSC_TRUE; 37 38 PetscClassId PETSCQUADRATURE_CLASSID = 0; 39 40 /*@ 41 PetscQuadratureCreate - Create a PetscQuadrature object 42 43 Collective 44 45 Input Parameter: 46 . comm - The communicator for the PetscQuadrature object 47 48 Output Parameter: 49 . q - The PetscQuadrature object 50 51 Level: beginner 52 53 .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 54 @*/ 55 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 56 { 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 PetscValidPointer(q, 2); 61 ierr = DMInitializePackage();CHKERRQ(ierr); 62 ierr = PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 63 (*q)->dim = -1; 64 (*q)->Nc = 1; 65 (*q)->order = -1; 66 (*q)->numPoints = 0; 67 (*q)->points = NULL; 68 (*q)->weights = NULL; 69 PetscFunctionReturn(0); 70 } 71 72 /*@ 73 PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 74 75 Collective on q 76 77 Input Parameter: 78 . q - The PetscQuadrature object 79 80 Output Parameter: 81 . r - The new PetscQuadrature object 82 83 Level: beginner 84 85 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 86 @*/ 87 PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 88 { 89 PetscInt order, dim, Nc, Nq; 90 const PetscReal *points, *weights; 91 PetscReal *p, *w; 92 PetscErrorCode ierr; 93 94 PetscFunctionBegin; 95 PetscValidPointer(q, 1); 96 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 97 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 98 ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 99 ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 100 ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 101 ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr); 102 ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr); 103 ierr = PetscArraycpy(w, weights, Nc * Nq);CHKERRQ(ierr); 104 ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr); 105 PetscFunctionReturn(0); 106 } 107 108 /*@ 109 PetscQuadratureDestroy - Destroys a PetscQuadrature object 110 111 Collective on q 112 113 Input Parameter: 114 . q - The PetscQuadrature object 115 116 Level: beginner 117 118 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 119 @*/ 120 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 121 { 122 PetscErrorCode ierr; 123 124 PetscFunctionBegin; 125 if (!*q) PetscFunctionReturn(0); 126 PetscValidHeaderSpecific((*q),PETSCQUADRATURE_CLASSID,1); 127 if (--((PetscObject)(*q))->refct > 0) { 128 *q = NULL; 129 PetscFunctionReturn(0); 130 } 131 ierr = PetscFree((*q)->points);CHKERRQ(ierr); 132 ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 133 ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 134 PetscFunctionReturn(0); 135 } 136 137 /*@ 138 PetscQuadratureGetOrder - Return the order of the method 139 140 Not collective 141 142 Input Parameter: 143 . q - The PetscQuadrature object 144 145 Output Parameter: 146 . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 147 148 Level: intermediate 149 150 .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 151 @*/ 152 PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 153 { 154 PetscFunctionBegin; 155 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 156 PetscValidPointer(order, 2); 157 *order = q->order; 158 PetscFunctionReturn(0); 159 } 160 161 /*@ 162 PetscQuadratureSetOrder - Return the order of the method 163 164 Not collective 165 166 Input Parameters: 167 + q - The PetscQuadrature object 168 - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 169 170 Level: intermediate 171 172 .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 173 @*/ 174 PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 175 { 176 PetscFunctionBegin; 177 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 178 q->order = order; 179 PetscFunctionReturn(0); 180 } 181 182 /*@ 183 PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated 184 185 Not collective 186 187 Input Parameter: 188 . q - The PetscQuadrature object 189 190 Output Parameter: 191 . Nc - The number of components 192 193 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 194 195 Level: intermediate 196 197 .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 198 @*/ 199 PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc) 200 { 201 PetscFunctionBegin; 202 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 203 PetscValidPointer(Nc, 2); 204 *Nc = q->Nc; 205 PetscFunctionReturn(0); 206 } 207 208 /*@ 209 PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated 210 211 Not collective 212 213 Input Parameters: 214 + q - The PetscQuadrature object 215 - Nc - The number of components 216 217 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 218 219 Level: intermediate 220 221 .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 222 @*/ 223 PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc) 224 { 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 227 q->Nc = Nc; 228 PetscFunctionReturn(0); 229 } 230 231 /*@C 232 PetscQuadratureGetData - Returns the data defining the quadrature 233 234 Not collective 235 236 Input Parameter: 237 . q - The PetscQuadrature object 238 239 Output Parameters: 240 + dim - The spatial dimension 241 . Nc - The number of components 242 . npoints - The number of quadrature points 243 . points - The coordinates of each quadrature point 244 - weights - The weight of each quadrature point 245 246 Level: intermediate 247 248 Fortran Notes: 249 From Fortran you must call PetscQuadratureRestoreData() when you are done with the data 250 251 .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 252 @*/ 253 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 254 { 255 PetscFunctionBegin; 256 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 257 if (dim) { 258 PetscValidPointer(dim, 2); 259 *dim = q->dim; 260 } 261 if (Nc) { 262 PetscValidPointer(Nc, 3); 263 *Nc = q->Nc; 264 } 265 if (npoints) { 266 PetscValidPointer(npoints, 4); 267 *npoints = q->numPoints; 268 } 269 if (points) { 270 PetscValidPointer(points, 5); 271 *points = q->points; 272 } 273 if (weights) { 274 PetscValidPointer(weights, 6); 275 *weights = q->weights; 276 } 277 PetscFunctionReturn(0); 278 } 279 280 static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[]) 281 { 282 PetscScalar *Js, *Jinvs; 283 PetscInt i, j, k; 284 PetscBLASInt bm, bn, info; 285 PetscErrorCode ierr; 286 287 PetscFunctionBegin; 288 if (!m || !n) PetscFunctionReturn(0); 289 ierr = PetscBLASIntCast(m, &bm);CHKERRQ(ierr); 290 ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr); 291 #if defined(PETSC_USE_COMPLEX) 292 ierr = PetscMalloc2(m*n, &Js, m*n, &Jinvs);CHKERRQ(ierr); 293 for (i = 0; i < m*n; i++) Js[i] = J[i]; 294 #else 295 Js = (PetscReal *) J; 296 Jinvs = Jinv; 297 #endif 298 if (m == n) { 299 PetscBLASInt *pivots; 300 PetscScalar *W; 301 302 ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 303 304 ierr = PetscArraycpy(Jinvs, Js, m * m);CHKERRQ(ierr); 305 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info)); 306 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 307 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info)); 308 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 309 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 310 } else if (m < n) { 311 PetscScalar *JJT; 312 PetscBLASInt *pivots; 313 PetscScalar *W; 314 315 ierr = PetscMalloc1(m*m, &JJT);CHKERRQ(ierr); 316 ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 317 for (i = 0; i < m; i++) { 318 for (j = 0; j < m; j++) { 319 PetscScalar val = 0.; 320 321 for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k]; 322 JJT[i * m + j] = val; 323 } 324 } 325 326 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info)); 327 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 328 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info)); 329 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 330 for (i = 0; i < n; i++) { 331 for (j = 0; j < m; j++) { 332 PetscScalar val = 0.; 333 334 for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j]; 335 Jinvs[i * m + j] = val; 336 } 337 } 338 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 339 ierr = PetscFree(JJT);CHKERRQ(ierr); 340 } else { 341 PetscScalar *JTJ; 342 PetscBLASInt *pivots; 343 PetscScalar *W; 344 345 ierr = PetscMalloc1(n*n, &JTJ);CHKERRQ(ierr); 346 ierr = PetscMalloc2(n, &pivots, n, &W);CHKERRQ(ierr); 347 for (i = 0; i < n; i++) { 348 for (j = 0; j < n; j++) { 349 PetscScalar val = 0.; 350 351 for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j]; 352 JTJ[i * n + j] = val; 353 } 354 } 355 356 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bn, pivots, &info)); 357 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 358 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info)); 359 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 360 for (i = 0; i < n; i++) { 361 for (j = 0; j < m; j++) { 362 PetscScalar val = 0.; 363 364 for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k]; 365 Jinvs[i * m + j] = val; 366 } 367 } 368 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 369 ierr = PetscFree(JTJ);CHKERRQ(ierr); 370 } 371 #if defined(PETSC_USE_COMPLEX) 372 for (i = 0; i < m*n; i++) Jinv[i] = PetscRealPart(Jinvs[i]); 373 ierr = PetscFree2(Js, Jinvs);CHKERRQ(ierr); 374 #endif 375 PetscFunctionReturn(0); 376 } 377 378 /*@ 379 PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation. 380 381 Collecive on PetscQuadrature 382 383 Input Parameters: 384 + q - the quadrature functional 385 . imageDim - the dimension of the image of the transformation 386 . origin - a point in the original space 387 . originImage - the image of the origin under the transformation 388 . J - the Jacobian of the image: an [imageDim x dim] matrix in row major order 389 - formDegree - transform the quadrature weights as k-forms of this form degree (if the number of components is a multiple of (dim choose formDegree), it is assumed that they represent multiple k-forms) [see PetscDTAltVPullback() for interpretation of formDegree] 390 391 Output Parameters: 392 . Jinvstarq - a quadrature rule where each point is the image of a point in the original quadrature rule, and where the k-form weights have been pulled-back by the pseudoinverse of J to the k-form weights in the image space. 393 394 Note: the new quadrature rule will have a different number of components if spaces have different dimensions. For example, pushing a 2-form forward from a two dimensional space to a three dimensional space changes the number of components from 1 to 3. 395 396 Level: intermediate 397 398 .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 399 @*/ 400 PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq) 401 { 402 PetscInt dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c; 403 const PetscReal *points; 404 const PetscReal *weights; 405 PetscReal *imagePoints, *imageWeights; 406 PetscReal *Jinv; 407 PetscReal *Jinvstar; 408 PetscErrorCode ierr; 409 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 412 PetscCheckFalse(imageDim < PetscAbsInt(formDegree),PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %D-form in %D dimensions", PetscAbsInt(formDegree), imageDim); 413 ierr = PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights);CHKERRQ(ierr); 414 ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize);CHKERRQ(ierr); 415 PetscCheckFalse(Nc % formSize,PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of formSize %D", Nc, formSize); 416 Ncopies = Nc / formSize; 417 ierr = PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize);CHKERRQ(ierr); 418 imageNc = Ncopies * imageFormSize; 419 ierr = PetscMalloc1(Npoints * imageDim, &imagePoints);CHKERRQ(ierr); 420 ierr = PetscMalloc1(Npoints * imageNc, &imageWeights);CHKERRQ(ierr); 421 ierr = PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar);CHKERRQ(ierr); 422 ierr = PetscDTJacobianInverse_Internal(imageDim, dim, J, Jinv);CHKERRQ(ierr); 423 ierr = PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar);CHKERRQ(ierr); 424 for (pt = 0; pt < Npoints; pt++) { 425 const PetscReal *point = &points[pt * dim]; 426 PetscReal *imagePoint = &imagePoints[pt * imageDim]; 427 428 for (i = 0; i < imageDim; i++) { 429 PetscReal val = originImage[i]; 430 431 for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]); 432 imagePoint[i] = val; 433 } 434 for (c = 0; c < Ncopies; c++) { 435 const PetscReal *form = &weights[pt * Nc + c * formSize]; 436 PetscReal *imageForm = &imageWeights[pt * imageNc + c * imageFormSize]; 437 438 for (i = 0; i < imageFormSize; i++) { 439 PetscReal val = 0.; 440 441 for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j]; 442 imageForm[i] = val; 443 } 444 } 445 } 446 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq);CHKERRQ(ierr); 447 ierr = PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights);CHKERRQ(ierr); 448 ierr = PetscFree2(Jinv, Jinvstar);CHKERRQ(ierr); 449 PetscFunctionReturn(0); 450 } 451 452 /*@C 453 PetscQuadratureSetData - Sets the data defining the quadrature 454 455 Not collective 456 457 Input Parameters: 458 + q - The PetscQuadrature object 459 . dim - The spatial dimension 460 . Nc - The number of components 461 . npoints - The number of quadrature points 462 . points - The coordinates of each quadrature point 463 - weights - The weight of each quadrature point 464 465 Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them. 466 467 Level: intermediate 468 469 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 470 @*/ 471 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 472 { 473 PetscFunctionBegin; 474 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 475 if (dim >= 0) q->dim = dim; 476 if (Nc >= 0) q->Nc = Nc; 477 if (npoints >= 0) q->numPoints = npoints; 478 if (points) { 479 PetscValidPointer(points, 5); 480 q->points = points; 481 } 482 if (weights) { 483 PetscValidPointer(weights, 6); 484 q->weights = weights; 485 } 486 PetscFunctionReturn(0); 487 } 488 489 static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v) 490 { 491 PetscInt q, d, c; 492 PetscViewerFormat format; 493 PetscErrorCode ierr; 494 495 PetscFunctionBegin; 496 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);CHKERRQ(ierr);} 497 else {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);} 498 ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 499 if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 500 for (q = 0; q < quad->numPoints; ++q) { 501 ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr); 502 ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr); 503 for (d = 0; d < quad->dim; ++d) { 504 if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 505 ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 506 } 507 ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr); 508 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);} 509 for (c = 0; c < quad->Nc; ++c) { 510 if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 511 ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 512 } 513 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);} 514 ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr); 515 ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr); 516 } 517 PetscFunctionReturn(0); 518 } 519 520 /*@C 521 PetscQuadratureView - Views a PetscQuadrature object 522 523 Collective on quad 524 525 Input Parameters: 526 + quad - The PetscQuadrature object 527 - viewer - The PetscViewer object 528 529 Level: beginner 530 531 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 532 @*/ 533 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 534 { 535 PetscBool iascii; 536 PetscErrorCode ierr; 537 538 PetscFunctionBegin; 539 PetscValidHeader(quad, 1); 540 if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 541 if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);} 542 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 543 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 544 if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);} 545 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 546 PetscFunctionReturn(0); 547 } 548 549 /*@C 550 PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 551 552 Not collective 553 554 Input Parameters: 555 + q - The original PetscQuadrature 556 . numSubelements - The number of subelements the original element is divided into 557 . v0 - An array of the initial points for each subelement 558 - jac - An array of the Jacobian mappings from the reference to each subelement 559 560 Output Parameters: 561 . dim - The dimension 562 563 Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 564 565 Not available from Fortran 566 567 Level: intermediate 568 569 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 570 @*/ 571 PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 572 { 573 const PetscReal *points, *weights; 574 PetscReal *pointsRef, *weightsRef; 575 PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 576 PetscErrorCode ierr; 577 578 PetscFunctionBegin; 579 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 580 PetscValidPointer(v0, 3); 581 PetscValidPointer(jac, 4); 582 PetscValidPointer(qref, 5); 583 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 584 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 585 ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 586 npointsRef = npoints*numSubelements; 587 ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 588 ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 589 for (c = 0; c < numSubelements; ++c) { 590 for (p = 0; p < npoints; ++p) { 591 for (d = 0; d < dim; ++d) { 592 pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 593 for (e = 0; e < dim; ++e) { 594 pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 595 } 596 } 597 /* Could also use detJ here */ 598 for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 599 } 600 } 601 ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 602 ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 603 PetscFunctionReturn(0); 604 } 605 606 /* Compute the coefficients for the Jacobi polynomial recurrence, 607 * 608 * J^{a,b}_n(x) = (cnm1 + cnm1x * x) * J^{a,b}_{n-1}(x) - cnm2 * J^{a,b}_{n-2}(x). 609 */ 610 #define PetscDTJacobiRecurrence_Internal(n,a,b,cnm1,cnm1x,cnm2) \ 611 do { \ 612 PetscReal _a = (a); \ 613 PetscReal _b = (b); \ 614 PetscReal _n = (n); \ 615 if (n == 1) { \ 616 (cnm1) = (_a-_b) * 0.5; \ 617 (cnm1x) = (_a+_b+2.)*0.5; \ 618 (cnm2) = 0.; \ 619 } else { \ 620 PetscReal _2n = _n+_n; \ 621 PetscReal _d = (_2n*(_n+_a+_b)*(_2n+_a+_b-2)); \ 622 PetscReal _n1 = (_2n+_a+_b-1.)*(_a*_a-_b*_b); \ 623 PetscReal _n1x = (_2n+_a+_b-1.)*(_2n+_a+_b)*(_2n+_a+_b-2); \ 624 PetscReal _n2 = 2.*((_n+_a-1.)*(_n+_b-1.)*(_2n+_a+_b)); \ 625 (cnm1) = _n1 / _d; \ 626 (cnm1x) = _n1x / _d; \ 627 (cnm2) = _n2 / _d; \ 628 } \ 629 } while (0) 630 631 /*@ 632 PetscDTJacobiNorm - Compute the weighted L2 norm of a Jacobi polynomial. 633 634 $\| P^{\alpha,\beta}_n \|_{\alpha,\beta}^2 = \int_{-1}^1 (1 + x)^{\alpha} (1 - x)^{\beta} P^{\alpha,\beta}_n (x)^2 dx.$ 635 636 Input Parameters: 637 - alpha - the left exponent > -1 638 . beta - the right exponent > -1 639 + n - the polynomial degree 640 641 Output Parameter: 642 . norm - the weighted L2 norm 643 644 Level: beginner 645 646 .seealso: PetscDTJacobiEval() 647 @*/ 648 PetscErrorCode PetscDTJacobiNorm(PetscReal alpha, PetscReal beta, PetscInt n, PetscReal *norm) 649 { 650 PetscReal twoab1; 651 PetscReal gr; 652 653 PetscFunctionBegin; 654 PetscCheckFalse(alpha <= -1.,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exponent alpha %g <= -1. invalid", (double) alpha); 655 PetscCheckFalse(beta <= -1.,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Exponent beta %g <= -1. invalid", (double) beta); 656 PetscCheckFalse(n < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "n %D < 0 invalid", n); 657 twoab1 = PetscPowReal(2., alpha + beta + 1.); 658 #if defined(PETSC_HAVE_LGAMMA) 659 if (!n) { 660 gr = PetscExpReal(PetscLGamma(alpha+1.) + PetscLGamma(beta+1.) - PetscLGamma(alpha+beta+2.)); 661 } else { 662 gr = PetscExpReal(PetscLGamma(n+alpha+1.) + PetscLGamma(n+beta+1.) - (PetscLGamma(n+1.) + PetscLGamma(n+alpha+beta+1.))) / (n+n+alpha+beta+1.); 663 } 664 #else 665 { 666 PetscInt alphai = (PetscInt) alpha; 667 PetscInt betai = (PetscInt) beta; 668 PetscInt i; 669 670 gr = n ? (1. / (n+n+alpha+beta+1.)) : 1.; 671 if ((PetscReal) alphai == alpha) { 672 if (!n) { 673 for (i = 0; i < alphai; i++) gr *= (i+1.) / (beta+i+1.); 674 gr /= (alpha+beta+1.); 675 } else { 676 for (i = 0; i < alphai; i++) gr *= (n+i+1.) / (n+beta+i+1.); 677 } 678 } else if ((PetscReal) betai == beta) { 679 if (!n) { 680 for (i = 0; i < betai; i++) gr *= (i+1.) / (alpha+i+2.); 681 gr /= (alpha+beta+1.); 682 } else { 683 for (i = 0; i < betai; i++) gr *= (n+i+1.) / (n+alpha+i+1.); 684 } 685 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable."); 686 } 687 #endif 688 *norm = PetscSqrtReal(twoab1 * gr); 689 PetscFunctionReturn(0); 690 } 691 692 static PetscErrorCode PetscDTJacobiEval_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscInt k, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *p) 693 { 694 PetscReal ak, bk; 695 PetscReal abk1; 696 PetscInt i,l,maxdegree; 697 698 PetscFunctionBegin; 699 maxdegree = degrees[ndegree-1] - k; 700 ak = a + k; 701 bk = b + k; 702 abk1 = a + b + k + 1.; 703 if (maxdegree < 0) { 704 for (i = 0; i < npoints; i++) for (l = 0; l < ndegree; l++) p[i*ndegree+l] = 0.; 705 PetscFunctionReturn(0); 706 } 707 for (i=0; i<npoints; i++) { 708 PetscReal pm1,pm2,x; 709 PetscReal cnm1, cnm1x, cnm2; 710 PetscInt j,m; 711 712 x = points[i]; 713 pm2 = 1.; 714 PetscDTJacobiRecurrence_Internal(1,ak,bk,cnm1,cnm1x,cnm2); 715 pm1 = (cnm1 + cnm1x*x); 716 l = 0; 717 while (l < ndegree && degrees[l] - k < 0) { 718 p[l++] = 0.; 719 } 720 while (l < ndegree && degrees[l] - k == 0) { 721 p[l] = pm2; 722 for (m = 0; m < k; m++) p[l] *= (abk1 + m) * 0.5; 723 l++; 724 } 725 while (l < ndegree && degrees[l] - k == 1) { 726 p[l] = pm1; 727 for (m = 0; m < k; m++) p[l] *= (abk1 + 1 + m) * 0.5; 728 l++; 729 } 730 for (j=2; j<=maxdegree; j++) { 731 PetscReal pp; 732 733 PetscDTJacobiRecurrence_Internal(j,ak,bk,cnm1,cnm1x,cnm2); 734 pp = (cnm1 + cnm1x*x)*pm1 - cnm2*pm2; 735 pm2 = pm1; 736 pm1 = pp; 737 while (l < ndegree && degrees[l] - k == j) { 738 p[l] = pp; 739 for (m = 0; m < k; m++) p[l] *= (abk1 + j + m) * 0.5; 740 l++; 741 } 742 } 743 p += ndegree; 744 } 745 PetscFunctionReturn(0); 746 } 747 748 /*@ 749 PetscDTJacobiEvalJet - Evaluate the jet (function and derivatives) of the Jacobi polynomials basis up to a given degree. The Jacobi polynomials with indices $\alpha$ and $\beta$ are orthogonal with respect to the weighted inner product $\langle f, g \rangle = \int_{-1}^1 (1+x)^{\alpha} (1-x)^{\beta) f(x) g(x) dx$. 750 751 Input Parameters: 752 + alpha - the left exponent of the weight 753 . beta - the right exponetn of the weight 754 . npoints - the number of points to evaluate the polynomials at 755 . points - [npoints] array of point coordinates 756 . degree - the maximm degree polynomial space to evaluate, (degree + 1) will be evaluated total. 757 - k - the maximum derivative to evaluate in the jet, (k + 1) will be evaluated total. 758 759 Output Argments: 760 - p - an array containing the evaluations of the Jacobi polynomials's jets on the points. the size is (degree + 1) x 761 (k + 1) x npoints, which also describes the order of the dimensions of this three-dimensional array: the first 762 (slowest varying) dimension is polynomial degree; the second dimension is derivative order; the third (fastest 763 varying) dimension is the index of the evaluation point. 764 765 Level: advanced 766 767 .seealso: PetscDTJacobiEval(), PetscDTPKDEvalJet() 768 @*/ 769 PetscErrorCode PetscDTJacobiEvalJet(PetscReal alpha, PetscReal beta, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[]) 770 { 771 PetscInt i, j, l; 772 PetscInt *degrees; 773 PetscReal *psingle; 774 PetscErrorCode ierr; 775 776 PetscFunctionBegin; 777 if (degree == 0) { 778 PetscInt zero = 0; 779 780 for (i = 0; i <= k; i++) { 781 ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, i, points, 1, &zero, &p[i*npoints]);CHKERRQ(ierr); 782 } 783 PetscFunctionReturn(0); 784 } 785 ierr = PetscMalloc1(degree + 1, °rees);CHKERRQ(ierr); 786 ierr = PetscMalloc1((degree + 1) * npoints, &psingle);CHKERRQ(ierr); 787 for (i = 0; i <= degree; i++) degrees[i] = i; 788 for (i = 0; i <= k; i++) { 789 ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, i, points, degree + 1, degrees, psingle);CHKERRQ(ierr); 790 for (j = 0; j <= degree; j++) { 791 for (l = 0; l < npoints; l++) { 792 p[(j * (k + 1) + i) * npoints + l] = psingle[l * (degree + 1) + j]; 793 } 794 } 795 } 796 ierr = PetscFree(psingle);CHKERRQ(ierr); 797 ierr = PetscFree(degrees);CHKERRQ(ierr); 798 PetscFunctionReturn(0); 799 } 800 801 /*@ 802 PetscDTJacobiEval - evaluate Jacobi polynomials for the weight function $(1.+x)^{\alpha} (1.-x)^{\beta}$ 803 at points 804 805 Not Collective 806 807 Input Parameters: 808 + npoints - number of spatial points to evaluate at 809 . alpha - the left exponent > -1 810 . beta - the right exponent > -1 811 . points - array of locations to evaluate at 812 . ndegree - number of basis degrees to evaluate 813 - degrees - sorted array of degrees to evaluate 814 815 Output Parameters: 816 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 817 . D - row-oriented derivative evaluation matrix (or NULL) 818 - D2 - row-oriented second derivative evaluation matrix (or NULL) 819 820 Level: intermediate 821 822 .seealso: PetscDTGaussQuadrature() 823 @*/ 824 PetscErrorCode PetscDTJacobiEval(PetscInt npoints,PetscReal alpha, PetscReal beta, const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 825 { 826 PetscErrorCode ierr; 827 828 PetscFunctionBegin; 829 PetscCheckFalse(alpha <= -1.,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 830 PetscCheckFalse(beta <= -1.,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 831 if (!npoints || !ndegree) PetscFunctionReturn(0); 832 if (B) {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 0, points, ndegree, degrees, B);CHKERRQ(ierr);} 833 if (D) {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 1, points, ndegree, degrees, D);CHKERRQ(ierr);} 834 if (D2) {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 2, points, ndegree, degrees, D2);CHKERRQ(ierr);} 835 PetscFunctionReturn(0); 836 } 837 838 /*@ 839 PetscDTLegendreEval - evaluate Legendre polynomials at points 840 841 Not Collective 842 843 Input Parameters: 844 + npoints - number of spatial points to evaluate at 845 . points - array of locations to evaluate at 846 . ndegree - number of basis degrees to evaluate 847 - degrees - sorted array of degrees to evaluate 848 849 Output Parameters: 850 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 851 . D - row-oriented derivative evaluation matrix (or NULL) 852 - D2 - row-oriented second derivative evaluation matrix (or NULL) 853 854 Level: intermediate 855 856 .seealso: PetscDTGaussQuadrature() 857 @*/ 858 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 859 { 860 PetscErrorCode ierr; 861 862 PetscFunctionBegin; 863 ierr = PetscDTJacobiEval(npoints, 0., 0., points, ndegree, degrees, B, D, D2);CHKERRQ(ierr); 864 PetscFunctionReturn(0); 865 } 866 867 /*@ 868 PetscDTIndexToGradedOrder - convert an index into a tuple of monomial degrees in a graded order (that is, if the degree sum of tuple x is less than the degree sum of tuple y, then the index of x is smaller than the index of y) 869 870 Input Parameters: 871 + len - the desired length of the degree tuple 872 - index - the index to convert: should be >= 0 873 874 Output Parameter: 875 . degtup - will be filled with a tuple of degrees 876 877 Level: beginner 878 879 Note: for two tuples x and y with the same degree sum, partial degree sums over the final elements of the tuples 880 acts as a tiebreaker. For example, (2, 1, 1) and (1, 2, 1) have the same degree sum, but the degree sum over the 881 last two elements is smaller for the former, so (2, 1, 1) < (1, 2, 1). 882 883 .seealso: PetscDTGradedOrderToIndex() 884 @*/ 885 PetscErrorCode PetscDTIndexToGradedOrder(PetscInt len, PetscInt index, PetscInt degtup[]) 886 { 887 PetscInt i, total; 888 PetscInt sum; 889 890 PetscFunctionBeginHot; 891 PetscCheckFalse(len < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative"); 892 PetscCheckFalse(index < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index must be non-negative"); 893 total = 1; 894 sum = 0; 895 while (index >= total) { 896 index -= total; 897 total = (total * (len + sum)) / (sum + 1); 898 sum++; 899 } 900 for (i = 0; i < len; i++) { 901 PetscInt c; 902 903 degtup[i] = sum; 904 for (c = 0, total = 1; c < sum; c++) { 905 /* going into the loop, total is the number of way to have a tuple of sum exactly c with length len - 1 - i */ 906 if (index < total) break; 907 index -= total; 908 total = (total * (len - 1 - i + c)) / (c + 1); 909 degtup[i]--; 910 } 911 sum -= degtup[i]; 912 } 913 PetscFunctionReturn(0); 914 } 915 916 /*@ 917 PetscDTGradedOrderToIndex - convert a tuple into an index in a graded order, the inverse of PetscDTIndexToGradedOrder(). 918 919 Input Parameters: 920 + len - the length of the degree tuple 921 - degtup - tuple with this length 922 923 Output Parameter: 924 . index - index in graded order: >= 0 925 926 Level: Beginner 927 928 Note: for two tuples x and y with the same degree sum, partial degree sums over the final elements of the tuples 929 acts as a tiebreaker. For example, (2, 1, 1) and (1, 2, 1) have the same degree sum, but the degree sum over the 930 last two elements is smaller for the former, so (2, 1, 1) < (1, 2, 1). 931 932 .seealso: PetscDTIndexToGradedOrder() 933 @*/ 934 PetscErrorCode PetscDTGradedOrderToIndex(PetscInt len, const PetscInt degtup[], PetscInt *index) 935 { 936 PetscInt i, idx, sum, total; 937 938 PetscFunctionBeginHot; 939 PetscCheckFalse(len < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative"); 940 for (i = 0, sum = 0; i < len; i++) sum += degtup[i]; 941 idx = 0; 942 total = 1; 943 for (i = 0; i < sum; i++) { 944 idx += total; 945 total = (total * (len + i)) / (i + 1); 946 } 947 for (i = 0; i < len - 1; i++) { 948 PetscInt c; 949 950 total = 1; 951 sum -= degtup[i]; 952 for (c = 0; c < sum; c++) { 953 idx += total; 954 total = (total * (len - 1 - i + c)) / (c + 1); 955 } 956 } 957 *index = idx; 958 PetscFunctionReturn(0); 959 } 960 961 static PetscBool PKDCite = PETSC_FALSE; 962 const char PKDCitation[] = "@article{Kirby2010,\n" 963 " title={Singularity-free evaluation of collapsed-coordinate orthogonal polynomials},\n" 964 " author={Kirby, Robert C},\n" 965 " journal={ACM Transactions on Mathematical Software (TOMS)},\n" 966 " volume={37},\n" 967 " number={1},\n" 968 " pages={1--16},\n" 969 " year={2010},\n" 970 " publisher={ACM New York, NY, USA}\n}\n"; 971 972 /*@ 973 PetscDTPKDEvalJet - Evaluate the jet (function and derivatives) of the Proriol-Koornwinder-Dubiner (PKD) basis for 974 the space of polynomials up to a given degree. The PKD basis is L2-orthonormal on the biunit simplex (which is used 975 as the reference element for finite elements in PETSc), which makes it a stable basis to use for evaluating 976 polynomials in that domain. 977 978 Input Parameters: 979 + dim - the number of variables in the multivariate polynomials 980 . npoints - the number of points to evaluate the polynomials at 981 . points - [npoints x dim] array of point coordinates 982 . degree - the degree (sum of degrees on the variables in a monomial) of the polynomial space to evaluate. There are ((dim + degree) choose dim) polynomials in this space. 983 - k - the maximum order partial derivative to evaluate in the jet. There are (dim + k choose dim) partial derivatives 984 in the jet. Choosing k = 0 means to evaluate just the function and no derivatives 985 986 Output Argments: 987 - p - an array containing the evaluations of the PKD polynomials' jets on the points. The size is ((dim + degree) 988 choose dim) x ((dim + k) choose dim) x npoints, which also describes the order of the dimensions of this 989 three-dimensional array: the first (slowest varying) dimension is basis function index; the second dimension is jet 990 index; the third (fastest varying) dimension is the index of the evaluation point. 991 992 Level: advanced 993 994 Note: The ordering of the basis functions, and the ordering of the derivatives in the jet, both follow the graded 995 ordering of PetscDTIndexToGradedOrder() and PetscDTGradedOrderToIndex(). For example, in 3D, the polynomial with 996 leading monomial x^2,y^0,z^1, which has degree tuple (2,0,1), which by PetscDTGradedOrderToIndex() has index 12 (it is the 13th basis function in the space); 997 the partial derivative $\partial_x \partial_z$ has order tuple (1,0,1), appears at index 6 in the jet (it is the 7th partial derivative in the jet). 998 999 The implementation uses Kirby's singularity-free evaluation algorithm, https://doi.org/10.1145/1644001.1644006. 1000 1001 .seealso: PetscDTGradedOrderToIndex(), PetscDTIndexToGradedOrder(), PetscDTJacobiEvalJet() 1002 @*/ 1003 PetscErrorCode PetscDTPKDEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt k, PetscReal p[]) 1004 { 1005 PetscInt degidx, kidx, d, pt; 1006 PetscInt Nk, Ndeg; 1007 PetscInt *ktup, *degtup; 1008 PetscReal *scales, initscale, scaleexp; 1009 PetscErrorCode ierr; 1010 1011 PetscFunctionBegin; 1012 ierr = PetscCitationsRegister(PKDCitation, &PKDCite);CHKERRQ(ierr); 1013 ierr = PetscDTBinomialInt(dim + k, k, &Nk);CHKERRQ(ierr); 1014 ierr = PetscDTBinomialInt(degree + dim, degree, &Ndeg);CHKERRQ(ierr); 1015 ierr = PetscMalloc2(dim, °tup, dim, &ktup);CHKERRQ(ierr); 1016 ierr = PetscMalloc1(Ndeg, &scales);CHKERRQ(ierr); 1017 initscale = 1.; 1018 if (dim > 1) { 1019 ierr = PetscDTBinomial(dim,2,&scaleexp);CHKERRQ(ierr); 1020 initscale = PetscPowReal(2.,scaleexp*0.5); 1021 } 1022 for (degidx = 0; degidx < Ndeg; degidx++) { 1023 PetscInt e, i; 1024 PetscInt m1idx = -1, m2idx = -1; 1025 PetscInt n; 1026 PetscInt degsum; 1027 PetscReal alpha; 1028 PetscReal cnm1, cnm1x, cnm2; 1029 PetscReal norm; 1030 1031 ierr = PetscDTIndexToGradedOrder(dim, degidx, degtup);CHKERRQ(ierr); 1032 for (d = dim - 1; d >= 0; d--) if (degtup[d]) break; 1033 if (d < 0) { /* constant is 1 everywhere, all derivatives are zero */ 1034 scales[degidx] = initscale; 1035 for (e = 0; e < dim; e++) { 1036 ierr = PetscDTJacobiNorm(e,0.,0,&norm);CHKERRQ(ierr); 1037 scales[degidx] /= norm; 1038 } 1039 for (i = 0; i < npoints; i++) p[degidx * Nk * npoints + i] = 1.; 1040 for (i = 0; i < (Nk - 1) * npoints; i++) p[(degidx * Nk + 1) * npoints + i] = 0.; 1041 continue; 1042 } 1043 n = degtup[d]; 1044 degtup[d]--; 1045 ierr = PetscDTGradedOrderToIndex(dim, degtup, &m1idx);CHKERRQ(ierr); 1046 if (degtup[d] > 0) { 1047 degtup[d]--; 1048 ierr = PetscDTGradedOrderToIndex(dim, degtup, &m2idx);CHKERRQ(ierr); 1049 degtup[d]++; 1050 } 1051 degtup[d]++; 1052 for (e = 0, degsum = 0; e < d; e++) degsum += degtup[e]; 1053 alpha = 2 * degsum + d; 1054 PetscDTJacobiRecurrence_Internal(n,alpha,0.,cnm1,cnm1x,cnm2); 1055 1056 scales[degidx] = initscale; 1057 for (e = 0, degsum = 0; e < dim; e++) { 1058 PetscInt f; 1059 PetscReal ealpha; 1060 PetscReal enorm; 1061 1062 ealpha = 2 * degsum + e; 1063 for (f = 0; f < degsum; f++) scales[degidx] *= 2.; 1064 ierr = PetscDTJacobiNorm(ealpha,0.,degtup[e],&enorm);CHKERRQ(ierr); 1065 scales[degidx] /= enorm; 1066 degsum += degtup[e]; 1067 } 1068 1069 for (pt = 0; pt < npoints; pt++) { 1070 /* compute the multipliers */ 1071 PetscReal thetanm1, thetanm1x, thetanm2; 1072 1073 thetanm1x = dim - (d+1) + 2.*points[pt * dim + d]; 1074 for (e = d+1; e < dim; e++) thetanm1x += points[pt * dim + e]; 1075 thetanm1x *= 0.5; 1076 thetanm1 = (2. - (dim-(d+1))); 1077 for (e = d+1; e < dim; e++) thetanm1 -= points[pt * dim + e]; 1078 thetanm1 *= 0.5; 1079 thetanm2 = thetanm1 * thetanm1; 1080 1081 for (kidx = 0; kidx < Nk; kidx++) { 1082 PetscInt f; 1083 1084 ierr = PetscDTIndexToGradedOrder(dim, kidx, ktup);CHKERRQ(ierr); 1085 /* first sum in the same derivative terms */ 1086 p[(degidx * Nk + kidx) * npoints + pt] = (cnm1 * thetanm1 + cnm1x * thetanm1x) * p[(m1idx * Nk + kidx) * npoints + pt]; 1087 if (m2idx >= 0) { 1088 p[(degidx * Nk + kidx) * npoints + pt] -= cnm2 * thetanm2 * p[(m2idx * Nk + kidx) * npoints + pt]; 1089 } 1090 1091 for (f = d; f < dim; f++) { 1092 PetscInt km1idx, mplty = ktup[f]; 1093 1094 if (!mplty) continue; 1095 ktup[f]--; 1096 ierr = PetscDTGradedOrderToIndex(dim, ktup, &km1idx);CHKERRQ(ierr); 1097 1098 /* the derivative of cnm1x * thetanm1x wrt x variable f is 0.5 * cnm1x if f > d otherwise it is cnm1x */ 1099 /* the derivative of cnm1 * thetanm1 wrt x variable f is 0 if f == d, otherwise it is -0.5 * cnm1 */ 1100 /* the derivative of -cnm2 * thetanm2 wrt x variable f is 0 if f == d, otherwise it is cnm2 * thetanm1 */ 1101 if (f > d) { 1102 PetscInt f2; 1103 1104 p[(degidx * Nk + kidx) * npoints + pt] += mplty * 0.5 * (cnm1x - cnm1) * p[(m1idx * Nk + km1idx) * npoints + pt]; 1105 if (m2idx >= 0) { 1106 p[(degidx * Nk + kidx) * npoints + pt] += mplty * cnm2 * thetanm1 * p[(m2idx * Nk + km1idx) * npoints + pt]; 1107 /* second derivatives of -cnm2 * thetanm2 wrt x variable f,f2 is like - 0.5 * cnm2 */ 1108 for (f2 = f; f2 < dim; f2++) { 1109 PetscInt km2idx, mplty2 = ktup[f2]; 1110 PetscInt factor; 1111 1112 if (!mplty2) continue; 1113 ktup[f2]--; 1114 ierr = PetscDTGradedOrderToIndex(dim, ktup, &km2idx);CHKERRQ(ierr); 1115 1116 factor = mplty * mplty2; 1117 if (f == f2) factor /= 2; 1118 p[(degidx * Nk + kidx) * npoints + pt] -= 0.5 * factor * cnm2 * p[(m2idx * Nk + km2idx) * npoints + pt]; 1119 ktup[f2]++; 1120 } 1121 } 1122 } else { 1123 p[(degidx * Nk + kidx) * npoints + pt] += mplty * cnm1x * p[(m1idx * Nk + km1idx) * npoints + pt]; 1124 } 1125 ktup[f]++; 1126 } 1127 } 1128 } 1129 } 1130 for (degidx = 0; degidx < Ndeg; degidx++) { 1131 PetscReal scale = scales[degidx]; 1132 PetscInt i; 1133 1134 for (i = 0; i < Nk * npoints; i++) p[degidx*Nk*npoints + i] *= scale; 1135 } 1136 ierr = PetscFree(scales);CHKERRQ(ierr); 1137 ierr = PetscFree2(degtup, ktup);CHKERRQ(ierr); 1138 PetscFunctionReturn(0); 1139 } 1140 1141 /*@ 1142 PetscDTPTrimmedSize - The size of the trimmed polynomial space of k-forms with a given degree and form degree, 1143 which can be evaluated in PetscDTPTrimmedEvalJet(). 1144 1145 Input Parameters: 1146 + dim - the number of variables in the multivariate polynomials 1147 . degree - the degree (sum of degrees on the variables in a monomial) of the trimmed polynomial space. 1148 - formDegree - the degree of the form 1149 1150 Output Argments: 1151 - size - The number ((dim + degree) choose (dim + formDegree)) x ((degree + formDegree - 1) choose (formDegree)) 1152 1153 Level: advanced 1154 1155 .seealso: PetscDTPTrimmedEvalJet() 1156 @*/ 1157 PetscErrorCode PetscDTPTrimmedSize(PetscInt dim, PetscInt degree, PetscInt formDegree, PetscInt *size) 1158 { 1159 PetscInt Nrk, Nbpt; // number of trimmed polynomials 1160 PetscErrorCode ierr; 1161 1162 PetscFunctionBegin; 1163 formDegree = PetscAbsInt(formDegree); 1164 ierr = PetscDTBinomialInt(degree + dim, degree + formDegree, &Nbpt);CHKERRQ(ierr); 1165 ierr = PetscDTBinomialInt(degree + formDegree - 1, formDegree, &Nrk);CHKERRQ(ierr); 1166 Nbpt *= Nrk; 1167 *size = Nbpt; 1168 PetscFunctionReturn(0); 1169 } 1170 1171 /* there was a reference implementation based on section 4.4 of Arnold, Falk & Winther (acta numerica, 2006), but it 1172 * was inferior to this implementation */ 1173 static PetscErrorCode PetscDTPTrimmedEvalJet_Internal(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt formDegree, PetscInt jetDegree, PetscReal p[]) 1174 { 1175 PetscInt formDegreeOrig = formDegree; 1176 PetscBool formNegative = (formDegreeOrig < 0) ? PETSC_TRUE : PETSC_FALSE; 1177 PetscErrorCode ierr; 1178 1179 PetscFunctionBegin; 1180 formDegree = PetscAbsInt(formDegreeOrig); 1181 if (formDegree == 0) { 1182 ierr = PetscDTPKDEvalJet(dim, npoints, points, degree, jetDegree, p);CHKERRQ(ierr); 1183 PetscFunctionReturn(0); 1184 } 1185 if (formDegree == dim) { 1186 ierr = PetscDTPKDEvalJet(dim, npoints, points, degree - 1, jetDegree, p);CHKERRQ(ierr); 1187 PetscFunctionReturn(0); 1188 } 1189 PetscInt Nbpt; 1190 ierr = PetscDTPTrimmedSize(dim, degree, formDegree, &Nbpt);CHKERRQ(ierr); 1191 PetscInt Nf; 1192 ierr = PetscDTBinomialInt(dim, formDegree, &Nf);CHKERRQ(ierr); 1193 PetscInt Nk; 1194 ierr = PetscDTBinomialInt(dim + jetDegree, dim, &Nk);CHKERRQ(ierr); 1195 ierr = PetscArrayzero(p, Nbpt * Nf * Nk * npoints);CHKERRQ(ierr); 1196 1197 PetscInt Nbpm1; // number of scalar polynomials up to degree - 1; 1198 ierr = PetscDTBinomialInt(dim + degree - 1, dim, &Nbpm1);CHKERRQ(ierr); 1199 PetscReal *p_scalar; 1200 ierr = PetscMalloc1(Nbpm1 * Nk * npoints, &p_scalar);CHKERRQ(ierr); 1201 ierr = PetscDTPKDEvalJet(dim, npoints, points, degree - 1, jetDegree, p_scalar);CHKERRQ(ierr); 1202 PetscInt total = 0; 1203 // First add the full polynomials up to degree - 1 into the basis: take the scalar 1204 // and copy one for each form component 1205 for (PetscInt i = 0; i < Nbpm1; i++) { 1206 const PetscReal *src = &p_scalar[i * Nk * npoints]; 1207 for (PetscInt f = 0; f < Nf; f++) { 1208 PetscReal *dest = &p[(total++ * Nf + f) * Nk * npoints]; 1209 ierr = PetscArraycpy(dest, src, Nk * npoints);CHKERRQ(ierr); 1210 } 1211 } 1212 PetscInt *form_atoms; 1213 ierr = PetscMalloc1(formDegree + 1, &form_atoms);CHKERRQ(ierr); 1214 // construct the interior product pattern 1215 PetscInt (*pattern)[3]; 1216 PetscInt Nf1; // number of formDegree + 1 forms 1217 ierr = PetscDTBinomialInt(dim, formDegree + 1, &Nf1);CHKERRQ(ierr); 1218 PetscInt nnz = Nf1 * (formDegree+1); 1219 ierr = PetscMalloc1(Nf1 * (formDegree+1), &pattern);CHKERRQ(ierr); 1220 ierr = PetscDTAltVInteriorPattern(dim, formDegree+1, pattern);CHKERRQ(ierr); 1221 PetscReal centroid = (1. - dim) / (dim + 1.); 1222 PetscInt *deriv; 1223 ierr = PetscMalloc1(dim, &deriv);CHKERRQ(ierr); 1224 for (PetscInt d = dim; d >= formDegree + 1; d--) { 1225 PetscInt Nfd1; // number of formDegree + 1 forms in dimension d that include dx_0 1226 // (equal to the number of formDegree forms in dimension d-1) 1227 ierr = PetscDTBinomialInt(d - 1, formDegree, &Nfd1);CHKERRQ(ierr); 1228 // The number of homogeneous (degree-1) scalar polynomials in d variables 1229 PetscInt Nh; 1230 ierr = PetscDTBinomialInt(d - 1 + degree - 1, d - 1, &Nh);CHKERRQ(ierr); 1231 const PetscReal *h_scalar = &p_scalar[(Nbpm1 - Nh) * Nk * npoints]; 1232 for (PetscInt b = 0; b < Nh; b++) { 1233 const PetscReal *h_s = &h_scalar[b * Nk * npoints]; 1234 for (PetscInt f = 0; f < Nfd1; f++) { 1235 // construct all formDegree+1 forms that start with dx_(dim - d) /\ ... 1236 form_atoms[0] = dim - d; 1237 ierr = PetscDTEnumSubset(d-1, formDegree, f, &form_atoms[1]);CHKERRQ(ierr); 1238 for (PetscInt i = 0; i < formDegree; i++) { 1239 form_atoms[1+i] += form_atoms[0] + 1; 1240 } 1241 PetscInt f_ind; // index of the resulting form 1242 ierr = PetscDTSubsetIndex(dim, formDegree + 1, form_atoms, &f_ind);CHKERRQ(ierr); 1243 PetscReal *p_f = &p[total++ * Nf * Nk * npoints]; 1244 for (PetscInt nz = 0; nz < nnz; nz++) { 1245 PetscInt i = pattern[nz][0]; // formDegree component 1246 PetscInt j = pattern[nz][1]; // (formDegree + 1) component 1247 PetscInt v = pattern[nz][2]; // coordinate component 1248 PetscReal scale = v < 0 ? -1. : 1.; 1249 1250 i = formNegative ? (Nf - 1 - i) : i; 1251 scale = (formNegative && (i & 1)) ? -scale : scale; 1252 v = v < 0 ? -(v + 1) : v; 1253 if (j != f_ind) { 1254 continue; 1255 } 1256 PetscReal *p_i = &p_f[i * Nk * npoints]; 1257 for (PetscInt jet = 0; jet < Nk; jet++) { 1258 const PetscReal *h_jet = &h_s[jet * npoints]; 1259 PetscReal *p_jet = &p_i[jet * npoints]; 1260 1261 for (PetscInt pt = 0; pt < npoints; pt++) { 1262 p_jet[pt] += scale * h_jet[pt] * (points[pt * dim + v] - centroid); 1263 } 1264 ierr = PetscDTIndexToGradedOrder(dim, jet, deriv);CHKERRQ(ierr); 1265 deriv[v]++; 1266 PetscReal mult = deriv[v]; 1267 PetscInt l; 1268 ierr = PetscDTGradedOrderToIndex(dim, deriv, &l);CHKERRQ(ierr); 1269 if (l >= Nk) { 1270 continue; 1271 } 1272 p_jet = &p_i[l * npoints]; 1273 for (PetscInt pt = 0; pt < npoints; pt++) { 1274 p_jet[pt] += scale * mult * h_jet[pt]; 1275 } 1276 deriv[v]--; 1277 } 1278 } 1279 } 1280 } 1281 } 1282 PetscCheckFalse(total != Nbpt,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrectly counted P trimmed polynomials"); 1283 ierr = PetscFree(deriv);CHKERRQ(ierr); 1284 ierr = PetscFree(pattern);CHKERRQ(ierr); 1285 ierr = PetscFree(form_atoms);CHKERRQ(ierr); 1286 ierr = PetscFree(p_scalar);CHKERRQ(ierr); 1287 PetscFunctionReturn(0); 1288 } 1289 1290 /*@ 1291 PetscDTPTrimmedEvalJet - Evaluate the jet (function and derivatives) of a basis of the trimmed polynomial k-forms up to 1292 a given degree. 1293 1294 Input Parameters: 1295 + dim - the number of variables in the multivariate polynomials 1296 . npoints - the number of points to evaluate the polynomials at 1297 . points - [npoints x dim] array of point coordinates 1298 . degree - the degree (sum of degrees on the variables in a monomial) of the trimmed polynomial space to evaluate. 1299 There are ((dim + degree) choose (dim + formDegree)) x ((degree + formDegree - 1) choose (formDegree)) polynomials in this space. 1300 (You can use PetscDTPTrimmedSize() to compute this size.) 1301 . formDegree - the degree of the form 1302 - jetDegree - the maximum order partial derivative to evaluate in the jet. There are ((dim + jetDegree) choose dim) partial derivatives 1303 in the jet. Choosing jetDegree = 0 means to evaluate just the function and no derivatives 1304 1305 Output Argments: 1306 - p - an array containing the evaluations of the PKD polynomials' jets on the points. The size is 1307 PetscDTPTrimmedSize() x ((dim + formDegree) choose dim) x ((dim + k) choose dim) x npoints, 1308 which also describes the order of the dimensions of this 1309 four-dimensional array: 1310 the first (slowest varying) dimension is basis function index; 1311 the second dimension is component of the form; 1312 the third dimension is jet index; 1313 the fourth (fastest varying) dimension is the index of the evaluation point. 1314 1315 Level: advanced 1316 1317 Note: The ordering of the basis functions is not graded, so the basis functions are not nested by degree like PetscDTPKDEvalJet(). 1318 The basis functions are not an L2-orthonormal basis on any particular domain. 1319 1320 The implementation is based on the description of the trimmed polynomials up to degree r as 1321 the direct sum of polynomials up to degree (r-1) and the Koszul differential applied to 1322 homogeneous polynomials of degree (r-1). 1323 1324 .seealso: PetscDTPKDEvalJet(), PetscDTPTrimmedSize() 1325 @*/ 1326 PetscErrorCode PetscDTPTrimmedEvalJet(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt degree, PetscInt formDegree, PetscInt jetDegree, PetscReal p[]) 1327 { 1328 PetscErrorCode ierr; 1329 1330 PetscFunctionBegin; 1331 ierr = PetscDTPTrimmedEvalJet_Internal(dim, npoints, points, degree, formDegree, jetDegree, p);CHKERRQ(ierr); 1332 PetscFunctionReturn(0); 1333 } 1334 1335 /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V 1336 * with lds n; diag and subdiag are overwritten */ 1337 static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[], 1338 PetscReal eigs[], PetscScalar V[]) 1339 { 1340 char jobz = 'V'; /* eigenvalues and eigenvectors */ 1341 char range = 'A'; /* all eigenvalues will be found */ 1342 PetscReal VL = 0.; /* ignored because range is 'A' */ 1343 PetscReal VU = 0.; /* ignored because range is 'A' */ 1344 PetscBLASInt IL = 0; /* ignored because range is 'A' */ 1345 PetscBLASInt IU = 0; /* ignored because range is 'A' */ 1346 PetscReal abstol = 0.; /* unused */ 1347 PetscBLASInt bn, bm, ldz; /* bm will equal bn on exit */ 1348 PetscBLASInt *isuppz; 1349 PetscBLASInt lwork, liwork; 1350 PetscReal workquery; 1351 PetscBLASInt iworkquery; 1352 PetscBLASInt *iwork; 1353 PetscBLASInt info; 1354 PetscReal *work = NULL; 1355 PetscErrorCode ierr; 1356 1357 PetscFunctionBegin; 1358 #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG) 1359 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found"); 1360 #endif 1361 ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr); 1362 ierr = PetscBLASIntCast(n, &ldz);CHKERRQ(ierr); 1363 #if !defined(PETSC_MISSING_LAPACK_STEGR) 1364 ierr = PetscMalloc1(2 * n, &isuppz);CHKERRQ(ierr); 1365 lwork = -1; 1366 liwork = -1; 1367 PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,&workquery,&lwork,&iworkquery,&liwork,&info)); 1368 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error"); 1369 lwork = (PetscBLASInt) workquery; 1370 liwork = (PetscBLASInt) iworkquery; 1371 ierr = PetscMalloc2(lwork, &work, liwork, &iwork);CHKERRQ(ierr); 1372 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1373 PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,work,&lwork,iwork,&liwork,&info)); 1374 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1375 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error"); 1376 ierr = PetscFree2(work, iwork);CHKERRQ(ierr); 1377 ierr = PetscFree(isuppz);CHKERRQ(ierr); 1378 #elif !defined(PETSC_MISSING_LAPACK_STEQR) 1379 jobz = 'I'; /* Compute eigenvalues and eigenvectors of the 1380 tridiagonal matrix. Z is initialized to the identity 1381 matrix. */ 1382 ierr = PetscMalloc1(PetscMax(1,2*n-2),&work);CHKERRQ(ierr); 1383 PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&bn,diag,subdiag,V,&ldz,work,&info)); 1384 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1385 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 1386 ierr = PetscFree(work);CHKERRQ(ierr); 1387 ierr = PetscArraycpy(eigs,diag,n);CHKERRQ(ierr); 1388 #endif 1389 PetscFunctionReturn(0); 1390 } 1391 1392 /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi 1393 * quadrature rules on the interval [-1, 1] */ 1394 static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw) 1395 { 1396 PetscReal twoab1; 1397 PetscInt m = n - 2; 1398 PetscReal a = alpha + 1.; 1399 PetscReal b = beta + 1.; 1400 PetscReal gra, grb; 1401 1402 PetscFunctionBegin; 1403 twoab1 = PetscPowReal(2., a + b - 1.); 1404 #if defined(PETSC_HAVE_LGAMMA) 1405 grb = PetscExpReal(2. * PetscLGamma(b+1.) + PetscLGamma(m+1.) + PetscLGamma(m+a+1.) - 1406 (PetscLGamma(m+b+1) + PetscLGamma(m+a+b+1.))); 1407 gra = PetscExpReal(2. * PetscLGamma(a+1.) + PetscLGamma(m+1.) + PetscLGamma(m+b+1.) - 1408 (PetscLGamma(m+a+1) + PetscLGamma(m+a+b+1.))); 1409 #else 1410 { 1411 PetscInt alphai = (PetscInt) alpha; 1412 PetscInt betai = (PetscInt) beta; 1413 PetscErrorCode ierr; 1414 1415 if ((PetscReal) alphai == alpha && (PetscReal) betai == beta) { 1416 PetscReal binom1, binom2; 1417 1418 ierr = PetscDTBinomial(m+b, b, &binom1);CHKERRQ(ierr); 1419 ierr = PetscDTBinomial(m+a+b, b, &binom2);CHKERRQ(ierr); 1420 grb = 1./ (binom1 * binom2); 1421 ierr = PetscDTBinomial(m+a, a, &binom1);CHKERRQ(ierr); 1422 ierr = PetscDTBinomial(m+a+b, a, &binom2);CHKERRQ(ierr); 1423 gra = 1./ (binom1 * binom2); 1424 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable."); 1425 } 1426 #endif 1427 *leftw = twoab1 * grb / b; 1428 *rightw = twoab1 * gra / a; 1429 PetscFunctionReturn(0); 1430 } 1431 1432 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 1433 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 1434 static inline PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 1435 { 1436 PetscReal pn1, pn2; 1437 PetscReal cnm1, cnm1x, cnm2; 1438 PetscInt k; 1439 1440 PetscFunctionBegin; 1441 if (!n) {*P = 1.0; PetscFunctionReturn(0);} 1442 PetscDTJacobiRecurrence_Internal(1,a,b,cnm1,cnm1x,cnm2); 1443 pn2 = 1.; 1444 pn1 = cnm1 + cnm1x*x; 1445 if (n == 1) {*P = pn1; PetscFunctionReturn(0);} 1446 *P = 0.0; 1447 for (k = 2; k < n+1; ++k) { 1448 PetscDTJacobiRecurrence_Internal(k,a,b,cnm1,cnm1x,cnm2); 1449 1450 *P = (cnm1 + cnm1x*x)*pn1 - cnm2*pn2; 1451 pn2 = pn1; 1452 pn1 = *P; 1453 } 1454 PetscFunctionReturn(0); 1455 } 1456 1457 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 1458 static inline PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P) 1459 { 1460 PetscReal nP; 1461 PetscInt i; 1462 PetscErrorCode ierr; 1463 1464 PetscFunctionBegin; 1465 *P = 0.0; 1466 if (k > n) PetscFunctionReturn(0); 1467 ierr = PetscDTComputeJacobi(a+k, b+k, n-k, x, &nP);CHKERRQ(ierr); 1468 for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5; 1469 *P = nP; 1470 PetscFunctionReturn(0); 1471 } 1472 1473 static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[]) 1474 { 1475 PetscInt maxIter = 100; 1476 PetscReal eps = PetscExpReal(0.75 * PetscLogReal(PETSC_MACHINE_EPSILON)); 1477 PetscReal a1, a6, gf; 1478 PetscInt k; 1479 PetscErrorCode ierr; 1480 1481 PetscFunctionBegin; 1482 1483 a1 = PetscPowReal(2.0, a+b+1); 1484 #if defined(PETSC_HAVE_LGAMMA) 1485 { 1486 PetscReal a2, a3, a4, a5; 1487 a2 = PetscLGamma(a + npoints + 1); 1488 a3 = PetscLGamma(b + npoints + 1); 1489 a4 = PetscLGamma(a + b + npoints + 1); 1490 a5 = PetscLGamma(npoints + 1); 1491 gf = PetscExpReal(a2 + a3 - (a4 + a5)); 1492 } 1493 #else 1494 { 1495 PetscInt ia, ib; 1496 1497 ia = (PetscInt) a; 1498 ib = (PetscInt) b; 1499 gf = 1.; 1500 if (ia == a && ia >= 0) { /* compute ratio of rising factorals wrt a */ 1501 for (k = 0; k < ia; k++) gf *= (npoints + 1. + k) / (npoints + b + 1. + k); 1502 } else if (b == b && ib >= 0) { /* compute ratio of rising factorials wrt b */ 1503 for (k = 0; k < ib; k++) gf *= (npoints + 1. + k) / (npoints + a + 1. + k); 1504 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable."); 1505 } 1506 #endif 1507 1508 a6 = a1 * gf; 1509 /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 1510 Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 1511 for (k = 0; k < npoints; ++k) { 1512 PetscReal r = PetscCosReal(PETSC_PI * (1. - (4.*k + 3. + 2.*b) / (4.*npoints + 2.*(a + b + 1.)))), dP; 1513 PetscInt j; 1514 1515 if (k > 0) r = 0.5 * (r + x[k-1]); 1516 for (j = 0; j < maxIter; ++j) { 1517 PetscReal s = 0.0, delta, f, fp; 1518 PetscInt i; 1519 1520 for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 1521 ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 1522 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp);CHKERRQ(ierr); 1523 delta = f / (fp - f * s); 1524 r = r - delta; 1525 if (PetscAbsReal(delta) < eps) break; 1526 } 1527 x[k] = r; 1528 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP);CHKERRQ(ierr); 1529 w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 1530 } 1531 PetscFunctionReturn(0); 1532 } 1533 1534 /* Compute the diagonals of the Jacobi matrix used in Golub & Welsch algorithms for Gauss-Jacobi 1535 * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */ 1536 static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s) 1537 { 1538 PetscInt i; 1539 1540 PetscFunctionBegin; 1541 for (i = 0; i < nPoints; i++) { 1542 PetscReal A, B, C; 1543 1544 PetscDTJacobiRecurrence_Internal(i+1,a,b,A,B,C); 1545 d[i] = -A / B; 1546 if (i) s[i-1] *= C / B; 1547 if (i < nPoints - 1) s[i] = 1. / B; 1548 } 1549 PetscFunctionReturn(0); 1550 } 1551 1552 static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 1553 { 1554 PetscReal mu0; 1555 PetscReal ga, gb, gab; 1556 PetscInt i; 1557 PetscErrorCode ierr; 1558 1559 PetscFunctionBegin; 1560 ierr = PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite);CHKERRQ(ierr); 1561 1562 #if defined(PETSC_HAVE_TGAMMA) 1563 ga = PetscTGamma(a + 1); 1564 gb = PetscTGamma(b + 1); 1565 gab = PetscTGamma(a + b + 2); 1566 #else 1567 { 1568 PetscInt ia, ib; 1569 1570 ia = (PetscInt) a; 1571 ib = (PetscInt) b; 1572 if (ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0) { /* All gamma(x) terms are (x-1)! terms */ 1573 ierr = PetscDTFactorial(ia, &ga);CHKERRQ(ierr); 1574 ierr = PetscDTFactorial(ib, &gb);CHKERRQ(ierr); 1575 ierr = PetscDTFactorial(ia + ib + 1, &gb);CHKERRQ(ierr); 1576 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 1577 } 1578 #endif 1579 mu0 = PetscPowReal(2.,a + b + 1.) * ga * gb / gab; 1580 1581 #if defined(PETSCDTGAUSSIANQUADRATURE_EIG) 1582 { 1583 PetscReal *diag, *subdiag; 1584 PetscScalar *V; 1585 1586 ierr = PetscMalloc2(npoints, &diag, npoints, &subdiag);CHKERRQ(ierr); 1587 ierr = PetscMalloc1(npoints*npoints, &V);CHKERRQ(ierr); 1588 ierr = PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag);CHKERRQ(ierr); 1589 for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]); 1590 ierr = PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V);CHKERRQ(ierr); 1591 for (i = 0; i < npoints; i++) w[i] = PetscSqr(PetscRealPart(V[i * npoints])) * mu0; 1592 ierr = PetscFree(V);CHKERRQ(ierr); 1593 ierr = PetscFree2(diag, subdiag);CHKERRQ(ierr); 1594 } 1595 #else 1596 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found"); 1597 #endif 1598 { /* As of March 2, 2020, The Sun Performance Library breaks the LAPACK contract for xstegr and xsteqr: the 1599 eigenvalues are not guaranteed to be in ascending order. So we heave a passive aggressive sigh and check that 1600 the eigenvalues are sorted */ 1601 PetscBool sorted; 1602 1603 ierr = PetscSortedReal(npoints, x, &sorted);CHKERRQ(ierr); 1604 if (!sorted) { 1605 PetscInt *order, i; 1606 PetscReal *tmp; 1607 1608 ierr = PetscMalloc2(npoints, &order, npoints, &tmp);CHKERRQ(ierr); 1609 for (i = 0; i < npoints; i++) order[i] = i; 1610 ierr = PetscSortRealWithPermutation(npoints, x, order);CHKERRQ(ierr); 1611 ierr = PetscArraycpy(tmp, x, npoints);CHKERRQ(ierr); 1612 for (i = 0; i < npoints; i++) x[i] = tmp[order[i]]; 1613 ierr = PetscArraycpy(tmp, w, npoints);CHKERRQ(ierr); 1614 for (i = 0; i < npoints; i++) w[i] = tmp[order[i]]; 1615 ierr = PetscFree2(order, tmp);CHKERRQ(ierr); 1616 } 1617 } 1618 PetscFunctionReturn(0); 1619 } 1620 1621 static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton) 1622 { 1623 PetscErrorCode ierr; 1624 1625 PetscFunctionBegin; 1626 PetscCheckFalse(npoints < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive"); 1627 /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */ 1628 PetscCheckFalse(alpha <= -1.,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 1629 PetscCheckFalse(beta <= -1.,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 1630 1631 if (newton) { 1632 ierr = PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr); 1633 } else { 1634 ierr = PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr); 1635 } 1636 if (alpha == beta) { /* symmetrize */ 1637 PetscInt i; 1638 for (i = 0; i < (npoints + 1) / 2; i++) { 1639 PetscInt j = npoints - 1 - i; 1640 PetscReal xi = x[i]; 1641 PetscReal xj = x[j]; 1642 PetscReal wi = w[i]; 1643 PetscReal wj = w[j]; 1644 1645 x[i] = (xi - xj) / 2.; 1646 x[j] = (xj - xi) / 2.; 1647 w[i] = w[j] = (wi + wj) / 2.; 1648 } 1649 } 1650 PetscFunctionReturn(0); 1651 } 1652 1653 /*@ 1654 PetscDTGaussJacobiQuadrature - quadrature for the interval [a, b] with the weight function 1655 $(x-a)^\alpha (x-b)^\beta$. 1656 1657 Not collective 1658 1659 Input Parameters: 1660 + npoints - the number of points in the quadrature rule 1661 . a - the left endpoint of the interval 1662 . b - the right endpoint of the interval 1663 . alpha - the left exponent 1664 - beta - the right exponent 1665 1666 Output Parameters: 1667 + x - array of length npoints, the locations of the quadrature points 1668 - w - array of length npoints, the weights of the quadrature points 1669 1670 Level: intermediate 1671 1672 Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 1. 1673 @*/ 1674 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[]) 1675 { 1676 PetscInt i; 1677 PetscErrorCode ierr; 1678 1679 PetscFunctionBegin; 1680 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr); 1681 if (a != -1. || b != 1.) { /* shift */ 1682 for (i = 0; i < npoints; i++) { 1683 x[i] = (x[i] + 1.) * ((b - a) / 2.) + a; 1684 w[i] *= (b - a) / 2.; 1685 } 1686 } 1687 PetscFunctionReturn(0); 1688 } 1689 1690 static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton) 1691 { 1692 PetscInt i; 1693 PetscErrorCode ierr; 1694 1695 PetscFunctionBegin; 1696 PetscCheckFalse(npoints < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive"); 1697 /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */ 1698 PetscCheckFalse(alpha <= -1.,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 1699 PetscCheckFalse(beta <= -1.,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 1700 1701 x[0] = -1.; 1702 x[npoints-1] = 1.; 1703 if (npoints > 2) { 1704 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints-2, alpha+1., beta+1., &x[1], &w[1], newton);CHKERRQ(ierr); 1705 } 1706 for (i = 1; i < npoints - 1; i++) { 1707 w[i] /= (1. - x[i]*x[i]); 1708 } 1709 ierr = PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints-1]);CHKERRQ(ierr); 1710 PetscFunctionReturn(0); 1711 } 1712 1713 /*@ 1714 PetscDTGaussLobattoJacobiQuadrature - quadrature for the interval [a, b] with the weight function 1715 $(x-a)^\alpha (x-b)^\beta$, with endpoints a and b included as quadrature points. 1716 1717 Not collective 1718 1719 Input Parameters: 1720 + npoints - the number of points in the quadrature rule 1721 . a - the left endpoint of the interval 1722 . b - the right endpoint of the interval 1723 . alpha - the left exponent 1724 - beta - the right exponent 1725 1726 Output Parameters: 1727 + x - array of length npoints, the locations of the quadrature points 1728 - w - array of length npoints, the weights of the quadrature points 1729 1730 Level: intermediate 1731 1732 Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 3. 1733 @*/ 1734 PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[]) 1735 { 1736 PetscInt i; 1737 PetscErrorCode ierr; 1738 1739 PetscFunctionBegin; 1740 ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr); 1741 if (a != -1. || b != 1.) { /* shift */ 1742 for (i = 0; i < npoints; i++) { 1743 x[i] = (x[i] + 1.) * ((b - a) / 2.) + a; 1744 w[i] *= (b - a) / 2.; 1745 } 1746 } 1747 PetscFunctionReturn(0); 1748 } 1749 1750 /*@ 1751 PetscDTGaussQuadrature - create Gauss-Legendre quadrature 1752 1753 Not Collective 1754 1755 Input Parameters: 1756 + npoints - number of points 1757 . a - left end of interval (often-1) 1758 - b - right end of interval (often +1) 1759 1760 Output Parameters: 1761 + x - quadrature points 1762 - w - quadrature weights 1763 1764 Level: intermediate 1765 1766 References: 1767 . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 1768 1769 .seealso: PetscDTLegendreEval() 1770 @*/ 1771 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 1772 { 1773 PetscInt i; 1774 PetscErrorCode ierr; 1775 1776 PetscFunctionBegin; 1777 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr); 1778 if (a != -1. || b != 1.) { /* shift */ 1779 for (i = 0; i < npoints; i++) { 1780 x[i] = (x[i] + 1.) * ((b - a) / 2.) + a; 1781 w[i] *= (b - a) / 2.; 1782 } 1783 } 1784 PetscFunctionReturn(0); 1785 } 1786 1787 /*@C 1788 PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre 1789 nodes of a given size on the domain [-1,1] 1790 1791 Not Collective 1792 1793 Input Parameters: 1794 + n - number of grid nodes 1795 - type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 1796 1797 Output Parameters: 1798 + x - quadrature points 1799 - w - quadrature weights 1800 1801 Notes: 1802 For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not 1803 close enough to the desired solution 1804 1805 These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes 1806 1807 See https://epubs.siam.org/doi/abs/10.1137/110855442 https://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes 1808 1809 Level: intermediate 1810 1811 .seealso: PetscDTGaussQuadrature() 1812 1813 @*/ 1814 PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w) 1815 { 1816 PetscBool newton; 1817 PetscErrorCode ierr; 1818 1819 PetscFunctionBegin; 1820 PetscCheckFalse(npoints < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element"); 1821 newton = (PetscBool) (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON); 1822 ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton);CHKERRQ(ierr); 1823 PetscFunctionReturn(0); 1824 } 1825 1826 /*@ 1827 PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 1828 1829 Not Collective 1830 1831 Input Parameters: 1832 + dim - The spatial dimension 1833 . Nc - The number of components 1834 . npoints - number of points in one dimension 1835 . a - left end of interval (often-1) 1836 - b - right end of interval (often +1) 1837 1838 Output Parameter: 1839 . q - A PetscQuadrature object 1840 1841 Level: intermediate 1842 1843 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 1844 @*/ 1845 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 1846 { 1847 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 1848 PetscReal *x, *w, *xw, *ww; 1849 PetscErrorCode ierr; 1850 1851 PetscFunctionBegin; 1852 ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 1853 ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 1854 /* Set up the Golub-Welsch system */ 1855 switch (dim) { 1856 case 0: 1857 ierr = PetscFree(x);CHKERRQ(ierr); 1858 ierr = PetscFree(w);CHKERRQ(ierr); 1859 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 1860 ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 1861 x[0] = 0.0; 1862 for (c = 0; c < Nc; ++c) w[c] = 1.0; 1863 break; 1864 case 1: 1865 ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 1866 ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 1867 for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 1868 ierr = PetscFree(ww);CHKERRQ(ierr); 1869 break; 1870 case 2: 1871 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 1872 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 1873 for (i = 0; i < npoints; ++i) { 1874 for (j = 0; j < npoints; ++j) { 1875 x[(i*npoints+j)*dim+0] = xw[i]; 1876 x[(i*npoints+j)*dim+1] = xw[j]; 1877 for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 1878 } 1879 } 1880 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 1881 break; 1882 case 3: 1883 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 1884 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 1885 for (i = 0; i < npoints; ++i) { 1886 for (j = 0; j < npoints; ++j) { 1887 for (k = 0; k < npoints; ++k) { 1888 x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 1889 x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 1890 x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 1891 for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 1892 } 1893 } 1894 } 1895 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 1896 break; 1897 default: 1898 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 1899 } 1900 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1901 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 1902 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 1903 ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr); 1904 PetscFunctionReturn(0); 1905 } 1906 1907 /*@ 1908 PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex 1909 1910 Not Collective 1911 1912 Input Parameters: 1913 + dim - The simplex dimension 1914 . Nc - The number of components 1915 . npoints - The number of points in one dimension 1916 . a - left end of interval (often-1) 1917 - b - right end of interval (often +1) 1918 1919 Output Parameter: 1920 . q - A PetscQuadrature object 1921 1922 Level: intermediate 1923 1924 References: 1925 . 1. - Karniadakis and Sherwin. FIAT 1926 1927 Note: For dim == 1, this is Gauss-Legendre quadrature 1928 1929 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 1930 @*/ 1931 PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 1932 { 1933 PetscInt totprev, totrem; 1934 PetscInt totpoints; 1935 PetscReal *p1, *w1; 1936 PetscReal *x, *w; 1937 PetscInt i, j, k, l, m, pt, c; 1938 PetscErrorCode ierr; 1939 1940 PetscFunctionBegin; 1941 PetscCheckFalse((a != -1.0) || (b != 1.0),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 1942 totpoints = 1; 1943 for (i = 0, totpoints = 1; i < dim; i++) totpoints *= npoints; 1944 ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 1945 ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 1946 ierr = PetscMalloc2(npoints, &p1, npoints, &w1);CHKERRQ(ierr); 1947 for (i = 0; i < totpoints*Nc; i++) w[i] = 1.; 1948 for (i = 0, totprev = 1, totrem = totpoints / npoints; i < dim; i++) { 1949 PetscReal mul; 1950 1951 mul = PetscPowReal(2.,-i); 1952 ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., i, 0.0, p1, w1);CHKERRQ(ierr); 1953 for (pt = 0, l = 0; l < totprev; l++) { 1954 for (j = 0; j < npoints; j++) { 1955 for (m = 0; m < totrem; m++, pt++) { 1956 for (k = 0; k < i; k++) x[pt*dim+k] = (x[pt*dim+k]+1.)*(1.-p1[j])*0.5 - 1.; 1957 x[pt * dim + i] = p1[j]; 1958 for (c = 0; c < Nc; c++) w[pt*Nc + c] *= mul * w1[j]; 1959 } 1960 } 1961 } 1962 totprev *= npoints; 1963 totrem /= npoints; 1964 } 1965 ierr = PetscFree2(p1, w1);CHKERRQ(ierr); 1966 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1967 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 1968 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 1969 ierr = PetscObjectChangeTypeName((PetscObject)*q,"StroudConical");CHKERRQ(ierr); 1970 PetscFunctionReturn(0); 1971 } 1972 1973 /*@ 1974 PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 1975 1976 Not Collective 1977 1978 Input Parameters: 1979 + dim - The cell dimension 1980 . level - The number of points in one dimension, 2^l 1981 . a - left end of interval (often-1) 1982 - b - right end of interval (often +1) 1983 1984 Output Parameter: 1985 . q - A PetscQuadrature object 1986 1987 Level: intermediate 1988 1989 .seealso: PetscDTGaussTensorQuadrature() 1990 @*/ 1991 PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 1992 { 1993 const PetscInt p = 16; /* Digits of precision in the evaluation */ 1994 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1995 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1996 const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 1997 PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 1998 PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 1999 PetscReal *x, *w; 2000 PetscInt K, k, npoints; 2001 PetscErrorCode ierr; 2002 2003 PetscFunctionBegin; 2004 PetscCheckFalse(dim > 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 2005 PetscCheckFalse(!level,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 2006 /* Find K such that the weights are < 32 digits of precision */ 2007 for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 2008 wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 2009 } 2010 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 2011 ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 2012 npoints = 2*K-1; 2013 ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 2014 ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 2015 /* Center term */ 2016 x[0] = beta; 2017 w[0] = 0.5*alpha*PETSC_PI; 2018 for (k = 1; k < K; ++k) { 2019 wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 2020 xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 2021 x[2*k-1] = -alpha*xk+beta; 2022 w[2*k-1] = wk; 2023 x[2*k+0] = alpha*xk+beta; 2024 w[2*k+0] = wk; 2025 } 2026 ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 2027 PetscFunctionReturn(0); 2028 } 2029 2030 PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 2031 { 2032 const PetscInt p = 16; /* Digits of precision in the evaluation */ 2033 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 2034 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 2035 PetscReal h = 1.0; /* Step size, length between x_k */ 2036 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 2037 PetscReal osum = 0.0; /* Integral on last level */ 2038 PetscReal psum = 0.0; /* Integral on the level before the last level */ 2039 PetscReal sum; /* Integral on current level */ 2040 PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 2041 PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 2042 PetscReal wk; /* Quadrature weight at x_k */ 2043 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 2044 PetscInt d; /* Digits of precision in the integral */ 2045 2046 PetscFunctionBegin; 2047 PetscCheckFalse(digits <= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 2048 /* Center term */ 2049 func(beta, &lval); 2050 sum = 0.5*alpha*PETSC_PI*lval; 2051 /* */ 2052 do { 2053 PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 2054 PetscInt k = 1; 2055 2056 ++l; 2057 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 2058 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 2059 psum = osum; 2060 osum = sum; 2061 h *= 0.5; 2062 sum *= 0.5; 2063 do { 2064 wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 2065 yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 2066 lx = -alpha*(1.0 - yk)+beta; 2067 rx = alpha*(1.0 - yk)+beta; 2068 func(lx, &lval); 2069 func(rx, &rval); 2070 lterm = alpha*wk*lval; 2071 maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 2072 sum += lterm; 2073 rterm = alpha*wk*rval; 2074 maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 2075 sum += rterm; 2076 ++k; 2077 /* Only need to evaluate every other point on refined levels */ 2078 if (l != 1) ++k; 2079 } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 2080 2081 d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 2082 d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 2083 d3 = PetscLog10Real(maxTerm) - p; 2084 if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 2085 else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 2086 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 2087 } while (d < digits && l < 12); 2088 *sol = sum; 2089 2090 PetscFunctionReturn(0); 2091 } 2092 2093 #if defined(PETSC_HAVE_MPFR) 2094 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 2095 { 2096 const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 2097 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 2098 mpfr_t alpha; /* Half-width of the integration interval */ 2099 mpfr_t beta; /* Center of the integration interval */ 2100 mpfr_t h; /* Step size, length between x_k */ 2101 mpfr_t osum; /* Integral on last level */ 2102 mpfr_t psum; /* Integral on the level before the last level */ 2103 mpfr_t sum; /* Integral on current level */ 2104 mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 2105 mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 2106 mpfr_t wk; /* Quadrature weight at x_k */ 2107 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 2108 PetscInt d; /* Digits of precision in the integral */ 2109 mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 2110 2111 PetscFunctionBegin; 2112 PetscCheckFalse(digits <= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 2113 /* Create high precision storage */ 2114 mpfr_inits2(PetscCeilReal(safetyFactor*digits*PetscLogReal(10.)/PetscLogReal(2.)), alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 2115 /* Initialization */ 2116 mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 2117 mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 2118 mpfr_set_d(osum, 0.0, MPFR_RNDN); 2119 mpfr_set_d(psum, 0.0, MPFR_RNDN); 2120 mpfr_set_d(h, 1.0, MPFR_RNDN); 2121 mpfr_const_pi(pi2, MPFR_RNDN); 2122 mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 2123 /* Center term */ 2124 func(0.5*(b+a), &lval); 2125 mpfr_set(sum, pi2, MPFR_RNDN); 2126 mpfr_mul(sum, sum, alpha, MPFR_RNDN); 2127 mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 2128 /* */ 2129 do { 2130 PetscReal d1, d2, d3, d4; 2131 PetscInt k = 1; 2132 2133 ++l; 2134 mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 2135 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 2136 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 2137 mpfr_set(psum, osum, MPFR_RNDN); 2138 mpfr_set(osum, sum, MPFR_RNDN); 2139 mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 2140 mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 2141 do { 2142 mpfr_set_si(kh, k, MPFR_RNDN); 2143 mpfr_mul(kh, kh, h, MPFR_RNDN); 2144 /* Weight */ 2145 mpfr_set(wk, h, MPFR_RNDN); 2146 mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 2147 mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 2148 mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 2149 mpfr_cosh(tmp, msinh, MPFR_RNDN); 2150 mpfr_sqr(tmp, tmp, MPFR_RNDN); 2151 mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 2152 mpfr_div(wk, wk, tmp, MPFR_RNDN); 2153 /* Abscissa */ 2154 mpfr_set_d(yk, 1.0, MPFR_RNDZ); 2155 mpfr_cosh(tmp, msinh, MPFR_RNDN); 2156 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 2157 mpfr_exp(tmp, msinh, MPFR_RNDN); 2158 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 2159 /* Quadrature points */ 2160 mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 2161 mpfr_mul(lx, lx, alpha, MPFR_RNDU); 2162 mpfr_add(lx, lx, beta, MPFR_RNDU); 2163 mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 2164 mpfr_mul(rx, rx, alpha, MPFR_RNDD); 2165 mpfr_add(rx, rx, beta, MPFR_RNDD); 2166 /* Evaluation */ 2167 func(mpfr_get_d(lx, MPFR_RNDU), &lval); 2168 func(mpfr_get_d(rx, MPFR_RNDD), &rval); 2169 /* Update */ 2170 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 2171 mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 2172 mpfr_add(sum, sum, tmp, MPFR_RNDN); 2173 mpfr_abs(tmp, tmp, MPFR_RNDN); 2174 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 2175 mpfr_set(curTerm, tmp, MPFR_RNDN); 2176 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 2177 mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 2178 mpfr_add(sum, sum, tmp, MPFR_RNDN); 2179 mpfr_abs(tmp, tmp, MPFR_RNDN); 2180 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 2181 mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 2182 ++k; 2183 /* Only need to evaluate every other point on refined levels */ 2184 if (l != 1) ++k; 2185 mpfr_log10(tmp, wk, MPFR_RNDN); 2186 mpfr_abs(tmp, tmp, MPFR_RNDN); 2187 } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 2188 mpfr_sub(tmp, sum, osum, MPFR_RNDN); 2189 mpfr_abs(tmp, tmp, MPFR_RNDN); 2190 mpfr_log10(tmp, tmp, MPFR_RNDN); 2191 d1 = mpfr_get_d(tmp, MPFR_RNDN); 2192 mpfr_sub(tmp, sum, psum, MPFR_RNDN); 2193 mpfr_abs(tmp, tmp, MPFR_RNDN); 2194 mpfr_log10(tmp, tmp, MPFR_RNDN); 2195 d2 = mpfr_get_d(tmp, MPFR_RNDN); 2196 mpfr_log10(tmp, maxTerm, MPFR_RNDN); 2197 d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 2198 mpfr_log10(tmp, curTerm, MPFR_RNDN); 2199 d4 = mpfr_get_d(tmp, MPFR_RNDN); 2200 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 2201 } while (d < digits && l < 8); 2202 *sol = mpfr_get_d(sum, MPFR_RNDN); 2203 /* Cleanup */ 2204 mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 2205 PetscFunctionReturn(0); 2206 } 2207 #else 2208 2209 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 2210 { 2211 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 2212 } 2213 #endif 2214 2215 /* Overwrites A. Can only handle full-rank problems with m>=n 2216 * A in column-major format 2217 * Ainv in row-major format 2218 * tau has length m 2219 * worksize must be >= max(1,n) 2220 */ 2221 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 2222 { 2223 PetscErrorCode ierr; 2224 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 2225 PetscScalar *A,*Ainv,*R,*Q,Alpha; 2226 2227 PetscFunctionBegin; 2228 #if defined(PETSC_USE_COMPLEX) 2229 { 2230 PetscInt i,j; 2231 ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 2232 for (j=0; j<n; j++) { 2233 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 2234 } 2235 mstride = m; 2236 } 2237 #else 2238 A = A_in; 2239 Ainv = Ainv_out; 2240 #endif 2241 2242 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 2243 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 2244 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 2245 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 2246 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 2247 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 2248 ierr = PetscFPTrapPop();CHKERRQ(ierr); 2249 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 2250 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 2251 2252 /* Extract an explicit representation of Q */ 2253 Q = Ainv; 2254 ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr); 2255 K = N; /* full rank */ 2256 PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 2257 PetscCheckFalse(info,PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 2258 2259 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 2260 Alpha = 1.0; 2261 ldb = lda; 2262 PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 2263 /* Ainv is Q, overwritten with inverse */ 2264 2265 #if defined(PETSC_USE_COMPLEX) 2266 { 2267 PetscInt i; 2268 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 2269 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 2270 } 2271 #endif 2272 PetscFunctionReturn(0); 2273 } 2274 2275 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 2276 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 2277 { 2278 PetscErrorCode ierr; 2279 PetscReal *Bv; 2280 PetscInt i,j; 2281 2282 PetscFunctionBegin; 2283 ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 2284 /* Point evaluation of L_p on all the source vertices */ 2285 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 2286 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 2287 for (i=0; i<ninterval; i++) { 2288 for (j=0; j<ndegree; j++) { 2289 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 2290 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 2291 } 2292 } 2293 ierr = PetscFree(Bv);CHKERRQ(ierr); 2294 PetscFunctionReturn(0); 2295 } 2296 2297 /*@ 2298 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 2299 2300 Not Collective 2301 2302 Input Parameters: 2303 + degree - degree of reconstruction polynomial 2304 . nsource - number of source intervals 2305 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 2306 . ntarget - number of target intervals 2307 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 2308 2309 Output Parameter: 2310 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 2311 2312 Level: advanced 2313 2314 .seealso: PetscDTLegendreEval() 2315 @*/ 2316 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 2317 { 2318 PetscErrorCode ierr; 2319 PetscInt i,j,k,*bdegrees,worksize; 2320 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 2321 PetscScalar *tau,*work; 2322 2323 PetscFunctionBegin; 2324 PetscValidRealPointer(sourcex,3); 2325 PetscValidRealPointer(targetx,5); 2326 PetscValidRealPointer(R,6); 2327 PetscCheckFalse(degree >= nsource,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Reconstruction degree %D must be less than number of source intervals %D",degree,nsource); 2328 if (PetscDefined(USE_DEBUG)) { 2329 for (i=0; i<nsource; i++) { 2330 PetscCheckFalse(sourcex[i] >= sourcex[i+1],PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Source interval %D has negative orientation (%g,%g)",i,(double)sourcex[i],(double)sourcex[i+1]); 2331 } 2332 for (i=0; i<ntarget; i++) { 2333 PetscCheckFalse(targetx[i] >= targetx[i+1],PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Target interval %D has negative orientation (%g,%g)",i,(double)targetx[i],(double)targetx[i+1]); 2334 } 2335 } 2336 xmin = PetscMin(sourcex[0],targetx[0]); 2337 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 2338 center = (xmin + xmax)/2; 2339 hscale = (xmax - xmin)/2; 2340 worksize = nsource; 2341 ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 2342 ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 2343 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 2344 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 2345 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 2346 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 2347 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 2348 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 2349 for (i=0; i<ntarget; i++) { 2350 PetscReal rowsum = 0; 2351 for (j=0; j<nsource; j++) { 2352 PetscReal sum = 0; 2353 for (k=0; k<degree+1; k++) { 2354 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 2355 } 2356 R[i*nsource+j] = sum; 2357 rowsum += sum; 2358 } 2359 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 2360 } 2361 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 2362 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 2363 PetscFunctionReturn(0); 2364 } 2365 2366 /*@C 2367 PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points 2368 2369 Not Collective 2370 2371 Input Parameters: 2372 + n - the number of GLL nodes 2373 . nodes - the GLL nodes 2374 . weights - the GLL weights 2375 - f - the function values at the nodes 2376 2377 Output Parameter: 2378 . in - the value of the integral 2379 2380 Level: beginner 2381 2382 .seealso: PetscDTGaussLobattoLegendreQuadrature() 2383 2384 @*/ 2385 PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in) 2386 { 2387 PetscInt i; 2388 2389 PetscFunctionBegin; 2390 *in = 0.; 2391 for (i=0; i<n; i++) { 2392 *in += f[i]*f[i]*weights[i]; 2393 } 2394 PetscFunctionReturn(0); 2395 } 2396 2397 /*@C 2398 PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element 2399 2400 Not Collective 2401 2402 Input Parameters: 2403 + n - the number of GLL nodes 2404 . nodes - the GLL nodes 2405 - weights - the GLL weights 2406 2407 Output Parameter: 2408 . A - the stiffness element 2409 2410 Level: beginner 2411 2412 Notes: 2413 Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy() 2414 2415 You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented (the array is symmetric) 2416 2417 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 2418 2419 @*/ 2420 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2421 { 2422 PetscReal **A; 2423 PetscErrorCode ierr; 2424 const PetscReal *gllnodes = nodes; 2425 const PetscInt p = n-1; 2426 PetscReal z0,z1,z2 = -1,x,Lpj,Lpr; 2427 PetscInt i,j,nn,r; 2428 2429 PetscFunctionBegin; 2430 ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 2431 ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 2432 for (i=1; i<n; i++) A[i] = A[i-1]+n; 2433 2434 for (j=1; j<p; j++) { 2435 x = gllnodes[j]; 2436 z0 = 1.; 2437 z1 = x; 2438 for (nn=1; nn<p; nn++) { 2439 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 2440 z0 = z1; 2441 z1 = z2; 2442 } 2443 Lpj=z2; 2444 for (r=1; r<p; r++) { 2445 if (r == j) { 2446 A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj); 2447 } else { 2448 x = gllnodes[r]; 2449 z0 = 1.; 2450 z1 = x; 2451 for (nn=1; nn<p; nn++) { 2452 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 2453 z0 = z1; 2454 z1 = z2; 2455 } 2456 Lpr = z2; 2457 A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r])); 2458 } 2459 } 2460 } 2461 for (j=1; j<p+1; j++) { 2462 x = gllnodes[j]; 2463 z0 = 1.; 2464 z1 = x; 2465 for (nn=1; nn<p; nn++) { 2466 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 2467 z0 = z1; 2468 z1 = z2; 2469 } 2470 Lpj = z2; 2471 A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j])); 2472 A[0][j] = A[j][0]; 2473 } 2474 for (j=0; j<p; j++) { 2475 x = gllnodes[j]; 2476 z0 = 1.; 2477 z1 = x; 2478 for (nn=1; nn<p; nn++) { 2479 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 2480 z0 = z1; 2481 z1 = z2; 2482 } 2483 Lpj=z2; 2484 2485 A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j])); 2486 A[j][p] = A[p][j]; 2487 } 2488 A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.; 2489 A[p][p]=A[0][0]; 2490 *AA = A; 2491 PetscFunctionReturn(0); 2492 } 2493 2494 /*@C 2495 PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element 2496 2497 Not Collective 2498 2499 Input Parameters: 2500 + n - the number of GLL nodes 2501 . nodes - the GLL nodes 2502 . weights - the GLL weightss 2503 - A - the stiffness element 2504 2505 Level: beginner 2506 2507 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate() 2508 2509 @*/ 2510 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2511 { 2512 PetscErrorCode ierr; 2513 2514 PetscFunctionBegin; 2515 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2516 ierr = PetscFree(*AA);CHKERRQ(ierr); 2517 *AA = NULL; 2518 PetscFunctionReturn(0); 2519 } 2520 2521 /*@C 2522 PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element 2523 2524 Not Collective 2525 2526 Input Parameter: 2527 + n - the number of GLL nodes 2528 . nodes - the GLL nodes 2529 . weights - the GLL weights 2530 2531 Output Parameters: 2532 . AA - the stiffness element 2533 - AAT - the transpose of AA (pass in NULL if you do not need this array) 2534 2535 Level: beginner 2536 2537 Notes: 2538 Destroy this with PetscGaussLobattoLegendreElementGradientDestroy() 2539 2540 You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented 2541 2542 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 2543 2544 @*/ 2545 PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 2546 { 2547 PetscReal **A, **AT = NULL; 2548 PetscErrorCode ierr; 2549 const PetscReal *gllnodes = nodes; 2550 const PetscInt p = n-1; 2551 PetscReal Li, Lj,d0; 2552 PetscInt i,j; 2553 2554 PetscFunctionBegin; 2555 ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 2556 ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 2557 for (i=1; i<n; i++) A[i] = A[i-1]+n; 2558 2559 if (AAT) { 2560 ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr); 2561 ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr); 2562 for (i=1; i<n; i++) AT[i] = AT[i-1]+n; 2563 } 2564 2565 if (n==1) {A[0][0] = 0.;} 2566 d0 = (PetscReal)p*((PetscReal)p+1.)/4.; 2567 for (i=0; i<n; i++) { 2568 for (j=0; j<n; j++) { 2569 A[i][j] = 0.; 2570 ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li);CHKERRQ(ierr); 2571 ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj);CHKERRQ(ierr); 2572 if (i!=j) A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j])); 2573 if ((j==i) && (i==0)) A[i][j] = -d0; 2574 if (j==i && i==p) A[i][j] = d0; 2575 if (AT) AT[j][i] = A[i][j]; 2576 } 2577 } 2578 if (AAT) *AAT = AT; 2579 *AA = A; 2580 PetscFunctionReturn(0); 2581 } 2582 2583 /*@C 2584 PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate() 2585 2586 Not Collective 2587 2588 Input Parameters: 2589 + n - the number of GLL nodes 2590 . nodes - the GLL nodes 2591 . weights - the GLL weights 2592 . AA - the stiffness element 2593 - AAT - the transpose of the element 2594 2595 Level: beginner 2596 2597 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate() 2598 2599 @*/ 2600 PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 2601 { 2602 PetscErrorCode ierr; 2603 2604 PetscFunctionBegin; 2605 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2606 ierr = PetscFree(*AA);CHKERRQ(ierr); 2607 *AA = NULL; 2608 if (*AAT) { 2609 ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr); 2610 ierr = PetscFree(*AAT);CHKERRQ(ierr); 2611 *AAT = NULL; 2612 } 2613 PetscFunctionReturn(0); 2614 } 2615 2616 /*@C 2617 PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element 2618 2619 Not Collective 2620 2621 Input Parameters: 2622 + n - the number of GLL nodes 2623 . nodes - the GLL nodes 2624 - weights - the GLL weightss 2625 2626 Output Parameter: 2627 . AA - the stiffness element 2628 2629 Level: beginner 2630 2631 Notes: 2632 Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy() 2633 2634 This is the same as the Gradient operator multiplied by the diagonal mass matrix 2635 2636 You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented 2637 2638 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy() 2639 2640 @*/ 2641 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2642 { 2643 PetscReal **D; 2644 PetscErrorCode ierr; 2645 const PetscReal *gllweights = weights; 2646 const PetscInt glln = n; 2647 PetscInt i,j; 2648 2649 PetscFunctionBegin; 2650 ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr); 2651 for (i=0; i<glln; i++) { 2652 for (j=0; j<glln; j++) { 2653 D[i][j] = gllweights[i]*D[i][j]; 2654 } 2655 } 2656 *AA = D; 2657 PetscFunctionReturn(0); 2658 } 2659 2660 /*@C 2661 PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element 2662 2663 Not Collective 2664 2665 Input Parameters: 2666 + n - the number of GLL nodes 2667 . nodes - the GLL nodes 2668 . weights - the GLL weights 2669 - A - advection 2670 2671 Level: beginner 2672 2673 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate() 2674 2675 @*/ 2676 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2677 { 2678 PetscErrorCode ierr; 2679 2680 PetscFunctionBegin; 2681 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2682 ierr = PetscFree(*AA);CHKERRQ(ierr); 2683 *AA = NULL; 2684 PetscFunctionReturn(0); 2685 } 2686 2687 PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2688 { 2689 PetscReal **A; 2690 PetscErrorCode ierr; 2691 const PetscReal *gllweights = weights; 2692 const PetscInt glln = n; 2693 PetscInt i,j; 2694 2695 PetscFunctionBegin; 2696 ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr); 2697 ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr); 2698 for (i=1; i<glln; i++) A[i] = A[i-1]+glln; 2699 if (glln==1) {A[0][0] = 0.;} 2700 for (i=0; i<glln; i++) { 2701 for (j=0; j<glln; j++) { 2702 A[i][j] = 0.; 2703 if (j==i) A[i][j] = gllweights[i]; 2704 } 2705 } 2706 *AA = A; 2707 PetscFunctionReturn(0); 2708 } 2709 2710 PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2711 { 2712 PetscErrorCode ierr; 2713 2714 PetscFunctionBegin; 2715 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2716 ierr = PetscFree(*AA);CHKERRQ(ierr); 2717 *AA = NULL; 2718 PetscFunctionReturn(0); 2719 } 2720 2721 /*@ 2722 PetscDTIndexToBary - convert an index into a barycentric coordinate. 2723 2724 Input Parameters: 2725 + len - the desired length of the barycentric tuple (usually 1 more than the dimension it represents, so a barycentric coordinate in a triangle has length 3) 2726 . sum - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to 2727 - index - the index to convert: should be >= 0 and < Binomial(len - 1 + sum, sum) 2728 2729 Output Parameter: 2730 . coord - will be filled with the barycentric coordinate 2731 2732 Level: beginner 2733 2734 Note: the indices map to barycentric coordinates in lexicographic order, where the first index is the 2735 least significant and the last index is the most significant. 2736 2737 .seealso: PetscDTBaryToIndex() 2738 @*/ 2739 PetscErrorCode PetscDTIndexToBary(PetscInt len, PetscInt sum, PetscInt index, PetscInt coord[]) 2740 { 2741 PetscInt c, d, s, total, subtotal, nexttotal; 2742 2743 PetscFunctionBeginHot; 2744 PetscCheckFalse(len < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative"); 2745 PetscCheckFalse(index < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index must be non-negative"); 2746 if (!len) { 2747 if (!sum && !index) PetscFunctionReturn(0); 2748 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate"); 2749 } 2750 for (c = 1, total = 1; c <= len; c++) { 2751 /* total is the number of ways to have a tuple of length c with sum */ 2752 if (index < total) break; 2753 total = (total * (sum + c)) / c; 2754 } 2755 PetscCheckFalse(c > len,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "index out of range"); 2756 for (d = c; d < len; d++) coord[d] = 0; 2757 for (s = 0, subtotal = 1, nexttotal = 1; c > 0;) { 2758 /* subtotal is the number of ways to have a tuple of length c with sum s */ 2759 /* nexttotal is the number of ways to have a tuple of length c-1 with sum s */ 2760 if ((index + subtotal) >= total) { 2761 coord[--c] = sum - s; 2762 index -= (total - subtotal); 2763 sum = s; 2764 total = nexttotal; 2765 subtotal = 1; 2766 nexttotal = 1; 2767 s = 0; 2768 } else { 2769 subtotal = (subtotal * (c + s)) / (s + 1); 2770 nexttotal = (nexttotal * (c - 1 + s)) / (s + 1); 2771 s++; 2772 } 2773 } 2774 PetscFunctionReturn(0); 2775 } 2776 2777 /*@ 2778 PetscDTBaryToIndex - convert a barycentric coordinate to an index 2779 2780 Input Parameters: 2781 + len - the desired length of the barycentric tuple (usually 1 more than the dimension it represents, so a barycentric coordinate in a triangle has length 3) 2782 . sum - the value that the sum of the barycentric coordinates (which will be non-negative integers) should sum to 2783 - coord - a barycentric coordinate with the given length and sum 2784 2785 Output Parameter: 2786 . index - the unique index for the coordinate, >= 0 and < Binomial(len - 1 + sum, sum) 2787 2788 Level: beginner 2789 2790 Note: the indices map to barycentric coordinates in lexicographic order, where the first index is the 2791 least significant and the last index is the most significant. 2792 2793 .seealso: PetscDTIndexToBary 2794 @*/ 2795 PetscErrorCode PetscDTBaryToIndex(PetscInt len, PetscInt sum, const PetscInt coord[], PetscInt *index) 2796 { 2797 PetscInt c; 2798 PetscInt i; 2799 PetscInt total; 2800 2801 PetscFunctionBeginHot; 2802 PetscCheckFalse(len < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "length must be non-negative"); 2803 if (!len) { 2804 if (!sum) { 2805 *index = 0; 2806 PetscFunctionReturn(0); 2807 } 2808 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid index or sum for length 0 barycentric coordinate"); 2809 } 2810 for (c = 1, total = 1; c < len; c++) total = (total * (sum + c)) / c; 2811 i = total - 1; 2812 c = len - 1; 2813 sum -= coord[c]; 2814 while (sum > 0) { 2815 PetscInt subtotal; 2816 PetscInt s; 2817 2818 for (s = 1, subtotal = 1; s < sum; s++) subtotal = (subtotal * (c + s)) / s; 2819 i -= subtotal; 2820 sum -= coord[--c]; 2821 } 2822 *index = i; 2823 PetscFunctionReturn(0); 2824 } 2825