1dba47a55SKris Buschelman #define PETSCKSP_DLL 26218e575SSatish Balay /* 36218e575SSatish Balay Provides an interface to the Tufo-Fischer parallel direct solver 46218e575SSatish Balay */ 56218e575SSatish Balay 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 77c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 87c4f633dSBarry Smith #include "../src/ksp/pc/impls/tfs/tfs.h" 96218e575SSatish Balay 106218e575SSatish Balay typedef struct { 116218e575SSatish Balay xxt_ADT xxt; 126218e575SSatish Balay xyt_ADT xyt; 136218e575SSatish Balay Vec b,xd,xo; 146218e575SSatish Balay PetscInt nd; 156218e575SSatish Balay } PC_TFS; 166218e575SSatish Balay 176218e575SSatish Balay #undef __FUNCT__ 186218e575SSatish Balay #define __FUNCT__ "PCDestroy_TFS" 196218e575SSatish Balay PetscErrorCode PCDestroy_TFS(PC pc) 206218e575SSatish Balay { 216218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 226218e575SSatish Balay PetscErrorCode ierr; 236218e575SSatish Balay 246218e575SSatish Balay PetscFunctionBegin; 256218e575SSatish Balay /* free the XXT datastructures */ 266218e575SSatish Balay if (tfs->xxt) { 276218e575SSatish Balay ierr = XXT_free(tfs->xxt);CHKERRQ(ierr); 286218e575SSatish Balay } 296218e575SSatish Balay if (tfs->xyt) { 306218e575SSatish Balay ierr = XYT_free(tfs->xyt);CHKERRQ(ierr); 316218e575SSatish Balay } 326218e575SSatish Balay if (tfs->b) { 336218e575SSatish Balay ierr = VecDestroy(tfs->b);CHKERRQ(ierr); 346218e575SSatish Balay } 356218e575SSatish Balay if (tfs->xd) { 366218e575SSatish Balay ierr = VecDestroy(tfs->xd);CHKERRQ(ierr); 376218e575SSatish Balay } 386218e575SSatish Balay if (tfs->xo) { 396218e575SSatish Balay ierr = VecDestroy(tfs->xo);CHKERRQ(ierr); 406218e575SSatish Balay } 416218e575SSatish Balay ierr = PetscFree(tfs);CHKERRQ(ierr); 426218e575SSatish Balay PetscFunctionReturn(0); 436218e575SSatish Balay } 446218e575SSatish Balay 456218e575SSatish Balay #undef __FUNCT__ 466218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XXT" 476218e575SSatish Balay static PetscErrorCode PCApply_TFS_XXT(PC pc,Vec x,Vec y) 486218e575SSatish Balay { 496218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 506218e575SSatish Balay PetscScalar *xx,*yy; 516218e575SSatish Balay PetscErrorCode ierr; 526218e575SSatish Balay 536218e575SSatish Balay PetscFunctionBegin; 546218e575SSatish Balay ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 556218e575SSatish Balay ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 566218e575SSatish Balay ierr = XXT_solve(tfs->xxt,yy,xx);CHKERRQ(ierr); 576218e575SSatish Balay ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 586218e575SSatish Balay ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 596218e575SSatish Balay PetscFunctionReturn(0); 606218e575SSatish Balay } 616218e575SSatish Balay 626218e575SSatish Balay #undef __FUNCT__ 636218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XYT" 646218e575SSatish Balay static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y) 656218e575SSatish Balay { 666218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 676218e575SSatish Balay PetscScalar *xx,*yy; 686218e575SSatish Balay PetscErrorCode ierr; 696218e575SSatish Balay 706218e575SSatish Balay PetscFunctionBegin; 716218e575SSatish Balay ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 726218e575SSatish Balay ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 736218e575SSatish Balay ierr = XYT_solve(tfs->xyt,yy,xx);CHKERRQ(ierr); 746218e575SSatish Balay ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 756218e575SSatish Balay ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 766218e575SSatish Balay PetscFunctionReturn(0); 776218e575SSatish Balay } 786218e575SSatish Balay 796218e575SSatish Balay #undef __FUNCT__ 806218e575SSatish Balay #define __FUNCT__ "LocalMult_TFS" 816218e575SSatish Balay static PetscErrorCode LocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout) 826218e575SSatish Balay { 836218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 846218e575SSatish Balay Mat A = pc->pmat; 856218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 866218e575SSatish Balay PetscErrorCode ierr; 876218e575SSatish Balay 886218e575SSatish Balay PetscFunctionBegin; 896218e575SSatish Balay ierr = VecPlaceArray(tfs->b,xout);CHKERRQ(ierr); 906218e575SSatish Balay ierr = VecPlaceArray(tfs->xd,xin);CHKERRQ(ierr); 916218e575SSatish Balay ierr = VecPlaceArray(tfs->xo,xin+tfs->nd);CHKERRQ(ierr); 926218e575SSatish Balay ierr = MatMult(a->A,tfs->xd,tfs->b);CHKERRQ(ierr); 936218e575SSatish Balay ierr = MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b);CHKERRQ(ierr); 94026a3f60SBarry Smith ierr = VecResetArray(tfs->b);CHKERRQ(ierr); 95026a3f60SBarry Smith ierr = VecResetArray(tfs->xd);CHKERRQ(ierr); 96026a3f60SBarry Smith ierr = VecResetArray(tfs->xo);CHKERRQ(ierr); 976218e575SSatish Balay PetscFunctionReturn(0); 986218e575SSatish Balay } 996218e575SSatish Balay 1006218e575SSatish Balay #undef __FUNCT__ 1016218e575SSatish Balay #define __FUNCT__ "PCSetUp_TFS" 1026218e575SSatish Balay static PetscErrorCode PCSetUp_TFS(PC pc) 1036218e575SSatish Balay { 1046218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 1056218e575SSatish Balay Mat A = pc->pmat; 1066218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 1076218e575SSatish Balay PetscErrorCode ierr; 1086218e575SSatish Balay PetscInt *localtoglobal,ncol,i; 1096218e575SSatish Balay PetscTruth ismpiaij; 1106218e575SSatish Balay 1116218e575SSatish Balay /* 1126218e575SSatish Balay PetscTruth issymmetric; 1136218e575SSatish Balay Petsc Real tol = 0.0; 1146218e575SSatish Balay */ 1156218e575SSatish Balay 1166218e575SSatish Balay PetscFunctionBegin; 117*e7e72b3dSBarry Smith if (A->cmap->N != A->rmap->N) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_SIZ,"matrix must be square"); 1186218e575SSatish Balay ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 119*e7e72b3dSBarry Smith if (!ismpiaij) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices"); 1206218e575SSatish Balay 1216218e575SSatish Balay /* generate the local to global mapping */ 122d0f46423SBarry Smith ncol = a->A->cmap->n + a->B->cmap->n; 12352f87cdaSBarry Smith ierr = PetscMalloc((ncol)*sizeof(PetscInt),&localtoglobal);CHKERRQ(ierr); 124d0f46423SBarry Smith for (i=0; i<a->A->cmap->n; i++) { 125d0f46423SBarry Smith localtoglobal[i] = A->cmap->rstart + i + 1; 1266218e575SSatish Balay } 127d0f46423SBarry Smith for (i=0; i<a->B->cmap->n; i++) { 128d0f46423SBarry Smith localtoglobal[i+a->A->cmap->n] = a->garray[i] + 1; 1296218e575SSatish Balay } 1306218e575SSatish Balay /* generate the vectors needed for the local solves */ 131d0f46423SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->rmap->n,PETSC_NULL,&tfs->b);CHKERRQ(ierr); 132d0f46423SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->cmap->n,PETSC_NULL,&tfs->xd);CHKERRQ(ierr); 133d0f46423SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->cmap->n,PETSC_NULL,&tfs->xo);CHKERRQ(ierr); 134d0f46423SBarry Smith tfs->nd = a->A->cmap->n; 1356218e575SSatish Balay 1366218e575SSatish Balay 1376218e575SSatish Balay /* ierr = MatIsSymmetric(A,tol,&issymmetric); */ 1386218e575SSatish Balay /* if (issymmetric) { */ 139585bbf96SBarry Smith ierr = PetscBarrier((PetscObject)pc);CHKERRQ(ierr); 1406218e575SSatish Balay if (A->symmetric) { 1416218e575SSatish Balay tfs->xxt = XXT_new(); 142d0f46423SBarry Smith ierr = XXT_factor(tfs->xxt,localtoglobal,A->rmap->n,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr); 1436218e575SSatish Balay pc->ops->apply = PCApply_TFS_XXT; 1446218e575SSatish Balay } else { 1456218e575SSatish Balay tfs->xyt = XYT_new(); 146d0f46423SBarry Smith ierr = XYT_factor(tfs->xyt,localtoglobal,A->rmap->n,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" 17294c0dd61SBarry Smith /*MC 17394c0dd61SBarry Smith PCTFS - A parallel direct solver intended for problems with very few unknowns (like the 17494c0dd61SBarry Smith coarse grid in multigrid). 17594c0dd61SBarry Smith 17694c0dd61SBarry Smith Implemented by Henry M. Tufo III and Paul Fischer 17794c0dd61SBarry Smith 17894c0dd61SBarry Smith Level: beginner 17994c0dd61SBarry Smith 18094c0dd61SBarry Smith Notes: Only implemented for the MPIAIJ matrices 18194c0dd61SBarry Smith 18294c0dd61SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 18394c0dd61SBarry Smith M*/ 184dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_TFS(PC pc) 1856218e575SSatish Balay { 1866218e575SSatish Balay PetscErrorCode ierr; 1876218e575SSatish Balay PC_TFS *tfs; 1886218e575SSatish Balay 1896218e575SSatish Balay PetscFunctionBegin; 19038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_TFS,&tfs);CHKERRQ(ierr); 1916218e575SSatish Balay 1926218e575SSatish Balay tfs->xxt = 0; 1936218e575SSatish Balay tfs->xyt = 0; 1946218e575SSatish Balay tfs->b = 0; 1956218e575SSatish Balay tfs->xd = 0; 1966218e575SSatish Balay tfs->xo = 0; 1976218e575SSatish Balay tfs->nd = 0; 1986218e575SSatish Balay 1996218e575SSatish Balay pc->ops->apply = 0; 2006218e575SSatish Balay pc->ops->applytranspose = 0; 2016218e575SSatish Balay pc->ops->setup = PCSetUp_TFS; 2026218e575SSatish Balay pc->ops->destroy = PCDestroy_TFS; 2036218e575SSatish Balay pc->ops->setfromoptions = PCSetFromOptions_TFS; 2046218e575SSatish Balay pc->ops->view = PCView_TFS; 2056218e575SSatish Balay pc->ops->applyrichardson = 0; 2066218e575SSatish Balay pc->ops->applysymmetricleft = 0; 2076218e575SSatish Balay pc->ops->applysymmetricright = 0; 2086218e575SSatish Balay pc->data = (void*)tfs; 2096218e575SSatish Balay PetscFunctionReturn(0); 2106218e575SSatish Balay } 2116218e575SSatish Balay EXTERN_C_END 2126218e575SSatish Balay 213