Dear Sean,
I think you can find an example for what you want to achieve in the code
of a recent paper https://epubs.siam.org/doi/10.1137/23M158262X
The code is available at
https://zenodo.org/records/10459309
where you can check the file
InverseHomogenization/geometry.h
Best regards
Winni
On 05.09.24 22:13, Sean Carney wrote:
Hi all,
I want to create a 2D domain in deal.II with an internal, curve interface
$\Gamma$. A simple example of what I want to create is in the attached
screenshot (credit to Bartels, Bonito & Tscherner
<https://doi.org/10.1093/imanum/drad004>).
In short, does anyone have any suggestions for how to do this?
The critical thing is that upon mesh refinement, deal.II knows where my
interface is and adds points accordingly (of course, this issue is already
discussed in Step-1). Here, however, the 1D manifold is internal, so maybe
this changes things?
What I tried so far is the following:
I have an analytic expression for my curve (i.e. it can be written in the
form $(\Gamma(y), y)$), so I created two different triangulation objects
using the
GridGenerator::subdivided_hyper_rectangle
and
GridTools::transform
functions. Then, I tried to merge them together. However, I am receiving an
error:
"An error occurred in line <1235> of file
<dealii/source/grid/grid_tools.cc> in function
void dealii::GridTools::{anonymous}::AdjacentCells<2>::push_back(const
dealii::GridTools::{anonymous}::AdjacentCell&)
The violated condition was:
adjacent_cells[1].cell_index == numbers::invalid_unsigned_int"
Of course, I know that the vertices along the common interface, in this
case $\Gamma$, need to be aligned. I believe that this is indeed the case.
My code is simple enough:
void grid_10()
{
const int N = 14;
// create 1st triang
Triangulation<2> tri_1;
std::vector<unsigned int> repetitions(2);
repetitions[0] = N;
repetitions[1] = N;
GridGenerator::subdivided_hyper_rectangle(tri_1,
repetitions,
Point<2>(0.0, 0.0),
Point<2>(2./3., 1.));
GridTools::transform(
[](const Point<2> &in) {
if (in[0] < 0.5)
{
return in;
}
else
{
return Point<2>(2./3. - 1./6.*std::sin(numbers::PI * in[1]), in[1]);
}
},
tri_1);
// create 2nd triang
Triangulation<2> tri_2;
repetitions[0] = N;
repetitions[1] = N;
GridGenerator::subdivided_hyper_rectangle(tri_2,
repetitions,
Point<2>(2./3., 0),
Point<2>(1., 1.));
GridTools::transform(
[](const Point<2> &in) {
if (in[0] > 0.7)
{
return in;
}
else
{
return Point<2>(2./3. - 1./6.*std::sin(numbers::PI * in[1]), in[1]);
}
},
tri_2);
// merge them
Triangulation<2> triangulation;
GridGenerator::merge_triangulations(tri_1,tri_2,triangulation);
print_mesh_info(triangulation, "grid.vtu");
}
Is there an obvious error that I am missing?
Even more critically, will this even accomplish what I want? I see in the
"merge_triangulations" documentation that " In particular, manifolds
attached to the two input triangulations will be lost in the result
triangulation."
Any help is greatly appreciated. Thank you!
Best regards,
-Sean
--
--------------------------------------------------------
Prof. Dr. Winnifried Wollner
Universität Hamburg
MIN Faculty
Department of Mathematics
Bundesstr. 55; 20146 Hamburg; Germany
--------------------------------------------------------
Room : 115
Phone: +49 (0)40 42838-4079
Web : https://www.math.uni-hamburg.de/personen/wollner/
eMail: [email protected]
--------------------------------------------------------
--
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/5caa9d86-1348-4db5-9859-1dd0c51d8461%40uni-hamburg.de.