1be1d678aSKris Buschelman 270f19b1fSKris Buschelman /* 37fb60732SBarry Smith Defines transpose routines for SeqAIJ matrices. 470f19b1fSKris Buschelman */ 570f19b1fSKris Buschelman 6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 770f19b1fSKris Buschelman 87fb60732SBarry Smith /* 97fb60732SBarry Smith The symbolic and full transpose versions share several similar code blocks but the macros to reuse the code would be confusing and ugly 107fb60732SBarry Smith */ 11d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat A, Mat *B) 12d71ae5a4SJacob Faibussowitsch { 132e111b49SBarry Smith PetscInt i, j, anzj; 147fb60732SBarry Smith Mat At; 157fb60732SBarry 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. */ 267fb60732SBarry Smith for (i = 0; i < ai[am]; i++) ati[aj[i] + 1] += 1; 2770f19b1fSKris Buschelman /* Form ati for csr format of A^T. */ 287fb60732SBarry 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 */ 317fb60732SBarry 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 447fb60732SBarry Smith PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A), an, am, ati, atj, NULL, &At)); 457fb60732SBarry Smith PetscCall(MatSetBlockSizes(At, PetscAbs(A->cmap->bs), PetscAbs(A->rmap->bs))); 467fb60732SBarry Smith PetscCall(MatSetType(At, ((PetscObject)A)->type_name)); 477fb60732SBarry Smith at = (Mat_SeqAIJ *)At->data; 487fb60732SBarry Smith at->free_a = PETSC_FALSE; 497fb60732SBarry Smith at->free_ij = PETSC_TRUE; 507fb60732SBarry Smith at->nonew = 0; 517fb60732SBarry Smith at->maxnz = ati[an]; 527fb60732SBarry Smith *B = At; 53*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 540390132cSHong Zhang } 5570f19b1fSKris Buschelman 56d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTranspose_SeqAIJ(Mat A, MatReuse reuse, Mat *B) 57d71ae5a4SJacob Faibussowitsch { 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; 657fb60732SBarry Smith PetscContainer rB; 667fb60732SBarry Smith MatParentState *rb; 677fb60732SBarry Smith PetscBool nonzerochange = PETSC_FALSE; 682e111b49SBarry Smith 6970f19b1fSKris Buschelman PetscFunctionBegin; 707fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 717fb60732SBarry Smith PetscCall(PetscObjectQuery((PetscObject)*B, "MatTransposeParent", (PetscObject *)&rB)); 727fb60732SBarry Smith PetscCheck(rB, PetscObjectComm((PetscObject)*B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from call to MatTranspose()"); 737fb60732SBarry Smith PetscCall(PetscContainerGetPointer(rB, (void **)&rb)); 747fb60732SBarry Smith if (rb->nonzerostate != A->nonzerostate) nonzerochange = PETSC_TRUE; 757fb60732SBarry Smith } 767fb60732SBarry Smith 779566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(A, &av)); 78ce496241SStefano Zampini aa = av; 797fb60732SBarry 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. */ 857fb60732SBarry Smith for (i = 0; i < ai[am]; i++) ati[aj[i] + 1] += 1; 8670f19b1fSKris Buschelman /* Form ati for csr format of A^T. */ 877fb60732SBarry Smith for (i = 0; i < an; i++) ati[i + 1] += ati[i]; 887fb60732SBarry Smith PetscCall(PetscMalloc1(ai[am], &ata)); 897fb60732SBarry 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)); 1117fb60732SBarry Smith PetscCall(MatSeqAIJRestoreArrayRead(A, &av)); 1127fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(PetscObjectStateIncrease((PetscObject)(*B))); 1137fb60732SBarry Smith 1147fb60732SBarry 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))); 1177fb60732SBarry 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 1257fb60732SBarry Smith if (reuse == MAT_INITIAL_MATRIX || (reuse == MAT_REUSE_MATRIX && !nonzerochange)) { 12670f19b1fSKris Buschelman *B = At; 1277fb60732SBarry Smith } else if (nonzerochange) { 1287fb60732SBarry Smith PetscCall(MatHeaderMerge(*B, &At)); 1297fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(A, *B)); 1307fb60732SBarry Smith } else if (reuse == MAT_INPLACE_MATRIX) { 1319566063dSJacob Faibussowitsch PetscCall(MatHeaderMerge(A, &At)); 13270f19b1fSKris Buschelman } 133*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13470f19b1fSKris Buschelman } 13570f19b1fSKris Buschelman 1367fb60732SBarry Smith /* 1377fb60732SBarry Smith Get symbolic matrix structure of a submatrix of A, A[rstart:rend,:], 1387fb60732SBarry Smith */ 139d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat A, PetscInt rstart, PetscInt rend, PetscInt *Ati[], PetscInt *Atj[]) 140d71ae5a4SJacob Faibussowitsch { 1417fb60732SBarry Smith PetscInt i, j, anzj; 1427fb60732SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 1437fb60732SBarry Smith PetscInt an = A->cmap->N; 1447fb60732SBarry Smith PetscInt *ati, *atj, *atfill, *ai = a->i, *aj = a->j, am = ai[rend] - ai[rstart]; 1457fb60732SBarry Smith 1467fb60732SBarry Smith PetscFunctionBegin; 1477fb60732SBarry Smith PetscCall(PetscLogEventBegin(MAT_Getsymtransreduced, A, 0, 0, 0)); 1487fb60732SBarry Smith 1497fb60732SBarry Smith /* Allocate space for symbolic transpose info and work array */ 1507fb60732SBarry Smith PetscCall(PetscCalloc1(an + 1, &ati)); 1517fb60732SBarry Smith PetscCall(PetscMalloc1(am + 1, &atj)); 1527fb60732SBarry Smith 1537fb60732SBarry Smith /* Walk through aj and count ## of non-zeros in each row of A^T. */ 1547fb60732SBarry Smith /* Note: offset by 1 for fast conversion into csr format. */ 1557fb60732SBarry Smith for (i = ai[rstart]; i < ai[rend]; i++) ati[aj[i] + 1] += 1; 1567fb60732SBarry Smith /* Form ati for csr format of A^T. */ 1577fb60732SBarry Smith for (i = 0; i < an; i++) ati[i + 1] += ati[i]; 1587fb60732SBarry Smith 1597fb60732SBarry Smith /* Copy ati into atfill so we have locations of the next free space in atj */ 1607fb60732SBarry Smith PetscCall(PetscMalloc1(an + 1, &atfill)); 1617fb60732SBarry Smith PetscCall(PetscArraycpy(atfill, ati, an)); 1627fb60732SBarry Smith 1637fb60732SBarry Smith /* Walk through A row-wise and mark nonzero entries of A^T. */ 1647fb60732SBarry Smith aj = aj + ai[rstart]; 1657fb60732SBarry Smith for (i = rstart; i < rend; i++) { 1667fb60732SBarry Smith anzj = ai[i + 1] - ai[i]; 1677fb60732SBarry Smith for (j = 0; j < anzj; j++) { 1687fb60732SBarry Smith atj[atfill[*aj]] = i - rstart; 1697fb60732SBarry Smith atfill[*aj++] += 1; 1707fb60732SBarry Smith } 1717fb60732SBarry Smith } 1727fb60732SBarry Smith PetscCall(PetscFree(atfill)); 1737fb60732SBarry Smith *Ati = ati; 1747fb60732SBarry Smith *Atj = atj; 1757fb60732SBarry Smith 1767fb60732SBarry Smith PetscCall(PetscLogEventEnd(MAT_Getsymtransreduced, A, 0, 0, 0)); 177*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1787fb60732SBarry Smith } 1797fb60732SBarry Smith 1807fb60732SBarry Smith /* 1817fb60732SBarry Smith Returns the i and j arrays for a symbolic transpose, this is used internally within SeqAIJ code when the full 1827fb60732SBarry Smith symbolic matrix (which can be obtained with MatTransposeSymbolic() is not needed. MatRestoreSymbolicTranspose_SeqAIJ() should be used to free the arrays. 1837fb60732SBarry Smith */ 184d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat A, PetscInt *Ati[], PetscInt *Atj[]) 185d71ae5a4SJacob Faibussowitsch { 1867fb60732SBarry Smith PetscFunctionBegin; 1877fb60732SBarry Smith PetscCall(MatGetSymbolicTransposeReduced_SeqAIJ(A, 0, A->rmap->N, Ati, Atj)); 188*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1897fb60732SBarry Smith } 1907fb60732SBarry Smith 191d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat A, PetscInt *ati[], PetscInt *atj[]) 192d71ae5a4SJacob Faibussowitsch { 19370f19b1fSKris Buschelman PetscFunctionBegin; 1949566063dSJacob Faibussowitsch PetscCall(PetscFree(*ati)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFree(*atj)); 196*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19770f19b1fSKris Buschelman } 198