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