[petsc-users] Advice on creating vectors defined on lower dimensional manifolds of a DMPlex

Aldo Bonfiglioli aldo.bonfiglioli at unibas.it
Tue Oct 7 08:24:06 CDT 2025


On 10/6/25 20:10, Matthew Knepley wrote:
> On Mon, Oct 6, 2025 at 1:11 PM Aldo Bonfiglioli 
> <aldo.bonfiglioli at unibas.it> wrote:
>
>     Dear all,
>
>     what is the best approach for defining vectors that "sit" on the
>     (vertices and/or faces) of a given stratum of the "Face Sets" of a
>     DMPlex?
>
>>     DM Object: 3D plex 1 MPI process
>>      type: plex
>>     3D plex in 3 dimensions:
>>      Number of 0-cells per rank: 9261
>>      Number of 1-cells per rank: 59660
>>      Number of 2-cells per rank: 98400
>>      Number of 3-cells per rank: 48000
>>     Labels:
>>      marker: 1 strata with value/size (1 (14402))
>>      celltype: 4 strata with value/size (0 (9261), 1 (59660), 3
>>     (98400), 6 (48000))
>>      depth: 4 strata with value/size (0 (9261), 1 (59660), 2 (98400),
>>     3 (48000))
>>      Face Sets: 6 strata with value/size (1 (800), 2 (800), 3 (800),
>>     4 (800), 5 (800), 6 (800))
>>
>     These vectors are going to be used (for example) to store stresses
>     and heat flux on solid surfaces.
>
>     To be more specific: suppose stratum 3 of the "Face Sets" is a
>     solid wall.
>
>     I want to create a vector that that stores quantities computed on
>     the (800) faces of that wall OR the vertices of that wall.
>
>
> It should be simple to just create such vectors. You request a submesh 
> using that label
>
>   DM subdm;
>   DMLabel label, sublabel;
>
>   DMGetLabel(dm, "Face Sets", &label);
>   DMLabelDuplicate(label, &sublabel);
>   DMPlexLabelComplete(dm, sublabel);
>   DMPlexCreateSubmesh(dm, sublabel, 3, PETSC_TRUE, &subdm)
>   DMLabelDestroy(&sublabel);
>
> Now you can define a PetscFE over this submesh in the same way as any 
> other mesh. Moreover, the subdm contains a mapping back to the 
> original DM, from which you can create a mapping of dofs, so that you 
> can inject the subvector into a larger field if you wish.
>
> If you want to use fields on submeshes inside a PetscDS, so that the 
> Plex manages the solve, the procedure is slightly different, but I can 
> detail it if you want.
>
>   Thanks,
>
>      Matt
>
>     Thanks,
>
>     Aldo
>
>     -- 
>     Dr. Aldo Bonfiglioli
>     Associate professor of Fluid Mechanics
>     Dipartimento di Ingegneria
>     Universita' della Basilicata
>     V.le dell'Ateneo Lucano, 10 85100 Potenza ITALY
>     tel:+39.0971.205203 fax:+39.0971.205215
>     web:https://urldefense.us/v3/__http://docenti.unibas.it/site/home/docente.html?m=002423__;!!G_uCfscf7eWS!Yk-uGKnWGfd4b87jjp07cII-ta6Rld0p_kmImii1ZuYnG2DCFwjpUQkeQ8z3mrWOHPd6-oKoBkZsW8rkoKSWFblHNxPyamZ0CHI$  <https://urldefense.us/v3/__http://docenti.unibas.it/site/home/docente.html?m=002423__;!!G_uCfscf7eWS!dx5g28NqJW34oxLLKP1Fjtp65c0KkvUjelPzjza0lBJtf6uu5ROFqpa2GTX5Cle8L7S_YjHssSDqe6szXd2PEYvYVHHq5mtW8EU$>
>
>
>
> -- 
> 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://urldefense.us/v3/__https://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!Yk-uGKnWGfd4b87jjp07cII-ta6Rld0p_kmImii1ZuYnG2DCFwjpUQkeQ8z3mrWOHPd6-oKoBkZsW8rkoKSWFblHNxPykk2-vZU$  
> <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!Yk-uGKnWGfd4b87jjp07cII-ta6Rld0p_kmImii1ZuYnG2DCFwjpUQkeQ8z3mrWOHPd6-oKoBkZsW8rkoKSWFblHNxPy-jdNsBo$ >

Matthew,

I followed your suggestions, but I face a deadlock when the enclosed 
demonstrator is run in parallel (for the attached dotfile, deadlock 
occurs when nproc = 3).

I suspect it might be due to the fact that in a parallel environment the 
various strata of the "Face Sets" are not necessarily available to all 
processes.

Therefore, not all processes are going to call DMPlexCreateSubmesh with 
a given "value".

The attached piece of code links with petsc-3.24.0.

Thanks,

Aldo


-- 
Dr. Aldo Bonfiglioli
Associate professor of Fluid Mechanics
Dipartimento di Ingegneria
Universita' della Basilicata
V.le dell'Ateneo Lucano, 10 85100 Potenza ITALY
tel:+39.0971.205203 fax:+39.0971.205215
web:https://urldefense.us/v3/__http://docenti.unibas.it/site/home/docente.html?m=002423__;!!G_uCfscf7eWS!Yk-uGKnWGfd4b87jjp07cII-ta6Rld0p_kmImii1ZuYnG2DCFwjpUQkeQ8z3mrWOHPd6-oKoBkZsW8rkoKSWFblHNxPyamZ0CHI$ 
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20251007/43ed4f50/attachment.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: rgmsh.F90
Type: text/x-fortran
Size: 8893 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20251007/43ed4f50/attachment.bin>
-------------- next part --------------
-dm_plex_dim 3
-dm_plex_shape box 
-dm_plex_box_faces 20,20,20
-dm_plex_box_lower 0.,0.,0. 
-dm_plex_box_upper 1.,1.,1.
##-dm_plex_filename cube6.msh
-dm_plex_simplex true
-dm_plex_interpolate 
-dm_plex_check_all
#-dm_view
##-dm_plex_view_labels "marker"
##-dm_plex_view_labels "Face Sets"
-petscpartitioner_view
####-dm_petscsection_view
-options_left


More information about the petsc-users mailing list