I've been thinking about this some more and there's something that I've
left out - my shape functions come from a MappingQEulerian object.
In particular, this is what the top of my assemble_system() looks like:
MappingQEulerian<dim,Vector<double>,spacedim> mapping(2, dof_handler,
euler_vector);
const QGauss<dim> quadrature_formula (2*fe.degree);
FEValues<dim,spacedim> fe_values (mapping, fe, quadrature_formula,
update_values |
update_normal_vectors |
update_gradients |
update_quadrature_points |
update_JxW_values);
I figure I either: (1) want to give DataPostProcessor the mapping, or (2)
give it all_normal_vectors from
std::vector<Tensor<1,spacedim> > all_normal_vectors = fe_values.
get_all_normal_vectors()
Do you suppose there's a better way to proceed from here? Looking ahead, I
can see that it would be useful to have H*n available for computations, not
just for output.
PS. In my previous message I meant to provide a link into the 8.4.1
documentation to what appears to be an inconsistency, here's what I meant
when I said:
> 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)
The documentation
https://www.dealii.org/8.4.1/doxygen/deal.II/classDataPostprocessor.html#af8f60ada7b7ae5417923613d87e16e0d
points to this code:
33 compute_derived_quantities_scalar (const std::vector<double>
&/*uh*/,
34 const std::vector<Tensor<1,dim> >
&/*duh*/,
35 const std::vector<Tensor<2,dim> >
&/*dduh*/,
36 *const std::vector<Point<dim> > *
&/*normals*/,
37 const std::vector<Point<dim> >
&/*evaluation_points*/,
38 std::vector<Vector<double> >
&computed_quantities) const
39 {
40 computed_quantities.clear();
41 AssertThrow(false,ExcPureFunctionCalled());
42 }
43
but I think normals is supposed to be a const
std::vector<Tensor<1,spacedim> >, according
to
https://www.dealii.org/8.4.1/doxygen/deal.II/classFEValuesBase.html#a22dae8b659a0b0d3a9c8ab6965bbe47f
If this passes tests then I guess there's no problem, it just appears to be
inconsistent. If this indeed needs a correction, then maybe we can walk
through it together, I am interested in learning how to help.
Tom
On Wednesday, October 26, 2016 at 6:08:00 PM UTC-4, Wolfgang Bangerth wrote:
>
> 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]
> <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.