Thank you so much. I will have a look!

On Wednesday, August 1, 2018 at 1:58:27 PM UTC+1, Jean-Paul Pelteret wrote:
>
> Dear Jane,
>
> The method described by Claire is exactly the approach used by my 
> colleagues to model additive manufacturing processes. See, for example, 
> this paper 
> 10.1007/s00170-016-8819-6 <https://dx.doi.org/10.1007/s00170-016-8819-6> and 
> maybe others by D. Riedlbauer 
> <https://www.ltm.tf.fau.eu/chair-of-applied-mechanics/publications/> that 
> perhaps explain the approach in more detail.
>
> Best,
> Jean-Paul
>
> On 01 Aug 2018, at 14:20, Claire <claire....@gmail.com <javascript:>> 
> wrote:
>
> Dear Jane,
>
> If I got well what you want to do I would suggest I a different approach, 
> namely an "element birth" technique.
> It can be implemented in deal.ii using FE_Nothing elements and the hp:: 
> namespace
>
> https://dealii.org/9.0.0/doxygen/deal.II/classFE__Nothing.html
> https://dealii.org/9.0.0/doxygen/deal.II/namespacehp.html
>
> More practically, instead of adding cells to your triangulation, create 
> from the beginning the final domain you want to simulate.
> Then you set the fe_index of the cells that belong to the domain you don't 
> want to solve yet the one of FE_Nothing in your FECollection.
> This finite element has no degree of freedom, so nothing will happen there.
>
> When you want to activate one cell, just turn its fe_index to the one 
> corresponding to "FE_Something" in your FECollection.
>
> I hope my brief explanation was clear enough, and that this suggestion is 
> relevant for your problem.
> If you ever need further information, don't hesitate to ask.
>
> Best,
>
> Claire
>
> Le lundi 30 juillet 2018 12:39:34 UTC+2, Jane Lee a écrit :
>>
>>
>> Hi all,
>>
>> I'm a little stuck on a problem and hope that someone can suggest 
>> something. 
>>
>> I am solving a problem where I need to add a thin layer onto a domain 
>> which I move after solving each solution;
>> 1. I solve for the velocity
>> 2. I then move the mesh by using a vertex displacement method like in 
>> step-18 (it works fine til up here)
>> 3. I create a domain using the top boundary vertices of the resulting 
>> mesh in 2. and I want to add it to the same mesh and solve again. 
>> (The background is that I want to deposit a layer of material on top of 
>> the existing one to represent sedimentation).
>>
>> I have tried using both create_union_triangulation and 
>> merge_triangulations and problem is as follows:
>> I am using a subdivided hyper rectangle as the original mesh then 
>> refining before I do 1. of above. This means that when i create the domain 
>> in 3. with the vertex points of the top boundary, the created mesh has a 
>> different coarse grid than the first, so I am getting the error that they 
>> are not derived from the coarse mesh. 
>>
>> Can anyone suggest a way to coarsen both meshes so that the merging 
>> works? I also seem to have the problem that my original grid doesn't 
>> coarsen all the way back to the initial subdivided rectangle and not sure 
>> why.....
>>
>> And after that, would i have to interpolate the solution on the old mesh 
>> back first then some other solution on the created mesh (which I know) 
>> BEFORE i merge? 
>>
>> Many thanks
>>
>> Sample codes at each step:
>>
>> 2. move mesh
>>
>> for (typename DoFHandler<dim>::active_cell_iterator
>>
>>              cell = dof_handler.begin_active ();
>>
>>              cell != dof_handler.end(); ++cell)
>>
>>             for (unsigned int v=0; 
>> v<GeometryInfo<dim>::vertices_per_cell; ++v)
>>
>>                 if (vertex_touched[cell->vertex_index(v)] == false)
>>
>>                 {
>>
>>                     vertex_touched[cell->vertex_index(v)] = true;
>>
>>                     Point<dim> vertex_displacement;
>>
>>                     for (unsigned int d=0; d<dim; ++d)
>>
>>                     {vertex_displacement[d]
>>
>>                         = solution.block(0)(cell->vertex_dof_index(v,d));
>>
>>                     }
>>
>>                     cell->vertex(v) += vertex_displacement*timestep_size;
>>
>>                     //                    std::cout << cell->diameter() 
>> << std::endl;
>>
>>                 }
>>
>>
>> 3. creating top domain (then the last bit doesn't work)
>>
>>         print_mesh_info(triangulation, "after_solving_grid.eps");
>>
>>
>>         std::vector<Point<2>> vertex_points;
>>
>>         vertex_points.push_back (Point<dim> (0,1.));
>>
>>         unsigned int p = 0;
>>
>>
>>         for (typename Triangulation<dim>::active_cell_iterator
>>
>>              cell = triangulation.begin_active();
>>
>>              cell != triangulation.end(); ++cell)
>>
>>         {
>>
>>             for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; 
>> ++f)
>>
>>             {
>>
>>                 if (cell->face(f)->boundary_id() == 1)
>>
>>                 {
>>
>>                     Point<2> points = cell->face(f)->vertex(1);
>>
>> //                    std::cout << points << std::endl;
>>
>>                     vertex_points.push_back (points);
>>
>>                     ++p;
>>
>>                 }
>>
>>             }
>>
>>         }
>>         
>>
>>         const double end_z_point = vertex_points[vertex_points.size()-1][
>> 1];
>>
>>
>>         // Output vertices of top of old domain.
>>
>>         for (unsigned int i=0; i<=vertex_points.size()-1; ++i)
>>
>>             std::cout << vertex_points[i] << std::endl;
>>
>>
>>         // Number of vertices along boundary
>>
>>         const int no_vertex_bdry = vertex_points.size();
>>
>>         std::cout << "Number of vertex points: " << no_vertex_bdry 
>> <<std::endl;
>>         
>>
>>         for (unsigned int i=0; i<no_vertex_bdry; ++i)
>>
>>         {
>>
>>             vertex_points.push_back (Point<dim> (vertex_points[i][0], 
>> vertex_points[i][1]+sed_height));
>>
>>         }
>>
>>
>>         // output vertex points of new deposited layer
>>
>>         const int no_vertex_total = no_vertex_bdry*2;
>>
>>         for (unsigned int i=0; i<no_vertex_total; ++i)
>>
>>             std::cout << vertex_points[i] << std::endl;
>>
>>
>> //
>>
>> //        // numbering will go from 0 to vertex_length-1
>>         
>>
>>         const int n = no_vertex_bdry-1;
>>
>>         const int m = GeometryInfo<dim>::vertices_per_cell;
>>
>>         unsigned int cell_vertices[n][m];// = {0};
>>         
>>
>>         for (unsigned int i=0; i<n; ++i)
>>
>>         {
>>
>>             cell_vertices[i][0] = i;
>>
>>             cell_vertices[i][1] = i+1;
>>
>>             cell_vertices[i][2] = no_vertex_bdry + i;
>>
>>             cell_vertices[i][3] = no_vertex_bdry + i+1;
>>
>>             std::cout << cell_vertices[i][3] << std::endl;
>>
>>         }
>>
>>
>>         const unsigned int
>>
>>         n_cells = sizeof(cell_vertices) / sizeof(cell_vertices[0]);
>>
>>         std::vector<CellData<dim> > cells (n_cells, CellData<dim>());
>>
>>         for (unsigned int i=0; i<n_cells; ++i)
>>
>>         {
>>
>>             for (unsigned int j=0;
>>
>>                  j<GeometryInfo<dim>::vertices_per_cell;
>>
>>                  ++j)
>>
>>                 cells[i].vertices[j] = cell_vertices[i][j];
>>
>>             cells[i].material_id = 0;
>>
>>         }
>>
>>
>>         sediment_tria.create_triangulation (vertex_points, cells, 
>> SubCellData());
>>
>>
>>         print_mesh_info(sediment_tria, "sediment_grid_refined.eps");
>>         
>>
>>             for (typename Triangulation<dim>::active_cell_iterator
>>
>>                  cell = sediment_tria.begin_active();
>>
>>                  cell != sediment_tria.end(); ++cell)
>>
>>             {
>>
>>                 cell->set_coarsen_flag();
>>
>>             }
>>     
>>
>>             sediment_tria.execute_coarsening_and_refinement ();
>>     
>>
>>             print_mesh_info(sediment_tria, "sediment_grid_coarse.eps");
>>
>>         //
>>
>> //        for (typename Triangulation<dim>::active_cell_iterator
>>
>> //             cell = triangulation.begin_active();
>>
>> //             cell != triangulation.end(); ++cell)
>>
>> //        {
>>
>> //            cell->set_coarsen_flag();
>>
>> //        }
>>
>> //
>>
>> //        triangulation.execute_coarsening_and_refinement ();
>>
>> //        print_mesh_info(triangulation, "initial_grid.eps");
>>
>> //
>>
>> //        GridGenerator::create_union_triangulation (sediment_tria, 
>> triangulation, final_tria);
>>
>> //                                             print_mesh_info(final_tria, 
>> "final_grid.eps");
>>         
>>
>
> -- 
> 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 <javascript:>.
> For more options, visit https://groups.google.com/d/optout.
>
>
>

-- 
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.
For more options, visit https://groups.google.com/d/optout.

Reply via email to