Dear Wolfgang,
Thanks very much for your reply. I'm happy to hear that you think my
problems can be resolved without too much difficulty. I'd like to try
writing a patch for the mixed 3D meshes. I have done some development work
before, but I have never sent it out of our group; your advice on writing
the steps to write a patch would be most welcome!
I will re-review the relevant code and see if I can put together a general
plan on how to approach this.
Best regards,
Alex
On Wednesday, November 9, 2022 at 4:44:25 PM UTC-5 Wolfgang Bangerth wrote:
> Alex:
>
> > I am coming from a background of solid mechanics FEA and I've been using
> > CalculiX as an openSource program for some time. My group has decided
> > to switch to deal.ii
>
> Good choice -- we applaud you :-)
>
>
> > Part 1: Abaqus .inp format
> > I had success importing hex meshes from the abaqus <.inp> format, but it
> > seems that the UCD reader does not accept tet meshes. I am basing this
> > on the following bit of code in grid_in::read_ucd around Line 948
> >
> > if (((cell_type == "line") && (dim == 1)) ||
> > ((cell_type == "quad") && (dim == 2)) ||
> > ((cell_type == "hex") && (dim == 3)))
> > // found a cell
> > {
> >
> > The absence of
> > ((cell_type == "tet") && (dim == 3))
> > suggests that simplex elements are not supported when importing from
> > .inp. Am I correct in this conclusion. How difficult would it be for
> > me to add this functionality in?
>
> Not very difficult in all likelihood. The key difficulty in adding
> support for other cells in this function is that it uses the
> GeometryInfo<dim> class, which only provides information about hypercube
> cells. For example, we read vertex indices like this in that function:
> for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
> in >> cells.back().vertices[GeometryInfo<dim>::ucd_to_deal[i]];
>
> This needs to be changed to account for the number of vertices a cell
> has, given the cell type, not the dimension. The right approach will be
> to use the ReferenceCell objects defined in namespace ReferenceCells
> (plural), and the other readers in that function will serve as guidance.
> You might want to take a look at the gmsh reader, for example.
>
>
> > Part 2: GMSH format
> > While looking at grid_in.cc, I saw that the MSH reader supported tet
> > elements. I found the FAQ page supplied these examples of simplex
> > meshes:
> https://www.dealii.org/developer/doxygen/deal.II/group__simplex.html
> >
> > I first ran the 2D simplex mesh example, and then I modified it to work
> > for 3D without much problem. Then, I switched to the mixed mesh, and
> > had success running the 2D example. However, I started having problems
> > when I switched to 3D. The mesh seems to import, but then I have issues
> > operating with it. I found several places within fe_simplex_p.cc (L714,
> > L759) where the dimension is asserted to be 2:
> >
> > AssertDimension(dim, 2);
> > This seems like a problem for a 3D analysis, so I am wondering if I need
> > to do something different so that this function doesn't get called.
>
> No, I think this is just not implemented. The function is called on
> mixed meshes when the library wants to know how elements of type
> FE_SimplexP on tetrahedra connect to another element types on other cell
> shapes (namely, on wedges and pyramids). The comparison is simply not
> implemented, but if you want to compute on mixed 3d meshes, you will
> have to implement the missing parts of this function. As long as you
> want to use the same polynomial degree on all cells, this should not
> even be particularly difficult and I would be happy to walk you through
> the necessary steps if you were interested in writing up a patch!
>
> Best
> Wolfgang
>
> --
> ------------------------------------------------------------------------
> 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/885a1c8c-a301-4467-8d16-9a5d156d6f95n%40googlegroups.com.