Dear community

in case of use, I attempted to write a small code that seems to solve the 
issue of linking manifold cells to volume cells. 

As first, two maps are defined:

  *// maps*

  *// the map surface_to_volumefaces_mapping contains the outcome of the 
method extract_boundary_mesh and connects*

  *// panels from the manifold triangulation to the faces on the boundary 
of the volume triangulation.*

  *// the map manifoldcells_to_volumecells_mapping pairs surface panels 
from the manifold triangulation*

  *// to the cells of the volume triangulation. It is used in the method 
pair_the_triangulations().*

  *//*


  std::map<

            *typename* parallel::shared::Triangulation<manifold_dim
,dim>::cell_iterator,

            *typename* parallel::shared::Triangulation<dim>::face_iterator

          >

          surface_to_volumefaces_mapping;


  std::map<

            *typename* parallel::shared::Triangulation<manifold_dim
,dim>::cell_iterator,

            *typename* parallel::shared::Triangulation<dim>::cell_iterator

          >

          manifoldcells_to_volumecells_mapping;




then the manifold grid is generated


    *//*

    *// Surface triangulation definition*

    *//*

    

    *this*->pcout << "\n   Generation of surface discretization ... " << std
::flush;

    

    surface_to_volumefaces_mapping = GridGenerator::extract_boundary_mesh ( 
*this*->triangulation, *this*->manifold_triangulation );

    *this*->manifold_triangulation.set_all_manifold_ids(0);

    

    *this*->manifold_dof_handler.distribute_dofs ( *this*->manifold_fe );

    *this*->manifold_accumulated_solution.reinit ( 
*this*->manifold_dof_handler.n_dofs() 
);

    *this*->manifold_former_step_solution.reinit ( 
*this*->manifold_dof_handler.n_dofs() 
);

    *this*->manifold_initial_guess.reinit ( 
*this*->manifold_dof_handler.n_dofs() 
);

    

    *this*->pcout << " done \n" << std::flush;

    

    *//*

    *// Surface triangulation and volume triangulation pairing*

    *//*

    

    pair_the_triangulations();


In the method pair_the_triangulations() a manifold cell is linked  to the 
corresponding volume cell for parallel::shared::triangulations as follows.

pair_the_triangulations ( *bool* write_map_report )


*//*

*// this function associates each*

*// surface panel to the volume cell from which it emanates.*

*// the map manifoldcells_to_volumecells_mapping contains the pair 
surface/volume panels*

*//*


{

  *// unsigned face_counter = 0;*

  

  *this*->pcout << "pairing the two triangulations ... " << std::flush;

  

  *// loop through the volumetric triangulation*

  

  *typename* parallel::shared::Triangulation<dim>::active_cell_iterator

  cell = *this*->triangulation.begin_active (),

  endc = *this*->triangulation.end();

  

  *for* (; cell!=endc; ++cell)

    

  *// We do not seek for locally_owned on the volumetric triangulation*

    

  {

    

    *// debug*

    *// this->pcout << " new cell " << cell->id() << " " << std::flush;*

    

    *// select a face at boundary*

    

    *for* (*unsigned* *int* face_number=0; 
face_number<GeometryInfo<dim>::faces_per_cell; 
++face_number)

    {

      *const* *typename* parallel::shared::Triangulation<dim>::face_iterator

      face = cell->face( face_number );

      

      *if* ( face->at_boundary() )

      {

        

        *// scan the surface_to_volumefaces_mapping and*

        *// associate maifold panels to volume panels*

        

        *typename*

        std::map< *typename* parallel::shared::Triangulation<manifold_dim
,dim>::cell_iterator,

        *typename* parallel::shared::Triangulation<dim>::face_iterator >

        ::iterator

        map_it = surface_to_volumefaces_mapping.begin(),

        map_end = surface_to_volumefaces_mapping.end() ;

        

        *for* (; map_it!=map_end; ++map_it)

          *if* ( face == map_it->second )

          {

            

            *// debug*

            *// this->pcout <<  cell->face_index(face_number) << ", " << 
std::flush;*

            *// face_counter++;*

            

            manifoldcells_to_volumecells_mapping[ map_it->first ] = cell ;

          }

      }

    }

  }

  

  *// debug*

  *// this->pcout <<  "\n I found " << face_counter << " faces. \n bye. \n 
" << std::flush;*

  

  *// map report*

  

  *if* ( write_map_report )

  {

    

    *typename* std::map< *typename* parallel::shared::Triangulation<
manifold_dim,dim>::cell_iterator,

    *typename* parallel::shared::Triangulation<dim>::cell_iterator >

    ::iterator

    manifoldcells_to_volumecells_it = 
manifoldcells_to_volumecells_mapping.begin() 
,

    manifoldcells_to_volumecells_end = 
manifoldcells_to_volumecells_mapping.end() 
;

    

    *this*->pcout <<  "\n map report: "  << std::flush;

    *for* (; 
manifoldcells_to_volumecells_it!=manifoldcells_to_volumecells_end; 
++manifoldcells_to_volumecells_it)

      *this*->pcout <<  "[ " << manifoldcells_to_volumecells_it->first->*id*() 
<< ", " << manifoldcells_to_volumecells_it->second->*id*() << " ], " << std
::flush;

    *this*->pcout <<  "\n done. \n"  << std::flush;

    

  }

  

  *// debug*

  *// std::exit( 0 );*

  

}

 
Comments and suggestions are very welcome.
Thank you in advance

Alberto


Il giorno sabato 7 settembre 2019 12:09:30 UTC+2, Alberto Salvadori ha 
scritto:
>
> Dear community
>
> I apologize if the question was too generic. I attempt here to be more 
> precise, asking for some suggestions.
>
> In solving a Laplace-Beltrami problem on an advecting surface one should 
> identify the Gauss points on the manifold dim-1 cell and retrieve at such 
> locations relevant information from the solution of the advection problem, 
> using the dim dof_handler of the volume. I wonder how to connect a manifold 
> dim-1 cell to the volumetric cell it was extracted from. The 
> GridGenerator::extract_boundary_mesh method seem to provide some 
> information, since "it returns A map that for each cell of the surface mesh 
> (key) returns an iterator to the corresponding face of a cell of the volume 
> mesh (value). " . I am not sure how I can use this map to link a manifold 
> cell to the corresponding volume cell. 
>
> Any help is very much appreciated. Thanks!
>
> Alberto
>
>
>
> Il giorno martedì 3 settembre 2019 18:36:59 UTC+2, Alberto Salvadori ha 
> scritto:
>>
>> Dear community
>>
>> I am trying to find a good way to deal with the following problem.
>>
>> I am interested in a Laplace-Beltrami simulation (say with unknown scalar 
>> field c) on a surface that advects (with displacement field u) enclosesing 
>> a volume of say elastic material.
>> The two problems are per se uncoupled, but for the geometrical evolution. 
>> In this regard, one may solve separately for displacements in the volume 
>> and for c on the surface, once the latter is informed on the displacement 
>> field (and related gradients) on the surface nodes (and Gauss points). 
>> Extracting the surface mesh from the volume is quite simple. Which is 
>> though the best way to extract (project) the displacement solution on the 
>> boundary? Any suggestion/ hint?
>>
>> Many thanks,
>> Alberto
>>
>
-- 


Informativa sulla Privacy: http://www.unibs.it/node/8155 
<http://www.unibs.it/node/8155>

-- 
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.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/57f2de0e-b449-4cbd-9332-3c66a1b078e7%40googlegroups.com.

Reply via email to