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