Dear Hamed

Yes it is the correct way. You can add each data of history_field into 
output using history_dof_handler. For example:

  DataOut<dim> data_out;

  std::vector<DataComponentInterpretation::DataComponentInterpretation>
              data_component_interpretation(1, 
DataComponentInterpretation::component_is_scalar);

  data_out.add_data_vector(history_dof_handler, history_field[0][0], 
"Name_of_out",
                data_component_interpretation);

and for the main field use same the DataOut object but with the main 
dofhandler object, for example for 2D:

  {
    std::vector<std::string> solution_names(dim, "dis");
    std::vector<DataComponentInterpretation::DataComponentInterpretation> 
data_component_interpretation(
      dim, DataComponentInterpretation::component_is_part_of_vector);

    data_out.add_data_vector(dof_handler, solution,
                             solution_names, data_component_interpretation);
  }

  {

    std::vector<std::string> solution_names;
    solution_names.push_back("displacement_x");
    solution_names.push_back("displacement_y");


    std::vector<DataComponentInterpretation::DataComponentInterpretation>
    data_component_interpretation( n_components,
            DataComponentInterpretation::component_is_scalar);

    data_out.add_data_vector(dof_handler, solution,
                             solution_names, data_component_interpretation);
  }

  data_out.build_patches();



On Sunday, July 24, 2016 at 8:01:15 PM UTC+4:30, Hamed Babaei wrote:
>
> Hi Gunnar,
>
> I have written a code based on step 18 for a elastic problem. 
> I need to output individual stress tensor components at some specific 
> points of the domain. I think I can use the following piece of 
> pseudo-code presented at the end of step 18 tutorial to make a conversion 
> from individual quadrature points to a global finite element field. 
>
> FE_DGQ<dim> <http://dealii.org/8.4.1/doxygen/deal.II/classFE__DGQ.html> 
> history_fe (1);
> DoFHandler<dim> 
> <http://dealii.org/8.4.1/doxygen/deal.II/classDoFHandler.html> 
> history_dof_handler (triangulation);
> history_dof_handler.distribute_dofs (history_fe);
> std::vector< std::vector< Vector<double> > >
> history_field (dim, std::vector< Vector<double> 
> <http://dealii.org/8.4.1/doxygen/deal.II/classVector.html> >(dim)),
> local_history_values_at_qpoints (dim, std::vector< Vector<double> 
> <http://dealii.org/8.4.1/doxygen/deal.II/classVector.html> >(dim)),
> local_history_fe_values (dim, std::vector< Vector<double> 
> <http://dealii.org/8.4.1/doxygen/deal.II/classVector.html> >(dim));
> for (unsigned int i=0; i<dim; i++)
> for (unsigned int j=0; j<dim; j++)
> {
> history_field[i][j].reinit(history_dof_handler.n_dofs());
> local_history_values_at_qpoints[i][j].reinit(quadrature.size());
> local_history_fe_values[i][j].reinit(history_fe.dofs_per_cell);
> }
> FullMatrix<double> 
> <http://dealii.org/8.4.1/doxygen/deal.II/classFullMatrix.html> 
> qpoint_to_dof_matrix (history_fe.dofs_per_cell,
> quadrature.size());
> FETools::compute_projection_from_quadrature_points_matrix 
> <http://dealii.org/8.4.1/doxygen/deal.II/namespaceFETools.html#ab1a0545c897ee022029f8c5f2c252735>
> (history_fe,
> quadrature, quadrature,
> qpoint_to_dof_matrix);
> typename DoFHandler<dim>::active_cell_iterator 
> <http://dealii.org/8.4.1/doxygen/deal.II/classDoFHandler.html> cell = 
> dof_handler.begin_active 
> <http://dealii.org/8.4.1/doxygen/deal.II/classDoFHandler.html#ad2df51a32781906218910731e063ac62>
> (),
> endc = dof_handler.end 
> <http://dealii.org/8.4.1/doxygen/deal.II/classDoFHandler.html#a05c70e1862a2ad145f91f9da1f44cc28>
> (),
> dg_cell = history_dof_handler.begin_active();
> for (; cell!=endc; ++cell, ++dg_cell)
> {
> PointHistory<dim> *local_quadrature_points_history
> = reinterpret_cast<PointHistory<dim> *>(cell->user_pointer());
> Assert 
> <http://dealii.org/8.4.1/doxygen/deal.II/group__Exceptions.html#ga70a0bb353656e704acf927945277bbc6>
>  
> (local_quadrature_points_history >=
> &quadrature_point_history.front(),
> ExcInternalError());
> Assert 
> <http://dealii.org/8.4.1/doxygen/deal.II/group__Exceptions.html#ga70a0bb353656e704acf927945277bbc6>
>  
> (local_quadrature_points_history <
> &quadrature_point_history.back(),
> ExcInternalError());
> for (unsigned int i=0; i<dim; i++)
> for (unsigned int j=0; j<dim; j++)
> {
> for (unsigned int q=0; q<quadrature.size(); ++q)
> local_history_values_at_qpoints[i][j](q)
> = local_quadrature_points_history[q].old_stress[i][j];
> qpoint_to_dof_matrix.vmult (local_history_fe_values[i][j],
> local_history_values_at_qpoints[i][j]);
> dg_cell->set_dof_values (local_history_fe_values[i][j],
> history_field[i][j]);
> }
> }
>
> would you believe this lines of code would be what I am looking for?
> In addition, I was wondering how to output data for a particular point in 
> the domain from history_field finite element field.
>
> Thanks,
> Hamed
>

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

Reply via email to