Actual source code: mmbdiag.c
1: /*$Id: mmbdiag.c,v 1.40 2001/03/23 23:22:05 balay Exp $*/
3: /*
4: Support for the MPIBDIAG matrix-vector multiply
5: */
6: #include src/mat/impls/bdiag/mpi/mpibdiag.h
10: int MatSetUpMultiply_MPIBDiag(Mat mat)
11: {
12: Mat_MPIBDiag *mbd = (Mat_MPIBDiag*)mat->data;
13: Mat_SeqBDiag *lmbd = (Mat_SeqBDiag*)mbd->A->data;
14: int ierr,N = mat->N,*indices,*garray,ec=0;
15: int bs = lmbd->bs,d,i,j,diag;
16: IS to,from;
17: Vec gvec;
20: /* We make an array as long as the number of columns */
21: /* mark those columns that are in mbd->A */
22: PetscMalloc((N+1)*sizeof(int),&indices);
23: PetscMemzero(indices,N*sizeof(int));
25: if (bs == 1) {
26: for (d=0; d<lmbd->nd; d++) {
27: diag = lmbd->diag[d];
28: if (diag > 0) { /* col = loc */
29: for (j=0; j<lmbd->bdlen[d]; j++) {
30: if (!indices[j]) ec++;
31: indices[j] = 1;
32: }
33: } else { /* col = loc-diag */
34: for (j=0; j<lmbd->bdlen[d]; j++) {
35: if (!indices[j-diag]) ec++;
36: indices[j-diag] = 1;
37: }
38: }
39: }
40: } else {
41: for (d=0; d<lmbd->nd; d++) {
42: diag = lmbd->diag[d];
43: if (diag > 0) { /* col = loc */
44: for (j=0; j<lmbd->bdlen[d]; j++) {
45: if (!indices[bs*j]) ec += bs;
46: for (i=0; i<bs; i++) indices[bs*j+i] = 1;
47: }
48: } else { /* col = loc-diag */
49: for (j=0; j<lmbd->bdlen[d]; j++) {
50: if (!indices[bs*(j-diag)]) ec += bs;
51: for (i=0; i<bs; i++) indices[bs*(j-diag)+i] = 1;
52: }
53: }
54: }
55: }
57: /* form array of columns we need */
58: PetscMalloc((ec+1)*sizeof(int),&garray);
59: ec = 0;
60: for (i=0; i<N; i++) {
61: if (indices[i]) garray[ec++] = i;
62: }
63: PetscFree(indices);
65: /* create local vector that is used to scatter into */
66: VecCreateSeq(PETSC_COMM_SELF,N,&mbd->lvec);
68: /* create temporary index set for building scatter-gather */
69: ISCreateGeneral(mat->comm,ec,garray,&from);
70: ISCreateGeneral(PETSC_COMM_SELF,ec,garray,&to);
71: PetscFree(garray);
73: /* create temporary global vector to generate scatter context */
74: /* this is inefficient, but otherwise we must do either
75: 1) save garray until the first actual scatter when the vector is known or
76: 2) have another way of generating a scatter context without a vector.*/
77: /*
78: This is not correct for a rectangular matrix mbd->m?
79: */
80: VecCreateMPI(mat->comm,mat->m,mat->N,&gvec);
82: /* generate the scatter context */
83: VecScatterCreate(gvec,from,mbd->lvec,to,&mbd->Mvctx);
84: PetscLogObjectParent(mat,mbd->Mvctx);
85: PetscLogObjectParent(mat,mbd->lvec);
86: PetscLogObjectParent(mat,to);
87: PetscLogObjectParent(mat,from);
88: PetscLogObjectParent(mat,gvec);
90: ISDestroy(to);
91: ISDestroy(from);
92: VecDestroy(gvec);
93: return(0);
94: }