xref: /petsc/src/ksp/pc/impls/tfs/tfs.c (revision a21bf2d209ea1c1aab2dc2bb0877a3a01be8875a)
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   PetscErrorCode ierr;
206218e575SSatish Balay 
216218e575SSatish Balay   PetscFunctionBegin;
226218e575SSatish Balay   /* free the XXT datastructures */
236218e575SSatish Balay   if (tfs->xxt) {
246218e575SSatish Balay     ierr = XXT_free(tfs->xxt);CHKERRQ(ierr);
256218e575SSatish Balay   }
266218e575SSatish Balay   if (tfs->xyt) {
276218e575SSatish Balay     ierr = XYT_free(tfs->xyt);CHKERRQ(ierr);
286218e575SSatish Balay   }
296bf464f9SBarry Smith   ierr = VecDestroy(&tfs->b);CHKERRQ(ierr);
306bf464f9SBarry Smith   ierr = VecDestroy(&tfs->xd);CHKERRQ(ierr);
316bf464f9SBarry Smith   ierr = VecDestroy(&tfs->xo);CHKERRQ(ierr);
32c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
336218e575SSatish Balay   PetscFunctionReturn(0);
346218e575SSatish Balay }
356218e575SSatish Balay 
366218e575SSatish Balay static PetscErrorCode PCApply_TFS_XXT(PC pc,Vec x,Vec y)
376218e575SSatish Balay {
386218e575SSatish Balay   PC_TFS            *tfs = (PC_TFS*)pc->data;
393468409bSBarry Smith   PetscScalar       *yy;
403468409bSBarry Smith   const PetscScalar *xx;
416218e575SSatish Balay   PetscErrorCode    ierr;
426218e575SSatish Balay 
436218e575SSatish Balay   PetscFunctionBegin;
443468409bSBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
456218e575SSatish Balay   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
463468409bSBarry Smith   ierr = XXT_solve(tfs->xxt,yy,(PetscScalar*)xx);CHKERRQ(ierr);
473468409bSBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
486218e575SSatish Balay   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
496218e575SSatish Balay   PetscFunctionReturn(0);
506218e575SSatish Balay }
516218e575SSatish Balay 
526218e575SSatish Balay static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y)
536218e575SSatish Balay {
546218e575SSatish Balay   PC_TFS            *tfs = (PC_TFS*)pc->data;
553468409bSBarry Smith   PetscScalar       *yy;
563468409bSBarry Smith   const PetscScalar *xx;
576218e575SSatish Balay   PetscErrorCode    ierr;
586218e575SSatish Balay 
596218e575SSatish Balay   PetscFunctionBegin;
603468409bSBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
616218e575SSatish Balay   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
623468409bSBarry Smith   ierr = XYT_solve(tfs->xyt,yy,(PetscScalar*)xx);CHKERRQ(ierr);
633468409bSBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
646218e575SSatish Balay   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
656218e575SSatish Balay   PetscFunctionReturn(0);
666218e575SSatish Balay }
676218e575SSatish Balay 
68fe05c926SBarry Smith static PetscErrorCode PCTFSLocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout)
696218e575SSatish Balay {
706218e575SSatish Balay   PC_TFS         *tfs = (PC_TFS*)pc->data;
716218e575SSatish Balay   Mat            A    = pc->pmat;
726218e575SSatish Balay   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data;
736218e575SSatish Balay   PetscErrorCode ierr;
746218e575SSatish Balay 
756218e575SSatish Balay   PetscFunctionBegin;
766218e575SSatish Balay   ierr = VecPlaceArray(tfs->b,xout);CHKERRQ(ierr);
776218e575SSatish Balay   ierr = VecPlaceArray(tfs->xd,xin);CHKERRQ(ierr);
786218e575SSatish Balay   ierr = VecPlaceArray(tfs->xo,xin+tfs->nd);CHKERRQ(ierr);
796218e575SSatish Balay   ierr = MatMult(a->A,tfs->xd,tfs->b);CHKERRQ(ierr);
806218e575SSatish Balay   ierr = MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b);CHKERRQ(ierr);
81026a3f60SBarry Smith   ierr = VecResetArray(tfs->b);CHKERRQ(ierr);
82026a3f60SBarry Smith   ierr = VecResetArray(tfs->xd);CHKERRQ(ierr);
83026a3f60SBarry Smith   ierr = VecResetArray(tfs->xo);CHKERRQ(ierr);
846218e575SSatish Balay   PetscFunctionReturn(0);
856218e575SSatish Balay }
866218e575SSatish Balay 
876218e575SSatish Balay static PetscErrorCode PCSetUp_TFS(PC pc)
886218e575SSatish Balay {
896218e575SSatish Balay   PC_TFS         *tfs = (PC_TFS*)pc->data;
906218e575SSatish Balay   Mat            A    = pc->pmat;
916218e575SSatish Balay   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data;
926218e575SSatish Balay   PetscErrorCode ierr;
936218e575SSatish Balay   PetscInt       *localtoglobal,ncol,i;
94ace3abfcSBarry Smith   PetscBool      ismpiaij;
956218e575SSatish Balay 
966218e575SSatish Balay   /*
97ace3abfcSBarry Smith   PetscBool      issymmetric;
986218e575SSatish Balay   Petsc Real tol = 0.0;
996218e575SSatish Balay   */
1006218e575SSatish Balay 
1016218e575SSatish Balay   PetscFunctionBegin;
102ce94432eSBarry Smith   if (A->cmap->N != A->rmap->N) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"matrix must be square");
103251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
104ce94432eSBarry Smith   if (!ismpiaij) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices");
1056218e575SSatish Balay 
1066218e575SSatish Balay   /* generate the local to global mapping */
107d0f46423SBarry Smith   ncol = a->A->cmap->n + a->B->cmap->n;
108854ce69bSBarry Smith   ierr = PetscMalloc1(ncol,&localtoglobal);CHKERRQ(ierr);
1092fa5cd67SKarl Rupp   for (i=0; i<a->A->cmap->n; i++) localtoglobal[i] = A->cmap->rstart + i + 1;
1102fa5cd67SKarl Rupp   for (i=0; i<a->B->cmap->n; i++) localtoglobal[i+a->A->cmap->n] = a->garray[i] + 1;
1112fa5cd67SKarl Rupp 
1126218e575SSatish Balay   /* generate the vectors needed for the local solves */
1130298fd71SBarry Smith   ierr    = VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->A->rmap->n,NULL,&tfs->b);CHKERRQ(ierr);
1140298fd71SBarry Smith   ierr    = VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->A->cmap->n,NULL,&tfs->xd);CHKERRQ(ierr);
1150298fd71SBarry Smith   ierr    = VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->B->cmap->n,NULL,&tfs->xo);CHKERRQ(ierr);
116d0f46423SBarry Smith   tfs->nd = a->A->cmap->n;
1176218e575SSatish Balay 
1186218e575SSatish Balay 
1196218e575SSatish Balay   /*  ierr =  MatIsSymmetric(A,tol,&issymmetric); */
1206218e575SSatish Balay   /*  if (issymmetric) { */
121585bbf96SBarry Smith   ierr = PetscBarrier((PetscObject)pc);CHKERRQ(ierr);
1226218e575SSatish Balay   if (A->symmetric) {
1236218e575SSatish Balay     tfs->xxt       = XXT_new();
1245c8f6a95SKarl Rupp     ierr           = XXT_factor(tfs->xxt,localtoglobal,A->rmap->n,ncol,(PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*))PCTFSLocalMult_TFS,pc);CHKERRQ(ierr);
1256218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XXT;
1266218e575SSatish Balay   } else {
1276218e575SSatish Balay     tfs->xyt       = XYT_new();
1285c8f6a95SKarl Rupp     ierr           = XYT_factor(tfs->xyt,localtoglobal,A->rmap->n,ncol,(PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*))PCTFSLocalMult_TFS,pc);CHKERRQ(ierr);
1296218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XYT;
1306218e575SSatish Balay   }
1316218e575SSatish Balay 
1326218e575SSatish Balay   ierr = PetscFree(localtoglobal);CHKERRQ(ierr);
1336218e575SSatish Balay   PetscFunctionReturn(0);
1346218e575SSatish Balay }
1356218e575SSatish Balay 
1364416b707SBarry Smith static PetscErrorCode PCSetFromOptions_TFS(PetscOptionItems *PetscOptionsObject,PC pc)
1376218e575SSatish Balay {
1386218e575SSatish Balay   PetscFunctionBegin;
1396218e575SSatish Balay   PetscFunctionReturn(0);
1406218e575SSatish Balay }
1416218e575SSatish Balay static PetscErrorCode PCView_TFS(PC pc,PetscViewer viewer)
1426218e575SSatish Balay {
1436218e575SSatish Balay   PetscFunctionBegin;
1446218e575SSatish Balay   PetscFunctionReturn(0);
1456218e575SSatish Balay }
1466218e575SSatish Balay 
14794c0dd61SBarry Smith /*MC
14894c0dd61SBarry Smith      PCTFS - A parallel direct solver intended for problems with very few unknowns (like the
149*a21bf2d2SBarry Smith          coarse grid in multigrid). Performs a Cholesky or LU factorization of a matrix defined by
150*a21bf2d2SBarry Smith          its local matrix vector product.
15194c0dd61SBarry Smith 
152*a21bf2d2SBarry Smith    Implemented by  Henry M. Tufo III and Paul Fischer originally for Nek5000 and called XXT or XYT
15394c0dd61SBarry Smith 
15494c0dd61SBarry Smith    Level: beginner
15594c0dd61SBarry Smith 
15695452b02SPatrick Sanan    Notes:
15795452b02SPatrick Sanan     Only implemented for the MPIAIJ matrices
15894c0dd61SBarry Smith 
1597773e740SBarry Smith     Only works on a solver object that lives on all of PETSC_COMM_WORLD!
1607773e740SBarry Smith 
161*a21bf2d2SBarry Smith     Only works for real numbers (is not built if PetscScalar is complex)
162*a21bf2d2SBarry Smith 
16394c0dd61SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
16494c0dd61SBarry Smith M*/
1658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_TFS(PC pc)
1666218e575SSatish Balay {
1676218e575SSatish Balay   PetscErrorCode ierr;
1686218e575SSatish Balay   PC_TFS         *tfs;
1693cd2be91SBarry Smith   PetscMPIInt    cmp;
1706218e575SSatish Balay 
1716218e575SSatish Balay   PetscFunctionBegin;
172ce94432eSBarry Smith   ierr = MPI_Comm_compare(PETSC_COMM_WORLD,PetscObjectComm((PetscObject)pc),&cmp);CHKERRQ(ierr);
173ce94432eSBarry Smith   if (cmp != MPI_IDENT && cmp != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"TFS only works with PETSC_COMM_WORLD objects");
174b00a9115SJed Brown   ierr = PetscNewLog(pc,&tfs);CHKERRQ(ierr);
1756218e575SSatish Balay 
1766218e575SSatish Balay   tfs->xxt = 0;
1776218e575SSatish Balay   tfs->xyt = 0;
1786218e575SSatish Balay   tfs->b   = 0;
1796218e575SSatish Balay   tfs->xd  = 0;
1806218e575SSatish Balay   tfs->xo  = 0;
1816218e575SSatish Balay   tfs->nd  = 0;
1826218e575SSatish Balay 
1836218e575SSatish Balay   pc->ops->apply               = 0;
1846218e575SSatish Balay   pc->ops->applytranspose      = 0;
1856218e575SSatish Balay   pc->ops->setup               = PCSetUp_TFS;
1866218e575SSatish Balay   pc->ops->destroy             = PCDestroy_TFS;
1876218e575SSatish Balay   pc->ops->setfromoptions      = PCSetFromOptions_TFS;
1886218e575SSatish Balay   pc->ops->view                = PCView_TFS;
1896218e575SSatish Balay   pc->ops->applyrichardson     = 0;
1906218e575SSatish Balay   pc->ops->applysymmetricleft  = 0;
1916218e575SSatish Balay   pc->ops->applysymmetricright = 0;
1926218e575SSatish Balay   pc->data                     = (void*)tfs;
1936218e575SSatish Balay   PetscFunctionReturn(0);
1946218e575SSatish Balay }
1956218e575SSatish Balay 
196