xref: /petsc/src/mat/impls/aij/seq/seqhashmatsetvalues.h (revision 26cec32642a02ddeaa9481e1a5bd50b8500ffeea)
1*26cec326SBarry Smith /*
2*26cec326SBarry Smith    used by SEQAIJ, BAIJ and SBAIJ to reduce code duplication
3*26cec326SBarry Smith 
4*26cec326SBarry Smith      define TYPE to AIJ BAIJ or SBAIJ
5*26cec326SBarry Smith             TYPE_BS_ON for BAIJ and SBAIJ and bs > 1
6*26cec326SBarry Smith             TYPE_SBAIJ for SBAIJ
7*26cec326SBarry Smith 
8*26cec326SBarry Smith */
9*26cec326SBarry Smith static PetscErrorCode PetscConcat(MatSetValues_Seq_Hash, TYPE_BS)(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode addv)
10*26cec326SBarry Smith {
11*26cec326SBarry Smith   PetscConcat(Mat_Seq, TYPE) *a = (PetscConcat(Mat_Seq, TYPE) *)A->data;
12*26cec326SBarry Smith #if defined(TYPE_BS_ON)
13*26cec326SBarry Smith   PetscInt bs;
14*26cec326SBarry Smith #endif
15*26cec326SBarry Smith 
16*26cec326SBarry Smith   PetscFunctionBegin;
17*26cec326SBarry Smith #if defined(TYPE_BS_ON)
18*26cec326SBarry Smith   PetscCall(MatGetBlockSize(A, &bs));
19*26cec326SBarry Smith #endif
20*26cec326SBarry Smith   for (PetscInt r = 0; r < m; ++r) {
21*26cec326SBarry Smith     PetscHashIJKey key;
22*26cec326SBarry Smith     PetscBool      missing;
23*26cec326SBarry Smith     PetscScalar    value;
24*26cec326SBarry Smith #if defined(TYPE_BS_ON)
25*26cec326SBarry Smith     PetscHashIJKey bkey;
26*26cec326SBarry Smith #endif
27*26cec326SBarry Smith 
28*26cec326SBarry Smith     key.i = rows[r];
29*26cec326SBarry Smith #if defined(TYPE_BS_ON)
30*26cec326SBarry Smith     bkey.i = key.i / bs;
31*26cec326SBarry Smith #endif
32*26cec326SBarry Smith     if (key.i < 0) continue;
33*26cec326SBarry Smith     for (PetscInt c = 0; c < n; ++c) {
34*26cec326SBarry Smith       key.j = cols[c];
35*26cec326SBarry Smith #if defined(TYPE_BS_ON)
36*26cec326SBarry Smith       bkey.j = key.j / bs;
37*26cec326SBarry Smith   #if defined(TYPE_SBAIJ)
38*26cec326SBarry Smith       if (bkey.j < bkey.i) continue;
39*26cec326SBarry Smith   #else
40*26cec326SBarry Smith       if (key.j < 0) continue;
41*26cec326SBarry Smith   #endif
42*26cec326SBarry Smith #else
43*26cec326SBarry Smith   #if defined(TYPE_SBAIJ)
44*26cec326SBarry Smith       if (key.j < key.i) continue;
45*26cec326SBarry Smith   #else
46*26cec326SBarry Smith       if (key.j < 0) continue;
47*26cec326SBarry Smith   #endif
48*26cec326SBarry Smith #endif
49*26cec326SBarry Smith       value = values ? ((a->roworiented) ? values[r * n + c] : values[r + m * c]) : 0;
50*26cec326SBarry Smith       switch (addv) {
51*26cec326SBarry Smith       case INSERT_VALUES:
52*26cec326SBarry Smith         PetscCall(PetscHMapIJVQuerySet(a->ht, key, value, &missing));
53*26cec326SBarry Smith         break;
54*26cec326SBarry Smith       case ADD_VALUES:
55*26cec326SBarry Smith         PetscCall(PetscHMapIJVQueryAdd(a->ht, key, value, &missing));
56*26cec326SBarry Smith         break;
57*26cec326SBarry Smith       default:
58*26cec326SBarry Smith         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "InsertMode not supported");
59*26cec326SBarry Smith       }
60*26cec326SBarry Smith       if (missing) ++a->dnz[key.i];
61*26cec326SBarry Smith #if defined(TYPE_BS_ON)
62*26cec326SBarry Smith       PetscCall(PetscHSetIJQueryAdd(a->bht, bkey, &missing));
63*26cec326SBarry Smith       if (missing) ++a->bdnz[bkey.i];
64*26cec326SBarry Smith #endif
65*26cec326SBarry Smith     }
66*26cec326SBarry Smith   }
67*26cec326SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
68*26cec326SBarry Smith }
69