<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, Apr 11, 2014 at 12:11 PM, Andrés Alessandro León Baldelli <span dir="ltr"><<a href="mailto:a.leon.baldelli@gmail.com" target="_blank">a.leon.baldelli@gmail.com</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Hi,<br>
I am trying to use the new Plex to HDF5 interface. From what I understand, it only works in sequential mode.<br>
I would like to write a sequential DM, created e.g. from an exodusII mesh via DMPlexCreateExodusFromFile.<br>
What puzzles me is that I only can do it if I run on a single processor, since it crashes as soon its  communicator contains more than one proc.<br>
<br>
My understanding is that, if numproc>1, the crash is due to the call to ISView in src/dm/impls/plex/plex.c:785 which fails because its two arguments do not share the same communicator.<br>
<br>
Could you help me understand a little more what is going on, what I am doing wrong and what could help solving this issue?<br></blockquote><div><br></div><div>Use 'next', not 'master'</div><div><br></div>
<div>   Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
I attach the error message, a minimal code to reproduce the crash and a test mesh.<br>
<br>
Best,<br>
<br>
————— begin output—————<br>
geothermal:dmcomplex kumiori$ mpirun -n 2 ./testHDF5Plex -i data/TwoTri.gen<br>
cubit(/Users/blaise/Development/DMComplex/data/junkTri.gen): 08/31/2013: 19:01:0 in 2 dimensions:<br>
  0-cells: 6 0<br>
  2-cells: 4 0<br>
Labels:<br>
  depth: 2 strata of sizes (6, 4)<br>
  Cell Sets: 1 strata of sizes (4)<br>
[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>
[0]PETSC ERROR: Arguments must have same communicators<br>
[0]PETSC ERROR: Different communicators in the two objects: Argument # 1 and 2 flag 3<br>
[0]PETSC ERROR: See http://<a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">http://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br>
[0]PETSC ERROR: Petsc Development GIT revision: v3.4.4-3957-g24a0f54  GIT Date: 2014-04-02 15:35:24 -0500<br>
[0]PETSC ERROR: ./testHDF5Plex on a Darwin-llvm50-g named <a href="http://geothermal.lsu.edu" target="_blank">geothermal.lsu.edu</a> by kumiori Fri Apr 11 12:06:44 2014<br>
[0]PETSC ERROR: Configure options --CFLAGS=-Wno-unused --download-chaco=1 --download-exodusii=1 --download-hdf5=1 --download-metis=1 --download-mpich=1 --download-netcdf=1 --download-parmetis=1 --download-sowing=1 --download-triangle=1 --download-yaml=1 --with-cmake=cmake --with-debugging=1 --with-pic --with-shared-libraries=1 --with-x11=1 PETSC_ARCH=Darwin-llvm50-g<br>

[0]PETSC ERROR: #1 ISView() line 756 in /opt/HPC/petsc/src/vec/is/is/interface/index.c<br>
[0]PETSC ERROR: #2 DMPlexView_HDF5() line 785 in /opt/HPC/petsc/src/dm/impls/plex/plex.c<br>
[0]PETSC ERROR: #3 DMView_Plex() line 812 in /opt/HPC/petsc/src/dm/impls/plex/plex.c<br>
[0]PETSC ERROR: #4 main() line 61 in /Users/kumiori/Code/dmcomplex/testHDF5Plex.c<br>
[0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint@mcs.anl.gov----------<br>
application called MPI_Abort(MPI_COMM_WORLD, 80) - process 0<br>
[cli_0]: aborting job:<br>
application called MPI_Abort(MPI_COMM_WORLD, 80) - process 0<br>
HDF5-DIAG: Error detected in HDF5 (1.8.10-patch1) MPI-process 1:<br>
  #000: H5F.c line 2058 in H5Fclose(): decrementing file ID failed<br>
    major: Object atom<br>
    minor: Unable to close file<br>
  #001: H5I.c line 1479 in H5I_dec_app_ref(): can't decrement ID ref count<br>
    major: Object atom<br>
    minor: Unable to decrement reference count<br>
  #002: H5F.c line 1835 in H5F_close(): can't close file<br>
    major: File accessability<br>
    minor: Unable to close file<br>
  #003: H5F.c line 1997 in H5F_try_close(): problems closing file<br>
    major: File accessability<br>
    minor: Unable to close file<br>
  #004: H5F.c line 1142 in H5F_dest(): low level truncate failed<br>
    major: File accessability<br>
    minor: Write failed<br>
  #005: H5FD.c line 1897 in H5FD_truncate(): driver truncate request failed<br>
    major: Virtual File Layer<br>
    minor: Can't update object<br>
  #006: H5FDmpio.c line 1984 in H5FD_mpio_truncate(): MPI_File_set_size failed<br>
    major: Internal error (too specific to document in detail)<br>
    minor: Some MPI function failed<br>
  #007: H5FDmpio.c line 1984 in H5FD_mpio_truncate(): Invalid argument, error stack:<br>
MPI_FILE_SET_SIZE(74): Inconsistent arguments to collective routine<br>
    major: Internal error (too specific to document in detail)<br>
    minor: MPI Error String<br>
[0]PETSC ERROR: ------------------------------------------------------------------------<br>
[0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range<br>
[0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger<br>
[0]PETSC ERROR: or see <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC" target="_blank">http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC</a> ERROR: or try <a href="http://valgrind.org" target="_blank">http://valgrind.org</a> on GNU/linux and Apple Mac OS X to find memory corruption errors<br>

[0]PETSC ERROR: likely location of problem given in stack below<br>
[0]PETSC ERROR: ---------------------  Stack Frames ------------------------------------<br>
[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,<br>
[0]PETSC ERROR:       INSTEAD the line number of the start of the function<br>
[0]PETSC ERROR:       is given.<br>
[0]PETSC ERROR: [0] PetscError line 363 /opt/HPC/petsc/src/sys/error/err.c<br>
[0]PETSC ERROR: [0] ISView line 750 /opt/HPC/petsc/src/vec/is/is/interface/index.c<br>
[0]PETSC ERROR: [0] DMPlexView_HDF5 line 656 /opt/HPC/petsc/src/dm/impls/plex/plex.c<br>
[0]PETSC ERROR: [0] DMView_Plex line 803 /opt/HPC/petsc/src/dm/impls/plex/plex.c<br>
[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------<br>
[0]PETSC ERROR: Signal received<br>
[0]PETSC ERROR: See http://<a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank">http://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br>
[0]PETSC ERROR: Petsc Development GIT revision: v3.4.4-3957-g24a0f54  GIT Date: 2014-04-02 15:35:24 -0500<br>
[0]PETSC ERROR: ./testHDF5Plex on a Darwin-llvm50-g named <a href="http://geothermal.lsu.edu" target="_blank">geothermal.lsu.edu</a> by kumiori Fri Apr 11 12:06:44 2014<br>
[0]PETSC ERROR: Configure options --CFLAGS=-Wno-unused --download-chaco=1 --download-exodusii=1 --download-hdf5=1 --download-metis=1 --download-mpich=1 --download-netcdf=1 --download-parmetis=1 --download-sowing=1 --download-triangle=1 --download-yaml=1 --with-cmake=cmake --with-debugging=1 --with-pic --with-shared-libraries=1 --with-x11=1 PETSC_ARCH=Darwin-llvm50-g<br>

[0]PETSC ERROR: #5 User provided function() line 0 in  unknown file<br>
application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0<br>
[cli_0]: aborting job:<br>
application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0<br>
<br>
===================================================================================<br>
=   BAD TERMINATION OF ONE OF YOUR APPLICATION PROCESSES<br>
=   EXIT CODE: 59<br>
=   CLEANING UP REMAINING PROCESSES<br>
=   YOU CAN IGNORE THE BELOW CLEANUP MESSAGES<br>
===================================================================================<br>
geothermal:dmcomplex kumiori$<br>
<br>
————— end output—————<br>
<br>
—————begin C source code testHDF5plex.c ———————<br>
/* usage: testHDF5plex -i meshfilename<br>
*/<br>
#include <petscdmplex.h><br>
#include <petscsf.h><br>
#include <exodusII.h><br>
#include <petsc-private/dmimpl.h>     /*I      "petscdm.h"     I*/<br>
#include <petscsf.h><br>
#include <petscviewerhdf5.h><br>
#include <petscdmplex.h><br>
<br>
#undef __FUNCT__<br>
#define __FUNCT__ "main"<br>
int main(int argc,char **argv)<br>
{<br>
        PetscErrorCode      ierr;<br>
        DM                  dm;<br>
        char                ifilename[PETSC_MAX_PATH_LEN];<br>
        PetscBool           flg;<br>
        int                 CPU_word_size = 0,IO_word_size = 0,exoid;<br>
        float               version;<br>
        PetscMPIInt         numproc,rank;<br>
        PetscViewer         hdf5Viewer;<br>
        PetscInt            dim;<br>
        PetscInt            numFields = 2;<br>
        PetscInt            numComp[2] = {1,1};<br>
        PetscInt            numDof[6] = {1, 0, 0,<br>
                                         0, 0, 1}; /*{Vertex, Edge, Cell} */<br>
        PetscInt            bcFields[1] = {0}, numBC=0;<br>
        IS                  bcPoints[1] = {NULL};<br>
        PetscSection        section;<br>
        MPI_Comm            comm;<br>
<br>
        ierr = PetscInitialize(&argc, &argv, (char*)0, help);CHKERRQ(ierr);<br>
        comm = PETSC_COMM_WORLD;<br>
<br>
        ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","testHDF5Plex Test Options","none");CHKERRQ(ierr);<br>
        ierr = PetscOptionsEnd();<br>
        ierr = PetscOptionsGetString(PETSC_NULL,"-i",ifilename,sizeof ifilename,&flg);CHKERRQ(ierr);<br>
<br>
        ierr = MPI_Comm_size(PETSC_COMM_WORLD,&numproc);<br>
        ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);<br>
<br>
        ierr = DMPlexCreateExodusFromFile(PETSC_COMM_WORLD, ifilename, PETSC_FALSE, &dm);CHKERRQ(ierr);<br>
        ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);<br>
        ierr = DMPlexCreateSection(dm, dim, numFields, numComp, numDof, numBC, bcFields, bcPoints, NULL, &section);CHKERRQ(ierr);<br>
        ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr);<br>
        ierr = PetscViewerHDF5Open(PetscObjectComm((PetscObject) dm), "out/seqDM.h5", FILE_MODE_WRITE, &hdf5Viewer); CHKERRQ(ierr);<br>
// Crashes if numproc > 1,<br>
        ierr = DMView_Plex(dm, hdf5Viewer); CHKERRQ(ierr);<br>
        PetscViewerDestroy(&hdf5Viewer);<br>
<br>
        ierr = DMDestroy(&dm);CHKERRQ(ierr);<br>
        ierr = PetscFinalize();<br>
        return 0;<br>
}<br>
<br>
—————end C source code testHDF5plex.c ———————<br>
<br>
————— test mesh<br>
<br>Andrés A. León Baldelli<br>
—<br>
LSU Center for Computation & Technology<br>
2099 Louisiana Digital Media Center<br>
<br>
<br></blockquote></div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
-- Norbert Wiener
</div></div>