[petsc-users] Advice on identifying boundary vertices
Matthew Knepley
knepley at gmail.com
Sat Feb 15 20:03:24 CST 2025
On Sat, Feb 15, 2025 at 11:40 AM Aldo Bonfiglioli <
aldo.bonfiglioli at unibas.it> wrote:
> Dear all,
>
> I am trying to identify the boundary vertices that belong to a given
> stratum of the Face sets. This is going to be used to prescribe
> Dirichlet-type bcs.
>
The label "Face Sets" is designed to only contain faces.
> I can select the boundary faces (edges in 2D) of a given stratum using
> PetscCall(DMLabelGetStratumIS(label, stratum, user%bndryfaces(i), ierr))
> and that looks ok;
>
Yes, you get the faces marked with some BC value.
> I then try to identify the boundary vertices by looping over the points
> (edges in 2D/faces in 3D) of the face set of a given stratum and
> retrieve the vertices that make up each individual edge/face.
>
You can do that. You can also do something like
DMLabelDuplicate(faceSets, &newLabel);
DMPlexLabelComplete(dm, newLabel);
which will put in all the points in the transitive closure (such as
vertices).
Then you can just loop over the points in the label, and check for vertices
using DMPlexGetPointDepth().
> For reasons I fail to understand, the above procedure fails to identify
> certain vertices (those circled in the enclosed pdf where different
> colours mark different ranks) in a parallel environment.
>
I do this all the time, so it should not happen. If the above fails, can you
send a small reproducer?
Thanks,
Matt
> Questions:
>
> 1. is there an available function that does what I am trying to do? I
> know that the boundary points can be found in the "marker" label, but I
> need to discriminate among Face Sets of different strata.
>
> 2. what might be wrong in the aforementioned approach?
>
> Thanks,
>
> Aldo
>
> --
> Dr. Aldo Bonfiglioli
> Associate professor of Fluid Machines
> Scuola 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!YLL7mzBlK9H_sEdtXob1AVklf8cVLhT3NTDdExIZsdI3xMfNHBoLXr92BmYzCXOuJqtm2L6OU4DJKV_81AkjvQf-Kra-Ym8P0tU$
>
--
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!cwtigEBvRpcpoLeV3aRNqP6fHwoAHN_2SW4Rjbh_ZqKelJ54Nvgncprg0QimeBjNfY5Ox-8qG6AzP5ZqCg-M$ <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!cwtigEBvRpcpoLeV3aRNqP6fHwoAHN_2SW4Rjbh_ZqKelJ54Nvgncprg0QimeBjNfY5Ox-8qG6AzP4FxXi_H$ >
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20250215/30bbf548/attachment.html>
More information about the petsc-users
mailing list