xref: /petsc/src/mat/impls/baij/seq/aijbaij.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
1c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
261e22dc2SBarry Smith 
MatConvert_SeqBAIJ_SeqAIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)3d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
4d71ae5a4SJacob Faibussowitsch {
5676c34cdSKris Buschelman   Mat          B;
60f6d62edSLisandro Dalcin   Mat_SeqAIJ  *b;
70f6d62edSLisandro Dalcin   PetscBool    roworiented;
861e22dc2SBarry Smith   Mat_SeqBAIJ *a  = (Mat_SeqBAIJ *)A->data;
9d0f46423SBarry Smith   PetscInt     bs = A->rmap->bs, *ai = a->i, *aj = a->j, n = A->rmap->N / bs, i, j, k;
10690b6cddSBarry Smith   PetscInt    *rowlengths, *rows, *cols, maxlen            = 0, ncols;
11dd6ea824SBarry Smith   MatScalar   *aa = a->a;
1261e22dc2SBarry Smith 
1361e22dc2SBarry Smith   PetscFunctionBegin;
140f6d62edSLisandro Dalcin   if (reuse == MAT_REUSE_MATRIX) {
150f6d62edSLisandro Dalcin     B = *newmat;
1657508eceSPierre Jolivet     for (i = 0; i < n; i++) maxlen = PetscMax(maxlen, ai[i + 1] - ai[i]);
170f6d62edSLisandro Dalcin   } else {
189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n * bs, &rowlengths));
1961e22dc2SBarry Smith     for (i = 0; i < n; i++) {
2057508eceSPierre Jolivet       maxlen = PetscMax(maxlen, ai[i + 1] - ai[i]);
21ad540459SPierre Jolivet       for (j = 0; j < bs; j++) rowlengths[i * bs + j] = bs * (ai[i + 1] - ai[i]);
2261e22dc2SBarry Smith     }
239566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
249566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQAIJ));
259566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
269566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
279566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, rowlengths));
289566063dSJacob Faibussowitsch     PetscCall(PetscFree(rowlengths));
290f6d62edSLisandro Dalcin   }
300f6d62edSLisandro Dalcin   b           = (Mat_SeqAIJ *)B->data;
310f6d62edSLisandro Dalcin   roworiented = b->roworiented;
3261e22dc2SBarry Smith 
339566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, PETSC_FALSE));
349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &rows));
359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * maxlen, &cols));
36b9b97703SBarry Smith   for (i = 0; i < n; i++) {
37ad540459SPierre Jolivet     for (j = 0; j < bs; j++) rows[j] = i * bs + j;
38b9b97703SBarry Smith     ncols = ai[i + 1] - ai[i];
39b9b97703SBarry Smith     for (k = 0; k < ncols; k++) {
40ad540459SPierre Jolivet       for (j = 0; j < bs; j++) cols[k * bs + j] = bs * (*aj) + j;
41b9b97703SBarry Smith       aj++;
42b9b97703SBarry Smith     }
439566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, bs, rows, bs * ncols, cols, aa, INSERT_VALUES));
448e3a54c0SPierre Jolivet     aa = PetscSafePointerPlusOffset(aa, ncols * bs * bs);
45b9b97703SBarry Smith   }
469566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
479566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_ROW_ORIENTED, roworiented));
51521d7252SBarry Smith 
52*ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
53*ac530a7eSPierre Jolivet   else *newmat = B;
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5561e22dc2SBarry Smith }
5661e22dc2SBarry Smith 
57c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
5885fc7724SBarry Smith 
MatConvert_SeqAIJ_SeqBAIJ_Preallocate(Mat A,PetscInt ** nnz)59247fa90bSBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqBAIJ_Preallocate(Mat A, PetscInt **nnz)
60247fa90bSBarry Smith {
61247fa90bSBarry Smith   Mat_SeqAIJ     *Aa = (Mat_SeqAIJ *)A->data;
6258b7e2c1SStefano Zampini   PetscInt        m, n, bs = A->rmap->bs;
63f85a0629SBarry Smith   const PetscInt *ai = Aa->i, *aj = Aa->j;
64247fa90bSBarry Smith 
65247fa90bSBarry Smith   PetscFunctionBegin;
66247fa90bSBarry Smith   PetscCall(MatGetSize(A, &m, &n));
67247fa90bSBarry Smith 
68f85a0629SBarry Smith   if (bs == 1) {
69f85a0629SBarry Smith     PetscCall(PetscMalloc1(m, nnz));
70f85a0629SBarry Smith     for (PetscInt i = 0; i < m; i++) (*nnz)[i] = ai[i + 1] - ai[i];
71f85a0629SBarry Smith   } else {
72f85a0629SBarry Smith     PetscHSetIJ    ht;
73f85a0629SBarry Smith     PetscInt       k;
74f85a0629SBarry Smith     PetscHashIJKey key;
75f85a0629SBarry Smith     PetscBool      missing;
76f85a0629SBarry Smith 
77f85a0629SBarry Smith     PetscCall(PetscHSetIJCreate(&ht));
78f85a0629SBarry Smith     PetscCall(PetscCalloc1(m / bs, nnz));
79f85a0629SBarry Smith     for (PetscInt i = 0; i < m; i++) {
80f85a0629SBarry Smith       key.i = i / bs;
81f85a0629SBarry Smith       for (k = ai[i]; k < ai[i + 1]; k++) {
82f85a0629SBarry Smith         key.j = aj[k] / bs;
83f85a0629SBarry Smith         PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
84f85a0629SBarry Smith         if (missing) (*nnz)[key.i]++;
85247fa90bSBarry Smith       }
86247fa90bSBarry Smith     }
87f85a0629SBarry Smith     PetscCall(PetscHSetIJDestroy(&ht));
88247fa90bSBarry Smith   }
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90247fa90bSBarry Smith }
91247fa90bSBarry Smith 
MatConvert_SeqAIJ_SeqBAIJ(Mat A,MatType newtype,MatReuse reuse,Mat * newmat)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;
9758b7e2c1SStefano Zampini   PetscInt     m = A->rmap->N, n = A->cmap->N, *rowlengths, bs = A->rmap->bs;
9885fc7724SBarry Smith 
9985fc7724SBarry Smith   PetscFunctionBegin;
100345dab1aSStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
101247fa90bSBarry 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) {
110f4f49eeaSPierre Jolivet     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 
127*ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &B));
128*ac530a7eSPierre Jolivet   else *newmat = B;
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13085fc7724SBarry Smith }
131