Actual source code: matelemimpl.h

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: #if !defined(_matelemimpl_h)
  2: #define _matelemimpl_h

  4: #include <El.hpp>
  5:  #include <petsc/private/matimpl.h>
  6:  #include <petscmatelemental.h>

  8: typedef struct {
  9:   PetscInt commsize;
 10:   PetscInt m[2];       /* Number of entries in a local block of the row (column) space */
 11:   PetscInt mr[2];      /* First incomplete/ragged rank of (row) column space.
 12:                           We expose a blocked ordering to the user because that is what all other PETSc infrastructure uses.
 13:                           With the blocked ordering when the number of processes do not evenly divide the vector size,
 14:                           we still need to be able to convert from PETSc/blocked ordering to VC/VR ordering. */
 15:   El::Grid                         *grid;
 16:   El::DistMatrix<PetscElemScalar>  *emat;
 17:   PetscInt                         pivoting;     /* 0: no pivoting; 1: partial pivoting; 2: full pivoting */
 18:   El::DistPermutation              *P,*Q;
 19:   PetscBool                        roworiented;  /* if true, row oriented input (default) */
 20: } Mat_Elemental;

 22: typedef struct {
 23:   El::Grid *grid;
 24:   PetscInt   grid_refct;
 25: } Mat_Elemental_Grid;

 27: /*
 28:  P2RO, RO2P, E2RO and RO2E convert indices between PETSc <-> (Rank,Offset) <-> Elemental 
 29:  (PETSc parallel vectors can be used with Elemental matries without changes)
 30: */
 31: PETSC_STATIC_INLINE void P2RO(Mat A,PetscInt rc,PetscInt p,PetscInt *rank,PetscInt *offset)
 32: {
 33:   Mat_Elemental *a       = (Mat_Elemental*)A->data;
 34:   PetscInt      critical = a->m[rc]*a->mr[rc];
 35:   if (p < critical) {
 36:     *rank   = p / a->m[rc];
 37:     *offset = p % a->m[rc];
 38:   } else {
 39:     *rank   = a->mr[rc] + (p - critical) / (a->m[rc] - 1);
 40:     *offset = (p - critical) % (a->m[rc] - 1);
 41:   }
 42: }
 43: PETSC_STATIC_INLINE void RO2P(Mat A,PetscInt rc,PetscInt rank,PetscInt offset,PetscInt *p)
 44: {
 45:   Mat_Elemental *a = (Mat_Elemental*)A->data;
 46:   if (rank < a->mr[rc]) {
 47:     *p = rank*a->m[rc] + offset;
 48:   } else {
 49:     *p = a->mr[rc]*a->m[rc] + (rank - a->mr[rc])*(a->m[rc]-1) + offset;
 50:   }
 51: }

 53: PETSC_STATIC_INLINE void E2RO(Mat A,PetscInt rc,PetscInt p,PetscInt *rank,PetscInt *offset)
 54: {
 55:   Mat_Elemental *a = (Mat_Elemental*)A->data;
 56:   *rank   = p % a->commsize;
 57:   *offset = p / a->commsize;
 58: }
 59: PETSC_STATIC_INLINE void RO2E(Mat A,PetscInt rc,PetscInt rank,PetscInt offset,PetscInt *e)
 60: {
 61:   Mat_Elemental *a = (Mat_Elemental*)A->data;
 62:   *e = offset * a->commsize + rank;
 63: }

 65: #endif