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