Selectively coarsen after a refinement step

39 views
Skip to first unread message

Kaushik Das

unread,
Nov 30, 2020, 10:33:16 AM11/30/20
to deal.II User Group
Hi:
I am new to Deal.ii and reading through the examples to learn more. But I have a question, that I did not find a suitable answer in the example problems. 

I am trying to solve a transient heat conduction problem over a domain due to laser heating in an additive manufacturing process. At each time step, I want to refine the mesh at the laser spot and then coarsen the end of the time step and then in the next time step laser spot moves and I can again refine the mesh where the laser is then. I know the position of the laser at all times.

Something like this:

Triangulation<3> triangulation;
GridGenerator::subdivided_hyper_rectangle(triangulation, repetitions, left, right, true);

// Refine the mesh based on some geometric criteria, like material interfaces, near edges, etc.

Loop over time step starts

// Refine at the laser spot

// do calculations. 

// coarsen to get back to the state where the mesh was at the beginning of the time step. 

Loop ends

I tried two different ways. 

1. At the end of the time step, call set_coarsen_flag on those cells that were marked with set_refine_flag at the beginning of the step.  But, cells that were refined due to smoothing are left refined which is not desirable.
2. Add a pre_refinement signal and then save the save_refine_flags during refinement and then load_coarsen_flags that vector<bool> from the save_refine_flags call. But here I am getting assert failure on mismatched size of that vector and the active number of cells whining calling execute_coarsening_and_refinement to do the coarsening. 

I also thought of using the persistenttrainregulation class, but I am not sure if that covers my use-case. 

I may be able to write some geometric criteria to figure out which cells need coarsening. But I am hoping there is a faster method to do this because there must be some way to know which cells got refined and by how much and just coarsen them again. 

Thank you very much.
-Kaushik 


 

Wolfgang Bangerth

unread,
Nov 30, 2020, 11:27:02 AM11/30/20
to dea...@googlegroups.com
On 11/30/20 8:33 AM, Kaushik Das wrote:
>
> I may be able to write some geometric criteria to figure out which cells need
> coarsening. But I am hoping there is a faster method to do this because there
> must be some way to know which cells got refined and by how much and just
> coarsen them again.

Not really. You'd have to save this information yourself somehow. The way to
do that would be like you do it now, but *after* you call
triangulation.prepare_coarsening_and_refinement() because then also those
cells that will be refined for mesh smoothing purposes will be correctly flagged.

But in the end, you should refine based on where the solution allows it, not
because the laser spot has moved on. Do it automatically, in a way like we do
in step-26, for example.

Best
W.

--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/

Kaushik Das

unread,
Nov 30, 2020, 2:37:33 PM11/30/20
to dea...@googlegroups.com
Dear Prof. Bangerth,
Thank you very much for your quick response. 

>The way to do that would be like you do it now, but *after* you call
>triangulation.prepare_coarsening_and_refinement() because then also those
>cells that will be refined for mesh smoothing purposes will be correctly flagged.

I thought the pre_refinement signal is called after prepare_coarsening_and_refinement call in execute_coarsening_and_refinement. That's why I tried: 

template <int dim>
class saveRefinement
{
...
  std::vector<bool> cellRefined;
...
  void operator()()
  {
    triangulation.save_refine_flags(cellRefined);
  } 
}
...
saveRefinement<3> SaveRefinement(triangulation);
triangulation.signals.pre_refinement.connect([&SaveRefinement] { SaveRefinement(); })

/// refine the mesh
....
triangulation.load_coarsen_flags(SaveRefinement.cellRefined);
triangulation.execute_coarsening_and_refinement();

At this point, I am hitting Assert(v.size() == n_active_cells(), ExcGridReadError())   in  load_coarsen_flags.

Which also makes sense to me now, because those flags saved from  save_refine_flags are now parent cells, and I have to mark their child cells for coarsening. 
I will greatly appreciate any hint on how that can be done. 

I looked into the Step-26 example. There is another reason I wanted to refine the mesh based on geometric criteria. In addition to the laser spot, I would like to refine the mesh near the layer that is being printed at the given time step. And coarsen the layers underneath it. I am trying to simulation a layer upon layer printing process following a paper by Claire Bruna-Rosso et al., who also used deal.ii. They used the FE_Nothing element to exclude the layers that are yet to be printed from the solution. 

Thank you,
Kaushik 

                                            

--
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 a topic in the Google Groups "deal.II User Group" group.
To unsubscribe from this topic, visit https://groups.google.com/d/topic/dealii/6hhi5Xf0_mQ/unsubscribe.
To unsubscribe from this group and all its topics, send an email to dealii+un...@googlegroups.com.
To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/b0e798b9-bf73-a3a3-d8d0-7f7d9783d274%40colostate.edu.

Wolfgang Bangerth

unread,
Nov 30, 2020, 2:52:36 PM11/30/20
to dea...@googlegroups.com
On 11/30/20 12:36 PM, Kaushik Das wrote:
> saveRefinement<3> SaveRefinement(triangulation);

Unrelated to your question, but this convention is very hard to read for most
programmers: You name a class with lower-case first letter and a variable with
upper-case first letter. Both use CamelCase. It's quite difficult to see from
the name what's a variable and what's a class.


> triangulation.signals.pre_refinement.connect([&SaveRefinement] {
> SaveRefinement(); })
>
> /// refine the mesh
> ....
> triangulation.load_coarsen_flags(SaveRefinement.cellRefined);
> triangulation.execute_coarsening_and_refinement();
>
> At this point, I am hitting Assert(v.size() == n_active_cells(),
> ExcGridReadError())   in load_coarsen_flags.
>
> Which also makes sense to me now, because those flags saved from
> save_refine_flags are now parent cells, and I have to mark their child cells
> for coarsening.

That's correct.


> I will greatly appreciate any hint on how that can be done.

You need to store a map from the old cell to the flags. You could use
std::map<Triangulation::cell_iterator,bool>
for example.


> I looked into the Step-26 example. There is another reason I wanted to refine
> the mesh based on geometric criteria. In addition to the laser spot, I would
> like to refine the mesh near the layer that is being printed at the given time
> step. And coarsen the layers underneath it. I am trying to simulation a layer
> upon layer printing process following a paper by Claire Bruna-Rosso et al.,
> who also used deal.ii. They used the FE_Nothing element to exclude the layers
> that are yet to be printed from the solution.

You can refine by your geometric argument, but it might be easiest to coarsen
based on the general indicators.

Kaushik Das

unread,
Dec 1, 2020, 1:48:25 PM12/1/20
to dea...@googlegroups.com
Thank you very much for your help. After some trial and error, I was able to do it. I have attached the code if that can be of any help to others searching for similar solutions. I was modifying step-1.cc and did not modify the file name. Thank you for your suggestion of using the error estimation for corsening, I will try that too.
Kaushik 

--
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 a topic in the Google Groups "deal.II User Group" group.
To unsubscribe from this topic, visit https://groups.google.com/d/topic/dealii/6hhi5Xf0_mQ/unsubscribe.
To unsubscribe from this group and all its topics, send an email to dealii+un...@googlegroups.com.
step-1.cc
Reply all
Reply to author
Forward
0 new messages