I think you are almost there. Once you have created your grid manually, you 
should be able to set the manifold ids (which you did not do in your case).

tria.set_all_manifold_ids_on_boundary(1);

and then, after creating the manifold with opencascade, you should attach it to 
the triangulation:

tria.set_manifold(1, manifold_object);

After this, the triangulation is aware of your manifold. 

The stp file you provide is not a three dimensional triangulation, but a two 
dimensional surface embedded on a three-dimensional space, and the opencascade 
utility won’t create a three-dimensional grid. Moreover, that utility is 
thought for single nurbs patches (that define internally a rectangular u-v 
space), and that’s how the grid is generated. In the case of a stanford bunny, 
the surface is triangulated already, so there is no way to define a unique u-v 
space and generate a Triangulation<2,3> out of it, and you have to proceed 
differently.

Ps: your procedure is ok, but remember that if your domain is topologically 
equivalent to a sphere, you should not start with a cube, but with a 
hypersphere (in the hypercube, the coarse cell has all faces on the boundary, 
and this will result after refinement in 8 cells that have three neighboring 
faces on a smooth boundary, with very high defromations, and give you distorted 
cells. If you refine instead the coarse mesh of a sphere attaching the manifold 
of the bunny to the its boundary, every cell only has one face on the boundary, 
and this will be preserved throughout refinement).

L.

> On 14 Jun 2023, at 15:23, Andrey Borzunov <andrey.borzu...@gmail.com> wrote:
> 
> Hi, all.
> 
> I'm trying to create a 3D mesh from CAD file (stanford bunny) to solve a 
> heating equation in this 3D volume.
> 
> I tried gmsh, but failed to produce coarse mesh, every try resulted in too 
> fine meshes (700k cells and more).
> 
> Than i tried manually create initial-mesh with one cell (because 
> `OpenCASCADE::create_triangulation` fails on 3D case). I used FreeCAD to 
> examine .stp file, selected 8 edges and obtained its vertex coordinates using 
> python console, and than used a code from Step14 Exercise_2_3 to create 
> triangulation.
> 
> Finally, i failed to set manifolds correctly.
> I have found similar discussions Bruno & Luca, yy.wayne & Wolfgang and wiki 
> page but failed to understand the way i should modify sources (see 
> attachment, it is slightly modified step-54).
> 
> Thanks in advance!
> 
> Best, Andrey.
>  
> 
> -- 
> 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/46a757a1-33be-42a6-898c-d025246564c5n%40googlegroups.com.
> <step-54.cc>

-- 
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/22DC976E-C1BC-4DBE-9FA3-9D3DA3381DA3%40gmail.com.

Reply via email to