xref: /petsc/src/ksp/pc/impls/tfs/tfs.c (revision 5f042095cd6794de959e41580753d6ff4dec299c)
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    *xx,*yy;
44   PetscErrorCode ierr;
45 
46   PetscFunctionBegin;
47   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
48   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
49   ierr = XXT_solve(tfs->xxt,yy,xx);CHKERRQ(ierr);
50   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
51   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
52   PetscFunctionReturn(0);
53 }
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "PCApply_TFS_XYT"
57 static PetscErrorCode PCApply_TFS_XYT(PC pc,Vec x,Vec y)
58 {
59   PC_TFS *tfs = (PC_TFS*)pc->data;
60   PetscScalar    *xx,*yy;
61   PetscErrorCode ierr;
62 
63   PetscFunctionBegin;
64   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
65   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
66   ierr = XYT_solve(tfs->xyt,yy,xx);CHKERRQ(ierr);
67   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
68   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
69   PetscFunctionReturn(0);
70 }
71 
72 #undef __FUNCT__
73 #define __FUNCT__ "LocalMult_TFS"
74 static PetscErrorCode LocalMult_TFS(PC pc,PetscScalar *xin,PetscScalar *xout)
75 {
76   PC_TFS        *tfs = (PC_TFS*)pc->data;
77   Mat           A = pc->pmat;
78   Mat_MPIAIJ    *a = (Mat_MPIAIJ*)A->data;
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   ierr = VecPlaceArray(tfs->b,xout);CHKERRQ(ierr);
83   ierr = VecPlaceArray(tfs->xd,xin);CHKERRQ(ierr);
84   ierr = VecPlaceArray(tfs->xo,xin+tfs->nd);CHKERRQ(ierr);
85   ierr = MatMult(a->A,tfs->xd,tfs->b);CHKERRQ(ierr);
86   ierr = MatMultAdd(a->B,tfs->xo,tfs->b,tfs->b);CHKERRQ(ierr);
87   ierr = VecResetArray(tfs->b);CHKERRQ(ierr);
88   ierr = VecResetArray(tfs->xd);CHKERRQ(ierr);
89   ierr = VecResetArray(tfs->xo);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "PCSetUp_TFS"
95 static PetscErrorCode PCSetUp_TFS(PC pc)
96 {
97   PC_TFS        *tfs = (PC_TFS*)pc->data;
98   Mat            A = pc->pmat;
99   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
100   PetscErrorCode ierr;
101   PetscInt      *localtoglobal,ncol,i;
102   PetscBool      ismpiaij;
103 
104   /*
105   PetscBool      issymmetric;
106   Petsc Real tol = 0.0;
107   */
108 
109   PetscFunctionBegin;
110   if (A->cmap->N != A->rmap->N) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_SIZ,"matrix must be square");
111   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
112   if (!ismpiaij) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Currently only supports MPIAIJ matrices");
113 
114   /* generate the local to global mapping */
115   ncol = a->A->cmap->n + a->B->cmap->n;
116   ierr = PetscMalloc((ncol)*sizeof(PetscInt),&localtoglobal);CHKERRQ(ierr);
117   for (i=0; i<a->A->cmap->n; i++) {
118     localtoglobal[i] = A->cmap->rstart + i + 1;
119   }
120   for (i=0; i<a->B->cmap->n; i++) {
121     localtoglobal[i+a->A->cmap->n] = a->garray[i] + 1;
122   }
123   /* generate the vectors needed for the local solves */
124   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->rmap->n,PETSC_NULL,&tfs->b);CHKERRQ(ierr);
125   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->cmap->n,PETSC_NULL,&tfs->xd);CHKERRQ(ierr);
126   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->cmap->n,PETSC_NULL,&tfs->xo);CHKERRQ(ierr);
127   tfs->nd = a->A->cmap->n;
128 
129 
130   /*  ierr =  MatIsSymmetric(A,tol,&issymmetric); */
131   /*  if (issymmetric) { */
132   ierr = PetscBarrier((PetscObject)pc);CHKERRQ(ierr);
133   if (A->symmetric) {
134     tfs->xxt       = XXT_new();
135     ierr           = XXT_factor(tfs->xxt,localtoglobal,A->rmap->n,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr);
136     pc->ops->apply = PCApply_TFS_XXT;
137   } else {
138     tfs->xyt       = XYT_new();
139     ierr           = XYT_factor(tfs->xyt,localtoglobal,A->rmap->n,ncol,(void*)LocalMult_TFS,pc);CHKERRQ(ierr);
140     pc->ops->apply = PCApply_TFS_XYT;
141   }
142 
143   ierr = PetscFree(localtoglobal);CHKERRQ(ierr);
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNCT__
148 #define __FUNCT__ "PCSetFromOptions_TFS"
149 static PetscErrorCode PCSetFromOptions_TFS(PC pc)
150 {
151   PetscFunctionBegin;
152   PetscFunctionReturn(0);
153 }
154 #undef __FUNCT__
155 #define __FUNCT__ "PCView_TFS"
156 static PetscErrorCode PCView_TFS(PC pc,PetscViewer viewer)
157 {
158   PetscFunctionBegin;
159   PetscFunctionReturn(0);
160 }
161 
162 EXTERN_C_BEGIN
163 #undef __FUNCT__
164 #define __FUNCT__ "PCCreate_TFS"
165 /*MC
166      PCTFS - A parallel direct solver intended for problems with very few unknowns (like the
167          coarse grid in multigrid).
168 
169    Implemented by  Henry M. Tufo III and Paul Fischer
170 
171    Level: beginner
172 
173    Notes: Only implemented for the MPIAIJ matrices
174 
175           Only works on a solver object that lives on all of PETSC_COMM_WORLD!
176 
177 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
178 M*/
179 PetscErrorCode  PCCreate_TFS(PC pc)
180 {
181   PetscErrorCode ierr;
182   PC_TFS         *tfs;
183   PetscMPIInt    cmp;
184 
185   PetscFunctionBegin;
186   ierr = MPI_Comm_compare(PETSC_COMM_WORLD,((PetscObject)pc)->comm,&cmp);CHKERRQ(ierr);
187   if (cmp != MPI_IDENT && cmp != MPI_CONGRUENT) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"TFS only works with PETSC_COMM_WORLD objects");
188   ierr = PetscNewLog(pc,PC_TFS,&tfs);CHKERRQ(ierr);
189 
190   tfs->xxt = 0;
191   tfs->xyt = 0;
192   tfs->b   = 0;
193   tfs->xd  = 0;
194   tfs->xo  = 0;
195   tfs->nd  = 0;
196 
197   pc->ops->apply               = 0;
198   pc->ops->applytranspose      = 0;
199   pc->ops->setup               = PCSetUp_TFS;
200   pc->ops->destroy             = PCDestroy_TFS;
201   pc->ops->setfromoptions      = PCSetFromOptions_TFS;
202   pc->ops->view                = PCView_TFS;
203   pc->ops->applyrichardson     = 0;
204   pc->ops->applysymmetricleft  = 0;
205   pc->ops->applysymmetricright = 0;
206   pc->data                     = (void*)tfs;
207   PetscFunctionReturn(0);
208 }
209 EXTERN_C_END
210 
211