1be1d678aSKris Buschelman 270f19b1fSKris Buschelman /* 3*7fb60732SBarry Smith Defines transpose routines for SeqAIJ matrices. 470f19b1fSKris Buschelman */ 570f19b1fSKris Buschelman 6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 770f19b1fSKris Buschelman 8*7fb60732SBarry Smith /* 9*7fb60732SBarry Smith The symbolic and full transpose versions share several similar code blocks but the macros to reuse the code would be confusing and ugly 10*7fb60732SBarry Smith */ 11*7fb60732SBarry Smith PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat A,Mat *B) 12dfbe8321SBarry Smith { 132e111b49SBarry Smith PetscInt i,j,anzj; 14*7fb60732SBarry Smith Mat At; 15*7fb60732SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*at; 16d0f46423SBarry Smith PetscInt an=A->cmap->N,am=A->rmap->N; 172e111b49SBarry Smith PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j; 1870f19b1fSKris Buschelman 1970f19b1fSKris Buschelman PetscFunctionBegin; 2070f19b1fSKris Buschelman /* Allocate space for symbolic transpose info and work array */ 219566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(an+1,&ati)); 229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ai[am],&atj)); 2370f19b1fSKris Buschelman 2470f19b1fSKris Buschelman /* Walk through aj and count ## of non-zeros in each row of A^T. */ 2570f19b1fSKris Buschelman /* Note: offset by 1 for fast conversion into csr format. */ 26*7fb60732SBarry Smith for (i=0;i<ai[am];i++) ati[aj[i]+1] += 1; 2770f19b1fSKris Buschelman /* Form ati for csr format of A^T. */ 28*7fb60732SBarry Smith for (i=0;i<an;i++) ati[i+1] += ati[i]; 2970f19b1fSKris Buschelman 3070f19b1fSKris Buschelman /* Copy ati into atfill so we have locations of the next free space in atj */ 31*7fb60732SBarry Smith PetscCall(PetscMalloc1(an,&atfill)); 329566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(atfill,ati,an)); 3370f19b1fSKris Buschelman 3470f19b1fSKris Buschelman /* Walk through A row-wise and mark nonzero entries of A^T. */ 3570f19b1fSKris Buschelman for (i=0;i<am;i++) { 3670f19b1fSKris Buschelman anzj = ai[i+1] - ai[i]; 3770f19b1fSKris Buschelman for (j=0;j<anzj;j++) { 3870f19b1fSKris Buschelman atj[atfill[*aj]] = i; 3970f19b1fSKris Buschelman atfill[*aj++] += 1; 4070f19b1fSKris Buschelman } 4170f19b1fSKris Buschelman } 429566063dSJacob Faibussowitsch PetscCall(PetscFree(atfill)); 4370f19b1fSKris Buschelman 44*7fb60732SBarry Smith PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),an,am,ati,atj,NULL,&At)); 45*7fb60732SBarry Smith PetscCall(MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs))); 46*7fb60732SBarry Smith PetscCall(MatSetType(At,((PetscObject)A)->type_name)); 47*7fb60732SBarry Smith at = (Mat_SeqAIJ*)At->data; 48*7fb60732SBarry Smith at->free_a = PETSC_FALSE; 49*7fb60732SBarry Smith at->free_ij = PETSC_TRUE; 50*7fb60732SBarry Smith at->nonew = 0; 51*7fb60732SBarry Smith at->maxnz = ati[an]; 52*7fb60732SBarry Smith *B = At; 530390132cSHong Zhang PetscFunctionReturn(0); 540390132cSHong Zhang } 5570f19b1fSKris Buschelman 5691e9d3e2SHong Zhang PetscErrorCode MatTranspose_SeqAIJ(Mat A,MatReuse reuse,Mat *B) 57dfbe8321SBarry Smith { 582e111b49SBarry Smith PetscInt i,j,anzj; 5970f19b1fSKris Buschelman Mat At; 6070f19b1fSKris Buschelman Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*at; 61d0f46423SBarry Smith PetscInt an=A->cmap->N,am=A->rmap->N; 622e111b49SBarry Smith PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j; 63ce496241SStefano Zampini MatScalar *ata; 64ce496241SStefano Zampini const MatScalar *aa,*av; 65*7fb60732SBarry Smith PetscContainer rB; 66*7fb60732SBarry Smith MatParentState *rb; 67*7fb60732SBarry Smith PetscBool nonzerochange = PETSC_FALSE; 682e111b49SBarry Smith 6970f19b1fSKris Buschelman PetscFunctionBegin; 70*7fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 71*7fb60732SBarry Smith PetscCall(PetscObjectQuery((PetscObject)*B,"MatTransposeParent",(PetscObject*)&rB)); 72*7fb60732SBarry Smith PetscCheck(rB,PetscObjectComm((PetscObject)*B),PETSC_ERR_ARG_WRONG,"Reuse matrix used was not generated from call to MatTranspose()"); 73*7fb60732SBarry Smith PetscCall(PetscContainerGetPointer(rB,(void**)&rb)); 74*7fb60732SBarry Smith if (rb->nonzerostate != A->nonzerostate) nonzerochange = PETSC_TRUE; 75*7fb60732SBarry Smith } 76*7fb60732SBarry Smith 779566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(A,&av)); 78ce496241SStefano Zampini aa = av; 79*7fb60732SBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX || nonzerochange) { 8070f19b1fSKris Buschelman /* Allocate space for symbolic transpose info and work array */ 819566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(an+1,&ati)); 829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ai[am],&atj)); 8370f19b1fSKris Buschelman /* Walk through aj and count ## of non-zeros in each row of A^T. */ 8470f19b1fSKris Buschelman /* Note: offset by 1 for fast conversion into csr format. */ 85*7fb60732SBarry Smith for (i=0;i<ai[am];i++) ati[aj[i]+1] += 1; 8670f19b1fSKris Buschelman /* Form ati for csr format of A^T. */ 87*7fb60732SBarry Smith for (i=0;i<an;i++) ati[i+1] += ati[i]; 88*7fb60732SBarry Smith PetscCall(PetscMalloc1(ai[am],&ata)); 89*7fb60732SBarry Smith } else { 90fc4dec0aSBarry Smith Mat_SeqAIJ *sub_B = (Mat_SeqAIJ*) (*B)->data; 91fc4dec0aSBarry Smith ati = sub_B->i; 92fc4dec0aSBarry Smith atj = sub_B->j; 93fc4dec0aSBarry Smith ata = sub_B->a; 94fc4dec0aSBarry Smith At = *B; 95fc4dec0aSBarry Smith } 9670f19b1fSKris Buschelman 9770f19b1fSKris Buschelman /* Copy ati into atfill so we have locations of the next free space in atj */ 989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(an,&atfill)); 999566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(atfill,ati,an)); 10070f19b1fSKris Buschelman 10170f19b1fSKris Buschelman /* Walk through A row-wise and mark nonzero entries of A^T. */ 10270f19b1fSKris Buschelman for (i=0;i<am;i++) { 10370f19b1fSKris Buschelman anzj = ai[i+1] - ai[i]; 10470f19b1fSKris Buschelman for (j=0;j<anzj;j++) { 10570f19b1fSKris Buschelman atj[atfill[*aj]] = i; 10670f19b1fSKris Buschelman ata[atfill[*aj]] = *aa++; 10770f19b1fSKris Buschelman atfill[*aj++] += 1; 10870f19b1fSKris Buschelman } 10970f19b1fSKris Buschelman } 1109566063dSJacob Faibussowitsch PetscCall(PetscFree(atfill)); 111*7fb60732SBarry Smith PetscCall(MatSeqAIJRestoreArrayRead(A,&av)); 112*7fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(PetscObjectStateIncrease((PetscObject)(*B))); 113*7fb60732SBarry Smith 114*7fb60732SBarry Smith if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX || nonzerochange) { 1159566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),an,am,ati,atj,ata,&At)); 1169566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs))); 117*7fb60732SBarry Smith PetscCall(MatSetType(At,((PetscObject)A)->type_name)); 11870f19b1fSKris Buschelman at = (Mat_SeqAIJ*)(At->data); 119e6b907acSBarry Smith at->free_a = PETSC_TRUE; 120e6b907acSBarry Smith at->free_ij = PETSC_TRUE; 12170f19b1fSKris Buschelman at->nonew = 0; 12265ac4ee0Sandi selinger at->maxnz = ati[an]; 123fc4dec0aSBarry Smith } 124fc4dec0aSBarry Smith 125*7fb60732SBarry Smith if (reuse == MAT_INITIAL_MATRIX || (reuse == MAT_REUSE_MATRIX && !nonzerochange)) { 12670f19b1fSKris Buschelman *B = At; 127*7fb60732SBarry Smith } else if (nonzerochange) { 128*7fb60732SBarry Smith PetscCall(MatHeaderMerge(*B,&At)); 129*7fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(A,*B)); 130*7fb60732SBarry Smith } else if (reuse == MAT_INPLACE_MATRIX) { 1319566063dSJacob Faibussowitsch PetscCall(MatHeaderMerge(A,&At)); 13270f19b1fSKris Buschelman } 13370f19b1fSKris Buschelman PetscFunctionReturn(0); 13470f19b1fSKris Buschelman } 13570f19b1fSKris Buschelman 136*7fb60732SBarry Smith /* 137*7fb60732SBarry Smith Get symbolic matrix structure of a submatrix of A, A[rstart:rend,:], 138*7fb60732SBarry Smith */ 139*7fb60732SBarry Smith PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A,PetscInt rstart,PetscInt rend,PetscInt *Ati[],PetscInt *Atj[]) 140*7fb60732SBarry Smith { 141*7fb60732SBarry Smith PetscInt i,j,anzj; 142*7fb60732SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 143*7fb60732SBarry Smith PetscInt an=A->cmap->N; 144*7fb60732SBarry Smith PetscInt *ati,*atj,*atfill,*ai=a->i,*aj=a->j, am=ai[rend] - ai[rstart]; 145*7fb60732SBarry Smith 146*7fb60732SBarry Smith PetscFunctionBegin; 147*7fb60732SBarry Smith PetscCall(PetscLogEventBegin(MAT_Getsymtransreduced,A,0,0,0)); 148*7fb60732SBarry Smith 149*7fb60732SBarry Smith /* Allocate space for symbolic transpose info and work array */ 150*7fb60732SBarry Smith PetscCall(PetscCalloc1(an+1,&ati)); 151*7fb60732SBarry Smith PetscCall(PetscMalloc1(am+1,&atj)); 152*7fb60732SBarry Smith 153*7fb60732SBarry Smith /* Walk through aj and count ## of non-zeros in each row of A^T. */ 154*7fb60732SBarry Smith /* Note: offset by 1 for fast conversion into csr format. */ 155*7fb60732SBarry Smith for (i=ai[rstart]; i<ai[rend]; i++) ati[aj[i]+1] += 1; 156*7fb60732SBarry Smith /* Form ati for csr format of A^T. */ 157*7fb60732SBarry Smith for (i=0;i<an;i++) ati[i+1] += ati[i]; 158*7fb60732SBarry Smith 159*7fb60732SBarry Smith /* Copy ati into atfill so we have locations of the next free space in atj */ 160*7fb60732SBarry Smith PetscCall(PetscMalloc1(an+1,&atfill)); 161*7fb60732SBarry Smith PetscCall(PetscArraycpy(atfill,ati,an)); 162*7fb60732SBarry Smith 163*7fb60732SBarry Smith /* Walk through A row-wise and mark nonzero entries of A^T. */ 164*7fb60732SBarry Smith aj = aj + ai[rstart]; 165*7fb60732SBarry Smith for (i=rstart; i<rend; i++) { 166*7fb60732SBarry Smith anzj = ai[i+1] - ai[i]; 167*7fb60732SBarry Smith for (j=0; j<anzj; j++) { 168*7fb60732SBarry Smith atj[atfill[*aj]] = i-rstart; 169*7fb60732SBarry Smith atfill[*aj++] += 1; 170*7fb60732SBarry Smith } 171*7fb60732SBarry Smith } 172*7fb60732SBarry Smith PetscCall(PetscFree(atfill)); 173*7fb60732SBarry Smith *Ati = ati; 174*7fb60732SBarry Smith *Atj = atj; 175*7fb60732SBarry Smith 176*7fb60732SBarry Smith PetscCall(PetscLogEventEnd(MAT_Getsymtransreduced,A,0,0,0)); 177*7fb60732SBarry Smith PetscFunctionReturn(0); 178*7fb60732SBarry Smith } 179*7fb60732SBarry Smith 180*7fb60732SBarry Smith /* 181*7fb60732SBarry Smith Returns the i and j arrays for a symbolic transpose, this is used internally within SeqAIJ code when the full 182*7fb60732SBarry Smith symbolic matrix (which can be obtained with MatTransposeSymbolic() is not needed. MatRestoreSymbolicTranspose_SeqAIJ() should be used to free the arrays. 183*7fb60732SBarry Smith */ 184*7fb60732SBarry Smith PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat A,PetscInt *Ati[],PetscInt *Atj[]) 185*7fb60732SBarry Smith { 186*7fb60732SBarry Smith PetscFunctionBegin; 187*7fb60732SBarry Smith PetscCall(MatGetSymbolicTransposeReduced_SeqAIJ(A,0,A->rmap->N,Ati,Atj)); 188*7fb60732SBarry Smith PetscFunctionReturn(0); 189*7fb60732SBarry Smith } 190*7fb60732SBarry Smith 1912e111b49SBarry Smith PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A,PetscInt *ati[],PetscInt *atj[]) 1922e111b49SBarry Smith { 19370f19b1fSKris Buschelman PetscFunctionBegin; 1949566063dSJacob Faibussowitsch PetscCall(PetscFree(*ati)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFree(*atj)); 19670f19b1fSKris Buschelman PetscFunctionReturn(0); 19770f19b1fSKris Buschelman } 198*7fb60732SBarry Smith 199