[petsc-users] DMPlex with spring elements

Abhyankar, Shrirang G. abhyshr at mcs.anl.gov
Tue Sep 30 11:22:43 CDT 2014


From:  Miguel Angel Salazar de Troya <salazardetroya at gmail.com>
Date:  Mon, 29 Sep 2014 16:55:14 -0500
To:  Shri <abhyshr at mcs.anl.gov>
Cc:  "petsc-users at mcs.anl.gov" <petsc-users at mcs.anl.gov>
Subject:  Re: [petsc-users] DMPlex with spring elements


>Hi all
>I'm bumping this post because I have more questions related to the same
>problem.
>
>I am looping over the edges of my DMNetwork, then I obtain the vertices
>that make up each edge with DMNetworkGetConnectedNode(). Each of these
>vertices have two variables (or actually, two degrees of freedom for my
>problem). My intentions are to modify the solution vector entries that
>are affected by these variables in each vertex. I would call the function
>DMNetworkGetVariableOffset() to do this. What happens if one of the
>vertices is a ghost vertex? Can I still modify the solution vector? My
>problem is that the edge has information to provide to these nodes.
>
>

Sorry for the delay. I think you would want to use
http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMNetworkIsG
hostVertex.html and not modify the value for the ghost vertex.

Shri

>
>Thanks
>Miguel
>
>
>
>On Fri, Sep 26, 2014 at 12:33 PM, Miguel Angel Salazar de Troya
><salazardetroya at gmail.com> wrote:
>
>I understand. Thanks a lot.
>Miguel
>
>
>On Fri, Sep 26, 2014 at 10:53 AM, Abhyankar, Shrirang G.
><abhyshr at mcs.anl.gov> wrote:
>
>What Matt is saying is that there are two interfaces in PETSc for setting
>the residual evaluation routine:
>
>i) SNESSetFunction takes in a function pointer for the residual evaluation
>routine that has the prototype
>                 PetscErrorCode xyzroutine(SNES snes, Vec X, Vec F, void*
>ctx);
>
>X and F are the "global" solution and residual vectors. To compute the
>global residual evaluation, typically one does -- (a) scattering X and F
>onto local vectors localX and localF (DMGlobalToLocal), (b) computing the
>local residual, and (c) gathering the localF in the global F
>(DMLocalToGlobal). This is what is done in the example.
>
>ii) DMSNESSetFunctionLocal takes in a function pointer for the residual
>evaluation routine that has the prototype
>                 PetscErrorCode xyzlocalroutine(DM, Vec localX, localF,
>void* ctx)
>
>In this case, the localX and localF get passed to the routine. So, you
>only have to do the local residual evaluation. PETSc does the
>LocalToGlobal gather to form the global residual.
>
>I chose to use SNESSetFunction in the example. You can use either of them.
>
>Shri
>
>From:  Matthew Knepley <knepley at gmail.com>
>Date:  Fri, 26 Sep 2014 10:28:26 -0500
>To:  Miguel Angel Salazar de Troya <salazardetroya at gmail.com>
>Cc:  Jed Brown <jed at jedbrown.org>, Shri <abhyshr at mcs.anl.gov>,
>"petsc-users at mcs.anl.gov" <petsc-users at mcs.anl.gov>
>Subject:  Re: [petsc-users] DMPlex with spring elements
>
>
>>On Fri, Sep 26, 2014 at 10:26 AM, Miguel Angel Salazar de Troya
>><salazardetroya at gmail.com> wrote:
>>
>>Yeah, but doesn't it only work with the local vectors localX and localF?
>>
>>
>>
>>I am telling you what the interface for the functions is. You can do
>>whatever you want inside.
>>
>>  Matt
>>
>>
>>Miguel
>>
>>On Fri, Sep 26, 2014 at 10:10 AM, Matthew Knepley <knepley at gmail.com>
>>wrote:
>>
>>On Fri, Sep 26, 2014 at 10:06 AM, Miguel Angel Salazar de Troya
>><salazardetroya at gmail.com> wrote:
>>
>>That means that if we call SNESSetFunction() we don't build the residual
>>vector in parallel? In the pflow example
>>(http://www.mcs.anl.gov/petsc/petsc-as/petsc-current/src/snes/examples/tu
>>t
>>orials/network/pflow/pf.c.html) the function FormFunction() (Input for
>>SNESSetFunction() works with the local vectors. I don't understand this.
>>
>>
>>
>>FormFunction() in that link clearly takes in a global vector X and
>>returns a global vector F. Inside, it
>>converts them to local vectors. This is exactly what you would do for a
>>function given to SNESSetFunction().
>>
>>  Matt
>>
>>
>>
>>Thanks
>>Miguel
>>
>>
>>On Fri, Sep 26, 2014 at 9:34 AM, Matthew Knepley <knepley at gmail.com>
>>wrote:
>>
>>On Fri, Sep 26, 2014 at 9:31 AM, Miguel Angel Salazar de Troya
>><salazardetroya at gmail.com> wrote:
>>
>>Thanks. I had another question about the DM and SNES and TS. There are
>>similar routines to assign the residual and jacobian evaluation to both
>>objects. For the SNES case are:
>>DMSNESSetFunctionLocal
>>DMSNESSetJacobianLocal
>>
>>What are the differences of these with:
>>
>>SNESSetFunction
>>SNESSetJacobian
>>
>>
>>
>>
>>SNESSetFunction() expects the user to construct the entire parallel
>>residual vector. DMSNESSetFunctionLocal()
>>expects the user to construct the local pieces of the residual, and then
>>it automatically calls DMLocalToGlobal()
>>to assembly the full residual. It also converts the input from global
>>vectors to local vectors, and in the case of
>>DMDA multidimensional arrays.
>>
>>  Thanks,
>>
>>    Matt
>>
>>
>>and when should we use each? With "Local", it is meant to evaluate the
>>function/jacobian for the elements in the local processor? I could get
>>the local edges in DMNetwork by calling DMNetworkGetEdgeRange?
>>
>>Miguel
>>
>>
>>On Thu, Sep 25, 2014 at 5:17 PM, Matthew Knepley <knepley at gmail.com>
>>wrote:
>>
>>
>>
>>On Thu, Sep 25, 2014 at 5:15 PM, Miguel Angel Salazar de Troya
>><salazardetroya at gmail.com> wrote:
>>
>>> If you need a symmetric Jacobian, you can use the BC facility in
>>> PetscSection, which eliminates the
>>> variables completely. This is how the FEM examples, like ex12, work.
>>Would that be with PetscSectionSetConstraintDof ? For that I will need
>>the PetscSection, DofSection, within DMNetwork, how can I obtain it? I
>>could cast it to DM_Network from the dm, networkdm,  declared in the main
>>program, maybe something like this:
>>DM_Network     *network = (DM_Network*) networkdm->data;Then I would loop
>>over the vertices and call PetscSectionSetConstraintDof if it's a
>>boundary node (by checking the corresponding component)
>>
>>
>>
>>
>>I admit to not completely understanding DMNetwork. However, it eventually
>>builds a PetscSection for data layout, which
>>you could get from DMGetDefaultSection(). The right thing to do is find
>>where it builds the Section, and put in your BC
>>there, but that sounds like it would entail coding.
>>
>>  Thanks,
>>
>>     Matt
>>
>>
>>
>>
>>Thanks for your responses.Miguel
>>
>>
>>
>>
>>On Thu, Sep 25, 2014 at 2:42 PM, Jed Brown <jed at jedbrown.org> wrote:
>>
>>Matthew Knepley <knepley at gmail.com> writes:
>>
>>> On Thu, Sep 25, 2014 at 1:46 PM, Abhyankar, Shrirang G.
>>><abhyshr at mcs.anl.gov
>>>> wrote:
>>>
>>>> You are right. The Jacobian for the power grid application is indeed
>>>> non-symmetric. Is that a problem for your application?
>>>>
>>>
>>> If you need a symmetric Jacobian, you can use the BC facility in
>>> PetscSection, which eliminates the
>>> variables completely. This is how the FEM examples, like ex12, work.
>>
>>You can also use MatZeroRowsColumns() or do the equivalent
>>transformation during assembly (my preference).
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>Miguel Angel Salazar de Troya
>>
>>
>>Graduate Research Assistant
>>Department of Mechanical Science and Engineering
>>University of Illinois at Urbana-Champaign
>
>
>>(217) 550-2360 <tel:%28217%29%20550-2360>
>>
>>
>>salaza11 at illinois.edu
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>What most experimenters take for granted before they begin their
>>experiments is infinitely more interesting than any results to which
>>their experiments lead.
>>-- Norbert Wiener
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>Miguel Angel Salazar de Troya
>>Graduate Research Assistant
>>Department of Mechanical Science and Engineering
>>University of Illinois at Urbana-Champaign
>>(217) 550-2360 <tel:%28217%29%20550-2360>
>>salaza11 at illinois.edu
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>What most experimenters take for granted before they begin their
>>experiments is infinitely more interesting than any results to which
>>their experiments lead.
>>-- Norbert Wiener
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>Miguel Angel Salazar de Troya
>>Graduate Research Assistant
>>Department of Mechanical Science and Engineering
>>University of Illinois at Urbana-Champaign
>>(217) 550-2360 <tel:%28217%29%20550-2360>
>>salaza11 at illinois.edu
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>What most experimenters take for granted before they begin their
>>experiments is infinitely more interesting than any results to which
>>their experiments lead.
>>-- Norbert Wiener
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>Miguel Angel Salazar de Troya
>>Graduate Research Assistant
>>Department of Mechanical Science and Engineering
>>University of Illinois at Urbana-Champaign
>>(217) 550-2360 <tel:%28217%29%20550-2360>
>>salaza11 at illinois.edu
>>
>>
>>
>>
>>
>>
>>
>>
>>
>>--
>>What most experimenters take for granted before they begin their
>>experiments is infinitely more interesting than any results to which
>>their experiments lead.
>>-- Norbert Wiener
>
>
>
>
>
>
>
>
>
>-- 
>
>
>Miguel Angel Salazar de Troya
>
>
>Graduate Research Assistant
>Department of Mechanical Science and Engineering
>University of Illinois at Urbana-Champaign
>(217) 550-2360
>
>
>salaza11 at illinois.edu
>
>
>
>
>
>
>
>
>-- 
>Miguel Angel Salazar de Troya
>Graduate Research Assistant
>Department of Mechanical Science and Engineering
>University of Illinois at Urbana-Champaign
>(217) 550-2360
>salaza11 at illinois.edu



More information about the petsc-users mailing list