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