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