Actual source code: symtranspose.c
petsc-3.3-p7 2013-05-11
2: /*
3: Defines symbolic transpose routines for SeqAIJ matrices.
5: Currently Get/Restore only allocates/frees memory for holding the
6: (i,j) info for the transpose. Someday, this info could be
7: maintained so successive calls to Get will not recompute the info.
9: Also defined is a "faster" implementation of MatTranspose for SeqAIJ
10: matrices which avoids calls to MatSetValues. This routine has not
11: been adopted as the standard yet as it is somewhat untested.
13: */
15: #include <../src/mat/impls/aij/seq/aij.h>
20: PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat A,PetscInt *Ati[],PetscInt *Atj[])
21: {
23: PetscInt i,j,anzj;
24: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data;
25: PetscInt an=A->cmap->N,am=A->rmap->N;
26: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
29: PetscInfo(A,"Getting Symbolic Transpose.\n");
31: /* Set up timers */
32: PetscLogEventBegin(MAT_Getsymtranspose,A,0,0,0);
34: /* Allocate space for symbolic transpose info and work array */
35: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
36: PetscMalloc(ai[am]*sizeof(PetscInt),&atj);
37: PetscMalloc(an*sizeof(PetscInt),&atfill);
38: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
40: /* Walk through aj and count ## of non-zeros in each row of A^T. */
41: /* Note: offset by 1 for fast conversion into csr format. */
42: for (i=0;i<ai[am];i++) {
43: ati[aj[i]+1] += 1;
44: }
45: /* Form ati for csr format of A^T. */
46: for (i=0;i<an;i++) {
47: ati[i+1] += ati[i];
48: }
50: /* Copy ati into atfill so we have locations of the next free space in atj */
51: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
53: /* Walk through A row-wise and mark nonzero entries of A^T. */
54: for (i=0;i<am;i++) {
55: anzj = ai[i+1] - ai[i];
56: for (j=0;j<anzj;j++) {
57: atj[atfill[*aj]] = i;
58: atfill[*aj++] += 1;
59: }
60: }
62: /* Clean up temporary space and complete requests. */
63: PetscFree(atfill);
64: *Ati = ati;
65: *Atj = atj;
67: PetscLogEventEnd(MAT_Getsymtranspose,A,0,0,0);
68: return(0);
69: }
70: /*
71: MatGetSymbolicTransposeReduced_SeqAIJ() - Get symbolic matrix structure of submatrix A[rstart:rend,:],
72: modified from MatGetSymbolicTranspose_SeqAIJ()
73: */
76: PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A,PetscInt rstart,PetscInt rend,PetscInt *Ati[],PetscInt *Atj[])
77: {
79: PetscInt i,j,anzj;
80: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data;
81: PetscInt an=A->cmap->N;
82: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
85: PetscInfo(A,"Getting Symbolic Transpose\n");
86: PetscLogEventBegin(MAT_Getsymtransreduced,A,0,0,0);
88: /* Allocate space for symbolic transpose info and work array */
89: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
90: anzj = ai[rend] - ai[rstart];
91: PetscMalloc((anzj+1)*sizeof(PetscInt),&atj);
92: PetscMalloc((an+1)*sizeof(PetscInt),&atfill);
93: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
95: /* Walk through aj and count ## of non-zeros in each row of A^T. */
96: /* Note: offset by 1 for fast conversion into csr format. */
97: for (i=ai[rstart]; i<ai[rend]; i++) {
98: ati[aj[i]+1] += 1;
99: }
100: /* Form ati for csr format of A^T. */
101: for (i=0;i<an;i++) {
102: ati[i+1] += ati[i];
103: }
105: /* Copy ati into atfill so we have locations of the next free space in atj */
106: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
108: /* Walk through A row-wise and mark nonzero entries of A^T. */
109: aj = aj + ai[rstart];
110: for (i=rstart; i<rend; i++) {
111: anzj = ai[i+1] - ai[i];
112: for (j=0;j<anzj;j++) {
113: atj[atfill[*aj]] = i-rstart;
114: atfill[*aj++] += 1;
115: }
116: }
118: /* Clean up temporary space and complete requests. */
119: PetscFree(atfill);
120: *Ati = ati;
121: *Atj = atj;
123: PetscLogEventEnd(MAT_Getsymtransreduced,A,0,0,0);
124: return(0);
125: }
129: PetscErrorCode MatTranspose_SeqAIJ_FAST(Mat A,MatReuse reuse,Mat *B)
130: {
132: PetscInt i,j,anzj;
133: Mat At;
134: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data,*at;
135: PetscInt an=A->cmap->N,am=A->rmap->N;
136: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
137: MatScalar *ata,*aa=a->a;
141: PetscLogEventBegin(MAT_Transpose_SeqAIJ,A,0,0,0);
143: if (reuse == MAT_INITIAL_MATRIX || *B == A) {
144: /* Allocate space for symbolic transpose info and work array */
145: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
146: PetscMalloc(ai[am]*sizeof(PetscInt),&atj);
147: PetscMalloc(ai[am]*sizeof(MatScalar),&ata);
148: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
149: /* Walk through aj and count ## of non-zeros in each row of A^T. */
150: /* Note: offset by 1 for fast conversion into csr format. */
151: for (i=0;i<ai[am];i++) {
152: ati[aj[i]+1] += 1;
153: }
154: /* Form ati for csr format of A^T. */
155: for (i=0;i<an;i++) {
156: ati[i+1] += ati[i];
157: }
158: } else {
159: Mat_SeqAIJ *sub_B = (Mat_SeqAIJ*) (*B)->data;
160: ati = sub_B->i;
161: atj = sub_B->j;
162: ata = sub_B->a;
163: At = *B;
164: }
166: /* Copy ati into atfill so we have locations of the next free space in atj */
167: PetscMalloc(an*sizeof(PetscInt),&atfill);
168: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
170: /* Walk through A row-wise and mark nonzero entries of A^T. */
171: for (i=0;i<am;i++) {
172: anzj = ai[i+1] - ai[i];
173: for (j=0;j<anzj;j++) {
174: atj[atfill[*aj]] = i;
175: ata[atfill[*aj]] = *aa++;
176: atfill[*aj++] += 1;
177: }
178: }
180: /* Clean up temporary space and complete requests. */
181: PetscFree(atfill);
182: if (reuse == MAT_INITIAL_MATRIX) {
183: MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,an,am,ati,atj,ata,&At);
184: at = (Mat_SeqAIJ *)(At->data);
185: at->free_a = PETSC_TRUE;
186: at->free_ij = PETSC_TRUE;
187: at->nonew = 0;
188: }
190: if (reuse == MAT_INITIAL_MATRIX || *B != A) {
191: *B = At;
192: } else {
193: MatHeaderMerge(A,At);
194: }
195: PetscLogEventEnd(MAT_Transpose_SeqAIJ,A,0,0,0);
196: return(0);
197: }
201: PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A,PetscInt *ati[],PetscInt *atj[])
202: {
206: PetscInfo(A,"Restoring Symbolic Transpose.\n");
207: PetscFree(*ati);
208: PetscFree(*atj);
209: return(0);
210: }