Actual source code: memc.c


  2: /*
  3:     We define the memory operations here. The reason we just do not use
  4:   the standard memory routines in the PETSc code is that on some machines
  5:   they are broken.

  7: */
  8: #include <petsc/private/petscimpl.h>
  9: #include <petscbt.h>
 10: #include <../src/sys/utils/ftn-kernels/fcopy.h>

 12: /*@
 13:    PetscMemcmp - Compares two byte streams in memory.

 15:    Not Collective

 17:    Input Parameters:
 18: +  str1 - Pointer to the first byte stream
 19: .  str2 - Pointer to the second byte stream
 20: -  len  - The length of the byte stream
 21:          (both str1 and str2 are assumed to be of length len)

 23:    Output Parameters:
 24: .   e - PETSC_TRUE if equal else PETSC_FALSE.

 26:    Level: intermediate

 28:    Note:
 29:    PetscArraycmp() is preferred
 30:    This routine is anologous to memcmp()

 32: .seealso: PetscMemcpy(), PetscMemcmp(), PetscArrayzero(), PetscMemzero(), PetscArraycmp(), PetscArraycpy(), PetscStrallocpy(),
 33:           PetscArraymove()
 34: @*/
 35: PetscErrorCode PetscMemcmp(const void *str1, const void *str2, size_t len, PetscBool *e)
 36: {
 37:   if (!len) {
 38:     // if e is a bad ptr I guess we just die here then?
 39:     *e = PETSC_TRUE;
 40:     return 0;
 41:   }

 46:   *e = memcmp((char*)str1,(char*)str2,len) ? PETSC_FALSE : PETSC_TRUE;
 47:   return 0;
 48: }

 50: #if defined(PETSC_HAVE_HWLOC)
 51: #include <petsc/private/petscimpl.h>
 52: #include <hwloc.h>

 54: /*@C
 55:      PetscProcessPlacementView - display the MPI process placement by core

 57:   Input Parameter:
 58: .   viewer - ASCII viewer to display the results on

 60:   Level: intermediate

 62:   Notes:
 63:     Requires that PETSc be installed with hwloc, for example using --download-hwloc
 64: @*/
 65: PetscErrorCode PetscProcessPlacementView(PetscViewer viewer)
 66: {
 67:   PetscBool        isascii;
 68:   PetscMPIInt      rank;
 69:   hwloc_bitmap_t   set;
 70:   hwloc_topology_t topology;

 73:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);

 76:   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
 77:   hwloc_topology_init ( &topology);
 78:   hwloc_topology_load ( topology);
 79:   set = hwloc_bitmap_alloc();

 81:   PetscStackCallStandard(hwloc_get_proc_cpubind,topology, getpid(), set, HWLOC_CPUBIND_PROCESS);
 82:   PetscViewerASCIIPushSynchronized(viewer);
 83:   PetscViewerASCIISynchronizedPrintf(viewer,"MPI rank %d Process id: %d coreid %d\n",rank,getpid(),hwloc_bitmap_first(set));
 84:   PetscViewerFlush(viewer);
 85:   hwloc_bitmap_free(set);
 86:   hwloc_topology_destroy(topology);
 87:   return 0;
 88: }
 89: #endif