Dear Jean-Paul,

Thank you for your response—I really appreciate your insights.

Regarding your suggestion: *"you should refactor the code (assembly, 
constraints creation + something before the Newton iterations start) such 
that you apply the constraints directly to the solution vector itself."*

When I apply my temperature boundary conditions differently using:


*MatrixTools::apply_boundary_values(boundary_values, tangent_matrix, 
solution_n, system_rhs); *

*(which modifies solution_n instead of newton_update)* before starting the 
Newton-Raphson iterations, everything works as expected, and I achieve the 
correct results.

However, when I instead use:










*VectorTools::interpolate_boundary_values(     dof_handler_ref,     
boundary_id,     
Functions::ConstantFunction<dim>(current_applied_temperature_0, dim + 1),   
  constraints,     fe.component_mask(Temperature) ); *

*  copy_local_to_global_ASM(const PerTaskData_ASM &data)*

*  {*

*  const AffineConstraints<double> 
<https://www.dealii.org/current/doxygen/deal.II/classAffineConstraints.html> 
&constraints 
= data.solid->constraints;*

*  BlockSparseMatrix<double> 
<https://www.dealii.org/current/doxygen/deal.II/classBlockSparseMatrix.html> 
&tangent_matrix 
= const_cast<Solid<dim,NumberType> *>(data.solid)->tangent_matrix;*

*  BlockVector<double> 
<https://www.dealii.org/current/doxygen/deal.II/classBlockVector.html> 
&system_rhs 
= const_cast<Solid<dim,NumberType> *>(data.solid)->system_rhs;*

 

*constraints.condense(tangent_matrix, system_rhs);*

*}*

 *constraints.distribute(newton_update);*

I notice a discrepancy—though the Newton-Raphson iterations still converge, 
the solution does not fully align with expectations. This seems quite odd, 
given that both approaches should theoretically produce the same results.

Would you have any thoughts on why this might be happening? Any guidance 
would be greatly appreciated.

 Regards,

Pradeep.

 
On Wednesday, March 19, 2025 at 7:51:37 AM UTC+1 Jean-Paul Pelteret wrote:

> Hi Alex,
>
> I haven't inspected your code in detail, but I would suggest that the 
> problem might be related to which system you are applying the constraints, 
> and into which vector you are distributing the constraints. It looks to me 
> that you are trying to specify the total displacement and temperature at 
> some locations, 
>
>       const double current_applied_Displacement = total_displacement 
>> *time_steps;
>>       const double current_applied_temperature_0 = 
>> total_temperature_0*time_steps;
>>       double current_applied_temperature_L= 
>> total_temperature_L*time_steps;
>
>
> while the constraints are being distributed to the solution increment
>
> constraints.distribute(newton_update);
>
>
> You might want to try either reformulating your boundary conditions as 
> updates to the BCs applied at the previous timestep (i.e. constrain the 
> increment) or you should refactor the code (assembly, constraints creation 
> + something before the Newton iterations start) such that you apply the 
> constraints directly to the solution vector itself.
>
> I hope that this helps.
>
> Best,
> Jean-Paul
>
> On Mon, 17 Mar 2025 at 18:50, Alex <alexander...@gmail.com> wrote:
>
>> Hi all,
>>
>> I would like to add a few statements to my above question, I have tried 
>> to boil down the issue and noticed that the issue arises with the *global 
>> right-hand side (RHS)* after applying AffineConstraints::condense(). 
>> Specifically, my RHS is initially zero, but after applying constraints, 
>> unexpected large negative values appear at certain degrees of freedom 
>> (DoFs).
>>
>> I have also attached my code below.
>>
>> *Setup*
>>
>>    - I have *2 elements* with *6 nodes*. 
>>    - The *far-left (200 kelvin) and far-right nodes  (0 Kelvin)* have 
>>    prescribed temperature Dirichlet boundary conditions, while the *middle 
>>    nodes remain unconstrained*. 
>>    - The *element-wise RHS* (cell_rhs) remains zero before assembling 
>>    the global system. 
>>    - I print the *global RHS before and after constraints* for debugging. 
>>
>> *Observed Behavior*
>>
>> *At time step = 1*, and *Newton iteration = 0*, I observe the following 
>> outputs:
>>
>>    1. *Before Applying Constraints:* 
>>
>> Global RHS BEFORE Constraints:
>>
>> 0.000e+00  0.000e+00  0.000e+00  0.000e+00  0.000e+00  0.000e+00  
>>
>> *2. After Applying Constraints:*
>>
>> Global RHS AFTER Constraints:
>>
>> 0.000e+00  0.000e+00  -1.900e+04  -1.900e+04  0.000e+00  0.000e+00  
>>
>>
>> Can someone guide me to where am I going wrong?
>>
>>
>> Thanks for your support.
>>
>>
>>
>>
>> On Mon, Mar 17, 2025 at 2:42 PM Alex <alexander...@gmail.com> wrote:
>>
>>> Hello everyone,
>>>
>>> I am modifying the Step-44 (Quasi-Static Finite-Strain Compressible 
>>> Elasticity) tutorial in deal.II to incorporate thermo-mechanical equations. 
>>> However, I am facing an issue with applying temperature Dirichlet boundary 
>>> conditions.
>>>
>>> Problem Description:
>>>
>>> I am applying the temperature boundary condition using:
>>>
>>> VectorTools::interpolate_boundary_values(
>>>
>>>     dof_handler_ref,
>>>
>>>     boundary_id,
>>>
>>>     Functions::ConstantFunction<dim>(current_applied_temperature_0, dim 
>>> + 1),  
>>>
>>>     constraints,
>>>
>>>     fe.component_mask(Temperature)
>>>
>>> );
>>>
>>> After that, before solving the system I apply:
>>>
>>> constraints.condense(tangent_matrix, system_rhs);
>>>
>>> Despite this, my solution converges, but the computed nodal temperatures 
>>> are negative and completely incorrect.
>>>
>>>  However, when I apply boundary conditions using:
>>>
>>>
>>> MatrixTools::apply_boundary_values(boundary_values,tangent_matrix,solution_n,system_rhs);
>>>
>>> I obtain the correct nodal temperatures.
>>>
>>> Could someone clarify what am I doing wrong in my first approach?  
>>>
>>> Thanks for your support !
>>>
>>>  
>>>
>>>  
>>>
>>>  
>>>
>>>  
>>>
>>>  
>>>
>>> -- 
>>> 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 visit 
>>> https://groups.google.com/d/msgid/dealii/0816bcff-fef9-40c9-972a-4ca7f626e75bn%40googlegroups.com
>>>  
>>> <https://groups.google.com/d/msgid/dealii/0816bcff-fef9-40c9-972a-4ca7f626e75bn%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+un...@googlegroups.com.
>>
> To view this discussion visit 
>> https://groups.google.com/d/msgid/dealii/CAK1JoDy69EbXrjK-aYx5DjBCbxdaFFEjWdeyQy_PO79601sd4Q%40mail.gmail.com
>>  
>> <https://groups.google.com/d/msgid/dealii/CAK1JoDy69EbXrjK-aYx5DjBCbxdaFFEjWdeyQy_PO79601sd4Q%40mail.gmail.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 visit 
https://groups.google.com/d/msgid/dealii/75e8854a-0c99-4e23-b98d-a8c05790606dn%40googlegroups.com.

Reply via email to