1 2 #include <petsc-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,Vt; 11 PetscInt nzero; 12 PetscReal zerosing; /* measure of smallest singular value treated as nonzero */ 13 PetscInt essrank; /* essential rank of operator */ 14 VecScatter left2red,right2red; 15 Vec leftred,rightred; 16 PetscViewer monitor; 17 } PC_SVD; 18 19 typedef enum {READ=1, WRITE=2, READ_WRITE=3} AccessMode; 20 21 /* -------------------------------------------------------------------------- */ 22 /* 23 PCSetUp_SVD - Prepares for the use of the SVD preconditioner 24 by setting data structures and options. 25 26 Input Parameter: 27 . pc - the preconditioner context 28 29 Application Interface Routine: PCSetUp() 30 31 Notes: 32 The interface routine PCSetUp() is not usually called directly by 33 the user, but instead is called by PCApply() if necessary. 34 */ 35 #undef __FUNCT__ 36 #define __FUNCT__ "PCSetUp_SVD" 37 static PetscErrorCode PCSetUp_SVD(PC pc) 38 { 39 #if defined(PETSC_MISSING_LAPACK_GESVD) 40 SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"GESVD - Lapack routine is unavailable\nNot able to provide singular value estimates."); 41 #else 42 PC_SVD *jac = (PC_SVD*)pc->data; 43 PetscErrorCode ierr; 44 PetscScalar *a,*u,*v,*d,*work; 45 PetscBLASInt nb,lwork; 46 PetscInt i,n; 47 PetscMPIInt size; 48 49 PetscFunctionBegin; 50 ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 51 ierr = MPI_Comm_size(((PetscObject)pc->pmat)->comm,&size);CHKERRQ(ierr); 52 if (size > 1) { 53 Mat redmat; 54 PetscInt M; 55 ierr = MatGetSize(pc->pmat,&M,PETSC_NULL);CHKERRQ(ierr); 56 ierr = MatGetRedundantMatrix(pc->pmat,size,PETSC_COMM_SELF,M,MAT_INITIAL_MATRIX,&redmat);CHKERRQ(ierr); 57 ierr = MatConvert(redmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr); 58 ierr = MatDestroy(&redmat);CHKERRQ(ierr); 59 } else { 60 ierr = MatConvert(pc->pmat,MATSEQDENSE,MAT_INITIAL_MATRIX,&jac->A);CHKERRQ(ierr); 61 } 62 if (!jac->diag) { /* assume square matrices */ 63 ierr = MatGetVecs(jac->A,&jac->diag,&jac->work);CHKERRQ(ierr); 64 } 65 if (!jac->U) { 66 ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->U);CHKERRQ(ierr); 67 ierr = MatDuplicate(jac->A,MAT_DO_NOT_COPY_VALUES,&jac->Vt);CHKERRQ(ierr); 68 } 69 ierr = MatGetSize(pc->pmat,&n,PETSC_NULL);CHKERRQ(ierr); 70 ierr = PetscBLASIntCast(n,&nb);CHKERRQ(ierr); 71 lwork = 5*nb; 72 ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr); 73 ierr = MatDenseGetArray(jac->A,&a);CHKERRQ(ierr); 74 ierr = MatDenseGetArray(jac->U,&u);CHKERRQ(ierr); 75 ierr = MatDenseGetArray(jac->Vt,&v);CHKERRQ(ierr); 76 ierr = VecGetArray(jac->diag,&d);CHKERRQ(ierr); 77 #if !defined(PETSC_USE_COMPLEX) 78 { 79 PetscBLASInt lierr; 80 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 81 LAPACKgesvd_("A","A",&nb,&nb,a,&nb,d,u,&nb,v,&nb,work,&lwork,&lierr); 82 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"gesv() error %d",lierr); 83 ierr = PetscFPTrapPop();CHKERRQ(ierr); 84 } 85 #else 86 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex"); 87 #endif 88 ierr = MatDenseRestoreArray(jac->A,&a);CHKERRQ(ierr); 89 ierr = MatDenseRestoreArray(jac->U,&u);CHKERRQ(ierr); 90 ierr = MatDenseRestoreArray(jac->Vt,&v);CHKERRQ(ierr); 91 for (i=n-1; i>=0; i--) if (PetscRealPart(d[i]) > jac->zerosing) break; 92 jac->nzero = n-1-i; 93 if (jac->monitor) { 94 ierr = PetscViewerASCIIAddTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 95 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); 96 if (n >= 10) { /* print 5 smallest and 5 largest */ 97 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); 98 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); 99 } else { /* print all singular values */ 100 char buf[256],*p; 101 size_t left = sizeof(buf),used; 102 PetscInt thisline; 103 for (p=buf,i=n-1,thisline=1; i>=0; i--,thisline++) { 104 ierr = PetscSNPrintfCount(p,left," %14.12e",&used,(double)PetscRealPart(d[i]));CHKERRQ(ierr); 105 left -= used; 106 p += used; 107 if (thisline > 4 || i==0) { 108 ierr = PetscViewerASCIIPrintf(jac->monitor," SVD: singular values:%s\n",buf);CHKERRQ(ierr); 109 p = buf; 110 thisline = 0; 111 } 112 } 113 } 114 ierr = PetscViewerASCIISubtractTab(jac->monitor,((PetscObject)pc)->tablevel);CHKERRQ(ierr); 115 } 116 ierr = PetscInfo2(pc,"Largest and smallest singular values %14.12e %14.12e\n",(double)PetscRealPart(d[0]),(double)PetscRealPart(d[n-1]));CHKERRQ(ierr); 117 for (i=0; i<n-jac->nzero; i++) d[i] = 1.0/d[i]; 118 for (; i<n; i++) d[i] = 0.0; 119 if (jac->essrank > 0) for (i=0; i<n-jac->nzero-jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */ 120 ierr = PetscInfo1(pc,"Number of zero or nearly singular values %D\n",jac->nzero);CHKERRQ(ierr); 121 ierr = VecRestoreArray(jac->diag,&d);CHKERRQ(ierr); 122 #if defined(foo) 123 { 124 PetscViewer viewer; 125 ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"joe",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 126 ierr = MatView(jac->A,viewer);CHKERRQ(ierr); 127 ierr = MatView(jac->U,viewer);CHKERRQ(ierr); 128 ierr = MatView(jac->Vt,viewer);CHKERRQ(ierr); 129 ierr = VecView(jac->diag,viewer);CHKERRQ(ierr); 130 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 131 } 132 #endif 133 ierr = PetscFree(work);CHKERRQ(ierr); 134 PetscFunctionReturn(0); 135 #endif 136 } 137 138 #undef __FUNCT__ 139 #define __FUNCT__ "PCSVDGetVec" 140 static PetscErrorCode PCSVDGetVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 141 { 142 PC_SVD *jac = (PC_SVD*)pc->data; 143 PetscErrorCode ierr; 144 PetscMPIInt size; 145 146 PetscFunctionBegin; 147 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 148 *xred = PETSC_NULL; 149 switch (side) { 150 case PC_LEFT: 151 if (size == 1) *xred = x; 152 else { 153 if (!jac->left2red) {ierr = VecScatterCreateToAll(x,&jac->left2red,&jac->leftred);CHKERRQ(ierr);} 154 if (amode & READ) { 155 ierr = VecScatterBegin(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 156 ierr = VecScatterEnd(jac->left2red,x,jac->leftred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 157 } 158 *xred = jac->leftred; 159 } 160 break; 161 case PC_RIGHT: 162 if (size == 1) *xred = x; 163 else { 164 if (!jac->right2red) {ierr = VecScatterCreateToAll(x,&jac->right2red,&jac->rightred);CHKERRQ(ierr);} 165 if (amode & READ) { 166 ierr = VecScatterBegin(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 167 ierr = VecScatterEnd(jac->right2red,x,jac->rightred,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 168 } 169 *xred = jac->rightred; 170 } 171 break; 172 default: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 173 } 174 PetscFunctionReturn(0); 175 } 176 177 #undef __FUNCT__ 178 #define __FUNCT__ "PCSVDRestoreVec" 179 static PetscErrorCode PCSVDRestoreVec(PC pc,PCSide side,AccessMode amode,Vec x,Vec *xred) 180 { 181 PC_SVD *jac = (PC_SVD*)pc->data; 182 PetscErrorCode ierr; 183 PetscMPIInt size; 184 185 PetscFunctionBegin; 186 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 187 switch (side) { 188 case PC_LEFT: 189 if (size != 1 && amode & WRITE) { 190 ierr = VecScatterBegin(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 191 ierr = VecScatterEnd(jac->left2red,jac->leftred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 192 } 193 break; 194 case PC_RIGHT: 195 if (size != 1 && amode & WRITE) { 196 ierr = VecScatterBegin(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 197 ierr = VecScatterEnd(jac->right2red,jac->rightred,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 198 } 199 break; 200 default: SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"Side must be LEFT or RIGHT"); 201 } 202 *xred = PETSC_NULL; 203 PetscFunctionReturn(0); 204 } 205 206 /* -------------------------------------------------------------------------- */ 207 /* 208 PCApply_SVD - Applies the SVD preconditioner to a vector. 209 210 Input Parameters: 211 . pc - the preconditioner context 212 . x - input vector 213 214 Output Parameter: 215 . y - output vector 216 217 Application Interface Routine: PCApply() 218 */ 219 #undef __FUNCT__ 220 #define __FUNCT__ "PCApply_SVD" 221 static PetscErrorCode PCApply_SVD(PC pc,Vec x,Vec y) 222 { 223 PC_SVD *jac = (PC_SVD*)pc->data; 224 Vec work = jac->work,xred,yred; 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 ierr = PCSVDGetVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr); 229 ierr = PCSVDGetVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr); 230 ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr); 231 ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 232 ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr); 233 ierr = PCSVDRestoreVec(pc,PC_RIGHT,READ,x,&xred);CHKERRQ(ierr); 234 ierr = PCSVDRestoreVec(pc,PC_LEFT,WRITE,y,&yred);CHKERRQ(ierr); 235 PetscFunctionReturn(0); 236 } 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "PCApplyTranspose_SVD" 240 static PetscErrorCode PCApplyTranspose_SVD(PC pc,Vec x,Vec y) 241 { 242 PC_SVD *jac = (PC_SVD*)pc->data; 243 Vec work = jac->work,xred,yred; 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 ierr = PCSVDGetVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr); 248 ierr = PCSVDGetVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr); 249 ierr = MatMultTranspose(jac->Vt,work,yred);CHKERRQ(ierr); 250 ierr = VecPointwiseMult(work,work,jac->diag);CHKERRQ(ierr); 251 ierr = MatMultTranspose(jac->U,xred,work);CHKERRQ(ierr); 252 ierr = PCSVDRestoreVec(pc,PC_LEFT,READ,x,&xred);CHKERRQ(ierr); 253 ierr = PCSVDRestoreVec(pc,PC_RIGHT,WRITE,y,&yred);CHKERRQ(ierr); 254 PetscFunctionReturn(0); 255 } 256 257 #undef __FUNCT__ 258 #define __FUNCT__ "PCReset_SVD" 259 static PetscErrorCode PCReset_SVD(PC pc) 260 { 261 PC_SVD *jac = (PC_SVD*)pc->data; 262 PetscErrorCode ierr; 263 264 PetscFunctionBegin; 265 ierr = MatDestroy(&jac->A);CHKERRQ(ierr); 266 ierr = MatDestroy(&jac->U);CHKERRQ(ierr); 267 ierr = MatDestroy(&jac->Vt);CHKERRQ(ierr); 268 ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 269 ierr = VecDestroy(&jac->work);CHKERRQ(ierr); 270 ierr = VecScatterDestroy(&jac->right2red);CHKERRQ(ierr); 271 ierr = VecScatterDestroy(&jac->left2red);CHKERRQ(ierr); 272 ierr = VecDestroy(&jac->rightred);CHKERRQ(ierr); 273 ierr = VecDestroy(&jac->leftred);CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 277 /* -------------------------------------------------------------------------- */ 278 /* 279 PCDestroy_SVD - Destroys the private context for the SVD preconditioner 280 that was created with PCCreate_SVD(). 281 282 Input Parameter: 283 . pc - the preconditioner context 284 285 Application Interface Routine: PCDestroy() 286 */ 287 #undef __FUNCT__ 288 #define __FUNCT__ "PCDestroy_SVD" 289 static PetscErrorCode PCDestroy_SVD(PC pc) 290 { 291 PC_SVD *jac = (PC_SVD*)pc->data; 292 PetscErrorCode ierr; 293 294 PetscFunctionBegin; 295 ierr = PCReset_SVD(pc);CHKERRQ(ierr); 296 ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 297 ierr = PetscFree(pc->data);CHKERRQ(ierr); 298 PetscFunctionReturn(0); 299 } 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "PCSetFromOptions_SVD" 303 static PetscErrorCode PCSetFromOptions_SVD(PC pc) 304 { 305 PetscErrorCode ierr; 306 PC_SVD *jac = (PC_SVD*)pc->data; 307 PetscBool flg,set; 308 309 PetscFunctionBegin; 310 ierr = PetscOptionsHead("SVD options");CHKERRQ(ierr); 311 ierr = PetscOptionsReal("-pc_svd_zero_sing","Singular values smaller than this treated as zero","None",jac->zerosing,&jac->zerosing,PETSC_NULL);CHKERRQ(ierr); 312 ierr = PetscOptionsInt("-pc_svd_ess_rank","Essential rank of operator (0 to use entire operator)","None",jac->essrank,&jac->essrank,PETSC_NULL);CHKERRQ(ierr); 313 ierr = PetscOptionsBool("-pc_svd_monitor","Monitor the conditioning, and extremal singular values","None",jac->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 314 if (set) { /* Should make PCSVDSetMonitor() */ 315 if (flg && !jac->monitor) { 316 ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,"stdout",&jac->monitor);CHKERRQ(ierr); 317 } else if (!flg) { 318 ierr = PetscViewerDestroy(&jac->monitor);CHKERRQ(ierr); 319 } 320 } 321 ierr = PetscOptionsTail();CHKERRQ(ierr); 322 PetscFunctionReturn(0); 323 } 324 325 /* -------------------------------------------------------------------------- */ 326 /* 327 PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 328 and sets this as the private data within the generic preconditioning 329 context, PC, that was created within PCCreate(). 330 331 Input Parameter: 332 . pc - the preconditioner context 333 334 Application Interface Routine: PCCreate() 335 */ 336 337 /*MC 338 PCSVD - Use pseudo inverse defined by SVD of operator 339 340 Level: advanced 341 342 Concepts: SVD 343 344 Options Database: 345 - -pc_svd_zero_sing <rtol> Singular values smaller than this are treated as zero 346 + -pc_svd_monitor Print information on the extreme singular values of the operator 347 348 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 349 M*/ 350 351 EXTERN_C_BEGIN 352 #undef __FUNCT__ 353 #define __FUNCT__ "PCCreate_SVD" 354 PetscErrorCode PCCreate_SVD(PC pc) 355 { 356 PC_SVD *jac; 357 PetscErrorCode ierr; 358 359 PetscFunctionBegin; 360 /* 361 Creates the private data structure for this preconditioner and 362 attach it to the PC object. 363 */ 364 ierr = PetscNewLog(pc,PC_SVD,&jac);CHKERRQ(ierr); 365 jac->zerosing = 1.e-12; 366 pc->data = (void*)jac; 367 368 /* 369 Set the pointers for the functions that are provided above. 370 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 371 are called, they will automatically call these functions. Note we 372 choose not to provide a couple of these functions since they are 373 not needed. 374 */ 375 pc->ops->apply = PCApply_SVD; 376 pc->ops->applytranspose = PCApplyTranspose_SVD; 377 pc->ops->setup = PCSetUp_SVD; 378 pc->ops->reset = PCReset_SVD; 379 pc->ops->destroy = PCDestroy_SVD; 380 pc->ops->setfromoptions = PCSetFromOptions_SVD; 381 pc->ops->view = 0; 382 pc->ops->applyrichardson = 0; 383 PetscFunctionReturn(0); 384 } 385 EXTERN_C_END 386 387