xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 7cc688d72b4457f28d1bfbdac8a417ba5f4246b6)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2a64fbb32SBarry Smith 
36eaac0f3SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h"
4a64fbb32SBarry Smith 
5dfbe8321SBarry Smith EXTERN PetscErrorCode CreateColmap_MPIAIJ_Private(Mat);
6b1d57f15SBarry Smith EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt*[],PetscInt*[],PetscTruth*);
7b1d57f15SBarry Smith EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt*[],PetscInt*[],PetscTruth*);
8a64fbb32SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
12a64fbb32SBarry Smith {
136eaac0f3SBarry Smith   Mat_MPIAIJ            *aij = (Mat_MPIAIJ*)mat->data;
146849ba73SBarry Smith   PetscErrorCode        ierr;
15b1d57f15SBarry Smith   PetscMPIInt           size,*ncolsonproc,*disp,nn;
16b1d57f15SBarry Smith   PetscInt              i,*is,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col;
17b1d57f15SBarry Smith   PetscInt              nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj;
18899cda47SBarry Smith   PetscInt              *rowhit,M = mat->rmap.n,cstart = mat->cmap.rstart,cend = mat->cmap.rend,colb;
19b1d57f15SBarry Smith   PetscInt              *columnsforrow,l;
20b9617806SBarry Smith   IS                    *isa;
21f1af5d2fSBarry Smith   PetscTruth             done,flg;
22b8f8c88eSHong Zhang   ISLocalToGlobalMapping map = mat->mapping;
23*7cc688d7SBarry Smith   PetscInt               *ltog = (map ? map->indices : 0) ,ctype=c->ctype;
24a64fbb32SBarry Smith 
253a40ed3dSBarry Smith   PetscFunctionBegin;
26522c5e43SBarry Smith   if (!mat->assembled) {
2729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
28522c5e43SBarry Smith   }
29*7cc688d7SBarry Smith   if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
30522c5e43SBarry Smith 
31b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
32899cda47SBarry Smith   c->M             = mat->rmap.N;  /* set the global rows and columns and local rows */
33899cda47SBarry Smith   c->N             = mat->cmap.N;
34899cda47SBarry Smith   c->m             = mat->rmap.n;
35899cda47SBarry Smith   c->rstart        = mat->rmap.rstart;
36005c665bSBarry Smith 
37a64fbb32SBarry Smith   c->ncolors       = nis;
38b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
39b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
40b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
41b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
42b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
4352e6d16bSBarry Smith   ierr = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
446eaac0f3SBarry Smith 
456eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
466eaac0f3SBarry Smith   if (!aij->colmap) {
476eaac0f3SBarry Smith     ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
486eaac0f3SBarry Smith   }
4943a90d84SBarry Smith   /*
5043a90d84SBarry Smith       Calls the _SeqAIJ() version of these routines to make sure it does not
5143a90d84SBarry Smith      get the reduced (by inodes) version of I and J
5243a90d84SBarry Smith   */
5343a90d84SBarry Smith   ierr = MatGetColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
5443a90d84SBarry Smith   ierr = MatGetColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
556eaac0f3SBarry Smith 
56b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
57b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
586eaac0f3SBarry Smith 
59a64fbb32SBarry Smith   for (i=0; i<nis; i++) {
60b9b97703SBarry Smith     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
61a64fbb32SBarry Smith     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
62a64fbb32SBarry Smith     c->ncolumns[i] = n;
63a64fbb32SBarry Smith     if (n) {
64b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
6552e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr);
66b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
67a64fbb32SBarry Smith     } else {
68a64fbb32SBarry Smith       c->columns[i]  = 0;
69a64fbb32SBarry Smith     }
70a64fbb32SBarry Smith 
71b8f8c88eSHong Zhang     if (ctype == IS_COLORING_LOCAL){
726eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
73b8f8c88eSHong Zhang       ierr = MPI_Comm_size(mat->comm,&size);CHKERRQ(ierr);
74b8f8c88eSHong Zhang       ierr = PetscMalloc(2*size*sizeof(PetscInt*),&ncolsonproc);CHKERRQ(ierr);
75b8f8c88eSHong Zhang       disp = ncolsonproc + size;
76b8f8c88eSHong Zhang 
77b1d57f15SBarry Smith       nn   = (PetscMPIInt)n;
78b1d57f15SBarry Smith       ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,mat->comm);CHKERRQ(ierr);
796eaac0f3SBarry Smith       nctot = 0; for (j=0; j<size; j++) {nctot += ncolsonproc[j];}
803a7fca6bSBarry Smith       if (!nctot) {
81ae15b995SBarry Smith         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
823a7fca6bSBarry Smith       }
836eaac0f3SBarry Smith 
846eaac0f3SBarry Smith       disp[0] = 0;
856eaac0f3SBarry Smith       for (j=1; j<size; j++) {
866eaac0f3SBarry Smith         disp[j] = disp[j-1] + ncolsonproc[j-1];
876eaac0f3SBarry Smith       }
886eaac0f3SBarry Smith 
896eaac0f3SBarry Smith       /* Get complete list of columns for color on each processor */
90b1d57f15SBarry Smith       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
91b1d57f15SBarry Smith       ierr = MPI_Allgatherv(is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,mat->comm);CHKERRQ(ierr);
92b8f8c88eSHong Zhang       ierr = PetscFree(ncolsonproc);CHKERRQ(ierr);
93b8f8c88eSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED){
94b8f8c88eSHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
95b8f8c88eSHong Zhang       nctot = n;
96b8f8c88eSHong Zhang       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
97b8f8c88eSHong Zhang       ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
98b8f8c88eSHong Zhang     } else {
99b8f8c88eSHong Zhang       SETERRQ(PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
100b8f8c88eSHong Zhang     }
1016eaac0f3SBarry Smith 
1026eaac0f3SBarry Smith     /*
1036eaac0f3SBarry Smith        Mark all rows affect by these columns
1046eaac0f3SBarry Smith     */
105b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
106b8f8c88eSHong Zhang     ierr = PetscOptionsHasName(PETSC_NULL,"-matfdcoloring_slow",&flg);CHKERRQ(ierr);
107f158e583SBarry Smith     if (!flg) {/*-----------------------------------------------------------------------------*/
108f158e583SBarry Smith       /* crude, fast version */
109b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
110a64fbb32SBarry Smith       /* loop over columns*/
1116eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
112b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
113b8f8c88eSHong Zhang           col = ltog[cols[j]];
114b8f8c88eSHong Zhang         } else {
1156eaac0f3SBarry Smith           col  = cols[j];
116b8f8c88eSHong Zhang         }
1176eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1186eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1196eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1206eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1216eaac0f3SBarry Smith         } else {
122aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
1230f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr)
124fa46199cSSatish Balay 	  colb --;
125b3d2dc96SSatish Balay #else
1266eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
127b3d2dc96SSatish Balay #endif
1286eaac0f3SBarry Smith           if (colb == -1) {
1296eaac0f3SBarry Smith             m = 0;
1306eaac0f3SBarry Smith           } else {
1316eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1326eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1336eaac0f3SBarry Smith           }
1346eaac0f3SBarry Smith         }
135a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
136a64fbb32SBarry Smith         for (k=0; k<m; k++) {
137a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
138a64fbb32SBarry Smith         }
139a64fbb32SBarry Smith       }
1406eaac0f3SBarry Smith 
141a64fbb32SBarry Smith       /* count the number of hits */
142a64fbb32SBarry Smith       nrows = 0;
1436eaac0f3SBarry Smith       for (j=0; j<M; j++) {
144a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
145a64fbb32SBarry Smith       }
146a64fbb32SBarry Smith       c->nrows[i]         = nrows;
147b1d57f15SBarry Smith       ierr                = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
148b1d57f15SBarry Smith       ierr                = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
14952e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
150a64fbb32SBarry Smith       nrows = 0;
1516eaac0f3SBarry Smith       for (j=0; j<M; j++) {
152a64fbb32SBarry Smith         if (rowhit[j]) {
153a64fbb32SBarry Smith           c->rows[i][nrows]           = j;
154a64fbb32SBarry Smith           c->columnsforrow[i][nrows] = rowhit[j] - 1;
155a64fbb32SBarry Smith           nrows++;
156a64fbb32SBarry Smith         }
157a64fbb32SBarry Smith       }
158a64fbb32SBarry Smith     } else {/*-------------------------------------------------------------------------------*/
159f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
160b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
1616eaac0f3SBarry Smith       rowhit[M] = M;
162a64fbb32SBarry Smith       nrows     = 0;
163a64fbb32SBarry Smith       /* loop over columns*/
1646eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
165b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
166b8f8c88eSHong Zhang           col = ltog[cols[j]];
167b8f8c88eSHong Zhang         } else {
1686eaac0f3SBarry Smith           col  = cols[j];
169b8f8c88eSHong Zhang         }
1706eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1716eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1726eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1736eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1746eaac0f3SBarry Smith         } else {
175aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
1760f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
177fa46199cSSatish Balay           colb --;
178b3d2dc96SSatish Balay #else
1796eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
180b3d2dc96SSatish Balay #endif
1816eaac0f3SBarry Smith           if (colb == -1) {
1826eaac0f3SBarry Smith             m = 0;
1836eaac0f3SBarry Smith           } else {
1846eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1856eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1866eaac0f3SBarry Smith           }
1876eaac0f3SBarry Smith         }
188b8f8c88eSHong Zhang 
189a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
1906eaac0f3SBarry Smith         fm    = M; /* fm points to first entry in linked list */
191a64fbb32SBarry Smith         for (k=0; k<m; k++) {
192a64fbb32SBarry Smith           currentcol = *rows++;
193a64fbb32SBarry Smith 	  /* is it already in the list? */
194a64fbb32SBarry Smith           do {
195a64fbb32SBarry Smith             mfm  = fm;
196a64fbb32SBarry Smith             fm   = rowhit[fm];
197a64fbb32SBarry Smith           } while (fm < currentcol);
198a64fbb32SBarry Smith           /* not in list so add it */
199a64fbb32SBarry Smith           if (fm != currentcol) {
200a64fbb32SBarry Smith             nrows++;
201a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
202a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
203a64fbb32SBarry Smith             rowhit[mfm]               = currentcol;
204a64fbb32SBarry Smith             rowhit[currentcol]        = fm;
205a64fbb32SBarry Smith             fm                        = currentcol;
206a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
207a64fbb32SBarry Smith           } else {
20829bbc08cSBarry Smith             SETERRQ(PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
209a64fbb32SBarry Smith           }
210a64fbb32SBarry Smith         }
211a64fbb32SBarry Smith       }
212a64fbb32SBarry Smith       c->nrows[i]         = nrows;
213b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
214b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
21552e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
216a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
217a64fbb32SBarry Smith       nrows = 0;
2186eaac0f3SBarry Smith       fm    = rowhit[M];
219a64fbb32SBarry Smith       do {
220a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
221a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
222a64fbb32SBarry Smith         fm                           = rowhit[fm];
2236eaac0f3SBarry Smith       } while (fm < M);
2246eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
225606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
2266eaac0f3SBarry Smith   }
22730b34957SBarry Smith 
22830b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
22930b34957SBarry Smith   /*
23030b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
23130b34957SBarry Smith   */
232b8f8c88eSHong Zhang   if (ctype == IS_COLORING_LOCAL) {
233b8f8c88eSHong Zhang     ierr = VecCreateGhost(mat->comm,aij->A->rmap.n,PETSC_DETERMINE,aij->B->cmap.n,aij->garray,&c->vscale);CHKERRQ(ierr);
234b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
23530b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
236b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
23730b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
23830b34957SBarry Smith         col = c->columnsforrow[k][l];
23930b34957SBarry Smith         if (col >= cstart && col < cend) {
24030b34957SBarry Smith           /* column is in diagonal block of matrix */
24130b34957SBarry Smith           colb = col - cstart;
24230b34957SBarry Smith         } else {
24330b34957SBarry Smith           /* column  is in "off-processor" part */
24430b34957SBarry Smith #if defined (PETSC_USE_CTABLE)
24530b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
24630b34957SBarry Smith           colb --;
24730b34957SBarry Smith #else
24830b34957SBarry Smith           colb = aij->colmap[col] - 1;
24930b34957SBarry Smith #endif
25030b34957SBarry Smith           colb += cend - cstart;
25130b34957SBarry Smith         }
25230b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
25330b34957SBarry Smith       }
25430b34957SBarry Smith     }
255b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
256b8f8c88eSHong Zhang     /* Get gtol mapping */
257b8f8c88eSHong Zhang     PetscInt N = mat->cmap.N, *gtol;
258b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
259b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
260b8f8c88eSHong Zhang     for (i=0; i<map->n; i++) gtol[ltog[i]] = i;
261b8f8c88eSHong Zhang 
262b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
263b8f8c88eSHong Zhang     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
264b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
265b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
266b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
267b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
268b8f8c88eSHong Zhang         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
269b8f8c88eSHong Zhang       }
270b8f8c88eSHong Zhang     }
271b8f8c88eSHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
272b8f8c88eSHong Zhang   }
273b9617806SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
27430b34957SBarry Smith 
275606d414cSSatish Balay   ierr = PetscFree(rowhit);CHKERRQ(ierr);
276606d414cSSatish Balay   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
27743a90d84SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
27843a90d84SBarry Smith   ierr = MatRestoreColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
2793a40ed3dSBarry Smith   PetscFunctionReturn(0);
280a64fbb32SBarry Smith }
281a64fbb32SBarry Smith 
282b9617806SBarry Smith 
283b9617806SBarry Smith 
284b9617806SBarry Smith 
285b9617806SBarry Smith 
286b9617806SBarry Smith 
287