xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision e5980f420c65071add1882455d9a7ac92e5ae3d9)
1 #ifndef lint
2 static char vcid[] = "$Id: mmaij.c,v 1.11 1995/06/08 03:09:29 bsmith Exp bsmith $";
3 #endif
4 
5 
6 /*
7    Support for the parallel AIJ matrix vector multiply
8 */
9 #include "mpiaij.h"
10 #include "vec/vecimpl.h"
11 #include "../seq/aij.h"
12 
13 int MatSetUpMultiply_MPIAIJ(Mat mat)
14 {
15   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) mat->data;
16   Mat_AIJ    *B = (Mat_AIJ *) (aij->B->data);
17   int        N = aij->N,i,j,*indices,*aj = B->j;
18   int        ierr,ec = 0,*garray;
19   IS         from,to;
20   Vec        gvec;
21 
22   /* For the first stab we make an array as long as the number of columns */
23   /* mark those columns that are in aij->B */
24   indices = (int *) PETSCMALLOC( N*sizeof(int) ); CHKPTRQ(indices);
25   PETSCMEMSET(indices,0,N*sizeof(int));
26   for ( i=0; i<B->m; i++ ) {
27     for ( j=0; j<B->ilen[i]; j++ ) {
28      if (!indices[aj[B->i[i] - 1 + j]-1]) ec++;
29      indices[aj[B->i[i] - 1 + j]-1] = 1;}
30   }
31 
32   /* form array of columns we need */
33   garray = (int *) PETSCMALLOC( (ec+1)*sizeof(int) ); CHKPTRQ(garray);
34   ec = 0;
35   for ( i=0; i<N; i++ ) {
36     if (indices[i]) garray[ec++] = i;
37   }
38 
39   /* make indices now point into garray */
40   for ( i=0; i<ec; i++ ) {
41     indices[garray[i]] = i+1;
42   }
43 
44   /* compact out the extra columns in B */
45   for ( i=0; i<B->m; i++ ) {
46     for ( j=0; j<B->ilen[i]; j++ ) {
47       aj[B->i[i] - 1 + j] = indices[aj[B->i[i] - 1 + j]-1];
48     }
49   }
50   B->n = ec;
51   PETSCFREE(indices);
52 
53   /* create local vector that is used to scatter into */
54   ierr = VecCreateSequential(MPI_COMM_SELF,ec,&aij->lvec); CHKERRQ(ierr);
55 
56   /* create two temporary Index sets for build scatter gather */
57   ierr = ISCreateSequential(MPI_COMM_SELF,ec,garray,&from); CHKERRQ(ierr);
58   ierr = ISCreateStrideSequential(MPI_COMM_SELF,ec,0,1,&to); CHKERRQ(ierr);
59 
60   /* create temporary global vector to generate scatter context */
61   /* this is inefficient, but otherwise we must do either
62      1) save garray until the first actual scatter when the vector is known or
63      2) have another way of generating a scatter context without a vector.*/
64   ierr = VecCreateMPI(mat->comm,aij->n,aij->N,&gvec); CHKERRQ(ierr);
65 
66   /* gnerate the scatter context */
67   ierr = VecScatterCtxCreate(gvec,from,aij->lvec,to,&aij->Mvctx); CHKERRQ(ierr);
68   PLogObjectParent(mat,aij->Mvctx);
69   PLogObjectParent(mat,aij->lvec);
70   aij->garray = garray;
71   ierr = ISDestroy(from); CHKERRQ(ierr);
72   ierr = ISDestroy(to); CHKERRQ(ierr);
73   ierr = VecDestroy(gvec);
74   return 0;
75 }
76 
77 
78 /*
79      Takes the local part of an already assembled MPIAIJ matrix
80    and disassembles it. This is to allow new nonzeros into the matrix
81    that require more communication in the matrix vector multiply.
82    Thus certain data-structures must be rebuilt.
83 
84    Kind of slow! But that's what application programmers get when
85    they are sloppy.
86 */
87 int DisAssemble_MPIAIJ(Mat A)
88 {
89   Mat_MPIAIJ *aij = (Mat_MPIAIJ *) A->data;
90   Mat        B = aij->B,Bnew;
91   Mat_AIJ    *Baij = (Mat_AIJ*)B->data;
92   int        ierr,i,j,m=Baij->m,n = aij->N,col,ct = 0,*garray = aij->garray;
93   Scalar     v;
94 
95 fprintf(stderr,"Warning: you are adding tricky new nonzeros\n");
96 
97   /* free stuff related to matrix-vec multiply */
98   ierr = VecDestroy(aij->lvec); CHKERRQ(ierr); aij->lvec = 0;
99   ierr = VecScatterCtxDestroy(aij->Mvctx); CHKERRQ(ierr); aij->Mvctx = 0;
100   PETSCFREE(aij->colmap); aij->colmap = 0;
101 
102   /* make sure that B is assembled so we can access its values */
103   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
104   MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
105 
106   /* invent new B and copy stuff over */
107   ierr = MatCreateSequentialAIJ(MPI_COMM_SELF,m,n,0,0,&Bnew); CHKERRQ(ierr);
108   for ( i=0; i<m; i++ ) {
109     for ( j=Baij->i[i]-1; j<Baij->i[i+1]-1; j++ ) {
110       col = garray[Baij->j[ct]-1];
111       v = Baij->a[ct++];
112       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,INSERTVALUES); CHKERRQ(ierr);
113     }
114   }
115   PETSCFREE(aij->garray); aij->garray = 0;
116   ierr = MatDestroy(B); CHKERRQ(ierr);
117   aij->B = Bnew;
118   aij->assembled = 0;
119   return 0;
120 }
121 
122