16218e575SSatish Balay /* 26218e575SSatish Balay Provides an interface to the Tufo-Fischer parallel direct solver 36218e575SSatish Balay */ 46218e575SSatish Balay 5af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 7c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h> 86218e575SSatish Balay 96218e575SSatish Balay typedef struct { 106218e575SSatish Balay xxt_ADT xxt; 116218e575SSatish Balay xyt_ADT xyt; 126218e575SSatish Balay Vec b, xd, xo; 136218e575SSatish Balay PetscInt nd; 146218e575SSatish Balay } PC_TFS; 156218e575SSatish Balay 169371c9d4SSatish Balay PetscErrorCode PCDestroy_TFS(PC pc) { 176218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data; 186218e575SSatish Balay 196218e575SSatish Balay PetscFunctionBegin; 206218e575SSatish Balay /* free the XXT datastructures */ 211baa6e33SBarry Smith if (tfs->xxt) PetscCall(XXT_free(tfs->xxt)); 221baa6e33SBarry Smith if (tfs->xyt) PetscCall(XYT_free(tfs->xyt)); 239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tfs->b)); 249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tfs->xd)); 259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tfs->xo)); 269566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 276218e575SSatish Balay PetscFunctionReturn(0); 286218e575SSatish Balay } 296218e575SSatish Balay 309371c9d4SSatish Balay static PetscErrorCode PCApply_TFS_XXT(PC pc, Vec x, Vec y) { 316218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data; 323468409bSBarry Smith PetscScalar *yy; 333468409bSBarry Smith const PetscScalar *xx; 346218e575SSatish Balay 356218e575SSatish Balay PetscFunctionBegin; 369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 379566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 389566063dSJacob Faibussowitsch PetscCall(XXT_solve(tfs->xxt, yy, (PetscScalar *)xx)); 399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 416218e575SSatish Balay PetscFunctionReturn(0); 426218e575SSatish Balay } 436218e575SSatish Balay 449371c9d4SSatish Balay static PetscErrorCode PCApply_TFS_XYT(PC pc, Vec x, Vec y) { 456218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data; 463468409bSBarry Smith PetscScalar *yy; 473468409bSBarry Smith const PetscScalar *xx; 486218e575SSatish Balay 496218e575SSatish Balay PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 519566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 529566063dSJacob Faibussowitsch PetscCall(XYT_solve(tfs->xyt, yy, (PetscScalar *)xx)); 539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 556218e575SSatish Balay PetscFunctionReturn(0); 566218e575SSatish Balay } 576218e575SSatish Balay 589371c9d4SSatish Balay static PetscErrorCode PCTFSLocalMult_TFS(PC pc, PetscScalar *xin, PetscScalar *xout) { 596218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data; 606218e575SSatish Balay Mat A = pc->pmat; 616218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data; 626218e575SSatish Balay 636218e575SSatish Balay PetscFunctionBegin; 649566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tfs->b, xout)); 659566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tfs->xd, xin)); 669566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tfs->xo, xin + tfs->nd)); 679566063dSJacob Faibussowitsch PetscCall(MatMult(a->A, tfs->xd, tfs->b)); 689566063dSJacob Faibussowitsch PetscCall(MatMultAdd(a->B, tfs->xo, tfs->b, tfs->b)); 699566063dSJacob Faibussowitsch PetscCall(VecResetArray(tfs->b)); 709566063dSJacob Faibussowitsch PetscCall(VecResetArray(tfs->xd)); 719566063dSJacob Faibussowitsch PetscCall(VecResetArray(tfs->xo)); 726218e575SSatish Balay PetscFunctionReturn(0); 736218e575SSatish Balay } 746218e575SSatish Balay 759371c9d4SSatish Balay static PetscErrorCode PCSetUp_TFS(PC pc) { 766218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data; 776218e575SSatish Balay Mat A = pc->pmat; 786218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data; 796218e575SSatish Balay PetscInt *localtoglobal, ncol, i; 80ace3abfcSBarry Smith PetscBool ismpiaij; 816218e575SSatish Balay 826218e575SSatish Balay /* 83ace3abfcSBarry Smith PetscBool issymmetric; 846218e575SSatish Balay Petsc Real tol = 0.0; 856218e575SSatish Balay */ 866218e575SSatish Balay 876218e575SSatish Balay PetscFunctionBegin; 8808401ef6SPierre Jolivet PetscCheck(A->cmap->N == A->rmap->N, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "matrix must be square"); 899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATMPIAIJ, &ismpiaij)); 9028b400f6SJacob Faibussowitsch PetscCheck(ismpiaij, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only supports MPIAIJ matrices"); 916218e575SSatish Balay 926218e575SSatish Balay /* generate the local to global mapping */ 93d0f46423SBarry Smith ncol = a->A->cmap->n + a->B->cmap->n; 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncol, &localtoglobal)); 952fa5cd67SKarl Rupp for (i = 0; i < a->A->cmap->n; i++) localtoglobal[i] = A->cmap->rstart + i + 1; 962fa5cd67SKarl Rupp for (i = 0; i < a->B->cmap->n; i++) localtoglobal[i + a->A->cmap->n] = a->garray[i] + 1; 972fa5cd67SKarl Rupp 986218e575SSatish Balay /* generate the vectors needed for the local solves */ 999566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, a->A->rmap->n, NULL, &tfs->b)); 1009566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, a->A->cmap->n, NULL, &tfs->xd)); 1019566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, a->B->cmap->n, NULL, &tfs->xo)); 102d0f46423SBarry Smith tfs->nd = a->A->cmap->n; 1036218e575SSatish Balay 1046218e575SSatish Balay /* ierr = MatIsSymmetric(A,tol,&issymmetric); */ 1056218e575SSatish Balay /* if (issymmetric) { */ 1069566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)pc)); 107b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE) { 1086218e575SSatish Balay tfs->xxt = XXT_new(); 1099566063dSJacob Faibussowitsch PetscCall(XXT_factor(tfs->xxt, localtoglobal, A->rmap->n, ncol, (PetscErrorCode(*)(void *, PetscScalar *, PetscScalar *))PCTFSLocalMult_TFS, pc)); 1106218e575SSatish Balay pc->ops->apply = PCApply_TFS_XXT; 1116218e575SSatish Balay } else { 1126218e575SSatish Balay tfs->xyt = XYT_new(); 1139566063dSJacob Faibussowitsch PetscCall(XYT_factor(tfs->xyt, localtoglobal, A->rmap->n, ncol, (PetscErrorCode(*)(void *, PetscScalar *, PetscScalar *))PCTFSLocalMult_TFS, pc)); 1146218e575SSatish Balay pc->ops->apply = PCApply_TFS_XYT; 1156218e575SSatish Balay } 1166218e575SSatish Balay 1179566063dSJacob Faibussowitsch PetscCall(PetscFree(localtoglobal)); 1186218e575SSatish Balay PetscFunctionReturn(0); 1196218e575SSatish Balay } 1206218e575SSatish Balay 1219371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_TFS(PC pc, PetscOptionItems *PetscOptionsObject) { 1226218e575SSatish Balay PetscFunctionBegin; 1236218e575SSatish Balay PetscFunctionReturn(0); 1246218e575SSatish Balay } 1259371c9d4SSatish Balay static PetscErrorCode PCView_TFS(PC pc, PetscViewer viewer) { 1266218e575SSatish Balay PetscFunctionBegin; 1276218e575SSatish Balay PetscFunctionReturn(0); 1286218e575SSatish Balay } 1296218e575SSatish Balay 13094c0dd61SBarry Smith /*MC 13194c0dd61SBarry Smith PCTFS - A parallel direct solver intended for problems with very few unknowns (like the 132a21bf2d2SBarry Smith coarse grid in multigrid). Performs a Cholesky or LU factorization of a matrix defined by 133a21bf2d2SBarry Smith its local matrix vector product. 13494c0dd61SBarry Smith 13594c0dd61SBarry Smith Level: beginner 13694c0dd61SBarry Smith 13795452b02SPatrick Sanan Notes: 138*f1580f4eSBarry Smith Only implemented for the `MATMPIAIJ` matrices 13994c0dd61SBarry Smith 140*f1580f4eSBarry Smith Only works on a solver object that lives on all of `PETSC_COMM_WORLD`! 1417773e740SBarry Smith 142*f1580f4eSBarry Smith Only works for real numbers (is not built if `PetscScalar` is complex) 143*f1580f4eSBarry Smith 144*f1580f4eSBarry Smith Implemented by Henry M. Tufo III and Paul Fischer originally for Nek5000 and called XXT or XYT 145a21bf2d2SBarry Smith 146db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC` 14794c0dd61SBarry Smith M*/ 1489371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_TFS(PC pc) { 1496218e575SSatish Balay PC_TFS *tfs; 1503cd2be91SBarry Smith PetscMPIInt cmp; 1516218e575SSatish Balay 1526218e575SSatish Balay PetscFunctionBegin; 1539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PETSC_COMM_WORLD, PetscObjectComm((PetscObject)pc), &cmp)); 1542472a847SBarry Smith PetscCheck(cmp == MPI_IDENT || cmp == MPI_CONGRUENT, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "TFS only works with PETSC_COMM_WORLD objects"); 1559566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &tfs)); 1566218e575SSatish Balay 1570a545947SLisandro Dalcin tfs->xxt = NULL; 1580a545947SLisandro Dalcin tfs->xyt = NULL; 1590a545947SLisandro Dalcin tfs->b = NULL; 1600a545947SLisandro Dalcin tfs->xd = NULL; 1610a545947SLisandro Dalcin tfs->xo = NULL; 1626218e575SSatish Balay tfs->nd = 0; 1636218e575SSatish Balay 1640a545947SLisandro Dalcin pc->ops->apply = NULL; 1650a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 1666218e575SSatish Balay pc->ops->setup = PCSetUp_TFS; 1676218e575SSatish Balay pc->ops->destroy = PCDestroy_TFS; 1686218e575SSatish Balay pc->ops->setfromoptions = PCSetFromOptions_TFS; 1696218e575SSatish Balay pc->ops->view = PCView_TFS; 1700a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 1710a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 1720a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 1736218e575SSatish Balay pc->data = (void *)tfs; 1746218e575SSatish Balay PetscFunctionReturn(0); 1756218e575SSatish Balay } 176