Dear Hermes,

Did you try VectorTools::point_values, https://www.dealii.org/current/doxygen/deal.II/namespaceVectorTools.html#a1ad6eceb6cbeaa505baf7f938289bbde ?

As an alternative, you can also check if the point is in the locally owned cell of the respective cell, and run the evaluation only on that MPI rank. You can then use an MPI collective, like MPI_Bcast or deal.II's version Utilities::MPI::broadcast, to get the solution to all MPI ranks (typically with another collective to make all ranks agree on who sends the information).

Best,
Martin

This is a variant that allows the evaluation even on remote processes.

On 24.08.21 14:05, Hermes Sampedro wrote:
Dear all,

I am trying to extend my implementation with parallel computing. I would like to get the solution in a single point when having a distributed system. For example, in step-40 I have the /locally_relevant_solution/ at the output_result(...) function.
In a non distributed implementation I used:

Vector<*double*> vecSol (fe.n_components());

VectorTools::point_value(dof_handler, solution,Point<2>(0.2, 0.2),vecSol);


What is the easiest way to do that? Is it possible to check if the point is inside the /locally_relevant_solution? /Or would be more appropriate to create a copy of the solution in a global vector and get the point from there?


Thank you

H


El mié, 11 ago 2021 a las 10:12, Hermes Sampedro (<[email protected] <mailto:[email protected]>>) escribió:

    Thank you very much, I could make it work.

    Regards,
    Hermes

    El mar, 10 ago 2021 a las 19:46, Jean-Paul Pelteret
    (<[email protected] <mailto:[email protected]>>) escribió:

        Another thing: You probably also need to initialise with the
        right number of components. So something like
        Vector<double> vecSol (fe.n_components());


        On 10. Aug 2021, at 19:40, Jean-Paul Pelteret
        <[email protected] <mailto:[email protected]>> wrote:

        Hi Hermes,

        You don’t say what errors you’re seeing, but my guess is that
        it now doesn’t compile. This variant of the function (the one
        that Daniel linked to) returns void, so you should call it
        before outputting the result:

        Vector<double> vecSol;
        VectorTools::point_value(dof_handler,
        solution,Point<2>(0.2, 0.2),vecSol)
        std::cout << "Solution at (0.2,0.2): "<< vecSol << std::endl;

        This will hopefully get you what you want.

        Best,
        Jean-Paul


        On 10. Aug 2021, at 16:09, Hermes Sampedro
        <[email protected] <mailto:[email protected]>>
        wrote:

        Thank you for your answer. I am not sure if I fully
        understand your suggestion. Do you mean something like that:

        Vector<*double*> vecSol;
        std::cout <<"Solution at (0.2,0.2): "<<
        VectorTools::point_value(dof_handler,
        solution,Point<2>(0.2,0.2),vecSol)<< std::endl;


        I still get some errors. Is there not any way to get for
        example the real part of/solution/ easily and use it
        directly on the point_value as in step-3?

        Thank you
        H

        El mar, 10 ago 2021 a las 15:49, Daniel Arndt
        (<[email protected] <mailto:[email protected]>>)
        escribió:

            Hermes,

            Use another overload. The one returning the solution as
            a parameter should
            
work:https://www.dealii.org/current/doxygen/deal.II/namespaceVectorTools.html#acd358e9b110ccbf4a7f76796d206b9c7
            
<https://www.dealii.org/current/doxygen/deal.II/namespaceVectorTools.html#acd358e9b110ccbf4a7f76796d206b9c7>

            Best,
            Daniel

            Am Di., 10. Aug. 2021 um 09:41 Uhr schrieb Hermes
            Sampedro <[email protected]
            <mailto:[email protected]>>:

                Dear all,

                It is explained in Step-3 how to evaluate the
                solution in a point. I am trying to do the same for
                Step-29, to evaluate the real and imaginary parts
                separately in a single point:

                /std::cout << "Solution at (0.2,0.2): "<<
                VectorTools::point_value(dof_handler,
                solution,Point<2>(0.2, 0.2))<< std::endl;/

                I do not have any problems compiling however, an
                error occurs when running:

                /The violated condition
                was:  dof.get_fe(0).n_components() == 1/

                What is the proper way to call the real and
                imaginary parts of the solution at a particular
                point here?


                Thank you very much!

                H.


                --
                The deal.II project is located
                athttp://www.dealii.org/ <http://www.dealii.org/>
                For mailing list/forum options,
                seehttps://groups.google.com/d/forum/dealii?hl=en
                <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
                [email protected]
                <mailto:[email protected]>.
                To view this discussion on the web
                
visithttps://groups.google.com/d/msgid/dealii/54a2d44e-194a-43c0-aa7f-9fddd5bd0dean%40googlegroups.com
                
<https://groups.google.com/d/msgid/dealii/54a2d44e-194a-43c0-aa7f-9fddd5bd0dean%40googlegroups.com?utm_medium=email&utm_source=footer>.


            --
            The deal.II project is located athttp://www.dealii.org/
            <http://www.dealii.org/>
            For mailing list/forum options,
            seehttps://groups.google.com/d/forum/dealii?hl=en
            <https://groups.google.com/d/forum/dealii?hl=en>
            ---
            You received this message because you are subscribed to
            a topic in the Google Groups "deal.II User Group" group.
            To unsubscribe from this topic,
            
visithttps://groups.google.com/d/topic/dealii/Ru1_uMbix30/unsubscribe
            <https://groups.google.com/d/topic/dealii/Ru1_uMbix30/unsubscribe>.
            To unsubscribe from this group and all its topics, send
            an email [email protected]
            <mailto:[email protected]>.
            To view this discussion on the web
            
visithttps://groups.google.com/d/msgid/dealii/CAOYDWb%2BR0bRTVxB-F22Hwk_ZexN4%3DVwCbZ7ZKqs3JjEaHFhncA%40mail.gmail.com
            
<https://groups.google.com/d/msgid/dealii/CAOYDWb%2BR0bRTVxB-F22Hwk_ZexN4%3DVwCbZ7ZKqs3JjEaHFhncA%40mail.gmail.com?utm_medium=email&utm_source=footer>.


        --
        The deal.II project is located athttp://www.dealii.org/
        <http://www.dealii.org/>
        For mailing list/forum options,
        seehttps://groups.google.com/d/forum/dealii?hl=en
        <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 [email protected]
        <mailto:[email protected]>.
        To view this discussion on the web
        
visithttps://groups.google.com/d/msgid/dealii/CAB%3DnHhZrnXUeYt5FrKa66oTwoV8WJ8b%3DFqLFJ9pfP%3DUHRTyV%2BQ%40mail.gmail.com
        
<https://groups.google.com/d/msgid/dealii/CAB%3DnHhZrnXUeYt5FrKa66oTwoV8WJ8b%3DFqLFJ9pfP%3DUHRTyV%2BQ%40mail.gmail.com?utm_medium=email&utm_source=footer>.


-- The deal.II project is located at http://www.dealii.org/
        <http://www.dealii.org/>
        For mailing list/forum options, see
        https://groups.google.com/d/forum/dealii?hl=en
        <https://groups.google.com/d/forum/dealii?hl=en>
        ---
        You received this message because you are subscribed to a
        topic in the Google Groups "deal.II User Group" group.
        To unsubscribe from this topic, visit
        https://groups.google.com/d/topic/dealii/Ru1_uMbix30/unsubscribe
        <https://groups.google.com/d/topic/dealii/Ru1_uMbix30/unsubscribe>.
        To unsubscribe from this group and all its topics, send an
        email to [email protected]
        <mailto:[email protected]>.
        To view this discussion on the web visit
        
https://groups.google.com/d/msgid/dealii/5DDE5428-30E3-48CC-9924-333605B3DAAA%40gmail.com
        
<https://groups.google.com/d/msgid/dealii/5DDE5428-30E3-48CC-9924-333605B3DAAA%40gmail.com?utm_medium=email&utm_source=footer>.

--
The deal.II project is located at http://www.dealii.org/ <http://www.dealii.org/> For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en <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]>. To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/CAB%3DnHhbh-ir_njBtYjvVafJ-Hyuhfwck5XK1yeHaxY%2BrNya03g%40mail.gmail.com <https://groups.google.com/d/msgid/dealii/CAB%3DnHhbh-ir_njBtYjvVafJ-Hyuhfwck5XK1yeHaxY%2BrNya03g%40mail.gmail.com?utm_medium=email&utm_source=footer>.

--
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].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/031efb5b-8628-4baa-da19-8869ced63165%40gmail.com.

Reply via email to