xref: /petsc/src/ksp/pc/impls/svd/svd.c (revision 14ce09a132c060283b0f1bcfa84acdede6679f6d)
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;
1727c67122SBarry Smith } PC_SVD;
1827c67122SBarry Smith 
193ed27f31SJed Brown typedef enum {READ=1, WRITE=2, READ_WRITE=3} AccessMode;
2027c67122SBarry Smith 
2127c67122SBarry Smith /* -------------------------------------------------------------------------- */
2227c67122SBarry Smith /*
2327c67122SBarry Smith    PCSetUp_SVD - Prepares for the use of the SVD preconditioner
2427c67122SBarry Smith                     by setting data structures and options.
2527c67122SBarry Smith 
2627c67122SBarry Smith    Input Parameter:
2727c67122SBarry Smith .  pc - the preconditioner context
2827c67122SBarry Smith 
2927c67122SBarry Smith    Application Interface Routine: PCSetUp()
3027c67122SBarry Smith 
3127c67122SBarry Smith    Notes:
3227c67122SBarry Smith    The interface routine PCSetUp() is not usually called directly by
3327c67122SBarry Smith    the user, but instead is called by PCApply() if necessary.
3427c67122SBarry Smith */
3527c67122SBarry Smith #undef __FUNCT__
3627c67122SBarry Smith #define __FUNCT__ "PCSetUp_SVD"
3727c67122SBarry Smith static PetscErrorCode PCSetUp_SVD(PC pc)
3827c67122SBarry Smith {
39a3c4e3ecSJed Brown #if defined(PETSC_MISSING_LAPACK_GESVD)
40ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates.");
41a3c4e3ecSJed Brown #else
4227c67122SBarry Smith   PC_SVD         *jac = (PC_SVD*)pc->data;
4327c67122SBarry Smith   PetscErrorCode ierr;
4427c67122SBarry Smith   PetscScalar    *a,*u,*v,*d,*work;
453f83f0d9SMatthew G Knepley   PetscBLASInt   nb,lwork;
4627c67122SBarry Smith   PetscInt       i,n;
473ed27f31SJed Brown   PetscMPIInt    size;
4827c67122SBarry Smith 
4927c67122SBarry Smith   PetscFunctionBegin;
506bf464f9SBarry Smith   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
513ed27f31SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc->pmat)->comm,&size);CHKERRQ(ierr);
523ed27f31SJed Brown   if (size > 1) {
533ed27f31SJed Brown     Mat          redmat;
5453cd1579SHong Zhang     ierr = MatCreateRedundantMatrix(pc->pmat,size,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&redmat);CHKERRQ(ierr);
553ed27f31SJed Brown     ierr = MatConvert(redmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
563ed27f31SJed Brown     ierr = MatDestroy(&redmat);CHKERRQ(ierr);
573ed27f31SJed Brown   } else {
588f1a2a5eSBarry Smith     ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr);
593ed27f31SJed Brown   }
603ed27f31SJed Brown   if (!jac->diag) {    /* assume square matrices */
612a7a6963SBarry Smith     ierr = MatCreateVecs(jac->A,&jac->diag,&jac->work);CHKERRQ(ierr);
623ed27f31SJed Brown   }
6327c67122SBarry Smith   if (!jac->U) {
6427c67122SBarry Smith     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr);
653ed27f31SJed Brown     ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->Vt);CHKERRQ(ierr);
6627c67122SBarry Smith   }
67*14ce09a1SBarry Smith   ierr  = MatGetSize(jac->A,&n,NULL);CHKERRQ(ierr);
68*14ce09a1SBarry Smith   if (!n) {
69*14ce09a1SBarry Smith     ierr = PetscInfo(pc,"Matrix has zero rows, skipping svd\n");
70*14ce09a1SBarry Smith     PetscFunctionReturn(0);
71*14ce09a1SBarry Smith   }
72c5df96a5SBarry Smith   ierr  = PetscBLASIntCast(n,&nb);CHKERRQ(ierr);
7327c67122SBarry Smith   lwork = 5*nb;
74785e854fSJed Brown   ierr  = PetscMalloc1(lwork,&work);CHKERRQ(ierr);
758c778c55SBarry Smith   ierr  = MatDenseGetArray(jac->A,&a);CHKERRQ(ierr);
768c778c55SBarry Smith   ierr  = MatDenseGetArray(jac->U,&u);CHKERRQ(ierr);
778c778c55SBarry Smith   ierr  = MatDenseGetArray(jac->Vt,&v);CHKERRQ(ierr);
7827c67122SBarry Smith   ierr  = VecGetArray(jac->diag,&d);CHKERRQ(ierr);
7927c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX)
803f83f0d9SMatthew G Knepley   {
813f83f0d9SMatthew G Knepley     PetscBLASInt lierr;
82670f3ff9SJed Brown     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
838b83055fSJed Brown     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr));
8498b909ebSSatish Balay     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
85670f3ff9SJed Brown     ierr = PetscFPTrapPop();CHKERRQ(ierr);
863f83f0d9SMatthew G Knepley   }
8727c67122SBarry Smith #else
88ef47b4b1SBarry Smith   {
89ef47b4b1SBarry Smith     PetscBLASInt lierr;
90ef47b4b1SBarry Smith     PetscReal    *rwork,*dd;
91785e854fSJed Brown     ierr = PetscMalloc1(5*nb,&rwork);CHKERRQ(ierr);
92785e854fSJed Brown     ierr = PetscMalloc1(nb,&dd);CHKERRQ(ierr);
93ef47b4b1SBarry Smith     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
94ef47b4b1SBarry Smith     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,dd,u,&nb,v,&nb,work,&lwork,rwork,&lierr));
95ef47b4b1SBarry Smith     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr);
96ef47b4b1SBarry Smith     ierr = PetscFree(rwork);CHKERRQ(ierr);
97*14ce09a1SBarry Smith     for (i=0; i<n; i++) d[i] = dd[i];
98ef47b4b1SBarry Smith     ierr = PetscFree(dd);CHKERRQ(ierr);
99ef47b4b1SBarry Smith     ierr = PetscFPTrapPop();CHKERRQ(ierr);
100ef47b4b1SBarry Smith   }
10127c67122SBarry Smith #endif
1028c778c55SBarry Smith   ierr = MatDenseRestoreArray(jac->A,&a);CHKERRQ(ierr);
1038c778c55SBarry Smith   ierr = MatDenseRestoreArray(jac->U,&u);CHKERRQ(ierr);
1048c778c55SBarry Smith   ierr = MatDenseRestoreArray(jac->Vt,&v);CHKERRQ(ierr);
105426160bdSJed Brown   for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break;
106426160bdSJed Brown   jac->nzero = n-1-i;
107426160bdSJed Brown   if (jac->monitor) {
108426160bdSJed Brown     ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
109426160bdSJed Brown     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);
110426160bdSJed Brown     if (n >= 10) {              /* print 5 smallest and 5 largest */
111426160bdSJed Brown       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);
112426160bdSJed Brown       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);
113426160bdSJed Brown     } else {                    /* print all singular values */
114426160bdSJed Brown       char     buf[256],*p;
1158caf3d72SBarry Smith       size_t   left = sizeof(buf),used;
116426160bdSJed Brown       PetscInt thisline;
117426160bdSJed Brown       for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) {
118426160bdSJed Brown         ierr  = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr);
119426160bdSJed Brown         left -= used;
120426160bdSJed Brown         p    += used;
121426160bdSJed Brown         if (thisline > 4 || i==0) {
122426160bdSJed Brown           ierr     = PetscViewerASCIIPrintf(jac->monitor,"    SVD: singular values:%s\n",buf);CHKERRQ(ierr);
123426160bdSJed Brown           p        = buf;
124426160bdSJed Brown           thisline = 0;
12527c67122SBarry Smith         }
126426160bdSJed Brown       }
127426160bdSJed Brown     }
128426160bdSJed Brown     ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr);
129426160bdSJed Brown   }
13022d28d08SBarry Smith   ierr = PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1]));CHKERRQ(ierr);
131426160bdSJed Brown   for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i];
132426160bdSJed Brown   for (; i<n; i++) d[i] = 0.0;
13385032590SJed Brown   if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */
13422d28d08SBarry Smith   ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);CHKERRQ(ierr);
13527c67122SBarry Smith   ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr);
13627c67122SBarry Smith #if defined(foo)
13727c67122SBarry Smith   {
13827c67122SBarry Smith     PetscViewer viewer;
13927c67122SBarry Smith     ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
14027c67122SBarry Smith     ierr = MatView(jac->A,viewer);CHKERRQ(ierr);
14127c67122SBarry Smith     ierr = MatView(jac->U,viewer);CHKERRQ(ierr);
1423ed27f31SJed Brown     ierr = MatView(jac->Vt,viewer);CHKERRQ(ierr);
14327c67122SBarry Smith     ierr = VecView(jac->diag,viewer);CHKERRQ(ierr);
14427c67122SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
14527c67122SBarry Smith   }
14627c67122SBarry Smith #endif
14722d28d08SBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
14827c67122SBarry Smith   PetscFunctionReturn(0);
149a3c4e3ecSJed Brown #endif
15027c67122SBarry Smith }
15127c67122SBarry Smith 
1523ed27f31SJed Brown #undef __FUNCT__
1533ed27f31SJed Brown #define __FUNCT__ "PCSVDGetVec"
1543ed27f31SJed Brown static PetscErrorCode PCSVDGetVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred)
1553ed27f31SJed Brown {
1563ed27f31SJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
1573ed27f31SJed Brown   PetscErrorCode ierr;
1583ed27f31SJed Brown   PetscMPIInt    size;
1593ed27f31SJed Brown 
1603ed27f31SJed Brown   PetscFunctionBegin;
161ce94432eSBarry Smith   ierr  = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
1620298fd71SBarry Smith   *xred = NULL;
1633ed27f31SJed Brown   switch (side) {
1643ed27f31SJed Brown   case PC_LEFT:
1653ed27f31SJed Brown     if (size == 1) *xred = x;
1663ed27f31SJed Brown     else {
1673ed27f31SJed Brown       if (!jac->left2red) {ierr = VecScatterCreateToAll(x,&jac->left2red,&jac->leftred);CHKERRQ(ierr);}
1683ed27f31SJed Brown       if (amode & READ) {
1693ed27f31SJed Brown         ierr = VecScatterBegin(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1703ed27f31SJed Brown         ierr = VecScatterEnd(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1713ed27f31SJed Brown       }
1723ed27f31SJed Brown       *xred = jac->leftred;
1733ed27f31SJed Brown     }
1743ed27f31SJed Brown     break;
1753ed27f31SJed Brown   case PC_RIGHT:
1763ed27f31SJed Brown     if (size == 1) *xred = x;
1773ed27f31SJed Brown     else {
1783ed27f31SJed Brown       if (!jac->right2red) {ierr = VecScatterCreateToAll(x,&jac->right2red,&jac->rightred);CHKERRQ(ierr);}
1793ed27f31SJed Brown       if (amode & READ) {
1803ed27f31SJed Brown         ierr = VecScatterBegin(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1813ed27f31SJed Brown         ierr = VecScatterEnd(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1823ed27f31SJed Brown       }
1833ed27f31SJed Brown       *xred = jac->rightred;
1843ed27f31SJed Brown     }
1853ed27f31SJed Brown     break;
186ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT");
1873ed27f31SJed Brown   }
1883ed27f31SJed Brown   PetscFunctionReturn(0);
1893ed27f31SJed Brown }
1903ed27f31SJed Brown 
1913ed27f31SJed Brown #undef __FUNCT__
1923ed27f31SJed Brown #define __FUNCT__ "PCSVDRestoreVec"
1933ed27f31SJed Brown static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred)
1943ed27f31SJed Brown {
1953ed27f31SJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
1963ed27f31SJed Brown   PetscErrorCode ierr;
1973ed27f31SJed Brown   PetscMPIInt    size;
1983ed27f31SJed Brown 
1993ed27f31SJed Brown   PetscFunctionBegin;
200ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
2013ed27f31SJed Brown   switch (side) {
2023ed27f31SJed Brown   case PC_LEFT:
2033ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
2043ed27f31SJed Brown       ierr = VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2053ed27f31SJed Brown       ierr = VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2063ed27f31SJed Brown     }
2073ed27f31SJed Brown     break;
2083ed27f31SJed Brown   case PC_RIGHT:
2093ed27f31SJed Brown     if (size != 1 && amode & WRITE) {
2103ed27f31SJed Brown       ierr = VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2113ed27f31SJed Brown       ierr = VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2123ed27f31SJed Brown     }
2133ed27f31SJed Brown     break;
214ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT");
2153ed27f31SJed Brown   }
2160298fd71SBarry Smith   *xred = NULL;
2173ed27f31SJed Brown   PetscFunctionReturn(0);
2183ed27f31SJed Brown }
2193ed27f31SJed Brown 
22027c67122SBarry Smith /* -------------------------------------------------------------------------- */
22127c67122SBarry Smith /*
22227c67122SBarry Smith    PCApply_SVD - Applies the SVD preconditioner to a vector.
22327c67122SBarry Smith 
22427c67122SBarry Smith    Input Parameters:
22527c67122SBarry Smith .  pc - the preconditioner context
22627c67122SBarry Smith .  x - input vector
22727c67122SBarry Smith 
22827c67122SBarry Smith    Output Parameter:
22927c67122SBarry Smith .  y - output vector
23027c67122SBarry Smith 
23127c67122SBarry Smith    Application Interface Routine: PCApply()
23227c67122SBarry Smith  */
23327c67122SBarry Smith #undef __FUNCT__
23427c67122SBarry Smith #define __FUNCT__ "PCApply_SVD"
23527c67122SBarry Smith static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y)
23627c67122SBarry Smith {
23727c67122SBarry Smith   PC_SVD         *jac = (PC_SVD*)pc->data;
2383ed27f31SJed Brown   Vec            work = jac->work,xred,yred;
23927c67122SBarry Smith   PetscErrorCode ierr;
24027c67122SBarry Smith 
24127c67122SBarry Smith   PetscFunctionBegin;
2423ed27f31SJed Brown   ierr = PCSVDGetVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr);
2433ed27f31SJed Brown   ierr = PCSVDGetVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr);
244ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2453ed27f31SJed Brown   ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr);
246ef47b4b1SBarry Smith #else
247ef47b4b1SBarry Smith   ierr = MatMultHermitianTranspose(jac->U,xred,work);CHKERRQ(ierr);
248ef47b4b1SBarry Smith #endif
24927c67122SBarry Smith   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
250ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2513ed27f31SJed Brown   ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
252ef47b4b1SBarry Smith #else
253ef47b4b1SBarry Smith   ierr = MatMultHermitianTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
254ef47b4b1SBarry Smith #endif
2553ed27f31SJed Brown   ierr = PCSVDRestoreVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr);
2563ed27f31SJed Brown   ierr = PCSVDRestoreVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr);
2573ed27f31SJed Brown   PetscFunctionReturn(0);
2583ed27f31SJed Brown }
2593ed27f31SJed Brown 
2603ed27f31SJed Brown #undef __FUNCT__
2613ed27f31SJed Brown #define __FUNCT__ "PCApplyTranspose_SVD"
2623ed27f31SJed Brown static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y)
2633ed27f31SJed Brown {
2643ed27f31SJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
2653ed27f31SJed Brown   Vec            work = jac->work,xred,yred;
2663ed27f31SJed Brown   PetscErrorCode ierr;
2673ed27f31SJed Brown 
2683ed27f31SJed Brown   PetscFunctionBegin;
2693ed27f31SJed Brown   ierr = PCSVDGetVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr);
2703ed27f31SJed Brown   ierr = PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr);
271ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2723ed27f31SJed Brown   ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
273ef47b4b1SBarry Smith #else
274ef47b4b1SBarry Smith   ierr = MatMultHermitianTranspose(jac->Vt,work,yred);CHKERRQ(ierr);
275ef47b4b1SBarry Smith #endif
2763ed27f31SJed Brown   ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr);
277ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
2783ed27f31SJed Brown   ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr);
279ef47b4b1SBarry Smith #else
280ef47b4b1SBarry Smith   ierr = MatMultHermitianTranspose(jac->U,xred,work);CHKERRQ(ierr);
281ef47b4b1SBarry Smith #endif
2823ed27f31SJed Brown   ierr = PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr);
2833ed27f31SJed Brown   ierr = PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr);
28427c67122SBarry Smith   PetscFunctionReturn(0);
28527c67122SBarry Smith }
28627c67122SBarry Smith 
287a2d70de2SBarry Smith #undef __FUNCT__
288a2d70de2SBarry Smith #define __FUNCT__ "PCReset_SVD"
289a2d70de2SBarry Smith static PetscErrorCode PCReset_SVD(PC pc)
290a2d70de2SBarry Smith {
291a2d70de2SBarry Smith   PC_SVD         *jac = (PC_SVD*)pc->data;
292a2d70de2SBarry Smith   PetscErrorCode ierr;
293a2d70de2SBarry Smith 
294a2d70de2SBarry Smith   PetscFunctionBegin;
295a2d70de2SBarry Smith   ierr = MatDestroy(&jac->A);CHKERRQ(ierr);
296a2d70de2SBarry Smith   ierr = MatDestroy(&jac->U);CHKERRQ(ierr);
2973ed27f31SJed Brown   ierr = MatDestroy(&jac->Vt);CHKERRQ(ierr);
298a2d70de2SBarry Smith   ierr = VecDestroy(&jac->diag);CHKERRQ(ierr);
299a2d70de2SBarry Smith   ierr = VecDestroy(&jac->work);CHKERRQ(ierr);
3003ed27f31SJed Brown   ierr = VecScatterDestroy(&jac->right2red);CHKERRQ(ierr);
3013ed27f31SJed Brown   ierr = VecScatterDestroy(&jac->left2red);CHKERRQ(ierr);
3023ed27f31SJed Brown   ierr = VecDestroy(&jac->rightred);CHKERRQ(ierr);
3033ed27f31SJed Brown   ierr = VecDestroy(&jac->leftred);CHKERRQ(ierr);
304a2d70de2SBarry Smith   PetscFunctionReturn(0);
305a2d70de2SBarry Smith }
306a2d70de2SBarry Smith 
30727c67122SBarry Smith /* -------------------------------------------------------------------------- */
30827c67122SBarry Smith /*
30927c67122SBarry Smith    PCDestroy_SVD - Destroys the private context for the SVD preconditioner
31027c67122SBarry Smith    that was created with PCCreate_SVD().
31127c67122SBarry Smith 
31227c67122SBarry Smith    Input Parameter:
31327c67122SBarry Smith .  pc - the preconditioner context
31427c67122SBarry Smith 
31527c67122SBarry Smith    Application Interface Routine: PCDestroy()
31627c67122SBarry Smith */
31727c67122SBarry Smith #undef __FUNCT__
31827c67122SBarry Smith #define __FUNCT__ "PCDestroy_SVD"
31927c67122SBarry Smith static PetscErrorCode PCDestroy_SVD(PC pc)
32027c67122SBarry Smith {
321426160bdSJed Brown   PC_SVD         *jac = (PC_SVD*)pc->data;
32227c67122SBarry Smith   PetscErrorCode ierr;
32327c67122SBarry Smith 
32427c67122SBarry Smith   PetscFunctionBegin;
325a2d70de2SBarry Smith   ierr = PCReset_SVD(pc);CHKERRQ(ierr);
326426160bdSJed Brown   ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
327c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
32827c67122SBarry Smith   PetscFunctionReturn(0);
32927c67122SBarry Smith }
33027c67122SBarry Smith 
33127c67122SBarry Smith #undef __FUNCT__
33227c67122SBarry Smith #define __FUNCT__ "PCSetFromOptions_SVD"
3334416b707SBarry Smith static PetscErrorCode PCSetFromOptions_SVD(PetscOptionItems *PetscOptionsObject,PC pc)
33427c67122SBarry Smith {
33527c67122SBarry Smith   PetscErrorCode ierr;
3368f1a2a5eSBarry Smith   PC_SVD         *jac = (PC_SVD*)pc->data;
337426160bdSJed Brown   PetscBool      flg,set;
33827c67122SBarry Smith 
33927c67122SBarry Smith   PetscFunctionBegin;
340e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SVD options");CHKERRQ(ierr);
3410298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,NULL);CHKERRQ(ierr);
3420298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,NULL);CHKERRQ(ierr);
3436ba663aaSJed Brown   ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
344426160bdSJed Brown   if (set) {                    /* Should make PCSVDSetMonitor() */
345426160bdSJed Brown     if (flg && !jac->monitor) {
346ce94432eSBarry Smith       ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),"stdout",&jac->monitor);CHKERRQ(ierr);
347426160bdSJed Brown     } else if (!flg) {
348426160bdSJed Brown       ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr);
349426160bdSJed Brown     }
350426160bdSJed Brown   }
35127c67122SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
35227c67122SBarry Smith   PetscFunctionReturn(0);
35327c67122SBarry Smith }
35427c67122SBarry Smith 
35527c67122SBarry Smith /* -------------------------------------------------------------------------- */
35627c67122SBarry Smith /*
35727c67122SBarry Smith    PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD,
35827c67122SBarry Smith    and sets this as the private data within the generic preconditioning
35927c67122SBarry Smith    context, PC, that was created within PCCreate().
36027c67122SBarry Smith 
36127c67122SBarry Smith    Input Parameter:
36227c67122SBarry Smith .  pc - the preconditioner context
36327c67122SBarry Smith 
36427c67122SBarry Smith    Application Interface Routine: PCCreate()
36527c67122SBarry Smith */
36627c67122SBarry Smith 
36727c67122SBarry Smith /*MC
36827c67122SBarry Smith      PCSVD - Use pseudo inverse defined by SVD of operator
36927c67122SBarry Smith 
37027c67122SBarry Smith    Level: advanced
37127c67122SBarry Smith 
37227c67122SBarry Smith   Concepts: SVD
37327c67122SBarry Smith 
3748f1a2a5eSBarry Smith   Options Database:
3758b0a5094SBarry Smith -  -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero
3768b0a5094SBarry Smith +  -pc_svd_monitor  Print information on the extreme singular values of the operator
37727c67122SBarry Smith 
3788997ae2eSBarry Smith   Developer Note: This implementation automatically creates a redundant copy of the
3798997ae2eSBarry Smith    matrix on each process and uses a sequential SVD solve. Why does it do this instead
3808997ae2eSBarry Smith    of using the composable PCREDUNDANT object?
3818997ae2eSBarry Smith 
38227c67122SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
38327c67122SBarry Smith M*/
38427c67122SBarry Smith 
38527c67122SBarry Smith #undef __FUNCT__
38627c67122SBarry Smith #define __FUNCT__ "PCCreate_SVD"
3878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc)
38827c67122SBarry Smith {
38927c67122SBarry Smith   PC_SVD         *jac;
39027c67122SBarry Smith   PetscErrorCode ierr;
39127c67122SBarry Smith 
39227c67122SBarry Smith   PetscFunctionBegin;
39327c67122SBarry Smith   /*
39427c67122SBarry Smith      Creates the private data structure for this preconditioner and
39527c67122SBarry Smith      attach it to the PC object.
39627c67122SBarry Smith   */
397b00a9115SJed Brown   ierr          = PetscNewLog(pc,&jac);CHKERRQ(ierr);
3988f1a2a5eSBarry Smith   jac->zerosing = 1.e-12;
39927c67122SBarry Smith   pc->data      = (void*)jac;
40027c67122SBarry Smith 
40127c67122SBarry Smith   /*
40227c67122SBarry Smith       Set the pointers for the functions that are provided above.
40327c67122SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
40427c67122SBarry Smith       are called, they will automatically call these functions.  Note we
40527c67122SBarry Smith       choose not to provide a couple of these functions since they are
40627c67122SBarry Smith       not needed.
40727c67122SBarry Smith   */
40827c67122SBarry Smith   pc->ops->apply           = PCApply_SVD;
4093ed27f31SJed Brown   pc->ops->applytranspose  = PCApplyTranspose_SVD;
41027c67122SBarry Smith   pc->ops->setup           = PCSetUp_SVD;
411a2d70de2SBarry Smith   pc->ops->reset           = PCReset_SVD;
41227c67122SBarry Smith   pc->ops->destroy         = PCDestroy_SVD;
41327c67122SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SVD;
41427c67122SBarry Smith   pc->ops->view            = 0;
41527c67122SBarry Smith   pc->ops->applyrichardson = 0;
41627c67122SBarry Smith   PetscFunctionReturn(0);
41727c67122SBarry Smith }
41827c67122SBarry Smith 
419