1 /* 2 Provides an interface to the Tufo-Fischer parallel direct solver 3 */ 4 5 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 6 #include <../src/mat/impls/aij/mpi/mpiaij.h> 7 #include <../src/ksp/pc/impls/tfs/tfs.h> 8 9 typedef struct { 10 xxt_ADT xxt; 11 xyt_ADT xyt; 12 Vec b,xd,xo; 13 PetscInt nd; 14 } PC_TFS; 15 16 #undef __FUNCT__ 17 #define __FUNCT__ "PCDestroy_TFS" 18 PetscErrorCode PCDestroy_TFS(PC pc) 19 { 20 PC_TFS *tfs = (PC_TFS*)pc->data; 21 PetscErrorCode ierr; 22 23 PetscFunctionBegin; 24 /* free the XXT datastructures */ 25 if (tfs->xxt) { 26 ierr = XXT_free(tfs->xxt);CHKERRQ(ierr); 27 } 28 if (tfs->xyt) { 29 ierr = XYT_free(tfs->xyt);CHKERRQ(ierr); 30 } 31 ierr = VecDestroy(&tfs->b);CHKERRQ(ierr); 32 ierr = VecDestroy(&tfs->xd);CHKERRQ(ierr); 33 ierr = VecDestroy(&tfs->xo);CHKERRQ(ierr); 34 ierr = PetscFree(pc->data);CHKERRQ(ierr); 35 PetscFunctionReturn(0); 36 } 37 38 #undef __FUNCT__ 39 #define __FUNCT__ "PCApply_TFS_XXT" 40 static PetscErrorCode PCApply_TFS_XXT(PC pc,Vec x,Vec y) 41 { 42 PC_TFS *tfs = (PC_TFS*)pc->data; 43 PetscScalar *yy; 44 const PetscScalar *xx; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 49 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 50 ierr = XXT_solve(tfs->xxt,yy,(PetscScalar*)xx);CHKERRQ(ierr); 51 ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 52 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 53 PetscFunctionReturn(0); 54 } 55 56 #undef __FUNCT__ 57 #define __FUNCT__ "PCApply_TFS_XYT" 58 static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y) 59 { 60 PC_TFS *tfs = (PC_TFS*)pc->data; 61 PetscScalar *yy; 62 const PetscScalar *xx; 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 67 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 68 ierr = XYT_solve(tfs->xyt,yy,(PetscScalar*)xx);CHKERRQ(ierr); 69 ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 70 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "PCTFSLocalMult_TFS" 76 static PetscErrorCode PCTFSLocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout) 77 { 78 PC_TFS *tfs = (PC_TFS*)pc->data; 79 Mat A = pc->pmat; 80 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 81 PetscErrorCode ierr; 82 83 PetscFunctionBegin; 84 ierr = VecPlaceArray(tfs->b,xout);CHKERRQ(ierr); 85 ierr = VecPlaceArray(tfs->xd,xin);CHKERRQ(ierr); 86 ierr = VecPlaceArray(tfs->xo,xin+tfs->nd);CHKERRQ(ierr); 87 ierr = MatMult(a->A,tfs->xd,tfs->b);CHKERRQ(ierr); 88 ierr = MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b);CHKERRQ(ierr); 89 ierr = VecResetArray(tfs->b);CHKERRQ(ierr); 90 ierr = VecResetArray(tfs->xd);CHKERRQ(ierr); 91 ierr = VecResetArray(tfs->xo);CHKERRQ(ierr); 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "PCSetUp_TFS" 97 static PetscErrorCode PCSetUp_TFS(PC pc) 98 { 99 PC_TFS *tfs = (PC_TFS*)pc->data; 100 Mat A = pc->pmat; 101 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 102 PetscErrorCode ierr; 103 PetscInt *localtoglobal,ncol,i; 104 PetscBool ismpiaij; 105 106 /* 107 PetscBool issymmetric; 108 Petsc Real tol = 0.0; 109 */ 110 111 PetscFunctionBegin; 112 if (A->cmap->N != A->rmap->N) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_SIZ,"matrix must be square"); 113 ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 114 if (!ismpiaij) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices"); 115 116 /* generate the local to global mapping */ 117 ncol = a->A->cmap->n + a->B->cmap->n; 118 ierr = PetscMalloc1(ncol,&localtoglobal);CHKERRQ(ierr); 119 for (i=0; i<a->A->cmap->n; i++) localtoglobal[i] = A->cmap->rstart + i + 1; 120 for (i=0; i<a->B->cmap->n; i++) localtoglobal[i+a->A->cmap->n] = a->garray[i] + 1; 121 122 /* generate the vectors needed for the local solves */ 123 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->A->rmap->n,NULL,&tfs->b);CHKERRQ(ierr); 124 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->A->cmap->n,NULL,&tfs->xd);CHKERRQ(ierr); 125 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,a->B->cmap->n,NULL,&tfs->xo);CHKERRQ(ierr); 126 tfs->nd = a->A->cmap->n; 127 128 129 /* ierr = MatIsSymmetric(A,tol,&issymmetric); */ 130 /* if (issymmetric) { */ 131 ierr = PetscBarrier((PetscObject)pc);CHKERRQ(ierr); 132 if (A->symmetric) { 133 tfs->xxt = XXT_new(); 134 ierr = XXT_factor(tfs->xxt,localtoglobal,A->rmap->n,ncol,(PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*))PCTFSLocalMult_TFS,pc);CHKERRQ(ierr); 135 pc->ops->apply = PCApply_TFS_XXT; 136 } else { 137 tfs->xyt = XYT_new(); 138 ierr = XYT_factor(tfs->xyt,localtoglobal,A->rmap->n,ncol,(PetscErrorCode (*)(void*,PetscScalar*,PetscScalar*))PCTFSLocalMult_TFS,pc);CHKERRQ(ierr); 139 pc->ops->apply = PCApply_TFS_XYT; 140 } 141 142 ierr = PetscFree(localtoglobal);CHKERRQ(ierr); 143 PetscFunctionReturn(0); 144 } 145 146 #undef __FUNCT__ 147 #define __FUNCT__ "PCSetFromOptions_TFS" 148 static PetscErrorCode PCSetFromOptions_TFS(PetscOptionItems *PetscOptionsObject,PC pc) 149 { 150 PetscFunctionBegin; 151 PetscFunctionReturn(0); 152 } 153 #undef __FUNCT__ 154 #define __FUNCT__ "PCView_TFS" 155 static PetscErrorCode PCView_TFS(PC pc,PetscViewer viewer) 156 { 157 PetscFunctionBegin; 158 PetscFunctionReturn(0); 159 } 160 161 #undef __FUNCT__ 162 #define __FUNCT__ "PCCreate_TFS" 163 /*MC 164 PCTFS - A parallel direct solver intended for problems with very few unknowns (like the 165 coarse grid in multigrid). 166 167 Implemented by Henry M. Tufo III and Paul Fischer 168 169 Level: beginner 170 171 Notes: Only implemented for the MPIAIJ matrices 172 173 Only works on a solver object that lives on all of PETSC_COMM_WORLD! 174 175 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 176 M*/ 177 PETSC_EXTERN PetscErrorCode PCCreate_TFS(PC pc) 178 { 179 PetscErrorCode ierr; 180 PC_TFS *tfs; 181 PetscMPIInt cmp; 182 183 PetscFunctionBegin; 184 ierr = MPI_Comm_compare(PETSC_COMM_WORLD,PetscObjectComm((PetscObject)pc),&cmp);CHKERRQ(ierr); 185 if (cmp != MPI_IDENT && cmp != MPI_CONGRUENT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"TFS only works with PETSC_COMM_WORLD objects"); 186 ierr = PetscNewLog(pc,&tfs);CHKERRQ(ierr); 187 188 tfs->xxt = 0; 189 tfs->xyt = 0; 190 tfs->b = 0; 191 tfs->xd = 0; 192 tfs->xo = 0; 193 tfs->nd = 0; 194 195 pc->ops->apply = 0; 196 pc->ops->applytranspose = 0; 197 pc->ops->setup = PCSetUp_TFS; 198 pc->ops->destroy = PCDestroy_TFS; 199 pc->ops->setfromoptions = PCSetFromOptions_TFS; 200 pc->ops->view = PCView_TFS; 201 pc->ops->applyrichardson = 0; 202 pc->ops->applysymmetricleft = 0; 203 pc->ops->applysymmetricright = 0; 204 pc->data = (void*)tfs; 205 PetscFunctionReturn(0); 206 } 207 208