[petsc-dev] DMPlexInterpolate after DMPlexDistribute

Pierre Jolivet pierre.jolivet at enseeiht.fr
Wed Apr 29 04:12:48 CDT 2020



> On 29 Apr 2020, at 10:14 AM, Hapla Vaclav <vaclav.hapla at erdw.ethz.ch> wrote:
> 
> 
> 
>> On 28 Apr 2020, at 10:30, Pierre Jolivet <pierre.jolivet at enseeiht.fr <mailto:pierre.jolivet at enseeiht.fr>> wrote:
>> 
>> 
>> 
>>> On 14 Apr 2020, at 2:36 PM, Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>> 
>>> On Tue, Apr 14, 2020 at 6:36 AM Pierre Jolivet <pierre.jolivet at enseeiht.fr <mailto:pierre.jolivet at enseeiht.fr>> wrote:
>>> Hello,
>>> I’d like to call DMPlexInterpolate after DMPlexDistribute and not the other way around for performance reasons (please stop me here if this is equivalent).
>>> When there is no overlap in DMPlexDistribute, it goes through fine.
>>> If there is overlap, I run into an error.
>>> Is this the expected behavior?
>> 
>> Sorry for taking so long to get back at this.
>> I thought I got everything setup, but when looking at the performance, I’m quite surprised.
>> I rewound and looked at src/dm/impls/plex/tutorials/ex2.c by just adding the DMPlexDistribute step. I guess this is equivalent to your steps #1 and #2.
>> diff --git a/src/dm/impls/plex/tutorials/ex2.c b/src/dm/impls/plex/tutorials/ex2.c
>> index a069d922b2..5e5c4fb584 100644
>> --- a/src/dm/impls/plex/tutorials/ex2.c
>> +++ b/src/dm/impls/plex/tutorials/ex2.c
>> @@ -65,2 +65,5 @@ static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
>>      ierr = DMPlexCreateFromFile(comm, user->filename, user->interpolate, dm);CHKERRQ(ierr);
>> +    DM dmParallel;
>> +    ierr = DMPlexDistribute(*dm, 0, NULL, &dmParallel);CHKERRQ(ierr);
>> +    ierr = DMDestroy(&dmParallel);CHKERRQ(ierr);
>>    }
>> 
>> With a cube of 741000 nodes 4574068 elements.
>> $ cat untitled.geo && ~/gmsh-4.5.4-Linux64/bin/gmsh untitled.geo -bin -3
>> //+
>> SetFactory("OpenCASCADE");
>> Box(1) = {-0.5, -0.5, 0, 1, 1, 1};
>> Characteristic Length {:} = 0.01;
>> 
>> I get the following timings (optimized build, scalar-type=real, 64-bit-indices=false, SKL cluster).
>> $ mpirun  -n 120 ./ex2 -filename untitled.msh -log_view -interpolate true
>> DMPlexInterp           1 1.0 1.0444e+02 1.0 0.00e+00 0.0 0.0e+00 0.0e+00 3.0e+00 73  0  0  0  5  73  0  0  0  6     0
>> DMPlexDistribute       1 1.0 2.8793e+01 1.0 0.00e+00 0.0 3.6e+03 6.2e+05 3.3e+01 21  0100100 60  21  0100100 69     0
>> $ mpirun  -n 120 ./ex2 -filename untitled.msh -log_view -interpolate false
>> DMPlexDistribute       1 1.0 7.0265e+02 1.0 0.00e+00 0.0 3.1e+03 2.0e+05 2.6e+01 99  0100100 58  99  0100100 68     0
>> 
>> Do you think I messed up something else?
>> Our legacy implementation, on the same mesh, takes around 20 seconds (accounting for the ParMETIS call) to partition and distribute the mesh and generate the underlying structures for our FEM kernels (kd tree, distributed numbering, so on and so forth).
>> Of course, DMPlex is much more versatile and generic, so I’d understand if it’s a little bit slower, but that’s quite a lot slower right now.
>> Are there any option I could try to play with to speed things up?
> 
> 
> Hi Pierre,
> 
> sorry, I got this thread delayed because our spam filter caught it all.

No problem!

> From just looking at your excerpt, I think the difference could be that we are loading HDF5/XDMF in parallel so DMPlexInterpolate() within DMPlexCreateFromFile() is also done in parallel.
> 
> However, if you load GMSH, I think it's loaded in serial, and then DMPlexInterpolate() is also serial if you just pass interpolate=PETSC_TRUE to DMPlexCreateFromFile().
> 
> So I would try
> ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_FALSE, dm);CHKERRQ(ierr);
> ierr = DMPlexDistribute(*dm, 0, NULL, &dmParallel);CHKERRQ(ierr);
> ierr = DMDestroy(dm);CHKERRQ(ierr);
> ierr = DMPlexInterpolate(dmParallel,dm);CHKERRQ(ierr);

I tried that with with -interpolate false in ex2.c, and I got terrible performance (DMPlexDistribute       1 1.0 7.0265e+02).
But Stefano helped me figure this out.

> Obviously, the API should be improved somehow to avoid these issues.
> https://gitlab.com/petsc/petsc/-/issues/615 <https://gitlab.com/petsc/petsc/-/issues/615>
I’ll track this closely.

Thank you,
Pierre

> Cheers,
> Vaclav
> 
>> 
>> Thanks in advance for your help,
>> Pierre
>> 
>>> Yes. What we want you to do is:
>>> 
>>> 1) Load/generate mesh
>>> 
>>> 2) Distribute (this can be done at the same time as load with parallel load)
>>> 
>>> 3) Interpolate (this is also an option from parallel load)
>>> 
>>> 4) If necessary, redistribute for load balance
>>> 
>>> 5) Construct overlap
>>> 
>>>     When you pass '1' below to DMDistribute(), it distributes as normal and then calls
>>> 
>>>       https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexDistributeOverlap.html <https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexDistributeOverlap.html>
>>> 
>>>     at the end. So you just postpone calling that until you have interpolated.
>>> 
>>>   Thanks,
>>> 
>>>     Matt
>>>  
>>> Here is a MWE.
>>> $ patch -p1 < patch.txt
>>> $ cd src/dm/impls/plex/tests/
>>> $ make ex18
>>> $ mpirun -n 2 ./ex18 -distribute -interpolate after_distribute
>>> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
>>> [0]PETSC ERROR: Petsc has generated inconsistent data
>>> [0]PETSC ERROR: Point SF contains 1 which is a cell
>>> 
>>> Thanks,
>>> Pierre
>>> 
>>> diff --git a/src/dm/impls/plex/tests/ex18.c b/src/dm/impls/plex/tests/ex18.c
>>> index 07421b3522..dd62be58e5 100644
>>> --- a/src/dm/impls/plex/tests/ex18.c
>>> +++ b/src/dm/impls/plex/tests/ex18.c
>>> @@ -806 +806 @@ static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
>>> -    ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
>>> +    ierr = DMPlexDistribute(*dm, 1, NULL, &pdm);CHKERRQ(ierr);
>>> 
>>> 
>>> -- 
>>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>>> -- Norbert Wiener
>>> 
>>> https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20200429/00ce27a0/attachment-0001.html>


More information about the petsc-dev mailing list