HI Karthik,

Glad we could help :-)

To Question 1:
So you estimate the error for the second component of each of your 
solutions, and then you add up the estimates on each cell. What is the 
reasoning behind summing up the error of multiple solutions?
Assume that the error estimate of one of your solutions is larger than all 
others by magnitudes, then this would be the dominant part of your estimate 
sum, and the other solutions would have no influence at all.
Maybe it would be reasonable to only pick one of your estimates as a 
measure for grid refinement.
First, you should find a detailed answer to the question why you want to 
refine your grid, and then find a way how you want to do it.

To Question 2:
`GridRefinement::refine_and_coarsen_fixed_number` always refines and 
coarsens a the specified fraction of ALL cells. With this you can control 
the growth of the mesh size, rather than the reduction of the error. For 
the latter, you can use `GridRefinement::refine_and_coarsen_fixed_fraction`. 
Otherwise, I would suggest to revise the fractions you are using for grid 
refinement.

Marc

On Wednesday, January 27, 2021 at 12:40:37 PM UTC-7 Karthi wrote:

> Hi Marc,
>
>
> Sorry for my delayed response, I was away for a couple of days. Thank you, 
> your suggestions were very helpful and it worked.
>
> I have a couple of follow up questions in regard to mesh refinement.
>
>
> Question (1)
>
>
> I employ FE_Q<dim>(degree, 2) for all my sub-problems (and the same 
> dof_handler).
>
>
> Hence, I have a std::vector<solution>, I created a container of 
> estimated_error_per_cell 
> for each solution and summed it all up in sum_estimater_error_per_cell, 
> which is then used in the Kelly Error Estimator. I used fe.componet_mask to 
> account for only the second component of my solution while calculating the 
> error estimate. 
>
>
> My objective is calculate a gradient error estimator as follow, where N 
> represents solution.size() (i.e the number of sub-problems)
>
>
> \sum_{i}^{N} |\nabla \alpha_i|
>
>
>
> I don't know if the below code achieves this?
>
>
> std::vector<Vector<float>> estimated_error_per_cell(num_index,Vector<float
> >(triangulation.n_active_cells()));
>
>
> Vector<float> sum_estimated_error_per_cell(triangulation.n_active_cells
> ());
>
> FEValuesExtractors::Scalar alpha(1);
>
>
> for(unsigned int i=0; i < num_index; i++)
>
>    KellyErrorEstimator<dim>::estimate(dof_handler,
>
> QGauss<dim-1>(fe.degree + 1),
>
> {},
>
> solution[i],
>
> estimated_error_per_cell[i],
>
> fe.component_mask(alpha));
>
>
>  for(unsigned int i=0; i < num_index; i++)
>
>     sum_estimated_error_per_cell += estimated_error_per_cell[i];
>
>
>
>  GridRefinement::refine_and_coarsen_fixed_number(triangulation,
>
>                                                   
> sum_estimated_error_per_cell,
>
>                                                     0.60,
>
>                                                     0.40);
>
>
> Question (2)
>
>
> Please see the attachment of mesh refinement plots of the transient 
> solution. As you can see, the initial mesh refinement at time zero makes 
> sense but as the solution decays I was hoping the mesh would coarsen (but 
> it refines further). I am clearly doing something wrong. I need some help 
> in fixing this issue.
>
>
> Thank you!
>
>
> Karthi.
>
> On Mon, Jan 25, 2021 at 12:19 AM Marc Fehling <mafe...@gmail.com> wrote:
>
>> Hi Karthi,
>>
>> if you work on the same DoFHandler, one SolutionTransfer object is 
>> sufficient.
>>
>> There are already member functions that take a container of solutions 
>> like yours as a parameter. Have a look at 
>> SolutionTransfer::prepare_for_coarsening_and_refinement 
>> <https://www.dealii.org/developer/doxygen/deal.II/classSolutionTransfer.html#ae6dc5e5a74b166b0dea35f5a64694e69>
>>  
>> and SolutionTransfer::interpolate 
>> <https://www.dealii.org/developer/doxygen/deal.II/classSolutionTransfer.html#ae067f9b520ed50c86a9ff4c7776d16cb>
>> .
>>
>> Best,
>> Marc
>>
>> On Saturday, January 23, 2021 at 6:30:34 AM UTC-7 Karthi wrote:
>>
>>> Dear All,
>>>
>>> I have a fourth order parabolic equation, which is split into two second 
>>> order equations. Hence I solve using components by 
>>> declaring FESystem<dim>  fe(FE_Q<dim>(degree,2). I have in total three such 
>>> sub-problems, which are coupled to each other in a semi-implicit manner. 
>>> Therefore I have a std::vector of solution for the entire system.
>>>
>>> std::vector<Vector<double>> solution.
>>>
>>> In addition, I am employing mesh adaptivity. After estimating the error 
>>> using Kelly, I would like to perform a solution transfer from old to new 
>>> mesh. Do I need to create a std::vector of SolutionTransfer objects; one 
>>> for each solution?
>>>
>>> The below code copied from step-26 seems to work. Is this the correct 
>>> approach?
>>>
>>> std::vector<SolutionTransfer<dim>> solution_trans(3,dof_handler);
>>>
>>> std::vector<Vector<double>> previous_solution(num_index);
>>>
>>> for(unsigned int i=0; i < num_index; i++){previous_solution[i] = 
>>> solution[i];}
>>>
>>> triangulation.prepare_coarsening_and_refinement();
>>>
>>> for(unsigned int i=0; i < num_index; 
>>> i++){solution_trans[i].prepare_for_coarsening_and_refinement(previous_solution[i]);}
>>>
>>> triangulation.execute_coarsening_and_refinement();
>>>
>>> setup_system();
>>>
>>> for(unsigned int i=0; i < num_index; i++){
>>>
>>> solution_trans[i].interpolate(previous_solution[i], solution[i]);
>>>
>>> constraints.distribute(solution[i]);}
>>>
>>> I look forward to your response. 
>>>
>>> Best regards,
>>>
>>> Karthi.
>>>
>> -- 
>> 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/YA-eLEs43bc/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/26b50c8f-90ec-4b34-9f64-2dae588a1cdfn%40googlegroups.com
>>  
>> <https://groups.google.com/d/msgid/dealii/26b50c8f-90ec-4b34-9f64-2dae588a1cdfn%40googlegroups.com?utm_medium=email&utm_source=footer>
>> .
>>
>

-- 
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/a3054b82-b1bf-4b7f-86bf-3d07be067fc6n%40googlegroups.com.

Reply via email to