xref: /petsc/src/mat/impls/aij/seq/symtranspose.c (revision fc4dec0ab267e831b41a38fdd32730b7f8abdc58)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
370f19b1fSKris Buschelman /*
470f19b1fSKris Buschelman   Defines symbolic transpose routines for SeqAIJ matrices.
570f19b1fSKris Buschelman 
670f19b1fSKris Buschelman   Currently Get/Restore only allocates/frees memory for holding the
770f19b1fSKris Buschelman   (i,j) info for the transpose.  Someday, this info could be
870f19b1fSKris Buschelman   maintained so successive calls to Get will not recompute the info.
970f19b1fSKris Buschelman 
1070f19b1fSKris Buschelman   Also defined is a "faster" implementation of MatTranspose for SeqAIJ
1170f19b1fSKris Buschelman   matrices which avoids calls to MatSetValues.  This routine has not
1270f19b1fSKris Buschelman   been adopted as the standard yet as it is somewhat untested.
1370f19b1fSKris Buschelman 
1470f19b1fSKris Buschelman */
1570f19b1fSKris Buschelman 
1670f19b1fSKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
1770f19b1fSKris Buschelman 
1870f19b1fSKris Buschelman 
1970f19b1fSKris Buschelman #undef __FUNCT__
2070f19b1fSKris Buschelman #define __FUNCT__ "MatGetSymbolicTranspose_SeqIJ"
212e111b49SBarry Smith PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat A,PetscInt *Ati[],PetscInt *Atj[])
22dfbe8321SBarry Smith {
23dfbe8321SBarry Smith   PetscErrorCode ierr;
242e111b49SBarry Smith   PetscInt       i,j,anzj;
2570f19b1fSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ *)A->data;
26899cda47SBarry Smith   PetscInt       an=A->cmap.N,am=A->rmap.N;
272e111b49SBarry Smith   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
2870f19b1fSKris Buschelman 
2970f19b1fSKris Buschelman   PetscFunctionBegin;
3070f19b1fSKris Buschelman 
31ae15b995SBarry Smith   ierr = PetscInfo(A,"Getting Symbolic Transpose.\n");CHKERRQ(ierr);
3270f19b1fSKris Buschelman 
3370f19b1fSKris Buschelman   /* Set up timers */
344ebed01fSBarry Smith   ierr = PetscLogEventBegin(MAT_Getsymtranspose,A,0,0,0);CHKERRQ(ierr);
3570f19b1fSKris Buschelman 
3670f19b1fSKris Buschelman   /* Allocate space for symbolic transpose info and work array */
372e111b49SBarry Smith   ierr = PetscMalloc((an+1)*sizeof(PetscInt),&ati);CHKERRQ(ierr);
382e111b49SBarry Smith   ierr = PetscMalloc(ai[am]*sizeof(PetscInt),&atj);CHKERRQ(ierr);
392e111b49SBarry Smith   ierr = PetscMalloc(an*sizeof(PetscInt),&atfill);CHKERRQ(ierr);
402e111b49SBarry Smith   ierr = PetscMemzero(ati,(an+1)*sizeof(PetscInt));CHKERRQ(ierr);
4170f19b1fSKris Buschelman 
4270f19b1fSKris Buschelman   /* Walk through aj and count ## of non-zeros in each row of A^T. */
4370f19b1fSKris Buschelman   /* Note: offset by 1 for fast conversion into csr format. */
4470f19b1fSKris Buschelman   for (i=0;i<ai[am];i++) {
4570f19b1fSKris Buschelman     ati[aj[i]+1] += 1;
4670f19b1fSKris Buschelman   }
4770f19b1fSKris Buschelman   /* Form ati for csr format of A^T. */
4870f19b1fSKris Buschelman   for (i=0;i<an;i++) {
4970f19b1fSKris Buschelman     ati[i+1] += ati[i];
5070f19b1fSKris Buschelman   }
5170f19b1fSKris Buschelman 
5270f19b1fSKris Buschelman   /* Copy ati into atfill so we have locations of the next free space in atj */
532e111b49SBarry Smith   ierr = PetscMemcpy(atfill,ati,an*sizeof(PetscInt));CHKERRQ(ierr);
5470f19b1fSKris Buschelman 
5570f19b1fSKris Buschelman   /* Walk through A row-wise and mark nonzero entries of A^T. */
5670f19b1fSKris Buschelman   for (i=0;i<am;i++) {
5770f19b1fSKris Buschelman     anzj = ai[i+1] - ai[i];
5870f19b1fSKris Buschelman     for (j=0;j<anzj;j++) {
5970f19b1fSKris Buschelman       atj[atfill[*aj]] = i;
6070f19b1fSKris Buschelman       atfill[*aj++]   += 1;
6170f19b1fSKris Buschelman     }
6270f19b1fSKris Buschelman   }
6370f19b1fSKris Buschelman 
6470f19b1fSKris Buschelman   /* Clean up temporary space and complete requests. */
6570f19b1fSKris Buschelman   ierr = PetscFree(atfill);CHKERRQ(ierr);
6670f19b1fSKris Buschelman   *Ati = ati;
6770f19b1fSKris Buschelman   *Atj = atj;
6870f19b1fSKris Buschelman 
694ebed01fSBarry Smith   ierr = PetscLogEventEnd(MAT_Getsymtranspose,A,0,0,0);CHKERRQ(ierr);
7070f19b1fSKris Buschelman   PetscFunctionReturn(0);
7170f19b1fSKris Buschelman }
720390132cSHong Zhang /*
730390132cSHong Zhang   MatGetSymbolicTransposeReduced_SeqAIJ() - Get symbolic matrix structure of submatrix A[rstart:rend,:],
740390132cSHong Zhang      modified from MatGetSymbolicTranspose_SeqAIJ()
750390132cSHong Zhang */
760390132cSHong Zhang #undef __FUNCT__
770390132cSHong Zhang #define __FUNCT__ "MatGetSymbolicTransposeReduced_SeqIJ"
780390132cSHong Zhang PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A,PetscInt rstart,PetscInt rend,PetscInt *Ati[],PetscInt *Atj[])
790390132cSHong Zhang {
800390132cSHong Zhang   PetscErrorCode ierr;
810390132cSHong Zhang   PetscInt       i,j,anzj;
820390132cSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ *)A->data;
83899cda47SBarry Smith   PetscInt       an=A->cmap.N;
840390132cSHong Zhang   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
850390132cSHong Zhang 
860390132cSHong Zhang   PetscFunctionBegin;
87ae15b995SBarry Smith   ierr = PetscInfo(A,"Getting Symbolic Transpose\n");CHKERRQ(ierr);
884ebed01fSBarry Smith   ierr = PetscLogEventBegin(MAT_Getsymtransreduced,A,0,0,0);CHKERRQ(ierr);
890390132cSHong Zhang 
900390132cSHong Zhang   /* Allocate space for symbolic transpose info and work array */
910390132cSHong Zhang   ierr = PetscMalloc((an+1)*sizeof(PetscInt),&ati);CHKERRQ(ierr);
920390132cSHong Zhang   anzj = ai[rend] - ai[rstart];
930390132cSHong Zhang   ierr = PetscMalloc((anzj+1)*sizeof(PetscInt),&atj);CHKERRQ(ierr);
940390132cSHong Zhang   ierr = PetscMalloc((an+1)*sizeof(PetscInt),&atfill);CHKERRQ(ierr);
950390132cSHong Zhang   ierr = PetscMemzero(ati,(an+1)*sizeof(PetscInt));CHKERRQ(ierr);
960390132cSHong Zhang 
970390132cSHong Zhang   /* Walk through aj and count ## of non-zeros in each row of A^T. */
980390132cSHong Zhang   /* Note: offset by 1 for fast conversion into csr format. */
990390132cSHong Zhang   for (i=ai[rstart]; i<ai[rend]; i++) {
1000390132cSHong Zhang     ati[aj[i]+1] += 1;
1010390132cSHong Zhang   }
1020390132cSHong Zhang   /* Form ati for csr format of A^T. */
1030390132cSHong Zhang   for (i=0;i<an;i++) {
1040390132cSHong Zhang     ati[i+1] += ati[i];
1050390132cSHong Zhang   }
1060390132cSHong Zhang 
1070390132cSHong Zhang   /* Copy ati into atfill so we have locations of the next free space in atj */
1080390132cSHong Zhang   ierr = PetscMemcpy(atfill,ati,an*sizeof(PetscInt));CHKERRQ(ierr);
1090390132cSHong Zhang 
1100390132cSHong Zhang   /* Walk through A row-wise and mark nonzero entries of A^T. */
1110390132cSHong Zhang   aj = aj + ai[rstart];
1120390132cSHong Zhang   for (i=rstart; i<rend; i++) {
1130390132cSHong Zhang     anzj = ai[i+1] - ai[i];
1140390132cSHong Zhang     for (j=0;j<anzj;j++) {
1150390132cSHong Zhang       atj[atfill[*aj]] = i-rstart;
1160390132cSHong Zhang       atfill[*aj++]   += 1;
1170390132cSHong Zhang     }
1180390132cSHong Zhang   }
1190390132cSHong Zhang 
1200390132cSHong Zhang   /* Clean up temporary space and complete requests. */
1210390132cSHong Zhang   ierr = PetscFree(atfill);CHKERRQ(ierr);
1220390132cSHong Zhang   *Ati = ati;
1230390132cSHong Zhang   *Atj = atj;
1240390132cSHong Zhang 
1254ebed01fSBarry Smith   ierr = PetscLogEventEnd(MAT_Getsymtransreduced,A,0,0,0);CHKERRQ(ierr);
1260390132cSHong Zhang   PetscFunctionReturn(0);
1270390132cSHong Zhang }
12870f19b1fSKris Buschelman 
12970f19b1fSKris Buschelman #undef __FUNCT__
13070f19b1fSKris Buschelman #define __FUNCT__ "MatTranspose_SeqIJ_FAST"
131*fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqAIJ_FAST(Mat A,MatReuse reuse,Mat *B)
132dfbe8321SBarry Smith {
133dfbe8321SBarry Smith   PetscErrorCode ierr;
1342e111b49SBarry Smith   PetscInt       i,j,anzj;
13570f19b1fSKris Buschelman   Mat            At;
13670f19b1fSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ *)A->data,*at;
137899cda47SBarry Smith   PetscInt       an=A->cmap.N,am=A->rmap.N;
1382e111b49SBarry Smith   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
13970f19b1fSKris Buschelman   MatScalar      *ata,*aa=a->a;
1402e111b49SBarry Smith 
14170f19b1fSKris Buschelman   PetscFunctionBegin;
14270f19b1fSKris Buschelman 
1434ebed01fSBarry Smith   ierr = PetscLogEventBegin(MAT_Transpose_SeqAIJ,A,0,0,0);CHKERRQ(ierr);
14470f19b1fSKris Buschelman 
145*fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || *B == A) {
14670f19b1fSKris Buschelman     /* Allocate space for symbolic transpose info and work array */
1472e111b49SBarry Smith     ierr = PetscMalloc((an+1)*sizeof(PetscInt),&ati);CHKERRQ(ierr);
1482e111b49SBarry Smith     ierr = PetscMalloc(ai[am]*sizeof(PetscInt),&atj);CHKERRQ(ierr);
14970f19b1fSKris Buschelman     ierr = PetscMalloc(ai[am]*sizeof(MatScalar),&ata);CHKERRQ(ierr);
1502e111b49SBarry Smith     ierr = PetscMemzero(ati,(an+1)*sizeof(PetscInt));CHKERRQ(ierr);
15170f19b1fSKris Buschelman     /* Walk through aj and count ## of non-zeros in each row of A^T. */
15270f19b1fSKris Buschelman     /* Note: offset by 1 for fast conversion into csr format. */
15370f19b1fSKris Buschelman     for (i=0;i<ai[am];i++) {
15470f19b1fSKris Buschelman       ati[aj[i]+1] += 1;
15570f19b1fSKris Buschelman     }
15670f19b1fSKris Buschelman     /* Form ati for csr format of A^T. */
15770f19b1fSKris Buschelman     for (i=0;i<an;i++) {
15870f19b1fSKris Buschelman       ati[i+1] += ati[i];
15970f19b1fSKris Buschelman     }
160*fc4dec0aSBarry Smith   } else {
161*fc4dec0aSBarry Smith     Mat_SeqAIJ *sub_B = (Mat_SeqAIJ*) (*B)->data;
162*fc4dec0aSBarry Smith     ati = sub_B->i;
163*fc4dec0aSBarry Smith     atj = sub_B->j;
164*fc4dec0aSBarry Smith     ata = sub_B->a;
165*fc4dec0aSBarry Smith     At = *B;
166*fc4dec0aSBarry Smith   }
16770f19b1fSKris Buschelman 
16870f19b1fSKris Buschelman   /* Copy ati into atfill so we have locations of the next free space in atj */
169*fc4dec0aSBarry Smith   ierr = PetscMalloc(an*sizeof(PetscInt),&atfill);CHKERRQ(ierr);
1702e111b49SBarry Smith   ierr = PetscMemcpy(atfill,ati,an*sizeof(PetscInt));CHKERRQ(ierr);
17170f19b1fSKris Buschelman 
17270f19b1fSKris Buschelman   /* Walk through A row-wise and mark nonzero entries of A^T. */
17370f19b1fSKris Buschelman   for (i=0;i<am;i++) {
17470f19b1fSKris Buschelman     anzj = ai[i+1] - ai[i];
17570f19b1fSKris Buschelman     for (j=0;j<anzj;j++) {
17670f19b1fSKris Buschelman       atj[atfill[*aj]] = i;
17770f19b1fSKris Buschelman       ata[atfill[*aj]] = *aa++;
17870f19b1fSKris Buschelman       atfill[*aj++]   += 1;
17970f19b1fSKris Buschelman     }
18070f19b1fSKris Buschelman   }
18170f19b1fSKris Buschelman 
18270f19b1fSKris Buschelman   /* Clean up temporary space and complete requests. */
18370f19b1fSKris Buschelman   ierr = PetscFree(atfill);CHKERRQ(ierr);
184*fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
1857adad957SLisandro Dalcin     ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,an,am,ati,atj,ata,&At);CHKERRQ(ierr);
18670f19b1fSKris Buschelman     at   = (Mat_SeqAIJ *)(At->data);
187e6b907acSBarry Smith     at->free_a  = PETSC_TRUE;
188e6b907acSBarry Smith     at->free_ij  = PETSC_TRUE;
18970f19b1fSKris Buschelman     at->nonew   = 0;
190*fc4dec0aSBarry Smith   }
191*fc4dec0aSBarry Smith 
192*fc4dec0aSBarry Smith   if (*B != A) {
19370f19b1fSKris Buschelman     *B = At;
19470f19b1fSKris Buschelman   } else {
19570f19b1fSKris Buschelman     ierr = MatHeaderCopy(A,At);
19670f19b1fSKris Buschelman   }
1974ebed01fSBarry Smith   ierr = PetscLogEventEnd(MAT_Transpose_SeqAIJ,A,0,0,0);CHKERRQ(ierr);
19870f19b1fSKris Buschelman   PetscFunctionReturn(0);
19970f19b1fSKris Buschelman }
20070f19b1fSKris Buschelman 
20170f19b1fSKris Buschelman #undef __FUNCT__
20270f19b1fSKris Buschelman #define __FUNCT__ "MatRestoreSymbolicTranspose_SeqAIJ"
2032e111b49SBarry Smith PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A,PetscInt *ati[],PetscInt *atj[])
2042e111b49SBarry Smith {
205dfbe8321SBarry Smith   PetscErrorCode ierr;
20670f19b1fSKris Buschelman 
20770f19b1fSKris Buschelman   PetscFunctionBegin;
208ae15b995SBarry Smith   ierr = PetscInfo(A,"Restoring Symbolic Transpose.\n");CHKERRQ(ierr);
20970f19b1fSKris Buschelman   ierr = PetscFree(*ati);CHKERRQ(ierr);
21070f19b1fSKris Buschelman   ati  = PETSC_NULL;
21170f19b1fSKris Buschelman   ierr = PetscFree(*atj);CHKERRQ(ierr);
21270f19b1fSKris Buschelman   atj  = PETSC_NULL;
21370f19b1fSKris Buschelman   PetscFunctionReturn(0);
21470f19b1fSKris Buschelman }
21570f19b1fSKris Buschelman 
216