137045ce4SJed Brown /* Discretization tools */ 237045ce4SJed Brown 337045ce4SJed Brown #include <petscdt.h> /*I "petscdt.h" I*/ 437045ce4SJed Brown #include <petscblaslapack.h> 5*194825f6SJed Brown #include <petsc-private/petscimpl.h> 637045ce4SJed Brown 737045ce4SJed Brown #undef __FUNCT__ 837045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval" 937045ce4SJed Brown /*@ 1037045ce4SJed Brown PetscDTLegendreEval - evaluate Legendre polynomial at points 1137045ce4SJed Brown 1237045ce4SJed Brown Not Collective 1337045ce4SJed Brown 1437045ce4SJed Brown Input Arguments: 1537045ce4SJed Brown + npoints - number of spatial points to evaluate at 1637045ce4SJed Brown . points - array of locations to evaluate at 1737045ce4SJed Brown . ndegree - number of basis degrees to evaluate 1837045ce4SJed Brown - degrees - sorted array of degrees to evaluate 1937045ce4SJed Brown 2037045ce4SJed Brown Output Arguments: 210298fd71SBarry Smith + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 220298fd71SBarry Smith . D - row-oriented derivative evaluation matrix (or NULL) 230298fd71SBarry Smith - D2 - row-oriented second derivative evaluation matrix (or NULL) 2437045ce4SJed Brown 2537045ce4SJed Brown Level: intermediate 2637045ce4SJed Brown 2737045ce4SJed Brown .seealso: PetscDTGaussQuadrature() 2837045ce4SJed Brown @*/ 2937045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 3037045ce4SJed Brown { 3137045ce4SJed Brown PetscInt i,maxdegree; 3237045ce4SJed Brown 3337045ce4SJed Brown PetscFunctionBegin; 3437045ce4SJed Brown if (!npoints || !ndegree) PetscFunctionReturn(0); 3537045ce4SJed Brown maxdegree = degrees[ndegree-1]; 3637045ce4SJed Brown for (i=0; i<npoints; i++) { 3737045ce4SJed Brown PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 3837045ce4SJed Brown PetscInt j,k; 3937045ce4SJed Brown x = points[i]; 4037045ce4SJed Brown pm2 = 0; 4137045ce4SJed Brown pm1 = 1; 4237045ce4SJed Brown pd2 = 0; 4337045ce4SJed Brown pd1 = 0; 4437045ce4SJed Brown pdd2 = 0; 4537045ce4SJed Brown pdd1 = 0; 4637045ce4SJed Brown k = 0; 4737045ce4SJed Brown if (degrees[k] == 0) { 4837045ce4SJed Brown if (B) B[i*ndegree+k] = pm1; 4937045ce4SJed Brown if (D) D[i*ndegree+k] = pd1; 5037045ce4SJed Brown if (D2) D2[i*ndegree+k] = pdd1; 5137045ce4SJed Brown k++; 5237045ce4SJed Brown } 5337045ce4SJed Brown for (j=1; j<=maxdegree; j++,k++) { 5437045ce4SJed Brown PetscReal p,d,dd; 5537045ce4SJed Brown p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 5637045ce4SJed Brown d = pd2 + (2*j-1)*pm1; 5737045ce4SJed Brown dd = pdd2 + (2*j-1)*pd1; 5837045ce4SJed Brown pm2 = pm1; 5937045ce4SJed Brown pm1 = p; 6037045ce4SJed Brown pd2 = pd1; 6137045ce4SJed Brown pd1 = d; 6237045ce4SJed Brown pdd2 = pdd1; 6337045ce4SJed Brown pdd1 = dd; 6437045ce4SJed Brown if (degrees[k] == j) { 6537045ce4SJed Brown if (B) B[i*ndegree+k] = p; 6637045ce4SJed Brown if (D) D[i*ndegree+k] = d; 6737045ce4SJed Brown if (D2) D2[i*ndegree+k] = dd; 6837045ce4SJed Brown } 6937045ce4SJed Brown } 7037045ce4SJed Brown } 7137045ce4SJed Brown PetscFunctionReturn(0); 7237045ce4SJed Brown } 7337045ce4SJed Brown 7437045ce4SJed Brown #undef __FUNCT__ 7537045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature" 7637045ce4SJed Brown /*@ 7737045ce4SJed Brown PetscDTGaussQuadrature - create Gauss quadrature 7837045ce4SJed Brown 7937045ce4SJed Brown Not Collective 8037045ce4SJed Brown 8137045ce4SJed Brown Input Arguments: 8237045ce4SJed Brown + npoints - number of points 8337045ce4SJed Brown . a - left end of interval (often-1) 8437045ce4SJed Brown - b - right end of interval (often +1) 8537045ce4SJed Brown 8637045ce4SJed Brown Output Arguments: 8737045ce4SJed Brown + x - quadrature points 8837045ce4SJed Brown - w - quadrature weights 8937045ce4SJed Brown 9037045ce4SJed Brown Level: intermediate 9137045ce4SJed Brown 9237045ce4SJed Brown References: 9337045ce4SJed Brown Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 9437045ce4SJed Brown 9537045ce4SJed Brown .seealso: PetscDTLegendreEval() 9637045ce4SJed Brown @*/ 9737045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 9837045ce4SJed Brown { 9937045ce4SJed Brown PetscErrorCode ierr; 10037045ce4SJed Brown PetscInt i; 10137045ce4SJed Brown PetscReal *work; 10237045ce4SJed Brown PetscScalar *Z; 10337045ce4SJed Brown PetscBLASInt N,LDZ,info; 10437045ce4SJed Brown 10537045ce4SJed Brown PetscFunctionBegin; 10637045ce4SJed Brown /* Set up the Golub-Welsch system */ 10737045ce4SJed Brown for (i=0; i<npoints; i++) { 10837045ce4SJed Brown x[i] = 0; /* diagonal is 0 */ 10937045ce4SJed Brown if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 11037045ce4SJed Brown } 11137045ce4SJed Brown ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 11237045ce4SJed Brown ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr); 113c5df96a5SBarry Smith ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 11437045ce4SJed Brown LDZ = N; 11537045ce4SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 116f75e95b9SBarry Smith PetscStackCall("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 11737045ce4SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 1181c3d6f74SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 11937045ce4SJed Brown 12037045ce4SJed Brown for (i=0; i<(npoints+1)/2; i++) { 12137045ce4SJed Brown PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 12237045ce4SJed Brown x[i] = (a+b)/2 - y*(b-a)/2; 12337045ce4SJed Brown x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 1240d644c17SKarl Rupp 12537045ce4SJed Brown w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints])); 12637045ce4SJed Brown } 12737045ce4SJed Brown ierr = PetscFree2(Z,work);CHKERRQ(ierr); 12837045ce4SJed Brown PetscFunctionReturn(0); 12937045ce4SJed Brown } 130*194825f6SJed Brown 131*194825f6SJed Brown #undef __FUNCT__ 132*194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR" 133*194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n 134*194825f6SJed Brown * A in column-major format 135*194825f6SJed Brown * Ainv in row-major format 136*194825f6SJed Brown * tau has length m 137*194825f6SJed Brown * worksize must be >= max(1,n) 138*194825f6SJed Brown */ 139*194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 140*194825f6SJed Brown { 141*194825f6SJed Brown PetscErrorCode ierr; 142*194825f6SJed Brown PetscBLASInt M,N,K,lda,ldb,ldwork,info; 143*194825f6SJed Brown PetscScalar *A,*Ainv,*R,*Q,Alpha; 144*194825f6SJed Brown 145*194825f6SJed Brown PetscFunctionBegin; 146*194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 147*194825f6SJed Brown { 148*194825f6SJed Brown PetscInt i,j; 149*194825f6SJed Brown ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr); 150*194825f6SJed Brown for (j=0; j<n; j++) { 151*194825f6SJed Brown for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 152*194825f6SJed Brown } 153*194825f6SJed Brown mstride = m; 154*194825f6SJed Brown } 155*194825f6SJed Brown #else 156*194825f6SJed Brown A = A_in; 157*194825f6SJed Brown Ainv = Ainv_out; 158*194825f6SJed Brown #endif 159*194825f6SJed Brown 160*194825f6SJed Brown ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 161*194825f6SJed Brown ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 162*194825f6SJed Brown ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 163*194825f6SJed Brown ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 164*194825f6SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 165*194825f6SJed Brown LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info); 166*194825f6SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 167*194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 168*194825f6SJed Brown R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 169*194825f6SJed Brown 170*194825f6SJed Brown /* Extract an explicit representation of Q */ 171*194825f6SJed Brown Q = Ainv; 172*194825f6SJed Brown ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 173*194825f6SJed Brown K = N; /* full rank */ 174*194825f6SJed Brown LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info); 175*194825f6SJed Brown if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 176*194825f6SJed Brown 177*194825f6SJed Brown /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 178*194825f6SJed Brown Alpha = 1.0; 179*194825f6SJed Brown ldb = lda; 180*194825f6SJed Brown BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 181*194825f6SJed Brown /* Ainv is Q, overwritten with inverse */ 182*194825f6SJed Brown 183*194825f6SJed Brown #if defined(PETSC_USE_COMPLEX) 184*194825f6SJed Brown { 185*194825f6SJed Brown PetscInt i; 186*194825f6SJed Brown for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 187*194825f6SJed Brown ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 188*194825f6SJed Brown } 189*194825f6SJed Brown #endif 190*194825f6SJed Brown PetscFunctionReturn(0); 191*194825f6SJed Brown } 192*194825f6SJed Brown 193*194825f6SJed Brown #undef __FUNCT__ 194*194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate" 195*194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 196*194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 197*194825f6SJed Brown { 198*194825f6SJed Brown PetscErrorCode ierr; 199*194825f6SJed Brown PetscReal *Bv; 200*194825f6SJed Brown PetscInt i,j; 201*194825f6SJed Brown 202*194825f6SJed Brown PetscFunctionBegin; 203*194825f6SJed Brown ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr); 204*194825f6SJed Brown /* Point evaluation of L_p on all the source vertices */ 205*194825f6SJed Brown ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 206*194825f6SJed Brown /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 207*194825f6SJed Brown for (i=0; i<ninterval; i++) { 208*194825f6SJed Brown for (j=0; j<ndegree; j++) { 209*194825f6SJed Brown if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 210*194825f6SJed Brown else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 211*194825f6SJed Brown } 212*194825f6SJed Brown } 213*194825f6SJed Brown ierr = PetscFree(Bv);CHKERRQ(ierr); 214*194825f6SJed Brown PetscFunctionReturn(0); 215*194825f6SJed Brown } 216*194825f6SJed Brown 217*194825f6SJed Brown #undef __FUNCT__ 218*194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly" 219*194825f6SJed Brown /*@ 220*194825f6SJed Brown PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 221*194825f6SJed Brown 222*194825f6SJed Brown Not Collective 223*194825f6SJed Brown 224*194825f6SJed Brown Input Arguments: 225*194825f6SJed Brown + degree - degree of reconstruction polynomial 226*194825f6SJed Brown . nsource - number of source intervals 227*194825f6SJed Brown . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 228*194825f6SJed Brown . ntarget - number of target intervals 229*194825f6SJed Brown - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 230*194825f6SJed Brown 231*194825f6SJed Brown Output Arguments: 232*194825f6SJed Brown . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 233*194825f6SJed Brown 234*194825f6SJed Brown Level: advanced 235*194825f6SJed Brown 236*194825f6SJed Brown .seealso: PetscDTLegendreEval() 237*194825f6SJed Brown @*/ 238*194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 239*194825f6SJed Brown { 240*194825f6SJed Brown PetscErrorCode ierr; 241*194825f6SJed Brown PetscInt i,j,k,*bdegrees,worksize; 242*194825f6SJed Brown PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 243*194825f6SJed Brown PetscScalar *tau,*work; 244*194825f6SJed Brown 245*194825f6SJed Brown PetscFunctionBegin; 246*194825f6SJed Brown PetscValidRealPointer(sourcex,3); 247*194825f6SJed Brown PetscValidRealPointer(targetx,5); 248*194825f6SJed Brown PetscValidRealPointer(R,6); 249*194825f6SJed Brown 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); 250*194825f6SJed Brown #if defined(PETSC_USE_DEBUG) 251*194825f6SJed Brown for (i=0; i<nsource; i++) { 252*194825f6SJed Brown if (sourcex[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]); 253*194825f6SJed Brown } 254*194825f6SJed Brown for (i=0; i<ntarget; i++) { 255*194825f6SJed Brown if (targetx[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]); 256*194825f6SJed Brown } 257*194825f6SJed Brown #endif 258*194825f6SJed Brown xmin = PetscMin(sourcex[0],targetx[0]); 259*194825f6SJed Brown xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 260*194825f6SJed Brown center = (xmin + xmax)/2; 261*194825f6SJed Brown hscale = (xmax - xmin)/2; 262*194825f6SJed Brown worksize = nsource; 263*194825f6SJed Brown ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr); 264*194825f6SJed Brown ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscScalar,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr); 265*194825f6SJed Brown for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 266*194825f6SJed Brown for (i=0; i<=degree; i++) bdegrees[i] = i+1; 267*194825f6SJed Brown ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 268*194825f6SJed Brown ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 269*194825f6SJed Brown for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 270*194825f6SJed Brown ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 271*194825f6SJed Brown for (i=0; i<ntarget; i++) { 272*194825f6SJed Brown PetscReal rowsum = 0; 273*194825f6SJed Brown for (j=0; j<nsource; j++) { 274*194825f6SJed Brown PetscReal sum = 0; 275*194825f6SJed Brown for (k=0; k<degree+1; k++) { 276*194825f6SJed Brown sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 277*194825f6SJed Brown } 278*194825f6SJed Brown R[i*nsource+j] = sum; 279*194825f6SJed Brown rowsum += sum; 280*194825f6SJed Brown } 281*194825f6SJed Brown for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 282*194825f6SJed Brown } 283*194825f6SJed Brown ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 284*194825f6SJed Brown ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 285*194825f6SJed Brown PetscFunctionReturn(0); 286*194825f6SJed Brown } 287