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