/* Discretization tools */ #include #if defined(PETSC_HAVE_MATHIMF_H) #include /* this needs to be included before math.h */ #endif #include /*I "petscdt.h" I*/ #include #include #include #include #include #include #undef __FUNCT__ #define __FUNCT__ "PetscQuadratureCreate" PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidPointer(q, 2); ierr = DMInitializePackage();CHKERRQ(ierr); ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr); (*q)->dim = -1; (*q)->numPoints = 0; (*q)->points = NULL; (*q)->weights = NULL; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscQuadratureDestroy" PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q) { PetscErrorCode ierr; PetscFunctionBegin; if (!*q) PetscFunctionReturn(0); PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1); if (--((PetscObject)(*q))->refct > 0) { *q = NULL; PetscFunctionReturn(0); } ierr = PetscFree((*q)->points);CHKERRQ(ierr); ierr = PetscFree((*q)->weights);CHKERRQ(ierr); ierr = PetscHeaderDestroy(q);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscQuadratureGetData" PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[]) { PetscFunctionBegin; PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); if (dim) { PetscValidPointer(dim, 2); *dim = q->dim; } if (npoints) { PetscValidPointer(npoints, 3); *npoints = q->numPoints; } if (points) { PetscValidPointer(points, 4); *points = q->points; } if (weights) { PetscValidPointer(weights, 5); *weights = q->weights; } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscQuadratureSetData" PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[]) { PetscFunctionBegin; PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1); if (dim >= 0) q->dim = dim; if (npoints >= 0) q->numPoints = npoints; if (points) { PetscValidPointer(points, 4); q->points = points; } if (weights) { PetscValidPointer(weights, 5); q->weights = weights; } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscQuadratureView" PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer) { PetscInt q, d; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n (", quad->numPoints);CHKERRQ(ierr); for (q = 0; q < quad->numPoints; ++q) { for (d = 0; d < quad->dim; ++d) { if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr); } ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscDTLegendreEval" /*@ PetscDTLegendreEval - evaluate Legendre polynomial at points Not Collective Input Arguments: + npoints - number of spatial points to evaluate at . points - array of locations to evaluate at . ndegree - number of basis degrees to evaluate - degrees - sorted array of degrees to evaluate Output Arguments: + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) . D - row-oriented derivative evaluation matrix (or NULL) - D2 - row-oriented second derivative evaluation matrix (or NULL) Level: intermediate .seealso: PetscDTGaussQuadrature() @*/ PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) { PetscInt i,maxdegree; PetscFunctionBegin; if (!npoints || !ndegree) PetscFunctionReturn(0); maxdegree = degrees[ndegree-1]; for (i=0; i 0) r = 0.5 * (r + x[k-1]); for (j = 0; j < maxIter; ++j) { PetscReal s = 0.0, delta, f, fp; PetscInt i; for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]); ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr); ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr); delta = f / (fp - f * s); r = r - delta; if (PetscAbs(delta) < eps) break; } x[k] = r; ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr); w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscDTGaussJacobiQuadrature" /*@C PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex Not Collective Input Arguments: + dim - The simplex dimension . order - The quadrature order . a - left end of interval (often-1) - b - right end of interval (often +1) Output Arguments: . q - A PetscQuadrature object Level: intermediate References: Karniadakis and Sherwin. FIAT .seealso: PetscDTGaussQuadrature() @*/ PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q) { PetscInt npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order; PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w; PetscInt i, j, k; PetscErrorCode ierr; PetscFunctionBegin; if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now"); ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr); ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr); switch (dim) { case 0: ierr = PetscFree(x);CHKERRQ(ierr); ierr = PetscFree(w);CHKERRQ(ierr); ierr = PetscMalloc1(1, &x);CHKERRQ(ierr); ierr = PetscMalloc1(1, &w);CHKERRQ(ierr); x[0] = 0.0; w[0] = 1.0; break; case 1: ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr); break; case 2: ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); for (i = 0; i < order; ++i) { for (j = 0; j < order; ++j) { ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr); w[i*order+j] = 0.5 * wx[i] * wy[j]; } } ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); break; case 3: ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr); for (i = 0; i < order; ++i) { for (j = 0; j < order; ++j) { for (k = 0; k < order; ++k) { ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*order+j)*order+k)*3+0], &x[((i*order+j)*order+k)*3+1], &x[((i*order+j)*order+k)*3+2]);CHKERRQ(ierr); w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k]; } } } ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr); break; default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim); } ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr); ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PetscDTPseudoInverseQR" /* Overwrites A. Can only handle full-rank problems with m>=n * A in column-major format * Ainv in row-major format * tau has length m * worksize must be >= max(1,n) */ static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) { PetscErrorCode ierr; PetscBLASInt M,N,K,lda,ldb,ldwork,info; PetscScalar *A,*Ainv,*R,*Q,Alpha; PetscFunctionBegin; #if defined(PETSC_USE_COMPLEX) { PetscInt i,j; ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr); for (j=0; j= nsource) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Reconstruction degree %D must be less than number of source intervals %D",degree,nsource); #if defined(PETSC_USE_DEBUG) for (i=0; 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]); } for (i=0; 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]); } #endif xmin = PetscMin(sourcex[0],targetx[0]); xmax = PetscMax(sourcex[nsource],targetx[ntarget]); center = (xmin + xmax)/2; hscale = (xmax - xmin)/2; worksize = nsource; ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr); ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr); for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; for (i=0; i<=degree; i++) bdegrees[i] = i+1; ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); for (i=0; i