Actual source code: pxxt.c

  1: /*$Id: pxxt.c,v 1.6 2001/04/10 19:35:25 bsmith 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)
 10:  #include src/contrib/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);

 21: int MatDestroy_MPIAIJ_XXT(Mat A)
 22: {
 23:   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data;
 24:   Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)a->spptr;
 25:   int            ierr;

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

 35: int MatSolve_MPIAIJ_XXT(Mat A,Vec b,Vec x)
 36: {
 37:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
 38:   Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)a->spptr;
 39:   Scalar         *bb,*xx;
 40:   int            ierr;

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

 51: int MatLUFactorNumeric_MPIAIJ_XXT(Mat A,Mat *F)
 52: {
 54:   return(0);
 55: }

 57: static int LocalMult_XXT(Mat_MPIAIJ *a,Scalar *xin,Scalar *xout)
 58: {
 59:   int            ierr;
 60:   Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)a->spptr;

 63:   VecPlaceArray(xxt->b,xout);
 64:   VecPlaceArray(xxt->xd,xin);
 65:   VecPlaceArray(xxt->xo,xin+xxt->nd);
 66:   MatMult(a->A,xxt->xd,xxt->b);
 67:   MatMultAdd(a->B,xxt->xo,xxt->b,xxt->b);
 68:   /*
 69:   PetscDoubleView(a->A->n+a->B->n,xin,PETSC_VIEWER_STDOUT_WORLD);
 70:   PetscDoubleView(a->A->m,xout,PETSC_VIEWER_STDOUT_WORLD);
 71:   */
 72:   return(0);
 73: }

 75: int MatLUFactorSymbolic_MPIAIJ_XXT(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
 76: {
 77:   Mat            B;
 78:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data,*b;
 79:   int            ierr,*localtoglobal,ncol,i;
 80:   Mat_MPIAIJ_XXT *xxt;

 83:   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
 84:   MatCreateMPIAIJ(A->comm,A->m,A->n,A->M,A->N,0,PETSC_NULL,0,PETSC_NULL,F);
 85:   B                       = *F;
 86:   B->ops->solve           = MatSolve_MPIAIJ_XXT;
 87:   B->ops->destroy         = MatDestroy_MPIAIJ_XXT;
 88:   B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XXT;
 89:   B->factor               = FACTOR_LU;
 90:   b                       = (Mat_MPIAIJ*)B->data;
 91:   ierr                    = PetscNew(Mat_MPIAIJ_XXT,&xxt);
 92:   b->spptr = a->spptr     = (void*)xxt;

 94:   xxt->xxt = XXT_new();
 95:   /* generate the local to global mapping */
 96:   ncol = a->A->n + a->B->n;
 97:   PetscMalloc((ncol)*sizeof(int),&localtoglobal);
 98:   for (i=0; i<a->A->n; i++) {
 99:     localtoglobal[i] = a->cstart + i + 1;
100:   }
101:   for (i=0; i<a->B->n; i++) {
102:     localtoglobal[i+a->A->n] = a->garray[i] + 1;
103:   }
104:   /*
105:   PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
106:   */

108:   /* generate the vectors needed for the local solves */
109:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xxt->b);
110:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xxt->xd);
111:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xxt->xo);
112:   xxt->nd = a->A->n;

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

117:   VecDestroy(xxt->b);
118:   VecDestroy(xxt->xd);
119:   VecDestroy(xxt->xo);
120:   PetscFree(localtoglobal);

122:   return(0);
123: }

125:  #include src/contrib/libtfs/xyt.h

127: typedef struct {
128:   xyt_ADT xyt;
129:   Vec     b,xd,xo;
130:   int     nd;
131: } Mat_MPIAIJ_XYT;


134: EXTERN int MatDestroy_MPIAIJ(Mat);

136: int MatDestroy_MPIAIJ_XYT(Mat A)
137: {
138:   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data;
139:   Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)a->spptr;
140:   int            ierr;

143:   /* free the XYT datastructures */
144:   XYT_free(xyt->xyt);
145:   PetscFree(xyt);
146:   MatDestroy_MPIAIJ(A);
147:   return(0);
148: }

150: int MatSolve_MPIAIJ_XYT(Mat A,Vec b,Vec x)
151: {
152:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data;
153:   Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)a->spptr;
154:   Scalar         *bb,*xx;
155:   int            ierr;

158:   VecGetArray(b,&bb);
159:   VecGetArray(x,&xx);
160:   XYT_solve(xyt->xyt,xx,bb);
161:   VecRestoreArray(b,&bb);
162:   VecRestoreArray(x,&xx);
163:   return(0);
164: }

166: int MatLUFactorNumeric_MPIAIJ_XYT(Mat A,Mat *F)
167: {
169:   return(0);
170: }

172: static int LocalMult_XYT(Mat_MPIAIJ *a,Scalar *xin,Scalar *xout)
173: {
174:   int            ierr;
175:   Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)a->spptr;

178:   VecPlaceArray(xyt->b,xout);
179:   VecPlaceArray(xyt->xd,xin);
180:   VecPlaceArray(xyt->xo,xin+xyt->nd);
181:   MatMult(a->A,xyt->xd,xyt->b);
182:   MatMultAdd(a->B,xyt->xo,xyt->b,xyt->b);
183:   /*
184:   PetscDoubleView(a->A->n+a->B->n,xin,PETSC_VIEWER_STDOUT_WORLD);
185:   PetscDoubleView(a->A->m,xout,PETSC_VIEWER_STDOUT_WORLD);
186:   */
187:   return(0);
188: }

190: int MatLUFactorSymbolic_MPIAIJ_XYT(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
191: {
192:   Mat            B;
193:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*)A->data,*b;
194:   int            ierr,*localtoglobal,ncol,i;
195:   Mat_MPIAIJ_XYT *xyt;

198:   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
199:   MatCreateMPIAIJ(A->comm,A->m,A->n,A->M,A->N,0,PETSC_NULL,0,PETSC_NULL,F);
200:   B                       = *F;
201:   B->ops->solve           = MatSolve_MPIAIJ_XYT;
202:   B->ops->destroy         = MatDestroy_MPIAIJ_XYT;
203:   B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XYT;
204:   B->factor               = FACTOR_LU;
205:   b                       = (Mat_MPIAIJ*)B->data;
206:   ierr                    = PetscNew(Mat_MPIAIJ_XYT,&xyt);
207:   b->spptr = a->spptr     = (void*)xyt;

209:   xyt->xyt = XYT_new();
210:   /* generate the local to global mapping */
211:   ncol = a->A->n + a->B->n;
212:   PetscMalloc((ncol)*sizeof(int),&localtoglobal);
213:   for (i=0; i<a->A->n; i++) {
214:     localtoglobal[i] = a->cstart + i + 1;
215:   }
216:   for (i=0; i<a->B->n; i++) {
217:     localtoglobal[i+a->A->n] = a->garray[i] + 1;
218:   }
219:   /*
220:   PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
221:   */

223:   /* generate the vectors needed for the local solves */
224:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xyt->b);
225:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xyt->xd);
226:   VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xyt->xo);
227:   xyt->nd = a->A->n;

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

232:   VecDestroy(xyt->b);
233:   VecDestroy(xyt->xd);
234:   VecDestroy(xyt->xo);
235:   PetscFree(localtoglobal);

237:   return(0);
238: }

240: int MatLUFactorSymbolic_MPIAIJ_TFS(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
241: {

245:   PetscLogInfo(0,"Using TFS for MPIAIJ LU factorization and solves");
246:   if (A->symmetric) {
247:     MatLUFactorSymbolic_MPIAIJ_XXT(A,r,c,info,F);
248:   } else {
249:     MatLUFactorSymbolic_MPIAIJ_XYT(A,r,c,info,F);
250:   }
251:   return(0);
252: }

254: #else

256: int dummy83(int a)
257: {
259:   return(0);
260: }

262: #endif