On 10/26/2016 03:30 PM, thomas stephens wrote:
I have already obtained H*H using a class derived from
DataPostProcessorScalar, and that makes sense to me - the function
signature of compute_derived_quantities_vector() looks like this:
template <int spacedim>
void
VectorValuedSolutionSquared<spacedim>::compute_derived_quantities_vector
(const std::vector<Vector<double> > &uh,
const
std::vector<std::vector<Tensor<1, spacedim> > > & /*duh*/,
const
std::vector<std::vector<Tensor<2, spacedim> > > & /*dduh*/,
const
std::vector<Point<spacedim> > & /*normals*/,
const
std::vector<Point<spacedim> > & /*evaluation_points*/,
std::vector<Vector<double> >
&computed_quantities) const
since H*H doesn't require anything other than uh, but now I would like
to take advantage of that vector of normals (which, for dealii 8.44.1,
should be of type const std::vector<Tensor<1,spacedim> >, see here - but
I digress)
At any rate, I don't know how to access this normals vector, and if I
were to compute it ahead of time using fe_values.get_all_normals() then
I am not sure how to pass it to compute_derived_quantities(), since my
code never actually calls that function.
At the end of the day I think this comes down to me fighting with a
solution H of type Vector<double> of size dof_handler.n_dofs(), where
that dof_handler was built with vector-valued finite elements, and now I
want to look at essentially a scalar function on my geometry.
Tom -- that function is called internally from DataOut for all
evaluation points on a cell. If you derive a class from
DataOutPostprocessorVector (as you already do) and tell the constructor
of the DataOutPostprocessorVector base class that it should also provide
you with the normal vectors, then your implementation of that function
(in your derived class) should have access to the normal vectors.
You then only need to return H*n in the last argument.
(I'm going to say that I'm not 100% sure that the normal vectors are
actually used for the *cell normals* when using DataOut, or if they are
only used for *face normals* by DataOutFaces. If the latter, it should
not be terrible difficult to add this functionality, and we'd be happy
to walk you through it.)
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].
For more options, visit https://groups.google.com/d/optout.