Here is perhaps a minimaler (non)working example:
#include <deal.II/grid/tria.h> #include <deal.II/grid/grid_generator.h> #
include <deal.II/grid/grid_refinement.h> #include 
<deal.II/dofs/dof_handler.h> #include <deal.II/dofs/dof_tools.h> #include 
<deal.II/fe/fe_dgq.h> #include <deal.II/lac/vector.h> #include 
<deal.II/numerics/vector_tools.h> using namespace dealii; int main() { 
FE_DGQ<2> fe(0); Triangulation<2> trig; Triangulation<2> trig_child; 
GridGenerator::hyper_cube(trig); trig.refine_global(1); 
trig_child.copy_triangulation(trig); trig_child.refine_global(1); 
DoFHandler<2> dof_handler(trig); dof_handler.distribute_dofs(fe); Vector<int> 
ref_level({1,2,3,4}); // ref_level is the vector on // the coarser mesh trig 
DoFHandler<2> dof_handler_child(trig_child); 
dof_handler_child.distribute_dofs(fe); Vector<int> 
ref_level_child(dof_handler_child.n_dofs()); // ref_level_child is the 
vector // on the finer mesh trig_child 
VectorTools::interpolate_to_different_mesh(dof_handler, ref_level, 
dof_handler_child, ref_level_child); } 

And I just found that the bug comes from feeding Vector<int> to
interpolate_to_different_mesh() while int is not a valid
VectorType::value_type therein. Changing int to double fixes it.

Best,
Yuan
​
On Friday, November 10, 2023 at 7:48:11 PM UTC Yuan Wang wrote:

> Hi Peter,
>
> Thanks a lot for your help! I was trying to implement a MWE with the first 
> method you suggested. However, the function 
> VectorTools::interpolate_to_different_mesh() seems to be causing a linker 
> error:
>
> /home/usrname/Projects/hdg_dealii/adr/ip/multilevel_refinement/MultilevelRefinement.cc:
> 66: error: undefined reference to 'void 
> dealii::VectorTools::interpolate_to_different_mesh<2, 2, 
> dealii::Vector<int> >(dealii::DoFHandler<2, 2> const&, dealii::Vector<int> 
> const&, dealii::DoFHandler<2, 2> const&, dealii::Vector<int>&)' collect2: 
> error: ld returned 1 exit status make[2]: *** 
> [CMakeFiles/MultilevelRefinement.dir/build.make:133: 
> MultilevelRefinement] Error 1 make[1]: *** [CMakeFiles/Makefile2:90: 
> CMakeFiles/MultilevelRefinement.dir/all] Error 2 make: *** [Makefile:91: 
> all] Error 2 
>
> Here is, I think, the relevant chunk of code (original code file is also 
> attached):
> // ... #include <deal.II/numerics/vector_tools.h> // ... 
> Triangulation<dim> triangulation_old; 
> triangulation_old.copy_triangulation(triangulation); DoFHandler<dim> 
> dof_handler_old(triangulation_old); dof_handler_old.distribute_dofs(fe); 
> Vector<int> ref_level_old(ref_level); 
> triangulation.execute_coarsening_and_refinement(); 
> dof_handler.reinit(triangulation); dof_handler.distribute_dofs(fe); 
> ref_level.reinit(dof_handler.n_dofs()); 
> VectorTools::interpolate_to_different_mesh(dof_handler_old, ref_level_old, 
> dof_handler, ref_level); // ... 
>
> Can you help shed some light on this?
>
> Best regards,
>
> Yuan
> ​
> On Thursday, November 9, 2023 at 6:25:45 AM UTC [email protected] wrote:
>
>> Hi Greg,
>>
>> I guess what you could do is to treat the vector as vector associated 
>> with a DoFHandler with FE_DGQ(0). During each refinement you "interpolate" 
>> the vector to the new mesh and decrease the value by one. This approach 
>> naturally extends to a parallel setting.
>>
>> Alternately, you could create a single vector associated the coarse cells 
>> and while looping over active cells of the refined mesh you would 
>> recursively call the parent cell until you reach the coarsest cell whose id 
>> you can use to access the vector. In the parallel setting, the vector 
>> associated to the coarse cells would need to replicated between processes.
>>
>> Hope this help,
>> PM
>>
>> On Thursday, 9 November 2023 at 03:32:42 UTC+1 [email protected] wrote:
>>
>>> Dear all,
>>>
>>> This might be a somewhat odd feature request but I was wondering whether 
>>> the functionality described below can be readily achieved in deal.II 
>>> already:
>>>
>>> Say we are given a coarse mesh with every cell never refined, and a 
>>> vector of non-negative integers associated with each cell. The integer is 
>>> the level of isotropic refinement to be applied to each cell. The number 
>>> can be 0, which means not doing anything here; it can be 1, which 
>>> corresponds to what set_refine_flag() does. But it can also be any 
>>> number bigger than 1. Say we have a quadrilateral cell and the number 
>>> associated with it is 2, then the cell is to be refined isotropically once 
>>> into four children and then every child is refined again, resulting in 16 
>>> grandchildren.
>>>
>>> If we were to hard code this, would it be sound to simply rely on 
>>> set_refine_flag() and execute the refinement one level after another 
>>> for M levels, with M being the maximum number of the given vector? And 
>>> perhaps to do this we would need iterators triangulation.begin(n), 
>>> triangulation.end(n) and isotropic_child() involved somehow?
>>>
>>> Any help and/or pro tips would be greatly appreciated!
>>>
>>> Best regards,
>>> Yuan
>>> ​
>>>
>>

-- 
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/4733aa4f-63b5-4a37-892a-cfa58b9e98aen%40googlegroups.com.

Reply via email to