Actual source code: fdaij.c
1: /*$Id: fdaij.c,v 1.40 2001/06/21 21:16:21 bsmith Exp $*/
3: #include src/mat/impls/aij/seq/aij.h
5: EXTERN int MatGetColumnIJ_SeqAIJ(Mat,int,PetscTruth,int*,int*[],int*[],PetscTruth*);
6: EXTERN int MatRestoreColumnIJ_SeqAIJ(Mat,int,PetscTruth,int*,int*[],int*[],PetscTruth*);
10: int MatFDColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
11: {
12: int i,*is,n,nrows,N = mat->N,j,k,m,*rows,ierr,*ci,*cj,ncols,col;
13: int nis = iscoloring->n,*rowhit,*columnsforrow,l;
14: IS *isa;
15: PetscTruth done,flg;
18: if (!mat->assembled) {
19: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
20: }
22: ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);
23: c->M = mat->M; /* set total rows, columns and local rows */
24: c->N = mat->N;
25: c->m = mat->M;
26: c->rstart = 0;
28: c->ncolors = nis;
29: PetscMalloc(nis*sizeof(int),&c->ncolumns);
30: PetscMalloc(nis*sizeof(int*),&c->columns);
31: PetscMalloc(nis*sizeof(int),&c->nrows);
32: PetscMalloc(nis*sizeof(int*),&c->rows);
33: PetscMalloc(nis*sizeof(int*),&c->columnsforrow);
35: /*
36: Calls the _SeqAIJ() version of these routines to make sure it does not
37: get the reduced (by inodes) version of I and J
38: */
39: MatGetColumnIJ_SeqAIJ(mat,0,PETSC_FALSE,&ncols,&ci,&cj,&done);
41: /*
42: Temporary option to allow for debugging/testing
43: */
44: PetscOptionsHasName(PETSC_NULL,"-matfdcoloring_slow",&flg);
46: PetscMalloc((N+1)*sizeof(int),&rowhit);
47: PetscMalloc((N+1)*sizeof(int),&columnsforrow);
49: for (i=0; i<nis; i++) {
50: ISGetLocalSize(isa[i],&n);
51: ISGetIndices(isa[i],&is);
52: c->ncolumns[i] = n;
53: if (n) {
54: PetscMalloc(n*sizeof(int),&c->columns[i]);
55: PetscMemcpy(c->columns[i],is,n*sizeof(int));
56: } else {
57: c->columns[i] = 0;
58: }
60: if (!flg) { /* ------------------------------------------------------------------------------*/
61: /* fast, crude version requires O(N*N) work */
62: PetscMemzero(rowhit,N*sizeof(int));
63: /* loop over columns*/
64: for (j=0; j<n; j++) {
65: col = is[j];
66: rows = cj + ci[col];
67: m = ci[col+1] - ci[col];
68: /* loop over columns marking them in rowhit */
69: for (k=0; k<m; k++) {
70: rowhit[*rows++] = col + 1;
71: }
72: }
73: /* count the number of hits */
74: nrows = 0;
75: for (j=0; j<N; j++) {
76: if (rowhit[j]) nrows++;
77: }
78: c->nrows[i] = nrows;
79: PetscMalloc((nrows+1)*sizeof(int),&c->rows[i]);
80: PetscMalloc((nrows+1)*sizeof(int),&c->columnsforrow[i]);
81: nrows = 0;
82: for (j=0; j<N; j++) {
83: if (rowhit[j]) {
84: c->rows[i][nrows] = j;
85: c->columnsforrow[i][nrows] = rowhit[j] - 1;
86: nrows++;
87: }
88: }
89: } else { /*-------------------------------------------------------------------------------*/
90: /* slow version, using rowhit as a linked list */
91: int currentcol,fm,mfm;
92: rowhit[N] = N;
93: nrows = 0;
94: /* loop over columns */
95: for (j=0; j<n; j++) {
96: col = is[j];
97: rows = cj + ci[col];
98: m = ci[col+1] - ci[col];
99: /* loop over columns marking them in rowhit */
100: fm = N; /* fm points to first entry in linked list */
101: for (k=0; k<m; k++) {
102: currentcol = *rows++;
103: /* is it already in the list? */
104: do {
105: mfm = fm;
106: fm = rowhit[fm];
107: } while (fm < currentcol);
108: /* not in list so add it */
109: if (fm != currentcol) {
110: nrows++;
111: columnsforrow[currentcol] = col;
112: /* next three lines insert new entry into linked list */
113: rowhit[mfm] = currentcol;
114: rowhit[currentcol] = fm;
115: fm = currentcol;
116: /* fm points to present position in list since we know the columns are sorted */
117: } else {
118: SETERRQ(PETSC_ERR_PLIB,"Detected invalid coloring");
119: }
121: }
122: }
123: c->nrows[i] = nrows;
124: PetscMalloc((nrows+1)*sizeof(int),&c->rows[i]);
125: PetscMalloc((nrows+1)*sizeof(int),&c->columnsforrow[i]);
126: /* now store the linked list of rows into c->rows[i] */
127: nrows = 0;
128: fm = rowhit[N];
129: do {
130: c->rows[i][nrows] = fm;
131: c->columnsforrow[i][nrows++] = columnsforrow[fm];
132: fm = rowhit[fm];
133: } while (fm < N);
134: } /* ---------------------------------------------------------------------------------------*/
135: ISRestoreIndices(isa[i],&is);
136: }
137: MatRestoreColumnIJ_SeqAIJ(mat,0,PETSC_FALSE,&ncols,&ci,&cj,&done);
139: PetscFree(rowhit);
140: PetscFree(columnsforrow);
142: /* Optimize by adding the vscale, and scaleforrow[][] fields */
143: /*
144: see the version for MPIAIJ
145: */
146: VecCreateGhost(mat->comm,mat->m,PETSC_DETERMINE,0,PETSC_NULL,&c->vscale);CHKERRQ(ierr)
147: PetscMalloc(c->ncolors*sizeof(int*),&c->vscaleforrow);
148: for (k=0; k<c->ncolors; k++) {
149: PetscMalloc((c->nrows[k]+1)*sizeof(int),&c->vscaleforrow[k]);
150: for (l=0; l<c->nrows[k]; l++) {
151: col = c->columnsforrow[k][l];
152: c->vscaleforrow[k][l] = col;
153: }
154: }
155: ISColoringRestoreIS(iscoloring,&isa);
156: return(0);
157: }
159: /*
160: Makes a longer coloring[] array and calls the usual code with that
161: */
164: int MatColoringPatch_SeqAIJ_Inode(Mat mat,int nin,int ncolors,const ISColoringValue coloring[],ISColoring *iscoloring)
165: {
166: Mat_SeqAIJ *a = (Mat_SeqAIJ*)mat->data;
167: int n = mat->n,ierr,m = a->inode.node_count,j,*ns = a->inode.size,row;
168: int *colorused,i;
169: ISColoringValue *newcolor;
172: PetscMalloc((n+1)*sizeof(int),&newcolor);
173: /* loop over inodes, marking a color for each column*/
174: row = 0;
175: for (i=0; i<m; i++){
176: for (j=0; j<ns[i]; j++) {
177: newcolor[row++] = coloring[i] + j*ncolors;
178: }
179: }
181: /* eliminate unneeded colors */
182: PetscMalloc(5*ncolors*sizeof(int),&colorused);
183: PetscMemzero(colorused,5*ncolors*sizeof(int));
184: for (i=0; i<n; i++) {
185: colorused[newcolor[i]] = 1;
186: }
188: for (i=1; i<5*ncolors; i++) {
189: colorused[i] += colorused[i-1];
190: }
191: ncolors = colorused[5*ncolors-1];
192: for (i=0; i<n; i++) {
193: newcolor[i] = colorused[newcolor[i]];
194: }
195: PetscFree(colorused);
196: ISColoringCreate(mat->comm,n,newcolor,iscoloring);
197: PetscFree((void*)coloring);
198: return(0);
199: }