Actual source code: sftype.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: #include <petsc/private/sfimpl.h>

  3: #if !defined(PETSC_HAVE_MPI_TYPE_GET_ENVELOPE)
  4: #define MPI_Type_get_envelope(datatype,num_ints,num_addrs,num_dtypes,combiner) (*(num_ints)=0,*(num_addrs)=0,*(num_dtypes)=0,*(combiner)=0,1);SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP_SYS,"Need an MPI-2 implementation")
  5: #define MPI_Type_get_contents(datatype,num_ints,num_addrs,num_dtypes,ints,addrs,dtypes) (*(ints)=0,*(addrs)=0,*(dtypes)=0,1);SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP_SYS,"Need an MPI-2 implementation")
  6: #endif
  7: #if !defined(PETSC_HAVE_MPI_COMBINER_DUP)  /* We have no way to interpret output of MPI_Type_get_envelope without this. */
  8: #  define MPI_COMBINER_DUP   0
  9: #endif
 10: #if !defined(PETSC_HAVE_MPI_COMBINER_CONTIGUOUS) && MPI_VERSION < 2
 11: #  define MPI_COMBINER_CONTIGUOUS -1
 12: #endif

 16: PetscErrorCode MPIPetsc_Type_unwrap(MPI_Datatype a,MPI_Datatype *atype)
 17: {
 18:   PetscMPIInt    nints,naddrs,ntypes,combiner;

 22:   MPI_Type_get_envelope(a,&nints,&naddrs,&ntypes,&combiner);
 23:   if (combiner == MPI_COMBINER_DUP) {
 24:     PetscMPIInt  ints[1];
 25:     MPI_Aint     addrs[1];
 26:     MPI_Datatype types[1];
 27:     if (nints != 0 || naddrs != 0 || ntypes != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unexpected returns from MPI_Type_get_envelope()");
 28:     MPI_Type_get_contents(a,0,0,1,ints,addrs,types);
 29:     *atype = types[0];
 30:   } else *atype = a;
 31:   return(0);
 32: }

 36: PetscErrorCode MPIPetsc_Type_compare(MPI_Datatype a,MPI_Datatype b,PetscBool *match)
 37: {
 39:   MPI_Datatype   atype,btype;
 40:   PetscMPIInt    aintcount,aaddrcount,atypecount,acombiner;
 41:   PetscMPIInt    bintcount,baddrcount,btypecount,bcombiner;

 44:   MPIPetsc_Type_unwrap(a,&atype);
 45:   MPIPetsc_Type_unwrap(b,&btype);
 46:   *match = PETSC_FALSE;
 47:   if (atype == btype) {
 48:     *match = PETSC_TRUE;
 49:     return(0);
 50:   }
 51:   MPI_Type_get_envelope(atype,&aintcount,&aaddrcount,&atypecount,&acombiner);
 52:   MPI_Type_get_envelope(btype,&bintcount,&baddrcount,&btypecount,&bcombiner);
 53:   if (acombiner == bcombiner && aintcount == bintcount && aaddrcount == baddrcount && atypecount == btypecount && (aintcount > 0 || aaddrcount > 0 || atypecount > 0)) {
 54:     PetscMPIInt  *aints,*bints;
 55:     MPI_Aint     *aaddrs,*baddrs;
 56:     MPI_Datatype *atypes,*btypes;
 57:     PetscBool    same;
 58:     PetscMalloc6(aintcount,&aints,bintcount,&bints,aaddrcount,&aaddrs,baddrcount,&baddrs,atypecount,&atypes,btypecount,&btypes);
 59:     MPI_Type_get_contents(atype,aintcount,aaddrcount,atypecount,aints,aaddrs,atypes);
 60:     MPI_Type_get_contents(btype,bintcount,baddrcount,btypecount,bints,baddrs,btypes);
 61:     PetscMemcmp(aints,bints,aintcount*sizeof(aints[0]),&same);
 62:     if (same) {
 63:       PetscMemcmp(aaddrs,baddrs,aaddrcount*sizeof(aaddrs[0]),&same);
 64:       if (same) {
 65:         /* This comparison should be recursive */
 66:         PetscMemcmp(atypes,btypes,atypecount*sizeof(atypes[0]),&same);
 67:       }
 68:     }
 69:     PetscFree6(aints,bints,aaddrs,baddrs,atypes,btypes);
 70:     if (same) *match = PETSC_TRUE;
 71:     return(0);
 72:   }
 73:   return(0);
 74: }

 78: /* Check whether a was created via MPI_Type_contiguous from b
 79:  *
 80:  */
 81: PetscErrorCode MPIPetsc_Type_compare_contig(MPI_Datatype a,MPI_Datatype b,PetscInt *n)
 82: {
 84:   MPI_Datatype   atype,btype;
 85:   PetscMPIInt    aintcount,aaddrcount,atypecount,acombiner;

 88:   MPIPetsc_Type_unwrap(a,&atype);
 89:   MPIPetsc_Type_unwrap(b,&btype);
 90:   *n = PETSC_FALSE;
 91:   if (atype == btype) {
 92:     *n = 1;
 93:     return(0);
 94:   }
 95:   MPI_Type_get_envelope(atype,&aintcount,&aaddrcount,&atypecount,&acombiner);
 96:   if (acombiner == MPI_COMBINER_CONTIGUOUS && aintcount >= 1) {
 97:     PetscMPIInt  *aints;
 98:     MPI_Aint     *aaddrs;
 99:     MPI_Datatype *atypes;
100:     PetscMalloc3(aintcount,&aints,aaddrcount,&aaddrs,atypecount,&atypes);
101:     MPI_Type_get_contents(atype,aintcount,aaddrcount,atypecount,aints,aaddrs,atypes);
102:     if (atypes[0] == btype) *n = aints[0];
103:     PetscFree3(aints,aaddrs,atypes);
104:     return(0);
105:   }
106:   return(0);
107: }