xref: /petsc/src/ksp/pc/impls/tfs/tfs.c (revision 585bbf96a41398e6a604d929a8899a97201ec332)
16218e575SSatish Balay /*
26218e575SSatish Balay         Provides an interface to the Tufo-Fischer parallel direct solver
36218e575SSatish Balay */
46218e575SSatish Balay 
56218e575SSatish Balay #include "src/ksp/pc/pcimpl.h"   /*I "petscpc.h" I*/
66218e575SSatish Balay #include "src/mat/impls/aij/mpi/mpiaij.h"
76218e575SSatish Balay 
86218e575SSatish Balay #include "src/ksp/pc/impls/tfs/xxt.h"
96218e575SSatish Balay #include "src/ksp/pc/impls/tfs/xyt.h"
106218e575SSatish Balay 
116218e575SSatish Balay typedef struct {
126218e575SSatish Balay   xxt_ADT  xxt;
136218e575SSatish Balay   xyt_ADT  xyt;
146218e575SSatish Balay   Vec      b,xd,xo;
156218e575SSatish Balay   PetscInt nd;
166218e575SSatish Balay } PC_TFS;
176218e575SSatish Balay 
186218e575SSatish Balay #undef __FUNCT__
196218e575SSatish Balay #define __FUNCT__ "PCDestroy_TFS"
206218e575SSatish Balay PetscErrorCode PCDestroy_TFS(PC pc)
216218e575SSatish Balay {
226218e575SSatish Balay   PC_TFS *tfs = (PC_TFS*)pc->data;
236218e575SSatish Balay   PetscErrorCode ierr;
246218e575SSatish Balay 
256218e575SSatish Balay   PetscFunctionBegin;
266218e575SSatish Balay   /* free the XXT datastructures */
276218e575SSatish Balay   if (tfs->xxt) {
286218e575SSatish Balay     ierr = XXT_free(tfs->xxt);CHKERRQ(ierr);
296218e575SSatish Balay   }
306218e575SSatish Balay   if (tfs->xyt) {
316218e575SSatish Balay     ierr = XYT_free(tfs->xyt);CHKERRQ(ierr);
326218e575SSatish Balay   }
336218e575SSatish Balay   if (tfs->b) {
346218e575SSatish Balay   ierr = VecDestroy(tfs->b);CHKERRQ(ierr);
356218e575SSatish Balay   }
366218e575SSatish Balay   if (tfs->xd) {
376218e575SSatish Balay   ierr = VecDestroy(tfs->xd);CHKERRQ(ierr);
386218e575SSatish Balay   }
396218e575SSatish Balay   if (tfs->xo) {
406218e575SSatish Balay   ierr = VecDestroy(tfs->xo);CHKERRQ(ierr);
416218e575SSatish Balay   }
426218e575SSatish Balay   ierr = PetscFree(tfs);CHKERRQ(ierr);
436218e575SSatish Balay   PetscFunctionReturn(0);
446218e575SSatish Balay }
456218e575SSatish Balay 
466218e575SSatish Balay #undef __FUNCT__
476218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XXT"
486218e575SSatish Balay static PetscErrorCode PCApply_TFS_XXT(PC pc,Vec x,Vec y)
496218e575SSatish Balay {
506218e575SSatish Balay   PC_TFS *tfs = (PC_TFS*)pc->data;
516218e575SSatish Balay   PetscScalar    *xx,*yy;
526218e575SSatish Balay   PetscErrorCode ierr;
536218e575SSatish Balay 
546218e575SSatish Balay   PetscFunctionBegin;
556218e575SSatish Balay   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
566218e575SSatish Balay   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
576218e575SSatish Balay   ierr = XXT_solve(tfs->xxt,yy,xx);CHKERRQ(ierr);
586218e575SSatish Balay   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
596218e575SSatish Balay   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
606218e575SSatish Balay   PetscFunctionReturn(0);
616218e575SSatish Balay }
626218e575SSatish Balay 
636218e575SSatish Balay #undef __FUNCT__
646218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XYT"
656218e575SSatish Balay static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y)
666218e575SSatish Balay {
676218e575SSatish Balay   PC_TFS *tfs = (PC_TFS*)pc->data;
686218e575SSatish Balay   PetscScalar    *xx,*yy;
696218e575SSatish Balay   PetscErrorCode ierr;
706218e575SSatish Balay 
716218e575SSatish Balay   PetscFunctionBegin;
726218e575SSatish Balay   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
736218e575SSatish Balay   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
746218e575SSatish Balay   ierr = XYT_solve(tfs->xyt,yy,xx);CHKERRQ(ierr);
756218e575SSatish Balay   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
766218e575SSatish Balay   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
776218e575SSatish Balay   PetscFunctionReturn(0);
786218e575SSatish Balay }
796218e575SSatish Balay 
806218e575SSatish Balay #undef __FUNCT__
816218e575SSatish Balay #define __FUNCT__ "LocalMult_TFS"
826218e575SSatish Balay static PetscErrorCode LocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout)
836218e575SSatish Balay {
846218e575SSatish Balay   PC_TFS        *tfs = (PC_TFS*)pc->data;
856218e575SSatish Balay   Mat           A = pc->pmat;
866218e575SSatish Balay   Mat_MPIAIJ    *a = (Mat_MPIAIJ*)A->data;
876218e575SSatish Balay   PetscErrorCode ierr;
886218e575SSatish Balay 
896218e575SSatish Balay   PetscFunctionBegin;
906218e575SSatish Balay   ierr = VecPlaceArray(tfs->b,xout);CHKERRQ(ierr);
916218e575SSatish Balay   ierr = VecPlaceArray(tfs->xd,xin);CHKERRQ(ierr);
926218e575SSatish Balay   ierr = VecPlaceArray(tfs->xo,xin+tfs->nd);CHKERRQ(ierr);
936218e575SSatish Balay   ierr = MatMult(a->A,tfs->xd,tfs->b);CHKERRQ(ierr);
946218e575SSatish Balay   ierr = MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b);CHKERRQ(ierr);
956218e575SSatish Balay   PetscFunctionReturn(0);
966218e575SSatish Balay }
976218e575SSatish Balay 
986218e575SSatish Balay #undef __FUNCT__
996218e575SSatish Balay #define __FUNCT__ "PCSetUp_TFS"
1006218e575SSatish Balay static PetscErrorCode PCSetUp_TFS(PC pc)
1016218e575SSatish Balay {
1026218e575SSatish Balay   PC_TFS        *tfs = (PC_TFS*)pc->data;
1036218e575SSatish Balay   Mat            A = pc->pmat;
1046218e575SSatish Balay   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
1056218e575SSatish Balay   PetscErrorCode ierr;
1066218e575SSatish Balay   PetscInt      *localtoglobal,ncol,i;
1076218e575SSatish Balay   PetscTruth     ismpiaij;
1086218e575SSatish Balay 
1096218e575SSatish Balay   /*
1106218e575SSatish Balay   PetscTruth     issymmetric;
1116218e575SSatish Balay   PetscReal tol = 0.0;
1126218e575SSatish Balay   */
1136218e575SSatish Balay 
1146218e575SSatish Balay   PetscFunctionBegin;
1156218e575SSatish Balay   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
1166218e575SSatish Balay   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1176218e575SSatish Balay   if (!ismpiaij) {
1186218e575SSatish Balay     SETERRQ(PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices");
1196218e575SSatish Balay   }
1206218e575SSatish Balay 
1216218e575SSatish Balay   /* generate the local to global mapping */
1226218e575SSatish Balay   ncol = a->A->n + a->B->n;
1236218e575SSatish Balay   ierr = PetscMalloc((ncol)*sizeof(int),&localtoglobal);CHKERRQ(ierr);
1246218e575SSatish Balay   for (i=0; i<a->A->n; i++) {
1256218e575SSatish Balay     localtoglobal[i] = a->cstart + i + 1;
1266218e575SSatish Balay   }
1276218e575SSatish Balay   for (i=0; i<a->B->n; i++) {
1286218e575SSatish Balay     localtoglobal[i+a->A->n] = a->garray[i] + 1;
1296218e575SSatish Balay   }
1306218e575SSatish Balay   /* generate the vectors needed for the local solves */
1316218e575SSatish Balay   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&tfs->b);CHKERRQ(ierr);
1326218e575SSatish Balay   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&tfs->xd);CHKERRQ(ierr);
1336218e575SSatish Balay   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&tfs->xo);CHKERRQ(ierr);
1346218e575SSatish Balay   tfs->nd = a->A->n;
1356218e575SSatish Balay 
1366218e575SSatish Balay 
1376218e575SSatish Balay   /*  ierr =  MatIsSymmetric(A,tol,&issymmetric); */
1386218e575SSatish Balay   /*  if (issymmetric) { */
139*585bbf96SBarry Smith   ierr = PetscBarrier((PetscObject)pc);CHKERRQ(ierr);
1406218e575SSatish Balay   if (A->symmetric) {
1416218e575SSatish Balay     tfs->xxt       = XXT_new();
1426218e575SSatish Balay     ierr           = XXT_factor(tfs->xxt,localtoglobal,A->m,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr);
1436218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XXT;
1446218e575SSatish Balay   } else {
1456218e575SSatish Balay     tfs->xyt       = XYT_new();
1466218e575SSatish Balay     ierr           = XYT_factor(tfs->xyt,localtoglobal,A->m,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr);
1476218e575SSatish Balay     pc->ops->apply = PCApply_TFS_XYT;
1486218e575SSatish Balay   }
1496218e575SSatish Balay 
1506218e575SSatish Balay   ierr = PetscFree(localtoglobal);CHKERRQ(ierr);
1516218e575SSatish Balay   PetscFunctionReturn(0);
1526218e575SSatish Balay }
1536218e575SSatish Balay 
1546218e575SSatish Balay #undef __FUNCT__
1556218e575SSatish Balay #define __FUNCT__ "PCSetFromOptions_TFS"
1566218e575SSatish Balay static PetscErrorCode PCSetFromOptions_TFS(PC pc)
1576218e575SSatish Balay {
1586218e575SSatish Balay   PetscFunctionBegin;
1596218e575SSatish Balay   PetscFunctionReturn(0);
1606218e575SSatish Balay }
1616218e575SSatish Balay #undef __FUNCT__
1626218e575SSatish Balay #define __FUNCT__ "PCView_TFS"
1636218e575SSatish Balay static PetscErrorCode PCView_TFS(PC pc,PetscViewer viewer)
1646218e575SSatish Balay {
1656218e575SSatish Balay   PetscFunctionBegin;
1666218e575SSatish Balay   PetscFunctionReturn(0);
1676218e575SSatish Balay }
1686218e575SSatish Balay 
1696218e575SSatish Balay EXTERN_C_BEGIN
1706218e575SSatish Balay #undef __FUNCT__
1716218e575SSatish Balay #define __FUNCT__ "PCCreate_TFS"
1726218e575SSatish Balay PetscErrorCode PCCreate_TFS(PC pc)
1736218e575SSatish Balay {
1746218e575SSatish Balay   PetscErrorCode ierr;
1756218e575SSatish Balay   PC_TFS         *tfs;
1766218e575SSatish Balay 
1776218e575SSatish Balay   PetscFunctionBegin;
1786218e575SSatish Balay   ierr = PetscNew(PC_TFS,&tfs);CHKERRQ(ierr);
1796218e575SSatish Balay   PetscLogObjectMemory(pc,sizeof(PC_TFS));
1806218e575SSatish Balay 
1816218e575SSatish Balay   tfs->xxt = 0;
1826218e575SSatish Balay   tfs->xyt = 0;
1836218e575SSatish Balay   tfs->b   = 0;
1846218e575SSatish Balay   tfs->xd  = 0;
1856218e575SSatish Balay   tfs->xo  = 0;
1866218e575SSatish Balay   tfs->nd  = 0;
1876218e575SSatish Balay 
1886218e575SSatish Balay   pc->ops->apply               = 0;
1896218e575SSatish Balay   pc->ops->applytranspose      = 0;
1906218e575SSatish Balay   pc->ops->setup               = PCSetUp_TFS;
1916218e575SSatish Balay   pc->ops->destroy             = PCDestroy_TFS;
1926218e575SSatish Balay   pc->ops->setfromoptions      = PCSetFromOptions_TFS;
1936218e575SSatish Balay   pc->ops->view                = PCView_TFS;
1946218e575SSatish Balay   pc->ops->applyrichardson     = 0;
1956218e575SSatish Balay   pc->ops->applysymmetricleft  = 0;
1966218e575SSatish Balay   pc->ops->applysymmetricright = 0;
1976218e575SSatish Balay   pc->data                     = (void*)tfs;
1986218e575SSatish Balay   PetscFunctionReturn(0);
1996218e575SSatish Balay }
2006218e575SSatish Balay EXTERN_C_END
2016218e575SSatish Balay 
202