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 PetscScalar *a,*u,*v,*d,*work; 393f83f0d9SMatthew G Knepley PetscBLASInt nb,lwork; 4027c67122SBarry Smith PetscInt i,n; 413ed27f31SJed Brown PetscMPIInt size; 4227c67122SBarry Smith 4327c67122SBarry Smith PetscFunctionBegin; 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&jac->A)); 45*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(((PetscObject)pc->pmat)->comm,&size)); 463ed27f31SJed Brown if (size > 1) { 473ed27f31SJed Brown Mat redmat; 4828d58a37SPierre Jolivet 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(pc->pmat,size,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&redmat)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(redmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&redmat)); 523ed27f31SJed Brown } else { 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A)); 543ed27f31SJed Brown } 553ed27f31SJed Brown if (!jac->diag) { /* assume square matrices */ 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(jac->A,&jac->diag,&jac->work)); 573ed27f31SJed Brown } 5827c67122SBarry Smith if (!jac->U) { 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->Vt)); 6127c67122SBarry Smith } 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(jac->A,&n,NULL)); 6314ce09a1SBarry Smith if (!n) { 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(pc,"Matrix has zero rows, skipping svd\n")); 6514ce09a1SBarry Smith PetscFunctionReturn(0); 6614ce09a1SBarry Smith } 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBLASIntCast(n,&nb)); 6827c67122SBarry Smith lwork = 5*nb; 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(lwork,&work)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(jac->A,&a)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(jac->U,&u)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(jac->Vt,&v)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(jac->diag,&d)); 7427c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX) 753f83f0d9SMatthew G Knepley { 763f83f0d9SMatthew G Knepley PetscBLASInt lierr; 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 788b83055fSJed Brown PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr)); 792c71b3e2SJacob Faibussowitsch PetscCheckFalse(lierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPTrapPop()); 813f83f0d9SMatthew G Knepley } 8227c67122SBarry Smith #else 83ef47b4b1SBarry Smith { 84ef47b4b1SBarry Smith PetscBLASInt lierr; 85ef47b4b1SBarry Smith PetscReal *rwork,*dd; 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(5*nb,&rwork)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nb,&dd)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 89ef47b4b1SBarry Smith PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("A","A",&nb,&nb,a,&nb,dd,u,&nb,v,&nb,work,&lwork,rwork,&lierr)); 902c71b3e2SJacob Faibussowitsch PetscCheckFalse(lierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rwork)); 9214ce09a1SBarry Smith for (i=0; i<n; i++) d[i] = dd[i]; 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(dd)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPTrapPop()); 95ef47b4b1SBarry Smith } 9627c67122SBarry Smith #endif 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(jac->A,&a)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(jac->U,&u)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(jac->Vt,&v)); 100426160bdSJed Brown for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break; 101426160bdSJed Brown jac->nzero = n-1-i; 102426160bdSJed Brown if (jac->monitor) { 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 105426160bdSJed Brown if (n >= 10) { /* print 5 smallest and 5 largest */ 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(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]))); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(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]))); 108426160bdSJed Brown } else { /* print all singular values */ 109426160bdSJed Brown char buf[256],*p; 1108caf3d72SBarry Smith size_t left = sizeof(buf),used; 111426160bdSJed Brown PetscInt thisline; 112426160bdSJed Brown for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) { 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]))); 114426160bdSJed Brown left -= used; 115426160bdSJed Brown p += used; 116426160bdSJed Brown if (thisline > 4 || i==0) { 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(jac->monitor," SVD: singular values:%s\n",buf)); 118426160bdSJed Brown p = buf; 119426160bdSJed Brown thisline = 0; 12027c67122SBarry Smith } 121426160bdSJed Brown } 122426160bdSJed Brown } 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel)); 124426160bdSJed Brown } 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 12885032590SJed Brown if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */ 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(pc,"Number of zero or nearly singular values %D\n",jac->nzero)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(jac->diag,&d)); 13127c67122SBarry Smith #if defined(foo) 13227c67122SBarry Smith { 13327c67122SBarry Smith PetscViewer viewer; 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(jac->A,viewer)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(jac->U,viewer)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(jac->Vt,viewer)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(jac->diag,viewer)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(viewer)); 14027c67122SBarry Smith } 14127c67122SBarry Smith #endif 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(work)); 14327c67122SBarry Smith PetscFunctionReturn(0); 14427c67122SBarry Smith } 14527c67122SBarry Smith 1463ed27f31SJed Brown static PetscErrorCode PCSVDGetVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 1473ed27f31SJed Brown { 1483ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 1493ed27f31SJed Brown PetscMPIInt size; 1503ed27f31SJed Brown 1513ed27f31SJed Brown PetscFunctionBegin; 152*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 1530298fd71SBarry Smith *xred = NULL; 1543ed27f31SJed Brown switch (side) { 1553ed27f31SJed Brown case PC_LEFT: 1563ed27f31SJed Brown if (size == 1) *xred = x; 1573ed27f31SJed Brown else { 158*5f80ce2aSJacob Faibussowitsch if (!jac->left2red) CHKERRQ(VecScatterCreateToAll(x,&jac->left2red,&jac->leftred)); 1593ed27f31SJed Brown if (amode & READ) { 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD)); 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD)); 1623ed27f31SJed Brown } 1633ed27f31SJed Brown *xred = jac->leftred; 1643ed27f31SJed Brown } 1653ed27f31SJed Brown break; 1663ed27f31SJed Brown case PC_RIGHT: 1673ed27f31SJed Brown if (size == 1) *xred = x; 1683ed27f31SJed Brown else { 169*5f80ce2aSJacob Faibussowitsch if (!jac->right2red) CHKERRQ(VecScatterCreateToAll(x,&jac->right2red,&jac->rightred)); 1703ed27f31SJed Brown if (amode & READ) { 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD)); 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD)); 1733ed27f31SJed Brown } 1743ed27f31SJed Brown *xred = jac->rightred; 1753ed27f31SJed Brown } 1763ed27f31SJed Brown break; 177ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 1783ed27f31SJed Brown } 1793ed27f31SJed Brown PetscFunctionReturn(0); 1803ed27f31SJed Brown } 1813ed27f31SJed Brown 1823ed27f31SJed Brown static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 1833ed27f31SJed Brown { 1843ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 1853ed27f31SJed Brown PetscMPIInt size; 1863ed27f31SJed Brown 1873ed27f31SJed Brown PetscFunctionBegin; 188*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 1893ed27f31SJed Brown switch (side) { 1903ed27f31SJed Brown case PC_LEFT: 1913ed27f31SJed Brown if (size != 1 && amode & WRITE) { 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE)); 1943ed27f31SJed Brown } 1953ed27f31SJed Brown break; 1963ed27f31SJed Brown case PC_RIGHT: 1973ed27f31SJed Brown if (size != 1 && amode & WRITE) { 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE)); 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE)); 2003ed27f31SJed Brown } 2013ed27f31SJed Brown break; 202ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 2033ed27f31SJed Brown } 2040298fd71SBarry Smith *xred = NULL; 2053ed27f31SJed Brown PetscFunctionReturn(0); 2063ed27f31SJed Brown } 2073ed27f31SJed Brown 20827c67122SBarry Smith /* -------------------------------------------------------------------------- */ 20927c67122SBarry Smith /* 21027c67122SBarry Smith PCApply_SVD - Applies the SVD preconditioner to a vector. 21127c67122SBarry Smith 21227c67122SBarry Smith Input Parameters: 21327c67122SBarry Smith . pc - the preconditioner context 21427c67122SBarry Smith . x - input vector 21527c67122SBarry Smith 21627c67122SBarry Smith Output Parameter: 21727c67122SBarry Smith . y - output vector 21827c67122SBarry Smith 21927c67122SBarry Smith Application Interface Routine: PCApply() 22027c67122SBarry Smith */ 22127c67122SBarry Smith static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y) 22227c67122SBarry Smith { 22327c67122SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 2243ed27f31SJed Brown Vec work = jac->work,xred,yred; 22527c67122SBarry Smith 22627c67122SBarry Smith PetscFunctionBegin; 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDGetVec(pc,PC_RIGHT,READ,x,&xred)); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDGetVec(pc,PC_LEFT,WRITE,y,&yred)); 229ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 230*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(jac->U,xred,work)); 231ef47b4b1SBarry Smith #else 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultHermitianTranspose(jac->U,xred,work)); 233ef47b4b1SBarry Smith #endif 234*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(work,work,jac->diag)); 235ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(jac->Vt,work,yred)); 237ef47b4b1SBarry Smith #else 238*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultHermitianTranspose(jac->Vt,work,yred)); 239ef47b4b1SBarry Smith #endif 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDRestoreVec(pc,PC_RIGHT,READ,x,&xred)); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDRestoreVec(pc,PC_LEFT,WRITE,y,&yred)); 2423ed27f31SJed Brown PetscFunctionReturn(0); 2433ed27f31SJed Brown } 2443ed27f31SJed Brown 2453ed27f31SJed Brown static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y) 2463ed27f31SJed Brown { 2473ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 2483ed27f31SJed Brown Vec work = jac->work,xred,yred; 2493ed27f31SJed Brown 2503ed27f31SJed Brown PetscFunctionBegin; 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDGetVec(pc,PC_LEFT,READ,x,&xred)); 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred)); 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(jac->Vt,xred,work)); 254*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(work,work,jac->diag)); 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(jac->U,work,yred)); 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred)); 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred)); 25827c67122SBarry Smith PetscFunctionReturn(0); 25927c67122SBarry Smith } 26027c67122SBarry Smith 261a2d70de2SBarry Smith static PetscErrorCode PCReset_SVD(PC pc) 262a2d70de2SBarry Smith { 263a2d70de2SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 264a2d70de2SBarry Smith 265a2d70de2SBarry Smith PetscFunctionBegin; 266*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&jac->A)); 267*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&jac->U)); 268*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&jac->Vt)); 269*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&jac->diag)); 270*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&jac->work)); 271*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&jac->right2red)); 272*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&jac->left2red)); 273*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&jac->rightred)); 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&jac->leftred)); 275a2d70de2SBarry Smith PetscFunctionReturn(0); 276a2d70de2SBarry Smith } 277a2d70de2SBarry Smith 27827c67122SBarry Smith /* -------------------------------------------------------------------------- */ 27927c67122SBarry Smith /* 28027c67122SBarry Smith PCDestroy_SVD - Destroys the private context for the SVD preconditioner 28127c67122SBarry Smith that was created with PCCreate_SVD(). 28227c67122SBarry Smith 28327c67122SBarry Smith Input Parameter: 28427c67122SBarry Smith . pc - the preconditioner context 28527c67122SBarry Smith 28627c67122SBarry Smith Application Interface Routine: PCDestroy() 28727c67122SBarry Smith */ 28827c67122SBarry Smith static PetscErrorCode PCDestroy_SVD(PC pc) 28927c67122SBarry Smith { 290426160bdSJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 29127c67122SBarry Smith 29227c67122SBarry Smith PetscFunctionBegin; 293*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_SVD(pc)); 294*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&jac->monitor)); 295*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc->data)); 29627c67122SBarry Smith PetscFunctionReturn(0); 29727c67122SBarry Smith } 29827c67122SBarry Smith 2994416b707SBarry Smith static PetscErrorCode PCSetFromOptions_SVD(PetscOptionItems *PetscOptionsObject,PC pc) 30027c67122SBarry Smith { 3018f1a2a5eSBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 302426160bdSJed Brown PetscBool flg,set; 30327c67122SBarry Smith 30427c67122SBarry Smith PetscFunctionBegin; 305*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHead(PetscOptionsObject,"SVD options")); 306*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,NULL)); 307*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,NULL)); 308*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set)); 309426160bdSJed Brown if (set) { /* Should make PCSVDSetMonitor() */ 310426160bdSJed Brown if (flg && !jac->monitor) { 311*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),"stdout",&jac->monitor)); 312426160bdSJed Brown } else if (!flg) { 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&jac->monitor)); 314426160bdSJed Brown } 315426160bdSJed Brown } 316*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsTail()); 31727c67122SBarry Smith PetscFunctionReturn(0); 31827c67122SBarry Smith } 31927c67122SBarry Smith 32033761216SBarry Smith static PetscErrorCode PCView_SVD(PC pc,PetscViewer viewer) 32133761216SBarry Smith { 32233761216SBarry Smith PC_SVD *svd = (PC_SVD*)pc->data; 32333761216SBarry Smith PetscBool iascii; 32433761216SBarry Smith 32533761216SBarry Smith PetscFunctionBegin; 326*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 32733761216SBarry Smith if (iascii) { 328*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," All singular values smaller than %g treated as zero\n",(double)svd->zerosing)); 329*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," Provided essential rank of the matrix %D (all other eigenvalues are zeroed)\n",svd->essrank)); 33033761216SBarry Smith } 33133761216SBarry Smith PetscFunctionReturn(0); 33233761216SBarry Smith } 33327c67122SBarry Smith /* -------------------------------------------------------------------------- */ 33427c67122SBarry Smith /* 33527c67122SBarry Smith PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 33627c67122SBarry Smith and sets this as the private data within the generic preconditioning 33727c67122SBarry Smith context, PC, that was created within PCCreate(). 33827c67122SBarry Smith 33927c67122SBarry Smith Input Parameter: 34027c67122SBarry Smith . pc - the preconditioner context 34127c67122SBarry Smith 34227c67122SBarry Smith Application Interface Routine: PCCreate() 34327c67122SBarry Smith */ 34427c67122SBarry Smith 34527c67122SBarry Smith /*MC 34627c67122SBarry Smith PCSVD - Use pseudo inverse defined by SVD of operator 34727c67122SBarry Smith 34827c67122SBarry Smith Level: advanced 34927c67122SBarry Smith 3508f1a2a5eSBarry Smith Options Database: 351147403d9SBarry Smith + -pc_svd_zero_sing <rtol> - Singular values smaller than this are treated as zero 352147403d9SBarry Smith - -pc_svd_monitor - Print information on the extreme singular values of the operator 35327c67122SBarry Smith 354a2b725a8SWilliam Gropp Developer Note: 355a2b725a8SWilliam Gropp This implementation automatically creates a redundant copy of the 3568997ae2eSBarry Smith matrix on each process and uses a sequential SVD solve. Why does it do this instead 3578997ae2eSBarry Smith of using the composable PCREDUNDANT object? 3588997ae2eSBarry Smith 35927c67122SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 36027c67122SBarry Smith M*/ 36127c67122SBarry Smith 3628cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc) 36327c67122SBarry Smith { 36427c67122SBarry Smith PC_SVD *jac; 36527c67122SBarry Smith 36627c67122SBarry Smith PetscFunctionBegin; 36727c67122SBarry Smith /* 36827c67122SBarry Smith Creates the private data structure for this preconditioner and 36927c67122SBarry Smith attach it to the PC object. 37027c67122SBarry Smith */ 371*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc,&jac)); 3728f1a2a5eSBarry Smith jac->zerosing = 1.e-12; 37327c67122SBarry Smith pc->data = (void*)jac; 37427c67122SBarry Smith 37527c67122SBarry Smith /* 37627c67122SBarry Smith Set the pointers for the functions that are provided above. 37727c67122SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 37827c67122SBarry Smith are called, they will automatically call these functions. Note we 37927c67122SBarry Smith choose not to provide a couple of these functions since they are 38027c67122SBarry Smith not needed. 38127c67122SBarry Smith */ 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