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