xref: /petsc/include/petsc/private/petscelemental.h (revision c56ceb406e7f1e4076bf88d4e7444df54c64fc7a)
1 #pragma once
2 
3 #include <petsc/private/matimpl.h>
4 #include <petscmatelemental.h>
5 
6 typedef struct {
7   PetscInt                         commsize;
8   PetscInt                         m[2];  /* Number of entries in a local block of the row (column) space */
9   PetscInt                         mr[2]; /* First incomplete/ragged rank of (row) column space.
10                           We expose a blocked ordering to the user because that is what all other PETSc infrastructure uses.
11                           With the blocked ordering when the number of processes do not evenly divide the vector size,
12                           we still need to be able to convert from PETSc/blocked ordering to VC/VR ordering. */
13   El::Grid                        *grid;
14   El::DistMatrix<PetscElemScalar> *emat;
15   PetscInt                         pivoting; /* 0: no pivoting; 1: partial pivoting; 2: full pivoting */
16   El::DistPermutation             *P, *Q;
17   PetscBool                        roworiented; /* if true, row-oriented input (default) */
18 } Mat_Elemental;
19 
20 typedef struct {
21   El::Grid *grid;
22   PetscInt  grid_refct;
23 } Mat_Elemental_Grid;
24 
25 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat, Mat, PetscReal, Mat);
26 PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat, Mat, Mat);
27 
28 /*
29  P2RO, RO2P, E2RO and RO2E convert indices between PETSc <-> (Rank,Offset) <-> Elemental
30  (PETSc parallel vectors can be used with Elemental matries without changes)
31 */
P2RO(Mat A,PetscInt rc,PetscInt p,PetscInt * rank,PetscInt * offset)32 static inline void P2RO(Mat A, PetscInt rc, PetscInt p, PetscInt *rank, PetscInt *offset)
33 {
34   Mat_Elemental *a        = (Mat_Elemental *)A->data;
35   PetscInt       critical = a->m[rc] * a->mr[rc];
36   if (p < critical) {
37     *rank   = p / a->m[rc];
38     *offset = p % a->m[rc];
39   } else {
40     *rank   = a->mr[rc] + (p - critical) / (a->m[rc] - 1);
41     *offset = (p - critical) % (a->m[rc] - 1);
42   }
43 }
RO2P(Mat A,PetscInt rc,PetscInt rank,PetscInt offset,PetscInt * p)44 static inline void RO2P(Mat A, PetscInt rc, PetscInt rank, PetscInt offset, PetscInt *p)
45 {
46   Mat_Elemental *a = (Mat_Elemental *)A->data;
47   if (rank < a->mr[rc]) {
48     *p = rank * a->m[rc] + offset;
49   } else {
50     *p = a->mr[rc] * a->m[rc] + (rank - a->mr[rc]) * (a->m[rc] - 1) + offset;
51   }
52 }
53 
E2RO(Mat A,PetscInt,PetscInt p,PetscInt * rank,PetscInt * offset)54 static inline void E2RO(Mat A, PetscInt, PetscInt p, PetscInt *rank, PetscInt *offset)
55 {
56   Mat_Elemental *a = (Mat_Elemental *)A->data;
57   *rank            = p % a->commsize;
58   *offset          = p / a->commsize;
59 }
RO2E(Mat A,PetscInt,PetscInt rank,PetscInt offset,PetscInt * e)60 static inline void RO2E(Mat A, PetscInt, PetscInt rank, PetscInt offset, PetscInt *e)
61 {
62   Mat_Elemental *a = (Mat_Elemental *)A->data;
63   *e               = offset * a->commsize + rank;
64 }
65