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: #undef __FUNCT__
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: }
36: #undef __FUNCT__
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: }
53: #undef __FUNCT__
55: int MatLUFactorNumeric_MPIAIJ_XXT(Mat A,Mat *F)
56: {
58: return(0);
59: }
61: #undef __FUNCT__
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: }
82: #undef __FUNCT__
84: int MatLUFactorSymbolic_MPIAIJ_XXT(Mat A,IS r,IS c,MatLUInfo *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: MatCreateMPIAIJ(A->comm,A->m,A->n,A->M,A->N,0,PETSC_NULL,0,PETSC_NULL,F);
94: B = *F;
95: B->ops->solve = MatSolve_MPIAIJ_XXT;
96: B->ops->destroy = MatDestroy_MPIAIJ_XXT;
97: B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XXT;
98: B->factor = FACTOR_LU;
99: B->assembled = PETSC_TRUE;
100: ierr = PetscNew(Mat_MPIAIJ_XXT,&xxt);
101: B->spptr = A->spptr = (void*)xxt;
103: xxt->xxt = XXT_new();
104: /* generate the local to global mapping */
105: ncol = a->A->n + a->B->n;
106: PetscMalloc((ncol)*sizeof(int),&localtoglobal);
107: for (i=0; i<a->A->n; i++) {
108: localtoglobal[i] = a->cstart + i + 1;
109: }
110: for (i=0; i<a->B->n; i++) {
111: localtoglobal[i+a->A->n] = a->garray[i] + 1;
112: }
113: /*
114: PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
115: */
117: /* generate the vectors needed for the local solves */
118: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xxt->b);
119: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xxt->xd);
120: VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xxt->xo);
121: xxt->nd = a->A->n;
123: /* factor the beast */
124: XXT_factor(xxt->xxt,localtoglobal,A->m,ncol,(void*)LocalMult_XXT,a);
126: VecDestroy(xxt->b);
127: VecDestroy(xxt->xd);
128: VecDestroy(xxt->xo);
129: PetscFree(localtoglobal);
131: return(0);
132: }
134: #include src/contrib/libtfs/xyt.h
136: typedef struct {
137: xyt_ADT xyt;
138: Vec b,xd,xo;
139: int nd;
140: } Mat_MPIAIJ_XYT;
143: EXTERN int MatDestroy_MPIAIJ(Mat);
145: #undef __FUNCT__
147: int MatDestroy_MPIAIJ_XYT(Mat A)
148: {
149: Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
150: int ierr;
153: /* free the XYT datastructures */
154: XYT_free(xyt->xyt);
155: PetscFree(xyt);
156: MatDestroy_MPIAIJ(A);
157: return(0);
158: }
160: #undef __FUNCT__
162: int MatSolve_MPIAIJ_XYT(Mat A,Vec b,Vec x)
163: {
164: Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
165: PetscScalar *bb,*xx;
166: int ierr;
169: VecGetArray(b,&bb);
170: VecGetArray(x,&xx);
171: XYT_solve(xyt->xyt,xx,bb);
172: VecRestoreArray(b,&bb);
173: VecRestoreArray(x,&xx);
174: return(0);
175: }
177: #undef __FUNCT__
179: int MatLUFactorNumeric_MPIAIJ_XYT(Mat A,Mat *F)
180: {
182: return(0);
183: }
185: #undef __FUNCT__
187: static int LocalMult_XYT(Mat A,PetscScalar *xin,PetscScalar *xout)
188: {
189: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
190: int ierr;
191: Mat_MPIAIJ_XYT *xyt = (Mat_MPIAIJ_XYT*)A->spptr;
194: VecPlaceArray(xyt->b,xout);
195: VecPlaceArray(xyt->xd,xin);
196: VecPlaceArray(xyt->xo,xin+xyt->nd);
197: MatMult(a->A,xyt->xd,xyt->b);
198: MatMultAdd(a->B,xyt->xo,xyt->b,xyt->b);
199: /*
200: PetscRealView(a->A->n+a->B->n,xin,PETSC_VIEWER_STDOUT_WORLD);
201: PetscRealView(a->A->m,xout,PETSC_VIEWER_STDOUT_WORLD);
202: */
203: return(0);
204: }
206: #undef __FUNCT__
208: int MatLUFactorSymbolic_MPIAIJ_XYT(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
209: {
210: Mat B;
211: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
212: int ierr,*localtoglobal,ncol,i;
213: Mat_MPIAIJ_XYT *xyt;
216: if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
217: MatCreateMPIAIJ(A->comm,A->m,A->n,A->M,A->N,0,PETSC_NULL,0,PETSC_NULL,F);
218: B = *F;
219: B->ops->solve = MatSolve_MPIAIJ_XYT;
220: B->ops->destroy = MatDestroy_MPIAIJ_XYT;
221: B->ops->lufactornumeric = MatLUFactorNumeric_MPIAIJ_XYT;
222: B->factor = FACTOR_LU;
223: B->assembled = PETSC_TRUE;
224: ierr = PetscNew(Mat_MPIAIJ_XYT,&xyt);
225: B->spptr = A->spptr = (void*)xyt;
227: xyt->xyt = XYT_new();
228: /* generate the local to global mapping */
229: ncol = a->A->n + a->B->n;
230: PetscMalloc((ncol)*sizeof(int),&localtoglobal);
231: for (i=0; i<a->A->n; i++) {
232: localtoglobal[i] = a->cstart + i + 1;
233: }
234: for (i=0; i<a->B->n; i++) {
235: localtoglobal[i+a->A->n] = a->garray[i] + 1;
236: }
237: /*
238: PetscIntView(ncol,localtoglobal,PETSC_VIEWER_STDOUT_WORLD);
239: */
241: /* generate the vectors needed for the local solves */
242: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->m,PETSC_NULL,&xyt->b);
243: VecCreateSeqWithArray(PETSC_COMM_SELF,a->A->n,PETSC_NULL,&xyt->xd);
244: VecCreateSeqWithArray(PETSC_COMM_SELF,a->B->n,PETSC_NULL,&xyt->xo);
245: xyt->nd = a->A->n;
247: /* factor the beast */
248: XYT_factor(xyt->xyt,localtoglobal,A->m,ncol,(void*)LocalMult_XYT,A);
250: VecDestroy(xyt->b);
251: VecDestroy(xyt->xd);
252: VecDestroy(xyt->xo);
253: PetscFree(localtoglobal);
255: return(0);
256: }
258: #undef __FUNCT__
260: int MatLUFactorSymbolic_MPIAIJ_TFS(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
261: {
265: PetscLogInfo(0,"Using TFS for MPIAIJ LU factorization and solves");
266: if (A->symmetric) {
267: MatLUFactorSymbolic_MPIAIJ_XXT(A,r,c,info,F);
268: } else {
269: MatLUFactorSymbolic_MPIAIJ_XYT(A,r,c,info,F);
270: }
271: return(0);
272: }
274: #else
276: #undef __FUNCT__
278: int dummy83(int a)
279: {
281: return(0);
282: }
284: #endif