Actual source code: symtranspose.c
1: #define PETSCMAT_DLL
3: /*
4: Defines symbolic transpose routines for SeqAIJ matrices.
6: Currently Get/Restore only allocates/frees memory for holding the
7: (i,j) info for the transpose. Someday, this info could be
8: maintained so successive calls to Get will not recompute the info.
10: Also defined is a "faster" implementation of MatTranspose for SeqAIJ
11: matrices which avoids calls to MatSetValues. This routine has not
12: been adopted as the standard yet as it is somewhat untested.
14: */
16: #include ../src/mat/impls/aij/seq/aij.h
21: PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat A,PetscInt *Ati[],PetscInt *Atj[])
22: {
24: PetscInt i,j,anzj;
25: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data;
26: PetscInt an=A->cmap->N,am=A->rmap->N;
27: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
31: PetscInfo(A,"Getting Symbolic Transpose.\n");
33: /* Set up timers */
34: PetscLogEventBegin(MAT_Getsymtranspose,A,0,0,0);
36: /* Allocate space for symbolic transpose info and work array */
37: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
38: PetscMalloc(ai[am]*sizeof(PetscInt),&atj);
39: PetscMalloc(an*sizeof(PetscInt),&atfill);
40: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
42: /* Walk through aj and count ## of non-zeros in each row of A^T. */
43: /* Note: offset by 1 for fast conversion into csr format. */
44: for (i=0;i<ai[am];i++) {
45: ati[aj[i]+1] += 1;
46: }
47: /* Form ati for csr format of A^T. */
48: for (i=0;i<an;i++) {
49: ati[i+1] += ati[i];
50: }
52: /* Copy ati into atfill so we have locations of the next free space in atj */
53: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
55: /* Walk through A row-wise and mark nonzero entries of A^T. */
56: for (i=0;i<am;i++) {
57: anzj = ai[i+1] - ai[i];
58: for (j=0;j<anzj;j++) {
59: atj[atfill[*aj]] = i;
60: atfill[*aj++] += 1;
61: }
62: }
64: /* Clean up temporary space and complete requests. */
65: PetscFree(atfill);
66: *Ati = ati;
67: *Atj = atj;
69: PetscLogEventEnd(MAT_Getsymtranspose,A,0,0,0);
70: return(0);
71: }
72: /*
73: MatGetSymbolicTransposeReduced_SeqAIJ() - Get symbolic matrix structure of submatrix A[rstart:rend,:],
74: modified from MatGetSymbolicTranspose_SeqAIJ()
75: */
78: PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A,PetscInt rstart,PetscInt rend,PetscInt *Ati[],PetscInt *Atj[])
79: {
81: PetscInt i,j,anzj;
82: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data;
83: PetscInt an=A->cmap->N;
84: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
87: PetscInfo(A,"Getting Symbolic Transpose\n");
88: PetscLogEventBegin(MAT_Getsymtransreduced,A,0,0,0);
90: /* Allocate space for symbolic transpose info and work array */
91: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
92: anzj = ai[rend] - ai[rstart];
93: PetscMalloc((anzj+1)*sizeof(PetscInt),&atj);
94: PetscMalloc((an+1)*sizeof(PetscInt),&atfill);
95: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
97: /* Walk through aj and count ## of non-zeros in each row of A^T. */
98: /* Note: offset by 1 for fast conversion into csr format. */
99: for (i=ai[rstart]; i<ai[rend]; i++) {
100: ati[aj[i]+1] += 1;
101: }
102: /* Form ati for csr format of A^T. */
103: for (i=0;i<an;i++) {
104: ati[i+1] += ati[i];
105: }
107: /* Copy ati into atfill so we have locations of the next free space in atj */
108: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
110: /* Walk through A row-wise and mark nonzero entries of A^T. */
111: aj = aj + ai[rstart];
112: for (i=rstart; i<rend; i++) {
113: anzj = ai[i+1] - ai[i];
114: for (j=0;j<anzj;j++) {
115: atj[atfill[*aj]] = i-rstart;
116: atfill[*aj++] += 1;
117: }
118: }
120: /* Clean up temporary space and complete requests. */
121: PetscFree(atfill);
122: *Ati = ati;
123: *Atj = atj;
125: PetscLogEventEnd(MAT_Getsymtransreduced,A,0,0,0);
126: return(0);
127: }
131: PetscErrorCode MatTranspose_SeqAIJ_FAST(Mat A,MatReuse reuse,Mat *B)
132: {
134: PetscInt i,j,anzj;
135: Mat At;
136: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data,*at;
137: PetscInt an=A->cmap->N,am=A->rmap->N;
138: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
139: MatScalar *ata,*aa=a->a;
143: PetscLogEventBegin(MAT_Transpose_SeqAIJ,A,0,0,0);
145: if (reuse == MAT_INITIAL_MATRIX || *B == A) {
146: /* Allocate space for symbolic transpose info and work array */
147: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
148: PetscMalloc(ai[am]*sizeof(PetscInt),&atj);
149: PetscMalloc(ai[am]*sizeof(MatScalar),&ata);
150: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
151: /* Walk through aj and count ## of non-zeros in each row of A^T. */
152: /* Note: offset by 1 for fast conversion into csr format. */
153: for (i=0;i<ai[am];i++) {
154: ati[aj[i]+1] += 1;
155: }
156: /* Form ati for csr format of A^T. */
157: for (i=0;i<an;i++) {
158: ati[i+1] += ati[i];
159: }
160: } else {
161: Mat_SeqAIJ *sub_B = (Mat_SeqAIJ*) (*B)->data;
162: ati = sub_B->i;
163: atj = sub_B->j;
164: ata = sub_B->a;
165: At = *B;
166: }
168: /* Copy ati into atfill so we have locations of the next free space in atj */
169: PetscMalloc(an*sizeof(PetscInt),&atfill);
170: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
172: /* Walk through A row-wise and mark nonzero entries of A^T. */
173: for (i=0;i<am;i++) {
174: anzj = ai[i+1] - ai[i];
175: for (j=0;j<anzj;j++) {
176: atj[atfill[*aj]] = i;
177: ata[atfill[*aj]] = *aa++;
178: atfill[*aj++] += 1;
179: }
180: }
182: /* Clean up temporary space and complete requests. */
183: PetscFree(atfill);
184: if (reuse == MAT_INITIAL_MATRIX) {
185: MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,an,am,ati,atj,ata,&At);
186: at = (Mat_SeqAIJ *)(At->data);
187: at->free_a = PETSC_TRUE;
188: at->free_ij = PETSC_TRUE;
189: at->nonew = 0;
190: }
192: if (reuse == MAT_INITIAL_MATRIX || *B != A) {
193: *B = At;
194: } else {
195: MatHeaderCopy(A,At);
196: }
197: PetscLogEventEnd(MAT_Transpose_SeqAIJ,A,0,0,0);
198: return(0);
199: }
203: PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A,PetscInt *ati[],PetscInt *atj[])
204: {
208: PetscInfo(A,"Restoring Symbolic Transpose.\n");
209: PetscFree(*ati);
210: ati = PETSC_NULL;
211: PetscFree(*atj);
212: atj = PETSC_NULL;
213: return(0);
214: }