xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
127c67122SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
3c6db04a5SJed Brown #include <petscblaslapack.h>
427c67122SBarry Smith 
527c67122SBarry Smith /*
627c67122SBarry Smith    Private context (data structure) for the SVD preconditioner.
727c67122SBarry Smith */
827c67122SBarry Smith typedef struct {
927c67122SBarry Smith   Vec               diag, work;
103ed27f31SJed Brown   Mat               A, U, Vt;
1127c67122SBarry Smith   PetscInt          nzero;
128f1a2a5eSBarry Smith   PetscReal         zerosing; /* measure of smallest singular value treated as nonzero */
1385032590SJed Brown   PetscInt          essrank;  /* essential rank of operator */
143ed27f31SJed Brown   VecScatter        left2red, right2red;
153ed27f31SJed Brown   Vec               leftred, rightred;
16426160bdSJed Brown   PetscViewer       monitor;
177962402dSFande Kong   PetscViewerFormat monitorformat;
1827c67122SBarry Smith } PC_SVD;
1927c67122SBarry Smith 
209371c9d4SSatish Balay typedef enum {
219371c9d4SSatish Balay   READ       = 1,
229371c9d4SSatish Balay   WRITE      = 2,
239371c9d4SSatish Balay   READ_WRITE = 3
249371c9d4SSatish Balay } AccessMode;
2527c67122SBarry Smith 
2627c67122SBarry Smith /* -------------------------------------------------------------------------- */
2727c67122SBarry Smith /*
2827c67122SBarry Smith    PCSetUp_SVD - Prepares for the use of the SVD preconditioner
2927c67122SBarry Smith                     by setting data structures and options.
3027c67122SBarry Smith 
3127c67122SBarry Smith    Input Parameter:
3227c67122SBarry Smith .  pc - the preconditioner context
3327c67122SBarry Smith 
3427c67122SBarry Smith    Application Interface Routine: PCSetUp()
3527c67122SBarry Smith 
3627c67122SBarry Smith    Notes:
3727c67122SBarry Smith    The interface routine PCSetUp() is not usually called directly by
3827c67122SBarry Smith    the user, but instead is called by PCApply() if necessary.
3927c67122SBarry Smith */
409371c9d4SSatish Balay static PetscErrorCode PCSetUp_SVD(PC pc) {
4127c67122SBarry Smith   PC_SVD      *jac = (PC_SVD *)pc->data;
4227c67122SBarry Smith   PetscScalar *a, *u, *v, *d, *work;
433f83f0d9SMatthew G Knepley   PetscBLASInt nb, lwork;
4427c67122SBarry Smith   PetscInt     i, n;
453ed27f31SJed Brown   PetscMPIInt  size;
4627c67122SBarry Smith 
4727c67122SBarry Smith   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->A));
499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(((PetscObject)pc->pmat)->comm, &size));
503ed27f31SJed Brown   if (size > 1) {
513ed27f31SJed Brown     Mat redmat;
5228d58a37SPierre Jolivet 
539566063dSJacob Faibussowitsch     PetscCall(MatCreateRedundantMatrix(pc->pmat, size, PETSC_COMM_SELF, MAT_INITIAL_MATRIX, &redmat));
549566063dSJacob Faibussowitsch     PetscCall(MatConvert(redmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A));
559566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&redmat));
563ed27f31SJed Brown   } else {
579566063dSJacob Faibussowitsch     PetscCall(MatConvert(pc->pmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A));
583ed27f31SJed Brown   }
593ed27f31SJed Brown   if (!jac->diag) { /* assume square matrices */
609566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(jac->A, &jac->diag, &jac->work));
613ed27f31SJed Brown   }
6227c67122SBarry Smith   if (!jac->U) {
639566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->U));
649566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->Vt));
6527c67122SBarry Smith   }
669566063dSJacob Faibussowitsch   PetscCall(MatGetSize(jac->A, &n, NULL));
6714ce09a1SBarry Smith   if (!n) {
689566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Matrix has zero rows, skipping svd\n"));
6914ce09a1SBarry Smith     PetscFunctionReturn(0);
7014ce09a1SBarry Smith   }
719566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &nb));
7227c67122SBarry Smith   lwork = 5 * nb;
739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lwork, &work));
749566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->A, &a));
759566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->U, &u));
769566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->Vt, &v));
779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(jac->diag, &d));
7827c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX)
793f83f0d9SMatthew G Knepley   {
803f83f0d9SMatthew G Knepley     PetscBLASInt lierr;
819566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
82792fecdfSBarry Smith     PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, d, u, &nb, v, &nb, work, &lwork, &lierr));
8363a3b9bcSJacob Faibussowitsch     PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr);
849566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
853f83f0d9SMatthew G Knepley   }
8627c67122SBarry Smith #else
87ef47b4b1SBarry Smith   {
88ef47b4b1SBarry Smith     PetscBLASInt lierr;
89ef47b4b1SBarry Smith     PetscReal   *rwork, *dd;
909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * nb, &rwork));
919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nb, &dd));
929566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
93792fecdfSBarry Smith     PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, dd, u, &nb, v, &nb, work, &lwork, rwork, &lierr));
9463a3b9bcSJacob Faibussowitsch     PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr);
959566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
9614ce09a1SBarry Smith     for (i = 0; i < n; i++) d[i] = dd[i];
979566063dSJacob Faibussowitsch     PetscCall(PetscFree(dd));
989566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
99ef47b4b1SBarry Smith   }
10027c67122SBarry Smith #endif
1019566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->A, &a));
1029566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->U, &u));
1039566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->Vt, &v));
1049371c9d4SSatish Balay   for (i = n - 1; i >= 0; i--)
1059371c9d4SSatish Balay     if (PetscRealPart(d[i]) > jac->zerosing) break;
106426160bdSJed Brown   jac->nzero = n - 1 - i;
107426160bdSJed Brown   if (jac->monitor) {
1089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(jac->monitor, ((PetscObject)pc)->tablevel));
10963a3b9bcSJacob Faibussowitsch     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));
1107962402dSFande Kong     if (n < 10 || jac->monitorformat == PETSC_VIEWER_ALL) {
1117962402dSFande Kong       PetscCall(PetscViewerASCIIPrintf(jac->monitor, "    SVD: singular values:\n"));
1127962402dSFande Kong       for (i = 0; i < n; i++) {
1137962402dSFande Kong         if (i % 5 == 0) {
114*48a46eb9SPierre Jolivet           if (i != 0) PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n"));
1157962402dSFande Kong           PetscCall(PetscViewerASCIIPrintf(jac->monitor, "        "));
1167962402dSFande Kong         }
1177962402dSFande Kong         PetscCall(PetscViewerASCIIPrintf(jac->monitor, " %14.12e", (double)PetscRealPart(d[i])));
1187962402dSFande Kong       }
1197962402dSFande Kong       PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n"));
1207962402dSFande Kong     } else { /* print 5 smallest and 5 largest */
1219566063dSJacob Faibussowitsch       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])));
1229566063dSJacob Faibussowitsch       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])));
123426160bdSJed Brown     }
1249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(jac->monitor, ((PetscObject)pc)->tablevel));
125426160bdSJed Brown   }
1269566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "Largest and smallest singular values %14.12e %14.12e\n", (double)PetscRealPart(d[0]), (double)PetscRealPart(d[n - 1])));
127426160bdSJed Brown   for (i = 0; i < n - jac->nzero; i++) d[i] = 1.0 / d[i];
128426160bdSJed Brown   for (; i < n; i++) d[i] = 0.0;
1299371c9d4SSatish Balay   if (jac->essrank > 0)
1309371c9d4SSatish Balay     for (i = 0; i < n - jac->nzero - jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */
13163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "Number of zero or nearly singular values %" PetscInt_FMT "\n", jac->nzero));
1329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(jac->diag, &d));
1339566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
13427c67122SBarry Smith   PetscFunctionReturn(0);
13527c67122SBarry Smith }
13627c67122SBarry Smith 
1379371c9d4SSatish Balay static PetscErrorCode PCSVDGetVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) {
1383ed27f31SJed Brown   PC_SVD     *jac = (PC_SVD *)pc->data;
1393ed27f31SJed Brown   PetscMPIInt size;
1403ed27f31SJed Brown 
1413ed27f31SJed Brown   PetscFunctionBegin;
1429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
1430298fd71SBarry Smith   *xred = NULL;
1443ed27f31SJed Brown   switch (side) {
1453ed27f31SJed Brown   case PC_LEFT:
1463ed27f31SJed Brown     if (size == 1) *xred = x;
1473ed27f31SJed Brown     else {
1489566063dSJacob Faibussowitsch       if (!jac->left2red) PetscCall(VecScatterCreateToAll(x, &jac->left2red, &jac->leftred));
1493ed27f31SJed Brown       if (amode & READ) {
1509566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD));
1519566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD));
1523ed27f31SJed Brown       }
1533ed27f31SJed Brown       *xred = jac->leftred;
1543ed27f31SJed Brown     }
1553ed27f31SJed Brown     break;
1563ed27f31SJed Brown   case PC_RIGHT:
1573ed27f31SJed Brown     if (size == 1) *xred = x;
1583ed27f31SJed Brown     else {
1599566063dSJacob Faibussowitsch       if (!jac->right2red) PetscCall(VecScatterCreateToAll(x, &jac->right2red, &jac->rightred));
1603ed27f31SJed Brown       if (amode & READ) {
1619566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD));
1629566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD));
1633ed27f31SJed Brown       }
1643ed27f31SJed Brown       *xred = jac->rightred;
1653ed27f31SJed Brown     }
1663ed27f31SJed Brown     break;
167ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT");
1683ed27f31SJed Brown   }
1693ed27f31SJed Brown   PetscFunctionReturn(0);
1703ed27f31SJed Brown }
1713ed27f31SJed Brown 
1729371c9d4SSatish Balay static PetscErrorCode PCSVDRestoreVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) {
1733ed27f31SJed Brown   PC_SVD     *jac = (PC_SVD *)pc->data;
1743ed27f31SJed Brown   PetscMPIInt size;
1753ed27f31SJed Brown 
1763ed27f31SJed Brown   PetscFunctionBegin;
1779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
1783ed27f31SJed Brown   switch (side) {
1793ed27f31SJed Brown   case PC_LEFT:
1803ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
1819566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE));
1829566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE));
1833ed27f31SJed Brown     }
1843ed27f31SJed Brown     break;
1853ed27f31SJed Brown   case PC_RIGHT:
1863ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
1879566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE));
1889566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE));
1893ed27f31SJed Brown     }
1903ed27f31SJed Brown     break;
191ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT");
1923ed27f31SJed Brown   }
1930298fd71SBarry Smith   *xred = NULL;
1943ed27f31SJed Brown   PetscFunctionReturn(0);
1953ed27f31SJed Brown }
1963ed27f31SJed Brown 
19727c67122SBarry Smith /* -------------------------------------------------------------------------- */
19827c67122SBarry Smith /*
19927c67122SBarry Smith    PCApply_SVD - Applies the SVD preconditioner to a vector.
20027c67122SBarry Smith 
20127c67122SBarry Smith    Input Parameters:
20227c67122SBarry Smith .  pc - the preconditioner context
20327c67122SBarry Smith .  x - input vector
20427c67122SBarry Smith 
20527c67122SBarry Smith    Output Parameter:
20627c67122SBarry Smith .  y - output vector
20727c67122SBarry Smith 
20827c67122SBarry Smith    Application Interface Routine: PCApply()
20927c67122SBarry Smith  */
2109371c9d4SSatish Balay static PetscErrorCode PCApply_SVD(PC pc, Vec x, Vec y) {
21127c67122SBarry Smith   PC_SVD *jac  = (PC_SVD *)pc->data;
2123ed27f31SJed Brown   Vec     work = jac->work, xred, yred;
21327c67122SBarry Smith 
21427c67122SBarry Smith   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_RIGHT, READ, x, &xred));
2169566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_LEFT, WRITE, y, &yred));
217ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2189566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(jac->U, xred, work));
219ef47b4b1SBarry Smith #else
2209566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(jac->U, xred, work));
221ef47b4b1SBarry Smith #endif
2229566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(work, work, jac->diag));
223ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2249566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(jac->Vt, work, yred));
225ef47b4b1SBarry Smith #else
2269566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(jac->Vt, work, yred));
227ef47b4b1SBarry Smith #endif
2289566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, READ, x, &xred));
2299566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_LEFT, WRITE, y, &yred));
2303ed27f31SJed Brown   PetscFunctionReturn(0);
2313ed27f31SJed Brown }
2323ed27f31SJed Brown 
2339371c9d4SSatish Balay static PetscErrorCode PCMatApply_SVD(PC pc, Mat X, Mat Y) {
234024db992SStefano Zampini   PC_SVD *jac = (PC_SVD *)pc->data;
235024db992SStefano Zampini   Mat     W;
236024db992SStefano Zampini 
237024db992SStefano Zampini   PetscFunctionBegin;
238024db992SStefano Zampini   PetscCall(MatTransposeMatMult(jac->U, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &W));
239024db992SStefano Zampini   PetscCall(MatDiagonalScale(W, jac->diag, NULL));
240024db992SStefano Zampini   PetscCall(MatTransposeMatMult(jac->Vt, W, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
241024db992SStefano Zampini   PetscCall(MatDestroy(&W));
242024db992SStefano Zampini   PetscFunctionReturn(0);
243024db992SStefano Zampini }
244024db992SStefano Zampini 
2459371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_SVD(PC pc, Vec x, Vec y) {
2463ed27f31SJed Brown   PC_SVD *jac  = (PC_SVD *)pc->data;
2473ed27f31SJed Brown   Vec     work = jac->work, xred, yred;
2483ed27f31SJed Brown 
2493ed27f31SJed Brown   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_LEFT, READ, x, &xred));
2519566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_RIGHT, WRITE, y, &yred));
2529566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->Vt, xred, work));
2539566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(work, work, jac->diag));
2549566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->U, work, yred));
2559566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_LEFT, READ, x, &xred));
2569566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, WRITE, y, &yred));
25727c67122SBarry Smith   PetscFunctionReturn(0);
25827c67122SBarry Smith }
25927c67122SBarry Smith 
2609371c9d4SSatish Balay static PetscErrorCode PCReset_SVD(PC pc) {
261a2d70de2SBarry Smith   PC_SVD *jac = (PC_SVD *)pc->data;
262a2d70de2SBarry Smith 
263a2d70de2SBarry Smith   PetscFunctionBegin;
2649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->A));
2659566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->U));
2669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->Vt));
2679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
2689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->work));
2699566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&jac->right2red));
2709566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&jac->left2red));
2719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->rightred));
2729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->leftred));
273a2d70de2SBarry Smith   PetscFunctionReturn(0);
274a2d70de2SBarry Smith }
275a2d70de2SBarry Smith 
27627c67122SBarry Smith /* -------------------------------------------------------------------------- */
27727c67122SBarry Smith /*
27827c67122SBarry Smith    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
27927c67122SBarry Smith    that was created with PCCreate_SVD().
28027c67122SBarry Smith 
28127c67122SBarry Smith    Input Parameter:
28227c67122SBarry Smith .  pc - the preconditioner context
28327c67122SBarry Smith 
28427c67122SBarry Smith    Application Interface Routine: PCDestroy()
28527c67122SBarry Smith */
2869371c9d4SSatish Balay static PetscErrorCode PCDestroy_SVD(PC pc) {
287426160bdSJed Brown   PC_SVD *jac = (PC_SVD *)pc->data;
28827c67122SBarry Smith 
28927c67122SBarry Smith   PetscFunctionBegin;
2909566063dSJacob Faibussowitsch   PetscCall(PCReset_SVD(pc));
2919566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&jac->monitor));
2929566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
29327c67122SBarry Smith   PetscFunctionReturn(0);
29427c67122SBarry Smith }
29527c67122SBarry Smith 
2969371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_SVD(PC pc, PetscOptionItems *PetscOptionsObject) {
2978f1a2a5eSBarry Smith   PC_SVD   *jac = (PC_SVD *)pc->data;
2987962402dSFande Kong   PetscBool flg;
29927c67122SBarry Smith 
30027c67122SBarry Smith   PetscFunctionBegin;
301d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SVD options");
3029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_svd_zero_sing", "Singular values smaller than this treated as zero", "None", jac->zerosing, &jac->zerosing, NULL));
3039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_svd_ess_rank", "Essential rank of operator (0 to use entire operator)", "None", jac->essrank, &jac->essrank, NULL));
3047962402dSFande Kong   PetscCall(PetscOptionsViewer("-pc_svd_monitor", "Monitor the conditioning, and extremal singular values", "None", &jac->monitor, &jac->monitorformat, &flg));
305d0609cedSBarry Smith   PetscOptionsHeadEnd();
30627c67122SBarry Smith   PetscFunctionReturn(0);
30727c67122SBarry Smith }
30827c67122SBarry Smith 
3099371c9d4SSatish Balay static PetscErrorCode PCView_SVD(PC pc, PetscViewer viewer) {
31033761216SBarry Smith   PC_SVD   *svd = (PC_SVD *)pc->data;
31133761216SBarry Smith   PetscBool iascii;
31233761216SBarry Smith 
31333761216SBarry Smith   PetscFunctionBegin;
3149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
31533761216SBarry Smith   if (iascii) {
3169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  All singular values smaller than %g treated as zero\n", (double)svd->zerosing));
31763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Provided essential rank of the matrix %" PetscInt_FMT " (all other eigenvalues are zeroed)\n", svd->essrank));
31833761216SBarry Smith   }
31933761216SBarry Smith   PetscFunctionReturn(0);
32033761216SBarry Smith }
32127c67122SBarry Smith /* -------------------------------------------------------------------------- */
32227c67122SBarry Smith /*
32327c67122SBarry Smith    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
32427c67122SBarry Smith    and sets this as the private data within the generic preconditioning
32527c67122SBarry Smith    context, PC, that was created within PCCreate().
32627c67122SBarry Smith 
32727c67122SBarry Smith    Input Parameter:
32827c67122SBarry Smith .  pc - the preconditioner context
32927c67122SBarry Smith 
33027c67122SBarry Smith    Application Interface Routine: PCCreate()
33127c67122SBarry Smith */
33227c67122SBarry Smith 
33327c67122SBarry Smith /*MC
33427c67122SBarry Smith      PCSVD - Use pseudo inverse defined by SVD of operator
33527c67122SBarry Smith 
33627c67122SBarry Smith    Level: advanced
33727c67122SBarry Smith 
3388f1a2a5eSBarry Smith   Options Database:
339147403d9SBarry Smith +  -pc_svd_zero_sing <rtol> - Singular values smaller than this are treated as zero
340147403d9SBarry Smith -  -pc_svd_monitor - Print information on the extreme singular values of the operator
34127c67122SBarry Smith 
342a2b725a8SWilliam Gropp   Developer Note:
343a2b725a8SWilliam Gropp   This implementation automatically creates a redundant copy of the
3448997ae2eSBarry Smith    matrix on each process and uses a sequential SVD solve. Why does it do this instead
3458997ae2eSBarry Smith    of using the composable PCREDUNDANT object?
3468997ae2eSBarry Smith 
347db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`
34827c67122SBarry Smith M*/
34927c67122SBarry Smith 
3509371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc) {
35127c67122SBarry Smith   PC_SVD     *jac;
352024db992SStefano Zampini   PetscMPIInt size = 0;
35327c67122SBarry Smith 
35427c67122SBarry Smith   PetscFunctionBegin;
35527c67122SBarry Smith   /*
35627c67122SBarry Smith      Creates the private data structure for this preconditioner and
35727c67122SBarry Smith      attach it to the PC object.
35827c67122SBarry Smith   */
3599566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &jac));
3608f1a2a5eSBarry Smith   jac->zerosing = 1.e-12;
36127c67122SBarry Smith   pc->data      = (void *)jac;
36227c67122SBarry Smith 
36327c67122SBarry Smith   /*
36427c67122SBarry Smith       Set the pointers for the functions that are provided above.
36527c67122SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
36627c67122SBarry Smith       are called, they will automatically call these functions.  Note we
36727c67122SBarry Smith       choose not to provide a couple of these functions since they are
36827c67122SBarry Smith       not needed.
36927c67122SBarry Smith   */
370024db992SStefano Zampini 
371024db992SStefano Zampini #if defined(PETSC_HAVE_COMPLEX)
372024db992SStefano Zampini   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
373024db992SStefano Zampini #endif
374024db992SStefano Zampini   if (size == 1) pc->ops->matapply = PCMatApply_SVD;
37527c67122SBarry Smith   pc->ops->apply           = PCApply_SVD;
3763ed27f31SJed Brown   pc->ops->applytranspose  = PCApplyTranspose_SVD;
37727c67122SBarry Smith   pc->ops->setup           = PCSetUp_SVD;
378a2d70de2SBarry Smith   pc->ops->reset           = PCReset_SVD;
37927c67122SBarry Smith   pc->ops->destroy         = PCDestroy_SVD;
38027c67122SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SVD;
38133761216SBarry Smith   pc->ops->view            = PCView_SVD;
3820a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
38327c67122SBarry Smith   PetscFunctionReturn(0);
38427c67122SBarry Smith }
385