xref: /petsc/src/mat/impls/baij/mpi/mpb_baij.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 #include <../src/mat/impls/baij/mpi/mpibaij.h>
2 
3 PetscErrorCode MatGetMultiProcBlock_MPIBAIJ(Mat mat, MPI_Comm subComm, MatReuse scall, Mat *subMat) {
4   Mat_MPIBAIJ *aij  = (Mat_MPIBAIJ *)mat->data;
5   Mat_SeqBAIJ *aijB = (Mat_SeqBAIJ *)aij->B->data;
6   PetscMPIInt  commRank, subCommSize, subCommRank;
7   PetscMPIInt *commRankMap, subRank, rank, commsize;
8   PetscInt    *garrayCMap, col, i, j, *nnz, newRow, newCol, *newbRow, *newbCol, k, k1;
9   PetscInt     bs = mat->rmap->bs;
10   PetscScalar *vals, *aijBvals;
11 
12   PetscFunctionBegin;
13   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &commsize));
14   PetscCallMPI(MPI_Comm_size(subComm, &subCommSize));
15 
16   /* create subMat object with the relevant layout */
17   if (scall == MAT_INITIAL_MATRIX) {
18     PetscCall(MatCreate(subComm, subMat));
19     PetscCall(MatSetType(*subMat, MATMPIBAIJ));
20     PetscCall(MatSetSizes(*subMat, mat->rmap->n, mat->cmap->n, PETSC_DECIDE, PETSC_DECIDE));
21     PetscCall(MatSetBlockSizes(*subMat, mat->rmap->bs, mat->cmap->bs));
22 
23     /* need to setup rmap and cmap before Preallocation */
24     PetscCall(PetscLayoutSetBlockSize((*subMat)->rmap, mat->rmap->bs));
25     PetscCall(PetscLayoutSetBlockSize((*subMat)->cmap, mat->cmap->bs));
26     PetscCall(PetscLayoutSetUp((*subMat)->rmap));
27     PetscCall(PetscLayoutSetUp((*subMat)->cmap));
28   }
29 
30   /* create a map of comm_rank from subComm to comm - should commRankMap and garrayCMap be kept for reused? */
31   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &commRank));
32   PetscCallMPI(MPI_Comm_rank(subComm, &subCommRank));
33   PetscCall(PetscMalloc1(subCommSize, &commRankMap));
34   PetscCallMPI(MPI_Allgather(&commRank, 1, MPI_INT, commRankMap, 1, MPI_INT, subComm));
35 
36   /* Traverse garray and identify blocked column indices [of offdiag mat] that
37    should be discarded. For the ones not discarded, store the newCol+1
38    value in garrayCMap */
39   PetscCall(PetscCalloc1(aij->B->cmap->n / bs, &garrayCMap));
40   for (i = 0; i < aij->B->cmap->n / bs; i++) {
41     col = aij->garray[i]; /* blocked column index */
42     for (subRank = 0; subRank < subCommSize; subRank++) {
43       rank = commRankMap[subRank];
44       if ((col >= mat->cmap->range[rank] / bs) && (col < mat->cmap->range[rank + 1] / bs)) {
45         garrayCMap[i] = (((*subMat)->cmap->range[subRank] - mat->cmap->range[rank]) / bs + col + 1);
46         break;
47       }
48     }
49   }
50 
51   if (scall == MAT_INITIAL_MATRIX) {
52     /* Now compute preallocation for the offdiag mat */
53     PetscCall(PetscCalloc1(aij->B->rmap->n / bs, &nnz));
54     for (i = 0; i < aij->B->rmap->n / bs; i++) {
55       for (j = aijB->i[i]; j < aijB->i[i + 1]; j++) {
56         if (garrayCMap[aijB->j[j]]) nnz[i]++;
57       }
58     }
59     PetscCall(MatMPIBAIJSetPreallocation(*(subMat), bs, 0, NULL, 0, nnz));
60 
61     /* reuse diag block with the new submat */
62     PetscCall(MatDestroy(&((Mat_MPIBAIJ *)((*subMat)->data))->A));
63 
64     ((Mat_MPIBAIJ *)((*subMat)->data))->A = aij->A;
65 
66     PetscCall(PetscObjectReference((PetscObject)aij->A));
67   } else if (((Mat_MPIBAIJ *)(*subMat)->data)->A != aij->A) {
68     PetscObject obj = (PetscObject)((Mat_MPIBAIJ *)((*subMat)->data))->A;
69 
70     PetscCall(PetscObjectReference((PetscObject)obj));
71 
72     ((Mat_MPIBAIJ *)((*subMat)->data))->A = aij->A;
73 
74     PetscCall(PetscObjectReference((PetscObject)aij->A));
75   }
76 
77   /* Now traverse aij->B and insert values into subMat */
78   PetscCall(PetscMalloc3(bs, &newbRow, bs, &newbCol, bs * bs, &vals));
79   for (i = 0; i < aij->B->rmap->n / bs; i++) {
80     newRow = (*subMat)->rmap->range[subCommRank] + i * bs;
81     for (j = aijB->i[i]; j < aijB->i[i + 1]; j++) {
82       newCol = garrayCMap[aijB->j[j]];
83       if (newCol) {
84         newCol--; /* remove the increment */
85         newCol *= bs;
86         for (k = 0; k < bs; k++) {
87           newbRow[k] = newRow + k;
88           newbCol[k] = newCol + k;
89         }
90         /* copy column-oriented aijB->a into row-oriented vals */
91         aijBvals = aijB->a + j * bs * bs;
92         for (k1 = 0; k1 < bs; k1++) {
93           for (k = 0; k < bs; k++) { vals[k1 + k * bs] = *aijBvals++; }
94         }
95         PetscCall(MatSetValues(*subMat, bs, newbRow, bs, newbCol, vals, INSERT_VALUES));
96       }
97     }
98   }
99   PetscCall(MatAssemblyBegin(*subMat, MAT_FINAL_ASSEMBLY));
100   PetscCall(MatAssemblyEnd(*subMat, MAT_FINAL_ASSEMBLY));
101 
102   /* deallocate temporary data */
103   PetscCall(PetscFree3(newbRow, newbCol, vals));
104   PetscCall(PetscFree(commRankMap));
105   PetscCall(PetscFree(garrayCMap));
106   if (scall == MAT_INITIAL_MATRIX) { PetscCall(PetscFree(nnz)); }
107   PetscFunctionReturn(0);
108 }
109