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
PCDestroy_TFS(PC pc)1666976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_TFS(PC pc)
17d71ae5a4SJacob Faibussowitsch {
186218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data;
196218e575SSatish Balay
206218e575SSatish Balay PetscFunctionBegin;
216218e575SSatish Balay /* free the XXT datastructures */
221baa6e33SBarry Smith if (tfs->xxt) PetscCall(XXT_free(tfs->xxt));
231baa6e33SBarry Smith if (tfs->xyt) PetscCall(XYT_free(tfs->xyt));
249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tfs->b));
259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tfs->xd));
269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tfs->xo));
279566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data));
283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
296218e575SSatish Balay }
306218e575SSatish Balay
PCApply_TFS_XXT(PC pc,Vec x,Vec y)31d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_TFS_XXT(PC pc, Vec x, Vec y)
32d71ae5a4SJacob Faibussowitsch {
336218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data;
343468409bSBarry Smith PetscScalar *yy;
353468409bSBarry Smith const PetscScalar *xx;
366218e575SSatish Balay
376218e575SSatish Balay PetscFunctionBegin;
389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx));
399566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy));
409566063dSJacob Faibussowitsch PetscCall(XXT_solve(tfs->xxt, yy, (PetscScalar *)xx));
419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx));
429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy));
433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
446218e575SSatish Balay }
456218e575SSatish Balay
PCApply_TFS_XYT(PC pc,Vec x,Vec y)46d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_TFS_XYT(PC pc, Vec x, Vec y)
47d71ae5a4SJacob Faibussowitsch {
486218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data;
493468409bSBarry Smith PetscScalar *yy;
503468409bSBarry Smith const PetscScalar *xx;
516218e575SSatish Balay
526218e575SSatish Balay PetscFunctionBegin;
539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx));
549566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy));
559566063dSJacob Faibussowitsch PetscCall(XYT_solve(tfs->xyt, yy, (PetscScalar *)xx));
569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx));
579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy));
583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
596218e575SSatish Balay }
606218e575SSatish Balay
PCTFSLocalMult_TFS(PC pc,PetscScalar * xin,PetscScalar * xout)61d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTFSLocalMult_TFS(PC pc, PetscScalar *xin, PetscScalar *xout)
62d71ae5a4SJacob Faibussowitsch {
636218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data;
646218e575SSatish Balay Mat A = pc->pmat;
656218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data;
666218e575SSatish Balay
676218e575SSatish Balay PetscFunctionBegin;
689566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tfs->b, xout));
699566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tfs->xd, xin));
709566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(tfs->xo, xin + tfs->nd));
719566063dSJacob Faibussowitsch PetscCall(MatMult(a->A, tfs->xd, tfs->b));
729566063dSJacob Faibussowitsch PetscCall(MatMultAdd(a->B, tfs->xo, tfs->b, tfs->b));
739566063dSJacob Faibussowitsch PetscCall(VecResetArray(tfs->b));
749566063dSJacob Faibussowitsch PetscCall(VecResetArray(tfs->xd));
759566063dSJacob Faibussowitsch PetscCall(VecResetArray(tfs->xo));
763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
776218e575SSatish Balay }
786218e575SSatish Balay
PCSetUp_TFS(PC pc)79d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_TFS(PC pc)
80d71ae5a4SJacob Faibussowitsch {
816218e575SSatish Balay PC_TFS *tfs = (PC_TFS *)pc->data;
826218e575SSatish Balay Mat A = pc->pmat;
836218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data;
846218e575SSatish Balay PetscInt *localtoglobal, ncol, i;
85ace3abfcSBarry Smith PetscBool ismpiaij;
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
1049566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)pc));
105b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE) {
1066218e575SSatish Balay tfs->xxt = XXT_new();
1079566063dSJacob Faibussowitsch PetscCall(XXT_factor(tfs->xxt, localtoglobal, A->rmap->n, ncol, (PetscErrorCode (*)(void *, PetscScalar *, PetscScalar *))PCTFSLocalMult_TFS, pc));
1086218e575SSatish Balay pc->ops->apply = PCApply_TFS_XXT;
1096218e575SSatish Balay } else {
1106218e575SSatish Balay tfs->xyt = XYT_new();
1119566063dSJacob Faibussowitsch PetscCall(XYT_factor(tfs->xyt, localtoglobal, A->rmap->n, ncol, (PetscErrorCode (*)(void *, PetscScalar *, PetscScalar *))PCTFSLocalMult_TFS, pc));
1126218e575SSatish Balay pc->ops->apply = PCApply_TFS_XYT;
1136218e575SSatish Balay }
1146218e575SSatish Balay
1159566063dSJacob Faibussowitsch PetscCall(PetscFree(localtoglobal));
1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1176218e575SSatish Balay }
1186218e575SSatish Balay
11994c0dd61SBarry Smith /*MC
12094c0dd61SBarry Smith PCTFS - A parallel direct solver intended for problems with very few unknowns (like the
121a21bf2d2SBarry Smith coarse grid in multigrid). Performs a Cholesky or LU factorization of a matrix defined by
122*562efe2eSBarry Smith its local matrix-vector product.
12394c0dd61SBarry Smith
12494c0dd61SBarry Smith Level: beginner
12594c0dd61SBarry Smith
12695452b02SPatrick Sanan Notes:
127f1580f4eSBarry Smith Only implemented for the `MATMPIAIJ` matrices
12894c0dd61SBarry Smith
129f1580f4eSBarry Smith Only works on a solver object that lives on all of `PETSC_COMM_WORLD`!
1307773e740SBarry Smith
131f1580f4eSBarry Smith Only works for real numbers (is not built if `PetscScalar` is complex)
132f1580f4eSBarry Smith
133f1580f4eSBarry Smith Implemented by Henry M. Tufo III and Paul Fischer originally for Nek5000 and called XXT or XYT
134a21bf2d2SBarry Smith
135*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`
13694c0dd61SBarry Smith M*/
PCCreate_TFS(PC pc)137d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_TFS(PC pc)
138d71ae5a4SJacob Faibussowitsch {
1396218e575SSatish Balay PC_TFS *tfs;
1403cd2be91SBarry Smith PetscMPIInt cmp;
1416218e575SSatish Balay
1426218e575SSatish Balay PetscFunctionBegin;
1439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PETSC_COMM_WORLD, PetscObjectComm((PetscObject)pc), &cmp));
1442472a847SBarry Smith PetscCheck(cmp == MPI_IDENT || cmp == MPI_CONGRUENT, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "TFS only works with PETSC_COMM_WORLD objects");
1454dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&tfs));
1466218e575SSatish Balay
1470a545947SLisandro Dalcin tfs->xxt = NULL;
1480a545947SLisandro Dalcin tfs->xyt = NULL;
1490a545947SLisandro Dalcin tfs->b = NULL;
1500a545947SLisandro Dalcin tfs->xd = NULL;
1510a545947SLisandro Dalcin tfs->xo = NULL;
1526218e575SSatish Balay tfs->nd = 0;
1536218e575SSatish Balay
1540a545947SLisandro Dalcin pc->ops->apply = NULL;
1550a545947SLisandro Dalcin pc->ops->applytranspose = NULL;
1566218e575SSatish Balay pc->ops->setup = PCSetUp_TFS;
1576218e575SSatish Balay pc->ops->destroy = PCDestroy_TFS;
1580a545947SLisandro Dalcin pc->ops->applyrichardson = NULL;
1590a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL;
1600a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL;
1616218e575SSatish Balay pc->data = (void *)tfs;
1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1636218e575SSatish Balay }
164