Hi,

Unlike a lot of FE software deal.II does not actually store, e.g., 27 points 
for a HEX27 element. They are always implicitly defined by the FiniteElement, 
Mapping, and Triangulation put together so you have to compute their values 
explicitly.

Hence, if you want globally unique nodes, the correct things to compute are the 
support points for all locally owned degrees of freedom. The function

DoFTools::map_dofs_to_support_points()

would be what you want but that function doesn't work in parallel. Fortunately 
this function is pretty short and you can write your own version that works on 
the locally owned dofs pretty easily by starting at

https://github.com/dealii/dealii/blob/master/source/dofs/dof_tools.cc#L2102

and doing something like

  1.  Set up a quadrature with fe.get_unit_support_points()
  2.  set up FEValues with your FE, mapping, and that quadrature. Do you have a 
finite element field that describes a position or displacement? If so take a 
look at MappingFEField.
  3.  do fe_values.get_quadrature_points() on each cell, but only save the 
results if the dof corresponding to the quadrature point (i.e., its support 
point) is locally owned.

Do you see what to do?

Best,
David Wells
________________________________
From: dealii@googlegroups.com <dealii@googlegroups.com> on behalf of Joss G. 
<jossgeo...@gmail.com>
Sent: Wednesday, November 17, 2021 7:24 AM
To: deal.II User Group <dealii@googlegroups.com>
Subject: [deal.II] Re: Mesh export, high order


Thank you for your answer.

Following the examples you mentioned, I got a vtu file where I can extract the 
points. As I am using distributed tirangulation, I am getting the (x,y) values 
of the points of each cell, meaning that the points are not unique, there are 
repeated.

DataOut<dim> data_out;
DataOutBase::VtkFlags flags;
flags.write_higher_order_cells = true;
data_out.set_flags(flags);
data_out.attach_dof_handler(dof_handler);
data_out.add_data_vector(solution, "sol");
Vector<float> subdomain(triangulation.n_active_cells());
for (unsigned int i = 0; i < subdomain.size(); ++i)
subdomain(i) = triangulation.locally_owned_subdomain();
data_out.add_data_vector(subdomain, "subdomain");
data_out.build_patches();
data_out.write_vtu_with_pvtu_record("./", "solution", 1, mpi_communicator, 2, 
8);

I would need to get the unique points corresponding to the domain in the same 
order as the operators are assembled. Similar to what I got with 
grid_out_write.msh but including the higher order points. COuld you please 
suggest me how can I approach that?

Thank you again.



El miércoles, 27 de octubre de 2021 a las 22:06:59 UTC+2, mafe...@gmail.com 
escribió:
Hello,

yes, you can use high-order output using the `DataOut` class after setting the 
corresponding `DataOutFlags` flag 
write_higher_order_cells<https://www.dealii.org/current/doxygen/deal.II/structDataOutBase_1_1VtkFlags.html#aa9dd1830c0ff35a2431704c4d45453eb>.

The tutorials step-11, step-44, step-65, step-67, and step-76 demonstrate that 
feature.

Best,
Marc

On Wednesday, October 27, 2021 at 5:05:22 AM UTC-6 jossg...@gmail.com wrote:
Hi all,

I am using high order methods ad I am trying to export the mesh points. I am 
using  grid_out.write_msh(triangulation, out).
Is it any other way to get the  mesh points including the ones corresponding to 
the high orders?

Thank you

--
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 
dealii+unsubscr...@googlegroups.com<mailto:dealii+unsubscr...@googlegroups.com>.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/7095bcb6-fccd-4a28-a2ba-5b6803f4f80cn%40googlegroups.com<https://groups.google.com/d/msgid/dealii/7095bcb6-fccd-4a28-a2ba-5b6803f4f80cn%40googlegroups.com?utm_medium=email&utm_source=footer>.

-- 
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 dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/BN7PR03MB43561ECAF512BA0062499316ED9A9%40BN7PR03MB4356.namprd03.prod.outlook.com.

Reply via email to