xref: /petsc/src/dm/dt/interface/dt.c (revision 8b83055f287c8deb1a16325ae3faee24d0648b7a)
137045ce4SJed Brown /* Discretization tools */
237045ce4SJed Brown 
337045ce4SJed Brown #include <petscdt.h>            /*I "petscdt.h" I*/
437045ce4SJed Brown #include <petscblaslapack.h>
5194825f6SJed Brown #include <petsc-private/petscimpl.h>
6665c2dedSJed Brown #include <petscviewer.h>
737045ce4SJed Brown 
837045ce4SJed Brown #undef __FUNCT__
937045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval"
1037045ce4SJed Brown /*@
1137045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
1237045ce4SJed Brown 
1337045ce4SJed Brown    Not Collective
1437045ce4SJed Brown 
1537045ce4SJed Brown    Input Arguments:
1637045ce4SJed Brown +  npoints - number of spatial points to evaluate at
1737045ce4SJed Brown .  points - array of locations to evaluate at
1837045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
1937045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
2037045ce4SJed Brown 
2137045ce4SJed Brown    Output Arguments:
220298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
230298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
240298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
2537045ce4SJed Brown 
2637045ce4SJed Brown    Level: intermediate
2737045ce4SJed Brown 
2837045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
2937045ce4SJed Brown @*/
3037045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
3137045ce4SJed Brown {
3237045ce4SJed Brown   PetscInt i,maxdegree;
3337045ce4SJed Brown 
3437045ce4SJed Brown   PetscFunctionBegin;
3537045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
3637045ce4SJed Brown   maxdegree = degrees[ndegree-1];
3737045ce4SJed Brown   for (i=0; i<npoints; i++) {
3837045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
3937045ce4SJed Brown     PetscInt  j,k;
4037045ce4SJed Brown     x    = points[i];
4137045ce4SJed Brown     pm2  = 0;
4237045ce4SJed Brown     pm1  = 1;
4337045ce4SJed Brown     pd2  = 0;
4437045ce4SJed Brown     pd1  = 0;
4537045ce4SJed Brown     pdd2 = 0;
4637045ce4SJed Brown     pdd1 = 0;
4737045ce4SJed Brown     k    = 0;
4837045ce4SJed Brown     if (degrees[k] == 0) {
4937045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
5037045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
5137045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
5237045ce4SJed Brown       k++;
5337045ce4SJed Brown     }
5437045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
5537045ce4SJed Brown       PetscReal p,d,dd;
5637045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
5737045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
5837045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
5937045ce4SJed Brown       pm2  = pm1;
6037045ce4SJed Brown       pm1  = p;
6137045ce4SJed Brown       pd2  = pd1;
6237045ce4SJed Brown       pd1  = d;
6337045ce4SJed Brown       pdd2 = pdd1;
6437045ce4SJed Brown       pdd1 = dd;
6537045ce4SJed Brown       if (degrees[k] == j) {
6637045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
6737045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
6837045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
6937045ce4SJed Brown       }
7037045ce4SJed Brown     }
7137045ce4SJed Brown   }
7237045ce4SJed Brown   PetscFunctionReturn(0);
7337045ce4SJed Brown }
7437045ce4SJed Brown 
7537045ce4SJed Brown #undef __FUNCT__
7637045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature"
7737045ce4SJed Brown /*@
7837045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
7937045ce4SJed Brown 
8037045ce4SJed Brown    Not Collective
8137045ce4SJed Brown 
8237045ce4SJed Brown    Input Arguments:
8337045ce4SJed Brown +  npoints - number of points
8437045ce4SJed Brown .  a - left end of interval (often-1)
8537045ce4SJed Brown -  b - right end of interval (often +1)
8637045ce4SJed Brown 
8737045ce4SJed Brown    Output Arguments:
8837045ce4SJed Brown +  x - quadrature points
8937045ce4SJed Brown -  w - quadrature weights
9037045ce4SJed Brown 
9137045ce4SJed Brown    Level: intermediate
9237045ce4SJed Brown 
9337045ce4SJed Brown    References:
9437045ce4SJed Brown    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
9537045ce4SJed Brown 
9637045ce4SJed Brown .seealso: PetscDTLegendreEval()
9737045ce4SJed Brown @*/
9837045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
9937045ce4SJed Brown {
10037045ce4SJed Brown   PetscErrorCode ierr;
10137045ce4SJed Brown   PetscInt       i;
10237045ce4SJed Brown   PetscReal      *work;
10337045ce4SJed Brown   PetscScalar    *Z;
10437045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
10537045ce4SJed Brown 
10637045ce4SJed Brown   PetscFunctionBegin;
10737045ce4SJed Brown   /* Set up the Golub-Welsch system */
10837045ce4SJed Brown   for (i=0; i<npoints; i++) {
10937045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
11037045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
11137045ce4SJed Brown   }
11237045ce4SJed Brown   ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
11337045ce4SJed Brown   ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr);
114c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
11537045ce4SJed Brown   LDZ  = N;
11637045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
117*8b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
11837045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1191c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
12037045ce4SJed Brown 
12137045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
12237045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
12337045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
12437045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
1250d644c17SKarl Rupp 
12637045ce4SJed Brown     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
12737045ce4SJed Brown   }
12837045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
12937045ce4SJed Brown   PetscFunctionReturn(0);
13037045ce4SJed Brown }
131194825f6SJed Brown 
132194825f6SJed Brown #undef __FUNCT__
133194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR"
134194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
135194825f6SJed Brown  * A in column-major format
136194825f6SJed Brown  * Ainv in row-major format
137194825f6SJed Brown  * tau has length m
138194825f6SJed Brown  * worksize must be >= max(1,n)
139194825f6SJed Brown  */
140194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
141194825f6SJed Brown {
142194825f6SJed Brown   PetscErrorCode ierr;
143194825f6SJed Brown   PetscBLASInt M,N,K,lda,ldb,ldwork,info;
144194825f6SJed Brown   PetscScalar *A,*Ainv,*R,*Q,Alpha;
145194825f6SJed Brown 
146194825f6SJed Brown   PetscFunctionBegin;
147194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
148194825f6SJed Brown   {
149194825f6SJed Brown     PetscInt i,j;
150194825f6SJed Brown     ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr);
151194825f6SJed Brown     for (j=0; j<n; j++) {
152194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
153194825f6SJed Brown     }
154194825f6SJed Brown     mstride = m;
155194825f6SJed Brown   }
156194825f6SJed Brown #else
157194825f6SJed Brown   A = A_in;
158194825f6SJed Brown   Ainv = Ainv_out;
159194825f6SJed Brown #endif
160194825f6SJed Brown 
161194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
162194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
163194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
164194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
165194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
166194825f6SJed Brown   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
167194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
168194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
169194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
170194825f6SJed Brown 
171194825f6SJed Brown   /* Extract an explicit representation of Q */
172194825f6SJed Brown   Q = Ainv;
173194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
174194825f6SJed Brown   K = N;                        /* full rank */
175194825f6SJed Brown   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
176194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
177194825f6SJed Brown 
178194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
179194825f6SJed Brown   Alpha = 1.0;
180194825f6SJed Brown   ldb = lda;
181194825f6SJed Brown   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
182194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
183194825f6SJed Brown 
184194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
185194825f6SJed Brown   {
186194825f6SJed Brown     PetscInt i;
187194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
188194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
189194825f6SJed Brown   }
190194825f6SJed Brown #endif
191194825f6SJed Brown   PetscFunctionReturn(0);
192194825f6SJed Brown }
193194825f6SJed Brown 
194194825f6SJed Brown #undef __FUNCT__
195194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate"
196194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
197194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
198194825f6SJed Brown {
199194825f6SJed Brown   PetscErrorCode ierr;
200194825f6SJed Brown   PetscReal *Bv;
201194825f6SJed Brown   PetscInt i,j;
202194825f6SJed Brown 
203194825f6SJed Brown   PetscFunctionBegin;
204194825f6SJed Brown   ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr);
205194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
206194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
207194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
208194825f6SJed Brown   for (i=0; i<ninterval; i++) {
209194825f6SJed Brown     for (j=0; j<ndegree; j++) {
210194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
211194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
212194825f6SJed Brown     }
213194825f6SJed Brown   }
214194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
215194825f6SJed Brown   PetscFunctionReturn(0);
216194825f6SJed Brown }
217194825f6SJed Brown 
218194825f6SJed Brown #undef __FUNCT__
219194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly"
220194825f6SJed Brown /*@
221194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
222194825f6SJed Brown 
223194825f6SJed Brown    Not Collective
224194825f6SJed Brown 
225194825f6SJed Brown    Input Arguments:
226194825f6SJed Brown +  degree - degree of reconstruction polynomial
227194825f6SJed Brown .  nsource - number of source intervals
228194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
229194825f6SJed Brown .  ntarget - number of target intervals
230194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
231194825f6SJed Brown 
232194825f6SJed Brown    Output Arguments:
233194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
234194825f6SJed Brown 
235194825f6SJed Brown    Level: advanced
236194825f6SJed Brown 
237194825f6SJed Brown .seealso: PetscDTLegendreEval()
238194825f6SJed Brown @*/
239194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
240194825f6SJed Brown {
241194825f6SJed Brown   PetscErrorCode ierr;
242194825f6SJed Brown   PetscInt i,j,k,*bdegrees,worksize;
243194825f6SJed Brown   PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
244194825f6SJed Brown   PetscScalar *tau,*work;
245194825f6SJed Brown 
246194825f6SJed Brown   PetscFunctionBegin;
247194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
248194825f6SJed Brown   PetscValidRealPointer(targetx,5);
249194825f6SJed Brown   PetscValidRealPointer(R,6);
250194825f6SJed 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);
251194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
252194825f6SJed Brown   for (i=0; i<nsource; i++) {
253194825f6SJed 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]);
254194825f6SJed Brown   }
255194825f6SJed Brown   for (i=0; i<ntarget; i++) {
256194825f6SJed 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]);
257194825f6SJed Brown   }
258194825f6SJed Brown #endif
259194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
260194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
261194825f6SJed Brown   center = (xmin + xmax)/2;
262194825f6SJed Brown   hscale = (xmax - xmin)/2;
263194825f6SJed Brown   worksize = nsource;
264194825f6SJed Brown   ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr);
26582772646SJed Brown   ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr);
266194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
267194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
268194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
269194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
270194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
271194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
272194825f6SJed Brown   for (i=0; i<ntarget; i++) {
273194825f6SJed Brown     PetscReal rowsum = 0;
274194825f6SJed Brown     for (j=0; j<nsource; j++) {
275194825f6SJed Brown       PetscReal sum = 0;
276194825f6SJed Brown       for (k=0; k<degree+1; k++) {
277194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
278194825f6SJed Brown       }
279194825f6SJed Brown       R[i*nsource+j] = sum;
280194825f6SJed Brown       rowsum += sum;
281194825f6SJed Brown     }
282194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
283194825f6SJed Brown   }
284194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
285194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
286194825f6SJed Brown   PetscFunctionReturn(0);
287194825f6SJed Brown }
288