1 2 /* 3 Defines symbolic transpose routines for SeqAIJ matrices. 4 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. 8 9 Also defined is a faster implementation of MatTranspose for SeqAIJ 10 matrices which avoids calls to MatSetValues. This routine is the new 11 standard since it is much faster than MatTranspose_AIJ. 12 13 */ 14 15 #include <../src/mat/impls/aij/seq/aij.h> 16 17 PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat A,PetscInt *Ati[],PetscInt *Atj[]) 18 { 19 PetscInt i,j,anzj; 20 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 21 PetscInt an=A->cmap->N,am=A->rmap->N; 22 PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j; 23 24 PetscFunctionBegin; 25 PetscCall(PetscInfo(A,"Getting Symbolic Transpose.\n")); 26 27 /* Set up timers */ 28 PetscCall(PetscLogEventBegin(MAT_Getsymtranspose,A,0,0,0)); 29 30 /* Allocate space for symbolic transpose info and work array */ 31 PetscCall(PetscCalloc1(an+1,&ati)); 32 PetscCall(PetscMalloc1(ai[am],&atj)); 33 PetscCall(PetscMalloc1(an,&atfill)); 34 35 /* Walk through aj and count ## of non-zeros in each row of A^T. */ 36 /* Note: offset by 1 for fast conversion into csr format. */ 37 for (i=0;i<ai[am];i++) { 38 ati[aj[i]+1] += 1; 39 } 40 /* Form ati for csr format of A^T. */ 41 for (i=0;i<an;i++) { 42 ati[i+1] += ati[i]; 43 } 44 45 /* Copy ati into atfill so we have locations of the next free space in atj */ 46 PetscCall(PetscArraycpy(atfill,ati,an)); 47 48 /* Walk through A row-wise and mark nonzero entries of A^T. */ 49 for (i=0; i<am; i++) { 50 anzj = ai[i+1] - ai[i]; 51 for (j=0; j<anzj; j++) { 52 atj[atfill[*aj]] = i; 53 atfill[*aj++] += 1; 54 } 55 } 56 57 /* Clean up temporary space and complete requests. */ 58 PetscCall(PetscFree(atfill)); 59 *Ati = ati; 60 *Atj = atj; 61 62 PetscCall(PetscLogEventEnd(MAT_Getsymtranspose,A,0,0,0)); 63 PetscFunctionReturn(0); 64 } 65 /* 66 MatGetSymbolicTransposeReduced_SeqAIJ() - Get symbolic matrix structure of submatrix A[rstart:rend,:], 67 modified from MatGetSymbolicTranspose_SeqAIJ() 68 */ 69 PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A,PetscInt rstart,PetscInt rend,PetscInt *Ati[],PetscInt *Atj[]) 70 { 71 PetscInt i,j,anzj; 72 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 73 PetscInt an=A->cmap->N; 74 PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j; 75 76 PetscFunctionBegin; 77 PetscCall(PetscInfo(A,"Getting Symbolic Transpose\n")); 78 PetscCall(PetscLogEventBegin(MAT_Getsymtransreduced,A,0,0,0)); 79 80 /* Allocate space for symbolic transpose info and work array */ 81 PetscCall(PetscCalloc1(an+1,&ati)); 82 anzj = ai[rend] - ai[rstart]; 83 PetscCall(PetscMalloc1(anzj+1,&atj)); 84 PetscCall(PetscMalloc1(an+1,&atfill)); 85 86 /* Walk through aj and count ## of non-zeros in each row of A^T. */ 87 /* Note: offset by 1 for fast conversion into csr format. */ 88 for (i=ai[rstart]; i<ai[rend]; i++) { 89 ati[aj[i]+1] += 1; 90 } 91 /* Form ati for csr format of A^T. */ 92 for (i=0;i<an;i++) { 93 ati[i+1] += ati[i]; 94 } 95 96 /* Copy ati into atfill so we have locations of the next free space in atj */ 97 PetscCall(PetscArraycpy(atfill,ati,an)); 98 99 /* Walk through A row-wise and mark nonzero entries of A^T. */ 100 aj = aj + ai[rstart]; 101 for (i=rstart; i<rend; i++) { 102 anzj = ai[i+1] - ai[i]; 103 for (j=0; j<anzj; j++) { 104 atj[atfill[*aj]] = i-rstart; 105 atfill[*aj++] += 1; 106 } 107 } 108 109 /* Clean up temporary space and complete requests. */ 110 PetscCall(PetscFree(atfill)); 111 *Ati = ati; 112 *Atj = atj; 113 114 PetscCall(PetscLogEventEnd(MAT_Getsymtransreduced,A,0,0,0)); 115 PetscFunctionReturn(0); 116 } 117 118 PetscErrorCode MatTranspose_SeqAIJ(Mat A,MatReuse reuse,Mat *B) 119 { 120 PetscInt i,j,anzj; 121 Mat At; 122 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*at; 123 PetscInt an=A->cmap->N,am=A->rmap->N; 124 PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j; 125 MatScalar *ata; 126 const MatScalar *aa,*av; 127 128 PetscFunctionBegin; 129 PetscCall(MatSeqAIJGetArrayRead(A,&av)); 130 aa = av; 131 if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 132 /* Allocate space for symbolic transpose info and work array */ 133 PetscCall(PetscCalloc1(an+1,&ati)); 134 PetscCall(PetscMalloc1(ai[am],&atj)); 135 PetscCall(PetscMalloc1(ai[am],&ata)); 136 /* Walk through aj and count ## of non-zeros in each row of A^T. */ 137 /* Note: offset by 1 for fast conversion into csr format. */ 138 for (i=0;i<ai[am];i++) { 139 ati[aj[i]+1] += 1; /* count ## of non-zeros for row aj[i] of A^T */ 140 } 141 /* Form ati for csr format of A^T. */ 142 for (i=0;i<an;i++) { 143 ati[i+1] += ati[i]; 144 } 145 } else { /* This segment is called by MatTranspose_MPIAIJ(...,MAT_INITIAL_MATRIX,..) directly! */ 146 Mat_SeqAIJ *sub_B = (Mat_SeqAIJ*) (*B)->data; 147 ati = sub_B->i; 148 atj = sub_B->j; 149 ata = sub_B->a; 150 At = *B; 151 } 152 153 /* Copy ati into atfill so we have locations of the next free space in atj */ 154 PetscCall(PetscMalloc1(an,&atfill)); 155 PetscCall(PetscArraycpy(atfill,ati,an)); 156 157 /* Walk through A row-wise and mark nonzero entries of A^T. */ 158 for (i=0;i<am;i++) { 159 anzj = ai[i+1] - ai[i]; 160 for (j=0;j<anzj;j++) { 161 atj[atfill[*aj]] = i; 162 ata[atfill[*aj]] = *aa++; 163 atfill[*aj++] += 1; 164 } 165 } 166 PetscCall(MatSeqAIJRestoreArrayRead(A,&av)); 167 168 /* Clean up temporary space and complete requests. */ 169 PetscCall(PetscFree(atfill)); 170 if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) { 171 PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),an,am,ati,atj,ata,&At)); 172 PetscCall(MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs))); 173 174 at = (Mat_SeqAIJ*)(At->data); 175 at->free_a = PETSC_TRUE; 176 at->free_ij = PETSC_TRUE; 177 at->nonew = 0; 178 at->maxnz = ati[an]; 179 180 PetscCall(MatSetType(At,((PetscObject)A)->type_name)); 181 } 182 183 if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) { 184 *B = At; 185 } else { 186 PetscCall(MatHeaderMerge(A,&At)); 187 } 188 PetscFunctionReturn(0); 189 } 190 191 PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A,PetscInt *ati[],PetscInt *atj[]) 192 { 193 PetscFunctionBegin; 194 PetscCall(PetscInfo(A,"Restoring Symbolic Transpose.\n")); 195 PetscCall(PetscFree(*ati)); 196 PetscCall(PetscFree(*atj)); 197 PetscFunctionReturn(0); 198 } 199