>
> [...] 
>
1.  When evaluating point is on the vertex of a cell, this means that there 
> should be several cells in 3D have this vertex. 
> Thus, my first question is the function :
>
> const std::pair<typename DoFHandler<dim>::active_cell_iterator, Point<dim>>
>             cell_point = GridTools::find_active_cell_around_point(mapping, 
> dof_handler, point); 
>
> I personally think that as I set it to be const, so the return cell 
> (cell_point.first) is only one cell, rather than a group of cells. Could you 
> please tell me about this?
>
>
This function gives you one cell that includes the given point, not all of 
them. In case, you are using continuous elements, the value you obtain 
should not depend on the choice of this cell. Otherwise, you have to be 
more specific which cell you want.
 

> 2. My second idea is with the codes:
> const Quadrature<dim> 
> quadrature(GeometryInfo<dim>::project_to_unit_cell(cell_point.second));
>     FEValues<dim> fe_values(mapping, dof_handler.get_fe(), quadrature,
>                             update_values | update_gradients | 
> update_quadrature_points);
>     fe_values.reinit(cell_point.first);
>
> I have only one quadrature defined on the cell, and the printing result show 
> me the same thing. 
> Can I use this way to get the solution, or it would be more efficient when 
> using codes like:
> QGaussLobatto<dim> quadrature_formula(2);
> FEValues<dim> fe_values(fe, quadrature_formula,
>                         update_values | update_gradients |
>                         update_quadrature_points | update_JxW_values);
> fe_values.reinit(cell_point.first);
>
> And then finding the quadrature point coincide with the evaluating point to 
> get the solution?
>
>
If you know beforehand that your evaluation coincides with a unit support 
point on a known cell, you can save the call to 
GridTools::find_active_cell_around_point().
Otherwise, you likely can't do much better than using its return value.
 

> 3. My third thought is in my function I have only 1 quadrature point so that 
> temp_solution_re[0] and temp_solution_im[0] are the variables containing the 
> solution. 
>
> The printing results showed that the solution from point_value() is identical 
> with those from temp_solution_re[0] and temp_solution_im[0].
>
> Do temp_solution_re[1], temp_solution_re[2], 
> temp_solution_im[1],temp_solution_im[2] make any sense for my solution? 
> Because I saw the code can run with these variables.
>
>
VectorTools::point_value() does more or less the same you are doing above 
manually, i.e. calling GridTools::find_active_cell_around_point() and then 
initializing a FEValues object for evaluating the given finite element 
vector. Hence, it is not surprising if the results coincide. Since the size 
of temp_solution_re and temp_solution_im is 1, it doesn't make sense to 
access the values in the second and third entry.

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