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 static PetscErrorCode PCSetUp_SVD(PC pc) 3627c67122SBarry Smith { 3727c67122SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 3827c67122SBarry Smith PetscErrorCode ierr; 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; 456bf464f9SBarry Smith ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 46ffc4695bSBarry Smith ierr = MPI_Comm_size(((PetscObject)pc->pmat)->comm,&size);CHKERRMPI(ierr); 473ed27f31SJed Brown if (size > 1) { 483ed27f31SJed Brown Mat redmat; 4928d58a37SPierre Jolivet 5053cd1579SHong Zhang ierr = MatCreateRedundantMatrix(pc->pmat,size,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&redmat);CHKERRQ(ierr); 513ed27f31SJed Brown ierr = MatConvert(redmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr); 523ed27f31SJed Brown ierr = MatDestroy(&redmat);CHKERRQ(ierr); 533ed27f31SJed Brown } else { 548f1a2a5eSBarry Smith ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr); 553ed27f31SJed Brown } 563ed27f31SJed Brown if (!jac->diag) { /* assume square matrices */ 572a7a6963SBarry Smith ierr = MatCreateVecs(jac->A,&jac->diag,&jac->work);CHKERRQ(ierr); 583ed27f31SJed Brown } 5927c67122SBarry Smith if (!jac->U) { 6027c67122SBarry Smith ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr); 613ed27f31SJed Brown ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->Vt);CHKERRQ(ierr); 6227c67122SBarry Smith } 6314ce09a1SBarry Smith ierr = MatGetSize(jac->A,&n,NULL);CHKERRQ(ierr); 6414ce09a1SBarry Smith if (!n) { 65459726d8SSatish Balay ierr = PetscInfo(pc,"Matrix has zero rows, skipping svd\n");CHKERRQ(ierr); 6614ce09a1SBarry Smith PetscFunctionReturn(0); 6714ce09a1SBarry Smith } 68c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&nb);CHKERRQ(ierr); 6927c67122SBarry Smith lwork = 5*nb; 70785e854fSJed Brown ierr = PetscMalloc1(lwork,&work);CHKERRQ(ierr); 718c778c55SBarry Smith ierr = MatDenseGetArray(jac->A,&a);CHKERRQ(ierr); 728c778c55SBarry Smith ierr = MatDenseGetArray(jac->U,&u);CHKERRQ(ierr); 738c778c55SBarry Smith ierr = MatDenseGetArray(jac->Vt,&v);CHKERRQ(ierr); 7427c67122SBarry Smith ierr = VecGetArray(jac->diag,&d);CHKERRQ(ierr); 7527c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX) 763f83f0d9SMatthew G Knepley { 773f83f0d9SMatthew G Knepley PetscBLASInt lierr; 78670f3ff9SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 798b83055fSJed Brown PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr)); 8098921bdaSJacob Faibussowitsch if (lierr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 81670f3ff9SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 823f83f0d9SMatthew G Knepley } 8327c67122SBarry Smith #else 84ef47b4b1SBarry Smith { 85ef47b4b1SBarry Smith PetscBLASInt lierr; 86ef47b4b1SBarry Smith PetscReal *rwork,*dd; 87785e854fSJed Brown ierr = PetscMalloc1(5*nb,&rwork);CHKERRQ(ierr); 88785e854fSJed Brown ierr = PetscMalloc1(nb,&dd);CHKERRQ(ierr); 89ef47b4b1SBarry Smith ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 90ef47b4b1SBarry Smith PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,dd,u,&nb,v,&nb,work,&lwork,rwork,&lierr)); 9198921bdaSJacob Faibussowitsch if (lierr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 92ef47b4b1SBarry Smith ierr = PetscFree(rwork);CHKERRQ(ierr); 9314ce09a1SBarry Smith for (i=0; i<n; i++) d[i] = dd[i]; 94ef47b4b1SBarry Smith ierr = PetscFree(dd);CHKERRQ(ierr); 95ef47b4b1SBarry Smith ierr = PetscFPTrapPop();CHKERRQ(ierr); 96ef47b4b1SBarry Smith } 9727c67122SBarry Smith #endif 988c778c55SBarry Smith ierr = MatDenseRestoreArray(jac->A,&a);CHKERRQ(ierr); 998c778c55SBarry Smith ierr = MatDenseRestoreArray(jac->U,&u);CHKERRQ(ierr); 1008c778c55SBarry Smith ierr = MatDenseRestoreArray(jac->Vt,&v);CHKERRQ(ierr); 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) { 104426160bdSJed Brown ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 105426160bdSJed 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); 106426160bdSJed Brown if (n >= 10) { /* print 5 smallest and 5 largest */ 107426160bdSJed 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); 108426160bdSJed 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); 109426160bdSJed Brown } else { /* print all singular values */ 110426160bdSJed Brown char buf[256],*p; 1118caf3d72SBarry Smith size_t left = sizeof(buf),used; 112426160bdSJed Brown PetscInt thisline; 113426160bdSJed Brown for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) { 114426160bdSJed Brown ierr = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr); 115426160bdSJed Brown left -= used; 116426160bdSJed Brown p += used; 117426160bdSJed Brown if (thisline > 4 || i==0) { 118426160bdSJed Brown ierr = PetscViewerASCIIPrintf(jac->monitor," SVD: singular values:%s\n",buf);CHKERRQ(ierr); 119426160bdSJed Brown p = buf; 120426160bdSJed Brown thisline = 0; 12127c67122SBarry Smith } 122426160bdSJed Brown } 123426160bdSJed Brown } 124426160bdSJed Brown ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 125426160bdSJed Brown } 126*7d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1]));CHKERRQ(ierr); 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; 12985032590SJed Brown if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */ 130*7d3de750SJacob Faibussowitsch ierr = PetscInfo(pc,"Number of zero or nearly singular values %D\n",jac->nzero);CHKERRQ(ierr); 13127c67122SBarry Smith ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr); 13227c67122SBarry Smith #if defined(foo) 13327c67122SBarry Smith { 13427c67122SBarry Smith PetscViewer viewer; 13527c67122SBarry Smith ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 13627c67122SBarry Smith ierr = MatView(jac->A,viewer);CHKERRQ(ierr); 13727c67122SBarry Smith ierr = MatView(jac->U,viewer);CHKERRQ(ierr); 1383ed27f31SJed Brown ierr = MatView(jac->Vt,viewer);CHKERRQ(ierr); 13927c67122SBarry Smith ierr = VecView(jac->diag,viewer);CHKERRQ(ierr); 14027c67122SBarry Smith ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 14127c67122SBarry Smith } 14227c67122SBarry Smith #endif 14322d28d08SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 14427c67122SBarry Smith PetscFunctionReturn(0); 14527c67122SBarry Smith } 14627c67122SBarry Smith 1473ed27f31SJed Brown static PetscErrorCode PCSVDGetVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 1483ed27f31SJed Brown { 1493ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 1503ed27f31SJed Brown PetscErrorCode ierr; 1513ed27f31SJed Brown PetscMPIInt size; 1523ed27f31SJed Brown 1533ed27f31SJed Brown PetscFunctionBegin; 15455b25c41SPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRMPI(ierr); 1550298fd71SBarry Smith *xred = NULL; 1563ed27f31SJed Brown switch (side) { 1573ed27f31SJed Brown case PC_LEFT: 1583ed27f31SJed Brown if (size == 1) *xred = x; 1593ed27f31SJed Brown else { 1603ed27f31SJed Brown if (!jac->left2red) {ierr = VecScatterCreateToAll(x,&jac->left2red,&jac->leftred);CHKERRQ(ierr);} 1613ed27f31SJed Brown if (amode & READ) { 1623ed27f31SJed Brown ierr = VecScatterBegin(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1633ed27f31SJed Brown ierr = VecScatterEnd(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1643ed27f31SJed Brown } 1653ed27f31SJed Brown *xred = jac->leftred; 1663ed27f31SJed Brown } 1673ed27f31SJed Brown break; 1683ed27f31SJed Brown case PC_RIGHT: 1693ed27f31SJed Brown if (size == 1) *xred = x; 1703ed27f31SJed Brown else { 1713ed27f31SJed Brown if (!jac->right2red) {ierr = VecScatterCreateToAll(x,&jac->right2red,&jac->rightred);CHKERRQ(ierr);} 1723ed27f31SJed Brown if (amode & READ) { 1733ed27f31SJed Brown ierr = VecScatterBegin(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1743ed27f31SJed Brown ierr = VecScatterEnd(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1753ed27f31SJed Brown } 1763ed27f31SJed Brown *xred = jac->rightred; 1773ed27f31SJed Brown } 1783ed27f31SJed Brown break; 179ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 1803ed27f31SJed Brown } 1813ed27f31SJed Brown PetscFunctionReturn(0); 1823ed27f31SJed Brown } 1833ed27f31SJed Brown 1843ed27f31SJed Brown static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 1853ed27f31SJed Brown { 1863ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 1873ed27f31SJed Brown PetscErrorCode ierr; 1883ed27f31SJed Brown PetscMPIInt size; 1893ed27f31SJed Brown 1903ed27f31SJed Brown PetscFunctionBegin; 191ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRMPI(ierr); 1923ed27f31SJed Brown switch (side) { 1933ed27f31SJed Brown case PC_LEFT: 1943ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1953ed27f31SJed Brown ierr = VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1963ed27f31SJed Brown ierr = VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1973ed27f31SJed Brown } 1983ed27f31SJed Brown break; 1993ed27f31SJed Brown case PC_RIGHT: 2003ed27f31SJed Brown if (size != 1 && amode & WRITE) { 2013ed27f31SJed Brown ierr = VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2023ed27f31SJed Brown ierr = VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2033ed27f31SJed Brown } 2043ed27f31SJed Brown break; 205ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 2063ed27f31SJed Brown } 2070298fd71SBarry Smith *xred = NULL; 2083ed27f31SJed Brown PetscFunctionReturn(0); 2093ed27f31SJed Brown } 2103ed27f31SJed Brown 21127c67122SBarry Smith /* -------------------------------------------------------------------------- */ 21227c67122SBarry Smith /* 21327c67122SBarry Smith PCApply_SVD - Applies the SVD preconditioner to a vector. 21427c67122SBarry Smith 21527c67122SBarry Smith Input Parameters: 21627c67122SBarry Smith . pc - the preconditioner context 21727c67122SBarry Smith . x - input vector 21827c67122SBarry Smith 21927c67122SBarry Smith Output Parameter: 22027c67122SBarry Smith . y - output vector 22127c67122SBarry Smith 22227c67122SBarry Smith Application Interface Routine: PCApply() 22327c67122SBarry Smith */ 22427c67122SBarry Smith static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y) 22527c67122SBarry Smith { 22627c67122SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 2273ed27f31SJed Brown Vec work = jac->work,xred,yred; 22827c67122SBarry Smith PetscErrorCode ierr; 22927c67122SBarry Smith 23027c67122SBarry Smith PetscFunctionBegin; 2313ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr); 2323ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr); 233ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2343ed27f31SJed Brown ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr); 235ef47b4b1SBarry Smith #else 236ef47b4b1SBarry Smith ierr = MatMultHermitianTranspose(jac->U,xred,work);CHKERRQ(ierr); 237ef47b4b1SBarry Smith #endif 23827c67122SBarry Smith ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 239ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2403ed27f31SJed Brown ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr); 241ef47b4b1SBarry Smith #else 242ef47b4b1SBarry Smith ierr = MatMultHermitianTranspose(jac->Vt,work,yred);CHKERRQ(ierr); 243ef47b4b1SBarry Smith #endif 2443ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr); 2453ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr); 2463ed27f31SJed Brown PetscFunctionReturn(0); 2473ed27f31SJed Brown } 2483ed27f31SJed Brown 2493ed27f31SJed Brown static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y) 2503ed27f31SJed Brown { 2513ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 2523ed27f31SJed Brown Vec work = jac->work,xred,yred; 2533ed27f31SJed Brown PetscErrorCode ierr; 2543ed27f31SJed Brown 2553ed27f31SJed Brown PetscFunctionBegin; 2563ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr); 2573ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr); 258eaf392e5SStefano Zampini ierr = MatMult(jac->Vt,xred,work);CHKERRQ(ierr); 2593ed27f31SJed Brown ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 260eaf392e5SStefano Zampini ierr = MatMult(jac->U,work,yred);CHKERRQ(ierr); 2613ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr); 2623ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr); 26327c67122SBarry Smith PetscFunctionReturn(0); 26427c67122SBarry Smith } 26527c67122SBarry Smith 266a2d70de2SBarry Smith static PetscErrorCode PCReset_SVD(PC pc) 267a2d70de2SBarry Smith { 268a2d70de2SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 269a2d70de2SBarry Smith PetscErrorCode ierr; 270a2d70de2SBarry Smith 271a2d70de2SBarry Smith PetscFunctionBegin; 272a2d70de2SBarry Smith ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 273a2d70de2SBarry Smith ierr = MatDestroy(&jac->U);CHKERRQ(ierr); 2743ed27f31SJed Brown ierr = MatDestroy(&jac->Vt);CHKERRQ(ierr); 275a2d70de2SBarry Smith ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 276a2d70de2SBarry Smith ierr = VecDestroy(&jac->work);CHKERRQ(ierr); 2773ed27f31SJed Brown ierr = VecScatterDestroy(&jac->right2red);CHKERRQ(ierr); 2783ed27f31SJed Brown ierr = VecScatterDestroy(&jac->left2red);CHKERRQ(ierr); 2793ed27f31SJed Brown ierr = VecDestroy(&jac->rightred);CHKERRQ(ierr); 2803ed27f31SJed Brown ierr = VecDestroy(&jac->leftred);CHKERRQ(ierr); 281a2d70de2SBarry Smith PetscFunctionReturn(0); 282a2d70de2SBarry Smith } 283a2d70de2SBarry Smith 28427c67122SBarry Smith /* -------------------------------------------------------------------------- */ 28527c67122SBarry Smith /* 28627c67122SBarry Smith PCDestroy_SVD - Destroys the private context for the SVD preconditioner 28727c67122SBarry Smith that was created with PCCreate_SVD(). 28827c67122SBarry Smith 28927c67122SBarry Smith Input Parameter: 29027c67122SBarry Smith . pc - the preconditioner context 29127c67122SBarry Smith 29227c67122SBarry Smith Application Interface Routine: PCDestroy() 29327c67122SBarry Smith */ 29427c67122SBarry Smith static PetscErrorCode PCDestroy_SVD(PC pc) 29527c67122SBarry Smith { 296426160bdSJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 29727c67122SBarry Smith PetscErrorCode ierr; 29827c67122SBarry Smith 29927c67122SBarry Smith PetscFunctionBegin; 300a2d70de2SBarry Smith ierr = PCReset_SVD(pc);CHKERRQ(ierr); 301426160bdSJed Brown ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 302c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 30327c67122SBarry Smith PetscFunctionReturn(0); 30427c67122SBarry Smith } 30527c67122SBarry Smith 3064416b707SBarry Smith static PetscErrorCode PCSetFromOptions_SVD(PetscOptionItems *PetscOptionsObject,PC pc) 30727c67122SBarry Smith { 30827c67122SBarry Smith PetscErrorCode ierr; 3098f1a2a5eSBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 310426160bdSJed Brown PetscBool flg,set; 31127c67122SBarry Smith 31227c67122SBarry Smith PetscFunctionBegin; 313e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SVD options");CHKERRQ(ierr); 3140298fd71SBarry Smith ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,NULL);CHKERRQ(ierr); 3150298fd71SBarry Smith ierr = PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,NULL);CHKERRQ(ierr); 3166ba663aaSJed Brown ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 317426160bdSJed Brown if (set) { /* Should make PCSVDSetMonitor() */ 318426160bdSJed Brown if (flg && !jac->monitor) { 319ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),"stdout",&jac->monitor);CHKERRQ(ierr); 320426160bdSJed Brown } else if (!flg) { 321426160bdSJed Brown ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 322426160bdSJed Brown } 323426160bdSJed Brown } 32427c67122SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 32527c67122SBarry Smith PetscFunctionReturn(0); 32627c67122SBarry Smith } 32727c67122SBarry Smith 32833761216SBarry Smith static PetscErrorCode PCView_SVD(PC pc,PetscViewer viewer) 32933761216SBarry Smith { 33033761216SBarry Smith PC_SVD *svd = (PC_SVD*)pc->data; 33133761216SBarry Smith PetscErrorCode ierr; 33233761216SBarry Smith PetscBool iascii; 33333761216SBarry Smith 33433761216SBarry Smith PetscFunctionBegin; 33533761216SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 33633761216SBarry Smith if (iascii) { 337efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," All singular values smaller than %g treated as zero\n",(double)svd->zerosing);CHKERRQ(ierr); 338efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Provided essential rank of the matrix %D (all other eigenvalues are zeroed)\n",svd->essrank);CHKERRQ(ierr); 33933761216SBarry Smith } 34033761216SBarry Smith PetscFunctionReturn(0); 34133761216SBarry Smith } 34227c67122SBarry Smith /* -------------------------------------------------------------------------- */ 34327c67122SBarry Smith /* 34427c67122SBarry Smith PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 34527c67122SBarry Smith and sets this as the private data within the generic preconditioning 34627c67122SBarry Smith context, PC, that was created within PCCreate(). 34727c67122SBarry Smith 34827c67122SBarry Smith Input Parameter: 34927c67122SBarry Smith . pc - the preconditioner context 35027c67122SBarry Smith 35127c67122SBarry Smith Application Interface Routine: PCCreate() 35227c67122SBarry Smith */ 35327c67122SBarry Smith 35427c67122SBarry Smith /*MC 35527c67122SBarry Smith PCSVD - Use pseudo inverse defined by SVD of operator 35627c67122SBarry Smith 35727c67122SBarry Smith Level: advanced 35827c67122SBarry Smith 3598f1a2a5eSBarry Smith Options Database: 360a2b725a8SWilliam Gropp + -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero 361a2b725a8SWilliam Gropp - -pc_svd_monitor Print information on the extreme singular values of the operator 36227c67122SBarry Smith 363a2b725a8SWilliam Gropp Developer Note: 364a2b725a8SWilliam Gropp This implementation automatically creates a redundant copy of the 3658997ae2eSBarry Smith matrix on each process and uses a sequential SVD solve. Why does it do this instead 3668997ae2eSBarry Smith of using the composable PCREDUNDANT object? 3678997ae2eSBarry Smith 36827c67122SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 36927c67122SBarry Smith M*/ 37027c67122SBarry Smith 3718cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc) 37227c67122SBarry Smith { 37327c67122SBarry Smith PC_SVD *jac; 37427c67122SBarry Smith PetscErrorCode ierr; 37527c67122SBarry Smith 37627c67122SBarry Smith PetscFunctionBegin; 37727c67122SBarry Smith /* 37827c67122SBarry Smith Creates the private data structure for this preconditioner and 37927c67122SBarry Smith attach it to the PC object. 38027c67122SBarry Smith */ 381b00a9115SJed Brown ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 3828f1a2a5eSBarry Smith jac->zerosing = 1.e-12; 38327c67122SBarry Smith pc->data = (void*)jac; 38427c67122SBarry Smith 38527c67122SBarry Smith /* 38627c67122SBarry Smith Set the pointers for the functions that are provided above. 38727c67122SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 38827c67122SBarry Smith are called, they will automatically call these functions. Note we 38927c67122SBarry Smith choose not to provide a couple of these functions since they are 39027c67122SBarry Smith not needed. 39127c67122SBarry Smith */ 39227c67122SBarry Smith pc->ops->apply = PCApply_SVD; 3933ed27f31SJed Brown pc->ops->applytranspose = PCApplyTranspose_SVD; 39427c67122SBarry Smith pc->ops->setup = PCSetUp_SVD; 395a2d70de2SBarry Smith pc->ops->reset = PCReset_SVD; 39627c67122SBarry Smith pc->ops->destroy = PCDestroy_SVD; 39727c67122SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SVD; 39833761216SBarry Smith pc->ops->view = PCView_SVD; 3990a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 40027c67122SBarry Smith PetscFunctionReturn(0); 40127c67122SBarry Smith } 40227c67122SBarry Smith 403