/* Discretization tools */ #include /*I "petscdt.h" I*/ #include #include #include #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 (fabs(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 . npoints - number of points . a - left end of interval (often-1) - b - right end of interval (often +1) Output Arguments: + points - quadrature points - weights - quadrature weights Level: intermediate References: Karniadakis and Sherwin. FIAT .seealso: PetscDTGaussQuadrature() @*/ PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscReal *points[], PetscReal *weights[]) { 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"); switch (dim) { case 1: ierr = PetscMalloc(npoints * sizeof(PetscReal), &x);CHKERRQ(ierr); ierr = PetscMalloc(npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, w);CHKERRQ(ierr); break; case 2: ierr = PetscMalloc(npoints*npoints*2 * sizeof(PetscReal), &x);CHKERRQ(ierr); ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); ierr = PetscMalloc4(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); for (i = 0; i < npoints; ++i) { for (j = 0; j < npoints; ++j) { ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr); w[i*npoints+j] = 0.5 * wx[i] * wy[j]; } } ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr); break; case 3: ierr = PetscMalloc(npoints*npoints*3 * sizeof(PetscReal), &x);CHKERRQ(ierr); ierr = PetscMalloc(npoints*npoints * sizeof(PetscReal), &w);CHKERRQ(ierr); ierr = PetscMalloc6(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy,npoints,PetscReal,&pz,npoints,PetscReal,&wz);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr); ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr); for (i = 0; i < npoints; ++i) { for (j = 0; j < npoints; ++j) { for (k = 0; k < npoints; ++k) { 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); w[(i*npoints+j)*npoints+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); } if (points) *points = x; if (weights) *weights = w; 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,PetscScalar,&A,m*n,PetscScalar,&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,sourcex[i],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,targetx[i],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,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr); ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&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