Muhammad,

Dear Colleagues,
>                            I am working on thermal analysis where deal.ii
> Step-26 is being considered. In my case the thermal diffusivity
> (Coefficient being multiplied with Laplace Matrix) is temperature dependent
> and hence varying in every cell (or node). That is why I want to make the
> current Step-26 code flexible to cater this feature.
> I have two questions regarding this:
>
> *1) *Instead of directly creating the global Mass and Laplace matrices, I
> am making the cell_matrix as well as the cell_rhs as per following simple
> code (where to keep it simple, so far source term and diffusivity values
> are not yet present) :
> [...]
> The issue is that it is taking very long time to run 
> *VectorTools::point_value(dof_handler_temperature,
> old_solution_thermal, cell_temperature->vertex(j),
> temperature_old_solution_spatial_point); *function for evaluating the 
> *temperature_old_solution_spatial_point.
> *But the results match with the original Step-26 code.
> However as an efficient alternative, I tried to use the
> *temperature_old_solution_qpoint* in cell_rhs instead of
> *temperature_old_solution_spatial_point* which is very fast to compute
> but gives me the correct solution only at first time step after that the
> difference between this *new *and the *old solution* (from original
> Step-26 code) starts increasing i.e. new modified code solution kind of
> lags behind the original code (old) solution as per shared in the result
> snap shots in attachment.
>
> Any suggestion to improve the efficiency of current approach or any
> correction (in case I am mistaking or missing something) would be more than
> welcomed.
>

Yes, VectorTools::point_value is very slow because it searches in the whole
mesh each time it is called. Assuming that you use the same DoFHandler for
assembling the system_matrix and the old_solution_thermal, i.e.
dof_handler_temperature,
you should use FEValues::get_function_values()
https://www.dealii.org/current/doxygen/deal.II/classFEValuesBase.html#a357b422e374f2f2207af3512093f3907
instead. This function gives you the function values evaluated at every
quadrature point.
Currently, you are asking for the value at the vertices of the cell
instead. This can only possibly work if the number of vertices and the
number of degrees of freedom per cell coincide, i.e. for FE_Q(1) elements.
Using FEValues::get_function_values() you also don't need the loop with
running index j for assembling the right-hand side term.


> 2) Rather than making and assembling the local cell matrix and cell rhs,
> Is it more efficient and flexible way in deal.ii to directly modify the
> global Laplace matrix and system rhs in Step-26 for variable diffusivity
> coefficients at different corresponding dof entries?    (Hopefully the
> dynamic sparsity pattern and preconditioning etc. might not disturb the
> indexing of dofs in this case)
>

No, just modify assembling the local matrix to take the coefficient into
account.

Best,
Daniel

-- 
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/CAOYDWb%2BjD9KGb116ggqYFYnNAm7OXhHK-h%3D%3DZ0fYCKxMCUN3QQ%40mail.gmail.com.

Reply via email to