xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision 5e71baeff2f3138f93cd4f5927dfd596eb8325cc)
1 
2 #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
3 #include <petscblaslapack.h>
4 
5 /*
6    Private context (data structure) for the SVD preconditioner.
7 */
8 typedef struct {
9   Vec         diag,work;
10   Mat         A,U,Vt;
11   PetscInt    nzero;
12   PetscReal   zerosing;         /* measure of smallest singular value treated as nonzero */
13   PetscInt    essrank;          /* essential rank of operator */
14   VecScatter  left2red,right2red;
15   Vec         leftred,rightred;
16   PetscViewer monitor;
17 } PC_SVD;
18 
19 typedef enum {READ=1, WRITE=2, READ_WRITE=3} AccessMode;
20 
21 /* -------------------------------------------------------------------------- */
22 /*
23    PCSetUp_SVD - Prepares for the use of the SVD preconditioner
24                     by setting data structures and options.
25 
26    Input Parameter:
27 .  pc - the preconditioner context
28 
29    Application Interface Routine: PCSetUp()
30 
31    Notes:
32    The interface routine PCSetUp() is not usually called directly by
33    the user, but instead is called by PCApply() if necessary.
34 */
35 static PetscErrorCode PCSetUp_SVD(PC pc)
36 {
37 #if defined(PETSC_MISSING_LAPACK_GESVD)
38   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates.");
39 #else
40   PC_SVD         *jac = (PC_SVD*)pc->data;
41   PetscErrorCode ierr;
42   PetscScalar    *a,*u,*v,*d,*work;
43   PetscBLASInt   nb,lwork;
44   PetscInt       i,n;
45   PetscMPIInt    size;
46 
47   PetscFunctionBegin;
48   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
49   ierr = MPI_Comm_size(((PetscObject)pc->pmat)->comm,&size);CHKERRQ(ierr);
50   if (size > 1) {
51     Mat          redmat;
52     ierr = MatCreateRedundantMatrix(pc->pmat,size,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&redmat);CHKERRQ(ierr);
53     ierr = MatConvert(redmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
54     ierr = MatDestroy(&redmat);CHKERRQ(ierr);
55   } else {
56     ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
57   }
58   if (!jac->diag) {    /* assume square matrices */
59     ierr = MatCreateVecs(jac->A,&jac->diag,&jac->work);CHKERRQ(ierr);
60   }
61   if (!jac->U) {
62     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr);
63     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->Vt);CHKERRQ(ierr);
64   }
65   ierr  = MatGetSize(jac->A,&n,NULL);CHKERRQ(ierr);
66   if (!n) {
67     ierr = PetscInfo(pc,"Matrix has zero rows, skipping svd\n");CHKERRQ(ierr);
68     PetscFunctionReturn(0);
69   }
70   ierr  = PetscBLASIntCast(n,&nb);CHKERRQ(ierr);
71   lwork = 5*nb;
72   ierr  = PetscMalloc1(lwork,&work);CHKERRQ(ierr);
73   ierr  = MatDenseGetArray(jac->A,&a);CHKERRQ(ierr);
74   ierr  = MatDenseGetArray(jac->U,&u);CHKERRQ(ierr);
75   ierr  = MatDenseGetArray(jac->Vt,&v);CHKERRQ(ierr);
76   ierr  = VecGetArray(jac->diag,&d);CHKERRQ(ierr);
77 #if !defined(PETSC_USE_COMPLEX)
78   {
79     PetscBLASInt lierr;
80     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
81     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr));
82     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
83     ierr = PetscFPTrapPop();CHKERRQ(ierr);
84   }
85 #else
86   {
87     PetscBLASInt lierr;
88     PetscReal    *rwork,*dd;
89     ierr = PetscMalloc1(5*nb,&rwork);CHKERRQ(ierr);
90     ierr = PetscMalloc1(nb,&dd);CHKERRQ(ierr);
91     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
92     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,dd,u,&nb,v,&nb,work,&lwork,rwork,&lierr));
93     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
94     ierr = PetscFree(rwork);CHKERRQ(ierr);
95     for (i=0; i<n; i++) d[i] = dd[i];
96     ierr = PetscFree(dd);CHKERRQ(ierr);
97     ierr = PetscFPTrapPop();CHKERRQ(ierr);
98   }
99 #endif
100   ierr = MatDenseRestoreArray(jac->A,&a);CHKERRQ(ierr);
101   ierr = MatDenseRestoreArray(jac->U,&u);CHKERRQ(ierr);
102   ierr = MatDenseRestoreArray(jac->Vt,&v);CHKERRQ(ierr);
103   for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break;
104   jac->nzero = n-1-i;
105   if (jac->monitor) {
106     ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
107     ierr = PetscViewerASCIIPrintf(jac->monitor,"    SVD: condition number %14.12e, %D of %D singular values are (nearly) zero\n",(double)PetscRealPart(d[0]/d[n-1]),jac->nzero,n);CHKERRQ(ierr);
108     if (n >= 10) {              /* print 5 smallest and 5 largest */
109       ierr = PetscViewerASCIIPrintf(jac->monitor,"    SVD: smallest singular values: %14.12e %14.12e %14.12e %14.12e %14.12e\n",(double)PetscRealPart(d[n-1]),(double)PetscRealPart(d[n-2]),(double)PetscRealPart(d[n-3]),(double)PetscRealPart(d[n-4]),(double)PetscRealPart(d[n-5]));CHKERRQ(ierr);
110       ierr = PetscViewerASCIIPrintf(jac->monitor,"    SVD: largest singular values : %14.12e %14.12e %14.12e %14.12e %14.12e\n",(double)PetscRealPart(d[4]),(double)PetscRealPart(d[3]),(double)PetscRealPart(d[2]),(double)PetscRealPart(d[1]),(double)PetscRealPart(d[0]));CHKERRQ(ierr);
111     } else {                    /* print all singular values */
112       char     buf[256],*p;
113       size_t   left = sizeof(buf),used;
114       PetscInt thisline;
115       for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) {
116         ierr  = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr);
117         left -= used;
118         p    += used;
119         if (thisline > 4 || i==0) {
120           ierr     = PetscViewerASCIIPrintf(jac->monitor,"    SVD: singular values:%s\n",buf);CHKERRQ(ierr);
121           p        = buf;
122           thisline = 0;
123         }
124       }
125     }
126     ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
127   }
128   ierr = PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1]));CHKERRQ(ierr);
129   for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i];
130   for (; i<n; i++) d[i] = 0.0;
131   if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */
132   ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);CHKERRQ(ierr);
133   ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr);
134 #if defined(foo)
135   {
136     PetscViewer viewer;
137     ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
138     ierr = MatView(jac->A,viewer);CHKERRQ(ierr);
139     ierr = MatView(jac->U,viewer);CHKERRQ(ierr);
140     ierr = MatView(jac->Vt,viewer);CHKERRQ(ierr);
141     ierr = VecView(jac->diag,viewer);CHKERRQ(ierr);
142     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
143   }
144 #endif
145   ierr = PetscFree(work);CHKERRQ(ierr);
146   PetscFunctionReturn(0);
147 #endif
148 }
149 
150 static PetscErrorCode PCSVDGetVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred)
151 {
152   PC_SVD         *jac = (PC_SVD*)pc->data;
153   PetscErrorCode ierr;
154   PetscMPIInt    size;
155 
156   PetscFunctionBegin;
157   ierr  = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
158   *xred = NULL;
159   switch (side) {
160   case PC_LEFT:
161     if (size == 1) *xred = x;
162     else {
163       if (!jac->left2red) {ierr = VecScatterCreateToAll(x,&jac->left2red,&jac->leftred);CHKERRQ(ierr);}
164       if (amode & READ) {
165         ierr = VecScatterBegin(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
166         ierr = VecScatterEnd(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
167       }
168       *xred = jac->leftred;
169     }
170     break;
171   case PC_RIGHT:
172     if (size == 1) *xred = x;
173     else {
174       if (!jac->right2red) {ierr = VecScatterCreateToAll(x,&jac->right2red,&jac->rightred);CHKERRQ(ierr);}
175       if (amode & READ) {
176         ierr = VecScatterBegin(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
177         ierr = VecScatterEnd(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178       }
179       *xred = jac->rightred;
180     }
181     break;
182   default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT");
183   }
184   PetscFunctionReturn(0);
185 }
186 
187 static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred)
188 {
189   PC_SVD         *jac = (PC_SVD*)pc->data;
190   PetscErrorCode ierr;
191   PetscMPIInt    size;
192 
193   PetscFunctionBegin;
194   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
195   switch (side) {
196   case PC_LEFT:
197     if (size != 1 && amode & WRITE) {
198       ierr = VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
199       ierr = VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
200     }
201     break;
202   case PC_RIGHT:
203     if (size != 1 && amode & WRITE) {
204       ierr = VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
205       ierr = VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
206     }
207     break;
208   default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT");
209   }
210   *xred = NULL;
211   PetscFunctionReturn(0);
212 }
213 
214 /* -------------------------------------------------------------------------- */
215 /*
216    PCApply_SVD - Applies the SVD preconditioner to a vector.
217 
218    Input Parameters:
219 .  pc - the preconditioner context
220 .  x - input vector
221 
222    Output Parameter:
223 .  y - output vector
224 
225    Application Interface Routine: PCApply()
226  */
227 static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y)
228 {
229   PC_SVD         *jac = (PC_SVD*)pc->data;
230   Vec            work = jac->work,xred,yred;
231   PetscErrorCode ierr;
232 
233   PetscFunctionBegin;
234   ierr = PCSVDGetVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr);
235   ierr = PCSVDGetVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr);
236 #if !defined(PETSC_USE_COMPLEX)
237   ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr);
238 #else
239   ierr = MatMultHermitianTranspose(jac->U,xred,work);CHKERRQ(ierr);
240 #endif
241   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
242 #if !defined(PETSC_USE_COMPLEX)
243   ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
244 #else
245   ierr = MatMultHermitianTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
246 #endif
247   ierr = PCSVDRestoreVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr);
248   ierr = PCSVDRestoreVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 
252 static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y)
253 {
254   PC_SVD         *jac = (PC_SVD*)pc->data;
255   Vec            work = jac->work,xred,yred;
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   ierr = PCSVDGetVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr);
260   ierr = PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr);
261 #if !defined(PETSC_USE_COMPLEX)
262   ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
263 #else
264   ierr = MatMultHermitianTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
265 #endif
266   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
267 #if !defined(PETSC_USE_COMPLEX)
268   ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr);
269 #else
270   ierr = MatMultHermitianTranspose(jac->U,xred,work);CHKERRQ(ierr);
271 #endif
272   ierr = PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr);
273   ierr = PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr);
274   PetscFunctionReturn(0);
275 }
276 
277 static PetscErrorCode PCReset_SVD(PC pc)
278 {
279   PC_SVD         *jac = (PC_SVD*)pc->data;
280   PetscErrorCode ierr;
281 
282   PetscFunctionBegin;
283   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
284   ierr = MatDestroy(&jac->U);CHKERRQ(ierr);
285   ierr = MatDestroy(&jac->Vt);CHKERRQ(ierr);
286   ierr = VecDestroy(&jac->diag);CHKERRQ(ierr);
287   ierr = VecDestroy(&jac->work);CHKERRQ(ierr);
288   ierr = VecScatterDestroy(&jac->right2red);CHKERRQ(ierr);
289   ierr = VecScatterDestroy(&jac->left2red);CHKERRQ(ierr);
290   ierr = VecDestroy(&jac->rightred);CHKERRQ(ierr);
291   ierr = VecDestroy(&jac->leftred);CHKERRQ(ierr);
292   PetscFunctionReturn(0);
293 }
294 
295 /* -------------------------------------------------------------------------- */
296 /*
297    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
298    that was created with PCCreate_SVD().
299 
300    Input Parameter:
301 .  pc - the preconditioner context
302 
303    Application Interface Routine: PCDestroy()
304 */
305 static PetscErrorCode PCDestroy_SVD(PC pc)
306 {
307   PC_SVD         *jac = (PC_SVD*)pc->data;
308   PetscErrorCode ierr;
309 
310   PetscFunctionBegin;
311   ierr = PCReset_SVD(pc);CHKERRQ(ierr);
312   ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
313   ierr = PetscFree(pc->data);CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 
317 static PetscErrorCode PCSetFromOptions_SVD(PetscOptionItems *PetscOptionsObject,PC pc)
318 {
319   PetscErrorCode ierr;
320   PC_SVD         *jac = (PC_SVD*)pc->data;
321   PetscBool      flg,set;
322 
323   PetscFunctionBegin;
324   ierr = PetscOptionsHead(PetscOptionsObject,"SVD options");CHKERRQ(ierr);
325   ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,NULL);CHKERRQ(ierr);
326   ierr = PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,NULL);CHKERRQ(ierr);
327   ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
328   if (set) {                    /* Should make PCSVDSetMonitor() */
329     if (flg && !jac->monitor) {
330       ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),"stdout",&jac->monitor);CHKERRQ(ierr);
331     } else if (!flg) {
332       ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
333     }
334   }
335   ierr = PetscOptionsTail();CHKERRQ(ierr);
336   PetscFunctionReturn(0);
337 }
338 
339 static PetscErrorCode PCView_SVD(PC pc,PetscViewer viewer)
340 {
341   PC_SVD         *svd = (PC_SVD*)pc->data;
342   PetscErrorCode ierr;
343   PetscBool      iascii;
344 
345   PetscFunctionBegin;
346   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
347   if (iascii) {
348     ierr = PetscViewerASCIIPrintf(viewer,"  All singular values smaller than %g treated as zero\n",(double)svd->zerosing);CHKERRQ(ierr);
349     ierr = PetscViewerASCIIPrintf(viewer,"  Provided essential rank of the matrix %D (all other eigenvalues are zeroed)\n",svd->essrank);CHKERRQ(ierr);
350   }
351   PetscFunctionReturn(0);
352 }
353 /* -------------------------------------------------------------------------- */
354 /*
355    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
356    and sets this as the private data within the generic preconditioning
357    context, PC, that was created within PCCreate().
358 
359    Input Parameter:
360 .  pc - the preconditioner context
361 
362    Application Interface Routine: PCCreate()
363 */
364 
365 /*MC
366      PCSVD - Use pseudo inverse defined by SVD of operator
367 
368    Level: advanced
369 
370   Concepts: SVD
371 
372   Options Database:
373 -  -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero
374 +  -pc_svd_monitor  Print information on the extreme singular values of the operator
375 
376   Developer Note: This implementation automatically creates a redundant copy of the
377    matrix on each process and uses a sequential SVD solve. Why does it do this instead
378    of using the composable PCREDUNDANT object?
379 
380 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
381 M*/
382 
383 PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc)
384 {
385   PC_SVD         *jac;
386   PetscErrorCode ierr;
387 
388   PetscFunctionBegin;
389   /*
390      Creates the private data structure for this preconditioner and
391      attach it to the PC object.
392   */
393   ierr          = PetscNewLog(pc,&jac);CHKERRQ(ierr);
394   jac->zerosing = 1.e-12;
395   pc->data      = (void*)jac;
396 
397   /*
398       Set the pointers for the functions that are provided above.
399       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
400       are called, they will automatically call these functions.  Note we
401       choose not to provide a couple of these functions since they are
402       not needed.
403   */
404   pc->ops->apply           = PCApply_SVD;
405   pc->ops->applytranspose  = PCApplyTranspose_SVD;
406   pc->ops->setup           = PCSetUp_SVD;
407   pc->ops->reset           = PCReset_SVD;
408   pc->ops->destroy         = PCDestroy_SVD;
409   pc->ops->setfromoptions  = PCSetFromOptions_SVD;
410   pc->ops->view            = PCView_SVD;
411   pc->ops->applyrichardson = 0;
412   PetscFunctionReturn(0);
413 }
414 
415