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/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_XXT *xxt = (Mat_MPIAIJ_XXT*)A->spptr;
24: int ierr;
27: /* free the XXT datastructures */
28: XXT_free(xxt->xxt);
29: PetscFree(xxt);
30: MatDestroy_MPIAIJ(A);
31: return(0);
32: }
34: int MatSolve_MPIAIJ_XXT(Mat A,Vec b,Vec x)
35: {
36: Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)A->spptr;
37: PetscScalar *bb,*xx;
38: int ierr;
41: VecGetArray(b,&bb);
42: VecGetArray(x,&xx);
43: XXT_solve(xxt->xxt,xx,bb);
44: VecRestoreArray(b,&bb);
45: VecRestoreArray(x,&xx);
46: return(0);
47: }
49: int MatLUFactorNumeric_MPIAIJ_XXT(Mat A,Mat *F)
50: {
52: return(0);
53: }
55: static int LocalMult_XXT(Mat A,PetscScalar *xin,PetscScalar *xout)
56: {
57: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
58: int ierr;
59: Mat_MPIAIJ_XXT *xxt = (Mat_MPIAIJ_XXT*)A->spptr;
62: VecPlaceArray(xxt->b,xout);
63: VecPlaceArray(xxt->xd,xin);
64: VecPlaceArray(xxt->xo,xin+xxt->nd);
65: MatMult(a->A,xxt->xd,xxt->b);
66: MatMultAdd(a->B,xxt->xo,xxt->b,xxt->b);
67: /*
68: PetscRealView(a->A->n+a->B->n,xin,PETSC_VIEWER_STDOUT_WORLD);
69: PetscRealView(a->A->m,xout,PETSC_VIEWER_STDOUT_WORLD);
70: */
71: return(0);
72: }
74: int MatLUFactorSymbolic_MPIAIJ_XXT(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
75: {
76: Mat B;
77: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
78: int ierr,*localtoglobal,ncol,i;
79: Mat_MPIAIJ_XXT *xxt;
82: if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
83: MatCreateMPIAIJ(A->comm,A->m,A->n,A->M,A->N,0,PETSC_NULL,0,PETSC_NULL,F);
84: B = *F;
85: B->ops->solve = MatSolve_MPIAIJ_XXT;
86: B->ops->destroy = MatDestroy_MPIAIJ_XXT;
87: B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XXT;
88: B->factor = FACTOR_LU;
89: ierr = PetscNew(Mat_MPIAIJ_XXT,&xxt);
90: B->spptr = A->spptr = (void*)xxt;
92: xxt->xxt = XXT_new();
93: /* generate the local to global mapping */
94: ncol = a->A->n + a->B->n;
95: PetscMalloc((ncol)*sizeof(int),&localtoglobal);
96: for (i=0; i<a->A->n; i++) {
97: localtoglobal[i] = a->cstart + i + 1;
98: }
99: for (i=0; i<a->B->n; i++) {
100: localtoglobal[i+a->A->n] = a->garray[i] + 1;
101: }
102: /*
103: PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
104: */
106: /* generate the vectors needed for the local solves */
107: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xxt->b);
108: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xxt->xd);
109: VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xxt->xo);
110: xxt->nd = a->A->n;
112: /* factor the beast */
113: XXT_factor(xxt->xxt,localtoglobal,A->m,ncol,(void*)LocalMult_XXT,a);
115: VecDestroy(xxt->b);
116: VecDestroy(xxt->xd);
117: VecDestroy(xxt->xo);
118: PetscFree(localtoglobal);
120: return(0);
121: }
123: #include src/contrib/libtfs/xyt.h
125: typedef struct {
126: xyt_ADT xyt;
127: Vec b,xd,xo;
128: int nd;
129: } Mat_MPIAIJ_XYT;
132: EXTERN int MatDestroy_MPIAIJ(Mat);
134: int MatDestroy_MPIAIJ_XYT(Mat A)
135: {
136: Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
137: int ierr;
140: /* free the XYT datastructures */
141: XYT_free(xyt->xyt);
142: PetscFree(xyt);
143: MatDestroy_MPIAIJ(A);
144: return(0);
145: }
147: int MatSolve_MPIAIJ_XYT(Mat A,Vec b,Vec x)
148: {
149: Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
150: PetscScalar *bb,*xx;
151: int ierr;
154: VecGetArray(b,&bb);
155: VecGetArray(x,&xx);
156: XYT_solve(xyt->xyt,xx,bb);
157: VecRestoreArray(b,&bb);
158: VecRestoreArray(x,&xx);
159: return(0);
160: }
162: int MatLUFactorNumeric_MPIAIJ_XYT(Mat A,Mat *F)
163: {
165: return(0);
166: }
168: static int LocalMult_XYT(Mat A,PetscScalar *xin,PetscScalar *xout)
169: {
170: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
171: int ierr;
172: Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
175: VecPlaceArray(xyt->b,xout);
176: VecPlaceArray(xyt->xd,xin);
177: VecPlaceArray(xyt->xo,xin+xyt->nd);
178: MatMult(a->A,xyt->xd,xyt->b);
179: MatMultAdd(a->B,xyt->xo,xyt->b,xyt->b);
180: /*
181: PetscRealView(a->A->n+a->B->n,xin,PETSC_VIEWER_STDOUT_WORLD);
182: PetscRealView(a->A->m,xout,PETSC_VIEWER_STDOUT_WORLD);
183: */
184: return(0);
185: }
187: int MatLUFactorSymbolic_MPIAIJ_XYT(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
188: {
189: Mat B;
190: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
191: int ierr,*localtoglobal,ncol,i;
192: Mat_MPIAIJ_XYT *xyt;
195: if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
196: MatCreateMPIAIJ(A->comm,A->m,A->n,A->M,A->N,0,PETSC_NULL,0,PETSC_NULL,F);
197: B = *F;
198: B->ops->solve = MatSolve_MPIAIJ_XYT;
199: B->ops->destroy = MatDestroy_MPIAIJ_XYT;
200: B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XYT;
201: B->factor = FACTOR_LU;
202: ierr = PetscNew(Mat_MPIAIJ_XYT,&xyt);
203: B->spptr = A->spptr = (void*)xyt;
205: xyt->xyt = XYT_new();
206: /* generate the local to global mapping */
207: ncol = a->A->n + a->B->n;
208: PetscMalloc((ncol)*sizeof(int),&localtoglobal);
209: for (i=0; i<a->A->n; i++) {
210: localtoglobal[i] = a->cstart + i + 1;
211: }
212: for (i=0; i<a->B->n; i++) {
213: localtoglobal[i+a->A->n] = a->garray[i] + 1;
214: }
215: /*
216: PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
217: */
219: /* generate the vectors needed for the local solves */
220: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xyt->b);
221: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xyt->xd);
222: VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xyt->xo);
223: xyt->nd = a->A->n;
225: /* factor the beast */
226: XYT_factor(xyt->xyt,localtoglobal,A->m,ncol,(void*)LocalMult_XYT,a);
228: VecDestroy(xyt->b);
229: VecDestroy(xyt->xd);
230: VecDestroy(xyt->xo);
231: PetscFree(localtoglobal);
233: return(0);
234: }
236: int MatLUFactorSymbolic_MPIAIJ_TFS(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
237: {
241: PetscLogInfo(0,"Using TFS for MPIAIJ LU factorization and solves");
242: if (A->symmetric) {
243: MatLUFactorSymbolic_MPIAIJ_XXT(A,r,c,info,F);
244: } else {
245: MatLUFactorSymbolic_MPIAIJ_XYT(A,r,c,info,F);
246: }
247: return(0);
248: }
250: #else
252: int dummy83(int a)
253: {
255: return(0);
256: }
258: #endif