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