Dear Wolfgang,

below is the code, where  the function
  template<int dim>
  void HeatEquation<dim>::tensor1_output()
  {
    ...
  }
performs the gmv-output of
  template <int dim>
  class AdvectionField : public TensorFunction<1,dim>
  {
     // define 2D velocity vector field
  }

     Best wishes, JD

**********************************************************

 // description of the advection field
  template <int dim>
  class AdvectionField : public TensorFunction<1,dim>
  {
  public:
    AdvectionField () : TensorFunction<1,dim> () {}

    virtual Tensor<1,dim> value (const Point<dim> &p) const;

    virtual void value_list (const std::vector<Point<dim> > &points,
                             std::vector<Tensor<1,dim> >    &values) const;

    // corresponding exception
    DeclException2 (ExcDimensionMismatch,
                    unsigned int, unsigned int,
                    << "The vector has size " << arg1 << " but should have "
                    << arg2 << " elements.");
  };


  template <int dim>
  Tensor<1,dim>
  AdvectionField<dim>::value (const Point<dim> &p) const
  {
    Point<dim> value;
    value[0] = p[0];
    for (unsigned int i=1; i<dim; ++i)
      value[1] = p[1];

    return value;
  }


  // all values of the advection field
  template <int dim>
  void
  AdvectionField<dim>::value_list (const std::vector<Point<dim> > &points,
                                   std::vector<Tensor<1,dim> >    &values)
const
  {
    Assert (values.size() == points.size(),
            ExcDimensionMismatch (values.size(), points.size()));

    for (unsigned int i=0; i<points.size(); ++i)
      values[i] = AdvectionField<dim>::value (points[i]);
  }

 // output of the tensor function
  template<int dim>
  void HeatEquation<dim>::tensor1_output()
  {

    DataOut<dim> data_out;
    Vector<double> vec_output;

    DoFHandler<dim>      local_dof_handler(triangulation);
    FESystem<dim>        my_fe(FE_Q<dim>(1), dim);
    local_dof_handler.distribute_dofs (my_fe);

    vec_output.reinit ( local_dof_handler.n_dofs() );

    VectorTools::interpolate( local_dof_handler,
                               VectorFunctionFromTensorFunction<dim>(
AdvectionField<dim>(), 0, dim),
                               vec_output );

    std::vector<std::string> solution_names(dim, "vec");
    for (int i=0; i<dim; i++)
        solution_names[i] += Utilities::int_to_string(i+1, 1);

    std::vector<DataComponentInterpretation::DataComponentInterpretation>
    data_component_interpretation
    (dim, DataComponentInterpretation::component_is_part_of_vector);

    data_out.attach_dof_handler( local_dof_handler );

    data_out.add_data_vector ( vec_output, solution_names,
                          DataOut<dim>::type_dof_data,
                          data_component_interpretation);

    data_out.build_patches();
    const std::string filename = "solution-"
                                 +
Utilities::int_to_string(timestep_number, 3) +
                                 ".gmv";
    std::ofstream output(filename.c_str());
    data_out.write_gmv(output);

    // release
    solution_names.clear();
    data_out.clear();
    local_dof_handler.clear();

  }



On Thu, Nov 3, 2016 at 3:44 PM, Wolfgang Bangerth <bange...@colostate.edu>
wrote:

> On 11/03/2016 07:54 AM, Julian Dorn wrote:
>
>>
>>
>>   template <int dim>
>>   class AdvectionField : public TensorFunction<1,dim>
>>   {
>>      // define 2D velocity vector field
>>   }
>>
>> How to do its Paraview/GMV-output?
>>
>
> The easiest way is probably to create an FESystem(FE_Q<dim>, dim) (i.e., a
> finite element with as many components as the vector field), then create a
> DoFHandler, use VectorTools::interpolate to interpolate your vector field
> onto this DoFHandler, and then use DataOut to write the result into a file.
>
> I think this would actually make for a nice self-contained example that
> one could put into the "Possibilities for extensions" of one of the
> tutorials. Want to give this a try and see whether you can come up with
> some code that does that?
>
> Cheers
>  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/fo
> rum/dealii?hl=en
> --- You received this message because you are subscribed to a topic in the
> Google Groups "deal.II User Group" group.
> To unsubscribe from this topic, visit https://groups.google.com/d/to
> pic/dealii/5yVcTyB8tC0/unsubscribe.
> To unsubscribe from this group and all its topics, send an email to
> dealii+unsubscr...@googlegroups.com.
> For more options, visit https://groups.google.com/d/optout.
>

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