xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
2c6db04a5SJed Brown #include <petscblaslapack.h>
327c67122SBarry Smith 
427c67122SBarry Smith /*
527c67122SBarry Smith    Private context (data structure) for the SVD preconditioner.
627c67122SBarry Smith */
727c67122SBarry Smith typedef struct {
827c67122SBarry Smith   Vec               diag, work;
93ed27f31SJed Brown   Mat               A, U, Vt;
1027c67122SBarry Smith   PetscInt          nzero;
118f1a2a5eSBarry Smith   PetscReal         zerosing; /* measure of smallest singular value treated as nonzero */
1285032590SJed Brown   PetscInt          essrank;  /* essential rank of operator */
133ed27f31SJed Brown   VecScatter        left2red, right2red;
143ed27f31SJed Brown   Vec               leftred, rightred;
15426160bdSJed Brown   PetscViewer       monitor;
167962402dSFande Kong   PetscViewerFormat monitorformat;
1727c67122SBarry Smith } PC_SVD;
1827c67122SBarry Smith 
199371c9d4SSatish Balay typedef enum {
209371c9d4SSatish Balay   READ       = 1,
219371c9d4SSatish Balay   WRITE      = 2,
229371c9d4SSatish Balay   READ_WRITE = 3
239371c9d4SSatish Balay } AccessMode;
2427c67122SBarry Smith 
2527c67122SBarry Smith /*
2627c67122SBarry Smith    PCSetUp_SVD - Prepares for the use of the SVD preconditioner
2727c67122SBarry Smith                     by setting data structures and options.
2827c67122SBarry Smith 
2927c67122SBarry Smith    Input Parameter:
3027c67122SBarry Smith .  pc - the preconditioner context
3127c67122SBarry Smith 
3227c67122SBarry Smith    Application Interface Routine: PCSetUp()
3327c67122SBarry Smith 
34f1580f4eSBarry Smith    Note:
3527c67122SBarry Smith    The interface routine PCSetUp() is not usually called directly by
3627c67122SBarry Smith    the user, but instead is called by PCApply() if necessary.
3727c67122SBarry Smith */
38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_SVD(PC pc)
39d71ae5a4SJacob Faibussowitsch {
4027c67122SBarry Smith   PC_SVD      *jac = (PC_SVD *)pc->data;
4127c67122SBarry Smith   PetscScalar *a, *u, *v, *d, *work;
423f83f0d9SMatthew G Knepley   PetscBLASInt nb, lwork;
4327c67122SBarry Smith   PetscInt     i, n;
443ed27f31SJed Brown   PetscMPIInt  size;
4527c67122SBarry Smith 
4627c67122SBarry Smith   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->A));
489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(((PetscObject)pc->pmat)->comm, &size));
493ed27f31SJed Brown   if (size > 1) {
503ed27f31SJed Brown     Mat redmat;
5128d58a37SPierre Jolivet 
529566063dSJacob Faibussowitsch     PetscCall(MatCreateRedundantMatrix(pc->pmat, size, PETSC_COMM_SELF, MAT_INITIAL_MATRIX, &redmat));
539566063dSJacob Faibussowitsch     PetscCall(MatConvert(redmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A));
549566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&redmat));
553ed27f31SJed Brown   } else {
569566063dSJacob Faibussowitsch     PetscCall(MatConvert(pc->pmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A));
573ed27f31SJed Brown   }
583ed27f31SJed Brown   if (!jac->diag) { /* assume square matrices */
599566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(jac->A, &jac->diag, &jac->work));
603ed27f31SJed Brown   }
6127c67122SBarry Smith   if (!jac->U) {
629566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->U));
639566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->Vt));
6427c67122SBarry Smith   }
659566063dSJacob Faibussowitsch   PetscCall(MatGetSize(jac->A, &n, NULL));
6614ce09a1SBarry Smith   if (!n) {
679566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "Matrix has zero rows, skipping svd\n"));
683ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6914ce09a1SBarry Smith   }
709566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &nb));
7127c67122SBarry Smith   lwork = 5 * nb;
729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lwork, &work));
739566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->A, &a));
749566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->U, &u));
759566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->Vt, &v));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(jac->diag, &d));
7727c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX)
783f83f0d9SMatthew G Knepley   {
793f83f0d9SMatthew G Knepley     PetscBLASInt lierr;
809566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
81792fecdfSBarry Smith     PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, d, u, &nb, v, &nb, work, &lwork, &lierr));
82c9c947b5SJed Brown     PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesvd() error %" PetscBLASInt_FMT, lierr);
839566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
843f83f0d9SMatthew G Knepley   }
8527c67122SBarry Smith #else
86ef47b4b1SBarry Smith   {
87ef47b4b1SBarry Smith     PetscBLASInt lierr;
88ef47b4b1SBarry Smith     PetscReal   *rwork, *dd;
899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5 * nb, &rwork));
909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nb, &dd));
919566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
92792fecdfSBarry Smith     PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, dd, u, &nb, v, &nb, work, &lwork, rwork, &lierr));
9363a3b9bcSJacob Faibussowitsch     PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr);
949566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
9514ce09a1SBarry Smith     for (i = 0; i < n; i++) d[i] = dd[i];
969566063dSJacob Faibussowitsch     PetscCall(PetscFree(dd));
979566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
98ef47b4b1SBarry Smith   }
9927c67122SBarry Smith #endif
1009566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->A, &a));
1019566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->U, &u));
1029566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->Vt, &v));
1039371c9d4SSatish Balay   for (i = n - 1; i >= 0; i--)
1049371c9d4SSatish Balay     if (PetscRealPart(d[i]) > jac->zerosing) break;
105426160bdSJed Brown   jac->nzero = n - 1 - i;
106426160bdSJed Brown   if (jac->monitor) {
1079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(jac->monitor, ((PetscObject)pc)->tablevel));
10863a3b9bcSJacob 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));
1097962402dSFande Kong     if (n < 10 || jac->monitorformat == PETSC_VIEWER_ALL) {
1107962402dSFande Kong       PetscCall(PetscViewerASCIIPrintf(jac->monitor, "    SVD: singular values:\n"));
1117962402dSFande Kong       for (i = 0; i < n; i++) {
1127962402dSFande Kong         if (i % 5 == 0) {
11348a46eb9SPierre Jolivet           if (i != 0) PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n"));
1147962402dSFande Kong           PetscCall(PetscViewerASCIIPrintf(jac->monitor, "        "));
1157962402dSFande Kong         }
1167962402dSFande Kong         PetscCall(PetscViewerASCIIPrintf(jac->monitor, " %14.12e", (double)PetscRealPart(d[i])));
1177962402dSFande Kong       }
1187962402dSFande Kong       PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n"));
1197962402dSFande Kong     } else { /* print 5 smallest and 5 largest */
1209566063dSJacob 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])));
1219566063dSJacob 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])));
122426160bdSJed Brown     }
1239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(jac->monitor, ((PetscObject)pc)->tablevel));
124426160bdSJed Brown   }
1259566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "Largest and smallest singular values %14.12e %14.12e\n", (double)PetscRealPart(d[0]), (double)PetscRealPart(d[n - 1])));
126426160bdSJed Brown   for (i = 0; i < n - jac->nzero; i++) d[i] = 1.0 / d[i];
127426160bdSJed Brown   for (; i < n; i++) d[i] = 0.0;
1289371c9d4SSatish Balay   if (jac->essrank > 0)
1299371c9d4SSatish Balay     for (i = 0; i < n - jac->nzero - jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */
13063a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "Number of zero or nearly singular values %" PetscInt_FMT "\n", jac->nzero));
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(jac->diag, &d));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13427c67122SBarry Smith }
13527c67122SBarry Smith 
136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSVDGetVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred)
137d71ae5a4SJacob Faibussowitsch {
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;
167d71ae5a4SJacob Faibussowitsch   default:
168d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT");
1693ed27f31SJed Brown   }
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1713ed27f31SJed Brown }
1723ed27f31SJed Brown 
173d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSVDRestoreVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred)
174d71ae5a4SJacob Faibussowitsch {
1753ed27f31SJed Brown   PC_SVD     *jac = (PC_SVD *)pc->data;
1763ed27f31SJed Brown   PetscMPIInt size;
1773ed27f31SJed Brown 
1783ed27f31SJed Brown   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
1803ed27f31SJed Brown   switch (side) {
1813ed27f31SJed Brown   case PC_LEFT:
1823ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
1839566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE));
1849566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE));
1853ed27f31SJed Brown     }
1863ed27f31SJed Brown     break;
1873ed27f31SJed Brown   case PC_RIGHT:
1883ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
1899566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE));
1909566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE));
1913ed27f31SJed Brown     }
1923ed27f31SJed Brown     break;
193d71ae5a4SJacob Faibussowitsch   default:
194d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT");
1953ed27f31SJed Brown   }
1960298fd71SBarry Smith   *xred = NULL;
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1983ed27f31SJed Brown }
1993ed27f31SJed Brown 
20027c67122SBarry Smith /*
20127c67122SBarry Smith    PCApply_SVD - Applies the SVD preconditioner to a vector.
20227c67122SBarry Smith 
20327c67122SBarry Smith    Input Parameters:
20427c67122SBarry Smith .  pc - the preconditioner context
20527c67122SBarry Smith .  x - input vector
20627c67122SBarry Smith 
20727c67122SBarry Smith    Output Parameter:
20827c67122SBarry Smith .  y - output vector
20927c67122SBarry Smith 
21027c67122SBarry Smith    Application Interface Routine: PCApply()
21127c67122SBarry Smith  */
212d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SVD(PC pc, Vec x, Vec y)
213d71ae5a4SJacob Faibussowitsch {
21427c67122SBarry Smith   PC_SVD *jac  = (PC_SVD *)pc->data;
2153ed27f31SJed Brown   Vec     work = jac->work, xred, yred;
21627c67122SBarry Smith 
21727c67122SBarry Smith   PetscFunctionBegin;
2189566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_RIGHT, READ, x, &xred));
2199566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_LEFT, WRITE, y, &yred));
220ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2219566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(jac->U, xred, work));
222ef47b4b1SBarry Smith #else
2239566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(jac->U, xred, work));
224ef47b4b1SBarry Smith #endif
2259566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(work, work, jac->diag));
226ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2279566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(jac->Vt, work, yred));
228ef47b4b1SBarry Smith #else
2299566063dSJacob Faibussowitsch   PetscCall(MatMultHermitianTranspose(jac->Vt, work, yred));
230ef47b4b1SBarry Smith #endif
2319566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, READ, x, &xred));
2329566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_LEFT, WRITE, y, &yred));
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2343ed27f31SJed Brown }
2353ed27f31SJed Brown 
236d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_SVD(PC pc, Mat X, Mat Y)
237d71ae5a4SJacob Faibussowitsch {
238024db992SStefano Zampini   PC_SVD *jac = (PC_SVD *)pc->data;
239024db992SStefano Zampini   Mat     W;
240024db992SStefano Zampini 
241024db992SStefano Zampini   PetscFunctionBegin;
242024db992SStefano Zampini   PetscCall(MatTransposeMatMult(jac->U, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &W));
243024db992SStefano Zampini   PetscCall(MatDiagonalScale(W, jac->diag, NULL));
244024db992SStefano Zampini   PetscCall(MatTransposeMatMult(jac->Vt, W, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
245024db992SStefano Zampini   PetscCall(MatDestroy(&W));
2463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
247024db992SStefano Zampini }
248024db992SStefano Zampini 
249d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_SVD(PC pc, Vec x, Vec y)
250d71ae5a4SJacob Faibussowitsch {
2513ed27f31SJed Brown   PC_SVD *jac  = (PC_SVD *)pc->data;
2523ed27f31SJed Brown   Vec     work = jac->work, xred, yred;
2533ed27f31SJed Brown 
2543ed27f31SJed Brown   PetscFunctionBegin;
2559566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_LEFT, READ, x, &xred));
2569566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc, PC_RIGHT, WRITE, y, &yred));
2579566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->Vt, xred, work));
2589566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(work, work, jac->diag));
2599566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->U, work, yred));
2609566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_LEFT, READ, x, &xred));
2619566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, WRITE, y, &yred));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26327c67122SBarry Smith }
26427c67122SBarry Smith 
265d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_SVD(PC pc)
266d71ae5a4SJacob Faibussowitsch {
267a2d70de2SBarry Smith   PC_SVD *jac = (PC_SVD *)pc->data;
268a2d70de2SBarry Smith 
269a2d70de2SBarry Smith   PetscFunctionBegin;
2709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->A));
2719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->U));
2729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->Vt));
2739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
2749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->work));
2759566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&jac->right2red));
2769566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&jac->left2red));
2779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->rightred));
2789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->leftred));
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
280a2d70de2SBarry Smith }
281a2d70de2SBarry Smith 
28227c67122SBarry Smith /*
28327c67122SBarry Smith    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
28427c67122SBarry Smith    that was created with PCCreate_SVD().
28527c67122SBarry Smith 
28627c67122SBarry Smith    Input Parameter:
28727c67122SBarry Smith .  pc - the preconditioner context
28827c67122SBarry Smith 
28927c67122SBarry Smith    Application Interface Routine: PCDestroy()
29027c67122SBarry Smith */
291d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SVD(PC pc)
292d71ae5a4SJacob Faibussowitsch {
293426160bdSJed Brown   PC_SVD *jac = (PC_SVD *)pc->data;
29427c67122SBarry Smith 
29527c67122SBarry Smith   PetscFunctionBegin;
2969566063dSJacob Faibussowitsch   PetscCall(PCReset_SVD(pc));
2979566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&jac->monitor));
2989566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30027c67122SBarry Smith }
30127c67122SBarry Smith 
302d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_SVD(PC pc, PetscOptionItems *PetscOptionsObject)
303d71ae5a4SJacob Faibussowitsch {
3048f1a2a5eSBarry Smith   PC_SVD   *jac = (PC_SVD *)pc->data;
3057962402dSFande Kong   PetscBool flg;
30627c67122SBarry Smith 
30727c67122SBarry Smith   PetscFunctionBegin;
308d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SVD options");
3099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_svd_zero_sing", "Singular values smaller than this treated as zero", "None", jac->zerosing, &jac->zerosing, NULL));
3109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_svd_ess_rank", "Essential rank of operator (0 to use entire operator)", "None", jac->essrank, &jac->essrank, NULL));
3117962402dSFande Kong   PetscCall(PetscOptionsViewer("-pc_svd_monitor", "Monitor the conditioning, and extremal singular values", "None", &jac->monitor, &jac->monitorformat, &flg));
312d0609cedSBarry Smith   PetscOptionsHeadEnd();
3133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31427c67122SBarry Smith }
31527c67122SBarry Smith 
316d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_SVD(PC pc, PetscViewer viewer)
317d71ae5a4SJacob Faibussowitsch {
31833761216SBarry Smith   PC_SVD   *svd = (PC_SVD *)pc->data;
31933761216SBarry Smith   PetscBool iascii;
32033761216SBarry Smith 
32133761216SBarry Smith   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
32333761216SBarry Smith   if (iascii) {
3249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  All singular values smaller than %g treated as zero\n", (double)svd->zerosing));
32563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Provided essential rank of the matrix %" PetscInt_FMT " (all other eigenvalues are zeroed)\n", svd->essrank));
32633761216SBarry Smith   }
3273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32833761216SBarry Smith }
329f1580f4eSBarry Smith 
33027c67122SBarry Smith /*
33127c67122SBarry Smith    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
33227c67122SBarry Smith    and sets this as the private data within the generic preconditioning
33327c67122SBarry Smith    context, PC, that was created within PCCreate().
33427c67122SBarry Smith 
33527c67122SBarry Smith    Input Parameter:
33627c67122SBarry Smith .  pc - the preconditioner context
33727c67122SBarry Smith 
33827c67122SBarry Smith    Application Interface Routine: PCCreate()
33927c67122SBarry Smith */
34027c67122SBarry Smith 
34127c67122SBarry Smith /*MC
34227c67122SBarry Smith      PCSVD - Use pseudo inverse defined by SVD of operator
34327c67122SBarry Smith 
34427c67122SBarry Smith    Level: advanced
34527c67122SBarry Smith 
346f1580f4eSBarry Smith   Options Database Keys:
347147403d9SBarry Smith +  -pc_svd_zero_sing <rtol> - Singular values smaller than this are treated as zero
348147403d9SBarry Smith -  -pc_svd_monitor - Print information on the extreme singular values of the operator
34927c67122SBarry Smith 
350a2b725a8SWilliam Gropp   Developer Note:
351a2b725a8SWilliam Gropp   This implementation automatically creates a redundant copy of the
3528997ae2eSBarry Smith    matrix on each process and uses a sequential SVD solve. Why does it do this instead
353f1580f4eSBarry Smith    of using the composable `PCREDUNDANT` object?
3548997ae2eSBarry Smith 
355*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCREDUNDANT`
35627c67122SBarry Smith M*/
35727c67122SBarry Smith 
358d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc)
359d71ae5a4SJacob Faibussowitsch {
36027c67122SBarry Smith   PC_SVD     *jac;
361024db992SStefano Zampini   PetscMPIInt size = 0;
36227c67122SBarry Smith 
36327c67122SBarry Smith   PetscFunctionBegin;
36427c67122SBarry Smith   /*
36527c67122SBarry Smith      Creates the private data structure for this preconditioner and
36627c67122SBarry Smith      attach it to the PC object.
36727c67122SBarry Smith   */
3684dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
3698f1a2a5eSBarry Smith   jac->zerosing = 1.e-12;
37027c67122SBarry Smith   pc->data      = (void *)jac;
37127c67122SBarry Smith 
37227c67122SBarry Smith   /*
37327c67122SBarry Smith       Set the pointers for the functions that are provided above.
37427c67122SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
37527c67122SBarry Smith       are called, they will automatically call these functions.  Note we
37627c67122SBarry Smith       choose not to provide a couple of these functions since they are
37727c67122SBarry Smith       not needed.
37827c67122SBarry Smith   */
379024db992SStefano Zampini 
380024db992SStefano Zampini #if defined(PETSC_HAVE_COMPLEX)
381024db992SStefano Zampini   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
382024db992SStefano Zampini #endif
383024db992SStefano Zampini   if (size == 1) pc->ops->matapply = PCMatApply_SVD;
38427c67122SBarry Smith   pc->ops->apply           = PCApply_SVD;
3853ed27f31SJed Brown   pc->ops->applytranspose  = PCApplyTranspose_SVD;
38627c67122SBarry Smith   pc->ops->setup           = PCSetUp_SVD;
387a2d70de2SBarry Smith   pc->ops->reset           = PCReset_SVD;
38827c67122SBarry Smith   pc->ops->destroy         = PCDestroy_SVD;
38927c67122SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SVD;
39033761216SBarry Smith   pc->ops->view            = PCView_SVD;
3910a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39327c67122SBarry Smith }
394