Actual source code: pxxt.c

  1: /*$Id: pxxt.c,v 1.8 2001/08/07 03:02:49 balay Exp $*/

  3: /* 
  4:         Provides an interface to the Tufo-Fischer parallel direct solver

  6: */
 7:  #include src/mat/impls/aij/mpi/mpiaij.h

  9: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
 10:  #include src/mat/impls/aij/mpi/libtfs/xxt.h

 12: typedef struct {
 13:   xxt_ADT xxt;
 14:   Vec     b,xd,xo;
 15:   int     nd;
 16: } Mat_MPIAIJ_XXT;


 19: EXTERN int MatDestroy_MPIAIJ(Mat);

 23: int MatDestroy_MPIAIJ_XXT(Mat A)
 24: {
 25:   Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)A->spptr;
 26:   int            ierr;

 29:   /* free the XXT datastructures */
 30:   XXT_free(xxt->xxt);
 31:   PetscFree(xxt);
 32:   MatDestroy_MPIAIJ(A);
 33:   return(0);
 34: }

 38: int MatSolve_MPIAIJ_XXT(Mat A,Vec b,Vec x)
 39: {
 40:   Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)A->spptr;
 41:   PetscScalar    *bb,*xx;
 42:   int            ierr;

 45:   VecGetArray(b,&bb);
 46:   VecGetArray(x,&xx);
 47:   XXT_solve(xxt->xxt,xx,bb);
 48:   VecRestoreArray(b,&bb);
 49:   VecRestoreArray(x,&xx);
 50:   return(0);
 51: }

 55: int MatLUFactorNumeric_MPIAIJ_XXT(Mat A,Mat *F)
 56: {
 58:   return(0);
 59: }

 63: static int LocalMult_XXT(Mat A,PetscScalar *xin,PetscScalar *xout)
 64: {
 65:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
 66:   int            ierr;
 67:   Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)A->spptr;

 70:   VecPlaceArray(xxt->b,xout);
 71:   VecPlaceArray(xxt->xd,xin);
 72:   VecPlaceArray(xxt->xo,xin+xxt->nd);
 73:   MatMult(a->A,xxt->xd,xxt->b);
 74:   MatMultAdd(a->B,xxt->xo,xxt->b,xxt->b);
 75:   /*
 76:   PetscRealView(a->A->n+a->B->n,xin,PETSC_VIEWER_STDOUT_WORLD);
 77:   PetscRealView(a->A->m,xout,PETSC_VIEWER_STDOUT_WORLD);
 78:   */
 79:   return(0);
 80: }

 84: int MatLUFactorSymbolic_MPIAIJ_XXT(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
 85: {
 86:   Mat            B;
 87:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
 88:   int            ierr,*localtoglobal,ncol,i;
 89:   Mat_MPIAIJ_XXT *xxt;

 92:   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
 93:   MatCreate(A->comm,A->m,A->n,A->M,A->N,F);
 94:   MatSetType(*F,A->type_name);
 95:   MatMPIAIJSetPreallocation(*F,0,PETSC_NULL,0,PETSC_NULL);
 96:   B                       = *F;
 97:   B->ops->solve           = MatSolve_MPIAIJ_XXT;
 98:   B->ops->destroy         = MatDestroy_MPIAIJ_XXT;
 99:   B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XXT;
100:   B->factor               = FACTOR_LU;
101:   B->assembled            = PETSC_TRUE;
102:   PetscNew(Mat_MPIAIJ_XXT,&xxt);
103:   B->spptr = A->spptr     = (void*)xxt;

105:   xxt->xxt = XXT_new();
106:   /* generate the local to global mapping */
107:   ncol = a->A->n + a->B->n;
108:   PetscMalloc((ncol)*sizeof(int),&localtoglobal);
109:   for (i=0; i<a->A->n; i++) {
110:     localtoglobal[i] = a->cstart + i + 1;
111:   }
112:   for (i=0; i<a->B->n; i++) {
113:     localtoglobal[i+a->A->n] = a->garray[i] + 1;
114:   }
115:   /*
116:   PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
117:   */

119:   /* generate the vectors needed for the local solves */
120:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xxt->b);
121:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xxt->xd);
122:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xxt->xo);
123:   xxt->nd = a->A->n;

125:   /* factor the beast */
126:   XXT_factor(xxt->xxt,localtoglobal,A->m,ncol,(void*)LocalMult_XXT,a);

128:   VecDestroy(xxt->b);
129:   VecDestroy(xxt->xd);
130:   VecDestroy(xxt->xo);
131:   PetscFree(localtoglobal);

133:   return(0);
134: }

136:  #include src/mat/impls/aij/mpi/libtfs/xyt.h

138: typedef struct {
139:   xyt_ADT xyt;
140:   Vec     b,xd,xo;
141:   int     nd;
142: } Mat_MPIAIJ_XYT;


145: EXTERN int MatDestroy_MPIAIJ(Mat);

149: int MatDestroy_MPIAIJ_XYT(Mat A)
150: {
151:   Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
152:   int            ierr;

155:   /* free the XYT datastructures */
156:   XYT_free(xyt->xyt);
157:   PetscFree(xyt);
158:   MatDestroy_MPIAIJ(A);
159:   return(0);
160: }

164: int MatSolve_MPIAIJ_XYT(Mat A,Vec b,Vec x)
165: {
166:   Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
167:   PetscScalar    *bb,*xx;
168:   int            ierr;

171:   VecGetArray(b,&bb);
172:   VecGetArray(x,&xx);
173:   XYT_solve(xyt->xyt,xx,bb);
174:   VecRestoreArray(b,&bb);
175:   VecRestoreArray(x,&xx);
176:   return(0);
177: }

181: int MatLUFactorNumeric_MPIAIJ_XYT(Mat A,Mat *F)
182: {
184:   return(0);
185: }

189: static int LocalMult_XYT(Mat A,PetscScalar *xin,PetscScalar *xout)
190: {
191:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
192:   int            ierr;
193:   Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;

196:   VecPlaceArray(xyt->b,xout);
197:   VecPlaceArray(xyt->xd,xin);
198:   VecPlaceArray(xyt->xo,xin+xyt->nd);
199:   MatMult(a->A,xyt->xd,xyt->b);
200:   MatMultAdd(a->B,xyt->xo,xyt->b,xyt->b);
201:   /*
202:   PetscRealView(a->A->n+a->B->n,xin,PETSC_VIEWER_STDOUT_WORLD);
203:   PetscRealView(a->A->m,xout,PETSC_VIEWER_STDOUT_WORLD);
204:   */
205:   return(0);
206: }

210: int MatLUFactorSymbolic_MPIAIJ_XYT(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
211: {
212:   Mat            B;
213:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
214:   int            ierr,*localtoglobal,ncol,i;
215:   Mat_MPIAIJ_XYT *xyt;

218:   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
219:   MatCreate(A->comm,A->m,A->n,A->M,A->N,F);
220:   MatSetType(*F,A->type_name);
221:   MatMPIAIJSetPreallocation(*F,0,PETSC_NULL,0,PETSC_NULL);
222:   B                       = *F;
223:   B->ops->solve           = MatSolve_MPIAIJ_XYT;
224:   B->ops->destroy         = MatDestroy_MPIAIJ_XYT;
225:   B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XYT;
226:   B->factor               = FACTOR_LU;
227:   B->assembled            = PETSC_TRUE;
228:   PetscNew(Mat_MPIAIJ_XYT,&xyt);
229:   B->spptr = A->spptr     = (void*)xyt;

231:   xyt->xyt = XYT_new();
232:   /* generate the local to global mapping */
233:   ncol = a->A->n + a->B->n;
234:   PetscMalloc((ncol)*sizeof(int),&localtoglobal);
235:   for (i=0; i<a->A->n; i++) {
236:     localtoglobal[i] = a->cstart + i + 1;
237:   }
238:   for (i=0; i<a->B->n; i++) {
239:     localtoglobal[i+a->A->n] = a->garray[i] + 1;
240:   }
241:   /*
242:   PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
243:   */

245:   /* generate the vectors needed for the local solves */
246:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xyt->b);
247:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xyt->xd);
248:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xyt->xo);
249:   xyt->nd = a->A->n;

251:   /* factor the beast */
252:   XYT_factor(xyt->xyt,localtoglobal,A->m,ncol,(void*)LocalMult_XYT,A);

254:   VecDestroy(xyt->b);
255:   VecDestroy(xyt->xd);
256:   VecDestroy(xyt->xo);
257:   PetscFree(localtoglobal);

259:   return(0);
260: }

264: int MatLUFactorSymbolic_MPIAIJ_TFS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
265: {

269:   PetscLogInfo(0,"Using TFS for MPIAIJ LU factorization and solves");
270:   if (A->symmetric) {
271:     MatLUFactorSymbolic_MPIAIJ_XXT(A,r,c,info,F);
272:   } else {
273:     MatLUFactorSymbolic_MPIAIJ_XYT(A,r,c,info,F);
274:   }
275:   return(0);
276: }

278: #else

282: int dummy83(int a)
283: {
285:   return(0);
286: }

288: #endif