[petsc-users] Inquiry regarding DMAdaptLabel function

Pierre Jolivet pierre at joliv.et
Mon Feb 27 09:44:49 CST 2023



> On 27 Feb 2023, at 4:42 PM, Matthew Knepley <knepley at gmail.com> wrote:
> 
> On Mon, Feb 27, 2023 at 10:26 AM Pierre Jolivet <pierre at joliv.et <mailto:pierre at joliv.et>> wrote:
>>> On 27 Feb 2023, at 4:16 PM, Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>> 
>>> On Mon, Feb 27, 2023 at 10:13 AM Pierre Jolivet <pierre at joliv.et <mailto:pierre at joliv.et>> wrote:
>>>>> On 27 Feb 2023, at 3:59 PM, Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>>>> 
>>>>> On Mon, Feb 27, 2023 at 9:53 AM Zongze Yang <yangzongze at gmail.com <mailto:yangzongze at gmail.com>> wrote:
>>>>>> Hi, Matt,
>>>>>> 
>>>>>> I tested coarsening a mesh by using ParMMg without firedrake, and found some issues:
>>>>>>  see the code and results here:  https://gitlab.com/petsc/petsc/-/issues/1331
>>>>>> 
>>>>>> Could you have a look and give some comments or suggestions?
>>>>> 
>>>>> I replied on the issue. More generally, the adaptive refinement software has not seen wide use
>>>> 
>>>> :)
>>>> Matt probably meant “the _DMPlex interface_ to adaptive refinement software has not seen wide use”, Mmg has been rather widely used for 10+ years (here is a 13-year old presentation https://www.ljll.math.upmc.fr/hecht/ftp/ff++days/2010/exposes/Morice-MeshMetric.pdf).
>>> 
>>> The interface is certainly new, but even ParMMG is only from Nov 2016, which is very new if you are an old person :)
>> 
>> Indeed. In fact, I do believe we should add a DMPlex mechanism to centralize (redistribute on a single process) a DMPlex and to call Mmg instead of ParMmg.
>> It would certainly not be scalable for large meshes but:
>> 1) there is no need for ParMmg on small-/medium-scale meshes
>> 2) Mmg is more robust than ParMmg at this point in time
>> 3) Mmg has more feature than ParMmg at this point in time, e.g., implicit remeshing using a level-set
>> 4) there is more industry money funnelled into Mmg than into ParMmg 
>> I think the mechanism I mentioned initially was in the TODO list of the Firedrake people (or yours?), maybe it’s already done, but in any case it’s not hooked in the Mmg adaptor code, though it should (erroring out in the case where the communicator is of size greater than one would then not happen anymore).
> 
> Yes, we used to do the same thing with partitioners. We can use DMPlexGather().
> 
> I thought MMG only did 2D and ParMMG only did 3D, but this must be wrong now. Can MMG do both?

Mmg does 2D, 3D, and 3D surfaces.
ParMmg only does 3D (with no short-term plan for 2D or 3D surfaces).

Thanks,
Pierre

>   Thanks,
> 
>      Matt
>  
>> Thanks,
>> Pierre
>> 
>>>   Thanks,
>>> 
>>>     Matt
>>>  
>>>> Thanks,
>>>> Pierre
>>>> 
>>>>> , and I expect
>>>>> more of these kinds of bugs until more people use it.
>>>>> 
>>>>>   Thanks,
>>>>> 
>>>>>      Matt
>>>>>  
>>>>>> Best wishes,
>>>>>> Zongze
>>>>>> 
>>>>>> 
>>>>>> On Mon, 27 Feb 2023 at 20:19, Matthew Knepley <knepley at gmail.com <mailto:knepley at gmail.com>> wrote:
>>>>>>> On Sat, Feb 18, 2023 at 6:41 AM Zongze Yang <yangzongze at gmail.com <mailto:yangzongze at gmail.com>> wrote:
>>>>>>>> Another question on mesh coarsening is about `DMCoarsen` which will fail when running in parallel.
>>>>>>>> 
>>>>>>>> I generate a mesh in Firedrake, and then create function space and functions, after that, I get the dmplex and coarsen it.
>>>>>>>> When running in serials, I get the mesh coarsened correctly. But it failed with errors in ParMMG when running parallel.
>>>>>>>> 
>>>>>>>> However, If I did not create function space and functions on the original mesh, everything works fine too.
>>>>>>>> 
>>>>>>>> The code and the error logs are attached.
>>>>>>> 
>>>>>>> I believe the problem is that Firedrake and PETSc currently have incompatible coordinate spaces. We are working
>>>>>>> to fix this, and I expect it to work by this summer.
>>>>>>> 
>>>>>>>   Thanks,
>>>>>>> 
>>>>>>>      Matt
>>>>>>>  
>>>>>>>> Thank you for your time and attention。
>>>>>>>> 
>>>>>>>> Best wishes,
>>>>>>>> Zongze
>>>>>>>> 
>>>>>>>> 
>>>>>>>> On Sat, 18 Feb 2023 at 15:24, Zongze Yang <yangzongze at gmail.com <mailto:yangzongze at gmail.com>> wrote:
>>>>>>>>> Dear PETSc Group,
>>>>>>>>> 
>>>>>>>>> I am writing to inquire about the function DMAdaptLabel in PETSc. 
>>>>>>>>> I am trying to use it coarse a mesh, but the resulting mesh is refined.
>>>>>>>>> 
>>>>>>>>> In the following code, all of the `adpat` label values were set to 2 (DM_ADAPT_COARSEN).
>>>>>>>>> There must be something wrong. Could you give some suggestions?
>>>>>>>>>  
>>>>>>>>> ```python
>>>>>>>>> from firedrake import *
>>>>>>>>> from firedrake.petsc import PETSc
>>>>>>>>> 
>>>>>>>>> def mark_all_cells(mesh):
>>>>>>>>>     plex = mesh.topology_dm
>>>>>>>>>     with PETSc.Log.Event("ADD_ADAPT_LABEL"):
>>>>>>>>>         plex.createLabel('adapt')
>>>>>>>>>         cs, ce = plex.getHeightStratum(0)
>>>>>>>>>         for i in range(cs, ce):
>>>>>>>>>             plex.setLabelValue('adapt', i, 2)
>>>>>>>>>             
>>>>>>>>>     return plex
>>>>>>>>> 
>>>>>>>>> mesh = RectangleMesh(10, 10, 1, 1)
>>>>>>>>> 
>>>>>>>>> x = SpatialCoordinate(mesh)
>>>>>>>>> V = FunctionSpace(mesh, 'CG', 1)
>>>>>>>>> f = Function(V).interpolate(10 + 10*sin(x[0]))
>>>>>>>>> triplot(mesh)
>>>>>>>>> 
>>>>>>>>> plex = mark_all_cells(mesh)
>>>>>>>>> new_plex = plex.adaptLabel('adapt')
>>>>>>>>> mesh = Mesh(new_plex)
>>>>>>>>> triplot(mesh)
>>>>>>>>> ```
>>>>>>>>> 
>>>>>>>>> Thank you very much for your time.
>>>>>>>>> 
>>>>>>>>> Best wishes,
>>>>>>>>> Zongze
>>>>>>> 
>>>>>>> 
>>>>>>> -- 
>>>>>>> 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/>
>>>>> 
>>>>> 
>>>>> -- 
>>>>> 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/>
>>> 
>>> 
>>> -- 
>>> 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/>
> 
> 
> -- 
> 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-users/attachments/20230227/00a3667e/attachment.html>


More information about the petsc-users mailing list