Dear Chucui,

Too add to all of the other comments, you can "sanity check” the result that 
you’re getting using this approach by performing the same calculations another 
way and comparing results. Specifically, you could easily calculate these norms 
by looping over all cells and integrating the quantities yourself. As 
pseudocode it would look like this:
1. Loop over all cells in DoFHandler
2. Reinit FEValues with the update_values and update_gradients flags enabled
3. Compute local solution and its gradient at the cell quadrature points using 
fe_values.get_solution_[values/gradients]()
4. Loop over all quadrature points
5. u_square +=  local_solution[q_point]*local_solution[q_point] * 
fe_values.JxW() ; grad_u_square +=  local_gradient[q_point]* 
local_gradient[q_point] * fe_values.JxW()

Maybe this approach could help you debug your issue.
Best,
Jean-Paul

> On 12 Feb 2019, at 06:49, Wolfgang Bangerth <[email protected]> wrote:
> 
> On 2/11/19 12:20 AM, [email protected] wrote:
>> Dear Prof. Bangerth:
>> 
>> Thank you very much for your quick reply!
>> 
>> When I apply hanging node constraints to my matrix, as the step-6 
>> says: https://www.dealii.org/developer/doxygen/deal.II/step_6.html
>> I make 4 steps:
>>   1.Create a Constraints Class:
>> |
>> ConstraintMatrix    constraints;
>> |
>> 
>>   2.Fill this object using the function 
>> DoFTools::make_hanging_node_constraints() to ensure continuity of the 
>> elements 
>> of the finite element space.
>> |
>> DoFTools::make_hanging_node_constraints (dof_handler,
>>                                          constraints);
>> |
>> 
>>   3.Copy the local contributions to the matrix and right hand side into the 
>> global objects:
>> |
>>     constraints.distribute_local_to_global (cell_matrix,
>>                                             cell_rhs,
>>                                             local_dof_indices,
>>                                             system_matrix,
>>                                             system_rhs);
>>     constraints.distribute_local_to_global (cell_volume_matrix,
>>                                             local_dof_indices,
>>                                             volume_matrix);
>>     constraints.distribute_local_to_global (cell_gradient_matrix,
>>                                             local_dof_indices,
>>                                             gradient_matrix);
>> |
>> 
>> 
>>   4.Make sure that the degrees of "freedom" located on hanging nodes get 
>> their correct (constrained) value:
>> |
>> constraints.distribute (solution);
>> |
> 
> Constraints are funny and sometimes require deep thought about what exactly 
> they mean. What happens if you don't apply constraints to the 
> cell_volume_matrix and cell_gradient_matrix -- i.e., you copy the elements 
> 1:1 
> into the matrix, without the 'constraints' object? (Like in step-4.)
> 
> Alternatively, you could do as you show above and after calling 
> 'constraint.distribute(solution)' you manually set all constrained degrees of 
> freedom in the solution vector to zero. That's not what you want to use for 
> visualization, error estimation, ... but it might work for the operation 
> you're trying here.
> 
> I'm pretty sure that both of these solutions should work individually (but 
> not 
> in combination). But explaining why this is so would take a page or two, I'm 
> afraid.
> 
> Best
> W.
> 
> 
> -- 
> ------------------------------------------------------------------------
> Wolfgang Bangerth          email:                 [email protected]
>                          www: http://www.math.colostate.edu/~bangerth/
> 
> -- 
> 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].
> 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 [email protected].
For more options, visit https://groups.google.com/d/optout.

Reply via email to