xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 247fa90bf5e10ed3512c8e8f64aafd98857b2846)
161e22dc2SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
361e22dc2SBarry Smith 
4d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
5d71ae5a4SJacob Faibussowitsch {
6676c34cdSKris Buschelman   Mat          B;
70f6d62edSLisandro Dalcin   Mat_SeqAIJ  *b;
80f6d62edSLisandro Dalcin   PetscBool    roworiented;
961e22dc2SBarry Smith   Mat_SeqBAIJ *a  = (Mat_SeqBAIJ *)A->data;
10d0f46423SBarry Smith   PetscInt     bs = A->rmap->bs, *ai = a->i, *aj = a->j, n = A->rmap->N / bs, i, j, k;
11690b6cddSBarry Smith   PetscInt    *rowlengths, *rows, *cols, maxlen            = 0, ncols;
12dd6ea824SBarry Smith   MatScalar   *aa = a->a;
1361e22dc2SBarry Smith 
1461e22dc2SBarry Smith   PetscFunctionBegin;
150f6d62edSLisandro Dalcin   if (reuse == MAT_REUSE_MATRIX) {
160f6d62edSLisandro Dalcin     B = *newmat;
17ad540459SPierre Jolivet     for (i = 0; i < n; i++) maxlen = PetscMax(maxlen, (ai[i + 1] - ai[i]));
180f6d62edSLisandro Dalcin   } else {
199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n * bs, &rowlengths));
2061e22dc2SBarry Smith     for (i = 0; i < n; i++) {
21b9b97703SBarry Smith       maxlen = PetscMax(maxlen, (ai[i + 1] - ai[i]));
22ad540459SPierre Jolivet       for (j = 0; j < bs; j++) rowlengths[i * bs + j] = bs * (ai[i + 1] - ai[i]);
2361e22dc2SBarry Smith     }
249566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
259566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQAIJ));
269566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
279566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
289566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, rowlengths));
299566063dSJacob Faibussowitsch     PetscCall(PetscFree(rowlengths));
300f6d62edSLisandro Dalcin   }
310f6d62edSLisandro Dalcin   b           = (Mat_SeqAIJ *)B->data;
320f6d62edSLisandro Dalcin   roworiented = b->roworiented;
3361e22dc2SBarry Smith 
349566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, PETSC_FALSE));
359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &rows));
369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * maxlen, &cols));
37b9b97703SBarry Smith   for (i = 0; i < n; i++) {
38ad540459SPierre Jolivet     for (j = 0; j < bs; j++) rows[j] = i * bs + j;
39b9b97703SBarry Smith     ncols = ai[i + 1] - ai[i];
40b9b97703SBarry Smith     for (k = 0; k < ncols; k++) {
41ad540459SPierre Jolivet       for (j = 0; j < bs; j++) cols[k * bs + j] = bs * (*aj) + j;
42b9b97703SBarry Smith       aj++;
43b9b97703SBarry Smith     }
449566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, bs, rows, bs * ncols, cols, aa, INSERT_VALUES));
45b9b97703SBarry Smith     aa += ncols * bs * bs;
46b9b97703SBarry Smith   }
479566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
489566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
519566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, roworiented));
52521d7252SBarry Smith 
53511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
549566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
55ae8d29abSPierre Jolivet   } else *newmat = B;
5661e22dc2SBarry Smith   PetscFunctionReturn(0);
5761e22dc2SBarry Smith }
5861e22dc2SBarry Smith 
59c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
6085fc7724SBarry Smith 
61*247fa90bSBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ_Preallocate(Mat A, PetscInt **nnz)
62*247fa90bSBarry Smith {
63*247fa90bSBarry Smith   PetscInt    m, n, bs = PetscAbs(A->rmap->bs);
64*247fa90bSBarry Smith   PetscInt   *ii;
65*247fa90bSBarry Smith   Mat_SeqAIJ *Aa = (Mat_SeqAIJ *)A->data;
66*247fa90bSBarry Smith 
67*247fa90bSBarry Smith   PetscFunctionBegin;
68*247fa90bSBarry Smith   PetscCall(MatGetSize(A, &m, &n));
69*247fa90bSBarry Smith   PetscCall(PetscCalloc1(m / bs, nnz));
70*247fa90bSBarry Smith   PetscCall(PetscMalloc1(bs, &ii));
71*247fa90bSBarry Smith 
72*247fa90bSBarry Smith   /* loop over all potential blocks in the matrix to see if the potential block has a nonzero */
73*247fa90bSBarry Smith   for (PetscInt i = 0; i < m / bs; i++) {
74*247fa90bSBarry Smith     for (PetscInt k = 0; k < bs; k++) ii[k] = Aa->i[i * bs + k];
75*247fa90bSBarry Smith     for (PetscInt j = 0; j < n / bs; j++) {
76*247fa90bSBarry Smith       for (PetscInt k = 0; k < bs; k++) {
77*247fa90bSBarry Smith         for (; ii[k] < Aa->i[i * bs + k + 1] && Aa->j[ii[k]] < (j + 1) * bs; ii[k]++) {
78*247fa90bSBarry Smith           if (Aa->j[ii[k]] >= j * bs) {
79*247fa90bSBarry Smith             /* found a nonzero in the potential block so allocate for that block and jump to check the next potential block */
80*247fa90bSBarry Smith             (*nnz)[i]++;
81*247fa90bSBarry Smith             goto theend;
82*247fa90bSBarry Smith           }
83*247fa90bSBarry Smith         }
84*247fa90bSBarry Smith       }
85*247fa90bSBarry Smith     theend:;
86*247fa90bSBarry Smith     }
87*247fa90bSBarry Smith   }
88*247fa90bSBarry Smith   PetscCall(PetscFree(ii));
89*247fa90bSBarry Smith   PetscFunctionReturn(0);
90*247fa90bSBarry Smith }
91*247fa90bSBarry Smith 
92d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
93d71ae5a4SJacob Faibussowitsch {
94676c34cdSKris Buschelman   Mat          B;
9585fc7724SBarry Smith   Mat_SeqAIJ  *a = (Mat_SeqAIJ *)A->data;
9685fc7724SBarry Smith   Mat_SeqBAIJ *b;
97*247fa90bSBarry Smith   PetscInt     m = A->rmap->N, n = A->cmap->N, *rowlengths, bs = PetscAbs(A->rmap->bs);
9885fc7724SBarry Smith 
9985fc7724SBarry Smith   PetscFunctionBegin;
100345dab1aSStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
101*247fa90bSBarry Smith     PetscCall(MatConvert_SeqAIJ_SeqBAIJ_Preallocate(A, &rowlengths));
1029566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1039566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
1049566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQBAIJ));
1059566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetPreallocation(B, bs, 0, rowlengths));
1069566063dSJacob Faibussowitsch     PetscCall(PetscFree(rowlengths));
107345dab1aSStefano Zampini   } else B = *newmat;
10885fc7724SBarry Smith 
109ae8d29abSPierre Jolivet   if (bs == 1) {
110676c34cdSKris Buschelman     b = (Mat_SeqBAIJ *)(B->data);
11185fc7724SBarry Smith 
1129566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->i, a->i, m + 1));
1139566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->ilen, a->ilen, m));
1149566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->j, a->j, a->nz));
1159566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->a, a->a, a->nz));
11685fc7724SBarry Smith 
1179566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, PETSC_TRUE));
1189566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1199566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
120ae8d29abSPierre Jolivet   } else {
121ae8d29abSPierre Jolivet     /* reuse may not be equal to MAT_REUSE_MATRIX, but the basic converter will reallocate or replace newmat if this value is not used */
122ae8d29abSPierre Jolivet     /* if reuse is equal to MAT_INITIAL_MATRIX, it has been appropriately preallocated before                                          */
123ae8d29abSPierre Jolivet     /*                      MAT_INPLACE_MATRIX, it will be replaced with MatHeaderReplace below                                        */
1249566063dSJacob Faibussowitsch     PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B));
125ae8d29abSPierre Jolivet   }
126676c34cdSKris Buschelman 
127511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
1289566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
129ae8d29abSPierre Jolivet   } else *newmat = B;
13085fc7724SBarry Smith   PetscFunctionReturn(0);
13185fc7724SBarry Smith }
132