1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 2c6db04a5SJed Brown #include <petscblaslapack.h> 327c67122SBarry Smith 427c67122SBarry Smith /* 527c67122SBarry Smith Private context (data structure) for the SVD preconditioner. 627c67122SBarry Smith */ 727c67122SBarry Smith typedef struct { 827c67122SBarry Smith Vec diag, work; 93ed27f31SJed Brown Mat A, U, Vt; 1027c67122SBarry Smith PetscInt nzero; 118f1a2a5eSBarry Smith PetscReal zerosing; /* measure of smallest singular value treated as nonzero */ 1285032590SJed Brown PetscInt essrank; /* essential rank of operator */ 133ed27f31SJed Brown VecScatter left2red, right2red; 143ed27f31SJed Brown Vec leftred, rightred; 15426160bdSJed Brown PetscViewer monitor; 167962402dSFande Kong PetscViewerFormat monitorformat; 1727c67122SBarry Smith } PC_SVD; 1827c67122SBarry Smith 199371c9d4SSatish Balay typedef enum { 209371c9d4SSatish Balay READ = 1, 219371c9d4SSatish Balay WRITE = 2, 229371c9d4SSatish Balay READ_WRITE = 3 239371c9d4SSatish Balay } AccessMode; 2427c67122SBarry Smith 2527c67122SBarry Smith /* 2627c67122SBarry Smith PCSetUp_SVD - Prepares for the use of the SVD preconditioner 2727c67122SBarry Smith by setting data structures and options. 2827c67122SBarry Smith 2927c67122SBarry Smith Input Parameter: 3027c67122SBarry Smith . pc - the preconditioner context 3127c67122SBarry Smith 3227c67122SBarry Smith Application Interface Routine: PCSetUp() 3327c67122SBarry Smith 34f1580f4eSBarry Smith Note: 3527c67122SBarry Smith The interface routine PCSetUp() is not usually called directly by 3627c67122SBarry Smith the user, but instead is called by PCApply() if necessary. 3727c67122SBarry Smith */ 38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_SVD(PC pc) 39d71ae5a4SJacob Faibussowitsch { 4027c67122SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 4127c67122SBarry Smith PetscScalar *a, *u, *v, *d, *work; 423f83f0d9SMatthew G Knepley PetscBLASInt nb, lwork; 4327c67122SBarry Smith PetscInt i, n; 443ed27f31SJed Brown PetscMPIInt size; 4527c67122SBarry Smith 4627c67122SBarry Smith PetscFunctionBegin; 479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->A)); 489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(((PetscObject)pc->pmat)->comm, &size)); 493ed27f31SJed Brown if (size > 1) { 503ed27f31SJed Brown Mat redmat; 5128d58a37SPierre Jolivet 529566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(pc->pmat, size, PETSC_COMM_SELF, MAT_INITIAL_MATRIX, &redmat)); 539566063dSJacob Faibussowitsch PetscCall(MatConvert(redmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A)); 549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&redmat)); 553ed27f31SJed Brown } else { 569566063dSJacob Faibussowitsch PetscCall(MatConvert(pc->pmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A)); 573ed27f31SJed Brown } 583ed27f31SJed Brown if (!jac->diag) { /* assume square matrices */ 599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(jac->A, &jac->diag, &jac->work)); 603ed27f31SJed Brown } 6127c67122SBarry Smith if (!jac->U) { 629566063dSJacob Faibussowitsch PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->U)); 639566063dSJacob Faibussowitsch PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->Vt)); 6427c67122SBarry Smith } 659566063dSJacob Faibussowitsch PetscCall(MatGetSize(jac->A, &n, NULL)); 66bf8ddafcSStefano Zampini if (!n) PetscFunctionReturn(PETSC_SUCCESS); 67bf8ddafcSStefano Zampini 689566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &nb)); 6927c67122SBarry Smith lwork = 5 * nb; 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lwork, &work)); 719566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->A, &a)); 729566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->U, &u)); 739566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->Vt, &v)); 749566063dSJacob Faibussowitsch PetscCall(VecGetArray(jac->diag, &d)); 7527c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX) 763f83f0d9SMatthew G Knepley { 773f83f0d9SMatthew G Knepley PetscBLASInt lierr; 789566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 79792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, d, u, &nb, v, &nb, work, &lwork, &lierr)); 80c9c947b5SJed Brown PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesvd() error %" PetscBLASInt_FMT, lierr); 819566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 823f83f0d9SMatthew G Knepley } 8327c67122SBarry Smith #else 84ef47b4b1SBarry Smith { 85ef47b4b1SBarry Smith PetscBLASInt lierr; 86ef47b4b1SBarry Smith PetscReal *rwork, *dd; 879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * nb, &rwork)); 889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nb, &dd)); 899566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 90792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, dd, u, &nb, v, &nb, work, &lwork, rwork, &lierr)); 9163a3b9bcSJacob Faibussowitsch PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr); 929566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 9314ce09a1SBarry Smith for (i = 0; i < n; i++) d[i] = dd[i]; 949566063dSJacob Faibussowitsch PetscCall(PetscFree(dd)); 959566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 96ef47b4b1SBarry Smith } 9727c67122SBarry Smith #endif 989566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->A, &a)); 999566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->U, &u)); 1009566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->Vt, &v)); 1019371c9d4SSatish Balay for (i = n - 1; i >= 0; i--) 1029371c9d4SSatish Balay if (PetscRealPart(d[i]) > jac->zerosing) break; 103426160bdSJed Brown jac->nzero = n - 1 - i; 104426160bdSJed Brown if (jac->monitor) { 1059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(jac->monitor, ((PetscObject)pc)->tablevel)); 10663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(jac->monitor, " SVD: condition number %14.12e, %" PetscInt_FMT " of %" PetscInt_FMT " singular values are (nearly) zero\n", (double)PetscRealPart(d[0] / d[n - 1]), jac->nzero, n)); 1077962402dSFande Kong if (n < 10 || jac->monitorformat == PETSC_VIEWER_ALL) { 1087962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " SVD: singular values:\n")); 1097962402dSFande Kong for (i = 0; i < n; i++) { 1107962402dSFande Kong if (i % 5 == 0) { 11148a46eb9SPierre Jolivet if (i != 0) PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n")); 1127962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " ")); 1137962402dSFande Kong } 1147962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " %14.12e", (double)PetscRealPart(d[i]))); 1157962402dSFande Kong } 1167962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n")); 1177962402dSFande Kong } else { /* print 5 smallest and 5 largest */ 1189566063dSJacob Faibussowitsch 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]))); 1199566063dSJacob Faibussowitsch 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]))); 120426160bdSJed Brown } 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(jac->monitor, ((PetscObject)pc)->tablevel)); 122426160bdSJed Brown } 1239566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Largest and smallest singular values %14.12e %14.12e\n", (double)PetscRealPart(d[0]), (double)PetscRealPart(d[n - 1]))); 124426160bdSJed Brown for (i = 0; i < n - jac->nzero; i++) d[i] = 1.0 / d[i]; 125426160bdSJed Brown for (; i < n; i++) d[i] = 0.0; 1269371c9d4SSatish Balay if (jac->essrank > 0) 1279371c9d4SSatish Balay for (i = 0; i < n - jac->nzero - jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */ 12863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Number of zero or nearly singular values %" PetscInt_FMT "\n", jac->nzero)); 1299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(jac->diag, &d)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13227c67122SBarry Smith } 13327c67122SBarry Smith 134d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSVDGetVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) 135d71ae5a4SJacob Faibussowitsch { 1363ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 1373ed27f31SJed Brown PetscMPIInt size; 1383ed27f31SJed Brown 1393ed27f31SJed Brown PetscFunctionBegin; 1409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1410298fd71SBarry Smith *xred = NULL; 1423ed27f31SJed Brown switch (side) { 1433ed27f31SJed Brown case PC_LEFT: 1443ed27f31SJed Brown if (size == 1) *xred = x; 1453ed27f31SJed Brown else { 1469566063dSJacob Faibussowitsch if (!jac->left2red) PetscCall(VecScatterCreateToAll(x, &jac->left2red, &jac->leftred)); 1473ed27f31SJed Brown if (amode & READ) { 1489566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD)); 1499566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD)); 1503ed27f31SJed Brown } 1513ed27f31SJed Brown *xred = jac->leftred; 1523ed27f31SJed Brown } 1533ed27f31SJed Brown break; 1543ed27f31SJed Brown case PC_RIGHT: 1553ed27f31SJed Brown if (size == 1) *xred = x; 1563ed27f31SJed Brown else { 1579566063dSJacob Faibussowitsch if (!jac->right2red) PetscCall(VecScatterCreateToAll(x, &jac->right2red, &jac->rightred)); 1583ed27f31SJed Brown if (amode & READ) { 1599566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD)); 1609566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD)); 1613ed27f31SJed Brown } 1623ed27f31SJed Brown *xred = jac->rightred; 1633ed27f31SJed Brown } 1643ed27f31SJed Brown break; 165d71ae5a4SJacob Faibussowitsch default: 166d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT"); 1673ed27f31SJed Brown } 1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1693ed27f31SJed Brown } 1703ed27f31SJed Brown 171d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSVDRestoreVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) 172d71ae5a4SJacob Faibussowitsch { 1733ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 1743ed27f31SJed Brown PetscMPIInt size; 1753ed27f31SJed Brown 1763ed27f31SJed Brown PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1783ed27f31SJed Brown switch (side) { 1793ed27f31SJed Brown case PC_LEFT: 1803ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1819566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE)); 1829566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE)); 1833ed27f31SJed Brown } 1843ed27f31SJed Brown break; 1853ed27f31SJed Brown case PC_RIGHT: 1863ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1879566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE)); 1889566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE)); 1893ed27f31SJed Brown } 1903ed27f31SJed Brown break; 191d71ae5a4SJacob Faibussowitsch default: 192d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT"); 1933ed27f31SJed Brown } 1940298fd71SBarry Smith *xred = NULL; 1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1963ed27f31SJed Brown } 1973ed27f31SJed Brown 19827c67122SBarry Smith /* 19927c67122SBarry Smith PCApply_SVD - Applies the SVD preconditioner to a vector. 20027c67122SBarry Smith 20127c67122SBarry Smith Input Parameters: 20227c67122SBarry Smith . pc - the preconditioner context 20327c67122SBarry Smith . x - input vector 20427c67122SBarry Smith 20527c67122SBarry Smith Output Parameter: 20627c67122SBarry Smith . y - output vector 20727c67122SBarry Smith 20827c67122SBarry Smith Application Interface Routine: PCApply() 20927c67122SBarry Smith */ 210d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SVD(PC pc, Vec x, Vec y) 211d71ae5a4SJacob Faibussowitsch { 21227c67122SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 2133ed27f31SJed Brown Vec work = jac->work, xred, yred; 21427c67122SBarry Smith 21527c67122SBarry Smith PetscFunctionBegin; 2169566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_RIGHT, READ, x, &xred)); 2179566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_LEFT, WRITE, y, &yred)); 218ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2199566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac->U, xred, work)); 220ef47b4b1SBarry Smith #else 2219566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->U, xred, work)); 222ef47b4b1SBarry Smith #endif 2239566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(work, work, jac->diag)); 224ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2259566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac->Vt, work, yred)); 226ef47b4b1SBarry Smith #else 2279566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->Vt, work, yred)); 228ef47b4b1SBarry Smith #endif 2299566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, READ, x, &xred)); 2309566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_LEFT, WRITE, y, &yred)); 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2323ed27f31SJed Brown } 2333ed27f31SJed Brown 234d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_SVD(PC pc, Mat X, Mat Y) 235d71ae5a4SJacob Faibussowitsch { 236024db992SStefano Zampini PC_SVD *jac = (PC_SVD *)pc->data; 237024db992SStefano Zampini Mat W; 238024db992SStefano Zampini 239024db992SStefano Zampini PetscFunctionBegin; 240fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(jac->U, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &W)); 241024db992SStefano Zampini PetscCall(MatDiagonalScale(W, jac->diag, NULL)); 242fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(jac->Vt, W, MAT_REUSE_MATRIX, PETSC_DETERMINE, &Y)); 243024db992SStefano Zampini PetscCall(MatDestroy(&W)); 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 245024db992SStefano Zampini } 246024db992SStefano Zampini 247d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_SVD(PC pc, Vec x, Vec y) 248d71ae5a4SJacob Faibussowitsch { 2493ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 2503ed27f31SJed Brown Vec work = jac->work, xred, yred; 2513ed27f31SJed Brown 2523ed27f31SJed Brown PetscFunctionBegin; 2539566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_LEFT, READ, x, &xred)); 2549566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_RIGHT, WRITE, y, &yred)); 2559566063dSJacob Faibussowitsch PetscCall(MatMult(jac->Vt, xred, work)); 2569566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(work, work, jac->diag)); 2579566063dSJacob Faibussowitsch PetscCall(MatMult(jac->U, work, yred)); 2589566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_LEFT, READ, x, &xred)); 2599566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, WRITE, y, &yred)); 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26127c67122SBarry Smith } 26227c67122SBarry Smith 263d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_SVD(PC pc) 264d71ae5a4SJacob Faibussowitsch { 265a2d70de2SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 266a2d70de2SBarry Smith 267a2d70de2SBarry Smith PetscFunctionBegin; 2689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->A)); 2699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->U)); 2709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->Vt)); 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diag)); 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->work)); 2739566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&jac->right2red)); 2749566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&jac->left2red)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->rightred)); 2769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->leftred)); 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 278a2d70de2SBarry Smith } 279a2d70de2SBarry Smith 28027c67122SBarry Smith /* 28127c67122SBarry Smith PCDestroy_SVD - Destroys the private context for the SVD preconditioner 28227c67122SBarry Smith that was created with PCCreate_SVD(). 28327c67122SBarry Smith 28427c67122SBarry Smith Input Parameter: 28527c67122SBarry Smith . pc - the preconditioner context 28627c67122SBarry Smith 28727c67122SBarry Smith Application Interface Routine: PCDestroy() 28827c67122SBarry Smith */ 289d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SVD(PC pc) 290d71ae5a4SJacob Faibussowitsch { 291426160bdSJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 29227c67122SBarry Smith 29327c67122SBarry Smith PetscFunctionBegin; 2949566063dSJacob Faibussowitsch PetscCall(PCReset_SVD(pc)); 295648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&jac->monitor)); 2969566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29827c67122SBarry Smith } 29927c67122SBarry Smith 300ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_SVD(PC pc, PetscOptionItems PetscOptionsObject) 301d71ae5a4SJacob Faibussowitsch { 3028f1a2a5eSBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 3037962402dSFande Kong PetscBool flg; 30427c67122SBarry Smith 30527c67122SBarry Smith PetscFunctionBegin; 306d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SVD options"); 3079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_svd_zero_sing", "Singular values smaller than this treated as zero", "None", jac->zerosing, &jac->zerosing, NULL)); 3089566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_svd_ess_rank", "Essential rank of operator (0 to use entire operator)", "None", jac->essrank, &jac->essrank, NULL)); 3097962402dSFande Kong PetscCall(PetscOptionsViewer("-pc_svd_monitor", "Monitor the conditioning, and extremal singular values", "None", &jac->monitor, &jac->monitorformat, &flg)); 310d0609cedSBarry Smith PetscOptionsHeadEnd(); 3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31227c67122SBarry Smith } 31327c67122SBarry Smith 314d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_SVD(PC pc, PetscViewer viewer) 315d71ae5a4SJacob Faibussowitsch { 31633761216SBarry Smith PC_SVD *svd = (PC_SVD *)pc->data; 317*9f196a02SMartin Diehl PetscBool isascii; 31833761216SBarry Smith 31933761216SBarry Smith PetscFunctionBegin; 320*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 321*9f196a02SMartin Diehl if (isascii) { 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " All singular values smaller than %g treated as zero\n", (double)svd->zerosing)); 32363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Provided essential rank of the matrix %" PetscInt_FMT " (all other eigenvalues are zeroed)\n", svd->essrank)); 32433761216SBarry Smith } 3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32633761216SBarry Smith } 327f1580f4eSBarry Smith 32827c67122SBarry Smith /* 32927c67122SBarry Smith PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 33027c67122SBarry Smith and sets this as the private data within the generic preconditioning 33127c67122SBarry Smith context, PC, that was created within PCCreate(). 33227c67122SBarry Smith 33327c67122SBarry Smith Input Parameter: 33427c67122SBarry Smith . pc - the preconditioner context 33527c67122SBarry Smith 33627c67122SBarry Smith Application Interface Routine: PCCreate() 33727c67122SBarry Smith */ 33827c67122SBarry Smith 33927c67122SBarry Smith /*MC 34027c67122SBarry Smith PCSVD - Use pseudo inverse defined by SVD of operator 34127c67122SBarry Smith 34227c67122SBarry Smith Level: advanced 34327c67122SBarry Smith 344f1580f4eSBarry Smith Options Database Keys: 345147403d9SBarry Smith + -pc_svd_zero_sing <rtol> - Singular values smaller than this are treated as zero 346147403d9SBarry Smith - -pc_svd_monitor - Print information on the extreme singular values of the operator 34727c67122SBarry Smith 348a2b725a8SWilliam Gropp Developer Note: 349a2b725a8SWilliam Gropp This implementation automatically creates a redundant copy of the 3508997ae2eSBarry Smith matrix on each process and uses a sequential SVD solve. Why does it do this instead 351f1580f4eSBarry Smith of using the composable `PCREDUNDANT` object? 3528997ae2eSBarry Smith 353562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCREDUNDANT` 35427c67122SBarry Smith M*/ 35527c67122SBarry Smith 356d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc) 357d71ae5a4SJacob Faibussowitsch { 35827c67122SBarry Smith PC_SVD *jac; 359024db992SStefano Zampini PetscMPIInt size = 0; 36027c67122SBarry Smith 36127c67122SBarry Smith PetscFunctionBegin; 36227c67122SBarry Smith /* 36327c67122SBarry Smith Creates the private data structure for this preconditioner and 36427c67122SBarry Smith attach it to the PC object. 36527c67122SBarry Smith */ 3664dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 3678f1a2a5eSBarry Smith jac->zerosing = 1.e-12; 36827c67122SBarry Smith pc->data = (void *)jac; 36927c67122SBarry Smith 37027c67122SBarry Smith /* 37127c67122SBarry Smith Set the pointers for the functions that are provided above. 37227c67122SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 37327c67122SBarry Smith are called, they will automatically call these functions. Note we 37427c67122SBarry Smith choose not to provide a couple of these functions since they are 37527c67122SBarry Smith not needed. 37627c67122SBarry Smith */ 377024db992SStefano Zampini 378024db992SStefano Zampini #if defined(PETSC_HAVE_COMPLEX) 379024db992SStefano Zampini PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 380024db992SStefano Zampini #endif 381024db992SStefano Zampini if (size == 1) pc->ops->matapply = PCMatApply_SVD; 38227c67122SBarry Smith pc->ops->apply = PCApply_SVD; 3833ed27f31SJed Brown pc->ops->applytranspose = PCApplyTranspose_SVD; 38427c67122SBarry Smith pc->ops->setup = PCSetUp_SVD; 385a2d70de2SBarry Smith pc->ops->reset = PCReset_SVD; 38627c67122SBarry Smith pc->ops->destroy = PCDestroy_SVD; 38727c67122SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SVD; 38833761216SBarry Smith pc->ops->view = PCView_SVD; 3890a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39127c67122SBarry Smith } 392