xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
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 #undef __FUNCT__
36 #define __FUNCT__ "PCSetUp_SVD"
37 static PetscErrorCode PCSetUp_SVD(PC pc)
38 {
39 #if defined(PETSC_MISSING_LAPACK_GESVD)
40   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates.");
41 #else
42   PC_SVD         *jac = (PC_SVD*)pc->data;
43   PetscErrorCode ierr;
44   PetscScalar    *a,*u,*v,*d,*work;
45   PetscBLASInt   nb,lwork;
46   PetscInt       i,n;
47   PetscMPIInt    size;
48 
49   PetscFunctionBegin;
50   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
51   ierr = MPI_Comm_size(((PetscObject)pc->pmat)->comm,&size);CHKERRQ(ierr);
52   if (size > 1) {
53     Mat          redmat;
54     ierr = MatCreateRedundantMatrix(pc->pmat,size,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&redmat);CHKERRQ(ierr);
55     ierr = MatConvert(redmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
56     ierr = MatDestroy(&redmat);CHKERRQ(ierr);
57   } else {
58     ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
59   }
60   if (!jac->diag) {    /* assume square matrices */
61     ierr = MatCreateVecs(jac->A,&jac->diag,&jac->work);CHKERRQ(ierr);
62   }
63   if (!jac->U) {
64     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr);
65     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->Vt);CHKERRQ(ierr);
66   }
67   ierr  = MatGetSize(pc->pmat,&n,NULL);CHKERRQ(ierr);
68   ierr  = PetscBLASIntCast(n,&nb);CHKERRQ(ierr);
69   lwork = 5*nb;
70   ierr  = PetscMalloc1(lwork,&work);CHKERRQ(ierr);
71   ierr  = MatDenseGetArray(jac->A,&a);CHKERRQ(ierr);
72   ierr  = MatDenseGetArray(jac->U,&u);CHKERRQ(ierr);
73   ierr  = MatDenseGetArray(jac->Vt,&v);CHKERRQ(ierr);
74   ierr  = VecGetArray(jac->diag,&d);CHKERRQ(ierr);
75 #if !defined(PETSC_USE_COMPLEX)
76   {
77     PetscBLASInt lierr;
78     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
79     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr));
80     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
81     ierr = PetscFPTrapPop();CHKERRQ(ierr);
82   }
83 #else
84   {
85     PetscBLASInt lierr;
86     PetscReal    *rwork,*dd;
87     ierr = PetscMalloc1(5*nb,&rwork);CHKERRQ(ierr);
88     ierr = PetscMalloc1(nb,&dd);CHKERRQ(ierr);
89     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
90     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,dd,u,&nb,v,&nb,work,&lwork,rwork,&lierr));
91     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
92     ierr = PetscFree(rwork);CHKERRQ(ierr);
93     for (i=0; i<nb; i++) d[i] = dd[i];
94     ierr = PetscFree(dd);CHKERRQ(ierr);
95     ierr = PetscFPTrapPop();CHKERRQ(ierr);
96   }
97 #endif
98   ierr = MatDenseRestoreArray(jac->A,&a);CHKERRQ(ierr);
99   ierr = MatDenseRestoreArray(jac->U,&u);CHKERRQ(ierr);
100   ierr = MatDenseRestoreArray(jac->Vt,&v);CHKERRQ(ierr);
101   for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break;
102   jac->nzero = n-1-i;
103   if (jac->monitor) {
104     ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
105     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);
106     if (n >= 10) {              /* print 5 smallest and 5 largest */
107       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);
108       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);
109     } else {                    /* print all singular values */
110       char     buf[256],*p;
111       size_t   left = sizeof(buf),used;
112       PetscInt thisline;
113       for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) {
114         ierr  = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr);
115         left -= used;
116         p    += used;
117         if (thisline > 4 || i==0) {
118           ierr     = PetscViewerASCIIPrintf(jac->monitor,"    SVD: singular values:%s\n",buf);CHKERRQ(ierr);
119           p        = buf;
120           thisline = 0;
121         }
122       }
123     }
124     ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
125   }
126   ierr = PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1]));CHKERRQ(ierr);
127   for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i];
128   for (; i<n; i++) d[i] = 0.0;
129   if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */
130   ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);CHKERRQ(ierr);
131   ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr);
132 #if defined(foo)
133   {
134     PetscViewer viewer;
135     ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
136     ierr = MatView(jac->A,viewer);CHKERRQ(ierr);
137     ierr = MatView(jac->U,viewer);CHKERRQ(ierr);
138     ierr = MatView(jac->Vt,viewer);CHKERRQ(ierr);
139     ierr = VecView(jac->diag,viewer);CHKERRQ(ierr);
140     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
141   }
142 #endif
143   ierr = PetscFree(work);CHKERRQ(ierr);
144   PetscFunctionReturn(0);
145 #endif
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "PCSVDGetVec"
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 #undef __FUNCT__
188 #define __FUNCT__ "PCSVDRestoreVec"
189 static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred)
190 {
191   PC_SVD         *jac = (PC_SVD*)pc->data;
192   PetscErrorCode ierr;
193   PetscMPIInt    size;
194 
195   PetscFunctionBegin;
196   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
197   switch (side) {
198   case PC_LEFT:
199     if (size != 1 && amode & WRITE) {
200       ierr = VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
201       ierr = VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
202     }
203     break;
204   case PC_RIGHT:
205     if (size != 1 && amode & WRITE) {
206       ierr = VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
207       ierr = VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
208     }
209     break;
210   default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT");
211   }
212   *xred = NULL;
213   PetscFunctionReturn(0);
214 }
215 
216 /* -------------------------------------------------------------------------- */
217 /*
218    PCApply_SVD - Applies the SVD preconditioner to a vector.
219 
220    Input Parameters:
221 .  pc - the preconditioner context
222 .  x - input vector
223 
224    Output Parameter:
225 .  y - output vector
226 
227    Application Interface Routine: PCApply()
228  */
229 #undef __FUNCT__
230 #define __FUNCT__ "PCApply_SVD"
231 static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y)
232 {
233   PC_SVD         *jac = (PC_SVD*)pc->data;
234   Vec            work = jac->work,xred,yred;
235   PetscErrorCode ierr;
236 
237   PetscFunctionBegin;
238   ierr = PCSVDGetVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr);
239   ierr = PCSVDGetVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr);
240 #if !defined(PETSC_USE_COMPLEX)
241   ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr);
242 #else
243   ierr = MatMultHermitianTranspose(jac->U,xred,work);CHKERRQ(ierr);
244 #endif
245   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
246 #if !defined(PETSC_USE_COMPLEX)
247   ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
248 #else
249   ierr = MatMultHermitianTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
250 #endif
251   ierr = PCSVDRestoreVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr);
252   ierr = PCSVDRestoreVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr);
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "PCApplyTranspose_SVD"
258 static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y)
259 {
260   PC_SVD         *jac = (PC_SVD*)pc->data;
261   Vec            work = jac->work,xred,yred;
262   PetscErrorCode ierr;
263 
264   PetscFunctionBegin;
265   ierr = PCSVDGetVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr);
266   ierr = PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr);
267 #if !defined(PETSC_USE_COMPLEX)
268   ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
269 #else
270   ierr = MatMultHermitianTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
271 #endif
272   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
273 #if !defined(PETSC_USE_COMPLEX)
274   ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr);
275 #else
276   ierr = MatMultHermitianTranspose(jac->U,xred,work);CHKERRQ(ierr);
277 #endif
278   ierr = PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr);
279   ierr = PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "PCReset_SVD"
285 static PetscErrorCode PCReset_SVD(PC pc)
286 {
287   PC_SVD         *jac = (PC_SVD*)pc->data;
288   PetscErrorCode ierr;
289 
290   PetscFunctionBegin;
291   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
292   ierr = MatDestroy(&jac->U);CHKERRQ(ierr);
293   ierr = MatDestroy(&jac->Vt);CHKERRQ(ierr);
294   ierr = VecDestroy(&jac->diag);CHKERRQ(ierr);
295   ierr = VecDestroy(&jac->work);CHKERRQ(ierr);
296   ierr = VecScatterDestroy(&jac->right2red);CHKERRQ(ierr);
297   ierr = VecScatterDestroy(&jac->left2red);CHKERRQ(ierr);
298   ierr = VecDestroy(&jac->rightred);CHKERRQ(ierr);
299   ierr = VecDestroy(&jac->leftred);CHKERRQ(ierr);
300   PetscFunctionReturn(0);
301 }
302 
303 /* -------------------------------------------------------------------------- */
304 /*
305    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
306    that was created with PCCreate_SVD().
307 
308    Input Parameter:
309 .  pc - the preconditioner context
310 
311    Application Interface Routine: PCDestroy()
312 */
313 #undef __FUNCT__
314 #define __FUNCT__ "PCDestroy_SVD"
315 static PetscErrorCode PCDestroy_SVD(PC pc)
316 {
317   PC_SVD         *jac = (PC_SVD*)pc->data;
318   PetscErrorCode ierr;
319 
320   PetscFunctionBegin;
321   ierr = PCReset_SVD(pc);CHKERRQ(ierr);
322   ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
323   ierr = PetscFree(pc->data);CHKERRQ(ierr);
324   PetscFunctionReturn(0);
325 }
326 
327 #undef __FUNCT__
328 #define __FUNCT__ "PCSetFromOptions_SVD"
329 static PetscErrorCode PCSetFromOptions_SVD(PetscOptions *PetscOptionsObject,PC pc)
330 {
331   PetscErrorCode ierr;
332   PC_SVD         *jac = (PC_SVD*)pc->data;
333   PetscBool      flg,set;
334 
335   PetscFunctionBegin;
336   ierr = PetscOptionsHead(PetscOptionsObject,"SVD options");CHKERRQ(ierr);
337   ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,NULL);CHKERRQ(ierr);
338   ierr = PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,NULL);CHKERRQ(ierr);
339   ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
340   if (set) {                    /* Should make PCSVDSetMonitor() */
341     if (flg && !jac->monitor) {
342       ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),"stdout",&jac->monitor);CHKERRQ(ierr);
343     } else if (!flg) {
344       ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
345     }
346   }
347   ierr = PetscOptionsTail();CHKERRQ(ierr);
348   PetscFunctionReturn(0);
349 }
350 
351 /* -------------------------------------------------------------------------- */
352 /*
353    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
354    and sets this as the private data within the generic preconditioning
355    context, PC, that was created within PCCreate().
356 
357    Input Parameter:
358 .  pc - the preconditioner context
359 
360    Application Interface Routine: PCCreate()
361 */
362 
363 /*MC
364      PCSVD - Use pseudo inverse defined by SVD of operator
365 
366    Level: advanced
367 
368   Concepts: SVD
369 
370   Options Database:
371 -  -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero
372 +  -pc_svd_monitor  Print information on the extreme singular values of the operator
373 
374   Developer Note: This implementation automatically creates a redundant copy of the
375    matrix on each process and uses a sequential SVD solve. Why does it do this instead
376    of using the composable PCREDUNDANT object?
377 
378 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
379 M*/
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "PCCreate_SVD"
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            = 0;
411   pc->ops->applyrichardson = 0;
412   PetscFunctionReturn(0);
413 }
414 
415