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 35f1580f4eSBarry 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 */ 39d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_SVD(PC pc) 40d71ae5a4SJacob Faibussowitsch { 4127c67122SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 4227c67122SBarry Smith PetscScalar *a, *u, *v, *d, *work; 433f83f0d9SMatthew G Knepley PetscBLASInt nb, lwork; 4427c67122SBarry Smith PetscInt i, n; 453ed27f31SJed Brown PetscMPIInt size; 4627c67122SBarry Smith 4727c67122SBarry Smith PetscFunctionBegin; 489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->A)); 499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(((PetscObject)pc->pmat)->comm, &size)); 503ed27f31SJed Brown if (size > 1) { 513ed27f31SJed Brown Mat redmat; 5228d58a37SPierre Jolivet 539566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(pc->pmat, size, PETSC_COMM_SELF, MAT_INITIAL_MATRIX, &redmat)); 549566063dSJacob Faibussowitsch PetscCall(MatConvert(redmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A)); 559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&redmat)); 563ed27f31SJed Brown } else { 579566063dSJacob Faibussowitsch PetscCall(MatConvert(pc->pmat, MATSEQDENSE, MAT_INITIAL_MATRIX, &jac->A)); 583ed27f31SJed Brown } 593ed27f31SJed Brown if (!jac->diag) { /* assume square matrices */ 609566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(jac->A, &jac->diag, &jac->work)); 613ed27f31SJed Brown } 6227c67122SBarry Smith if (!jac->U) { 639566063dSJacob Faibussowitsch PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->U)); 649566063dSJacob Faibussowitsch PetscCall(MatDuplicate(jac->A, MAT_DO_NOT_COPY_VALUES, &jac->Vt)); 6527c67122SBarry Smith } 669566063dSJacob Faibussowitsch PetscCall(MatGetSize(jac->A, &n, NULL)); 6714ce09a1SBarry Smith if (!n) { 689566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Matrix has zero rows, skipping svd\n")); 69*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7014ce09a1SBarry Smith } 719566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &nb)); 7227c67122SBarry Smith lwork = 5 * nb; 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lwork, &work)); 749566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->A, &a)); 759566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->U, &u)); 769566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(jac->Vt, &v)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArray(jac->diag, &d)); 7827c67122SBarry Smith #if !defined(PETSC_USE_COMPLEX) 793f83f0d9SMatthew G Knepley { 803f83f0d9SMatthew G Knepley PetscBLASInt lierr; 819566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 82792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, d, u, &nb, v, &nb, work, &lwork, &lierr)); 8363a3b9bcSJacob Faibussowitsch PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr); 849566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 853f83f0d9SMatthew G Knepley } 8627c67122SBarry Smith #else 87ef47b4b1SBarry Smith { 88ef47b4b1SBarry Smith PetscBLASInt lierr; 89ef47b4b1SBarry Smith PetscReal *rwork, *dd; 909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * nb, &rwork)); 919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nb, &dd)); 929566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 93792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("A", "A", &nb, &nb, a, &nb, dd, u, &nb, v, &nb, work, &lwork, rwork, &lierr)); 9463a3b9bcSJacob Faibussowitsch PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "gesv() error %" PetscBLASInt_FMT, lierr); 959566063dSJacob Faibussowitsch PetscCall(PetscFree(rwork)); 9614ce09a1SBarry Smith for (i = 0; i < n; i++) d[i] = dd[i]; 979566063dSJacob Faibussowitsch PetscCall(PetscFree(dd)); 989566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 99ef47b4b1SBarry Smith } 10027c67122SBarry Smith #endif 1019566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->A, &a)); 1029566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->U, &u)); 1039566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(jac->Vt, &v)); 1049371c9d4SSatish Balay for (i = n - 1; i >= 0; i--) 1059371c9d4SSatish Balay if (PetscRealPart(d[i]) > jac->zerosing) break; 106426160bdSJed Brown jac->nzero = n - 1 - i; 107426160bdSJed Brown if (jac->monitor) { 1089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(jac->monitor, ((PetscObject)pc)->tablevel)); 10963a3b9bcSJacob 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)); 1107962402dSFande Kong if (n < 10 || jac->monitorformat == PETSC_VIEWER_ALL) { 1117962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " SVD: singular values:\n")); 1127962402dSFande Kong for (i = 0; i < n; i++) { 1137962402dSFande Kong if (i % 5 == 0) { 11448a46eb9SPierre Jolivet if (i != 0) PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n")); 1157962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " ")); 1167962402dSFande Kong } 1177962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, " %14.12e", (double)PetscRealPart(d[i]))); 1187962402dSFande Kong } 1197962402dSFande Kong PetscCall(PetscViewerASCIIPrintf(jac->monitor, "\n")); 1207962402dSFande Kong } else { /* print 5 smallest and 5 largest */ 1219566063dSJacob 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]))); 1229566063dSJacob 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]))); 123426160bdSJed Brown } 1249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(jac->monitor, ((PetscObject)pc)->tablevel)); 125426160bdSJed Brown } 1269566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Largest and smallest singular values %14.12e %14.12e\n", (double)PetscRealPart(d[0]), (double)PetscRealPart(d[n - 1]))); 127426160bdSJed Brown for (i = 0; i < n - jac->nzero; i++) d[i] = 1.0 / d[i]; 128426160bdSJed Brown for (; i < n; i++) d[i] = 0.0; 1299371c9d4SSatish Balay if (jac->essrank > 0) 1309371c9d4SSatish Balay for (i = 0; i < n - jac->nzero - jac->essrank; i++) d[i] = 0.0; /* Skip all but essrank eigenvalues */ 13163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Number of zero or nearly singular values %" PetscInt_FMT "\n", jac->nzero)); 1329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(jac->diag, &d)); 1339566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 134*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13527c67122SBarry Smith } 13627c67122SBarry Smith 137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSVDGetVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) 138d71ae5a4SJacob Faibussowitsch { 1393ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 1403ed27f31SJed Brown PetscMPIInt size; 1413ed27f31SJed Brown 1423ed27f31SJed Brown PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1440298fd71SBarry Smith *xred = NULL; 1453ed27f31SJed Brown switch (side) { 1463ed27f31SJed Brown case PC_LEFT: 1473ed27f31SJed Brown if (size == 1) *xred = x; 1483ed27f31SJed Brown else { 1499566063dSJacob Faibussowitsch if (!jac->left2red) PetscCall(VecScatterCreateToAll(x, &jac->left2red, &jac->leftred)); 1503ed27f31SJed Brown if (amode & READ) { 1519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD)); 1529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->left2red, x, jac->leftred, INSERT_VALUES, SCATTER_FORWARD)); 1533ed27f31SJed Brown } 1543ed27f31SJed Brown *xred = jac->leftred; 1553ed27f31SJed Brown } 1563ed27f31SJed Brown break; 1573ed27f31SJed Brown case PC_RIGHT: 1583ed27f31SJed Brown if (size == 1) *xred = x; 1593ed27f31SJed Brown else { 1609566063dSJacob Faibussowitsch if (!jac->right2red) PetscCall(VecScatterCreateToAll(x, &jac->right2red, &jac->rightred)); 1613ed27f31SJed Brown if (amode & READ) { 1629566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD)); 1639566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->right2red, x, jac->rightred, INSERT_VALUES, SCATTER_FORWARD)); 1643ed27f31SJed Brown } 1653ed27f31SJed Brown *xred = jac->rightred; 1663ed27f31SJed Brown } 1673ed27f31SJed Brown break; 168d71ae5a4SJacob Faibussowitsch default: 169d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT"); 1703ed27f31SJed Brown } 171*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1723ed27f31SJed Brown } 1733ed27f31SJed Brown 174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSVDRestoreVec(PC pc, PCSide side, AccessMode amode, Vec x, Vec *xred) 175d71ae5a4SJacob Faibussowitsch { 1763ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 1773ed27f31SJed Brown PetscMPIInt size; 1783ed27f31SJed Brown 1793ed27f31SJed Brown PetscFunctionBegin; 1809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1813ed27f31SJed Brown switch (side) { 1823ed27f31SJed Brown case PC_LEFT: 1833ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1849566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE)); 1859566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->left2red, jac->leftred, x, INSERT_VALUES, SCATTER_REVERSE)); 1863ed27f31SJed Brown } 1873ed27f31SJed Brown break; 1883ed27f31SJed Brown case PC_RIGHT: 1893ed27f31SJed Brown if (size != 1 && amode & WRITE) { 1909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE)); 1919566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(jac->right2red, jac->rightred, x, INSERT_VALUES, SCATTER_REVERSE)); 1923ed27f31SJed Brown } 1933ed27f31SJed Brown break; 194d71ae5a4SJacob Faibussowitsch default: 195d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Side must be LEFT or RIGHT"); 1963ed27f31SJed Brown } 1970298fd71SBarry Smith *xred = NULL; 198*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1993ed27f31SJed Brown } 2003ed27f31SJed Brown 20127c67122SBarry Smith /* 20227c67122SBarry Smith PCApply_SVD - Applies the SVD preconditioner to a vector. 20327c67122SBarry Smith 20427c67122SBarry Smith Input Parameters: 20527c67122SBarry Smith . pc - the preconditioner context 20627c67122SBarry Smith . x - input vector 20727c67122SBarry Smith 20827c67122SBarry Smith Output Parameter: 20927c67122SBarry Smith . y - output vector 21027c67122SBarry Smith 21127c67122SBarry Smith Application Interface Routine: PCApply() 21227c67122SBarry Smith */ 213d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SVD(PC pc, Vec x, Vec y) 214d71ae5a4SJacob Faibussowitsch { 21527c67122SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 2163ed27f31SJed Brown Vec work = jac->work, xred, yred; 21727c67122SBarry Smith 21827c67122SBarry Smith PetscFunctionBegin; 2199566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_RIGHT, READ, x, &xred)); 2209566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_LEFT, WRITE, y, &yred)); 221ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2229566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac->U, xred, work)); 223ef47b4b1SBarry Smith #else 2249566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->U, xred, work)); 225ef47b4b1SBarry Smith #endif 2269566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(work, work, jac->diag)); 227ef47b4b1SBarry Smith #if !defined(PETSC_USE_COMPLEX) 2289566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac->Vt, work, yred)); 229ef47b4b1SBarry Smith #else 2309566063dSJacob Faibussowitsch PetscCall(MatMultHermitianTranspose(jac->Vt, work, yred)); 231ef47b4b1SBarry Smith #endif 2329566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, READ, x, &xred)); 2339566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_LEFT, WRITE, y, &yred)); 234*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2353ed27f31SJed Brown } 2363ed27f31SJed Brown 237d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_SVD(PC pc, Mat X, Mat Y) 238d71ae5a4SJacob Faibussowitsch { 239024db992SStefano Zampini PC_SVD *jac = (PC_SVD *)pc->data; 240024db992SStefano Zampini Mat W; 241024db992SStefano Zampini 242024db992SStefano Zampini PetscFunctionBegin; 243024db992SStefano Zampini PetscCall(MatTransposeMatMult(jac->U, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &W)); 244024db992SStefano Zampini PetscCall(MatDiagonalScale(W, jac->diag, NULL)); 245024db992SStefano Zampini PetscCall(MatTransposeMatMult(jac->Vt, W, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 246024db992SStefano Zampini PetscCall(MatDestroy(&W)); 247*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 248024db992SStefano Zampini } 249024db992SStefano Zampini 250d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_SVD(PC pc, Vec x, Vec y) 251d71ae5a4SJacob Faibussowitsch { 2523ed27f31SJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 2533ed27f31SJed Brown Vec work = jac->work, xred, yred; 2543ed27f31SJed Brown 2553ed27f31SJed Brown PetscFunctionBegin; 2569566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_LEFT, READ, x, &xred)); 2579566063dSJacob Faibussowitsch PetscCall(PCSVDGetVec(pc, PC_RIGHT, WRITE, y, &yred)); 2589566063dSJacob Faibussowitsch PetscCall(MatMult(jac->Vt, xred, work)); 2599566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(work, work, jac->diag)); 2609566063dSJacob Faibussowitsch PetscCall(MatMult(jac->U, work, yred)); 2619566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_LEFT, READ, x, &xred)); 2629566063dSJacob Faibussowitsch PetscCall(PCSVDRestoreVec(pc, PC_RIGHT, WRITE, y, &yred)); 263*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26427c67122SBarry Smith } 26527c67122SBarry Smith 266d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_SVD(PC pc) 267d71ae5a4SJacob Faibussowitsch { 268a2d70de2SBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 269a2d70de2SBarry Smith 270a2d70de2SBarry Smith PetscFunctionBegin; 2719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->A)); 2729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->U)); 2739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&jac->Vt)); 2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diag)); 2759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->work)); 2769566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&jac->right2red)); 2779566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&jac->left2red)); 2789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->rightred)); 2799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->leftred)); 280*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 281a2d70de2SBarry Smith } 282a2d70de2SBarry Smith 28327c67122SBarry Smith /* 28427c67122SBarry Smith PCDestroy_SVD - Destroys the private context for the SVD preconditioner 28527c67122SBarry Smith that was created with PCCreate_SVD(). 28627c67122SBarry Smith 28727c67122SBarry Smith Input Parameter: 28827c67122SBarry Smith . pc - the preconditioner context 28927c67122SBarry Smith 29027c67122SBarry Smith Application Interface Routine: PCDestroy() 29127c67122SBarry Smith */ 292d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SVD(PC pc) 293d71ae5a4SJacob Faibussowitsch { 294426160bdSJed Brown PC_SVD *jac = (PC_SVD *)pc->data; 29527c67122SBarry Smith 29627c67122SBarry Smith PetscFunctionBegin; 2979566063dSJacob Faibussowitsch PetscCall(PCReset_SVD(pc)); 2989566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&jac->monitor)); 2999566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 300*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30127c67122SBarry Smith } 30227c67122SBarry Smith 303d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_SVD(PC pc, PetscOptionItems *PetscOptionsObject) 304d71ae5a4SJacob Faibussowitsch { 3058f1a2a5eSBarry Smith PC_SVD *jac = (PC_SVD *)pc->data; 3067962402dSFande Kong PetscBool flg; 30727c67122SBarry Smith 30827c67122SBarry Smith PetscFunctionBegin; 309d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SVD options"); 3109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_svd_zero_sing", "Singular values smaller than this treated as zero", "None", jac->zerosing, &jac->zerosing, NULL)); 3119566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_svd_ess_rank", "Essential rank of operator (0 to use entire operator)", "None", jac->essrank, &jac->essrank, NULL)); 3127962402dSFande Kong PetscCall(PetscOptionsViewer("-pc_svd_monitor", "Monitor the conditioning, and extremal singular values", "None", &jac->monitor, &jac->monitorformat, &flg)); 313d0609cedSBarry Smith PetscOptionsHeadEnd(); 314*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31527c67122SBarry Smith } 31627c67122SBarry Smith 317d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_SVD(PC pc, PetscViewer viewer) 318d71ae5a4SJacob Faibussowitsch { 31933761216SBarry Smith PC_SVD *svd = (PC_SVD *)pc->data; 32033761216SBarry Smith PetscBool iascii; 32133761216SBarry Smith 32233761216SBarry Smith PetscFunctionBegin; 3239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 32433761216SBarry Smith if (iascii) { 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " All singular values smaller than %g treated as zero\n", (double)svd->zerosing)); 32663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Provided essential rank of the matrix %" PetscInt_FMT " (all other eigenvalues are zeroed)\n", svd->essrank)); 32733761216SBarry Smith } 328*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32933761216SBarry Smith } 330f1580f4eSBarry Smith 33127c67122SBarry Smith /* 33227c67122SBarry Smith PCCreate_SVD - Creates a SVD preconditioner context, PC_SVD, 33327c67122SBarry Smith and sets this as the private data within the generic preconditioning 33427c67122SBarry Smith context, PC, that was created within PCCreate(). 33527c67122SBarry Smith 33627c67122SBarry Smith Input Parameter: 33727c67122SBarry Smith . pc - the preconditioner context 33827c67122SBarry Smith 33927c67122SBarry Smith Application Interface Routine: PCCreate() 34027c67122SBarry Smith */ 34127c67122SBarry Smith 34227c67122SBarry Smith /*MC 34327c67122SBarry Smith PCSVD - Use pseudo inverse defined by SVD of operator 34427c67122SBarry Smith 34527c67122SBarry Smith Level: advanced 34627c67122SBarry Smith 347f1580f4eSBarry Smith Options Database Keys: 348147403d9SBarry Smith + -pc_svd_zero_sing <rtol> - Singular values smaller than this are treated as zero 349147403d9SBarry Smith - -pc_svd_monitor - Print information on the extreme singular values of the operator 35027c67122SBarry Smith 351a2b725a8SWilliam Gropp Developer Note: 352a2b725a8SWilliam Gropp This implementation automatically creates a redundant copy of the 3538997ae2eSBarry Smith matrix on each process and uses a sequential SVD solve. Why does it do this instead 354f1580f4eSBarry Smith of using the composable `PCREDUNDANT` object? 3558997ae2eSBarry Smith 356f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCREDUNDANT` 35727c67122SBarry Smith M*/ 35827c67122SBarry Smith 359d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SVD(PC pc) 360d71ae5a4SJacob Faibussowitsch { 36127c67122SBarry Smith PC_SVD *jac; 362024db992SStefano Zampini PetscMPIInt size = 0; 36327c67122SBarry Smith 36427c67122SBarry Smith PetscFunctionBegin; 36527c67122SBarry Smith /* 36627c67122SBarry Smith Creates the private data structure for this preconditioner and 36727c67122SBarry Smith attach it to the PC object. 36827c67122SBarry Smith */ 3694dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 3708f1a2a5eSBarry Smith jac->zerosing = 1.e-12; 37127c67122SBarry Smith pc->data = (void *)jac; 37227c67122SBarry Smith 37327c67122SBarry Smith /* 37427c67122SBarry Smith Set the pointers for the functions that are provided above. 37527c67122SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 37627c67122SBarry Smith are called, they will automatically call these functions. Note we 37727c67122SBarry Smith choose not to provide a couple of these functions since they are 37827c67122SBarry Smith not needed. 37927c67122SBarry Smith */ 380024db992SStefano Zampini 381024db992SStefano Zampini #if defined(PETSC_HAVE_COMPLEX) 382024db992SStefano Zampini PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 383024db992SStefano Zampini #endif 384024db992SStefano Zampini if (size == 1) pc->ops->matapply = PCMatApply_SVD; 38527c67122SBarry Smith pc->ops->apply = PCApply_SVD; 3863ed27f31SJed Brown pc->ops->applytranspose = PCApplyTranspose_SVD; 38727c67122SBarry Smith pc->ops->setup = PCSetUp_SVD; 388a2d70de2SBarry Smith pc->ops->reset = PCReset_SVD; 38927c67122SBarry Smith pc->ops->destroy = PCDestroy_SVD; 39027c67122SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SVD; 39133761216SBarry Smith pc->ops->view = PCView_SVD; 3920a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 393*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39427c67122SBarry Smith } 395