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" 77758a8cdSBarry Smith #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 } 316218e575SSatish Balay if (tfs->b) { 326218e575SSatish Balay ierr = VecDestroy(tfs->b);CHKERRQ(ierr); 336218e575SSatish Balay } 346218e575SSatish Balay if (tfs->xd) { 356218e575SSatish Balay ierr = VecDestroy(tfs->xd);CHKERRQ(ierr); 366218e575SSatish Balay } 376218e575SSatish Balay if (tfs->xo) { 386218e575SSatish Balay ierr = VecDestroy(tfs->xo);CHKERRQ(ierr); 396218e575SSatish Balay } 406218e575SSatish Balay ierr = PetscFree(tfs);CHKERRQ(ierr); 416218e575SSatish Balay PetscFunctionReturn(0); 426218e575SSatish Balay } 436218e575SSatish Balay 446218e575SSatish Balay #undef __FUNCT__ 456218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XXT" 466218e575SSatish Balay static PetscErrorCode PCApply_TFS_XXT(PC pc,Vec x,Vec y) 476218e575SSatish Balay { 486218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 496218e575SSatish Balay PetscScalar *xx,*yy; 506218e575SSatish Balay PetscErrorCode ierr; 516218e575SSatish Balay 526218e575SSatish Balay PetscFunctionBegin; 536218e575SSatish Balay ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 546218e575SSatish Balay ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 556218e575SSatish Balay ierr = XXT_solve(tfs->xxt,yy,xx);CHKERRQ(ierr); 566218e575SSatish Balay ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 576218e575SSatish Balay ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 586218e575SSatish Balay PetscFunctionReturn(0); 596218e575SSatish Balay } 606218e575SSatish Balay 616218e575SSatish Balay #undef __FUNCT__ 626218e575SSatish Balay #define __FUNCT__ "PCApply_TFS_XYT" 636218e575SSatish Balay static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y) 646218e575SSatish Balay { 656218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 666218e575SSatish Balay PetscScalar *xx,*yy; 676218e575SSatish Balay PetscErrorCode ierr; 686218e575SSatish Balay 696218e575SSatish Balay PetscFunctionBegin; 706218e575SSatish Balay ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 716218e575SSatish Balay ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 726218e575SSatish Balay ierr = XYT_solve(tfs->xyt,yy,xx);CHKERRQ(ierr); 736218e575SSatish Balay ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 746218e575SSatish Balay ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 756218e575SSatish Balay PetscFunctionReturn(0); 766218e575SSatish Balay } 776218e575SSatish Balay 786218e575SSatish Balay #undef __FUNCT__ 796218e575SSatish Balay #define __FUNCT__ "LocalMult_TFS" 806218e575SSatish Balay static PetscErrorCode LocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout) 816218e575SSatish Balay { 826218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 836218e575SSatish Balay Mat A = pc->pmat; 846218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 856218e575SSatish Balay PetscErrorCode ierr; 866218e575SSatish Balay 876218e575SSatish Balay PetscFunctionBegin; 886218e575SSatish Balay ierr = VecPlaceArray(tfs->b,xout);CHKERRQ(ierr); 896218e575SSatish Balay ierr = VecPlaceArray(tfs->xd,xin);CHKERRQ(ierr); 906218e575SSatish Balay ierr = VecPlaceArray(tfs->xo,xin+tfs->nd);CHKERRQ(ierr); 916218e575SSatish Balay ierr = MatMult(a->A,tfs->xd,tfs->b);CHKERRQ(ierr); 926218e575SSatish Balay ierr = MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b);CHKERRQ(ierr); 936218e575SSatish Balay PetscFunctionReturn(0); 946218e575SSatish Balay } 956218e575SSatish Balay 966218e575SSatish Balay #undef __FUNCT__ 976218e575SSatish Balay #define __FUNCT__ "PCSetUp_TFS" 986218e575SSatish Balay static PetscErrorCode PCSetUp_TFS(PC pc) 996218e575SSatish Balay { 1006218e575SSatish Balay PC_TFS *tfs = (PC_TFS*)pc->data; 1016218e575SSatish Balay Mat A = pc->pmat; 1026218e575SSatish Balay Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 1036218e575SSatish Balay PetscErrorCode ierr; 1046218e575SSatish Balay PetscInt *localtoglobal,ncol,i; 1056218e575SSatish Balay PetscTruth ismpiaij; 1066218e575SSatish Balay 1076218e575SSatish Balay /* 1086218e575SSatish Balay PetscTruth issymmetric; 1096218e575SSatish Balay Petsc Real tol = 0.0; 1106218e575SSatish Balay */ 1116218e575SSatish Balay 1126218e575SSatish Balay PetscFunctionBegin; 1136218e575SSatish Balay if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square"); 1146218e575SSatish Balay ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 1156218e575SSatish Balay if (!ismpiaij) { 1166218e575SSatish Balay SETERRQ(PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices"); 1176218e575SSatish Balay } 1186218e575SSatish Balay 1196218e575SSatish Balay /* generate the local to global mapping */ 1206218e575SSatish Balay ncol = a->A->n + a->B->n; 1216218e575SSatish Balay ierr = PetscMalloc((ncol)*sizeof(int),&localtoglobal);CHKERRQ(ierr); 1226218e575SSatish Balay for (i=0; i<a->A->n; i++) { 1236218e575SSatish Balay localtoglobal[i] = a->cstart + i + 1; 1246218e575SSatish Balay } 1256218e575SSatish Balay for (i=0; i<a->B->n; i++) { 1266218e575SSatish Balay localtoglobal[i+a->A->n] = a->garray[i] + 1; 1276218e575SSatish Balay } 1286218e575SSatish Balay /* generate the vectors needed for the local solves */ 1296218e575SSatish Balay ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&tfs->b);CHKERRQ(ierr); 1306218e575SSatish Balay ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&tfs->xd);CHKERRQ(ierr); 1316218e575SSatish Balay ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&tfs->xo);CHKERRQ(ierr); 1326218e575SSatish Balay tfs->nd = a->A->n; 1336218e575SSatish Balay 1346218e575SSatish Balay 1356218e575SSatish Balay /* ierr = MatIsSymmetric(A,tol,&issymmetric); */ 1366218e575SSatish Balay /* if (issymmetric) { */ 137585bbf96SBarry Smith ierr = PetscBarrier((PetscObject)pc);CHKERRQ(ierr); 1386218e575SSatish Balay if (A->symmetric) { 1396218e575SSatish Balay tfs->xxt = XXT_new(); 1406218e575SSatish Balay ierr = XXT_factor(tfs->xxt,localtoglobal,A->m,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr); 1416218e575SSatish Balay pc->ops->apply = PCApply_TFS_XXT; 1426218e575SSatish Balay } else { 1436218e575SSatish Balay tfs->xyt = XYT_new(); 1446218e575SSatish Balay ierr = XYT_factor(tfs->xyt,localtoglobal,A->m,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr); 1456218e575SSatish Balay pc->ops->apply = PCApply_TFS_XYT; 1466218e575SSatish Balay } 1476218e575SSatish Balay 1486218e575SSatish Balay ierr = PetscFree(localtoglobal);CHKERRQ(ierr); 1496218e575SSatish Balay PetscFunctionReturn(0); 1506218e575SSatish Balay } 1516218e575SSatish Balay 1526218e575SSatish Balay #undef __FUNCT__ 1536218e575SSatish Balay #define __FUNCT__ "PCSetFromOptions_TFS" 1546218e575SSatish Balay static PetscErrorCode PCSetFromOptions_TFS(PC pc) 1556218e575SSatish Balay { 1566218e575SSatish Balay PetscFunctionBegin; 1576218e575SSatish Balay PetscFunctionReturn(0); 1586218e575SSatish Balay } 1596218e575SSatish Balay #undef __FUNCT__ 1606218e575SSatish Balay #define __FUNCT__ "PCView_TFS" 1616218e575SSatish Balay static PetscErrorCode PCView_TFS(PC pc,PetscViewer viewer) 1626218e575SSatish Balay { 1636218e575SSatish Balay PetscFunctionBegin; 1646218e575SSatish Balay PetscFunctionReturn(0); 1656218e575SSatish Balay } 1666218e575SSatish Balay 1676218e575SSatish Balay EXTERN_C_BEGIN 1686218e575SSatish Balay #undef __FUNCT__ 1696218e575SSatish Balay #define __FUNCT__ "PCCreate_TFS" 1706218e575SSatish Balay PetscErrorCode PCCreate_TFS(PC pc) 1716218e575SSatish Balay { 1726218e575SSatish Balay PetscErrorCode ierr; 1736218e575SSatish Balay PC_TFS *tfs; 1746218e575SSatish Balay 1756218e575SSatish Balay PetscFunctionBegin; 1766218e575SSatish Balay ierr = PetscNew(PC_TFS,&tfs);CHKERRQ(ierr); 177*52e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_TFS));CHKERRQ(ierr); 1786218e575SSatish Balay 1796218e575SSatish Balay tfs->xxt = 0; 1806218e575SSatish Balay tfs->xyt = 0; 1816218e575SSatish Balay tfs->b = 0; 1826218e575SSatish Balay tfs->xd = 0; 1836218e575SSatish Balay tfs->xo = 0; 1846218e575SSatish Balay tfs->nd = 0; 1856218e575SSatish Balay 1866218e575SSatish Balay pc->ops->apply = 0; 1876218e575SSatish Balay pc->ops->applytranspose = 0; 1886218e575SSatish Balay pc->ops->setup = PCSetUp_TFS; 1896218e575SSatish Balay pc->ops->destroy = PCDestroy_TFS; 1906218e575SSatish Balay pc->ops->setfromoptions = PCSetFromOptions_TFS; 1916218e575SSatish Balay pc->ops->view = PCView_TFS; 1926218e575SSatish Balay pc->ops->applyrichardson = 0; 1936218e575SSatish Balay pc->ops->applysymmetricleft = 0; 1946218e575SSatish Balay pc->ops->applysymmetricright = 0; 1956218e575SSatish Balay pc->data = (void*)tfs; 1966218e575SSatish Balay PetscFunctionReturn(0); 1976218e575SSatish Balay } 1986218e575SSatish Balay EXTERN_C_END 1996218e575SSatish Balay 200