Hi Jean,

Thanks a lot for your detailed answer to my questions. Now I am able to run 
my program and the newton method is converging. I completely misunderstood 
the documentation. Also, thank you for opening up the issue. It would 
really be helpful if the complications involved in making constraints can 
be removed. 

Thanks again!

Animesh

On Monday, October 12, 2020 at 11:20:32 AM UTC+5:30 Jean-Paul Pelteret 
wrote:

> Hi Animesh,
>
> Here are some answers to your questions:
>
>
> Essentially, I will have to apply the homogeneous boundary condition at 
> the zeroth and first Newton's step, and Non-homogeneous BC only at the 
> zeroth step, as mentioned in the documentation 
> <https://www.dealii.org/current/doxygen/deal.II/step_44.html#Solidmake_constraints>
> .
>
> 1. Why would I need to give the homogeneous boundary conditions twice?
>
>
> No, I’m afraid you’ve misunderstood what’s written there. Here’s the 
> important part of the documentation:
> “”"
> However, since we are dealing with an iterative Newton method, it should 
> be noted that any displacement constraints should only be specified at the 
> zeroth iteration and subsequently no additional contributions are to be 
> made since the constraints are already exactly satisfied.
> “””
>
> You only apply (inhomogeneous) displacement constraints at the zeroth 
> Newton iteration. Because this tutorial uses an incremental formulation, 
> and the constraints are applied to the increment itself, you want to make 
> sure that you apply the correct implement only once per timestep. This is 
> because you effectively displace your boundary from its position at the 
> previous tilmestep to the one at the current tilmestep using this single 
> update. After this update is done, you want to make sure that the boundary 
> doesn’t move away from that location, so you apply homogeneous constraints 
> on the same displacement components for each subsequent Newton step. That 
> is the other part of the logic that’s laid out in 
> step-44’s Solid<dim>::make_constraints().
>
> 2. Is it correct the way I am applying the constraints? I am copying the 
> part of the code in red? The major concern that I have is that if I do 
> constraints.clear(), then it will also clear my non-homogeneous constraint 
> applied at the zeroth Newton step. How should I handle this?
>
>
> Not quite. You do want the inhomogeneities cleared, but you still only 
> want to fix the displacement components that you originally applied 
> inhomogeneous constraints to. You’ll want to duplicate the logic applied in 
> the other branches of the make_constraints() function. So something like 
> this:
>
>     // Apply Non-homogeneous B.C.
>     {
>      const int boundary_id = 5;
>       if (apply_dirichlet_bc == true)
>         VectorTools::interpolate_boundary_values(dof_handler_ref,
>                                                  boundary_id,
>                                                  ConstantFunction<dim>(dim 
> == 3 ? std::vector<double> { -0.0001, 0.0, 0.0 } : std::vector<double> { 
> -0.0001, 0.0 }), constraints, fe.component_mask(x_displacement));
>
>         else
>
>         VectorTools::interpolate_boundary_values(dof_handler_ref,
>                                                  boundary_id,
>                                                  ZeroFunction<dim>(), 
> constraints, fe.component_mask(x_displacement));
>
>     }
>
>
> I must admit that this code is overly complicated, and it can be 
> simplified so that extending it to the case of inhomogeneous constraints is 
> easier. I’ve opened up an issue for this task here 
> <https://github.com/dealii/dealii/issues/11035> and I’ll update the 
> tutorial in the coming weeks.
>
> I hope that this helps you in understanding the logic and resolving your 
> issue!
>
> Best,
> Jean-Paul
>
>
> On 10 Oct 2020, at 19:05, Animesh Rastogi IIT Gandhinagar <
> animesh...@alumni.iitgn.ac.in> wrote:
>
> Hi All, 
>
> I am working on hyper-elastic beam buckling problem. I wish to apply 
> non-homogeneous Dirichlet Boundary condition on the right hand side of the 
> beam (to compress the beam). I was reading the step-44 
> <https://www.dealii.org/current/doxygen/deal.II/step_44.html#Solidmake_constraints>
>  
> where they have mentioned about how to apply non-homogeneous boundary 
> condition.  Essentially, I will have to apply the homogeneous boundary 
> condition at the zeroth and first Newton's step, and Non-homogeneous BC 
> only at the zeroth step, as mentioned in the documentation 
> <https://www.dealii.org/current/doxygen/deal.II/step_44.html#Solidmake_constraints>
> .
>
> I have following two questions regarding this - 
>
> 1. Why would I need to give the homogeneous boundary conditions twice?
> 2. Is it correct the way I am applying the constraints? I am copying the 
> part of the code in red? The major concern that I have is that if I do 
> constraints.clear(), then it will also clear my non-homogeneous constraint 
> applied at the zeroth Newton step. How should I handle this?
>
> I have edited the following code from code gallery - 
>
>
> https://dealii.org/developer/doxygen/deal.II/code_gallery_Quasi_static_Finite_strain_Compressible_Elasticity.html
>
> Thanks!
>
> Animesh
>
>
>
>
> template <int dim,typename NumberType>
>   void Solid<dim,NumberType>::make_constraints(const int &it_nr)
>   {
>     std::cout << " CST " << std::flush;
>     if (it_nr > 1)
>       return;
>     constraints.clear();
>     const bool apply_dirichlet_bc = (it_nr == 0);
>     const FEValuesExtractors::Scalar x_displacement(0);
>     // Apply Non-homogeneous B.C.
>     {
>      const int boundary_id = 5;
>       if (apply_dirichlet_bc == true)
>         VectorTools::interpolate_boundary_values(dof_handler_ref,
>                                                  boundary_id,
>                                                  ConstantFunction<dim>(dim 
> == 3 ? std::vector<double> { -0.0001, 0.0, 0.0 } : std::vector<double> { 
> -0.0001, 0.0 }), constraints, fe.component_mask(x_displacement));
>     }
> //Fix the left hand side..
>     {
>       const int boundary_id = 1;
>
>       if (apply_dirichlet_bc == true)
>         VectorTools::interpolate_boundary_values(dof_handler_ref,
>                                                  boundary_id,
>                                                  
> ZeroFunction<dim>(n_components),
>                                                  constraints,
>                                                  fe.component_mask(u_fe));
>       else
>         VectorTools::interpolate_boundary_values(dof_handler_ref,
>                         boundary_id,
>                                                  
> ZeroFunction<dim>(n_components),
>                                                  constraints,
>                                                  fe.component_mask(u_fe));
>     }
>
>
>     // Zero Z-displacement through thickness direction
>     // This corresponds to a plane strain condition being imposed on the 
> beam
>     if (dim == 3)
>       {
>         const int boundary_id = 2;
>         const FEValuesExtractors::Scalar z_displacement(2);
>
>         if (apply_dirichlet_bc == true)
>           VectorTools::interpolate_boundary_values(dof_handler_ref,
>                                                    boundary_id,
>                                                    
> ZeroFunction<dim>(n_components),
>                                                    constraints,
>                                                    
> fe.component_mask(z_displacement));
>         else
>           VectorTools::interpolate_boundary_values(dof_handler_ref,
>                                                    boundary_id,
>                                                    
> ZeroFunction<dim>(n_components),
>                                                    constraints,
>                                                    
> fe.component_mask(z_displacement));
>       }
>
>     constraints.close();
>   }
>
> -- 
> 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.
> To view this discussion on the web visit 
> https://groups.google.com/d/msgid/dealii/69400803-2427-4b62-be28-7d164c1cb1aan%40googlegroups.com
>  
> <https://groups.google.com/d/msgid/dealii/69400803-2427-4b62-be28-7d164c1cb1aan%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/35915df0-25ee-43c4-b59a-6c4961e879e0n%40googlegroups.com.

Reply via email to