1d4002b98SHong Zhang #include <../src/mat/impls/sell/seq/sell.h> 2d4002b98SHong Zhang #include <../src/mat/impls/aij/seq/aij.h> 3d4002b98SHong Zhang #include <petsc/private/isimpl.h> 4d4002b98SHong Zhang 5d4002b98SHong Zhang /* 6d4002b98SHong Zhang MatGetColumnIJ_SeqSELL_Color() and MatRestoreColumnIJ_SeqSELL_Color() are customized from 7d4002b98SHong Zhang MatGetColumnIJ_SeqSELL() and MatRestoreColumnIJ_SeqSELL() by adding an output 8d4002b98SHong Zhang spidx[], index of a->a, to be used in MatTransposeColoringCreate_SeqSELL() and MatFDColoringCreate_SeqSELL() 9d4002b98SHong Zhang */ 10d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetColumnIJ_SeqSELL_Color(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *nn, const PetscInt *ia[], const PetscInt *ja[], PetscInt *spidx[], PetscBool *done) 11d71ae5a4SJacob Faibussowitsch { 12d4002b98SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 13d4002b98SHong Zhang PetscInt i, j, *collengths, *cia, *cja, n = A->cmap->n, totalslices; 14d4002b98SHong Zhang PetscInt row, col; 15d4002b98SHong Zhang PetscInt *cspidx; 16d4002b98SHong Zhang PetscBool isnonzero; 17d4002b98SHong Zhang 18d4002b98SHong Zhang PetscFunctionBegin; 19*a748edf9SJed Brown PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected cmap->n %" PetscInt_FMT " >= 0", n); 20d4002b98SHong Zhang *nn = n; 213ba16761SJacob Faibussowitsch if (!ia) PetscFunctionReturn(PETSC_SUCCESS); 22d4002b98SHong Zhang 23*a748edf9SJed Brown PetscCall(PetscCalloc1(n, &collengths)); 249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &cia)); 259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->nz + 1, &cja)); 269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->nz + 1, &cspidx)); 27d4002b98SHong Zhang 28*a748edf9SJed Brown totalslices = PetscCeilInt(A->rmap->n, 8); 29d4002b98SHong Zhang for (i = 0; i < totalslices; i++) { /* loop over slices */ 30d4002b98SHong Zhang for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) & 0x07)) { 31d4002b98SHong Zhang isnonzero = (PetscBool)((j - a->sliidx[i]) / 8 < a->rlen[8 * i + row]); 32d4002b98SHong Zhang if (isnonzero) collengths[a->colidx[j]]++; 33d4002b98SHong Zhang } 34d4002b98SHong Zhang } 35d4002b98SHong Zhang 36d4002b98SHong Zhang cia[0] = oshift; 37ad540459SPierre Jolivet for (i = 0; i < n; i++) cia[i + 1] = cia[i] + collengths[i]; 389566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(collengths, n)); 39d4002b98SHong Zhang 40d4002b98SHong Zhang for (i = 0; i < totalslices; i++) { /* loop over slices */ 41d4002b98SHong Zhang for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) & 0x07)) { 42d4002b98SHong Zhang isnonzero = (PetscBool)((j - a->sliidx[i]) / 8 < a->rlen[8 * i + row]); 43d4002b98SHong Zhang if (isnonzero) { 44d4002b98SHong Zhang col = a->colidx[j]; 45d4002b98SHong Zhang cspidx[cia[col] + collengths[col] - oshift] = j; /* index of a->colidx */ 46d4002b98SHong Zhang cja[cia[col] + collengths[col] - oshift] = 8 * i + row + oshift; /* row index */ 47d4002b98SHong Zhang collengths[col]++; 48d4002b98SHong Zhang } 49d4002b98SHong Zhang } 50d4002b98SHong Zhang } 51d4002b98SHong Zhang 529566063dSJacob Faibussowitsch PetscCall(PetscFree(collengths)); 539371c9d4SSatish Balay *ia = cia; 549371c9d4SSatish Balay *ja = cja; 55d4002b98SHong Zhang *spidx = cspidx; 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57d4002b98SHong Zhang } 58d4002b98SHong Zhang 59d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreColumnIJ_SeqSELL_Color(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscInt *spidx[], PetscBool *done) 60d71ae5a4SJacob Faibussowitsch { 61d4002b98SHong Zhang PetscFunctionBegin; 623ba16761SJacob Faibussowitsch if (!ia) PetscFunctionReturn(PETSC_SUCCESS); 639566063dSJacob Faibussowitsch PetscCall(PetscFree(*ia)); 649566063dSJacob Faibussowitsch PetscCall(PetscFree(*ja)); 659566063dSJacob Faibussowitsch PetscCall(PetscFree(*spidx)); 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67d4002b98SHong Zhang } 68