On 1/13/23 09:55, Alex Quinlan wrote:
Is there a good way to find the global_dof numbers for a vertex/node?
No, and it's not just because of mesh refinement, but also because you may
want to renumber degrees of freedom as implemented in the DoFRenumbering
namespace. For example, if you have a vector displacement problem, you may
want to number DoFs in such a way that the three displacements at each node
are numbered consecutively, or you may want all x-displacements to be numbered
before all y-displacements before all z-displacements. In other words, at
least *in general*, there is no way to tell what the mapping from vertex to
DoF index is.
There is of course also the issue that DoFs live on only vertices only if you
are using linear elements. For quadratic elements, they also live on edge and
face midpoints.
If not, can I iterate through all of the vertices directly (rather than going
through the cell)? Something like:
for (const auto &vertex : dof_handler.vertex_iterators())
{
if (vertex.distance(proot1) < 0.001)
{
constraints.add_line(vertex->vertex_dof_index(0) );
}}
I'm afraid there is no such way, though you can make your life a bit cheaper
if you do something along the lines of
std::vector<bool> vertex_already_handled (triangulation.n_vertices(),
false);
// loop through all cells
for (const auto &cell : dof_handler.active_cell_iterators())
{
// loop through all of the vertices in a given cell
for (unsigned int v = 0; v<cell->n_vertices() ; v++)
if (vertex_already_handled[cell->vertex_index(v)] == false)
{
if (condition)
register constraint;
vertex_already_handled[cell->vertex_index(v)] = true;
}
}
There is of course also the question of whether you have to do the code
snippet above for every imported node, i.e., whether all of the above is
enclosed in a loop of the form
for (unsigned int boundary_node=0; boundary_node<n_nodes; ++...)
It may be cheaper to move that loop into the loop over all cells, but you're
still ending up with a quadratic complexity in the number of cells/boundary
nodes. If that turns out to be too expensive, you need to sort your boundary
nodes into something like an rtree or similar data structure so that looking
up which boundary nodes are close to a cell's vertex becomes substantially
cheaper.
Best
W.
--
------------------------------------------------------------------------
Wolfgang Bangerth email: [email protected]
www: http://www.math.colostate.edu/~bangerth/
--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see
https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email
to [email protected].
To view this discussion on the web visit
https://groups.google.com/d/msgid/dealii/fecc2d5a-d0c6-25fb-9e00-11231d27b403%40colostate.edu.