xref: /petsc/src/ksp/pc/impls/tfs/tfs.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
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 
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 
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 
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 
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 
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   /*
88ace3abfcSBarry Smith   PetscBool      issymmetric;
896218e575SSatish Balay   Petsc Real tol = 0.0;
906218e575SSatish Balay   */
916218e575SSatish Balay 
926218e575SSatish Balay   PetscFunctionBegin;
9308401ef6SPierre Jolivet   PetscCheck(A->cmap->N == A->rmap->N, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_SIZ, "matrix must be square");
949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATMPIAIJ, &ismpiaij));
9528b400f6SJacob Faibussowitsch   PetscCheck(ismpiaij, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only supports MPIAIJ matrices");
966218e575SSatish Balay 
976218e575SSatish Balay   /* generate the local to global mapping */
98d0f46423SBarry Smith   ncol = a->A->cmap->n + a->B->cmap->n;
999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ncol, &localtoglobal));
1002fa5cd67SKarl Rupp   for (i = 0; i < a->A->cmap->n; i++) localtoglobal[i] = A->cmap->rstart + i + 1;
1012fa5cd67SKarl Rupp   for (i = 0; i < a->B->cmap->n; i++) localtoglobal[i + a->A->cmap->n] = a->garray[i] + 1;
1022fa5cd67SKarl Rupp 
1036218e575SSatish Balay   /* generate the vectors needed for the local solves */
1049566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, a->A->rmap->n, NULL, &tfs->b));
1059566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, a->A->cmap->n, NULL, &tfs->xd));
1069566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, a->B->cmap->n, NULL, &tfs->xo));
107d0f46423SBarry Smith   tfs->nd = a->A->cmap->n;
1086218e575SSatish Balay 
1096218e575SSatish Balay   /*  ierr =  MatIsSymmetric(A,tol,&issymmetric); */
1106218e575SSatish Balay   /*  if (issymmetric) { */
1119566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)pc));
112b94d7dedSBarry Smith   if (A->symmetric == PETSC_BOOL3_TRUE) {
1136218e575SSatish Balay     tfs->xxt = XXT_new();
1149566063dSJacob Faibussowitsch     PetscCall(XXT_factor(tfs->xxt, localtoglobal, A->rmap->n, ncol, (PetscErrorCode(*)(void *, PetscScalar *, PetscScalar *))PCTFSLocalMult_TFS, pc));
1156218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XXT;
1166218e575SSatish Balay   } else {
1176218e575SSatish Balay     tfs->xyt = XYT_new();
1189566063dSJacob Faibussowitsch     PetscCall(XYT_factor(tfs->xyt, localtoglobal, A->rmap->n, ncol, (PetscErrorCode(*)(void *, PetscScalar *, PetscScalar *))PCTFSLocalMult_TFS, pc));
1196218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XYT;
1206218e575SSatish Balay   }
1216218e575SSatish Balay 
1229566063dSJacob Faibussowitsch   PetscCall(PetscFree(localtoglobal));
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1246218e575SSatish Balay }
1256218e575SSatish Balay 
126d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_TFS(PC pc, PetscOptionItems *PetscOptionsObject)
127d71ae5a4SJacob Faibussowitsch {
1286218e575SSatish Balay   PetscFunctionBegin;
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1306218e575SSatish Balay }
131d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_TFS(PC pc, PetscViewer viewer)
132d71ae5a4SJacob Faibussowitsch {
1336218e575SSatish Balay   PetscFunctionBegin;
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1356218e575SSatish Balay }
1366218e575SSatish Balay 
13794c0dd61SBarry Smith /*MC
13894c0dd61SBarry Smith      PCTFS - A parallel direct solver intended for problems with very few unknowns (like the
139a21bf2d2SBarry Smith          coarse grid in multigrid). Performs a Cholesky or LU factorization of a matrix defined by
140*562efe2eSBarry Smith          its local matrix-vector product.
14194c0dd61SBarry Smith 
14294c0dd61SBarry Smith    Level: beginner
14394c0dd61SBarry Smith 
14495452b02SPatrick Sanan    Notes:
145f1580f4eSBarry Smith     Only implemented for the `MATMPIAIJ` matrices
14694c0dd61SBarry Smith 
147f1580f4eSBarry Smith     Only works on a solver object that lives on all of `PETSC_COMM_WORLD`!
1487773e740SBarry Smith 
149f1580f4eSBarry Smith     Only works for real numbers (is not built if `PetscScalar` is complex)
150f1580f4eSBarry Smith 
151f1580f4eSBarry Smith    Implemented by  Henry M. Tufo III and Paul Fischer originally for Nek5000 and called XXT or XYT
152a21bf2d2SBarry Smith 
153*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`
15494c0dd61SBarry Smith M*/
155d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_TFS(PC pc)
156d71ae5a4SJacob Faibussowitsch {
1576218e575SSatish Balay   PC_TFS     *tfs;
1583cd2be91SBarry Smith   PetscMPIInt cmp;
1596218e575SSatish Balay 
1606218e575SSatish Balay   PetscFunctionBegin;
1619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_compare(PETSC_COMM_WORLD, PetscObjectComm((PetscObject)pc), &cmp));
1622472a847SBarry Smith   PetscCheck(cmp == MPI_IDENT || cmp == MPI_CONGRUENT, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "TFS only works with PETSC_COMM_WORLD objects");
1634dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&tfs));
1646218e575SSatish Balay 
1650a545947SLisandro Dalcin   tfs->xxt = NULL;
1660a545947SLisandro Dalcin   tfs->xyt = NULL;
1670a545947SLisandro Dalcin   tfs->b   = NULL;
1680a545947SLisandro Dalcin   tfs->xd  = NULL;
1690a545947SLisandro Dalcin   tfs->xo  = NULL;
1706218e575SSatish Balay   tfs->nd  = 0;
1716218e575SSatish Balay 
1720a545947SLisandro Dalcin   pc->ops->apply               = NULL;
1730a545947SLisandro Dalcin   pc->ops->applytranspose      = NULL;
1746218e575SSatish Balay   pc->ops->setup               = PCSetUp_TFS;
1756218e575SSatish Balay   pc->ops->destroy             = PCDestroy_TFS;
1766218e575SSatish Balay   pc->ops->setfromoptions      = PCSetFromOptions_TFS;
1776218e575SSatish Balay   pc->ops->view                = PCView_TFS;
1780a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
1790a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
1800a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
1816218e575SSatish Balay   pc->data                     = (void *)tfs;
1823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1836218e575SSatish Balay }
184