Dear Michal,

thanks for the use case. I have now created a pull request to the
deal.II master branch that should add the appropriate assertion:

https://github.com/dealii/dealii/pull/5212

The problem was that we previously set all FE components to the single
scalar vector component, effectively duplicating the components and
duplicating the results. We now correctly identify the "single" vector
given to the evaluation routine and set nullptr for the others which
allows us to catch unintended use.

Best,
Martin


On 08.10.2017 13:19, Michał Wichrowski wrote:
>
>
> W dniu niedziela, 8 października 2017 12:55:02 UTC+2 użytkownik Martin
> Kronbichler napisał:
>
>     Dear Michal,
>
>     Could you please be a bit more specific? What was the error that
>     you obtained? What vector did you pass into
>     FEEvaluation::read_dof_values/distribute_local_to_global? I guess
>     that you handed in a simple vector, e.g.
>     LinearAlgebra::distributed::Vector and got a segmentation fault.
>     But to add the assertion I need more details. The issues that in
>     principle, we do support using multiple components out of a scalar
>     FE - but in that case the user must hand in a BlockVector or
>     std::vector<VectorType>. We could catch that case in the
>     machinery, though.
>
>     Thanks!
>
>     Best,
>     Martin
>
> Martin,
> There was no segmentation fault, program runned without any errors. I
> realised that something is wrong when I outputed results.  This is
> mine local apply,  VectorType =LinearAlgebra::distributed::Vector,
> problem (Stokes) consists of dim components of velocity and 1
> component of pressure. 
> Instead of using FESystem(FE_Q, dim) I used FE_Q for velocity.
> After fixing this it work OK, so it is just a problem with assertion.
> I'll try to make minimal code reproducing error  (or rather lack of
> error).
>   
> template <int dim,int degree_u,  typename Number >
> void
> SymGradMass<dim, degree_u,Number>
> ::local_apply (const MatrixFree<dim,Number> &data,
>              VectorType          &dst,
>              const VectorType    &src,
>              const std::pair<unsigned int,unsigned int> &cell_range) const
> {
>
>   typedef VectorizedArray<Number> vector_t;
>   FEEvaluation<dim,degree_u,degree_u+1+0,dim,Number> velocity (data,
> 0);  //only velocity is used here 0-th component.
>
>   for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
>     {
>       velocity.reinit (cell);
>       velocity.read_dof_values (src);
>         integrate_A(velocity);
>       velocity.distribute_local_to_global (dst);
>
>     }
> }
>
> Michał
> -- 
> 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]
> <mailto:[email protected]>.
> 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 [email protected].
For more options, visit https://groups.google.com/d/optout.

Reply via email to