127c67122SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 3c6db04a5SJed Brown #include <petscblaslapack.h> 427c67122SBarry Smith 527c67122SBarry Smith /* 627c67122SBarry Smith Private context (data structure) for the SVD preconditioner. 727c67122SBarry Smith */ 827c67122SBarry Smith typedef struct { 927c67122SBarry Smith Vec diag, work; 103ed27f31SJed Brown Mat A, U, Vt; 1127c67122SBarry Smith PetscInt nzero; 128f1a2a5eSBarry Smith PetscReal zerosing; /* measure of smallest singular value treated as nonzero */ 1385032590SJed Brown PetscInt essrank; /* essential rank of operator */ 143ed27f31SJed Brown VecScatter left2red, right2red; 153ed27f31SJed Brown Vec leftred, rightred; 16426160bdSJed Brown PetscViewer monitor; 177962402dSFande Kong PetscViewerFormat monitorformat; 1827c67122SBarry Smith } PC_SVD; 1927c67122SBarry Smith 209371c9d4SSatish Balay typedef enum { 219371c9d4SSatish Balay READ = 1, 229371c9d4SSatish Balay WRITE = 2, 239371c9d4SSatish Balay READ_WRITE = 3 249371c9d4SSatish Balay } AccessMode; 2527c67122SBarry Smith 2627c67122SBarry Smith /* 2727c67122SBarry Smith PCSetUp_SVD - Prepares for the use of the SVD preconditioner 2827c67122SBarry Smith by setting data structures and options. 2927c67122SBarry Smith 3027c67122SBarry Smith Input Parameter: 3127c67122SBarry Smith . pc - the preconditioner context 3227c67122SBarry Smith 3327c67122SBarry Smith Application Interface Routine: PCSetUp() 3427c67122SBarry Smith 35*f1580f4eSBarry Smith Note: 3627c67122SBarry Smith The interface routine PCSetUp() is not usually called directly by 3727c67122SBarry Smith the user, but instead is called by PCApply() if necessary. 3827c67122SBarry Smith */ 399371c9d4SSatish Balay static PetscErrorCode PCSetUp_SVD(PC pc) { 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)); 6614ce09a1SBarry Smith if (!n) { 679566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Matrix has zero rows, skipping svd\n")); 6814ce09a1SBarry Smith PetscFunctionReturn(0); 6914ce09a1SBarry Smith } 709566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &nb)); 7127c67122SBarry Smith lwork = 5 * nb; 729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lwork, &work)); 739566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->A, &a)); 749566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->U, &u)); 759566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->Vt, &v)); 769566063dSJacob Faibussowitsch PetscCall(VecGetArray(jac->diag, &d)); 7727c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX) 783f83f0d9SMatthew G Knepley { 793f83f0d9SMatthew G Knepley PetscBLASInt lierr; 809566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 81792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, d, u, &nb, v, &nb, work, &lwork, &lierr)); 8263a3b9bcSJacob Faibussowitsch PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr); 839566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 843f83f0d9SMatthew G Knepley } 8527c67122SBarry Smith #else 86ef47b4b1SBarry Smith { 87ef47b4b1SBarry Smith PetscBLASInt lierr; 88ef47b4b1SBarry Smith PetscReal *rwork, *dd; 899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * nb, &rwork)); 909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nb, &dd)); 919566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 92792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, dd, u, &nb, v, &nb, work, &lwork, rwork, &lierr)); 9363a3b9bcSJacob Faibussowitsch PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr); 949566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 9514ce09a1SBarry Smith for (i = 0; i < n; i++) d[i] = dd[i]; 969566063dSJacob Faibussowitsch PetscCall(PetscFree(dd)); 979566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 98ef47b4b1SBarry Smith } 9927c67122SBarry Smith #endif 1009566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->A, &a)); 1019566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->U, &u)); 1029566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->Vt, &v)); 1039371c9d4SSatish Balay for (i = n - 1; i >= 0; i--) 1049371c9d4SSatish Balay if (PetscRealPart(d[i]) > jac->zerosing) break; 105426160bdSJed Brown jac->nzero = n - 1 - i; 106426160bdSJed Brown if (jac->monitor) { 1079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(jac->monitor, ((PetscObject)pc)->tablevel)); 10863a3b9bcSJacob 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)); 1097962402dSFande Kong if (n < 10 || jac->monitorformat == PETSC_VIEWER_ALL) { 1107962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " SVD: singular values:\n")); 1117962402dSFande Kong for (i = 0; i < n; i++) { 1127962402dSFande Kong if (i % 5 == 0) { 11348a46eb9SPierre Jolivet if (i != 0) PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n")); 1147962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " ")); 1157962402dSFande Kong } 1167962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " %14.12e", (double)PetscRealPart(d[i]))); 1177962402dSFande Kong } 1187962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n")); 1197962402dSFande Kong } else { /* print 5 smallest and 5 largest */ 1209566063dSJacob 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]))); 1219566063dSJacob 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]))); 122426160bdSJed Brown } 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(jac->monitor, ((PetscObject)pc)->tablevel)); 124426160bdSJed Brown } 1259566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Largest and smallest singular values %14.12e %14.12e\n", (double)PetscRealPart(d[0]), (double)PetscRealPart(d[n - 1]))); 126426160bdSJed Brown for (i = 0; i < n - jac->nzero; i++) d[i] = 1.0 / d[i]; 127426160bdSJed Brown for (; i < n; i++) d[i] = 0.0; 1289371c9d4SSatish Balay if (jac->essrank > 0) 1299371c9d4SSatish Balay for (i = 0; i < n - jac->nzero - jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */ 13063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Number of zero or nearly singular values %" PetscInt_FMT "\n", jac->nzero)); 1319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(jac->diag, &d)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 13327c67122SBarry Smith PetscFunctionReturn(0); 13427c67122SBarry Smith } 13527c67122SBarry Smith 1369371c9d4SSatish Balay static PetscErrorCode PCSVDGetVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) { 1373ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 1383ed27f31SJed Brown PetscMPIInt size; 1393ed27f31SJed Brown 1403ed27f31SJed Brown PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1420298fd71SBarry Smith *xred = NULL; 1433ed27f31SJed Brown switch (side) { 1443ed27f31SJed Brown case PC_LEFT: 1453ed27f31SJed Brown if (size == 1) *xred = x; 1463ed27f31SJed Brown else { 1479566063dSJacob Faibussowitsch if (!jac->left2red) PetscCall(VecScatterCreateToAll(x, &jac->left2red, &jac->leftred)); 1483ed27f31SJed Brown if (amode & READ) { 1499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD)); 1509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD)); 1513ed27f31SJed Brown } 1523ed27f31SJed Brown *xred = jac->leftred; 1533ed27f31SJed Brown } 1543ed27f31SJed Brown break; 1553ed27f31SJed Brown case PC_RIGHT: 1563ed27f31SJed Brown if (size == 1) *xred = x; 1573ed27f31SJed Brown else { 1589566063dSJacob Faibussowitsch if (!jac->right2red) PetscCall(VecScatterCreateToAll(x, &jac->right2red, &jac->rightred)); 1593ed27f31SJed Brown if (amode & READ) { 1609566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD)); 1619566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD)); 1623ed27f31SJed Brown } 1633ed27f31SJed Brown *xred = jac->rightred; 1643ed27f31SJed Brown } 1653ed27f31SJed Brown break; 166ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT"); 1673ed27f31SJed Brown } 1683ed27f31SJed Brown PetscFunctionReturn(0); 1693ed27f31SJed Brown } 1703ed27f31SJed Brown 1719371c9d4SSatish Balay static PetscErrorCode PCSVDRestoreVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) { 1723ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 1733ed27f31SJed Brown PetscMPIInt size; 1743ed27f31SJed Brown 1753ed27f31SJed Brown PetscFunctionBegin; 1769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1773ed27f31SJed Brown switch (side) { 1783ed27f31SJed Brown case PC_LEFT: 1793ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1809566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE)); 1819566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE)); 1823ed27f31SJed Brown } 1833ed27f31SJed Brown break; 1843ed27f31SJed Brown case PC_RIGHT: 1853ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1869566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE)); 1879566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE)); 1883ed27f31SJed Brown } 1893ed27f31SJed Brown break; 190ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT"); 1913ed27f31SJed Brown } 1920298fd71SBarry Smith *xred = NULL; 1933ed27f31SJed Brown PetscFunctionReturn(0); 1943ed27f31SJed Brown } 1953ed27f31SJed Brown 19627c67122SBarry Smith /* 19727c67122SBarry Smith PCApply_SVD - Applies the SVD preconditioner to a vector. 19827c67122SBarry Smith 19927c67122SBarry Smith Input Parameters: 20027c67122SBarry Smith . pc - the preconditioner context 20127c67122SBarry Smith . x - input vector 20227c67122SBarry Smith 20327c67122SBarry Smith Output Parameter: 20427c67122SBarry Smith . y - output vector 20527c67122SBarry Smith 20627c67122SBarry Smith Application Interface Routine: PCApply() 20727c67122SBarry Smith */ 2089371c9d4SSatish Balay static PetscErrorCode PCApply_SVD(PC pc, Vec x, Vec y) { 20927c67122SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 2103ed27f31SJed Brown Vec work = jac->work, xred, yred; 21127c67122SBarry Smith 21227c67122SBarry Smith PetscFunctionBegin; 2139566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_RIGHT, READ, x, &xred)); 2149566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_LEFT, WRITE, y, &yred)); 215ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2169566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac->U, xred, work)); 217ef47b4b1SBarry Smith #else 2189566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->U, xred, work)); 219ef47b4b1SBarry Smith #endif 2209566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(work, work, jac->diag)); 221ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2229566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac->Vt, work, yred)); 223ef47b4b1SBarry Smith #else 2249566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->Vt, work, yred)); 225ef47b4b1SBarry Smith #endif 2269566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, READ, x, &xred)); 2279566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_LEFT, WRITE, y, &yred)); 2283ed27f31SJed Brown PetscFunctionReturn(0); 2293ed27f31SJed Brown } 2303ed27f31SJed Brown 2319371c9d4SSatish Balay static PetscErrorCode PCMatApply_SVD(PC pc, Mat X, Mat Y) { 232024db992SStefano Zampini PC_SVD *jac = (PC_SVD *)pc->data; 233024db992SStefano Zampini Mat W; 234024db992SStefano Zampini 235024db992SStefano Zampini PetscFunctionBegin; 236024db992SStefano Zampini PetscCall(MatTransposeMatMult(jac->U, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &W)); 237024db992SStefano Zampini PetscCall(MatDiagonalScale(W, jac->diag, NULL)); 238024db992SStefano Zampini PetscCall(MatTransposeMatMult(jac->Vt, W, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 239024db992SStefano Zampini PetscCall(MatDestroy(&W)); 240024db992SStefano Zampini PetscFunctionReturn(0); 241024db992SStefano Zampini } 242024db992SStefano Zampini 2439371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_SVD(PC pc, Vec x, Vec y) { 2443ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 2453ed27f31SJed Brown Vec work = jac->work, xred, yred; 2463ed27f31SJed Brown 2473ed27f31SJed Brown PetscFunctionBegin; 2489566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_LEFT, READ, x, &xred)); 2499566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_RIGHT, WRITE, y, &yred)); 2509566063dSJacob Faibussowitsch PetscCall(MatMult(jac->Vt, xred, work)); 2519566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(work, work, jac->diag)); 2529566063dSJacob Faibussowitsch PetscCall(MatMult(jac->U, work, yred)); 2539566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_LEFT, READ, x, &xred)); 2549566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, WRITE, y, &yred)); 25527c67122SBarry Smith PetscFunctionReturn(0); 25627c67122SBarry Smith } 25727c67122SBarry Smith 2589371c9d4SSatish Balay static PetscErrorCode PCReset_SVD(PC pc) { 259a2d70de2SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 260a2d70de2SBarry Smith 261a2d70de2SBarry Smith PetscFunctionBegin; 2629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->A)); 2639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->U)); 2649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->Vt)); 2659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diag)); 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->work)); 2679566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&jac->right2red)); 2689566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&jac->left2red)); 2699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->rightred)); 2709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->leftred)); 271a2d70de2SBarry Smith PetscFunctionReturn(0); 272a2d70de2SBarry Smith } 273a2d70de2SBarry Smith 27427c67122SBarry Smith /* 27527c67122SBarry Smith PCDestroy_SVD - Destroys the private context for the SVD preconditioner 27627c67122SBarry Smith that was created with PCCreate_SVD(). 27727c67122SBarry Smith 27827c67122SBarry Smith Input Parameter: 27927c67122SBarry Smith . pc - the preconditioner context 28027c67122SBarry Smith 28127c67122SBarry Smith Application Interface Routine: PCDestroy() 28227c67122SBarry Smith */ 2839371c9d4SSatish Balay static PetscErrorCode PCDestroy_SVD(PC pc) { 284426160bdSJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 28527c67122SBarry Smith 28627c67122SBarry Smith PetscFunctionBegin; 2879566063dSJacob Faibussowitsch PetscCall(PCReset_SVD(pc)); 2889566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&jac->monitor)); 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 29027c67122SBarry Smith PetscFunctionReturn(0); 29127c67122SBarry Smith } 29227c67122SBarry Smith 2939371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_SVD(PC pc, PetscOptionItems *PetscOptionsObject) { 2948f1a2a5eSBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 2957962402dSFande Kong PetscBool flg; 29627c67122SBarry Smith 29727c67122SBarry Smith PetscFunctionBegin; 298d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SVD options"); 2999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_svd_zero_sing", "Singular values smaller than this treated as zero", "None", jac->zerosing, &jac->zerosing, NULL)); 3009566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_svd_ess_rank", "Essential rank of operator (0 to use entire operator)", "None", jac->essrank, &jac->essrank, NULL)); 3017962402dSFande Kong PetscCall(PetscOptionsViewer("-pc_svd_monitor", "Monitor the conditioning, and extremal singular values", "None", &jac->monitor, &jac->monitorformat, &flg)); 302d0609cedSBarry Smith PetscOptionsHeadEnd(); 30327c67122SBarry Smith PetscFunctionReturn(0); 30427c67122SBarry Smith } 30527c67122SBarry Smith 3069371c9d4SSatish Balay static PetscErrorCode PCView_SVD(PC pc, PetscViewer viewer) { 30733761216SBarry Smith PC_SVD *svd = (PC_SVD *)pc->data; 30833761216SBarry Smith PetscBool iascii; 30933761216SBarry Smith 31033761216SBarry Smith PetscFunctionBegin; 3119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 31233761216SBarry Smith if (iascii) { 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " All singular values smaller than %g treated as zero\n", (double)svd->zerosing)); 31463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Provided essential rank of the matrix %" PetscInt_FMT " (all other eigenvalues are zeroed)\n", svd->essrank)); 31533761216SBarry Smith } 31633761216SBarry Smith PetscFunctionReturn(0); 31733761216SBarry Smith } 318*f1580f4eSBarry Smith 31927c67122SBarry Smith /* 32027c67122SBarry Smith PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 32127c67122SBarry Smith and sets this as the private data within the generic preconditioning 32227c67122SBarry Smith context, PC, that was created within PCCreate(). 32327c67122SBarry Smith 32427c67122SBarry Smith Input Parameter: 32527c67122SBarry Smith . pc - the preconditioner context 32627c67122SBarry Smith 32727c67122SBarry Smith Application Interface Routine: PCCreate() 32827c67122SBarry Smith */ 32927c67122SBarry Smith 33027c67122SBarry Smith /*MC 33127c67122SBarry Smith PCSVD - Use pseudo inverse defined by SVD of operator 33227c67122SBarry Smith 33327c67122SBarry Smith Level: advanced 33427c67122SBarry Smith 335*f1580f4eSBarry Smith Options Database Keys: 336147403d9SBarry Smith + -pc_svd_zero_sing <rtol> - Singular values smaller than this are treated as zero 337147403d9SBarry Smith - -pc_svd_monitor - Print information on the extreme singular values of the operator 33827c67122SBarry Smith 339a2b725a8SWilliam Gropp Developer Note: 340a2b725a8SWilliam Gropp This implementation automatically creates a redundant copy of the 3418997ae2eSBarry Smith matrix on each process and uses a sequential SVD solve. Why does it do this instead 342*f1580f4eSBarry Smith of using the composable `PCREDUNDANT` object? 3438997ae2eSBarry Smith 344*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCREDUNDANT` 34527c67122SBarry Smith M*/ 34627c67122SBarry Smith 3479371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc) { 34827c67122SBarry Smith PC_SVD *jac; 349024db992SStefano Zampini PetscMPIInt size = 0; 35027c67122SBarry Smith 35127c67122SBarry Smith PetscFunctionBegin; 35227c67122SBarry Smith /* 35327c67122SBarry Smith Creates the private data structure for this preconditioner and 35427c67122SBarry Smith attach it to the PC object. 35527c67122SBarry Smith */ 3569566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &jac)); 3578f1a2a5eSBarry Smith jac->zerosing = 1.e-12; 35827c67122SBarry Smith pc->data = (void *)jac; 35927c67122SBarry Smith 36027c67122SBarry Smith /* 36127c67122SBarry Smith Set the pointers for the functions that are provided above. 36227c67122SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 36327c67122SBarry Smith are called, they will automatically call these functions. Note we 36427c67122SBarry Smith choose not to provide a couple of these functions since they are 36527c67122SBarry Smith not needed. 36627c67122SBarry Smith */ 367024db992SStefano Zampini 368024db992SStefano Zampini #if defined(PETSC_HAVE_COMPLEX) 369024db992SStefano Zampini PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 370024db992SStefano Zampini #endif 371024db992SStefano Zampini if (size == 1) pc->ops->matapply = PCMatApply_SVD; 37227c67122SBarry Smith pc->ops->apply = PCApply_SVD; 3733ed27f31SJed Brown pc->ops->applytranspose = PCApplyTranspose_SVD; 37427c67122SBarry Smith pc->ops->setup = PCSetUp_SVD; 375a2d70de2SBarry Smith pc->ops->reset = PCReset_SVD; 37627c67122SBarry Smith pc->ops->destroy = PCDestroy_SVD; 37727c67122SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SVD; 37833761216SBarry Smith pc->ops->view = PCView_SVD; 3790a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 38027c67122SBarry Smith PetscFunctionReturn(0); 38127c67122SBarry Smith } 382