Hi Nik,

> My question is: for boundary_worker. What am I missing that the
calculated residual contributions
don’t show up in the global RHS of the system?

I think that I might have an idea as to what's going on here. I feel that
the problem likely lies with what these lines

Vector<double> &cell_rhs = copy_data.cell_rhs;
...
ad_helper.compute_residual(cell_rhs);

do in unison, versus the

copy_data.cell_rhs(i) +=

that you have in the hand written assembly code. Normally you'd be
accumulating into the cell RHS/matrix (as you've demonstrated), but the
compute_residual/linearization functions don't do this. Instead, they
overwrite the local matrix/vector (although its a little deep into the AD
framework, you can see this here
https://github.com/dealii/dealii/blob/master/source/differentiation/ad/ad_drivers.cc#L1836-L1854).
This is incompatible with the meshloop() framework as it reuses the same
local matrices and vectors on a cell and its faces before scattering their
contributions to the global linear system.

What you could do as a workaround is to simply create an intermediate local
matrix and RHS vector in each of the workers, have the ADHelper write into
those, and then accumulate the result into the data structures in copy_data.

I hope that my thoughts were on the right track, and that this helps. If it
does, then this would motivate the addition of a set of functions along the
lines of

ad_helper.accumulate_residual(copy_data.cell_rhs, -1.0 /* factor */);
ad_helper.accumulate_linearization(copy_data.cell_matrix);

that would work better with mesh_loop(). We'd be happy to accept any
patches that help to realize this functionality!

Best,
Jean-Paul


On Tue, 9 Jul 2024 at 05:17, Wolfgang Bangerth <[email protected]>
wrote:

> On 7/8/24 08:51, Nik Markus Leuenberger wrote:
> >
> > I'm opening a new thread because the question is no longer whether or
> how to
> > use deal.ii for DG/finite volume implementations, but what I hope is a
> more
> > specific question regarding a combination of step-74 and step-72.
> >
> > I want to use automatic differentiation to obtain the Jacobian for a
> > (non)-linear residual of a Poisson equation. (Non) is in brackets
> because I
> > would like to solve the linear case first in a single Newton-Raphson
> step.
> >
> > The question described in more detail in the attached .pdf is how to use
> > automatic differentiation for the boundary_worker function of step-74.
>
> Nik:
> I'm not sure anyone can help with this problem short of debugging it in
> detail
> -- which I think you'll have to do yourself. But it may be useful to see
> how
> others have done something similar, for example in this recent code
> gallery
> program:
>
>
> https://dealii.org/developer/doxygen/deal.II/code_gallery_nonlinear_heat_transfer_with_AD_NOX.html
>
> 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].
> To view this discussion on the web visit
> https://groups.google.com/d/msgid/dealii/15834e9c-a00c-4f1a-a30e-602c2aaca296%40colostate.edu
> .
>

-- 
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].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/CAMvdQMDaHYHCB%3DSHFdjxKJMwm9%2BTqJ1wtoZawJ0gw%3DTQ0DK_8A%40mail.gmail.com.

Reply via email to