xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
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 
203ed27f31SJed Brown typedef enum {READ=1, WRITE=2, READ_WRITE=3} AccessMode;
2127c67122SBarry Smith 
2227c67122SBarry Smith /* -------------------------------------------------------------------------- */
2327c67122SBarry Smith /*
2427c67122SBarry Smith    PCSetUp_SVD - Prepares for the use of the SVD preconditioner
2527c67122SBarry Smith                     by setting data structures and options.
2627c67122SBarry Smith 
2727c67122SBarry Smith    Input Parameter:
2827c67122SBarry Smith .  pc - the preconditioner context
2927c67122SBarry Smith 
3027c67122SBarry Smith    Application Interface Routine: PCSetUp()
3127c67122SBarry Smith 
3227c67122SBarry Smith    Notes:
3327c67122SBarry Smith    The interface routine PCSetUp() is not usually called directly by
3427c67122SBarry Smith    the user, but instead is called by PCApply() if necessary.
3527c67122SBarry Smith */
3627c67122SBarry Smith static PetscErrorCode PCSetUp_SVD(PC pc)
3727c67122SBarry Smith {
3827c67122SBarry Smith   PC_SVD         *jac = (PC_SVD*)pc->data;
3927c67122SBarry Smith   PetscScalar    *a,*u,*v,*d,*work;
403f83f0d9SMatthew G Knepley   PetscBLASInt   nb,lwork;
4127c67122SBarry Smith   PetscInt       i,n;
423ed27f31SJed Brown   PetscMPIInt    size;
4327c67122SBarry Smith 
4427c67122SBarry Smith   PetscFunctionBegin;
459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->A));
469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(((PetscObject)pc->pmat)->comm,&size));
473ed27f31SJed Brown   if (size > 1) {
483ed27f31SJed Brown     Mat redmat;
4928d58a37SPierre Jolivet 
509566063dSJacob Faibussowitsch     PetscCall(MatCreateRedundantMatrix(pc->pmat,size,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&redmat));
519566063dSJacob Faibussowitsch     PetscCall(MatConvert(redmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A));
529566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&redmat));
533ed27f31SJed Brown   } else {
549566063dSJacob Faibussowitsch     PetscCall(MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A));
553ed27f31SJed Brown   }
563ed27f31SJed Brown   if (!jac->diag) {    /* assume square matrices */
579566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(jac->A,&jac->diag,&jac->work));
583ed27f31SJed Brown   }
5927c67122SBarry Smith   if (!jac->U) {
609566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U));
619566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->Vt));
6227c67122SBarry Smith   }
639566063dSJacob Faibussowitsch   PetscCall(MatGetSize(jac->A,&n,NULL));
6414ce09a1SBarry Smith   if (!n) {
659566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"Matrix has zero rows, skipping svd\n"));
6614ce09a1SBarry Smith     PetscFunctionReturn(0);
6714ce09a1SBarry Smith   }
689566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n,&nb));
6927c67122SBarry Smith   lwork = 5*nb;
709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(lwork,&work));
719566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->A,&a));
729566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->U,&u));
739566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(jac->Vt,&v));
749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(jac->diag,&d));
7527c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX)
763f83f0d9SMatthew G Knepley   {
773f83f0d9SMatthew G Knepley     PetscBLASInt lierr;
789566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
79*792fecdfSBarry Smith     PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr));
8063a3b9bcSJacob Faibussowitsch     PetscCheck(!lierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %" PetscBLASInt_FMT,lierr);
819566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
823f83f0d9SMatthew G Knepley   }
8327c67122SBarry Smith #else
84ef47b4b1SBarry Smith   {
85ef47b4b1SBarry Smith     PetscBLASInt lierr;
86ef47b4b1SBarry Smith     PetscReal    *rwork,*dd;
879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(5*nb,&rwork));
889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nb,&dd));
899566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
90*792fecdfSBarry Smith     PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,dd,u,&nb,v,&nb,work,&lwork,rwork,&lierr));
9163a3b9bcSJacob Faibussowitsch     PetscCheck(!lierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %" PetscBLASInt_FMT,lierr);
929566063dSJacob Faibussowitsch     PetscCall(PetscFree(rwork));
9314ce09a1SBarry Smith     for (i=0; i<n; i++) d[i] = dd[i];
949566063dSJacob Faibussowitsch     PetscCall(PetscFree(dd));
959566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
96ef47b4b1SBarry Smith   }
9727c67122SBarry Smith #endif
989566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->A,&a));
999566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->U,&u));
1009566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(jac->Vt,&v));
101426160bdSJed Brown   for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break;
102426160bdSJed Brown   jac->nzero = n-1-i;
103426160bdSJed Brown   if (jac->monitor) {
1049566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel));
10563a3b9bcSJacob 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));
1067962402dSFande Kong     if (n < 10 || jac->monitorformat == PETSC_VIEWER_ALL) {
1077962402dSFande Kong       PetscCall(PetscViewerASCIIPrintf(jac->monitor,"    SVD: singular values:\n"));
1087962402dSFande Kong       for (i=0; i<n; i++) {
1097962402dSFande Kong         if (i%5 == 0) {
1107962402dSFande Kong             if (i != 0) {
1117962402dSFande Kong               PetscCall(PetscViewerASCIIPrintf(jac->monitor,"\n"));
1127962402dSFande Kong             }
1137962402dSFande Kong             PetscCall(PetscViewerASCIIPrintf(jac->monitor,"        "));
1147962402dSFande Kong           }
1157962402dSFande Kong         PetscCall(PetscViewerASCIIPrintf(jac->monitor," %14.12e",(double)PetscRealPart(d[i])));
1167962402dSFande Kong       }
1177962402dSFande Kong       PetscCall(PetscViewerASCIIPrintf(jac->monitor,"\n"));
1187962402dSFande Kong     } else {              /* print 5 smallest and 5 largest */
1199566063dSJacob 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])));
1209566063dSJacob 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])));
121426160bdSJed Brown     }
1229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel));
123426160bdSJed Brown   }
1249566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1])));
125426160bdSJed Brown   for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i];
126426160bdSJed Brown   for (; i<n; i++) d[i] = 0.0;
12785032590SJed Brown   if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */
12863a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"Number of zero or nearly singular values %" PetscInt_FMT "\n",jac->nzero));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(jac->diag,&d));
1309566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
13127c67122SBarry Smith   PetscFunctionReturn(0);
13227c67122SBarry Smith }
13327c67122SBarry Smith 
1343ed27f31SJed Brown static PetscErrorCode PCSVDGetVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred)
1353ed27f31SJed Brown {
1363ed27f31SJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
1373ed27f31SJed Brown   PetscMPIInt    size;
1383ed27f31SJed Brown 
1393ed27f31SJed Brown   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size));
1410298fd71SBarry Smith   *xred = NULL;
1423ed27f31SJed Brown   switch (side) {
1433ed27f31SJed Brown   case PC_LEFT:
1443ed27f31SJed Brown     if (size == 1) *xred = x;
1453ed27f31SJed Brown     else {
1469566063dSJacob Faibussowitsch       if (!jac->left2red) PetscCall(VecScatterCreateToAll(x,&jac->left2red,&jac->leftred));
1473ed27f31SJed Brown       if (amode & READ) {
1489566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD));
1499566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD));
1503ed27f31SJed Brown       }
1513ed27f31SJed Brown       *xred = jac->leftred;
1523ed27f31SJed Brown     }
1533ed27f31SJed Brown     break;
1543ed27f31SJed Brown   case PC_RIGHT:
1553ed27f31SJed Brown     if (size == 1) *xred = x;
1563ed27f31SJed Brown     else {
1579566063dSJacob Faibussowitsch       if (!jac->right2red) PetscCall(VecScatterCreateToAll(x,&jac->right2red,&jac->rightred));
1583ed27f31SJed Brown       if (amode & READ) {
1599566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD));
1609566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD));
1613ed27f31SJed Brown       }
1623ed27f31SJed Brown       *xred = jac->rightred;
1633ed27f31SJed Brown     }
1643ed27f31SJed Brown     break;
165ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT");
1663ed27f31SJed Brown   }
1673ed27f31SJed Brown   PetscFunctionReturn(0);
1683ed27f31SJed Brown }
1693ed27f31SJed Brown 
1703ed27f31SJed Brown static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred)
1713ed27f31SJed Brown {
1723ed27f31SJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
1733ed27f31SJed Brown   PetscMPIInt    size;
1743ed27f31SJed Brown 
1753ed27f31SJed Brown   PetscFunctionBegin;
1769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size));
1773ed27f31SJed Brown   switch (side) {
1783ed27f31SJed Brown   case PC_LEFT:
1793ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
1809566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE));
1819566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE));
1823ed27f31SJed Brown     }
1833ed27f31SJed Brown     break;
1843ed27f31SJed Brown   case PC_RIGHT:
1853ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
1869566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE));
1879566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE));
1883ed27f31SJed Brown     }
1893ed27f31SJed Brown     break;
190ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT");
1913ed27f31SJed Brown   }
1920298fd71SBarry Smith   *xred = NULL;
1933ed27f31SJed Brown   PetscFunctionReturn(0);
1943ed27f31SJed Brown }
1953ed27f31SJed Brown 
19627c67122SBarry Smith /* -------------------------------------------------------------------------- */
19727c67122SBarry Smith /*
19827c67122SBarry Smith    PCApply_SVD - Applies the SVD preconditioner to a vector.
19927c67122SBarry Smith 
20027c67122SBarry Smith    Input Parameters:
20127c67122SBarry Smith .  pc - the preconditioner context
20227c67122SBarry Smith .  x - input vector
20327c67122SBarry Smith 
20427c67122SBarry Smith    Output Parameter:
20527c67122SBarry Smith .  y - output vector
20627c67122SBarry Smith 
20727c67122SBarry Smith    Application Interface Routine: PCApply()
20827c67122SBarry Smith  */
20927c67122SBarry Smith static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y)
21027c67122SBarry Smith {
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 
233024db992SStefano Zampini static PetscErrorCode PCMatApply_SVD(PC pc,Mat X,Mat Y)
234024db992SStefano Zampini {
235024db992SStefano Zampini   PC_SVD *jac = (PC_SVD*)pc->data;
236024db992SStefano Zampini   Mat    W;
237024db992SStefano Zampini 
238024db992SStefano Zampini   PetscFunctionBegin;
239024db992SStefano Zampini   PetscCall(MatTransposeMatMult(jac->U,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&W));
240024db992SStefano Zampini   PetscCall(MatDiagonalScale(W,jac->diag,NULL));
241024db992SStefano Zampini   PetscCall(MatTransposeMatMult(jac->Vt,W,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
242024db992SStefano Zampini   PetscCall(MatDestroy(&W));
243024db992SStefano Zampini   PetscFunctionReturn(0);
244024db992SStefano Zampini }
245024db992SStefano Zampini 
2463ed27f31SJed Brown static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y)
2473ed27f31SJed Brown {
2483ed27f31SJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
2493ed27f31SJed Brown   Vec            work = jac->work,xred,yred;
2503ed27f31SJed Brown 
2513ed27f31SJed Brown   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc,PC_LEFT,READ,x,&xred));
2539566063dSJacob Faibussowitsch   PetscCall(PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred));
2549566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->Vt,xred,work));
2559566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(work,work,jac->diag));
2569566063dSJacob Faibussowitsch   PetscCall(MatMult(jac->U,work,yred));
2579566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred));
2589566063dSJacob Faibussowitsch   PetscCall(PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred));
25927c67122SBarry Smith   PetscFunctionReturn(0);
26027c67122SBarry Smith }
26127c67122SBarry Smith 
262a2d70de2SBarry Smith static PetscErrorCode PCReset_SVD(PC pc)
263a2d70de2SBarry Smith {
264a2d70de2SBarry Smith   PC_SVD         *jac = (PC_SVD*)pc->data;
265a2d70de2SBarry Smith 
266a2d70de2SBarry Smith   PetscFunctionBegin;
2679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->A));
2689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->U));
2699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->Vt));
2709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
2719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->work));
2729566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&jac->right2red));
2739566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&jac->left2red));
2749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->rightred));
2759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->leftred));
276a2d70de2SBarry Smith   PetscFunctionReturn(0);
277a2d70de2SBarry Smith }
278a2d70de2SBarry Smith 
27927c67122SBarry Smith /* -------------------------------------------------------------------------- */
28027c67122SBarry Smith /*
28127c67122SBarry Smith    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
28227c67122SBarry Smith    that was created with PCCreate_SVD().
28327c67122SBarry Smith 
28427c67122SBarry Smith    Input Parameter:
28527c67122SBarry Smith .  pc - the preconditioner context
28627c67122SBarry Smith 
28727c67122SBarry Smith    Application Interface Routine: PCDestroy()
28827c67122SBarry Smith */
28927c67122SBarry Smith static PetscErrorCode PCDestroy_SVD(PC pc)
29027c67122SBarry Smith {
291426160bdSJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
29227c67122SBarry Smith 
29327c67122SBarry Smith   PetscFunctionBegin;
2949566063dSJacob Faibussowitsch   PetscCall(PCReset_SVD(pc));
2959566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&jac->monitor));
2969566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
29727c67122SBarry Smith   PetscFunctionReturn(0);
29827c67122SBarry Smith }
29927c67122SBarry Smith 
3004416b707SBarry Smith static PetscErrorCode PCSetFromOptions_SVD(PetscOptionItems *PetscOptionsObject,PC pc)
30127c67122SBarry Smith {
3028f1a2a5eSBarry Smith   PC_SVD         *jac = (PC_SVD*)pc->data;
3037962402dSFande Kong   PetscBool      flg;
30427c67122SBarry Smith 
30527c67122SBarry Smith   PetscFunctionBegin;
306d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SVD options");
3079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,NULL));
3089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,NULL));
3097962402dSFande Kong   PetscCall(PetscOptionsViewer("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",&jac->monitor,&jac->monitorformat,&flg));
310d0609cedSBarry Smith   PetscOptionsHeadEnd();
31127c67122SBarry Smith   PetscFunctionReturn(0);
31227c67122SBarry Smith }
31327c67122SBarry Smith 
31433761216SBarry Smith static PetscErrorCode PCView_SVD(PC pc,PetscViewer viewer)
31533761216SBarry Smith {
31633761216SBarry Smith   PC_SVD         *svd = (PC_SVD*)pc->data;
31733761216SBarry Smith   PetscBool      iascii;
31833761216SBarry Smith 
31933761216SBarry Smith   PetscFunctionBegin;
3209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
32133761216SBarry Smith   if (iascii) {
3229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  All singular values smaller than %g treated as zero\n",(double)svd->zerosing));
32363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Provided essential rank of the matrix %" PetscInt_FMT " (all other eigenvalues are zeroed)\n",svd->essrank));
32433761216SBarry Smith   }
32533761216SBarry Smith   PetscFunctionReturn(0);
32633761216SBarry Smith }
32727c67122SBarry Smith /* -------------------------------------------------------------------------- */
32827c67122SBarry Smith /*
32927c67122SBarry Smith    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
33027c67122SBarry Smith    and sets this as the private data within the generic preconditioning
33127c67122SBarry Smith    context, PC, that was created within PCCreate().
33227c67122SBarry Smith 
33327c67122SBarry Smith    Input Parameter:
33427c67122SBarry Smith .  pc - the preconditioner context
33527c67122SBarry Smith 
33627c67122SBarry Smith    Application Interface Routine: PCCreate()
33727c67122SBarry Smith */
33827c67122SBarry Smith 
33927c67122SBarry Smith /*MC
34027c67122SBarry Smith      PCSVD - Use pseudo inverse defined by SVD of operator
34127c67122SBarry Smith 
34227c67122SBarry Smith    Level: advanced
34327c67122SBarry Smith 
3448f1a2a5eSBarry Smith   Options Database:
345147403d9SBarry Smith +  -pc_svd_zero_sing <rtol> - Singular values smaller than this are treated as zero
346147403d9SBarry Smith -  -pc_svd_monitor - Print information on the extreme singular values of the operator
34727c67122SBarry Smith 
348a2b725a8SWilliam Gropp   Developer Note:
349a2b725a8SWilliam Gropp   This implementation automatically creates a redundant copy of the
3508997ae2eSBarry Smith    matrix on each process and uses a sequential SVD solve. Why does it do this instead
3518997ae2eSBarry Smith    of using the composable PCREDUNDANT object?
3528997ae2eSBarry Smith 
353db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`
35427c67122SBarry Smith M*/
35527c67122SBarry Smith 
3568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc)
35727c67122SBarry Smith {
35827c67122SBarry Smith   PC_SVD      *jac;
359024db992SStefano Zampini   PetscMPIInt size = 0;
36027c67122SBarry Smith 
36127c67122SBarry Smith   PetscFunctionBegin;
36227c67122SBarry Smith   /*
36327c67122SBarry Smith      Creates the private data structure for this preconditioner and
36427c67122SBarry Smith      attach it to the PC object.
36527c67122SBarry Smith   */
3669566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&jac));
3678f1a2a5eSBarry Smith   jac->zerosing = 1.e-12;
36827c67122SBarry Smith   pc->data      = (void*)jac;
36927c67122SBarry Smith 
37027c67122SBarry Smith   /*
37127c67122SBarry Smith       Set the pointers for the functions that are provided above.
37227c67122SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
37327c67122SBarry Smith       are called, they will automatically call these functions.  Note we
37427c67122SBarry Smith       choose not to provide a couple of these functions since they are
37527c67122SBarry Smith       not needed.
37627c67122SBarry Smith   */
377024db992SStefano Zampini 
378024db992SStefano Zampini #if defined(PETSC_HAVE_COMPLEX)
379024db992SStefano Zampini   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size));
380024db992SStefano Zampini #endif
381024db992SStefano Zampini   if (size == 1) pc->ops->matapply = PCMatApply_SVD;
38227c67122SBarry Smith   pc->ops->apply           = PCApply_SVD;
3833ed27f31SJed Brown   pc->ops->applytranspose  = PCApplyTranspose_SVD;
38427c67122SBarry Smith   pc->ops->setup           = PCSetUp_SVD;
385a2d70de2SBarry Smith   pc->ops->reset           = PCReset_SVD;
38627c67122SBarry Smith   pc->ops->destroy         = PCDestroy_SVD;
38727c67122SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SVD;
38833761216SBarry Smith   pc->ops->view            = PCView_SVD;
3890a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
39027c67122SBarry Smith   PetscFunctionReturn(0);
39127c67122SBarry Smith }
392