On 10/27/2016 08:43 AM, thomas stephens wrote:
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()
|
You need to give the mapping to DataOut::build_patches. It knows how to
use it internally to evaluate things, including 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.
That's something for which you need to use your own FEValues then.
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
This is a bit of a sore point -- a design flaw:
https://github.com/dealii/dealii/issues/391
See the second comment there. It *should* be a vector of Tensor<1,dim>,
but it's not, unfortunately :-(
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.