127c67122SBarry Smith 2*b45d2f2cSJed 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; 1027c67122SBarry Smith Mat A,U,V; 1127c67122SBarry Smith PetscInt nzero; 128f1a2a5eSBarry Smith PetscReal zerosing; /* measure of smallest singular value treated as nonzero */ 13426160bdSJed Brown PetscViewer monitor; 1427c67122SBarry Smith } PC_SVD; 1527c67122SBarry Smith 1627c67122SBarry Smith 1727c67122SBarry Smith /* -------------------------------------------------------------------------- */ 1827c67122SBarry Smith /* 1927c67122SBarry Smith PCSetUp_SVD - Prepares for the use of the SVD preconditioner 2027c67122SBarry Smith by setting data structures and options. 2127c67122SBarry Smith 2227c67122SBarry Smith Input Parameter: 2327c67122SBarry Smith . pc - the preconditioner context 2427c67122SBarry Smith 2527c67122SBarry Smith Application Interface Routine: PCSetUp() 2627c67122SBarry Smith 2727c67122SBarry Smith Notes: 2827c67122SBarry Smith The interface routine PCSetUp() is not usually called directly by 2927c67122SBarry Smith the user, but instead is called by PCApply() if necessary. 3027c67122SBarry Smith */ 3127c67122SBarry Smith #undef __FUNCT__ 3227c67122SBarry Smith #define __FUNCT__ "PCSetUp_SVD" 3327c67122SBarry Smith static PetscErrorCode PCSetUp_SVD(PC pc) 3427c67122SBarry Smith { 35a3c4e3ecSJed Brown #if defined(PETSC_MISSING_LAPACK_GESVD) 36a3c4e3ecSJed Brown SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates."); 37a3c4e3ecSJed Brown #else 3827c67122SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 3927c67122SBarry Smith PetscErrorCode ierr; 4027c67122SBarry Smith PetscScalar *a,*u,*v,*d,*work; 413f83f0d9SMatthew G Knepley PetscBLASInt nb,lwork; 4227c67122SBarry Smith PetscInt i,n; 4327c67122SBarry Smith 4427c67122SBarry Smith PetscFunctionBegin; 4527c67122SBarry Smith if (!jac->diag) { 4627c67122SBarry Smith /* assume square matrices */ 478f1a2a5eSBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diag,&jac->work);CHKERRQ(ierr); 4827c67122SBarry Smith } 496bf464f9SBarry Smith ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 508f1a2a5eSBarry Smith ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr); 5127c67122SBarry Smith if (!jac->U) { 5227c67122SBarry Smith ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr); 5327c67122SBarry Smith ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->V);CHKERRQ(ierr); 5427c67122SBarry Smith } 558f1a2a5eSBarry Smith ierr = MatGetSize(pc->pmat,&n,PETSC_NULL);CHKERRQ(ierr); 5627c67122SBarry Smith nb = PetscBLASIntCast(n); 5727c67122SBarry Smith lwork = 5*nb; 5827c67122SBarry Smith ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr); 5927c67122SBarry Smith ierr = MatGetArray(jac->A,&a);CHKERRQ(ierr); 6027c67122SBarry Smith ierr = MatGetArray(jac->U,&u);CHKERRQ(ierr); 6127c67122SBarry Smith ierr = MatGetArray(jac->V,&v);CHKERRQ(ierr); 6227c67122SBarry Smith ierr = VecGetArray(jac->diag,&d);CHKERRQ(ierr); 6327c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX) 643f83f0d9SMatthew G Knepley { 653f83f0d9SMatthew G Knepley PetscBLASInt lierr; 66670f3ff9SJed Brown ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 6727c67122SBarry Smith LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr); 6898b909ebSSatish Balay if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 69670f3ff9SJed Brown ierr = PetscFPTrapPop();CHKERRQ(ierr); 703f83f0d9SMatthew G Knepley } 7127c67122SBarry Smith #else 7227c67122SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex"); 7327c67122SBarry Smith #endif 7427c67122SBarry Smith ierr = MatRestoreArray(jac->A,&a);CHKERRQ(ierr); 7527c67122SBarry Smith ierr = MatRestoreArray(jac->U,&u);CHKERRQ(ierr); 7627c67122SBarry Smith ierr = MatRestoreArray(jac->V,&v);CHKERRQ(ierr); 77426160bdSJed Brown for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break; 78426160bdSJed Brown jac->nzero = n-1-i; 79426160bdSJed Brown if (jac->monitor) { 80426160bdSJed Brown ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 81426160bdSJed 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); 82426160bdSJed Brown if (n >= 10) { /* print 5 smallest and 5 largest */ 83426160bdSJed 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); 84426160bdSJed 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); 85426160bdSJed Brown } else { /* print all singular values */ 86426160bdSJed Brown char buf[256],*p; 87426160bdSJed Brown size_t left = sizeof buf,used; 88426160bdSJed Brown PetscInt thisline; 89426160bdSJed Brown for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) { 90426160bdSJed Brown ierr = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr); 91426160bdSJed Brown left -= used; 92426160bdSJed Brown p += used; 93426160bdSJed Brown if (thisline > 4 || i==0) { 94426160bdSJed Brown ierr = PetscViewerASCIIPrintf(jac->monitor," SVD: singular values:%s\n",buf);CHKERRQ(ierr); 95426160bdSJed Brown p = buf; 96426160bdSJed Brown thisline = 0; 9727c67122SBarry Smith } 98426160bdSJed Brown } 99426160bdSJed Brown } 100426160bdSJed Brown ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 101426160bdSJed Brown } 102426160bdSJed Brown ierr = PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1])); 103426160bdSJed Brown for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i]; 104426160bdSJed Brown for (; i<n; i++) d[i] = 0.0; 10527c67122SBarry Smith ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero); 10627c67122SBarry Smith ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr); 10727c67122SBarry Smith #if defined(foo) 10827c67122SBarry Smith { 10927c67122SBarry Smith PetscViewer viewer; 11027c67122SBarry Smith ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 11127c67122SBarry Smith ierr = MatView(jac->A,viewer);CHKERRQ(ierr); 11227c67122SBarry Smith ierr = MatView(jac->U,viewer);CHKERRQ(ierr); 11327c67122SBarry Smith ierr = MatView(jac->V,viewer);CHKERRQ(ierr); 11427c67122SBarry Smith ierr = VecView(jac->diag,viewer);CHKERRQ(ierr); 11527c67122SBarry Smith ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 11627c67122SBarry Smith } 11727c67122SBarry Smith #endif 11827c67122SBarry Smith ierr = PetscFree(work); 11927c67122SBarry Smith PetscFunctionReturn(0); 120a3c4e3ecSJed Brown #endif 12127c67122SBarry Smith } 12227c67122SBarry Smith 12327c67122SBarry Smith /* -------------------------------------------------------------------------- */ 12427c67122SBarry Smith /* 12527c67122SBarry Smith PCApply_SVD - Applies the SVD preconditioner to a vector. 12627c67122SBarry Smith 12727c67122SBarry Smith Input Parameters: 12827c67122SBarry Smith . pc - the preconditioner context 12927c67122SBarry Smith . x - input vector 13027c67122SBarry Smith 13127c67122SBarry Smith Output Parameter: 13227c67122SBarry Smith . y - output vector 13327c67122SBarry Smith 13427c67122SBarry Smith Application Interface Routine: PCApply() 13527c67122SBarry Smith */ 13627c67122SBarry Smith #undef __FUNCT__ 13727c67122SBarry Smith #define __FUNCT__ "PCApply_SVD" 13827c67122SBarry Smith static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y) 13927c67122SBarry Smith { 14027c67122SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 14127c67122SBarry Smith Vec work = jac->work; 14227c67122SBarry Smith PetscErrorCode ierr; 14327c67122SBarry Smith 14427c67122SBarry Smith PetscFunctionBegin; 14527c67122SBarry Smith ierr = MatMultTranspose(jac->U,x,work);CHKERRQ(ierr); 14627c67122SBarry Smith ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 14727c67122SBarry Smith ierr = MatMultTranspose(jac->V,work,y);CHKERRQ(ierr); 14827c67122SBarry Smith PetscFunctionReturn(0); 14927c67122SBarry Smith } 15027c67122SBarry Smith 151a2d70de2SBarry Smith #undef __FUNCT__ 152a2d70de2SBarry Smith #define __FUNCT__ "PCReset_SVD" 153a2d70de2SBarry Smith static PetscErrorCode PCReset_SVD(PC pc) 154a2d70de2SBarry Smith { 155a2d70de2SBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 156a2d70de2SBarry Smith PetscErrorCode ierr; 157a2d70de2SBarry Smith 158a2d70de2SBarry Smith PetscFunctionBegin; 159a2d70de2SBarry Smith ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 160a2d70de2SBarry Smith ierr = MatDestroy(&jac->U);CHKERRQ(ierr); 161a2d70de2SBarry Smith ierr = MatDestroy(&jac->V);CHKERRQ(ierr); 162a2d70de2SBarry Smith ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 163a2d70de2SBarry Smith ierr = VecDestroy(&jac->work);CHKERRQ(ierr); 164a2d70de2SBarry Smith PetscFunctionReturn(0); 165a2d70de2SBarry Smith } 166a2d70de2SBarry Smith 16727c67122SBarry Smith /* -------------------------------------------------------------------------- */ 16827c67122SBarry Smith /* 16927c67122SBarry Smith PCDestroy_SVD - Destroys the private context for the SVD preconditioner 17027c67122SBarry Smith that was created with PCCreate_SVD(). 17127c67122SBarry Smith 17227c67122SBarry Smith Input Parameter: 17327c67122SBarry Smith . pc - the preconditioner context 17427c67122SBarry Smith 17527c67122SBarry Smith Application Interface Routine: PCDestroy() 17627c67122SBarry Smith */ 17727c67122SBarry Smith #undef __FUNCT__ 17827c67122SBarry Smith #define __FUNCT__ "PCDestroy_SVD" 17927c67122SBarry Smith static PetscErrorCode PCDestroy_SVD(PC pc) 18027c67122SBarry Smith { 181426160bdSJed Brown PC_SVD *jac = (PC_SVD*)pc->data; 18227c67122SBarry Smith PetscErrorCode ierr; 18327c67122SBarry Smith 18427c67122SBarry Smith PetscFunctionBegin; 185a2d70de2SBarry Smith ierr = PCReset_SVD(pc);CHKERRQ(ierr); 186426160bdSJed Brown ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 187c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 18827c67122SBarry Smith PetscFunctionReturn(0); 18927c67122SBarry Smith } 19027c67122SBarry Smith 19127c67122SBarry Smith #undef __FUNCT__ 19227c67122SBarry Smith #define __FUNCT__ "PCSetFromOptions_SVD" 19327c67122SBarry Smith static PetscErrorCode PCSetFromOptions_SVD(PC pc) 19427c67122SBarry Smith { 19527c67122SBarry Smith PetscErrorCode ierr; 1968f1a2a5eSBarry Smith PC_SVD *jac = (PC_SVD*)pc->data; 197426160bdSJed Brown PetscBool flg,set; 19827c67122SBarry Smith 19927c67122SBarry Smith PetscFunctionBegin; 20027c67122SBarry Smith ierr = PetscOptionsHead("SVD options");CHKERRQ(ierr); 2018f1a2a5eSBarry Smith ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,PETSC_NULL);CHKERRQ(ierr); 2026ba663aaSJed Brown ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor?PETSC_TRUE:PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 203426160bdSJed Brown if (set) { /* Should make PCSVDSetMonitor() */ 204426160bdSJed Brown if (flg && !jac->monitor) { 205426160bdSJed Brown ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,"stdout",&jac->monitor);CHKERRQ(ierr); 206426160bdSJed Brown } else if (!flg) { 207426160bdSJed Brown ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 208426160bdSJed Brown } 209426160bdSJed Brown } 21027c67122SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 21127c67122SBarry Smith PetscFunctionReturn(0); 21227c67122SBarry Smith } 21327c67122SBarry Smith 21427c67122SBarry Smith /* -------------------------------------------------------------------------- */ 21527c67122SBarry Smith /* 21627c67122SBarry Smith PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 21727c67122SBarry Smith and sets this as the private data within the generic preconditioning 21827c67122SBarry Smith context, PC, that was created within PCCreate(). 21927c67122SBarry Smith 22027c67122SBarry Smith Input Parameter: 22127c67122SBarry Smith . pc - the preconditioner context 22227c67122SBarry Smith 22327c67122SBarry Smith Application Interface Routine: PCCreate() 22427c67122SBarry Smith */ 22527c67122SBarry Smith 22627c67122SBarry Smith /*MC 22727c67122SBarry Smith PCSVD - Use pseudo inverse defined by SVD of operator 22827c67122SBarry Smith 22927c67122SBarry Smith Level: advanced 23027c67122SBarry Smith 23127c67122SBarry Smith Concepts: SVD 23227c67122SBarry Smith 2338f1a2a5eSBarry Smith Options Database: 2348b0a5094SBarry Smith - -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero 2358b0a5094SBarry Smith + -pc_svd_monitor Print information on the extreme singular values of the operator 23627c67122SBarry Smith 23727c67122SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 23827c67122SBarry Smith M*/ 23927c67122SBarry Smith 24027c67122SBarry Smith EXTERN_C_BEGIN 24127c67122SBarry Smith #undef __FUNCT__ 24227c67122SBarry Smith #define __FUNCT__ "PCCreate_SVD" 2437087cfbeSBarry Smith PetscErrorCode PCCreate_SVD(PC pc) 24427c67122SBarry Smith { 24527c67122SBarry Smith PC_SVD *jac; 24627c67122SBarry Smith PetscErrorCode ierr; 24727c67122SBarry Smith 24827c67122SBarry Smith PetscFunctionBegin; 24927c67122SBarry Smith /* 25027c67122SBarry Smith Creates the private data structure for this preconditioner and 25127c67122SBarry Smith attach it to the PC object. 25227c67122SBarry Smith */ 25327c67122SBarry Smith ierr = PetscNewLog(pc,PC_SVD,&jac);CHKERRQ(ierr); 2548f1a2a5eSBarry Smith jac->zerosing = 1.e-12; 25527c67122SBarry Smith pc->data = (void*)jac; 25627c67122SBarry Smith 25727c67122SBarry Smith /* 25827c67122SBarry Smith Set the pointers for the functions that are provided above. 25927c67122SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 26027c67122SBarry Smith are called, they will automatically call these functions. Note we 26127c67122SBarry Smith choose not to provide a couple of these functions since they are 26227c67122SBarry Smith not needed. 26327c67122SBarry Smith */ 26427c67122SBarry Smith pc->ops->apply = PCApply_SVD; 26527c67122SBarry Smith pc->ops->applytranspose = PCApply_SVD; 26627c67122SBarry Smith pc->ops->setup = PCSetUp_SVD; 267a2d70de2SBarry Smith pc->ops->reset = PCReset_SVD; 26827c67122SBarry Smith pc->ops->destroy = PCDestroy_SVD; 26927c67122SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SVD; 27027c67122SBarry Smith pc->ops->view = 0; 27127c67122SBarry Smith pc->ops->applyrichardson = 0; 27227c67122SBarry Smith PetscFunctionReturn(0); 27327c67122SBarry Smith } 27427c67122SBarry Smith EXTERN_C_END 27527c67122SBarry Smith 276