Re: [deal.II] FEValues pre-computed values

2019-02-11 Thread Wolfgang Bangerth
On 2/11/19 1:52 PM, Doug wrote:
> 
> 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?

Yes. We don't store this information for all cells -- in fact, FEValues 
doesn't actually know anything about the mesh as a whole, all it knows is 
about the current cell. So there is no place to look this information up: if 
the current cell is different from the previous cell, then this information 
needs to be recomputed by FEValues. Of course, that doesn't prevent you from 
making a run through all cells, storing this information for each cell, and in 
all following loops over all cells use what you cached rather than utilizing 
FEValues. That's the classic memory-vs-compute trade-off.

(I will note that FEValues has optimizations that make some of this cheaper. 
For example, if the cell for which you call reinit() is simply a translation 
of the previous cell, then it is not necessary to recompute the Jacobian, 
transpose, inverse, and determinant because it will be the same as for the 
previous cell. There is some trickery involved in this when using multiple 
threads that might defeat the optimization, but it can be enabled.)


> 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".

Well said. Implement the simplest version of the code first and optimize once 
you have hard evidence that something is an actual problem.

Best
  W.

-- 

Wolfgang Bangerth  email: bange...@colostate.edu
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 dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.


Re: [deal.II] FEValues pre-computed values

2019-02-11 Thread Doug
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: bang...@colostate.edu 
>  
> 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 dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.


Re: [deal.II] FEValues pre-computed values

2019-02-10 Thread Wolfgang Bangerth


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: bange...@colostate.edu
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 dealii+unsubscr...@googlegroups.com.
For more options, visit https://groups.google.com/d/optout.


[deal.II] FEValues pre-computed values

2019-02-08 Thread Doug
Hello again,

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.

My question is: What is being re-computed every time I call 
FEValues::reinit(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.

Doug

-- 
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.
For more options, visit https://groups.google.com/d/optout.