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