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