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 static PetscBool GolubWelschCite = PETSC_FALSE; 16 const char GolubWelschCitation[] = "@article{GolubWelsch1969,\n" 17 " author = {Golub and Welsch},\n" 18 " title = {Calculation of Quadrature Rules},\n" 19 " journal = {Math. Comp.},\n" 20 " volume = {23},\n" 21 " number = {106},\n" 22 " pages = {221--230},\n" 23 " year = {1969}\n}\n"; 24 25 /* Numerical tests in src/dm/dt/examples/tests/ex1.c show that when computing the nodes and weights of Gauss-Jacobi 26 quadrature rules: 27 28 - in double precision, Newton's method and Golub & Welsch both work for moderate degrees (< 100), 29 - in single precision, Newton's method starts producing incorrect roots around n = 15, but 30 the weights from Golub & Welsch become a problem before then: they produces errors 31 in computing the Jacobi-polynomial Gram matrix around n = 6. 32 33 So we default to Newton's method (required fewer dependencies) */ 34 PetscBool PetscDTGaussQuadratureNewton_Internal = PETSC_TRUE; 35 36 PetscClassId PETSCQUADRATURE_CLASSID = 0; 37 38 /*@ 39 PetscQuadratureCreate - Create a PetscQuadrature object 40 41 Collective 42 43 Input Parameter: 44 . comm - The communicator for the PetscQuadrature object 45 46 Output Parameter: 47 . q - The PetscQuadrature object 48 49 Level: beginner 50 51 .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData() 52 @*/ 53 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) 54 { 55 PetscErrorCode ierr; 56 57 PetscFunctionBegin; 58 PetscValidPointer(q, 2); 59 ierr = DMInitializePackage();CHKERRQ(ierr); 60 ierr = PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); 61 (*q)->dim = -1; 62 (*q)->Nc = 1; 63 (*q)->order = -1; 64 (*q)->numPoints = 0; 65 (*q)->points = NULL; 66 (*q)->weights = NULL; 67 PetscFunctionReturn(0); 68 } 69 70 /*@ 71 PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object 72 73 Collective on q 74 75 Input Parameter: 76 . q - The PetscQuadrature object 77 78 Output Parameter: 79 . r - The new PetscQuadrature object 80 81 Level: beginner 82 83 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData() 84 @*/ 85 PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r) 86 { 87 PetscInt order, dim, Nc, Nq; 88 const PetscReal *points, *weights; 89 PetscReal *p, *w; 90 PetscErrorCode ierr; 91 92 PetscFunctionBegin; 93 PetscValidPointer(q, 2); 94 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr); 95 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 96 ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr); 97 ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 98 ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr); 99 ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr); 100 ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr); 101 ierr = PetscArraycpy(w, weights, Nc * Nq);CHKERRQ(ierr); 102 ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr); 103 PetscFunctionReturn(0); 104 } 105 106 /*@ 107 PetscQuadratureDestroy - Destroys a PetscQuadrature object 108 109 Collective on q 110 111 Input Parameter: 112 . q - The PetscQuadrature object 113 114 Level: beginner 115 116 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 117 @*/ 118 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) 119 { 120 PetscErrorCode ierr; 121 122 PetscFunctionBegin; 123 if (!*q) PetscFunctionReturn(0); 124 PetscValidHeaderSpecific((*q),PETSCQUADRATURE_CLASSID,1); 125 if (--((PetscObject)(*q))->refct > 0) { 126 *q = NULL; 127 PetscFunctionReturn(0); 128 } 129 ierr = PetscFree((*q)->points);CHKERRQ(ierr); 130 ierr = PetscFree((*q)->weights);CHKERRQ(ierr); 131 ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); 132 PetscFunctionReturn(0); 133 } 134 135 /*@ 136 PetscQuadratureGetOrder - Return the order of the method 137 138 Not collective 139 140 Input Parameter: 141 . q - The PetscQuadrature object 142 143 Output Parameter: 144 . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 145 146 Level: intermediate 147 148 .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 149 @*/ 150 PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order) 151 { 152 PetscFunctionBegin; 153 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 154 PetscValidPointer(order, 2); 155 *order = q->order; 156 PetscFunctionReturn(0); 157 } 158 159 /*@ 160 PetscQuadratureSetOrder - Return the order of the method 161 162 Not collective 163 164 Input Parameters: 165 + q - The PetscQuadrature object 166 - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated 167 168 Level: intermediate 169 170 .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData() 171 @*/ 172 PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order) 173 { 174 PetscFunctionBegin; 175 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 176 q->order = order; 177 PetscFunctionReturn(0); 178 } 179 180 /*@ 181 PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated 182 183 Not collective 184 185 Input Parameter: 186 . q - The PetscQuadrature object 187 188 Output Parameter: 189 . Nc - The number of components 190 191 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 192 193 Level: intermediate 194 195 .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 196 @*/ 197 PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc) 198 { 199 PetscFunctionBegin; 200 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 201 PetscValidPointer(Nc, 2); 202 *Nc = q->Nc; 203 PetscFunctionReturn(0); 204 } 205 206 /*@ 207 PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated 208 209 Not collective 210 211 Input Parameters: 212 + q - The PetscQuadrature object 213 - Nc - The number of components 214 215 Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components. 216 217 Level: intermediate 218 219 .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData() 220 @*/ 221 PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc) 222 { 223 PetscFunctionBegin; 224 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 225 q->Nc = Nc; 226 PetscFunctionReturn(0); 227 } 228 229 /*@C 230 PetscQuadratureGetData - Returns the data defining the quadrature 231 232 Not collective 233 234 Input Parameter: 235 . q - The PetscQuadrature object 236 237 Output Parameters: 238 + dim - The spatial dimension 239 . Nc - The number of components 240 . npoints - The number of quadrature points 241 . points - The coordinates of each quadrature point 242 - weights - The weight of each quadrature point 243 244 Level: intermediate 245 246 Fortran Notes: 247 From Fortran you must call PetscQuadratureRestoreData() when you are done with the data 248 249 .seealso: PetscQuadratureCreate(), PetscQuadratureSetData() 250 @*/ 251 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) 252 { 253 PetscFunctionBegin; 254 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 255 if (dim) { 256 PetscValidPointer(dim, 2); 257 *dim = q->dim; 258 } 259 if (Nc) { 260 PetscValidPointer(Nc, 3); 261 *Nc = q->Nc; 262 } 263 if (npoints) { 264 PetscValidPointer(npoints, 4); 265 *npoints = q->numPoints; 266 } 267 if (points) { 268 PetscValidPointer(points, 5); 269 *points = q->points; 270 } 271 if (weights) { 272 PetscValidPointer(weights, 6); 273 *weights = q->weights; 274 } 275 PetscFunctionReturn(0); 276 } 277 278 static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[]) 279 { 280 PetscScalar *Js, *Jinvs; 281 PetscInt i, j, k; 282 PetscBLASInt bm, bn, info; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 ierr = PetscBLASIntCast(m, &bm);CHKERRQ(ierr); 287 ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr); 288 #if defined(PETSC_USE_COMPLEX) 289 ierr = PetscMalloc2(m*n, &Js, m*n, &Jinvs);CHKERRQ(ierr); 290 for (i = 0; i < m*n; i++) Js[i] = J[i]; 291 #else 292 Js = (PetscReal *) J; 293 Jinvs = Jinv; 294 #endif 295 if (m == n) { 296 PetscBLASInt *pivots; 297 PetscScalar *W; 298 299 ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 300 301 ierr = PetscArraycpy(Jinvs, Js, m * m);CHKERRQ(ierr); 302 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info)); 303 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 304 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info)); 305 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 306 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 307 } else if (m < n) { 308 PetscScalar *JJT; 309 PetscBLASInt *pivots; 310 PetscScalar *W; 311 312 ierr = PetscMalloc1(m*m, &JJT);CHKERRQ(ierr); 313 ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr); 314 for (i = 0; i < m; i++) { 315 for (j = 0; j < m; j++) { 316 PetscScalar val = 0.; 317 318 for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k]; 319 JJT[i * m + j] = val; 320 } 321 } 322 323 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info)); 324 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 325 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info)); 326 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 327 for (i = 0; i < n; i++) { 328 for (j = 0; j < m; j++) { 329 PetscScalar val = 0.; 330 331 for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j]; 332 Jinvs[i * m + j] = val; 333 } 334 } 335 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 336 ierr = PetscFree(JJT);CHKERRQ(ierr); 337 } else { 338 PetscScalar *JTJ; 339 PetscBLASInt *pivots; 340 PetscScalar *W; 341 342 ierr = PetscMalloc1(n*n, &JTJ);CHKERRQ(ierr); 343 ierr = PetscMalloc2(n, &pivots, n, &W);CHKERRQ(ierr); 344 for (i = 0; i < n; i++) { 345 for (j = 0; j < n; j++) { 346 PetscScalar val = 0.; 347 348 for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j]; 349 JTJ[i * n + j] = val; 350 } 351 } 352 353 PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bm, pivots, &info)); 354 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 355 PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info)); 356 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 357 for (i = 0; i < n; i++) { 358 for (j = 0; j < m; j++) { 359 PetscScalar val = 0.; 360 361 for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k]; 362 Jinvs[i * m + j] = val; 363 } 364 } 365 ierr = PetscFree2(pivots, W);CHKERRQ(ierr); 366 ierr = PetscFree(JTJ);CHKERRQ(ierr); 367 } 368 #if defined(PETSC_USE_COMPLEX) 369 for (i = 0; i < m*n; i++) Jinv[i] = PetscRealPart(Jinvs[i]); 370 ierr = PetscFree2(Js, Jinvs);CHKERRQ(ierr); 371 #endif 372 PetscFunctionReturn(0); 373 } 374 375 /*@ 376 PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation. 377 378 Collecive on PetscQuadrature 379 380 Input Arguments: 381 + q - the quadrature functional 382 . imageDim - the dimension of the image of the transformation 383 . origin - a point in the original space 384 . originImage - the image of the origin under the transformation 385 . J - the Jacobian of the image: an [imageDim x dim] matrix in row major order 386 - 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] 387 388 Output Arguments: 389 . 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. 390 391 Note: the new quadrature rule will have a different number of components if spaces have different dimensions. For example, pushing a 2-form forward from a two dimensional space to a three dimensional space changes the number of components from 1 to 3. 392 393 .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 394 @*/ 395 PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq) 396 { 397 PetscInt dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c; 398 const PetscReal *points; 399 const PetscReal *weights; 400 PetscReal *imagePoints, *imageWeights; 401 PetscReal *Jinv; 402 PetscReal *Jinvstar; 403 PetscErrorCode ierr; 404 405 PetscFunctionBegin; 406 PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); 407 if (imageDim < PetscAbsInt(formDegree)) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %D-form in %D dimensions", PetscAbsInt(formDegree), imageDim); 408 ierr = PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights);CHKERRQ(ierr); 409 ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize);CHKERRQ(ierr); 410 if (Nc % formSize) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of formSize %D\n", Nc, formSize); 411 Ncopies = Nc / formSize; 412 ierr = PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize);CHKERRQ(ierr); 413 imageNc = Ncopies * imageFormSize; 414 ierr = PetscMalloc1(Npoints * imageDim, &imagePoints);CHKERRQ(ierr); 415 ierr = PetscMalloc1(Npoints * imageNc, &imageWeights);CHKERRQ(ierr); 416 ierr = PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar);CHKERRQ(ierr); 417 ierr = PetscDTJacobianInverse_Internal(dim, imageDim, J, Jinv);CHKERRQ(ierr); 418 ierr = PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar);CHKERRQ(ierr); 419 for (pt = 0; pt < Npoints; pt++) { 420 const PetscReal *point = &points[pt * dim]; 421 PetscReal *imagePoint = &imagePoints[pt * imageDim]; 422 423 for (i = 0; i < imageDim; i++) { 424 PetscReal val = originImage[i]; 425 426 for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]); 427 imagePoint[i] = val; 428 } 429 for (c = 0; c < Ncopies; c++) { 430 const PetscReal *form = &weights[pt * Nc + c * formSize]; 431 PetscReal *imageForm = &imageWeights[pt * imageNc + c * imageFormSize]; 432 433 for (i = 0; i < imageFormSize; i++) { 434 PetscReal val = 0.; 435 436 for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j]; 437 imageForm[i] = val; 438 } 439 } 440 } 441 ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq);CHKERRQ(ierr); 442 ierr = PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights);CHKERRQ(ierr); 443 ierr = PetscFree2(Jinv, Jinvstar);CHKERRQ(ierr); 444 PetscFunctionReturn(0); 445 } 446 447 /*@C 448 PetscQuadratureSetData - Sets the data defining the quadrature 449 450 Not collective 451 452 Input Parameters: 453 + q - The PetscQuadrature object 454 . dim - The spatial dimension 455 . Nc - The number of components 456 . npoints - The number of quadrature points 457 . points - The coordinates of each quadrature point 458 - weights - The weight of each quadrature point 459 460 Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them. 461 462 Level: intermediate 463 464 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 465 @*/ 466 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) 467 { 468 PetscFunctionBegin; 469 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 470 if (dim >= 0) q->dim = dim; 471 if (Nc >= 0) q->Nc = Nc; 472 if (npoints >= 0) q->numPoints = npoints; 473 if (points) { 474 PetscValidPointer(points, 4); 475 q->points = points; 476 } 477 if (weights) { 478 PetscValidPointer(weights, 5); 479 q->weights = weights; 480 } 481 PetscFunctionReturn(0); 482 } 483 484 static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v) 485 { 486 PetscInt q, d, c; 487 PetscViewerFormat format; 488 PetscErrorCode ierr; 489 490 PetscFunctionBegin; 491 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);CHKERRQ(ierr);} 492 else {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);} 493 ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 494 if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0); 495 for (q = 0; q < quad->numPoints; ++q) { 496 ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr); 497 ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr); 498 for (d = 0; d < quad->dim; ++d) { 499 if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 500 ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); 501 } 502 ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr); 503 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);} 504 for (c = 0; c < quad->Nc; ++c) { 505 if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);} 506 ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr); 507 } 508 if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);} 509 ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr); 510 ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr); 511 } 512 PetscFunctionReturn(0); 513 } 514 515 /*@C 516 PetscQuadratureView - Views a PetscQuadrature object 517 518 Collective on quad 519 520 Input Parameters: 521 + quad - The PetscQuadrature object 522 - viewer - The PetscViewer object 523 524 Level: beginner 525 526 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData() 527 @*/ 528 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) 529 { 530 PetscBool iascii; 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 PetscValidHeader(quad, 1); 535 if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 536 if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);} 537 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 538 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 539 if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);} 540 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 541 PetscFunctionReturn(0); 542 } 543 544 /*@C 545 PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement 546 547 Not collective 548 549 Input Parameter: 550 + q - The original PetscQuadrature 551 . numSubelements - The number of subelements the original element is divided into 552 . v0 - An array of the initial points for each subelement 553 - jac - An array of the Jacobian mappings from the reference to each subelement 554 555 Output Parameters: 556 . dim - The dimension 557 558 Note: Together v0 and jac define an affine mapping from the original reference element to each subelement 559 560 Not available from Fortran 561 562 Level: intermediate 563 564 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 565 @*/ 566 PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref) 567 { 568 const PetscReal *points, *weights; 569 PetscReal *pointsRef, *weightsRef; 570 PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e; 571 PetscErrorCode ierr; 572 573 PetscFunctionBegin; 574 PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1); 575 PetscValidPointer(v0, 3); 576 PetscValidPointer(jac, 4); 577 PetscValidPointer(qref, 5); 578 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr); 579 ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr); 580 ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr); 581 npointsRef = npoints*numSubelements; 582 ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr); 583 ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr); 584 for (c = 0; c < numSubelements; ++c) { 585 for (p = 0; p < npoints; ++p) { 586 for (d = 0; d < dim; ++d) { 587 pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d]; 588 for (e = 0; e < dim; ++e) { 589 pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0); 590 } 591 } 592 /* Could also use detJ here */ 593 for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements; 594 } 595 } 596 ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr); 597 ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr); 598 PetscFunctionReturn(0); 599 } 600 601 /* Compute the coefficients for the Jacobi polynomial recurrence, 602 * 603 * J^{a,b}_n(x) = (cnm1 + cnm1x * x) * J^{a,b}_{n-1}(x) - cnm2 * J^{a,b}_{n-2}(x). 604 */ 605 #define PetscDTJacobiRecurrence_Internal(n,a,b,cnm1,cnm1x,cnm2) \ 606 do { \ 607 PetscReal _a = (a); \ 608 PetscReal _b = (b); \ 609 PetscReal _n = (n); \ 610 if (n == 1) { \ 611 (cnm1) = (_a-_b) * 0.5; \ 612 (cnm1x) = (_a+_b+2.)*0.5; \ 613 (cnm2) = 0.; \ 614 } else { \ 615 PetscReal _2n = _n+_n; \ 616 PetscReal _d = (_2n*(_n+_a+_b)*(_2n+_a+_b-2)); \ 617 PetscReal _n1 = (_2n+_a+_b-1.)*(_a*_a-_b*_b); \ 618 PetscReal _n1x = (_2n+_a+_b-1.)*(_2n+_a+_b)*(_2n+_a+_b-2); \ 619 PetscReal _n2 = 2.*((_n+_a-1.)*(_n+_b-1.)*(_2n+_a+_b)); \ 620 (cnm1) = _n1 / _d; \ 621 (cnm1x) = _n1x / _d; \ 622 (cnm2) = _n2 / _d; \ 623 } \ 624 } while (0) 625 626 static PetscErrorCode PetscDTJacobiEval_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscInt k, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *p) 627 { 628 PetscReal ak, bk; 629 PetscReal abk1; 630 PetscInt i,l,maxdegree; 631 632 PetscFunctionBegin; 633 maxdegree = degrees[ndegree-1] - k; 634 ak = a + k; 635 bk = b + k; 636 abk1 = a + b + k + 1.; 637 if (maxdegree < 0) { 638 for (i = 0; i < npoints; i++) for (l = 0; l < ndegree; l++) p[i*ndegree+l] = 0.; 639 PetscFunctionReturn(0); 640 } 641 for (i=0; i<npoints; i++) { 642 PetscReal pm1,pm2,x; 643 PetscReal cnm1, cnm1x, cnm2; 644 PetscInt j,m; 645 646 x = points[i]; 647 pm2 = 1.; 648 PetscDTJacobiRecurrence_Internal(1,ak,bk,cnm1,cnm1x,cnm2); 649 pm1 = (cnm1 + cnm1x*x); 650 l = 0; 651 while (l < ndegree && degrees[l] - k < 0) { 652 p[l++] = 0.; 653 } 654 while (l < ndegree && degrees[l] - k == 0) { 655 p[l] = pm2; 656 for (m = 0; m < k; m++) p[l] *= (abk1 + m) * 0.5; 657 l++; 658 } 659 while (l < ndegree && degrees[l] - k == 1) { 660 p[l] = pm1; 661 for (m = 0; m < k; m++) p[l] *= (abk1 + 1 + m) * 0.5; 662 l++; 663 } 664 for (j=2; j<=maxdegree; j++) { 665 PetscReal pp; 666 667 PetscDTJacobiRecurrence_Internal(j,ak,bk,cnm1,cnm1x,cnm2); 668 pp = (cnm1 + cnm1x*x)*pm1 - cnm2*pm2; 669 pm2 = pm1; 670 pm1 = pp; 671 while (l < ndegree && degrees[l] - k == j) { 672 p[l] = pp; 673 for (m = 0; m < k; m++) p[l] *= (abk1 + j + m) * 0.5; 674 l++; 675 } 676 } 677 p += ndegree; 678 } 679 PetscFunctionReturn(0); 680 } 681 682 /*@ 683 PetscDTJacobiEval - evaluate Jacobi polynomials for the weight function $(1.+x)^{\alpha} (1.-x)^{\beta}$ 684 at points 685 686 Not Collective 687 688 Input Arguments: 689 + npoints - number of spatial points to evaluate at 690 . alpha - the left exponent > -1 691 . beta - the right exponent > -1 692 . points - array of locations to evaluate at 693 . ndegree - number of basis degrees to evaluate 694 - degrees - sorted array of degrees to evaluate 695 696 Output Arguments: 697 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 698 . D - row-oriented derivative evaluation matrix (or NULL) 699 - D2 - row-oriented second derivative evaluation matrix (or NULL) 700 701 Level: intermediate 702 703 .seealso: PetscDTGaussQuadrature() 704 @*/ 705 PetscErrorCode PetscDTJacobiEval(PetscInt npoints,PetscReal alpha, PetscReal beta, const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 706 { 707 PetscErrorCode ierr; 708 709 PetscFunctionBegin; 710 if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 711 if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 712 if (!npoints || !ndegree) PetscFunctionReturn(0); 713 if (B) {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 0, points, ndegree, degrees, B);CHKERRQ(ierr);} 714 if (D) {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 1, points, ndegree, degrees, D);CHKERRQ(ierr);} 715 if (D2) {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 2, points, ndegree, degrees, D2);CHKERRQ(ierr);} 716 PetscFunctionReturn(0); 717 } 718 719 /*@ 720 PetscDTLegendreEval - evaluate Legendre polynomials at points 721 722 Not Collective 723 724 Input Arguments: 725 + npoints - number of spatial points to evaluate at 726 . points - array of locations to evaluate at 727 . ndegree - number of basis degrees to evaluate 728 - degrees - sorted array of degrees to evaluate 729 730 Output Arguments: 731 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 732 . D - row-oriented derivative evaluation matrix (or NULL) 733 - D2 - row-oriented second derivative evaluation matrix (or NULL) 734 735 Level: intermediate 736 737 .seealso: PetscDTGaussQuadrature() 738 @*/ 739 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 740 { 741 PetscErrorCode ierr; 742 743 PetscFunctionBegin; 744 ierr = PetscDTJacobiEval(npoints, 0., 0., points, ndegree, degrees, B, D, D2);CHKERRQ(ierr); 745 PetscFunctionReturn(0); 746 } 747 748 /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V 749 * with lds n; diag and subdiag are overwritten */ 750 static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[], 751 PetscReal eigs[], PetscScalar V[]) 752 { 753 char jobz = 'V'; /* eigenvalues and eigenvectors */ 754 char range = 'A'; /* all eigenvalues will be found */ 755 PetscReal VL = 0.; /* ignored because range is 'A' */ 756 PetscReal VU = 0.; /* ignored because range is 'A' */ 757 PetscBLASInt IL = 0; /* ignored because range is 'A' */ 758 PetscBLASInt IU = 0; /* ignored because range is 'A' */ 759 PetscReal abstol = 0.; /* unused */ 760 PetscBLASInt bn, bm, ldz; /* bm will equal bn on exit */ 761 PetscBLASInt *isuppz; 762 PetscBLASInt lwork, liwork; 763 PetscReal workquery; 764 PetscBLASInt iworkquery; 765 PetscBLASInt *iwork; 766 PetscBLASInt info; 767 PetscReal *work = NULL; 768 PetscErrorCode ierr; 769 770 PetscFunctionBegin; 771 #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG) 772 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found"); 773 #endif 774 ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr); 775 ierr = PetscBLASIntCast(n, &ldz);CHKERRQ(ierr); 776 #if !defined(PETSC_MISSING_LAPACK_STEGR) 777 ierr = PetscMalloc1(2 * n, &isuppz);CHKERRQ(ierr); 778 lwork = -1; 779 liwork = -1; 780 PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,&workquery,&lwork,&iworkquery,&liwork,&info)); 781 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error"); 782 lwork = (PetscBLASInt) workquery; 783 liwork = (PetscBLASInt) iworkquery; 784 ierr = PetscMalloc2(lwork, &work, liwork, &iwork);CHKERRQ(ierr); 785 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 786 PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,work,&lwork,iwork,&liwork,&info)); 787 ierr = PetscFPTrapPop();CHKERRQ(ierr); 788 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error"); 789 ierr = PetscFree2(work, iwork);CHKERRQ(ierr); 790 ierr = PetscFree(isuppz);CHKERRQ(ierr); 791 #elif !defined(PETSC_MISSING_LAPACK_STEQR) 792 jobz = 'I'; /* Compute eigenvalues and eigenvectors of the 793 tridiagonal matrix. Z is initialized to the identity 794 matrix. */ 795 ierr = PetscMalloc1(PetscMax(1,2*n-2),&work);CHKERRQ(ierr); 796 PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&bn,diag,subdiag,V,&ldz,work,&info)); 797 ierr = PetscFPTrapPop();CHKERRQ(ierr); 798 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 799 ierr = PetscFree(work);CHKERRQ(ierr); 800 ierr = PetscArraycpy(eigs,diag,n);CHKERRQ(ierr); 801 #endif 802 PetscFunctionReturn(0); 803 } 804 805 /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi 806 * quadrature rules on the interval [-1, 1] */ 807 static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw) 808 { 809 PetscReal twoab1; 810 PetscInt m = n - 2; 811 PetscReal a = alpha + 1.; 812 PetscReal b = beta + 1.; 813 PetscReal gra, grb; 814 815 PetscFunctionBegin; 816 twoab1 = PetscPowReal(2., a + b - 1.); 817 #if defined(PETSC_HAVE_LGAMMA) 818 grb = PetscExpReal(2. * PetscLGamma(b+1.) + PetscLGamma(m+1.) + PetscLGamma(m+a+1.) - 819 (PetscLGamma(m+b+1) + PetscLGamma(m+a+b+1.))); 820 gra = PetscExpReal(2. * PetscLGamma(a+1.) + PetscLGamma(m+1.) + PetscLGamma(m+b+1.) - 821 (PetscLGamma(m+a+1) + PetscLGamma(m+a+b+1.))); 822 #else 823 { 824 PetscInt alphai = (PetscInt) alpha; 825 PetscInt betai = (PetscInt) beta; 826 PetscErrorCode ierr; 827 828 if ((PetscReal) alphai == alpha && (PetscReal) betai == beta) { 829 PetscReal binom1, binom2; 830 831 ierr = PetscDTBinomial(m+b, b, &binom1);CHKERRQ(ierr); 832 ierr = PetscDTBinomial(m+a+b, b, &binom2);CHKERRQ(ierr); 833 grb = 1./ (binom1 * binom2); 834 ierr = PetscDTBinomial(m+a, a, &binom1);CHKERRQ(ierr); 835 ierr = PetscDTBinomial(m+a+b, a, &binom2);CHKERRQ(ierr); 836 gra = 1./ (binom1 * binom2); 837 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable."); 838 } 839 #endif 840 *leftw = twoab1 * grb / b; 841 *rightw = twoab1 * gra / a; 842 PetscFunctionReturn(0); 843 } 844 845 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x. 846 Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */ 847 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P) 848 { 849 PetscReal pn1, pn2; 850 PetscReal cnm1, cnm1x, cnm2; 851 PetscInt k; 852 853 PetscFunctionBegin; 854 if (!n) {*P = 1.0; PetscFunctionReturn(0);} 855 PetscDTJacobiRecurrence_Internal(1,a,b,cnm1,cnm1x,cnm2); 856 pn2 = 1.; 857 pn1 = cnm1 + cnm1x*x; 858 if (n == 1) {*P = pn1; PetscFunctionReturn(0);} 859 *P = 0.0; 860 for (k = 2; k < n+1; ++k) { 861 PetscDTJacobiRecurrence_Internal(k,a,b,cnm1,cnm1x,cnm2); 862 863 *P = (cnm1 + cnm1x*x)*pn1 - cnm2*pn2; 864 pn2 = pn1; 865 pn1 = *P; 866 } 867 PetscFunctionReturn(0); 868 } 869 870 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */ 871 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P) 872 { 873 PetscReal nP; 874 PetscInt i; 875 PetscErrorCode ierr; 876 877 PetscFunctionBegin; 878 if (k > n) {*P = 0.0; PetscFunctionReturn(0);} 879 ierr = PetscDTComputeJacobi(a+k, b+k, n-k, x, &nP);CHKERRQ(ierr); 880 for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5; 881 *P = nP; 882 PetscFunctionReturn(0); 883 } 884 885 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 886 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta) 887 { 888 PetscFunctionBegin; 889 *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0; 890 *eta = y; 891 PetscFunctionReturn(0); 892 } 893 894 static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[]) 895 { 896 PetscInt maxIter = 100; 897 PetscReal eps = PetscExpReal(0.75 * PetscLogReal(PETSC_MACHINE_EPSILON)); 898 PetscReal a1, a2, a3, a4, a5, a6, gf; 899 PetscInt k; 900 PetscErrorCode ierr; 901 902 PetscFunctionBegin; 903 904 a1 = PetscPowReal(2.0, a+b+1); 905 #if defined(PETSC_HAVE_LGAMMA) 906 a2 = PetscLGamma(a + npoints + 1); 907 a3 = PetscLGamma(b + npoints + 1); 908 a4 = PetscLGamma(a + b + npoints + 1); 909 a5 = PetscLGamma(npoints + 1); 910 gf = PetscExpReal(a2 + a3 - (a4 + a5)); 911 #else 912 { 913 PetscInt ia, ib; 914 915 ia = (PetscInt) a; 916 ib = (PetscInt) b; 917 gf = 1.; 918 if (ia == a && ia >= 0) { /* compute ratio of rising factorals wrt a */ 919 for (k = 0; k < ia; k++) gf *= (npoints + 1. + k) / (npoints + b + 1. + k); 920 } else if (b == b && ib >= 0) { /* compute ratio of rising factorials wrt b */ 921 for (k = 0; k < ib; k++) gf *= (npoints + 1. + k) / (npoints + a + 1. + k); 922 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable."); 923 } 924 #endif 925 926 a6 = a1 * gf; 927 /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses. 928 Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */ 929 for (k = 0; k < npoints; ++k) { 930 PetscReal r = PetscCosReal(PETSC_PI * (1. - (4.*k + 3. + 2.*b) / (4.*npoints + 2.*(a + b + 1.)))), dP; 931 PetscInt j; 932 933 if (k > 0) r = 0.5 * (r + x[k-1]); 934 for (j = 0; j < maxIter; ++j) { 935 PetscReal s = 0.0, delta, f, fp; 936 PetscInt i; 937 938 for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); 939 ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); 940 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp);CHKERRQ(ierr); 941 delta = f / (fp - f * s); 942 r = r - delta; 943 if (PetscAbsReal(delta) < eps) break; 944 } 945 x[k] = r; 946 ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP);CHKERRQ(ierr); 947 w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); 948 } 949 PetscFunctionReturn(0); 950 } 951 952 /* Compute the diagonals of the Jacobi matrix used in Golub & Welsch algorithms for Gauss-Jacobi 953 * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */ 954 static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s) 955 { 956 PetscInt i; 957 958 PetscFunctionBegin; 959 for (i = 0; i < nPoints; i++) { 960 PetscReal A, B, C; 961 962 PetscDTJacobiRecurrence_Internal(i+1,a,b,A,B,C); 963 d[i] = -A / B; 964 if (i) s[i-1] *= C / B; 965 if (i < nPoints - 1) s[i] = 1. / B; 966 } 967 PetscFunctionReturn(0); 968 } 969 970 static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w) 971 { 972 PetscReal mu0; 973 PetscReal ga, gb, gab; 974 PetscInt i; 975 PetscErrorCode ierr; 976 977 PetscFunctionBegin; 978 ierr = PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite);CHKERRQ(ierr); 979 980 #if defined(PETSC_HAVE_TGAMMA) 981 ga = PetscTGamma(a + 1); 982 gb = PetscTGamma(b + 1); 983 gab = PetscTGamma(a + b + 2); 984 #else 985 { 986 PetscInt ia, ib; 987 988 ia = (PetscInt) a; 989 ib = (PetscInt) b; 990 if (ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0) { /* All gamma(x) terms are (x-1)! terms */ 991 ierr = PetscDTFactorial(ia, &ga);CHKERRQ(ierr); 992 ierr = PetscDTFactorial(ib, &gb);CHKERRQ(ierr); 993 ierr = PetscDTFactorial(ia + ib + 1, &gb);CHKERRQ(ierr); 994 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable."); 995 } 996 #endif 997 mu0 = PetscPowReal(2.,a + b + 1.) * ga * gb / gab; 998 999 #if defined(PETSCDTGAUSSIANQUADRATURE_EIG) 1000 { 1001 PetscReal *diag, *subdiag; 1002 PetscScalar *V; 1003 1004 ierr = PetscMalloc2(npoints, &diag, npoints, &subdiag);CHKERRQ(ierr); 1005 ierr = PetscMalloc1(npoints*npoints, &V);CHKERRQ(ierr); 1006 ierr = PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag);CHKERRQ(ierr); 1007 for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]); 1008 ierr = PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V);CHKERRQ(ierr); 1009 for (i = 0; i < npoints; i++) w[i] = PetscSqr(PetscRealPart(V[i * npoints])) * mu0; 1010 ierr = PetscFree(V);CHKERRQ(ierr); 1011 ierr = PetscFree2(diag, subdiag);CHKERRQ(ierr); 1012 } 1013 #else 1014 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found"); 1015 #endif 1016 { /* As of March 2, 2020, The Sun Performance Library breaks the LAPACK contract for xstegr and xsteqr: the 1017 eigenvalues are not guaranteed to be in ascending order. So we heave a passive aggressive sigh and check that 1018 the eigenvalues are sorted */ 1019 PetscBool sorted; 1020 1021 ierr = PetscSortedReal(npoints, x, &sorted);CHKERRQ(ierr); 1022 if (!sorted) { 1023 PetscInt *order, i; 1024 PetscReal *tmp; 1025 1026 ierr = PetscMalloc2(npoints, &order, npoints, &tmp);CHKERRQ(ierr); 1027 for (i = 0; i < npoints; i++) order[i] = i; 1028 ierr = PetscSortRealWithPermutation(npoints, x, order);CHKERRQ(ierr); 1029 ierr = PetscArraycpy(tmp, x, npoints);CHKERRQ(ierr); 1030 for (i = 0; i < npoints; i++) x[i] = tmp[order[i]]; 1031 ierr = PetscArraycpy(tmp, w, npoints);CHKERRQ(ierr); 1032 for (i = 0; i < npoints; i++) w[i] = tmp[order[i]]; 1033 ierr = PetscFree2(order, tmp);CHKERRQ(ierr); 1034 } 1035 } 1036 PetscFunctionReturn(0); 1037 } 1038 1039 static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton) 1040 { 1041 PetscErrorCode ierr; 1042 1043 PetscFunctionBegin; 1044 if (npoints < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive"); 1045 /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */ 1046 if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 1047 if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 1048 1049 if (newton) { 1050 ierr = PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr); 1051 } else { 1052 ierr = PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr); 1053 } 1054 if (alpha == beta) { /* symmetrize */ 1055 PetscInt i; 1056 for (i = 0; i < (npoints + 1) / 2; i++) { 1057 PetscInt j = npoints - 1 - i; 1058 PetscReal xi = x[i]; 1059 PetscReal xj = x[j]; 1060 PetscReal wi = w[i]; 1061 PetscReal wj = w[j]; 1062 1063 x[i] = (xi - xj) / 2.; 1064 x[j] = (xj - xi) / 2.; 1065 w[i] = w[j] = (wi + wj) / 2.; 1066 } 1067 } 1068 PetscFunctionReturn(0); 1069 } 1070 1071 /*@ 1072 PetscDTGaussJacobiQuadrature - quadrature for the interval [a, b] with the weight function 1073 $(x-a)^\alpha (x-b)^\beta$. 1074 1075 Not collective 1076 1077 Input Parameters: 1078 + npoints - the number of points in the quadrature rule 1079 . a - the left endpoint of the interval 1080 . b - the right endpoint of the interval 1081 . alpha - the left exponent 1082 - beta - the right exponent 1083 1084 Output Parameters: 1085 + x - array of length npoints, the locations of the quadrature points 1086 - w - array of length npoints, the weights of the quadrature points 1087 1088 Level: intermediate 1089 1090 Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 1. 1091 @*/ 1092 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[]) 1093 { 1094 PetscInt i; 1095 PetscErrorCode ierr; 1096 1097 PetscFunctionBegin; 1098 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr); 1099 if (a != -1. || b != 1.) { /* shift */ 1100 for (i = 0; i < npoints; i++) { 1101 x[i] = (x[i] + 1.) * ((b - a) / 2.) + a; 1102 w[i] *= (b - a) / 2.; 1103 } 1104 } 1105 PetscFunctionReturn(0); 1106 } 1107 1108 static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton) 1109 { 1110 PetscInt i; 1111 PetscErrorCode ierr; 1112 1113 PetscFunctionBegin; 1114 if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive"); 1115 /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */ 1116 if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1."); 1117 if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1."); 1118 1119 x[0] = -1.; 1120 x[npoints-1] = 1.; 1121 if (npoints > 2) { 1122 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints-2, alpha+1., beta+1., &x[1], &w[1], newton);CHKERRQ(ierr); 1123 } 1124 for (i = 1; i < npoints - 1; i++) { 1125 w[i] /= (1. - x[i]*x[i]); 1126 } 1127 ierr = PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints-1]);CHKERRQ(ierr); 1128 PetscFunctionReturn(0); 1129 } 1130 1131 /*@ 1132 PetscDTGaussLobattoJacobiQuadrature - quadrature for the interval [a, b] with the weight function 1133 $(x-a)^\alpha (x-b)^\beta$, with endpoints a and b included as quadrature points. 1134 1135 Not collective 1136 1137 Input Parameters: 1138 + npoints - the number of points in the quadrature rule 1139 . a - the left endpoint of the interval 1140 . b - the right endpoint of the interval 1141 . alpha - the left exponent 1142 - beta - the right exponent 1143 1144 Output Parameters: 1145 + x - array of length npoints, the locations of the quadrature points 1146 - w - array of length npoints, the weights of the quadrature points 1147 1148 Level: intermediate 1149 1150 Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 3. 1151 @*/ 1152 PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[]) 1153 { 1154 PetscInt i; 1155 PetscErrorCode ierr; 1156 1157 PetscFunctionBegin; 1158 ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr); 1159 if (a != -1. || b != 1.) { /* shift */ 1160 for (i = 0; i < npoints; i++) { 1161 x[i] = (x[i] + 1.) * ((b - a) / 2.) + a; 1162 w[i] *= (b - a) / 2.; 1163 } 1164 } 1165 PetscFunctionReturn(0); 1166 } 1167 1168 /*@ 1169 PetscDTGaussQuadrature - create Gauss-Legendre quadrature 1170 1171 Not Collective 1172 1173 Input Arguments: 1174 + npoints - number of points 1175 . a - left end of interval (often-1) 1176 - b - right end of interval (often +1) 1177 1178 Output Arguments: 1179 + x - quadrature points 1180 - w - quadrature weights 1181 1182 Level: intermediate 1183 1184 References: 1185 . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969. 1186 1187 .seealso: PetscDTLegendreEval() 1188 @*/ 1189 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 1190 { 1191 PetscInt i; 1192 PetscErrorCode ierr; 1193 1194 PetscFunctionBegin; 1195 ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr); 1196 if (a != -1. || b != 1.) { /* shift */ 1197 for (i = 0; i < npoints; i++) { 1198 x[i] = (x[i] + 1.) * ((b - a) / 2.) + a; 1199 w[i] *= (b - a) / 2.; 1200 } 1201 } 1202 PetscFunctionReturn(0); 1203 } 1204 1205 /*@C 1206 PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre 1207 nodes of a given size on the domain [-1,1] 1208 1209 Not Collective 1210 1211 Input Parameter: 1212 + n - number of grid nodes 1213 - type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON 1214 1215 Output Arguments: 1216 + x - quadrature points 1217 - w - quadrature weights 1218 1219 Notes: 1220 For n > 30 the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not 1221 close enough to the desired solution 1222 1223 These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes 1224 1225 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 1226 1227 Level: intermediate 1228 1229 .seealso: PetscDTGaussQuadrature() 1230 1231 @*/ 1232 PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w) 1233 { 1234 PetscBool newton; 1235 PetscErrorCode ierr; 1236 1237 PetscFunctionBegin; 1238 if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element"); 1239 newton = (PetscBool) (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON); 1240 ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton);CHKERRQ(ierr); 1241 PetscFunctionReturn(0); 1242 } 1243 1244 /*@ 1245 PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature 1246 1247 Not Collective 1248 1249 Input Arguments: 1250 + dim - The spatial dimension 1251 . Nc - The number of components 1252 . npoints - number of points in one dimension 1253 . a - left end of interval (often-1) 1254 - b - right end of interval (often +1) 1255 1256 Output Argument: 1257 . q - A PetscQuadrature object 1258 1259 Level: intermediate 1260 1261 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval() 1262 @*/ 1263 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 1264 { 1265 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c; 1266 PetscReal *x, *w, *xw, *ww; 1267 PetscErrorCode ierr; 1268 1269 PetscFunctionBegin; 1270 ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr); 1271 ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr); 1272 /* Set up the Golub-Welsch system */ 1273 switch (dim) { 1274 case 0: 1275 ierr = PetscFree(x);CHKERRQ(ierr); 1276 ierr = PetscFree(w);CHKERRQ(ierr); 1277 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 1278 ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 1279 x[0] = 0.0; 1280 for (c = 0; c < Nc; ++c) w[c] = 1.0; 1281 break; 1282 case 1: 1283 ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr); 1284 ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr); 1285 for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i]; 1286 ierr = PetscFree(ww);CHKERRQ(ierr); 1287 break; 1288 case 2: 1289 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 1290 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 1291 for (i = 0; i < npoints; ++i) { 1292 for (j = 0; j < npoints; ++j) { 1293 x[(i*npoints+j)*dim+0] = xw[i]; 1294 x[(i*npoints+j)*dim+1] = xw[j]; 1295 for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j]; 1296 } 1297 } 1298 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 1299 break; 1300 case 3: 1301 ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr); 1302 ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr); 1303 for (i = 0; i < npoints; ++i) { 1304 for (j = 0; j < npoints; ++j) { 1305 for (k = 0; k < npoints; ++k) { 1306 x[((i*npoints+j)*npoints+k)*dim+0] = xw[i]; 1307 x[((i*npoints+j)*npoints+k)*dim+1] = xw[j]; 1308 x[((i*npoints+j)*npoints+k)*dim+2] = xw[k]; 1309 for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k]; 1310 } 1311 } 1312 } 1313 ierr = PetscFree2(xw,ww);CHKERRQ(ierr); 1314 break; 1315 default: 1316 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 1317 } 1318 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1319 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 1320 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 1321 ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr); 1322 PetscFunctionReturn(0); 1323 } 1324 1325 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */ 1326 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta) 1327 { 1328 PetscFunctionBegin; 1329 *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0; 1330 *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0; 1331 *zeta = z; 1332 PetscFunctionReturn(0); 1333 } 1334 1335 1336 /*@ 1337 PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex 1338 1339 Not Collective 1340 1341 Input Arguments: 1342 + dim - The simplex dimension 1343 . Nc - The number of components 1344 . npoints - The number of points in one dimension 1345 . a - left end of interval (often-1) 1346 - b - right end of interval (often +1) 1347 1348 Output Argument: 1349 . q - A PetscQuadrature object 1350 1351 Level: intermediate 1352 1353 References: 1354 . 1. - Karniadakis and Sherwin. FIAT 1355 1356 Note: For dim == 1, this is Gauss-Legendre quadrature 1357 1358 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature() 1359 @*/ 1360 PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q) 1361 { 1362 PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints; 1363 PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; 1364 PetscInt i, j, k, c; PetscErrorCode ierr; 1365 1366 PetscFunctionBegin; 1367 if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); 1368 ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr); 1369 ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr); 1370 switch (dim) { 1371 case 0: 1372 ierr = PetscFree(x);CHKERRQ(ierr); 1373 ierr = PetscFree(w);CHKERRQ(ierr); 1374 ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); 1375 ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr); 1376 x[0] = 0.0; 1377 for (c = 0; c < Nc; ++c) w[c] = 1.0; 1378 break; 1379 case 1: 1380 ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr); 1381 ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, x, wx);CHKERRQ(ierr); 1382 for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i]; 1383 ierr = PetscFree(wx);CHKERRQ(ierr); 1384 break; 1385 case 2: 1386 ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr); 1387 ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, px, wx);CHKERRQ(ierr); 1388 ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 1.0, 0.0, py, wy);CHKERRQ(ierr); 1389 for (i = 0; i < npoints; ++i) { 1390 for (j = 0; j < npoints; ++j) { 1391 ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); 1392 for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j]; 1393 } 1394 } 1395 ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); 1396 break; 1397 case 3: 1398 ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr); 1399 ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, px, wx);CHKERRQ(ierr); 1400 ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 1.0, 0.0, py, wy);CHKERRQ(ierr); 1401 ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 2.0, 0.0, pz, wz);CHKERRQ(ierr); 1402 for (i = 0; i < npoints; ++i) { 1403 for (j = 0; j < npoints; ++j) { 1404 for (k = 0; k < npoints; ++k) { 1405 ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*npoints+j)*npoints+k)*3+0], &x[((i*npoints+j)*npoints+k)*3+1], &x[((i*npoints+j)*npoints+k)*3+2]);CHKERRQ(ierr); 1406 for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k]; 1407 } 1408 } 1409 } 1410 ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); 1411 break; 1412 default: 1413 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); 1414 } 1415 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1416 ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr); 1417 ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr); 1418 ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr); 1419 PetscFunctionReturn(0); 1420 } 1421 1422 /*@ 1423 PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell 1424 1425 Not Collective 1426 1427 Input Arguments: 1428 + dim - The cell dimension 1429 . level - The number of points in one dimension, 2^l 1430 . a - left end of interval (often-1) 1431 - b - right end of interval (often +1) 1432 1433 Output Argument: 1434 . q - A PetscQuadrature object 1435 1436 Level: intermediate 1437 1438 .seealso: PetscDTGaussTensorQuadrature() 1439 @*/ 1440 PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q) 1441 { 1442 const PetscInt p = 16; /* Digits of precision in the evaluation */ 1443 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1444 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1445 const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */ 1446 PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */ 1447 PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */ 1448 PetscReal *x, *w; 1449 PetscInt K, k, npoints; 1450 PetscErrorCode ierr; 1451 1452 PetscFunctionBegin; 1453 if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim); 1454 if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits"); 1455 /* Find K such that the weights are < 32 digits of precision */ 1456 for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) { 1457 wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h))); 1458 } 1459 ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); 1460 ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr); 1461 npoints = 2*K-1; 1462 ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); 1463 ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); 1464 /* Center term */ 1465 x[0] = beta; 1466 w[0] = 0.5*alpha*PETSC_PI; 1467 for (k = 1; k < K; ++k) { 1468 wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1469 xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h)); 1470 x[2*k-1] = -alpha*xk+beta; 1471 w[2*k-1] = wk; 1472 x[2*k+0] = alpha*xk+beta; 1473 w[2*k+0] = wk; 1474 } 1475 ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr); 1476 PetscFunctionReturn(0); 1477 } 1478 1479 PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1480 { 1481 const PetscInt p = 16; /* Digits of precision in the evaluation */ 1482 const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */ 1483 const PetscReal beta = (b+a)/2.; /* Center of the integration interval */ 1484 PetscReal h = 1.0; /* Step size, length between x_k */ 1485 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 1486 PetscReal osum = 0.0; /* Integral on last level */ 1487 PetscReal psum = 0.0; /* Integral on the level before the last level */ 1488 PetscReal sum; /* Integral on current level */ 1489 PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 1490 PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 1491 PetscReal wk; /* Quadrature weight at x_k */ 1492 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 1493 PetscInt d; /* Digits of precision in the integral */ 1494 1495 PetscFunctionBegin; 1496 if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 1497 /* Center term */ 1498 func(beta, &lval); 1499 sum = 0.5*alpha*PETSC_PI*lval; 1500 /* */ 1501 do { 1502 PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4; 1503 PetscInt k = 1; 1504 1505 ++l; 1506 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 1507 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 1508 psum = osum; 1509 osum = sum; 1510 h *= 0.5; 1511 sum *= 0.5; 1512 do { 1513 wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1514 yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h))); 1515 lx = -alpha*(1.0 - yk)+beta; 1516 rx = alpha*(1.0 - yk)+beta; 1517 func(lx, &lval); 1518 func(rx, &rval); 1519 lterm = alpha*wk*lval; 1520 maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm); 1521 sum += lterm; 1522 rterm = alpha*wk*rval; 1523 maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm); 1524 sum += rterm; 1525 ++k; 1526 /* Only need to evaluate every other point on refined levels */ 1527 if (l != 1) ++k; 1528 } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1529 1530 d1 = PetscLog10Real(PetscAbsReal(sum - osum)); 1531 d2 = PetscLog10Real(PetscAbsReal(sum - psum)); 1532 d3 = PetscLog10Real(maxTerm) - p; 1533 if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0; 1534 else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm))); 1535 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1536 } while (d < digits && l < 12); 1537 *sol = sum; 1538 1539 PetscFunctionReturn(0); 1540 } 1541 1542 #if defined(PETSC_HAVE_MPFR) 1543 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1544 { 1545 const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */ 1546 PetscInt l = 0; /* Level of refinement, h = 2^{-l} */ 1547 mpfr_t alpha; /* Half-width of the integration interval */ 1548 mpfr_t beta; /* Center of the integration interval */ 1549 mpfr_t h; /* Step size, length between x_k */ 1550 mpfr_t osum; /* Integral on last level */ 1551 mpfr_t psum; /* Integral on the level before the last level */ 1552 mpfr_t sum; /* Integral on current level */ 1553 mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */ 1554 mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */ 1555 mpfr_t wk; /* Quadrature weight at x_k */ 1556 PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */ 1557 PetscInt d; /* Digits of precision in the integral */ 1558 mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp; 1559 1560 PetscFunctionBegin; 1561 if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits"); 1562 /* Create high precision storage */ 1563 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); 1564 /* Initialization */ 1565 mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN); 1566 mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN); 1567 mpfr_set_d(osum, 0.0, MPFR_RNDN); 1568 mpfr_set_d(psum, 0.0, MPFR_RNDN); 1569 mpfr_set_d(h, 1.0, MPFR_RNDN); 1570 mpfr_const_pi(pi2, MPFR_RNDN); 1571 mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN); 1572 /* Center term */ 1573 func(0.5*(b+a), &lval); 1574 mpfr_set(sum, pi2, MPFR_RNDN); 1575 mpfr_mul(sum, sum, alpha, MPFR_RNDN); 1576 mpfr_mul_d(sum, sum, lval, MPFR_RNDN); 1577 /* */ 1578 do { 1579 PetscReal d1, d2, d3, d4; 1580 PetscInt k = 1; 1581 1582 ++l; 1583 mpfr_set_d(maxTerm, 0.0, MPFR_RNDN); 1584 /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */ 1585 /* At each level of refinement, h --> h/2 and sum --> sum/2 */ 1586 mpfr_set(psum, osum, MPFR_RNDN); 1587 mpfr_set(osum, sum, MPFR_RNDN); 1588 mpfr_mul_d(h, h, 0.5, MPFR_RNDN); 1589 mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN); 1590 do { 1591 mpfr_set_si(kh, k, MPFR_RNDN); 1592 mpfr_mul(kh, kh, h, MPFR_RNDN); 1593 /* Weight */ 1594 mpfr_set(wk, h, MPFR_RNDN); 1595 mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN); 1596 mpfr_mul(msinh, msinh, pi2, MPFR_RNDN); 1597 mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN); 1598 mpfr_cosh(tmp, msinh, MPFR_RNDN); 1599 mpfr_sqr(tmp, tmp, MPFR_RNDN); 1600 mpfr_mul(wk, wk, mcosh, MPFR_RNDN); 1601 mpfr_div(wk, wk, tmp, MPFR_RNDN); 1602 /* Abscissa */ 1603 mpfr_set_d(yk, 1.0, MPFR_RNDZ); 1604 mpfr_cosh(tmp, msinh, MPFR_RNDN); 1605 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 1606 mpfr_exp(tmp, msinh, MPFR_RNDN); 1607 mpfr_div(yk, yk, tmp, MPFR_RNDZ); 1608 /* Quadrature points */ 1609 mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ); 1610 mpfr_mul(lx, lx, alpha, MPFR_RNDU); 1611 mpfr_add(lx, lx, beta, MPFR_RNDU); 1612 mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ); 1613 mpfr_mul(rx, rx, alpha, MPFR_RNDD); 1614 mpfr_add(rx, rx, beta, MPFR_RNDD); 1615 /* Evaluation */ 1616 func(mpfr_get_d(lx, MPFR_RNDU), &lval); 1617 func(mpfr_get_d(rx, MPFR_RNDD), &rval); 1618 /* Update */ 1619 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 1620 mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN); 1621 mpfr_add(sum, sum, tmp, MPFR_RNDN); 1622 mpfr_abs(tmp, tmp, MPFR_RNDN); 1623 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 1624 mpfr_set(curTerm, tmp, MPFR_RNDN); 1625 mpfr_mul(tmp, wk, alpha, MPFR_RNDN); 1626 mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN); 1627 mpfr_add(sum, sum, tmp, MPFR_RNDN); 1628 mpfr_abs(tmp, tmp, MPFR_RNDN); 1629 mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN); 1630 mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN); 1631 ++k; 1632 /* Only need to evaluate every other point on refined levels */ 1633 if (l != 1) ++k; 1634 mpfr_log10(tmp, wk, MPFR_RNDN); 1635 mpfr_abs(tmp, tmp, MPFR_RNDN); 1636 } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */ 1637 mpfr_sub(tmp, sum, osum, MPFR_RNDN); 1638 mpfr_abs(tmp, tmp, MPFR_RNDN); 1639 mpfr_log10(tmp, tmp, MPFR_RNDN); 1640 d1 = mpfr_get_d(tmp, MPFR_RNDN); 1641 mpfr_sub(tmp, sum, psum, MPFR_RNDN); 1642 mpfr_abs(tmp, tmp, MPFR_RNDN); 1643 mpfr_log10(tmp, tmp, MPFR_RNDN); 1644 d2 = mpfr_get_d(tmp, MPFR_RNDN); 1645 mpfr_log10(tmp, maxTerm, MPFR_RNDN); 1646 d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits; 1647 mpfr_log10(tmp, curTerm, MPFR_RNDN); 1648 d4 = mpfr_get_d(tmp, MPFR_RNDN); 1649 d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4))); 1650 } while (d < digits && l < 8); 1651 *sol = mpfr_get_d(sum, MPFR_RNDN); 1652 /* Cleanup */ 1653 mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL); 1654 PetscFunctionReturn(0); 1655 } 1656 #else 1657 1658 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol) 1659 { 1660 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp"); 1661 } 1662 #endif 1663 1664 /* Overwrites A. Can only handle full-rank problems with m>=n 1665 * A in column-major format 1666 * Ainv in row-major format 1667 * tau has length m 1668 * worksize must be >= max(1,n) 1669 */ 1670 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1671 { 1672 PetscErrorCode ierr; 1673 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1674 PetscScalar *A,*Ainv,*R,*Q,Alpha; 1675 1676 PetscFunctionBegin; 1677 #if defined(PETSC_USE_COMPLEX) 1678 { 1679 PetscInt i,j; 1680 ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); 1681 for (j=0; j<n; j++) { 1682 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 1683 } 1684 mstride = m; 1685 } 1686 #else 1687 A = A_in; 1688 Ainv = Ainv_out; 1689 #endif 1690 1691 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 1692 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 1693 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 1694 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 1695 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 1696 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1697 ierr = PetscFPTrapPop();CHKERRQ(ierr); 1698 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1699 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1700 1701 /* Extract an explicit representation of Q */ 1702 Q = Ainv; 1703 ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr); 1704 K = N; /* full rank */ 1705 PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1706 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1707 1708 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1709 Alpha = 1.0; 1710 ldb = lda; 1711 PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb)); 1712 /* Ainv is Q, overwritten with inverse */ 1713 1714 #if defined(PETSC_USE_COMPLEX) 1715 { 1716 PetscInt i; 1717 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 1718 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 1719 } 1720 #endif 1721 PetscFunctionReturn(0); 1722 } 1723 1724 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 1725 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 1726 { 1727 PetscErrorCode ierr; 1728 PetscReal *Bv; 1729 PetscInt i,j; 1730 1731 PetscFunctionBegin; 1732 ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr); 1733 /* Point evaluation of L_p on all the source vertices */ 1734 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 1735 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 1736 for (i=0; i<ninterval; i++) { 1737 for (j=0; j<ndegree; j++) { 1738 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1739 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 1740 } 1741 } 1742 ierr = PetscFree(Bv);CHKERRQ(ierr); 1743 PetscFunctionReturn(0); 1744 } 1745 1746 /*@ 1747 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 1748 1749 Not Collective 1750 1751 Input Arguments: 1752 + degree - degree of reconstruction polynomial 1753 . nsource - number of source intervals 1754 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 1755 . ntarget - number of target intervals 1756 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 1757 1758 Output Arguments: 1759 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 1760 1761 Level: advanced 1762 1763 .seealso: PetscDTLegendreEval() 1764 @*/ 1765 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 1766 { 1767 PetscErrorCode ierr; 1768 PetscInt i,j,k,*bdegrees,worksize; 1769 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 1770 PetscScalar *tau,*work; 1771 1772 PetscFunctionBegin; 1773 PetscValidRealPointer(sourcex,3); 1774 PetscValidRealPointer(targetx,5); 1775 PetscValidRealPointer(R,6); 1776 if (degree >= nsource) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Reconstruction degree %D must be less than number of source intervals %D",degree,nsource); 1777 #if defined(PETSC_USE_DEBUG) 1778 for (i=0; i<nsource; i++) { 1779 if (sourcex[i] >= sourcex[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Source interval %D has negative orientation (%g,%g)",i,(double)sourcex[i],(double)sourcex[i+1]); 1780 } 1781 for (i=0; i<ntarget; i++) { 1782 if (targetx[i] >= targetx[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Target interval %D has negative orientation (%g,%g)",i,(double)targetx[i],(double)targetx[i+1]); 1783 } 1784 #endif 1785 xmin = PetscMin(sourcex[0],targetx[0]); 1786 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 1787 center = (xmin + xmax)/2; 1788 hscale = (xmax - xmin)/2; 1789 worksize = nsource; 1790 ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); 1791 ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); 1792 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 1793 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 1794 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 1795 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 1796 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 1797 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 1798 for (i=0; i<ntarget; i++) { 1799 PetscReal rowsum = 0; 1800 for (j=0; j<nsource; j++) { 1801 PetscReal sum = 0; 1802 for (k=0; k<degree+1; k++) { 1803 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 1804 } 1805 R[i*nsource+j] = sum; 1806 rowsum += sum; 1807 } 1808 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 1809 } 1810 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 1811 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 1812 PetscFunctionReturn(0); 1813 } 1814 1815 /*@C 1816 PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points 1817 1818 Not Collective 1819 1820 Input Parameter: 1821 + n - the number of GLL nodes 1822 . nodes - the GLL nodes 1823 . weights - the GLL weights 1824 . f - the function values at the nodes 1825 1826 Output Parameter: 1827 . in - the value of the integral 1828 1829 Level: beginner 1830 1831 .seealso: PetscDTGaussLobattoLegendreQuadrature() 1832 1833 @*/ 1834 PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in) 1835 { 1836 PetscInt i; 1837 1838 PetscFunctionBegin; 1839 *in = 0.; 1840 for (i=0; i<n; i++) { 1841 *in += f[i]*f[i]*weights[i]; 1842 } 1843 PetscFunctionReturn(0); 1844 } 1845 1846 /*@C 1847 PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element 1848 1849 Not Collective 1850 1851 Input Parameter: 1852 + n - the number of GLL nodes 1853 . nodes - the GLL nodes 1854 . weights - the GLL weights 1855 1856 Output Parameter: 1857 . A - the stiffness element 1858 1859 Level: beginner 1860 1861 Notes: 1862 Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy() 1863 1864 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) 1865 1866 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1867 1868 @*/ 1869 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1870 { 1871 PetscReal **A; 1872 PetscErrorCode ierr; 1873 const PetscReal *gllnodes = nodes; 1874 const PetscInt p = n-1; 1875 PetscReal z0,z1,z2 = -1,x,Lpj,Lpr; 1876 PetscInt i,j,nn,r; 1877 1878 PetscFunctionBegin; 1879 ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 1880 ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 1881 for (i=1; i<n; i++) A[i] = A[i-1]+n; 1882 1883 for (j=1; j<p; j++) { 1884 x = gllnodes[j]; 1885 z0 = 1.; 1886 z1 = x; 1887 for (nn=1; nn<p; nn++) { 1888 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1889 z0 = z1; 1890 z1 = z2; 1891 } 1892 Lpj=z2; 1893 for (r=1; r<p; r++) { 1894 if (r == j) { 1895 A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj); 1896 } else { 1897 x = gllnodes[r]; 1898 z0 = 1.; 1899 z1 = x; 1900 for (nn=1; nn<p; nn++) { 1901 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1902 z0 = z1; 1903 z1 = z2; 1904 } 1905 Lpr = z2; 1906 A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r])); 1907 } 1908 } 1909 } 1910 for (j=1; j<p+1; j++) { 1911 x = gllnodes[j]; 1912 z0 = 1.; 1913 z1 = x; 1914 for (nn=1; nn<p; nn++) { 1915 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1916 z0 = z1; 1917 z1 = z2; 1918 } 1919 Lpj = z2; 1920 A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j])); 1921 A[0][j] = A[j][0]; 1922 } 1923 for (j=0; j<p; j++) { 1924 x = gllnodes[j]; 1925 z0 = 1.; 1926 z1 = x; 1927 for (nn=1; nn<p; nn++) { 1928 z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.)); 1929 z0 = z1; 1930 z1 = z2; 1931 } 1932 Lpj=z2; 1933 1934 A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j])); 1935 A[j][p] = A[p][j]; 1936 } 1937 A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.; 1938 A[p][p]=A[0][0]; 1939 *AA = A; 1940 PetscFunctionReturn(0); 1941 } 1942 1943 /*@C 1944 PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element 1945 1946 Not Collective 1947 1948 Input Parameter: 1949 + n - the number of GLL nodes 1950 . nodes - the GLL nodes 1951 . weights - the GLL weightss 1952 - A - the stiffness element 1953 1954 Level: beginner 1955 1956 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate() 1957 1958 @*/ 1959 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 1960 { 1961 PetscErrorCode ierr; 1962 1963 PetscFunctionBegin; 1964 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 1965 ierr = PetscFree(*AA);CHKERRQ(ierr); 1966 *AA = NULL; 1967 PetscFunctionReturn(0); 1968 } 1969 1970 /*@C 1971 PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element 1972 1973 Not Collective 1974 1975 Input Parameter: 1976 + n - the number of GLL nodes 1977 . nodes - the GLL nodes 1978 . weights - the GLL weights 1979 1980 Output Parameter: 1981 . AA - the stiffness element 1982 - AAT - the transpose of AA (pass in NULL if you do not need this array) 1983 1984 Level: beginner 1985 1986 Notes: 1987 Destroy this with PetscGaussLobattoLegendreElementGradientDestroy() 1988 1989 You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented 1990 1991 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy() 1992 1993 @*/ 1994 PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 1995 { 1996 PetscReal **A, **AT = NULL; 1997 PetscErrorCode ierr; 1998 const PetscReal *gllnodes = nodes; 1999 const PetscInt p = n-1; 2000 PetscReal Li, Lj,d0; 2001 PetscInt i,j; 2002 2003 PetscFunctionBegin; 2004 ierr = PetscMalloc1(n,&A);CHKERRQ(ierr); 2005 ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr); 2006 for (i=1; i<n; i++) A[i] = A[i-1]+n; 2007 2008 if (AAT) { 2009 ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr); 2010 ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr); 2011 for (i=1; i<n; i++) AT[i] = AT[i-1]+n; 2012 } 2013 2014 if (n==1) {A[0][0] = 0.;} 2015 d0 = (PetscReal)p*((PetscReal)p+1.)/4.; 2016 for (i=0; i<n; i++) { 2017 for (j=0; j<n; j++) { 2018 A[i][j] = 0.; 2019 ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li);CHKERRQ(ierr); 2020 ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj);CHKERRQ(ierr); 2021 if (i!=j) A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j])); 2022 if ((j==i) && (i==0)) A[i][j] = -d0; 2023 if (j==i && i==p) A[i][j] = d0; 2024 if (AT) AT[j][i] = A[i][j]; 2025 } 2026 } 2027 if (AAT) *AAT = AT; 2028 *AA = A; 2029 PetscFunctionReturn(0); 2030 } 2031 2032 /*@C 2033 PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate() 2034 2035 Not Collective 2036 2037 Input Parameter: 2038 + n - the number of GLL nodes 2039 . nodes - the GLL nodes 2040 . weights - the GLL weights 2041 . AA - the stiffness element 2042 - AAT - the transpose of the element 2043 2044 Level: beginner 2045 2046 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate() 2047 2048 @*/ 2049 PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT) 2050 { 2051 PetscErrorCode ierr; 2052 2053 PetscFunctionBegin; 2054 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2055 ierr = PetscFree(*AA);CHKERRQ(ierr); 2056 *AA = NULL; 2057 if (*AAT) { 2058 ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr); 2059 ierr = PetscFree(*AAT);CHKERRQ(ierr); 2060 *AAT = NULL; 2061 } 2062 PetscFunctionReturn(0); 2063 } 2064 2065 /*@C 2066 PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element 2067 2068 Not Collective 2069 2070 Input Parameter: 2071 + n - the number of GLL nodes 2072 . nodes - the GLL nodes 2073 . weights - the GLL weightss 2074 2075 Output Parameter: 2076 . AA - the stiffness element 2077 2078 Level: beginner 2079 2080 Notes: 2081 Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy() 2082 2083 This is the same as the Gradient operator multiplied by the diagonal mass matrix 2084 2085 You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented 2086 2087 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy() 2088 2089 @*/ 2090 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2091 { 2092 PetscReal **D; 2093 PetscErrorCode ierr; 2094 const PetscReal *gllweights = weights; 2095 const PetscInt glln = n; 2096 PetscInt i,j; 2097 2098 PetscFunctionBegin; 2099 ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr); 2100 for (i=0; i<glln; i++){ 2101 for (j=0; j<glln; j++) { 2102 D[i][j] = gllweights[i]*D[i][j]; 2103 } 2104 } 2105 *AA = D; 2106 PetscFunctionReturn(0); 2107 } 2108 2109 /*@C 2110 PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element 2111 2112 Not Collective 2113 2114 Input Parameter: 2115 + n - the number of GLL nodes 2116 . nodes - the GLL nodes 2117 . weights - the GLL weights 2118 - A - advection 2119 2120 Level: beginner 2121 2122 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate() 2123 2124 @*/ 2125 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2126 { 2127 PetscErrorCode ierr; 2128 2129 PetscFunctionBegin; 2130 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2131 ierr = PetscFree(*AA);CHKERRQ(ierr); 2132 *AA = NULL; 2133 PetscFunctionReturn(0); 2134 } 2135 2136 PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2137 { 2138 PetscReal **A; 2139 PetscErrorCode ierr; 2140 const PetscReal *gllweights = weights; 2141 const PetscInt glln = n; 2142 PetscInt i,j; 2143 2144 PetscFunctionBegin; 2145 ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr); 2146 ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr); 2147 for (i=1; i<glln; i++) A[i] = A[i-1]+glln; 2148 if (glln==1) {A[0][0] = 0.;} 2149 for (i=0; i<glln; i++) { 2150 for (j=0; j<glln; j++) { 2151 A[i][j] = 0.; 2152 if (j==i) A[i][j] = gllweights[i]; 2153 } 2154 } 2155 *AA = A; 2156 PetscFunctionReturn(0); 2157 } 2158 2159 PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA) 2160 { 2161 PetscErrorCode ierr; 2162 2163 PetscFunctionBegin; 2164 ierr = PetscFree((*AA)[0]);CHKERRQ(ierr); 2165 ierr = PetscFree(*AA);CHKERRQ(ierr); 2166 *AA = NULL; 2167 PetscFunctionReturn(0); 2168 } 2169