Prof. Bangerth,
Just to clarify, the Jacobian, its transpose inverse and its determinant
will be re-computed every time reinit() is called? Even if the mesh has not
changed from the previous time step?
If so, I might end up pre-computing and storing the Jacobian transpose
inverse and the Jacobian determinant to avoid re-computation. However, its
effect on the total time should be relatively small compared to flux
computations. "Premature optimization is the root of all evil".
Doug
On Sunday, February 10, 2019 at 2:15:04 PM UTC-5, Wolfgang Bangerth wrote:
>
>
> Doug,
>
> > I am currently looking into the implementation FEValues. Here is my
> understanding
> >
> > 1. FEValues takes some Mapping, FiniteElement, and Quadrature to
> evaluate and
> > store the values with the required UpdateFlags.
> > 2. FEValues::initialize through the FiniteElement::get_data will then
> > pre-compute all the required data that are possibly computable on
> the
> > reference cell.
> > 3. FEValues::reinit will then compute/fetch the data on a physical
> cell,
> > using the Jacobian information by calling the various
> fill_fe_values()
> > functions.
> >
> > I see that most of the examples use implicit methods to solve the
> various
> > problems. Therefore, the cost of initialize() and reinit() on each cell
> is
> > quite low compared to solving the actual system at every step. However,
> for
> > explicit methods such as RK, we wouldn't want to re-compute the Jacobian
> or
> > reference gradients (\nabla_{physical} \phi = J^{-T} \nabla _{reference}
> \phi)
> > at every time step. Therefore, the FEValues declaration would be outside
> the
> > "assemble_system()" seen in most tutorials, and passed as an argument as
> > "assemble_system(FEValues fe_v)". Unfortunately, I am getting lost in
> the code
> > where things values are pre-computed versus fetching them.
>
> I think that's not necessary. The creation of an FEValues object is
> expensive,
> but on meshes with a non-trivial number of cells, the creation of the
> FEValues
> object is still cheap compared to what you're going to do with the object
> on
> the cells. (I.e., it is expensive compared to the operations on one cell;
> but
> the creation is cheap compared to the operations on *all* cells.)
>
>
> > My question is: What is being re-computed every time I call
> > FEValues::reinit(cell)?
>
> initialize() computes everything that can be computed on the reference
> cell.
> This includes the values and gradients of the shape functions with regard
> to
> the reference coordinates.
>
> reinit() then only computes everything that depends on the actual cell.
> This
> means for example computing the Jacobian of the mapping on a given cell,
> and
> multiplying it by the (precomputed) gradients of the shape functions on
> the
> reference cell to obtain the gradient on the real cell.
>
>
> > I am mainly worried about it re-computing Jacobians (its determinant and
> its
> > inverse) and evaluating its product with the reference gradients
> > (\nabla_{reference} phi) every time reinit(cell) is called. Otherwise,
> the
> > solution will be to pre-compute and store those outside in some other
> > container and avoid using FEValues.
>
> Yes. This (the Jacobian of the mapping and its determinant) is a quantity
> that
> changes from each cell to the next, and so it is computed in reinit()
> because
> it simply can't be precomputed without knowledge of the actual cell.
>
> In explicit methods, you will likely visit the same mesh over and over,
> and so
> it does make sense to compute things such as the gradient of the shape
> functions on each cell (rather than the Jacobian or its determinant, which
> goes into the computation) once at the beginning of the program and then
> re-use it where necessary. You can of course do this pre-computation using
> FEValues :-)
>
> An alternative is the matrix-free framework that can do some of these
> things
> as well. You may want to look at the corresponding example programs.
>
> Best
> W.
>
>
> --
> ------------------------------------------------------------------------
> Wolfgang Bangerth email: [email protected]
> <javascript:>
> 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.