Dear Joel,

In general, no you cannot. Here you assume a specific ordering of the DoFs 
in the global solution vector, and this changes for each finite element. 
For example, using FE_Q's might result the ordering [x0 y0 z0 x1 y1 z1 ...] 
and FE_DGQ's might lead to [x0 x1 x2 ... y0 y1 y2 ... z0 z1 z2...]. To do 
this safely you should definitely do this using the built in tools, such as 
FEValueExtractors 
<https://www.dealii.org/developer/doxygen/deal.II/namespaceFEValuesExtractors.html>
 
or  fe.system_to_component_index 
<https://www.dealii.org/developer/doxygen/deal.II/classFiniteElement.html#a27220a135402b96c7e6eecbb04acda56>.
 
See the tutorials for examples on how to use them.

Regards,
J-P

On Tuesday, September 13, 2016 at 2:05:56 PM UTC+2, Joel Davidsson wrote:
>
> Dear all,
>
> If I have a Vector-valued solution in 3D and I would like to get the size 
> of each vector in another smaller Vector, can I do it like this:
>
> for (unsigned int i=0;i<solution.size();i=i+3)
> {
>  unsigned int x=i;
>  unsigned int y=i+1;
>  unsigned int z=i+2;
>  unsigned int j = i / 3.0;
>
>  size(j)=pow(solution(x),2.0)+pow(solution(y),2.0)+pow(solution(z),2.0);
> }
>
> Where solution is the vector valued Vector that is 3 times bigger that the 
> scalar Vector called size. 
>
> Best,
> Joel
>
> On Tuesday, August 30, 2016 at 6:48:46 PM UTC+2, Jean-Paul Pelteret wrote:
>>
>> Dear Joel,
>>
>> No, this functionality doesn't exist. However, you can create a 
>> quadrature rule with the quadrature points located at the support points 
>> (see the FAQ 
>> <https://github.com/dealii/dealii/wiki/Frequently-Asked-Questions#how-to-get-the-mapped-position-of-support-points-of-my-element>),
>>  
>> and could then extract the solution based on your custom quadrature rule.
>>
>> Regards,
>> J-P
>>
>> On Tuesday, August 30, 2016 at 5:32:13 PM UTC+2, Joel Davidsson wrote:
>>>
>>> Dear Daniel,
>>>
>>> I was just wondering if their exist a similar function 
>>> to DoFTools::map_dofs_to_support_points, but that maps a vector-valued 
>>> solution Vector into to a std::vector<Point<dim>> instead?
>>>
>>> Best,
>>> Joel
>>>
>>> On Monday, August 22, 2016 at 9:03:22 PM UTC+2, Daniel Arndt wrote:
>>>>
>>>> Joel,
>>>>
>>>> You should not rely on any particular sorting of the dofs in the 
>>>> solution vector. Instead you can ask the FiniteElement on each cell for 
>>>> the 
>>>> component of a local dof by
>>>> const unsigned int component = fe.system_to_component_index(i).first; 
>>>>
>>>> Apart from that the DataOut object knows about components and you can 
>>>> do some postprocessing on the solution vector using DataPostprocessor [1]
>>>>
>>>> Best,
>>>> Daniel
>>>>
>>>>  [1] 
>>>> https://dealii.org/8.4.1/doxygen/deal.II/classDataPostprocessor.html
>>>>
>>>> Am Montag, 22. August 2016 17:05:08 UTC+2 schrieb Joel Davidsson:
>>>>>
>>>>> Dear Daniel,
>>>>>
>>>>> Thank you for a very good answer, adding the fe_values_scalar and 
>>>>> cell_scalar fixed the problem.
>>>>>
>>>>> I have a follow-up question about the solution I get out, how is the 
>>>>> data organized in the solution vector? Say for example I want to loop 
>>>>> over 
>>>>> all the x components, how would I do that?
>>>>>
>>>>> Best,
>>>>> Joel
>>>>>
>>>>> On Friday, August 12, 2016 at 5:35:07 PM UTC+2, Daniel Arndt wrote:
>>>>>>
>>>>>> Joel,
>>>>>>
>>>>>> Yes the matrix you are assembling is a vector-valued mass matrix now.
>>>>>> For me your code is failing with
>>>>>>
>>>>>> void dealii::FEValuesBase<dim, 
>>>>>> spacedim>::get_function_gradients(const InputVector&, 
>>>>>> std::vector<dealii::Tensor<1, spacedim, typename 
>>>>>> InputVector::value_type> 
>>>>>> >&) const [with InputVector = dealii::Vector<double>; int dim = 3; int 
>>>>>> spacedim = 3; typename InputVector::value_type = double]
>>>>>> The violated condition was: 
>>>>>>     (fe->n_components()) == (1)
>>>>>> The name and call sequence of the exception was:
>>>>>>     dealii::ExcDimensionMismatch((fe->n_components()),(1))
>>>>>> Additional Information: 
>>>>>> Dimension 3 not equal to 1
>>>>>>
>>>>>> and this is not surprising. What you are doing is to extract the 
>>>>>> gradients of a vector-valued finite element solution. The object that 
>>>>>> should store this should therefore be a 
>>>>>> std::vector<std::vector<Tensor<1,dim>> as you want to store for each 
>>>>>> quadrature point and each component a Tensor<1,dim>.
>>>>>>
>>>>>> What you want to do is really that in "Do not work". As you have two 
>>>>>> DoFHandlers, you should also have two FEValues objects and two 
>>>>>> cell_iterator corresponding to the correct DoFHandler. Then you can 
>>>>>> extract 
>>>>>> the gradient of your scalar field and project it onto the ansatz space 
>>>>>> for 
>>>>>> the vector-valued field. This should look like:
>>>>>>
>>>>>> typename DoFHandler<dim>::active_cell_iterator cell_scalar = 
>>>>>> dof_handler_scalar.begin_active();
>>>>>> typename DoFHandler<dim>::active_cell_iterator cell = 
>>>>>> dof_handler.begin_active();
>>>>>> ...
>>>>>>
>>>>>>  for (; cell!=endc; ++cell, ++cell_vector)
>>>>>>     {
>>>>>>       fe_values.reinit (cell);
>>>>>>       fe_values_scalar.reinit (cell_scalar);
>>>>>>       fe_values_scalar.get_function_gradients(test,fe_function_grad);
>>>>>>       ...
>>>>>>     }
>>>>>>
>>>>>> Best,
>>>>>> Daniel
>>>>>>
>>>>>

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