Refining a cylinder only in x

37 views
Skip to first unread message

Lucas Campos

unread,
Sep 18, 2017, 3:50:20 AM9/18/17
to deal.II User Group
Dear all,

I am trying to refine a cylinder in the x direction only, but it seems that when I use 

> cell->set_refine_flag(RefinementCase<3>::cut_x);

the whole cylinder is refined. Indeed, if I run the above line for every active cell *once* I would expect the number of cells to double. However, they quintuple!

Any ideas on how to fix this?

I am using an adapted version of step-49. The offending function is outlined below, and the full program is attached, together with its output.

void grid_cylinder()
{
  Triangulation<3> triangulation;
  static const CylindricalManifold<3> boundary;

  GridGenerator::cylinder (triangulation, 0.05, 2);
  triangulation.set_manifold (0, boundary);

  print_mesh_info (triangulation, "pre_grid");
  for (int i=0; i<1; i++) {
      for(auto cell: triangulation.active_cell_iterators()) {
          cell->set_refine_flag(RefinementCase<3>::cut_x);
      }
      triangulation.execute_coarsening_and_refinement();
  }

  print_mesh_info (triangulation, "grid");
}
cylinder.cc
cylinder.out

Bruno Turcksin

unread,
Sep 18, 2017, 8:44:40 AM9/18/17
to deal.II User Group
Lucas,


On Monday, September 18, 2017 at 3:50:20 AM UTC-4, Lucas Campos wrote:
I am trying to refine a cylinder in the x direction only, but it seems that when I use 

> cell->set_refine_flag(RefinementCase<3>::cut_x);

the whole cylinder is refined. Indeed, if I run the above line for every active cell *once* I would expect the number of cells to double. However, they quintuple!
It looks good but how does the mesh look like? Here, you can see that cut_x is in the local coordinate system not the global one. So I would not be surprise if the mesh doesn't look like you think it should. Also if you do more than one refinement you should also use prepare_coarsening_and_refinement() before you refine your mesh.

Best,

Bruno

Lucas Campos

unread,
Sep 18, 2017, 10:29:09 AM9/18/17
to deal.II User Group
Dear Bruno,

You will find attached the resulting grids. While I originally expected to find new divisions only along the x-direction, your previous answer tells me it actually cuts in some local coordinate (whatever that might be). In this case, what is the best way to refine the cylinder in a single direction? I tried adding an if(cell->boundary_id() == 0) before setting the flag, but it did not really help. The mesh is still being refined in all directions.

I ask this because I plan to do some transformations on this cylinder, to turn it into a coil. In this case, the final result is much better if there are enough cuts in the x-axis.

Bests
grid_unrefined.png
grid_refined.png

Lucas Campos

unread,
Sep 18, 2017, 10:32:20 AM9/18/17
to deal.II User Group
Also, here is the new program, with the added if and prepare_coarsening_and_refinement(). You will find the interesting part between lines 85 and 104.
cylinder.cc

Bruno Turcksin

unread,
Sep 18, 2017, 10:40:41 AM9/18/17
to dea...@googlegroups.com
2017-09-18 10:32 GMT-04:00 Lucas Campos <rmk...@gmail.com>:
Also, here is the new program, with the added if and prepare_coarsening_and_refinement(). You will find the interesting part between lines 85 and 104.

You want to use prepare_coarsening_and_refinement after you have set the flags not before.

On Monday, 18 September 2017 16:29:09 UTC+2, Lucas Campos wrote:
Dear Bruno,

You will find attached the resulting grids. While I originally expected to find new divisions only along the x-direction, your previous answer tells me it actually cuts in some local coordinate (whatever that might be). In this case, what is the best way to refine the cylinder in a single direction? I tried adding an if(cell->boundary_id() == 0) before setting the flag, but it did not really help. The mesh is still being refined in all directions.
Have you looked at step-30? 

Best,

Bruno

luca.heltai

unread,
Sep 19, 2017, 5:00:57 AM9/19/17
to Deal.II Users
Dear Lucas,

you could use GridTools::remove_anisotropy

Best,
Luca.
> --
> 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+un...@googlegroups.com.
> For more options, visit https://groups.google.com/d/optout.
> <grid_unrefined.png><grid_refined.png>

Timo Heister

unread,
Sep 19, 2017, 1:48:29 PM9/19/17
to dea...@googlegroups.com
If you are just interested in a cylinder with more than a single
layer, it might be easier to make k copies of the cylinder, shift
them, and use merge_triangulations() to combine them into a single
mesh.

On Tue, Sep 19, 2017 at 5:00 AM, luca.heltai <luca....@gmail.com> wrote:
> Dear Lucas,
>
> you could use GridTools::remove_anisotropy
>
> Best,
> Luca.
>
>> On 18 Sep 2017, at 16:29, Lucas Campos <rmk...@gmail.com> wrote:
>>
>> Dear Bruno,
>>
>> You will find attached the resulting grids. While I originally expected to find new divisions only along the x-direction, your previous answer tells me it actually cuts in some local coordinate (whatever that might be). In this case, what is the best way to refine the cylinder in a single direction? I tried adding an if(cell->boundary_id() == 0) before setting the flag, but it did not really help. The mesh is still being refined in all directions.
>>
>> I ask this because I plan to do some transformations on this cylinder, to turn it into a coil. In this case, the final result is much better if there are enough cuts in the x-axis.
>>
>> Bests
>>
>> On Monday, 18 September 2017 14:44:40 UTC+2, Bruno Turcksin wrote:
>> Lucas,
>>
>> On Monday, September 18, 2017 at 3:50:20 AM UTC-4, Lucas Campos wrote:
>> I am trying to refine a cylinder in the x direction only, but it seems that when I use
>>
>> > cell->set_refine_flag(RefinementCase<3>::cut_x);
>>
>> the whole cylinder is refined. Indeed, if I run the above line for every active cell *once* I would expect the number of cells to double. However, they quintuple!
>> It looks good but how does the mesh look like? Here, you can see that cut_x is in the local coordinate system not the global one. So I would not be surprise if the mesh doesn't look like you think it should. Also if you do more than one refinement you should also use prepare_coarsening_and_refinement() before you refine your mesh.
>>
>> Best,
>>
>> Bruno
>>
>> --
>> The deal.II project is located at https://urldefense.proofpoint.com/v2/url?u=http-3A__www.dealii.org_&d=DwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=R5lvg9JC99XvuTgScgbY_QFS80R7PEA2q0EPwDy7VQw&m=Mtc_8RpRiOkRLKkBxN-XmY8kb7kcarfzigGzEoOBiPM&s=tE77erYniO6ZJQX0HXjseStuQ1pFyTwrB95kiqc4TyM&e=
>> For mailing list/forum options, see https://urldefense.proofpoint.com/v2/url?u=https-3A__groups.google.com_d_forum_dealii-3Fhl-3Den&d=DwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=R5lvg9JC99XvuTgScgbY_QFS80R7PEA2q0EPwDy7VQw&m=Mtc_8RpRiOkRLKkBxN-XmY8kb7kcarfzigGzEoOBiPM&s=3YruSttO5wH-jhuY8FkmIFxomIEok3SQ6VQJiVGSpHU&e=
>> ---
>> 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+un...@googlegroups.com.
>> For more options, visit https://urldefense.proofpoint.com/v2/url?u=https-3A__groups.google.com_d_optout&d=DwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=R5lvg9JC99XvuTgScgbY_QFS80R7PEA2q0EPwDy7VQw&m=Mtc_8RpRiOkRLKkBxN-XmY8kb7kcarfzigGzEoOBiPM&s=xt3MRzQwT7XF4s7ezHiQC0V6qwZtrYvic1H2vsymFZI&e= .
>> <grid_unrefined.png><grid_refined.png>
>
> --
> The deal.II project is located at https://urldefense.proofpoint.com/v2/url?u=http-3A__www.dealii.org_&d=DwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=R5lvg9JC99XvuTgScgbY_QFS80R7PEA2q0EPwDy7VQw&m=Mtc_8RpRiOkRLKkBxN-XmY8kb7kcarfzigGzEoOBiPM&s=tE77erYniO6ZJQX0HXjseStuQ1pFyTwrB95kiqc4TyM&e=
> For mailing list/forum options, see https://urldefense.proofpoint.com/v2/url?u=https-3A__groups.google.com_d_forum_dealii-3Fhl-3Den&d=DwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=R5lvg9JC99XvuTgScgbY_QFS80R7PEA2q0EPwDy7VQw&m=Mtc_8RpRiOkRLKkBxN-XmY8kb7kcarfzigGzEoOBiPM&s=3YruSttO5wH-jhuY8FkmIFxomIEok3SQ6VQJiVGSpHU&e=
> ---
> 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+un...@googlegroups.com.
> For more options, visit https://urldefense.proofpoint.com/v2/url?u=https-3A__groups.google.com_d_optout&d=DwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=R5lvg9JC99XvuTgScgbY_QFS80R7PEA2q0EPwDy7VQw&m=Mtc_8RpRiOkRLKkBxN-XmY8kb7kcarfzigGzEoOBiPM&s=xt3MRzQwT7XF4s7ezHiQC0V6qwZtrYvic1H2vsymFZI&e= .



--
Timo Heister
http://www.math.clemson.edu/~heister/

Timo Heister

unread,
Sep 19, 2017, 1:49:16 PM9/19/17
to dea...@googlegroups.com
Oh, even easier would be creating a 2d cylinder and then using
GridGenerator::extrude_triangulation()
Reply all
Reply to author
Forward
0 new messages