127c67122SBarry Smith 2b45d2f2cSJed Brown #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; 54*b2bf6370SHong Zhang ierr = MatGetRedundantMatrix(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 */ 613ed27f31SJed Brown ierr = MatGetVecs(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 } 670298fd71SBarry Smith ierr = MatGetSize(pc->pmat,&n,NULL);CHKERRQ(ierr); 68c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&nb);CHKERRQ(ierr); 6927c67122SBarry Smith lwork = 5*nb; 7027c67122SBarry Smith ierr = PetscMalloc(lwork*sizeof(PetscScalar),&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)); 8098b909ebSSatish Balay if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 81670f3ff9SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 823f83f0d9SMatthew G Knepley } 8327c67122SBarry Smith #else 8427c67122SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex"); 8527c67122SBarry Smith #endif 868c778c55SBarry Smith ierr = MatDenseRestoreArray(jac->A,&a);CHKERRQ(ierr); 878c778c55SBarry Smith ierr = MatDenseRestoreArray(jac->U,&u);CHKERRQ(ierr); 888c778c55SBarry Smith ierr = MatDenseRestoreArray(jac->Vt,&v);CHKERRQ(ierr); 89426160bdSJed Brown for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break; 90426160bdSJed Brown jac->nzero = n-1-i; 91426160bdSJed Brown if (jac->monitor) { 92426160bdSJed Brown ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 93426160bdSJed 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); 94426160bdSJed Brown if (n >= 10) { /* print 5 smallest and 5 largest */ 95426160bdSJed 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); 96426160bdSJed 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); 97426160bdSJed Brown } else { /* print all singular values */ 98426160bdSJed Brown char buf[256],*p; 998caf3d72SBarry Smith size_t left = sizeof(buf),used; 100426160bdSJed Brown PetscInt thisline; 101426160bdSJed Brown for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) { 102426160bdSJed Brown ierr = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr); 103426160bdSJed Brown left -= used; 104426160bdSJed Brown p += used; 105426160bdSJed Brown if (thisline > 4 || i==0) { 106426160bdSJed Brown ierr = PetscViewerASCIIPrintf(jac->monitor," SVD: singular values:%s\n",buf);CHKERRQ(ierr); 107426160bdSJed Brown p = buf; 108426160bdSJed Brown thisline = 0; 10927c67122SBarry Smith } 110426160bdSJed Brown } 111426160bdSJed Brown } 112426160bdSJed Brown ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 113426160bdSJed Brown } 11422d28d08SBarry 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); 115426160bdSJed Brown for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i]; 116426160bdSJed Brown for (; i<n; i++) d[i] = 0.0; 11785032590SJed Brown if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */ 11822d28d08SBarry Smith ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);CHKERRQ(ierr); 11927c67122SBarry Smith ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr); 12027c67122SBarry Smith #if defined(foo) 12127c67122SBarry Smith { 12227c67122SBarry Smith PetscViewer viewer; 12327c67122SBarry Smith ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 12427c67122SBarry Smith ierr = MatView(jac->A,viewer);CHKERRQ(ierr); 12527c67122SBarry Smith ierr = MatView(jac->U,viewer);CHKERRQ(ierr); 1263ed27f31SJed Brown ierr = MatView(jac->Vt,viewer);CHKERRQ(ierr); 12727c67122SBarry Smith ierr = VecView(jac->diag,viewer);CHKERRQ(ierr); 12827c67122SBarry Smith ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 12927c67122SBarry Smith } 13027c67122SBarry Smith #endif 13122d28d08SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 13227c67122SBarry Smith PetscFunctionReturn(0); 133a3c4e3ecSJed Brown #endif 13427c67122SBarry Smith } 13527c67122SBarry Smith 1363ed27f31SJed Brown #undef __FUNCT__ 1373ed27f31SJed Brown #define __FUNCT__ "PCSVDGetVec" 1383ed27f31SJed Brown static PetscErrorCode PCSVDGetVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 1393ed27f31SJed Brown { 1403ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 1413ed27f31SJed Brown PetscErrorCode ierr; 1423ed27f31SJed Brown PetscMPIInt size; 1433ed27f31SJed Brown 1443ed27f31SJed Brown PetscFunctionBegin; 145ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 1460298fd71SBarry Smith *xred = NULL; 1473ed27f31SJed Brown switch (side) { 1483ed27f31SJed Brown case PC_LEFT: 1493ed27f31SJed Brown if (size == 1) *xred = x; 1503ed27f31SJed Brown else { 1513ed27f31SJed Brown if (!jac->left2red) {ierr = VecScatterCreateToAll(x,&jac->left2red,&jac->leftred);CHKERRQ(ierr);} 1523ed27f31SJed Brown if (amode & READ) { 1533ed27f31SJed Brown ierr = VecScatterBegin(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1543ed27f31SJed Brown ierr = VecScatterEnd(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1553ed27f31SJed Brown } 1563ed27f31SJed Brown *xred = jac->leftred; 1573ed27f31SJed Brown } 1583ed27f31SJed Brown break; 1593ed27f31SJed Brown case PC_RIGHT: 1603ed27f31SJed Brown if (size == 1) *xred = x; 1613ed27f31SJed Brown else { 1623ed27f31SJed Brown if (!jac->right2red) {ierr = VecScatterCreateToAll(x,&jac->right2red,&jac->rightred);CHKERRQ(ierr);} 1633ed27f31SJed Brown if (amode & READ) { 1643ed27f31SJed Brown ierr = VecScatterBegin(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1653ed27f31SJed Brown ierr = VecScatterEnd(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1663ed27f31SJed Brown } 1673ed27f31SJed Brown *xred = jac->rightred; 1683ed27f31SJed Brown } 1693ed27f31SJed Brown break; 170ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 1713ed27f31SJed Brown } 1723ed27f31SJed Brown PetscFunctionReturn(0); 1733ed27f31SJed Brown } 1743ed27f31SJed Brown 1753ed27f31SJed Brown #undef __FUNCT__ 1763ed27f31SJed Brown #define __FUNCT__ "PCSVDRestoreVec" 1773ed27f31SJed Brown static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 1783ed27f31SJed Brown { 1793ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 1803ed27f31SJed Brown PetscErrorCode ierr; 1813ed27f31SJed Brown PetscMPIInt size; 1823ed27f31SJed Brown 1833ed27f31SJed Brown PetscFunctionBegin; 184ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 1853ed27f31SJed Brown switch (side) { 1863ed27f31SJed Brown case PC_LEFT: 1873ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1883ed27f31SJed Brown ierr = VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1893ed27f31SJed Brown ierr = VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1903ed27f31SJed Brown } 1913ed27f31SJed Brown break; 1923ed27f31SJed Brown case PC_RIGHT: 1933ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1943ed27f31SJed Brown ierr = VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1953ed27f31SJed Brown ierr = VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1963ed27f31SJed Brown } 1973ed27f31SJed Brown break; 198ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 1993ed27f31SJed Brown } 2000298fd71SBarry Smith *xred = NULL; 2013ed27f31SJed Brown PetscFunctionReturn(0); 2023ed27f31SJed Brown } 2033ed27f31SJed Brown 20427c67122SBarry Smith /* -------------------------------------------------------------------------- */ 20527c67122SBarry Smith /* 20627c67122SBarry Smith PCApply_SVD - Applies the SVD preconditioner to a vector. 20727c67122SBarry Smith 20827c67122SBarry Smith Input Parameters: 20927c67122SBarry Smith . pc - the preconditioner context 21027c67122SBarry Smith . x - input vector 21127c67122SBarry Smith 21227c67122SBarry Smith Output Parameter: 21327c67122SBarry Smith . y - output vector 21427c67122SBarry Smith 21527c67122SBarry Smith Application Interface Routine: PCApply() 21627c67122SBarry Smith */ 21727c67122SBarry Smith #undef __FUNCT__ 21827c67122SBarry Smith #define __FUNCT__ "PCApply_SVD" 21927c67122SBarry Smith static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y) 22027c67122SBarry Smith { 22127c67122SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 2223ed27f31SJed Brown Vec work = jac->work,xred,yred; 22327c67122SBarry Smith PetscErrorCode ierr; 22427c67122SBarry Smith 22527c67122SBarry Smith PetscFunctionBegin; 2263ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr); 2273ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr); 2283ed27f31SJed Brown ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr); 22927c67122SBarry Smith ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 2303ed27f31SJed Brown ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr); 2313ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr); 2323ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr); 2333ed27f31SJed Brown PetscFunctionReturn(0); 2343ed27f31SJed Brown } 2353ed27f31SJed Brown 2363ed27f31SJed Brown #undef __FUNCT__ 2373ed27f31SJed Brown #define __FUNCT__ "PCApplyTranspose_SVD" 2383ed27f31SJed Brown static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y) 2393ed27f31SJed Brown { 2403ed27f31SJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 2413ed27f31SJed Brown Vec work = jac->work,xred,yred; 2423ed27f31SJed Brown PetscErrorCode ierr; 2433ed27f31SJed Brown 2443ed27f31SJed Brown PetscFunctionBegin; 2453ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr); 2463ed27f31SJed Brown ierr = PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr); 2473ed27f31SJed Brown ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr); 2483ed27f31SJed Brown ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 2493ed27f31SJed Brown ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr); 2503ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr); 2513ed27f31SJed Brown ierr = PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr); 25227c67122SBarry Smith PetscFunctionReturn(0); 25327c67122SBarry Smith } 25427c67122SBarry Smith 255a2d70de2SBarry Smith #undef __FUNCT__ 256a2d70de2SBarry Smith #define __FUNCT__ "PCReset_SVD" 257a2d70de2SBarry Smith static PetscErrorCode PCReset_SVD(PC pc) 258a2d70de2SBarry Smith { 259a2d70de2SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 260a2d70de2SBarry Smith PetscErrorCode ierr; 261a2d70de2SBarry Smith 262a2d70de2SBarry Smith PetscFunctionBegin; 263a2d70de2SBarry Smith ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 264a2d70de2SBarry Smith ierr = MatDestroy(&jac->U);CHKERRQ(ierr); 2653ed27f31SJed Brown ierr = MatDestroy(&jac->Vt);CHKERRQ(ierr); 266a2d70de2SBarry Smith ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 267a2d70de2SBarry Smith ierr = VecDestroy(&jac->work);CHKERRQ(ierr); 2683ed27f31SJed Brown ierr = VecScatterDestroy(&jac->right2red);CHKERRQ(ierr); 2693ed27f31SJed Brown ierr = VecScatterDestroy(&jac->left2red);CHKERRQ(ierr); 2703ed27f31SJed Brown ierr = VecDestroy(&jac->rightred);CHKERRQ(ierr); 2713ed27f31SJed Brown ierr = VecDestroy(&jac->leftred);CHKERRQ(ierr); 272a2d70de2SBarry Smith PetscFunctionReturn(0); 273a2d70de2SBarry Smith } 274a2d70de2SBarry Smith 27527c67122SBarry Smith /* -------------------------------------------------------------------------- */ 27627c67122SBarry Smith /* 27727c67122SBarry Smith PCDestroy_SVD - Destroys the private context for the SVD preconditioner 27827c67122SBarry Smith that was created with PCCreate_SVD(). 27927c67122SBarry Smith 28027c67122SBarry Smith Input Parameter: 28127c67122SBarry Smith . pc - the preconditioner context 28227c67122SBarry Smith 28327c67122SBarry Smith Application Interface Routine: PCDestroy() 28427c67122SBarry Smith */ 28527c67122SBarry Smith #undef __FUNCT__ 28627c67122SBarry Smith #define __FUNCT__ "PCDestroy_SVD" 28727c67122SBarry Smith static PetscErrorCode PCDestroy_SVD(PC pc) 28827c67122SBarry Smith { 289426160bdSJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 29027c67122SBarry Smith PetscErrorCode ierr; 29127c67122SBarry Smith 29227c67122SBarry Smith PetscFunctionBegin; 293a2d70de2SBarry Smith ierr = PCReset_SVD(pc);CHKERRQ(ierr); 294426160bdSJed Brown ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 295c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 29627c67122SBarry Smith PetscFunctionReturn(0); 29727c67122SBarry Smith } 29827c67122SBarry Smith 29927c67122SBarry Smith #undef __FUNCT__ 30027c67122SBarry Smith #define __FUNCT__ "PCSetFromOptions_SVD" 30127c67122SBarry Smith static PetscErrorCode PCSetFromOptions_SVD(PC pc) 30227c67122SBarry Smith { 30327c67122SBarry Smith PetscErrorCode ierr; 3048f1a2a5eSBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 305426160bdSJed Brown PetscBool flg,set; 30627c67122SBarry Smith 30727c67122SBarry Smith PetscFunctionBegin; 30827c67122SBarry Smith ierr = PetscOptionsHead("SVD options");CHKERRQ(ierr); 3090298fd71SBarry Smith ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,NULL);CHKERRQ(ierr); 3100298fd71SBarry Smith ierr = PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,NULL);CHKERRQ(ierr); 3116ba663aaSJed Brown ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 312426160bdSJed Brown if (set) { /* Should make PCSVDSetMonitor() */ 313426160bdSJed Brown if (flg && !jac->monitor) { 314ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),"stdout",&jac->monitor);CHKERRQ(ierr); 315426160bdSJed Brown } else if (!flg) { 316426160bdSJed Brown ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 317426160bdSJed Brown } 318426160bdSJed Brown } 31927c67122SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 32027c67122SBarry Smith PetscFunctionReturn(0); 32127c67122SBarry Smith } 32227c67122SBarry Smith 32327c67122SBarry Smith /* -------------------------------------------------------------------------- */ 32427c67122SBarry Smith /* 32527c67122SBarry Smith PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 32627c67122SBarry Smith and sets this as the private data within the generic preconditioning 32727c67122SBarry Smith context, PC, that was created within PCCreate(). 32827c67122SBarry Smith 32927c67122SBarry Smith Input Parameter: 33027c67122SBarry Smith . pc - the preconditioner context 33127c67122SBarry Smith 33227c67122SBarry Smith Application Interface Routine: PCCreate() 33327c67122SBarry Smith */ 33427c67122SBarry Smith 33527c67122SBarry Smith /*MC 33627c67122SBarry Smith PCSVD - Use pseudo inverse defined by SVD of operator 33727c67122SBarry Smith 33827c67122SBarry Smith Level: advanced 33927c67122SBarry Smith 34027c67122SBarry Smith Concepts: SVD 34127c67122SBarry Smith 3428f1a2a5eSBarry Smith Options Database: 3438b0a5094SBarry Smith - -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero 3448b0a5094SBarry Smith + -pc_svd_monitor Print information on the extreme singular values of the operator 34527c67122SBarry Smith 34627c67122SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 34727c67122SBarry Smith M*/ 34827c67122SBarry Smith 34927c67122SBarry Smith #undef __FUNCT__ 35027c67122SBarry Smith #define __FUNCT__ "PCCreate_SVD" 3518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc) 35227c67122SBarry Smith { 35327c67122SBarry Smith PC_SVD *jac; 35427c67122SBarry Smith PetscErrorCode ierr; 35527c67122SBarry Smith 35627c67122SBarry Smith PetscFunctionBegin; 35727c67122SBarry Smith /* 35827c67122SBarry Smith Creates the private data structure for this preconditioner and 35927c67122SBarry Smith attach it to the PC object. 36027c67122SBarry Smith */ 36127c67122SBarry Smith ierr = PetscNewLog(pc,PC_SVD,&jac);CHKERRQ(ierr); 3628f1a2a5eSBarry Smith jac->zerosing = 1.e-12; 36327c67122SBarry Smith pc->data = (void*)jac; 36427c67122SBarry Smith 36527c67122SBarry Smith /* 36627c67122SBarry Smith Set the pointers for the functions that are provided above. 36727c67122SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 36827c67122SBarry Smith are called, they will automatically call these functions. Note we 36927c67122SBarry Smith choose not to provide a couple of these functions since they are 37027c67122SBarry Smith not needed. 37127c67122SBarry Smith */ 37227c67122SBarry Smith pc->ops->apply = PCApply_SVD; 3733ed27f31SJed Brown pc->ops->applytranspose = PCApplyTranspose_SVD; 37427c67122SBarry Smith pc->ops->setup = PCSetUp_SVD; 375a2d70de2SBarry Smith pc->ops->reset = PCReset_SVD; 37627c67122SBarry Smith pc->ops->destroy = PCDestroy_SVD; 37727c67122SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SVD; 37827c67122SBarry Smith pc->ops->view = 0; 37927c67122SBarry Smith pc->ops->applyrichardson = 0; 38027c67122SBarry Smith PetscFunctionReturn(0); 38127c67122SBarry Smith } 38227c67122SBarry Smith 383