<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Mon, Sep 29, 2014 at 4:55 PM, Miguel Angel Salazar de Troya <span dir="ltr"><<a href="mailto:salazardetroya@gmail.com" target="_blank">salazardetroya@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr">Hi all<div><br></div><div>I'm bumping this post because I have more questions related to the same problem.</div><div><div><br></div><div>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.</div></div></div></blockquote><div><br></div><div>This is unfortunately not explained well enough in our documentation. I will try to put a more substantial explanation in the manual.</div><div><br></div><div>PETSc operates in two different spaces: local space and global space. In the global space, every degree of freedom (dof, or unknown)</div><div>is owned by a single process and is a single entry which lives on that process in a global Vec. This is what solvers, like GMRES, deal with</div><div>and also what is plotted for visualization. In the local space, each process has all of the dofs it owns AND all of those that it shares with</div><div>other owning processes (ghosts). Thus in the local space a given dof can be entries on several processes in a local Vec. You usually use</div><div>the local space to compute residuals/Jacobians.</div><div><br></div><div>Our paradigm is: compute residual in the local space, so that you can set values for ghosts, and then call DMLocalToGlobal() to put those</div><div>values into a global Vec (using some combination operator, like +) which can be handed to the solver.</div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div><div>Thanks</div><div>Miguel</div></div></div><div class="gmail_extra"><br><div class="gmail_quote">On Fri, Sep 26, 2014 at 12:33 PM, Miguel Angel Salazar de Troya <span dir="ltr"><<a href="mailto:salazardetroya@gmail.com" target="_blank">salazardetroya@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr">I understand. Thanks a lot.<div><br></div><div>Miguel</div></div><div class="gmail_extra"><div><div><br><div class="gmail_quote">On Fri, Sep 26, 2014 at 10:53 AM, Abhyankar, Shrirang G. <span dir="ltr"><<a href="mailto:abhyshr@mcs.anl.gov" target="_blank">abhyshr@mcs.anl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">What Matt is saying is that there are two interfaces in PETSc for setting<br>
the residual evaluation routine:<br>
<br>
i) SNESSetFunction takes in a function pointer for the residual evaluation<br>
routine that has the prototype<br>
PetscErrorCode xyzroutine(SNES snes, Vec X, Vec F, void*<br>
ctx);<br>
<br>
X and F are the "global" solution and residual vectors. To compute the<br>
global residual evaluation, typically one does -- (a) scattering X and F<br>
onto local vectors localX and localF (DMGlobalToLocal), (b) computing the<br>
local residual, and (c) gathering the localF in the global F<br>
(DMLocalToGlobal). This is what is done in the example.<br>
<br>
ii) DMSNESSetFunctionLocal takes in a function pointer for the residual<br>
evaluation routine that has the prototype<br>
PetscErrorCode xyzlocalroutine(DM, Vec localX, localF,<br>
void* ctx)<br>
<br>
In this case, the localX and localF get passed to the routine. So, you<br>
only have to do the local residual evaluation. PETSc does the<br>
LocalToGlobal gather to form the global residual.<br>
<br>
I chose to use SNESSetFunction in the example. You can use either of them.<br>
<span><br>
Shri<br>
<br>
From: Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
</span>Date: Fri, 26 Sep 2014 10:28:26 -0500<br>
<span>To: Miguel Angel Salazar de Troya <<a href="mailto:salazardetroya@gmail.com" target="_blank">salazardetroya@gmail.com</a>><br>
Cc: Jed Brown <<a href="mailto:jed@jedbrown.org" target="_blank">jed@jedbrown.org</a>>, Shri <<a href="mailto:abhyshr@mcs.anl.gov" target="_blank">abhyshr@mcs.anl.gov</a>>,<br>
"<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>" <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>><br>
Subject: Re: [petsc-users] DMPlex with spring elements<br>
<br>
<br>
</span><div><div>>On Fri, Sep 26, 2014 at 10:26 AM, Miguel Angel Salazar de Troya<br>
><<a href="mailto:salazardetroya@gmail.com" target="_blank">salazardetroya@gmail.com</a>> wrote:<br>
><br>
>Yeah, but doesn't it only work with the local vectors localX and localF?<br>
><br>
><br>
><br>
>I am telling you what the interface for the functions is. You can do<br>
>whatever you want inside.<br>
><br>
> Matt<br>
><br>
><br>
>Miguel<br>
><br>
>On Fri, Sep 26, 2014 at 10:10 AM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
>wrote:<br>
><br>
>On Fri, Sep 26, 2014 at 10:06 AM, Miguel Angel Salazar de Troya<br>
><<a href="mailto:salazardetroya@gmail.com" target="_blank">salazardetroya@gmail.com</a>> wrote:<br>
><br>
>That means that if we call SNESSetFunction() we don't build the residual<br>
>vector in parallel? In the pflow example<br>
>(<a href="http://www.mcs.anl.gov/petsc/petsc-as/petsc-current/src/snes/examples/tut" target="_blank">http://www.mcs.anl.gov/petsc/petsc-as/petsc-current/src/snes/examples/tut</a><br>
>orials/network/pflow/pf.c.html) the function FormFunction() (Input for<br>
>SNESSetFunction() works with the local vectors. I don't understand this.<br>
><br>
><br>
><br>
>FormFunction() in that link clearly takes in a global vector X and<br>
>returns a global vector F. Inside, it<br>
>converts them to local vectors. This is exactly what you would do for a<br>
>function given to SNESSetFunction().<br>
><br>
> Matt<br>
><br>
><br>
><br>
>Thanks<br>
>Miguel<br>
><br>
><br>
>On Fri, Sep 26, 2014 at 9:34 AM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
>wrote:<br>
><br>
>On Fri, Sep 26, 2014 at 9:31 AM, Miguel Angel Salazar de Troya<br>
><<a href="mailto:salazardetroya@gmail.com" target="_blank">salazardetroya@gmail.com</a>> wrote:<br>
><br>
>Thanks. I had another question about the DM and SNES and TS. There are<br>
>similar routines to assign the residual and jacobian evaluation to both<br>
>objects. For the SNES case are:<br>
>DMSNESSetFunctionLocal<br>
>DMSNESSetJacobianLocal<br>
><br>
>What are the differences of these with:<br>
><br>
>SNESSetFunction<br>
>SNESSetJacobian<br>
><br>
><br>
><br>
><br>
>SNESSetFunction() expects the user to construct the entire parallel<br>
>residual vector. DMSNESSetFunctionLocal()<br>
>expects the user to construct the local pieces of the residual, and then<br>
>it automatically calls DMLocalToGlobal()<br>
>to assembly the full residual. It also converts the input from global<br>
>vectors to local vectors, and in the case of<br>
>DMDA multidimensional arrays.<br>
><br>
> Thanks,<br>
><br>
> Matt<br>
><br>
><br>
>and when should we use each? With "Local", it is meant to evaluate the<br>
>function/jacobian for the elements in the local processor? I could get<br>
>the local edges in DMNetwork by calling DMNetworkGetEdgeRange?<br>
><br>
>Miguel<br>
><br>
><br>
>On Thu, Sep 25, 2014 at 5:17 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
>wrote:<br>
><br>
><br>
><br>
>On Thu, Sep 25, 2014 at 5:15 PM, Miguel Angel Salazar de Troya<br>
><<a href="mailto:salazardetroya@gmail.com" target="_blank">salazardetroya@gmail.com</a>> wrote:<br>
><br>
>> If you need a symmetric Jacobian, you can use the BC facility in<br>
>> PetscSection, which eliminates the<br>
>> variables completely. This is how the FEM examples, like ex12, work.<br>
>Would that be with PetscSectionSetConstraintDof ? For that I will need<br>
>the PetscSection, DofSection, within DMNetwork, how can I obtain it? I<br>
>could cast it to DM_Network from the dm, networkdm, declared in the main<br>
>program, maybe something like this:<br>
>DM_Network *network = (DM_Network*) networkdm->data;Then I would loop<br>
>over the vertices and call PetscSectionSetConstraintDof if it's a<br>
>boundary node (by checking the corresponding component)<br>
><br>
><br>
><br>
><br>
>I admit to not completely understanding DMNetwork. However, it eventually<br>
>builds a PetscSection for data layout, which<br>
>you could get from DMGetDefaultSection(). The right thing to do is find<br>
>where it builds the Section, and put in your BC<br>
>there, but that sounds like it would entail coding.<br>
><br>
> Thanks,<br>
><br>
> Matt<br>
><br>
><br>
><br>
><br>
>Thanks for your responses.Miguel<br>
><br>
><br>
><br>
><br>
>On Thu, Sep 25, 2014 at 2:42 PM, Jed Brown <<a href="mailto:jed@jedbrown.org" target="_blank">jed@jedbrown.org</a>> wrote:<br>
><br>
>Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> writes:<br>
><br>
>> On Thu, Sep 25, 2014 at 1:46 PM, Abhyankar, Shrirang G.<br>
>><<a href="mailto:abhyshr@mcs.anl.gov" target="_blank">abhyshr@mcs.anl.gov</a><br>
>>> wrote:<br>
>><br>
>>> You are right. The Jacobian for the power grid application is indeed<br>
>>> non-symmetric. Is that a problem for your application?<br>
>>><br>
>><br>
>> If you need a symmetric Jacobian, you can use the BC facility in<br>
>> PetscSection, which eliminates the<br>
>> variables completely. This is how the FEM examples, like ex12, work.<br>
><br>
>You can also use MatZeroRowsColumns() or do the equivalent<br>
>transformation during assembly (my preference).<br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>Miguel Angel Salazar de Troya<br>
><br>
><br>
>Graduate Research Assistant<br>
>Department of Mechanical Science and Engineering<br>
>University of Illinois at Urbana-Champaign<br>
</div></div>><a href="tel:%28217%29%20550-2360" value="+12175502360" target="_blank">(217) 550-2360</a> <tel:%28217%29%20550-2360><br>
<span>><br>
><br>
><a href="mailto:salaza11@illinois.edu" target="_blank">salaza11@illinois.edu</a><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>What most experimenters take for granted before they begin their<br>
>experiments is infinitely more interesting than any results to which<br>
>their experiments lead.<br>
>-- Norbert Wiener<br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>Miguel Angel Salazar de Troya<br>
>Graduate Research Assistant<br>
>Department of Mechanical Science and Engineering<br>
>University of Illinois at Urbana-Champaign<br>
</span>><a href="tel:%28217%29%20550-2360" value="+12175502360" target="_blank">(217) 550-2360</a> <tel:%28217%29%20550-2360><br>
<span>><a href="mailto:salaza11@illinois.edu" target="_blank">salaza11@illinois.edu</a><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>What most experimenters take for granted before they begin their<br>
>experiments is infinitely more interesting than any results to which<br>
>their experiments lead.<br>
>-- Norbert Wiener<br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>Miguel Angel Salazar de Troya<br>
>Graduate Research Assistant<br>
>Department of Mechanical Science and Engineering<br>
>University of Illinois at Urbana-Champaign<br>
</span>><a href="tel:%28217%29%20550-2360" value="+12175502360" target="_blank">(217) 550-2360</a> <tel:%28217%29%20550-2360><br>
<span>><a href="mailto:salaza11@illinois.edu" target="_blank">salaza11@illinois.edu</a><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>What most experimenters take for granted before they begin their<br>
>experiments is infinitely more interesting than any results to which<br>
>their experiments lead.<br>
>-- Norbert Wiener<br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>Miguel Angel Salazar de Troya<br>
>Graduate Research Assistant<br>
>Department of Mechanical Science and Engineering<br>
>University of Illinois at Urbana-Champaign<br>
</span>><a href="tel:%28217%29%20550-2360" value="+12175502360" target="_blank">(217) 550-2360</a> <tel:%28217%29%20550-2360><br>
<div><div>><a href="mailto:salaza11@illinois.edu" target="_blank">salaza11@illinois.edu</a><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
><br>
>--<br>
>What most experimenters take for granted before they begin their<br>
>experiments is infinitely more interesting than any results to which<br>
>their experiments lead.<br>
>-- Norbert Wiener<br>
<br>
</div></div></blockquote></div><br><br clear="all"><span class="HOEnZb"><font color="#888888"><div><br></div>-- <br></font></span></div></div><span class="HOEnZb"><font color="#888888"><div dir="ltr"><div><div><font face="verdana, sans-serif"><b>Miguel Angel Salazar de Troya</b></font></div></div><span><font color="#888888"><div><div><br><font face="arial,helvetica,sans-serif">Graduate Research Assistant<br>Department of Mechanical Science and Engineering<br></font>University of Illinois at Urbana-Champaign<br><a href="tel:%28217%29%20550-2360" value="+12175502360" target="_blank">(217) 550-2360</a><br>
</div></div><a href="mailto:salaza11@illinois.edu" target="_blank">salaza11@illinois.edu</a></font></span><div><br></div></div>
</font></span></div><span class="HOEnZb"><font color="#888888">
</font></span></blockquote></div><span class="HOEnZb"><font color="#888888"><br><br clear="all"><div><br></div>-- <br><div dir="ltr"><font face="verdana, sans-serif"><b>Miguel Angel Salazar de Troya</b></font><span><font color="#888888"><br><font face="arial,helvetica,sans-serif">Graduate Research Assistant<br>Department of Mechanical Science and Engineering<br></font>University of Illinois at Urbana-Champaign<br><a href="tel:%28217%29%20550-2360" value="+12175502360" target="_blank">(217) 550-2360</a><br>
<a href="mailto:salaza11@illinois.edu" target="_blank">salaza11@illinois.edu</a></font></span><div><br></div></div>
</font></span></div>
</blockquote></div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener
</div></div>