Dear Alberto,

>From the link in your previous post I presume that you're working with 
deal.II 8.4.1? The feature that I've referred to was introduced in version 
8.5.0. In fact, I see now that the latest version of step-33 
<https://www.dealii.org/8.5.0/doxygen/deal.II/step_33.html#EulerEquationsPostprocessor>
 
has been updated to demonstrate this new functionality. Does that help?

Regards,
Jean-Paul

On Friday, August 18, 2017 at 2:30:12 PM UTC+2, Alberto Salvadori wrote:
>
> Thank you Jean-Paul.
>
> From the CommonInputs class 
> <https://www.dealii.org/8.5.0/doxygen/deal.II/structDataPostprocessorInputs_1_1CommonInputs.html>,
>  
> I understand that I can retrieve a pointer to the current cell by using 
> this instruction, provided that the right DoFHandler is passed:
>
> const typename hp::DoFHandler<dim>::cell_iterator 
> <https://www.dealii.org/8.5.0/doxygen/deal.II/group__Iterators.html#ga60277a8a3957ba4b41c1e76a87decd30>
> current_cell = input_data.template get_cell<hp::DoFHandler<dim> >();
>
> which is indeed what I am looking for. However, I have no clue about the 
> input_data  variable, which is passed to the evaluate_vector_field function 
> in the example that computes the fluid norm of the stress, and how it 
> connects to the Post-process class that is described in step-33, i.e.
>
> template <int dim>
> class BoussinesqFlowProblem<dim>::Postprocessor : public DataPostprocessor 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classDataPostprocessor.html>
> <dim>
> {
> public:
> Postprocessor (const unsigned int partition,
> const double minimal_pressure);
> virtual
> void
> compute_derived_quantities_vector 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classDataPostprocessor.html#adf1b8b97f2b2f9c65d4c7b143f61865f>
>  
> (const std::vector<Vector<double> 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classVector.html> > &uh,
> const std::vector<std::vector<Tensor<1,dim> 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classTensor.html> > > &duh,
> const std::vector<std::vector<Tensor<2,dim> 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classTensor.html> > > &dduh,
> const std::vector<Point<dim> 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classPoint.html> > &normals,
> const std::vector<Point<dim> 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classPoint.html> > 
> &evaluation_points,
> std::vector<Vector<double> 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classVector.html> > 
> &computed_quantities) const;
> virtual std::vector<std::string> get_names 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classDataPostprocessor.html#a254f38bcdf4bdb5aa94231b695da7d55>
>  
> () const;
> virtual
> std::vector<DataComponentInterpretation::DataComponentInterpretation>
> get_data_component_interpretation 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classDataPostprocessor.html#a4d808c941d39a0c41cc9ac4afa82d47e>
>  
> () const;
> virtual UpdateFlags 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/group__feaccess.html#gaa94b67d2fdcc390690c523f28019e52f>
>  
> get_needed_update_flags 
> <https://www.dealii.org/8.4.1/doxygen/deal.II/classDataPostprocessor.html#aadecdd040447b395164397ea1196f721>
>  
> () const;
> private:
> const unsigned int partition;
> const double minimal_pressure;
> };
>
> Can you perhaps provide some further information, to fill the gap?
> Thank you,
>
> Alberto
>
>
> *Alberto Salvadori* Dipartimento di Ingegneria Civile, Architettura, 
> Territorio, Ambiente e di Matematica (DICATAM)
>  Universita` di Brescia, via Branze 43, 25123 Brescia
>  Italy
>  tel 030 3711239
>  fax 030 3711312
>
> e-mail: 
>  alberto.salvad...@unibs.it
> web-pages:
>  http://m4lab.unibs.it/faculty.html
>  http://dicata.ing.unibs.it/salvadori
>
> On Fri, Aug 18, 2017 at 12:38 PM, Jean-Paul Pelteret <> wrote:
>
>> Dear Alberto,
>>
>> It looks like one can already do this, although its somewhat obscurely 
>> documented. Here's a link to the documentation of the CommonInputs class 
>> <https://www.dealii.org/8.5.0/doxygen/deal.II/structDataPostprocessorInputs_1_1CommonInputs.html>,
>>  
>> wherein you can find a documented example of computing some solution and 
>> cell dependent data. It seems to be a relatively new feature (circa 2016) 
>> and I didn't know about it, so I'm happy to have found out about it as well!
>>
>> Regards,
>> Jean-Paul
>>
>> On Friday, August 18, 2017 at 11:16:40 AM UTC+2, Alberto Salvadori wrote:
>>>
>>> Dear all,
>>>
>>> I have been studying step-32 and I have found the class Postprocessor nice 
>>> and effective. 
>>> It is my understanding that the method compute_derived_quantities_vector 
>>> within the class operates at a call level, i.e. that beyond the hood the 
>>> class Postprocessor
>>> implements loops over the triangulation. 
>>>
>>> I wonder if there is a way to use the cell->user_pointer in the method 
>>> compute_derived_quantities_vector, i.e. if there is a way to access a 
>>> reference to the cell itself.
>>>
>>> This would be very useful in case of post-processing of data that depend 
>>> on the history at each Gauss point, which is indeed my case of interest.
>>>
>>> Thanks!
>>> Alberto
>>>
>> -- 
>> 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.
>>
>
>
>
> Informativa sulla Privacy: http://www.unibs.it/node/8155
>

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

Reply via email to