Hi Vachan,

I don't think you need to use RPE plane. My guess is that you can use 
VectorTools::point_values(), which internally uses PRE (see 
https://github.com/dealii/dealii/blob/e413cbf8ac77606cf4a04e1e2fa75232c08533b4/include/deal.II/numerics/vector_tools_evaluate.h#L230-L343).
 
Here you could collect the support points from the new DoFHandler, let 
VT::point_values() do the work and finally set the values at the support 
points.

Further resources are:
- the test 
folder 
https://github.com/dealii/dealii/tree/master/tests/remote_point_evaluation 
(in 
particular 
https://github.com/dealii/dealii/blob/master/tests/remote_point_evaluation/vector_tools_evaluate_at_points_01.cc;
 
here, the solution is interpolated between non-matching grids)
- the usage of RPE and VT::point_values() in 
DataOutResample 
(https://github.com/dealii/dealii/blob/master/source/numerics/data_out_resample.cc)
 
- the usage in the two-phase solver adaflo 
(https://github.com/kronbichler/adaflo/blob/master/include/adaflo/sharp_interface_util.h)

PM

On Tuesday, 3 August 2021 at 09:03:02 UTC+2 vachanpo...@gmail.com wrote:

> Dr. Wolfgang,
>
> Thank you for still taking interest in this thread :)
>
> This becomes a very difficult data transfer problem because you want to
>> evaluate the solution at a point that the current process may know nothing
>> about. In essence, you will have to ask all of the other processes about 
>> who
>> owns a part of the mesh on which a given interpolation point is located, 
>> and
>> then that process has to send you the value of the function at that 
>> point. You
>> have to do that for all points. This is going to be an expensive 
>> operation.
>> You might want to check out the Utilities::MPI::RemotePointEvaluation 
>> class in
>> the latest (9.3) release for help with this.
>
> I understand that the approach I had mentioned is inefficient. I only need 
> to do this once, to start a simulation, so I thought it might be okay.
>
> I had a look at Utilities::MPI::RemotePointEvaluation. This is how I 
> think it can be used (correct me if wrong)
>
>    1. Call reinit() with the dof locations, along with the corresponding 
>    triangulation and mapping, where I want to evaluate FEFieldFunction (I 
>    suppose the triangulation and mapping here are not of the 
> FEFieldFunction?).
>    2. Call evaluate_and_process() by passing in the FEFieldFunction.
>
> I have a few questions.
>
>    1. What is the "buffer" argument in evaluate_and_process()?
>    2. The documentation for this function says get_point_ptrs() must be 
>    used to "process" the output in case the point-cell map is not one-one. I 
>    will surely encounter such cases. How can I use the data returned by 
>    get_point_ptrs() and how exactly should I "process" the output?
>
> I couldn't find this used in any examples. Any clarifications would be 
> greatly helpful.
>
> On Tue, 3 Aug 2021 at 04:47, Wolfgang Bangerth <bang...@colostate.edu> 
> wrote:
>
>> On 7/5/21 11:04 PM, vachan potluri wrote:
>> > 
>> >     So is your fine mesh a refinement of the coarse one? If not, you 
>> may want to
>> >     look at FEFieldFunction.
>> > 
>> > Yes, it is. But the "refinement" is done by the meshing 
>> software, outside 
>> > dealii. Is there any simplification possible in such a case?
>>
>> Not really. If deal.II has no knowledge about the relationship between 
>> cells 
>> on the two meshes, then it is in essence the interpolation from one 
>> unstructured grid to another unstructured grid.
>>
>>
>> > Otherwise, I think FEFieldFunction would be a safe choice. Since I want 
>> to use 
>> > it with p::d::Triangulation, will setting all dofs as relevant do the 
>> job? 
>> > This way the partitioning can also be different.
>>
>> This becomes a very difficult data transfer problem because you want to 
>> evaluate the solution at a point that the current process may know 
>> nothing 
>> about. In essence, you will have to ask all of the other processes about 
>> who 
>> owns a part of the mesh on which a given interpolation point is located, 
>> and 
>> then that process has to send you the value of the function at that 
>> point. You 
>> have to do that for all points. This is going to be an expensive 
>> operation.
>>
>> You might want to check out the Utilities::MPI::RemotePointEvaluation 
>> class in 
>> the latest (9.3) release for help with this.
>>
>> Best
>>   W.
>>
>>
>> -- 
>> ------------------------------------------------------------------------
>> Wolfgang Bangerth          email:                 bang...@colostate.edu
>>                             www: http://www.math.colostate.edu/~bangerth/
>>
>> -- 
>> 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+un...@googlegroups.com.
>> To view this discussion on the web visit 
>> https://groups.google.com/d/msgid/dealii/61601d83-856e-dd27-dc78-3653bd45ec2f%40colostate.edu
>> .
>>
>

-- 
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.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/5032c623-7ba0-4965-b62d-e55d06f045b4n%40googlegroups.com.

Reply via email to