1 /* Discretization tools */ 2 3 #include <petscdt.h> /*I "petscdt.h" I*/ 4 #include <petscblaslapack.h> 5 #include <petsc-private/petscimpl.h> 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "PetscDTLegendreEval" 9 /*@ 10 PetscDTLegendreEval - evaluate Legendre polynomial at points 11 12 Not Collective 13 14 Input Arguments: 15 + npoints - number of spatial points to evaluate at 16 . points - array of locations to evaluate at 17 . ndegree - number of basis degrees to evaluate 18 - degrees - sorted array of degrees to evaluate 19 20 Output Arguments: 21 + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL) 22 . D - row-oriented derivative evaluation matrix (or NULL) 23 - D2 - row-oriented second derivative evaluation matrix (or NULL) 24 25 Level: intermediate 26 27 .seealso: PetscDTGaussQuadrature() 28 @*/ 29 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2) 30 { 31 PetscInt i,maxdegree; 32 33 PetscFunctionBegin; 34 if (!npoints || !ndegree) PetscFunctionReturn(0); 35 maxdegree = degrees[ndegree-1]; 36 for (i=0; i<npoints; i++) { 37 PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x; 38 PetscInt j,k; 39 x = points[i]; 40 pm2 = 0; 41 pm1 = 1; 42 pd2 = 0; 43 pd1 = 0; 44 pdd2 = 0; 45 pdd1 = 0; 46 k = 0; 47 if (degrees[k] == 0) { 48 if (B) B[i*ndegree+k] = pm1; 49 if (D) D[i*ndegree+k] = pd1; 50 if (D2) D2[i*ndegree+k] = pdd1; 51 k++; 52 } 53 for (j=1; j<=maxdegree; j++,k++) { 54 PetscReal p,d,dd; 55 p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j; 56 d = pd2 + (2*j-1)*pm1; 57 dd = pdd2 + (2*j-1)*pd1; 58 pm2 = pm1; 59 pm1 = p; 60 pd2 = pd1; 61 pd1 = d; 62 pdd2 = pdd1; 63 pdd1 = dd; 64 if (degrees[k] == j) { 65 if (B) B[i*ndegree+k] = p; 66 if (D) D[i*ndegree+k] = d; 67 if (D2) D2[i*ndegree+k] = dd; 68 } 69 } 70 } 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "PetscDTGaussQuadrature" 76 /*@ 77 PetscDTGaussQuadrature - create Gauss quadrature 78 79 Not Collective 80 81 Input Arguments: 82 + npoints - number of points 83 . a - left end of interval (often-1) 84 - b - right end of interval (often +1) 85 86 Output Arguments: 87 + x - quadrature points 88 - w - quadrature weights 89 90 Level: intermediate 91 92 References: 93 Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969. 94 95 .seealso: PetscDTLegendreEval() 96 @*/ 97 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w) 98 { 99 PetscErrorCode ierr; 100 PetscInt i; 101 PetscReal *work; 102 PetscScalar *Z; 103 PetscBLASInt N,LDZ,info; 104 105 PetscFunctionBegin; 106 /* Set up the Golub-Welsch system */ 107 for (i=0; i<npoints; i++) { 108 x[i] = 0; /* diagonal is 0 */ 109 if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i)); 110 } 111 ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 112 ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr); 113 ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr); 114 LDZ = N; 115 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 116 PetscStackCall("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info)); 117 ierr = PetscFPTrapPop();CHKERRQ(ierr); 118 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error"); 119 120 for (i=0; i<(npoints+1)/2; i++) { 121 PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */ 122 x[i] = (a+b)/2 - y*(b-a)/2; 123 x[npoints-i-1] = (a+b)/2 + y*(b-a)/2; 124 125 w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints])); 126 } 127 ierr = PetscFree2(Z,work);CHKERRQ(ierr); 128 PetscFunctionReturn(0); 129 } 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "PetscDTPseudoInverseQR" 133 /* Overwrites A. Can only handle full-rank problems with m>=n 134 * A in column-major format 135 * Ainv in row-major format 136 * tau has length m 137 * worksize must be >= max(1,n) 138 */ 139 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 140 { 141 PetscErrorCode ierr; 142 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 143 PetscScalar *A,*Ainv,*R,*Q,Alpha; 144 145 PetscFunctionBegin; 146 #if defined(PETSC_USE_COMPLEX) 147 { 148 PetscInt i,j; 149 ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr); 150 for (j=0; j<n; j++) { 151 for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j]; 152 } 153 mstride = m; 154 } 155 #else 156 A = A_in; 157 Ainv = Ainv_out; 158 #endif 159 160 ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr); 161 ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr); 162 ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr); 163 ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr); 164 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 165 LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info); 166 ierr = PetscFPTrapPop();CHKERRQ(ierr); 167 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 168 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 169 170 /* Extract an explicit representation of Q */ 171 Q = Ainv; 172 ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr); 173 K = N; /* full rank */ 174 LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info); 175 if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 176 177 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 178 Alpha = 1.0; 179 ldb = lda; 180 BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 181 /* Ainv is Q, overwritten with inverse */ 182 183 #if defined(PETSC_USE_COMPLEX) 184 { 185 PetscInt i; 186 for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]); 187 ierr = PetscFree2(A,Ainv);CHKERRQ(ierr); 188 } 189 #endif 190 PetscFunctionReturn(0); 191 } 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "PetscDTLegendreIntegrate" 195 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */ 196 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B) 197 { 198 PetscErrorCode ierr; 199 PetscReal *Bv; 200 PetscInt i,j; 201 202 PetscFunctionBegin; 203 ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr); 204 /* Point evaluation of L_p on all the source vertices */ 205 ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr); 206 /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */ 207 for (i=0; i<ninterval; i++) { 208 for (j=0; j<ndegree; j++) { 209 if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 210 else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j]; 211 } 212 } 213 ierr = PetscFree(Bv);CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "PetscDTReconstructPoly" 219 /*@ 220 PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals 221 222 Not Collective 223 224 Input Arguments: 225 + degree - degree of reconstruction polynomial 226 . nsource - number of source intervals 227 . sourcex - sorted coordinates of source cell boundaries (length nsource+1) 228 . ntarget - number of target intervals 229 - targetx - sorted coordinates of target cell boundaries (length ntarget+1) 230 231 Output Arguments: 232 . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s] 233 234 Level: advanced 235 236 .seealso: PetscDTLegendreEval() 237 @*/ 238 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R) 239 { 240 PetscErrorCode ierr; 241 PetscInt i,j,k,*bdegrees,worksize; 242 PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget; 243 PetscScalar *tau,*work; 244 245 PetscFunctionBegin; 246 PetscValidRealPointer(sourcex,3); 247 PetscValidRealPointer(targetx,5); 248 PetscValidRealPointer(R,6); 249 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 #if defined(PETSC_USE_DEBUG) 251 for (i=0; i<nsource; i++) { 252 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 } 254 for (i=0; i<ntarget; i++) { 255 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 } 257 #endif 258 xmin = PetscMin(sourcex[0],targetx[0]); 259 xmax = PetscMax(sourcex[nsource],targetx[ntarget]); 260 center = (xmin + xmax)/2; 261 hscale = (xmax - xmin)/2; 262 worksize = nsource; 263 ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr); 264 ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscScalar,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr); 265 for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale; 266 for (i=0; i<=degree; i++) bdegrees[i] = i+1; 267 ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr); 268 ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr); 269 for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale; 270 ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr); 271 for (i=0; i<ntarget; i++) { 272 PetscReal rowsum = 0; 273 for (j=0; j<nsource; j++) { 274 PetscReal sum = 0; 275 for (k=0; k<degree+1; k++) { 276 sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j]; 277 } 278 R[i*nsource+j] = sum; 279 rowsum += sum; 280 } 281 for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */ 282 } 283 ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr); 284 ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr); 285 PetscFunctionReturn(0); 286 } 287