xref: /petsc/src/mat/impls/aij/seq/symtranspose.c (revision cf37664f2dbf3effe674ec5ac2718028c1621ae0)
1be1d678aSKris Buschelman 
270f19b1fSKris Buschelman /*
370f19b1fSKris Buschelman   Defines symbolic transpose routines for SeqAIJ matrices.
470f19b1fSKris Buschelman 
570f19b1fSKris Buschelman   Currently Get/Restore only allocates/frees memory for holding the
670f19b1fSKris Buschelman   (i,j) info for the transpose.  Someday, this info could be
770f19b1fSKris Buschelman   maintained so successive calls to Get will not recompute the info.
870f19b1fSKris Buschelman 
970f19b1fSKris Buschelman   Also defined is a "faster" implementation of MatTranspose for SeqAIJ
1070f19b1fSKris Buschelman   matrices which avoids calls to MatSetValues.  This routine has not
1170f19b1fSKris Buschelman   been adopted as the standard yet as it is somewhat untested.
1270f19b1fSKris Buschelman 
1370f19b1fSKris Buschelman */
1470f19b1fSKris Buschelman 
15c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
1670f19b1fSKris Buschelman 
1770f19b1fSKris Buschelman 
1870f19b1fSKris Buschelman #undef __FUNCT__
1953c77d0aSJed Brown #define __FUNCT__ "MatGetSymbolicTranspose_SeqAIJ"
202e111b49SBarry Smith PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat A,PetscInt *Ati[],PetscInt *Atj[])
21dfbe8321SBarry Smith {
22dfbe8321SBarry Smith   PetscErrorCode ierr;
232e111b49SBarry Smith   PetscInt       i,j,anzj;
2470f19b1fSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
25d0f46423SBarry Smith   PetscInt       an=A->cmap->N,am=A->rmap->N;
262e111b49SBarry Smith   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
2770f19b1fSKris Buschelman 
2870f19b1fSKris Buschelman   PetscFunctionBegin;
29ae15b995SBarry Smith   ierr = PetscInfo(A,"Getting Symbolic Transpose.\n");CHKERRQ(ierr);
3070f19b1fSKris Buschelman 
3170f19b1fSKris Buschelman   /* Set up timers */
324ebed01fSBarry Smith   ierr = PetscLogEventBegin(MAT_Getsymtranspose,A,0,0,0);CHKERRQ(ierr);
3370f19b1fSKris Buschelman 
3470f19b1fSKris Buschelman   /* Allocate space for symbolic transpose info and work array */
35854ce69bSBarry Smith   ierr = PetscCalloc1(an+1,&ati);CHKERRQ(ierr);
36785e854fSJed Brown   ierr = PetscMalloc1(ai[am],&atj);CHKERRQ(ierr);
37785e854fSJed Brown   ierr = PetscMalloc1(an,&atfill);CHKERRQ(ierr);
3870f19b1fSKris Buschelman 
3970f19b1fSKris Buschelman   /* Walk through aj and count ## of non-zeros in each row of A^T. */
4070f19b1fSKris Buschelman   /* Note: offset by 1 for fast conversion into csr format. */
4170f19b1fSKris Buschelman   for (i=0;i<ai[am];i++) {
4270f19b1fSKris Buschelman     ati[aj[i]+1] += 1;
4370f19b1fSKris Buschelman   }
4470f19b1fSKris Buschelman   /* Form ati for csr format of A^T. */
4570f19b1fSKris Buschelman   for (i=0;i<an;i++) {
4670f19b1fSKris Buschelman     ati[i+1] += ati[i];
4770f19b1fSKris Buschelman   }
4870f19b1fSKris Buschelman 
4970f19b1fSKris Buschelman   /* Copy ati into atfill so we have locations of the next free space in atj */
502e111b49SBarry Smith   ierr = PetscMemcpy(atfill,ati,an*sizeof(PetscInt));CHKERRQ(ierr);
5170f19b1fSKris Buschelman 
5270f19b1fSKris Buschelman   /* Walk through A row-wise and mark nonzero entries of A^T. */
5370f19b1fSKris Buschelman   for (i=0; i<am; i++) {
5470f19b1fSKris Buschelman     anzj = ai[i+1] - ai[i];
5570f19b1fSKris Buschelman     for (j=0; j<anzj; j++) {
5670f19b1fSKris Buschelman       atj[atfill[*aj]] = i;
5770f19b1fSKris Buschelman       atfill[*aj++]   += 1;
5870f19b1fSKris Buschelman     }
5970f19b1fSKris Buschelman   }
6070f19b1fSKris Buschelman 
6170f19b1fSKris Buschelman   /* Clean up temporary space and complete requests. */
6270f19b1fSKris Buschelman   ierr = PetscFree(atfill);CHKERRQ(ierr);
6370f19b1fSKris Buschelman   *Ati = ati;
6470f19b1fSKris Buschelman   *Atj = atj;
6570f19b1fSKris Buschelman 
664ebed01fSBarry Smith   ierr = PetscLogEventEnd(MAT_Getsymtranspose,A,0,0,0);CHKERRQ(ierr);
6770f19b1fSKris Buschelman   PetscFunctionReturn(0);
6870f19b1fSKris Buschelman }
690390132cSHong Zhang /*
700390132cSHong Zhang   MatGetSymbolicTransposeReduced_SeqAIJ() - Get symbolic matrix structure of submatrix A[rstart:rend,:],
710390132cSHong Zhang      modified from MatGetSymbolicTranspose_SeqAIJ()
720390132cSHong Zhang */
730390132cSHong Zhang #undef __FUNCT__
7453c77d0aSJed Brown #define __FUNCT__ "MatGetSymbolicTransposeReduced_SeqAIJ"
750390132cSHong Zhang PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A,PetscInt rstart,PetscInt rend,PetscInt *Ati[],PetscInt *Atj[])
760390132cSHong Zhang {
770390132cSHong Zhang   PetscErrorCode ierr;
780390132cSHong Zhang   PetscInt       i,j,anzj;
790390132cSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
80d0f46423SBarry Smith   PetscInt       an=A->cmap->N;
810390132cSHong Zhang   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
820390132cSHong Zhang 
830390132cSHong Zhang   PetscFunctionBegin;
84ae15b995SBarry Smith   ierr = PetscInfo(A,"Getting Symbolic Transpose\n");CHKERRQ(ierr);
854ebed01fSBarry Smith   ierr = PetscLogEventBegin(MAT_Getsymtransreduced,A,0,0,0);CHKERRQ(ierr);
860390132cSHong Zhang 
870390132cSHong Zhang   /* Allocate space for symbolic transpose info and work array */
88854ce69bSBarry Smith   ierr = PetscCalloc1(an+1,&ati);CHKERRQ(ierr);
890390132cSHong Zhang   anzj = ai[rend] - ai[rstart];
90854ce69bSBarry Smith   ierr = PetscMalloc1(anzj+1,&atj);CHKERRQ(ierr);
91854ce69bSBarry Smith   ierr = PetscMalloc1(an+1,&atfill);CHKERRQ(ierr);
920390132cSHong Zhang 
930390132cSHong Zhang   /* Walk through aj and count ## of non-zeros in each row of A^T. */
940390132cSHong Zhang   /* Note: offset by 1 for fast conversion into csr format. */
950390132cSHong Zhang   for (i=ai[rstart]; i<ai[rend]; i++) {
960390132cSHong Zhang     ati[aj[i]+1] += 1;
970390132cSHong Zhang   }
980390132cSHong Zhang   /* Form ati for csr format of A^T. */
990390132cSHong Zhang   for (i=0;i<an;i++) {
1000390132cSHong Zhang     ati[i+1] += ati[i];
1010390132cSHong Zhang   }
1020390132cSHong Zhang 
1030390132cSHong Zhang   /* Copy ati into atfill so we have locations of the next free space in atj */
1040390132cSHong Zhang   ierr = PetscMemcpy(atfill,ati,an*sizeof(PetscInt));CHKERRQ(ierr);
1050390132cSHong Zhang 
1060390132cSHong Zhang   /* Walk through A row-wise and mark nonzero entries of A^T. */
1070390132cSHong Zhang   aj = aj + ai[rstart];
1080390132cSHong Zhang   for (i=rstart; i<rend; i++) {
1090390132cSHong Zhang     anzj = ai[i+1] - ai[i];
1100390132cSHong Zhang     for (j=0; j<anzj; j++) {
1110390132cSHong Zhang       atj[atfill[*aj]] = i-rstart;
1120390132cSHong Zhang       atfill[*aj++]   += 1;
1130390132cSHong Zhang     }
1140390132cSHong Zhang   }
1150390132cSHong Zhang 
1160390132cSHong Zhang   /* Clean up temporary space and complete requests. */
1170390132cSHong Zhang   ierr = PetscFree(atfill);CHKERRQ(ierr);
1180390132cSHong Zhang   *Ati = ati;
1190390132cSHong Zhang   *Atj = atj;
1200390132cSHong Zhang 
1214ebed01fSBarry Smith   ierr = PetscLogEventEnd(MAT_Getsymtransreduced,A,0,0,0);CHKERRQ(ierr);
1220390132cSHong Zhang   PetscFunctionReturn(0);
1230390132cSHong Zhang }
12470f19b1fSKris Buschelman 
12570f19b1fSKris Buschelman #undef __FUNCT__
12653c77d0aSJed Brown #define __FUNCT__ "MatTranspose_SeqAIJ_FAST"
127fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqAIJ_FAST(Mat A,MatReuse reuse,Mat *B)
128dfbe8321SBarry Smith {
129dfbe8321SBarry Smith   PetscErrorCode ierr;
1302e111b49SBarry Smith   PetscInt       i,j,anzj;
13170f19b1fSKris Buschelman   Mat            At;
13270f19b1fSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*at;
133d0f46423SBarry Smith   PetscInt       an=A->cmap->N,am=A->rmap->N;
1342e111b49SBarry Smith   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;
13570f19b1fSKris Buschelman   MatScalar      *ata,*aa=a->a;
1362e111b49SBarry Smith 
13770f19b1fSKris Buschelman   PetscFunctionBegin;
1384ebed01fSBarry Smith   ierr = PetscLogEventBegin(MAT_Transpose_SeqAIJ,A,0,0,0);CHKERRQ(ierr);
13970f19b1fSKris Buschelman 
140*cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
14170f19b1fSKris Buschelman     /* Allocate space for symbolic transpose info and work array */
142854ce69bSBarry Smith     ierr = PetscCalloc1(an+1,&ati);CHKERRQ(ierr);
143785e854fSJed Brown     ierr = PetscMalloc1(ai[am],&atj);CHKERRQ(ierr);
144785e854fSJed Brown     ierr = PetscMalloc1(ai[am],&ata);CHKERRQ(ierr);
14570f19b1fSKris Buschelman     /* Walk through aj and count ## of non-zeros in each row of A^T. */
14670f19b1fSKris Buschelman     /* Note: offset by 1 for fast conversion into csr format. */
14770f19b1fSKris Buschelman     for (i=0;i<ai[am];i++) {
14870f19b1fSKris Buschelman       ati[aj[i]+1] += 1;
14970f19b1fSKris Buschelman     }
15070f19b1fSKris Buschelman     /* Form ati for csr format of A^T. */
15170f19b1fSKris Buschelman     for (i=0;i<an;i++) {
15270f19b1fSKris Buschelman       ati[i+1] += ati[i];
15370f19b1fSKris Buschelman     }
154fc4dec0aSBarry Smith   } else {
155fc4dec0aSBarry Smith     Mat_SeqAIJ *sub_B = (Mat_SeqAIJ*) (*B)->data;
156fc4dec0aSBarry Smith     ati = sub_B->i;
157fc4dec0aSBarry Smith     atj = sub_B->j;
158fc4dec0aSBarry Smith     ata = sub_B->a;
159fc4dec0aSBarry Smith     At  = *B;
160fc4dec0aSBarry Smith   }
16170f19b1fSKris Buschelman 
16270f19b1fSKris Buschelman   /* Copy ati into atfill so we have locations of the next free space in atj */
163785e854fSJed Brown   ierr = PetscMalloc1(an,&atfill);CHKERRQ(ierr);
1642e111b49SBarry Smith   ierr = PetscMemcpy(atfill,ati,an*sizeof(PetscInt));CHKERRQ(ierr);
16570f19b1fSKris Buschelman 
16670f19b1fSKris Buschelman   /* Walk through A row-wise and mark nonzero entries of A^T. */
16770f19b1fSKris Buschelman   for (i=0;i<am;i++) {
16870f19b1fSKris Buschelman     anzj = ai[i+1] - ai[i];
16970f19b1fSKris Buschelman     for (j=0;j<anzj;j++) {
17070f19b1fSKris Buschelman       atj[atfill[*aj]] = i;
17170f19b1fSKris Buschelman       ata[atfill[*aj]] = *aa++;
17270f19b1fSKris Buschelman       atfill[*aj++]   += 1;
17370f19b1fSKris Buschelman     }
17470f19b1fSKris Buschelman   }
17570f19b1fSKris Buschelman 
17670f19b1fSKris Buschelman   /* Clean up temporary space and complete requests. */
17770f19b1fSKris Buschelman   ierr = PetscFree(atfill);CHKERRQ(ierr);
178fc4dec0aSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
179ce94432eSBarry Smith     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),an,am,ati,atj,ata,&At);CHKERRQ(ierr);
1802205254eSKarl Rupp 
18170f19b1fSKris Buschelman     at          = (Mat_SeqAIJ*)(At->data);
182e6b907acSBarry Smith     at->free_a  = PETSC_TRUE;
183e6b907acSBarry Smith     at->free_ij = PETSC_TRUE;
18470f19b1fSKris Buschelman     at->nonew   = 0;
185fc4dec0aSBarry Smith   }
186fc4dec0aSBarry Smith 
187*cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
18870f19b1fSKris Buschelman     *B = At;
18970f19b1fSKris Buschelman   } else {
19028be2f97SBarry Smith     ierr = MatHeaderMerge(A,&At);CHKERRQ(ierr);
19170f19b1fSKris Buschelman   }
1924ebed01fSBarry Smith   ierr = PetscLogEventEnd(MAT_Transpose_SeqAIJ,A,0,0,0);CHKERRQ(ierr);
19370f19b1fSKris Buschelman   PetscFunctionReturn(0);
19470f19b1fSKris Buschelman }
19570f19b1fSKris Buschelman 
19670f19b1fSKris Buschelman #undef __FUNCT__
19770f19b1fSKris Buschelman #define __FUNCT__ "MatRestoreSymbolicTranspose_SeqAIJ"
1982e111b49SBarry Smith PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A,PetscInt *ati[],PetscInt *atj[])
1992e111b49SBarry Smith {
200dfbe8321SBarry Smith   PetscErrorCode ierr;
20170f19b1fSKris Buschelman 
20270f19b1fSKris Buschelman   PetscFunctionBegin;
203ae15b995SBarry Smith   ierr = PetscInfo(A,"Restoring Symbolic Transpose.\n");CHKERRQ(ierr);
20470f19b1fSKris Buschelman   ierr = PetscFree(*ati);CHKERRQ(ierr);
20570f19b1fSKris Buschelman   ierr = PetscFree(*atj);CHKERRQ(ierr);
20670f19b1fSKris Buschelman   PetscFunctionReturn(0);
20770f19b1fSKris Buschelman }
20870f19b1fSKris Buschelman 
209