Hi, all
I was practicing using post_processor, referring step-29 and step 31.
If I want to use post_processor for vector-valued problem, where my
solution is (u,v,p)
and I want to compute magnitude of velocity, |u^2+v^2| ignoring pressure
parts.
i wrote as step 29;
template <int dim>
class Magnitude : public DataPostprocessorScalar<dim>
{
public:
Magnitude ();
virtual
void
compute_derived_quantities_vector (const std::vector<Vector<double> >
&uh,
const std::vector<std::vector<Tensor<
1, dim> > > &duh,
const std::vector<std::vector<Tensor<
2, dim> > > &dduh,
const std::vector<Point<dim> >
&normals,
const std::vector<Point<dim> >
&evaluation_points,
std::vector<Vector<double> >
&computed_quantities) const;
};
template <int dim>
Magnitude<dim>::Magnitude ()
:
DataPostprocessorScalar<dim> ("ShearRate", update_values)
{}
// Remember that this function may only use data for which the
// respective update flag is specified by
template <int dim>
void
Magnitude<dim>::compute_derived_quantities_vector (
const
std::vector<Vector<double> > &uh,
const
std::vector<std::vector<Tensor<1, dim> > > & /*duh*/,
const
std::vector<std::vector<Tensor<2, dim> > > & /*dduh*/,
const
std::vector<Point<dim> > & /*normals*/,
const
std::vector<Point<dim> > & /*evaluation_points*/,
std::vector<Vector<double> > &computed_quantities
) const
{
for (unsigned int i=0; i<computed_quantities.size(); i++)
{
Assert(computed_quantities[i].size() == 1,
ExcDimensionMismatch (computed_quantities[i].size(), 1));
Assert(uh[i].size() == 2, ExcDimensionMismatch (uh[i].size(), 2));
computed_quantities[i](0) = std::sqrt(uh[i](0)*uh[i](0) + uh[i](1
)*uh[i](1));
}
}
and I used it as wrote as step 29;
template <int dim>
void
StokesProblem<dim>::output_results (const unsigned int
refinement_cycle) const
{
std::vector<std::string> solution_names (dim, "velocity");
solution_names.push_back ("pressure");
std::vector<DataComponentInterpretation::DataComponentInterpretation>
data_component_interpretation
(dim, DataComponentInterpretation::component_is_part_of_vector);
data_component_interpretation
.push_back (DataComponentInterpretation::component_is_scalar);
DataOut<dim> data_out;
data_out.attach_dof_handler (dof_handler);
data_out.add_data_vector (solution, solution_names,
DataOut<dim>::type_dof_data,
data_component_interpretation);
// added
Magnitude<dim> magnitude;
*data_out.add_data_vector (solution.block(**0), magnitude);
//because I only want to process velocity solution*
// added - end ;
data_out.add_data_vector (cellwise_shear_rate,"Shear_rate");
data_out.add_data_vector (cellwise_viscosity,"Viscosity");
data_out.build_patches ();
std::ostringstream filenameeps;
filenameeps << "solution-"<< Utilities::int_to_string
(refinement_cycle, 2)<< ".vtk";
std::ofstream output (filenameeps.str().c_str());
data_out.write_vtk (output);
}
but I get error like
*An error occurred in line <1022> of file
</Users/kimjaekwang/dealii-8.4.1/source/numerics/data_out_dof_data.cc> in
function*
* void dealii::DataOut_DoFData<dealii::DoFHandler<2, 2>, 2,
2>::add_data_vector(const VectorType &, const
DataPostprocessor<DoFHandlerType::space_dimension> &) [VectorType =
dealii::Vector<double>]*
*The violated condition was: *
* vec.size() == dofs->n_dofs()*
*The name and call sequence of the exception was:*
* Exceptions::DataOut::ExcInvalidVectorSize (vec.size(), dofs->n_dofs(),
dofs->get_triangulation().n_active_cells())*
*Additional Information: *
*The vector has size 238 but the DoFHandler object says that there are 274
degrees of freedom and there are 24 active cells. The size of your vector
needs to be either equal to the number of degrees of freedom, or equal to
the number of active cells.*
because my dof_handler has doc for pressure, size does not match...
how can I avoid this problem ?
--
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.