Actual source code: symtranspose.c
1: /*
2: Defines symbolic transpose routines for SeqAIJ matrices.
4: Currently Get/Restore only allocates/frees memory for holding the
5: (i,j) info for the transpose. Someday, this info could be
6: maintained so successive calls to Get will not recompute the info.
8: Also defined is a "faster" implementation of MatTranspose for SeqAIJ
9: matrices which avoids calls to MatSetValues. This routine has not
10: been adopted as the standard yet as it is somewhat untested.
12: */
14: #include src/mat/impls/aij/seq/aij.h
16: static PetscEvent logkey_matgetsymtranspose = 0;
17: static PetscEvent logkey_mattranspose = 0;
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->N,am=A->M;
27: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
31: PetscLogInfo(A,"Getting Symbolic Transpose.\n");
33: /* Set up timers */
34: if (!logkey_matgetsymtranspose) {
35: PetscLogEventRegister(&logkey_matgetsymtranspose,"MatGetSymbolicTranspose",MAT_COOKIE);
36: }
37: PetscLogEventBegin(logkey_matgetsymtranspose,A,0,0,0);
39: /* Allocate space for symbolic transpose info and work array */
40: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
41: PetscMalloc(ai[am]*sizeof(PetscInt),&atj);
42: PetscMalloc(an*sizeof(PetscInt),&atfill);
43: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
45: /* Walk through aj and count ## of non-zeros in each row of A^T. */
46: /* Note: offset by 1 for fast conversion into csr format. */
47: for (i=0;i<ai[am];i++) {
48: ati[aj[i]+1] += 1;
49: }
50: /* Form ati for csr format of A^T. */
51: for (i=0;i<an;i++) {
52: ati[i+1] += ati[i];
53: }
55: /* Copy ati into atfill so we have locations of the next free space in atj */
56: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
58: /* Walk through A row-wise and mark nonzero entries of A^T. */
59: for (i=0;i<am;i++) {
60: anzj = ai[i+1] - ai[i];
61: for (j=0;j<anzj;j++) {
62: atj[atfill[*aj]] = i;
63: atfill[*aj++] += 1;
64: }
65: }
67: /* Clean up temporary space and complete requests. */
68: PetscFree(atfill);
69: *Ati = ati;
70: *Atj = atj;
72: PetscLogEventEnd(logkey_matgetsymtranspose,A,0,0,0);
73: return(0);
74: }
75: /*
76: MatGetSymbolicTransposeReduced_SeqAIJ() - Get symbolic matrix structure of submatrix A[rstart:rend,:],
77: modified from MatGetSymbolicTranspose_SeqAIJ()
78: */
79: static PetscEvent logkey_matgetsymtransreduced = 0;
82: PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A,PetscInt rstart,PetscInt rend,PetscInt *Ati[],PetscInt *Atj[])
83: {
85: PetscInt i,j,anzj;
86: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data;
87: PetscInt an=A->N;
88: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
92: PetscLogInfo(A,"Getting Symbolic Transpose.\n");
94: /* Set up timers */
95: if (!logkey_matgetsymtransreduced) {
96: PetscLogEventRegister(&logkey_matgetsymtransreduced,"MatGetSymbolicTransposeReduced",MAT_COOKIE);
97: }
98: PetscLogEventBegin(logkey_matgetsymtransreduced,A,0,0,0);
100: /* Allocate space for symbolic transpose info and work array */
101: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
102: anzj = ai[rend] - ai[rstart];
103: PetscMalloc((anzj+1)*sizeof(PetscInt),&atj);
104: PetscMalloc((an+1)*sizeof(PetscInt),&atfill);
105: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
107: /* Walk through aj and count ## of non-zeros in each row of A^T. */
108: /* Note: offset by 1 for fast conversion into csr format. */
109: for (i=ai[rstart]; i<ai[rend]; i++) {
110: ati[aj[i]+1] += 1;
111: }
112: /* Form ati for csr format of A^T. */
113: for (i=0;i<an;i++) {
114: ati[i+1] += ati[i];
115: }
117: /* Copy ati into atfill so we have locations of the next free space in atj */
118: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
120: /* Walk through A row-wise and mark nonzero entries of A^T. */
121: aj = aj + ai[rstart];
122: for (i=rstart; i<rend; i++) {
123: anzj = ai[i+1] - ai[i];
124: for (j=0;j<anzj;j++) {
125: atj[atfill[*aj]] = i-rstart;
126: atfill[*aj++] += 1;
127: }
128: }
130: /* Clean up temporary space and complete requests. */
131: PetscFree(atfill);
132: *Ati = ati;
133: *Atj = atj;
135: PetscLogEventEnd(logkey_matgetsymtransreduced,A,0,0,0);
136: return(0);
137: }
141: PetscErrorCode MatTranspose_SeqAIJ_FAST(Mat A,Mat *B)
142: {
144: PetscInt i,j,anzj;
145: Mat At;
146: Mat_SeqAIJ *a=(Mat_SeqAIJ *)A->data,*at;
147: PetscInt an=A->N,am=A->M;
148: PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
149: MatScalar *ata,*aa=a->a;
153: /* Set up timers */
154: if (!logkey_mattranspose) {
155: PetscLogEventRegister(&logkey_mattranspose,"MatTranspose_SeqAIJ_FAST",MAT_COOKIE);
156: }
157: PetscLogEventBegin(logkey_mattranspose,A,0,0,0);
159: /* Allocate space for symbolic transpose info and work array */
160: PetscMalloc((an+1)*sizeof(PetscInt),&ati);
161: PetscMalloc(ai[am]*sizeof(PetscInt),&atj);
162: PetscMalloc(ai[am]*sizeof(MatScalar),&ata);
163: PetscMalloc(an*sizeof(PetscInt),&atfill);
164: PetscMemzero(ati,(an+1)*sizeof(PetscInt));
165: /* Walk through aj and count ## of non-zeros in each row of A^T. */
166: /* Note: offset by 1 for fast conversion into csr format. */
167: for (i=0;i<ai[am];i++) {
168: ati[aj[i]+1] += 1;
169: }
170: /* Form ati for csr format of A^T. */
171: for (i=0;i<an;i++) {
172: ati[i+1] += ati[i];
173: }
175: /* Copy ati into atfill so we have locations of the next free space in atj */
176: PetscMemcpy(atfill,ati,an*sizeof(PetscInt));
178: /* Walk through A row-wise and mark nonzero entries of A^T. */
179: for (i=0;i<am;i++) {
180: anzj = ai[i+1] - ai[i];
181: for (j=0;j<anzj;j++) {
182: atj[atfill[*aj]] = i;
183: ata[atfill[*aj]] = *aa++;
184: atfill[*aj++] += 1;
185: }
186: }
188: /* Clean up temporary space and complete requests. */
189: PetscFree(atfill);
190: MatCreateSeqAIJWithArrays(A->comm,an,am,ati,atj,ata,&At);
191: at = (Mat_SeqAIJ *)(At->data);
192: at->freedata = PETSC_TRUE;
193: at->nonew = 0;
194: if (B) {
195: *B = At;
196: } else {
197: MatHeaderCopy(A,At);
198: }
199: PetscLogEventEnd(logkey_mattranspose,A,0,0,0);
200: return(0);
201: }
205: PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A,PetscInt *ati[],PetscInt *atj[])
206: {
210: PetscLogInfo(A,"Restoring Symbolic Transpose.\n");
211: PetscFree(*ati);
212: ati = PETSC_NULL;
213: PetscFree(*atj);
214: atj = PETSC_NULL;
215: return(0);
216: }