xref: /petsc/src/ksp/pc/impls/tfs/tfs.c (revision 7773e740a59dc995f6fb3ee2f4d2bd5cd7cda562)
16218e575SSatish Balay /*
26218e575SSatish Balay         Provides an interface to the Tufo-Fischer parallel direct solver
36218e575SSatish Balay */
46218e575SSatish Balay 
5b45d2f2cSJed Brown #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 #undef __FUNCT__
176218e575SSatish Balay #define __FUNCT__ "PCDestroy_TFS"
186218e575SSatish Balay PetscErrorCode PCDestroy_TFS(PC pc)
196218e575SSatish Balay {
206218e575SSatish Balay   PC_TFS *tfs = (PC_TFS*)pc->data;
216218e575SSatish Balay   PetscErrorCode ierr;
226218e575SSatish Balay 
236218e575SSatish Balay   PetscFunctionBegin;
246218e575SSatish Balay   /* free the XXT datastructures */
256218e575SSatish Balay   if (tfs->xxt) {
266218e575SSatish Balay     ierr = XXT_free(tfs->xxt);CHKERRQ(ierr);
276218e575SSatish Balay   }
286218e575SSatish Balay   if (tfs->xyt) {
296218e575SSatish Balay     ierr = XYT_free(tfs->xyt);CHKERRQ(ierr);
306218e575SSatish Balay   }
316bf464f9SBarry Smith   ierr = VecDestroy(&tfs->b);CHKERRQ(ierr);
326bf464f9SBarry Smith   ierr = VecDestroy(&tfs->xd);CHKERRQ(ierr);
336bf464f9SBarry Smith   ierr = VecDestroy(&tfs->xo);CHKERRQ(ierr);
34c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
356218e575SSatish Balay   PetscFunctionReturn(0);
366218e575SSatish Balay }
376218e575SSatish Balay 
386218e575SSatish Balay #undef __FUNCT__
396218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XXT"
406218e575SSatish Balay static PetscErrorCode PCApply_TFS_XXT(PC pc,Vec x,Vec y)
416218e575SSatish Balay {
426218e575SSatish Balay   PC_TFS *tfs = (PC_TFS*)pc->data;
436218e575SSatish Balay   PetscScalar    *xx,*yy;
446218e575SSatish Balay   PetscErrorCode ierr;
456218e575SSatish Balay 
466218e575SSatish Balay   PetscFunctionBegin;
476218e575SSatish Balay   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
486218e575SSatish Balay   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
496218e575SSatish Balay   ierr = XXT_solve(tfs->xxt,yy,xx);CHKERRQ(ierr);
506218e575SSatish Balay   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
516218e575SSatish Balay   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
526218e575SSatish Balay   PetscFunctionReturn(0);
536218e575SSatish Balay }
546218e575SSatish Balay 
556218e575SSatish Balay #undef __FUNCT__
566218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XYT"
576218e575SSatish Balay static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y)
586218e575SSatish Balay {
596218e575SSatish Balay   PC_TFS *tfs = (PC_TFS*)pc->data;
606218e575SSatish Balay   PetscScalar    *xx,*yy;
616218e575SSatish Balay   PetscErrorCode ierr;
626218e575SSatish Balay 
636218e575SSatish Balay   PetscFunctionBegin;
646218e575SSatish Balay   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
656218e575SSatish Balay   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
666218e575SSatish Balay   ierr = XYT_solve(tfs->xyt,yy,xx);CHKERRQ(ierr);
676218e575SSatish Balay   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
686218e575SSatish Balay   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
696218e575SSatish Balay   PetscFunctionReturn(0);
706218e575SSatish Balay }
716218e575SSatish Balay 
726218e575SSatish Balay #undef __FUNCT__
736218e575SSatish Balay #define __FUNCT__ "LocalMult_TFS"
746218e575SSatish Balay static PetscErrorCode LocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout)
756218e575SSatish Balay {
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   PetscErrorCode ierr;
806218e575SSatish Balay 
816218e575SSatish Balay   PetscFunctionBegin;
826218e575SSatish Balay   ierr = VecPlaceArray(tfs->b,xout);CHKERRQ(ierr);
836218e575SSatish Balay   ierr = VecPlaceArray(tfs->xd,xin);CHKERRQ(ierr);
846218e575SSatish Balay   ierr = VecPlaceArray(tfs->xo,xin+tfs->nd);CHKERRQ(ierr);
856218e575SSatish Balay   ierr = MatMult(a->A,tfs->xd,tfs->b);CHKERRQ(ierr);
866218e575SSatish Balay   ierr = MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b);CHKERRQ(ierr);
87026a3f60SBarry Smith   ierr = VecResetArray(tfs->b);CHKERRQ(ierr);
88026a3f60SBarry Smith   ierr = VecResetArray(tfs->xd);CHKERRQ(ierr);
89026a3f60SBarry Smith   ierr = VecResetArray(tfs->xo);CHKERRQ(ierr);
906218e575SSatish Balay   PetscFunctionReturn(0);
916218e575SSatish Balay }
926218e575SSatish Balay 
936218e575SSatish Balay #undef __FUNCT__
946218e575SSatish Balay #define __FUNCT__ "PCSetUp_TFS"
956218e575SSatish Balay static PetscErrorCode PCSetUp_TFS(PC pc)
966218e575SSatish Balay {
976218e575SSatish Balay   PC_TFS        *tfs = (PC_TFS*)pc->data;
986218e575SSatish Balay   Mat            A = pc->pmat;
996218e575SSatish Balay   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
1006218e575SSatish Balay   PetscErrorCode ierr;
1016218e575SSatish Balay   PetscInt      *localtoglobal,ncol,i;
102ace3abfcSBarry Smith   PetscBool      ismpiaij;
1036218e575SSatish Balay 
1046218e575SSatish Balay   /*
105ace3abfcSBarry Smith   PetscBool      issymmetric;
1066218e575SSatish Balay   Petsc Real tol = 0.0;
1076218e575SSatish Balay   */
1086218e575SSatish Balay 
1096218e575SSatish Balay   PetscFunctionBegin;
110e7e72b3dSBarry Smith   if (A->cmap->N != A->rmap->N) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_SIZ,"matrix must be square");
1116218e575SSatish Balay   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
112e7e72b3dSBarry Smith   if (!ismpiaij) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices");
1136218e575SSatish Balay 
1146218e575SSatish Balay   /* generate the local to global mapping */
115d0f46423SBarry Smith   ncol = a->A->cmap->n + a->B->cmap->n;
11652f87cdaSBarry Smith   ierr = PetscMalloc((ncol)*sizeof(PetscInt),&localtoglobal);CHKERRQ(ierr);
117d0f46423SBarry Smith   for (i=0; i<a->A->cmap->n; i++) {
118d0f46423SBarry Smith     localtoglobal[i] = A->cmap->rstart + i + 1;
1196218e575SSatish Balay   }
120d0f46423SBarry Smith   for (i=0; i<a->B->cmap->n; i++) {
121d0f46423SBarry Smith     localtoglobal[i+a->A->cmap->n] = a->garray[i] + 1;
1226218e575SSatish Balay   }
1236218e575SSatish Balay   /* generate the vectors needed for the local solves */
124d0f46423SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->rmap->n,PETSC_NULL,&tfs->b);CHKERRQ(ierr);
125d0f46423SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->cmap->n,PETSC_NULL,&tfs->xd);CHKERRQ(ierr);
126d0f46423SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->cmap->n,PETSC_NULL,&tfs->xo);CHKERRQ(ierr);
127d0f46423SBarry Smith   tfs->nd = a->A->cmap->n;
1286218e575SSatish Balay 
1296218e575SSatish Balay 
1306218e575SSatish Balay   /*  ierr =  MatIsSymmetric(A,tol,&issymmetric); */
1316218e575SSatish Balay   /*  if (issymmetric) { */
132585bbf96SBarry Smith   ierr = PetscBarrier((PetscObject)pc);CHKERRQ(ierr);
1336218e575SSatish Balay   if (A->symmetric) {
1346218e575SSatish Balay     tfs->xxt       = XXT_new();
135d0f46423SBarry Smith     ierr           = XXT_factor(tfs->xxt,localtoglobal,A->rmap->n,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr);
1366218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XXT;
1376218e575SSatish Balay   } else {
1386218e575SSatish Balay     tfs->xyt       = XYT_new();
139d0f46423SBarry Smith     ierr           = XYT_factor(tfs->xyt,localtoglobal,A->rmap->n,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr);
1406218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XYT;
1416218e575SSatish Balay   }
1426218e575SSatish Balay 
1436218e575SSatish Balay   ierr = PetscFree(localtoglobal);CHKERRQ(ierr);
1446218e575SSatish Balay   PetscFunctionReturn(0);
1456218e575SSatish Balay }
1466218e575SSatish Balay 
1476218e575SSatish Balay #undef __FUNCT__
1486218e575SSatish Balay #define __FUNCT__ "PCSetFromOptions_TFS"
1496218e575SSatish Balay static PetscErrorCode PCSetFromOptions_TFS(PC pc)
1506218e575SSatish Balay {
1516218e575SSatish Balay   PetscFunctionBegin;
1526218e575SSatish Balay   PetscFunctionReturn(0);
1536218e575SSatish Balay }
1546218e575SSatish Balay #undef __FUNCT__
1556218e575SSatish Balay #define __FUNCT__ "PCView_TFS"
1566218e575SSatish Balay static PetscErrorCode PCView_TFS(PC pc,PetscViewer viewer)
1576218e575SSatish Balay {
1586218e575SSatish Balay   PetscFunctionBegin;
1596218e575SSatish Balay   PetscFunctionReturn(0);
1606218e575SSatish Balay }
1616218e575SSatish Balay 
1626218e575SSatish Balay EXTERN_C_BEGIN
1636218e575SSatish Balay #undef __FUNCT__
1646218e575SSatish Balay #define __FUNCT__ "PCCreate_TFS"
16594c0dd61SBarry Smith /*MC
16694c0dd61SBarry Smith      PCTFS - A parallel direct solver intended for problems with very few unknowns (like the
16794c0dd61SBarry Smith          coarse grid in multigrid).
16894c0dd61SBarry Smith 
16994c0dd61SBarry Smith    Implemented by  Henry M. Tufo III and Paul Fischer
17094c0dd61SBarry Smith 
17194c0dd61SBarry Smith    Level: beginner
17294c0dd61SBarry Smith 
17394c0dd61SBarry Smith    Notes: Only implemented for the MPIAIJ matrices
17494c0dd61SBarry Smith 
175*7773e740SBarry Smith           Only works on a solver object that lives on all of PETSC_COMM_WORLD!
176*7773e740SBarry Smith 
17794c0dd61SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
17894c0dd61SBarry Smith M*/
1797087cfbeSBarry Smith PetscErrorCode  PCCreate_TFS(PC pc)
1806218e575SSatish Balay {
1816218e575SSatish Balay   PetscErrorCode ierr;
1826218e575SSatish Balay   PC_TFS         *tfs;
1836218e575SSatish Balay 
1846218e575SSatish Balay   PetscFunctionBegin;
18538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_TFS,&tfs);CHKERRQ(ierr);
1866218e575SSatish Balay 
1876218e575SSatish Balay   tfs->xxt = 0;
1886218e575SSatish Balay   tfs->xyt = 0;
1896218e575SSatish Balay   tfs->b   = 0;
1906218e575SSatish Balay   tfs->xd  = 0;
1916218e575SSatish Balay   tfs->xo  = 0;
1926218e575SSatish Balay   tfs->nd  = 0;
1936218e575SSatish Balay 
1946218e575SSatish Balay   pc->ops->apply               = 0;
1956218e575SSatish Balay   pc->ops->applytranspose      = 0;
1966218e575SSatish Balay   pc->ops->setup               = PCSetUp_TFS;
1976218e575SSatish Balay   pc->ops->destroy             = PCDestroy_TFS;
1986218e575SSatish Balay   pc->ops->setfromoptions      = PCSetFromOptions_TFS;
1996218e575SSatish Balay   pc->ops->view                = PCView_TFS;
2006218e575SSatish Balay   pc->ops->applyrichardson     = 0;
2016218e575SSatish Balay   pc->ops->applysymmetricleft  = 0;
2026218e575SSatish Balay   pc->ops->applysymmetricright = 0;
2036218e575SSatish Balay   pc->data                     = (void*)tfs;
2046218e575SSatish Balay   PetscFunctionReturn(0);
2056218e575SSatish Balay }
2066218e575SSatish Balay EXTERN_C_END
2076218e575SSatish Balay 
208