xref: /petsc/src/ksp/pc/impls/tfs/tfs.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 
166218e575SSatish Balay PetscErrorCode PCDestroy_TFS(PC pc)
176218e575SSatish Balay {
186218e575SSatish Balay   PC_TFS         *tfs = (PC_TFS*)pc->data;
196218e575SSatish Balay 
206218e575SSatish Balay   PetscFunctionBegin;
216218e575SSatish Balay   /* free the XXT datastructures */
226218e575SSatish Balay   if (tfs->xxt) {
235f80ce2aSJacob Faibussowitsch     CHKERRQ(XXT_free(tfs->xxt));
246218e575SSatish Balay   }
256218e575SSatish Balay   if (tfs->xyt) {
265f80ce2aSJacob Faibussowitsch     CHKERRQ(XYT_free(tfs->xyt));
276218e575SSatish Balay   }
285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&tfs->b));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&tfs->xd));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&tfs->xo));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pc->data));
326218e575SSatish Balay   PetscFunctionReturn(0);
336218e575SSatish Balay }
346218e575SSatish Balay 
356218e575SSatish Balay static PetscErrorCode PCApply_TFS_XXT(PC pc,Vec x,Vec y)
366218e575SSatish Balay {
376218e575SSatish Balay   PC_TFS            *tfs = (PC_TFS*)pc->data;
383468409bSBarry Smith   PetscScalar       *yy;
393468409bSBarry Smith   const PetscScalar *xx;
406218e575SSatish Balay 
416218e575SSatish Balay   PetscFunctionBegin;
425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(y,&yy));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(XXT_solve(tfs->xxt,yy,(PetscScalar*)xx));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(y,&yy));
476218e575SSatish Balay   PetscFunctionReturn(0);
486218e575SSatish Balay }
496218e575SSatish Balay 
506218e575SSatish Balay static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y)
516218e575SSatish Balay {
526218e575SSatish Balay   PC_TFS            *tfs = (PC_TFS*)pc->data;
533468409bSBarry Smith   PetscScalar       *yy;
543468409bSBarry Smith   const PetscScalar *xx;
556218e575SSatish Balay 
566218e575SSatish Balay   PetscFunctionBegin;
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xx));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(y,&yy));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(XYT_solve(tfs->xyt,yy,(PetscScalar*)xx));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xx));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(y,&yy));
626218e575SSatish Balay   PetscFunctionReturn(0);
636218e575SSatish Balay }
646218e575SSatish Balay 
65fe05c926SBarry Smith static PetscErrorCode PCTFSLocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout)
666218e575SSatish Balay {
676218e575SSatish Balay   PC_TFS         *tfs = (PC_TFS*)pc->data;
686218e575SSatish Balay   Mat            A    = pc->pmat;
696218e575SSatish Balay   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data;
706218e575SSatish Balay 
716218e575SSatish Balay   PetscFunctionBegin;
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPlaceArray(tfs->b,xout));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPlaceArray(tfs->xd,xin));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPlaceArray(tfs->xo,xin+tfs->nd));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(a->A,tfs->xd,tfs->b));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecResetArray(tfs->b));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecResetArray(tfs->xd));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecResetArray(tfs->xo));
806218e575SSatish Balay   PetscFunctionReturn(0);
816218e575SSatish Balay }
826218e575SSatish Balay 
836218e575SSatish Balay static PetscErrorCode PCSetUp_TFS(PC pc)
846218e575SSatish Balay {
856218e575SSatish Balay   PC_TFS         *tfs = (PC_TFS*)pc->data;
866218e575SSatish Balay   Mat            A    = pc->pmat;
876218e575SSatish Balay   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data;
886218e575SSatish Balay   PetscInt       *localtoglobal,ncol,i;
89ace3abfcSBarry Smith   PetscBool      ismpiaij;
906218e575SSatish Balay 
916218e575SSatish Balay   /*
92ace3abfcSBarry Smith   PetscBool      issymmetric;
936218e575SSatish Balay   Petsc Real tol = 0.0;
946218e575SSatish Balay   */
956218e575SSatish Balay 
966218e575SSatish Balay   PetscFunctionBegin;
972c71b3e2SJacob Faibussowitsch   PetscCheckFalse(A->cmap->N != A->rmap->N,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"matrix must be square");
985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij));
99*28b400f6SJacob Faibussowitsch   PetscCheck(ismpiaij,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices");
1006218e575SSatish Balay 
1016218e575SSatish Balay   /* generate the local to global mapping */
102d0f46423SBarry Smith   ncol = a->A->cmap->n + a->B->cmap->n;
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ncol,&localtoglobal));
1042fa5cd67SKarl Rupp   for (i=0; i<a->A->cmap->n; i++) localtoglobal[i] = A->cmap->rstart + i + 1;
1052fa5cd67SKarl Rupp   for (i=0; i<a->B->cmap->n; i++) localtoglobal[i+a->A->cmap->n] = a->garray[i] + 1;
1062fa5cd67SKarl Rupp 
1076218e575SSatish Balay   /* generate the vectors needed for the local solves */
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->A->rmap->n,NULL,&tfs->b));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->A->cmap->n,NULL,&tfs->xd));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->B->cmap->n,NULL,&tfs->xo));
111d0f46423SBarry Smith   tfs->nd = a->A->cmap->n;
1126218e575SSatish Balay 
1136218e575SSatish Balay   /*  ierr =  MatIsSymmetric(A,tol,&issymmetric); */
1146218e575SSatish Balay   /*  if (issymmetric) { */
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBarrier((PetscObject)pc));
1166218e575SSatish Balay   if (A->symmetric) {
1176218e575SSatish Balay     tfs->xxt       = XXT_new();
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(XXT_factor(tfs->xxt,localtoglobal,A->rmap->n,ncol,(PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*))PCTFSLocalMult_TFS,pc));
1196218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XXT;
1206218e575SSatish Balay   } else {
1216218e575SSatish Balay     tfs->xyt       = XYT_new();
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(XYT_factor(tfs->xyt,localtoglobal,A->rmap->n,ncol,(PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*))PCTFSLocalMult_TFS,pc));
1236218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XYT;
1246218e575SSatish Balay   }
1256218e575SSatish Balay 
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(localtoglobal));
1276218e575SSatish Balay   PetscFunctionReturn(0);
1286218e575SSatish Balay }
1296218e575SSatish Balay 
1304416b707SBarry Smith static PetscErrorCode PCSetFromOptions_TFS(PetscOptionItems *PetscOptionsObject,PC pc)
1316218e575SSatish Balay {
1326218e575SSatish Balay   PetscFunctionBegin;
1336218e575SSatish Balay   PetscFunctionReturn(0);
1346218e575SSatish Balay }
1356218e575SSatish Balay static PetscErrorCode PCView_TFS(PC pc,PetscViewer viewer)
1366218e575SSatish Balay {
1376218e575SSatish Balay   PetscFunctionBegin;
1386218e575SSatish Balay   PetscFunctionReturn(0);
1396218e575SSatish Balay }
1406218e575SSatish Balay 
14194c0dd61SBarry Smith /*MC
14294c0dd61SBarry Smith      PCTFS - A parallel direct solver intended for problems with very few unknowns (like the
143a21bf2d2SBarry Smith          coarse grid in multigrid). Performs a Cholesky or LU factorization of a matrix defined by
144a21bf2d2SBarry Smith          its local matrix vector product.
14594c0dd61SBarry Smith 
146a21bf2d2SBarry Smith    Implemented by  Henry M. Tufo III and Paul Fischer originally for Nek5000 and called XXT or XYT
14794c0dd61SBarry Smith 
14894c0dd61SBarry Smith    Level: beginner
14994c0dd61SBarry Smith 
15095452b02SPatrick Sanan    Notes:
15195452b02SPatrick Sanan     Only implemented for the MPIAIJ matrices
15294c0dd61SBarry Smith 
1537773e740SBarry Smith     Only works on a solver object that lives on all of PETSC_COMM_WORLD!
1547773e740SBarry Smith 
155a21bf2d2SBarry Smith     Only works for real numbers (is not built if PetscScalar is complex)
156a21bf2d2SBarry Smith 
15794c0dd61SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
15894c0dd61SBarry Smith M*/
1598cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_TFS(PC pc)
1606218e575SSatish Balay {
1616218e575SSatish Balay   PC_TFS         *tfs;
1623cd2be91SBarry Smith   PetscMPIInt    cmp;
1636218e575SSatish Balay 
1646218e575SSatish Balay   PetscFunctionBegin;
1655f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_compare(PETSC_COMM_WORLD,PetscObjectComm((PetscObject)pc),&cmp));
1662c71b3e2SJacob Faibussowitsch   PetscCheckFalse(cmp != MPI_IDENT && cmp != MPI_CONGRUENT,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"TFS only works with PETSC_COMM_WORLD objects");
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(pc,&tfs));
1686218e575SSatish Balay 
1690a545947SLisandro Dalcin   tfs->xxt = NULL;
1700a545947SLisandro Dalcin   tfs->xyt = NULL;
1710a545947SLisandro Dalcin   tfs->b   = NULL;
1720a545947SLisandro Dalcin   tfs->xd  = NULL;
1730a545947SLisandro Dalcin   tfs->xo  = NULL;
1746218e575SSatish Balay   tfs->nd  = 0;
1756218e575SSatish Balay 
1760a545947SLisandro Dalcin   pc->ops->apply               = NULL;
1770a545947SLisandro Dalcin   pc->ops->applytranspose      = NULL;
1786218e575SSatish Balay   pc->ops->setup               = PCSetUp_TFS;
1796218e575SSatish Balay   pc->ops->destroy             = PCDestroy_TFS;
1806218e575SSatish Balay   pc->ops->setfromoptions      = PCSetFromOptions_TFS;
1816218e575SSatish Balay   pc->ops->view                = PCView_TFS;
1820a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
1830a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
1840a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
1856218e575SSatish Balay   pc->data                     = (void*)tfs;
1866218e575SSatish Balay   PetscFunctionReturn(0);
1876218e575SSatish Balay }
188