1 2 #include <private/pcimpl.h> /*I "petscpc.h" I*/ 3 #include <petscblaslapack.h> 4 5 /* 6 Private context (data structure) for the SVD preconditioner. 7 */ 8 typedef struct { 9 Vec diag,work; 10 Mat A,U,V; 11 PetscInt nzero; 12 PetscReal zerosing; /* measure of smallest singular value treated as nonzero */ 13 PetscViewer monitor; 14 } PC_SVD; 15 16 17 /* -------------------------------------------------------------------------- */ 18 /* 19 PCSetUp_SVD - Prepares for the use of the SVD preconditioner 20 by setting data structures and options. 21 22 Input Parameter: 23 . pc - the preconditioner context 24 25 Application Interface Routine: PCSetUp() 26 27 Notes: 28 The interface routine PCSetUp() is not usually called directly by 29 the user, but instead is called by PCApply() if necessary. 30 */ 31 #undef __FUNCT__ 32 #define __FUNCT__ "PCSetUp_SVD" 33 static PetscErrorCode PCSetUp_SVD(PC pc) 34 { 35 #if defined(PETSC_MISSING_LAPACK_GESVD) 36 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates."); 37 #else 38 PC_SVD *jac = (PC_SVD*)pc->data; 39 PetscErrorCode ierr; 40 PetscScalar *a,*u,*v,*d,*work; 41 PetscBLASInt nb,lwork; 42 PetscInt i,n; 43 44 PetscFunctionBegin; 45 if (!jac->diag) { 46 /* assume square matrices */ 47 ierr = MatGetVecs(pc->pmat,&jac->diag,&jac->work);CHKERRQ(ierr); 48 } 49 ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 50 ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr); 51 if (!jac->U) { 52 ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr); 53 ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->V);CHKERRQ(ierr); 54 } 55 ierr = MatGetSize(pc->pmat,&n,PETSC_NULL);CHKERRQ(ierr); 56 nb = PetscBLASIntCast(n); 57 lwork = 5*nb; 58 ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr); 59 ierr = MatGetArray(jac->A,&a);CHKERRQ(ierr); 60 ierr = MatGetArray(jac->U,&u);CHKERRQ(ierr); 61 ierr = MatGetArray(jac->V,&v);CHKERRQ(ierr); 62 ierr = VecGetArray(jac->diag,&d);CHKERRQ(ierr); 63 #if !defined(PETSC_USE_COMPLEX) 64 { 65 PetscBLASInt lierr; 66 LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr); 67 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 68 } 69 #else 70 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex"); 71 #endif 72 ierr = MatRestoreArray(jac->A,&a);CHKERRQ(ierr); 73 ierr = MatRestoreArray(jac->U,&u);CHKERRQ(ierr); 74 ierr = MatRestoreArray(jac->V,&v);CHKERRQ(ierr); 75 for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break; 76 jac->nzero = n-1-i; 77 if (jac->monitor) { 78 ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 79 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); 80 if (n >= 10) { /* print 5 smallest and 5 largest */ 81 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); 82 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); 83 } else { /* print all singular values */ 84 char buf[256],*p; 85 size_t left = sizeof buf,used; 86 PetscInt thisline; 87 for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) { 88 ierr = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr); 89 left -= used; 90 p += used; 91 if (thisline > 4 || i==0) { 92 ierr = PetscViewerASCIIPrintf(jac->monitor," SVD: singular values:%s\n",buf);CHKERRQ(ierr); 93 p = buf; 94 thisline = 0; 95 } 96 } 97 } 98 ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 99 } 100 ierr = PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1])); 101 for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i]; 102 for (; i<n; i++) d[i] = 0.0; 103 ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero); 104 ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr); 105 #if defined(foo) 106 { 107 PetscViewer viewer; 108 ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 109 ierr = MatView(jac->A,viewer);CHKERRQ(ierr); 110 ierr = MatView(jac->U,viewer);CHKERRQ(ierr); 111 ierr = MatView(jac->V,viewer);CHKERRQ(ierr); 112 ierr = VecView(jac->diag,viewer);CHKERRQ(ierr); 113 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 114 } 115 #endif 116 ierr = PetscFree(work); 117 PetscFunctionReturn(0); 118 #endif 119 } 120 121 /* -------------------------------------------------------------------------- */ 122 /* 123 PCApply_SVD - Applies the SVD preconditioner to a vector. 124 125 Input Parameters: 126 . pc - the preconditioner context 127 . x - input vector 128 129 Output Parameter: 130 . y - output vector 131 132 Application Interface Routine: PCApply() 133 */ 134 #undef __FUNCT__ 135 #define __FUNCT__ "PCApply_SVD" 136 static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y) 137 { 138 PC_SVD *jac = (PC_SVD*)pc->data; 139 Vec work = jac->work; 140 PetscErrorCode ierr; 141 142 PetscFunctionBegin; 143 ierr = MatMultTranspose(jac->U,x,work);CHKERRQ(ierr); 144 ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 145 ierr = MatMultTranspose(jac->V,work,y);CHKERRQ(ierr); 146 PetscFunctionReturn(0); 147 } 148 149 #undef __FUNCT__ 150 #define __FUNCT__ "PCReset_SVD" 151 static PetscErrorCode PCReset_SVD(PC pc) 152 { 153 PC_SVD *jac = (PC_SVD*)pc->data; 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 158 ierr = MatDestroy(&jac->U);CHKERRQ(ierr); 159 ierr = MatDestroy(&jac->V);CHKERRQ(ierr); 160 ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 161 ierr = VecDestroy(&jac->work);CHKERRQ(ierr); 162 PetscFunctionReturn(0); 163 } 164 165 /* -------------------------------------------------------------------------- */ 166 /* 167 PCDestroy_SVD - Destroys the private context for the SVD preconditioner 168 that was created with PCCreate_SVD(). 169 170 Input Parameter: 171 . pc - the preconditioner context 172 173 Application Interface Routine: PCDestroy() 174 */ 175 #undef __FUNCT__ 176 #define __FUNCT__ "PCDestroy_SVD" 177 static PetscErrorCode PCDestroy_SVD(PC pc) 178 { 179 PC_SVD *jac = (PC_SVD*)pc->data; 180 PetscErrorCode ierr; 181 182 PetscFunctionBegin; 183 ierr = PCReset_SVD(pc);CHKERRQ(ierr); 184 ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 185 ierr = PetscFree(pc->data);CHKERRQ(ierr); 186 PetscFunctionReturn(0); 187 } 188 189 #undef __FUNCT__ 190 #define __FUNCT__ "PCSetFromOptions_SVD" 191 static PetscErrorCode PCSetFromOptions_SVD(PC pc) 192 { 193 PetscErrorCode ierr; 194 PC_SVD *jac = (PC_SVD*)pc->data; 195 PetscBool flg,set; 196 197 PetscFunctionBegin; 198 ierr = PetscOptionsHead("SVD options");CHKERRQ(ierr); 199 ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,PETSC_NULL);CHKERRQ(ierr); 200 ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor?PETSC_TRUE:PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 201 if (set) { /* Should make PCSVDSetMonitor() */ 202 if (flg && !jac->monitor) { 203 ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,"stdout",&jac->monitor);CHKERRQ(ierr); 204 } else if (!flg) { 205 ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 206 } 207 } 208 ierr = PetscOptionsTail();CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 /* -------------------------------------------------------------------------- */ 213 /* 214 PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 215 and sets this as the private data within the generic preconditioning 216 context, PC, that was created within PCCreate(). 217 218 Input Parameter: 219 . pc - the preconditioner context 220 221 Application Interface Routine: PCCreate() 222 */ 223 224 /*MC 225 PCSVD - Use pseudo inverse defined by SVD of operator 226 227 Level: advanced 228 229 Concepts: SVD 230 231 Options Database: 232 - -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero 233 + -pc_svd_monitor Print information on the extreme singular values of the operator 234 235 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 236 M*/ 237 238 EXTERN_C_BEGIN 239 #undef __FUNCT__ 240 #define __FUNCT__ "PCCreate_SVD" 241 PetscErrorCode PCCreate_SVD(PC pc) 242 { 243 PC_SVD *jac; 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 /* 248 Creates the private data structure for this preconditioner and 249 attach it to the PC object. 250 */ 251 ierr = PetscNewLog(pc,PC_SVD,&jac);CHKERRQ(ierr); 252 jac->zerosing = 1.e-12; 253 pc->data = (void*)jac; 254 255 /* 256 Set the pointers for the functions that are provided above. 257 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 258 are called, they will automatically call these functions. Note we 259 choose not to provide a couple of these functions since they are 260 not needed. 261 */ 262 pc->ops->apply = PCApply_SVD; 263 pc->ops->applytranspose = PCApply_SVD; 264 pc->ops->setup = PCSetUp_SVD; 265 pc->ops->reset = PCReset_SVD; 266 pc->ops->destroy = PCDestroy_SVD; 267 pc->ops->setfromoptions = PCSetFromOptions_SVD; 268 pc->ops->view = 0; 269 pc->ops->applyrichardson = 0; 270 PetscFunctionReturn(0); 271 } 272 EXTERN_C_END 273 274