#include /*I "petscdmplex.h" I*/ /* TODO PetscArrayExchangeBegin/End */ /* TODO blocksize */ /* TODO move to API ? */ static PetscErrorCode ExchangeArrayByRank_Private(PetscObject obj, MPI_Datatype dt, PetscInt nsranks, const PetscMPIInt sranks[], PetscInt ssize[], const void *sarr[], PetscInt nrranks, const PetscMPIInt rranks[], PetscInt *rsize_out[], void **rarr_out[]) { PetscInt r; PetscInt *rsize; void **rarr; MPI_Request *sreq, *rreq; PetscMPIInt tag, unitsize; MPI_Comm comm; PetscFunctionBegin; CHKERRMPI(MPI_Type_size(dt, &unitsize)); CHKERRQ(PetscObjectGetComm(obj, &comm)); CHKERRQ(PetscMalloc2(nrranks, &rsize, nrranks, &rarr)); CHKERRQ(PetscMalloc2(nrranks, &rreq, nsranks, &sreq)); /* exchange array size */ CHKERRQ(PetscObjectGetNewTag(obj,&tag)); for (r=0; rcoordinates && !dm->coordinatesLocal,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM coordinates must be set"); CHKERRQ(PetscSFSetUp(sf)); CHKERRQ(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote)); /* Expand sent cones per rank */ CHKERRQ(SortByRemote_Private(sf, &rmine1, &rremote1)); CHKERRQ(GetRecursiveConeCoordinatesPerRank_Private(dm, sf, rmine1, &sntCoordinatesPerRank)); /* Create inverse SF */ CHKERRQ(PetscSFGetMultiSF(sf,&msf)); CHKERRQ(PetscSFCreateInverseSF(msf,&imsf)); CHKERRQ(PetscSFSetUp(imsf)); CHKERRQ(PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL)); CHKERRQ(PetscSFGetRootRanks(imsf, &niranks, &iranks, &iroffset, &irmine, &irremote)); /* Compute original numbering of multi-roots (referenced points) */ CHKERRQ(PetscSFComputeMultiRootOriginalNumberingByRank_Private(sf, imsf, &mine_orig_numbering)); /* Expand coordinates of the referred cones per rank */ CHKERRQ(GetRecursiveConeCoordinatesPerRank_Private(dm, imsf, mine_orig_numbering, &refCoordinatesPerRank)); /* Send the coordinates */ CHKERRQ(ExchangeVecByRank_Private((PetscObject)sf, nranks, ranks, sntCoordinatesPerRank, niranks, iranks, &recCoordinatesPerRank)); /* verbose output */ CHKERRQ(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_check_cones_conform_on_interfaces_verbose", &verbose, NULL)); if (verbose) { PetscViewer sv, v = PETSC_VIEWER_STDOUT_WORLD; CHKERRQ(PetscViewerASCIIPrintf(v, "============\nDMPlexCheckInterfaceCones output\n============\n")); CHKERRQ(PetscViewerASCIIPushSynchronized(v)); CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", myrank)); for (r=0; r