<div dir="ltr"><div dir="ltr">On Tue, Feb 18, 2020 at 8:20 AM Yuyun Yang <<a href="mailto:yyang85@stanford.edu">yyang85@stanford.edu</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">Thanks for the clarification.<br>
<br>
Got one more question: if I have variable coefficients, my stencil will be updated at every time step, so will the coefficients in myMatMult. In that case, is it necessary to destroy the shell matrix and create it all over again, or can I use it as it is, only calling the stencil update function, assuming the result will be passed into the matrix operation automatically?<br></blockquote><div><br></div><div>You update the information in the context associated with the shell matrix. No need to destroy it.</div><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
Thanks,<br>
Yuyun<br>
<br>
On 2/18/20, 7:34 AM, "Smith, Barry F." <<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>> wrote:<br>
<br>
<br>
<br>
    > On Feb 17, 2020, at 7:56 AM, Yuyun Yang <<a href="mailto:yyang85@stanford.edu" target="_blank">yyang85@stanford.edu</a>> wrote:<br>
    > <br>
    > Hello,<br>
    > <br>
    > I actually have a question about the usage of DMDA since I'm quite new to this. I wonder if the DMDA suite of functions can be directly called on vectors created from VecCreate?<br>
<br>
       Yes, but you have to make sure the ones you create have the same sizes and parallel layouts. Generally best to get them from the DMDA or VecDuplicate() than the hassle of figuring out sizes.<br>
<br>
    > Or the vectors have to be formed by DMDACreateGlobalVector? I'm also not sure about what the dof and stencil width arguments do.<br>
    > <br>
    > I'm still unsure about the usage of MatCreateShell and MatShellSetOperation, since it seems that MyMatMult should still have 3 inputs just like MatMult (the matrix and two vectors). Since I'm not forming the matrix, does that mean the matrix input is meaningless but still needs to exist for the sake of this format?<br>
<br>
        Well the matrix input is your shell matrix so it likely has information you need to do your multiply routine. MatShellGetContext() (No you do not want to put your information about the matrix stencil inside global variables!)<br>
<br>
<br>
    > <br>
    > After I create such a shell matrix, can I use it like a regular matrix in KSP and utilize preconditioners?<br>
    > <br>
    > Thanks!<br>
    > Yuyun<br>
    > From: petsc-users <<a href="mailto:petsc-users-bounces@mcs.anl.gov" target="_blank">petsc-users-bounces@mcs.anl.gov</a>> on behalf of Yuyun Yang <<a href="mailto:yyang85@stanford.edu" target="_blank">yyang85@stanford.edu</a>><br>
    > Sent: Sunday, February 16, 2020 3:12 AM<br>
    > To: Smith, Barry F. <<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>><br>
    > Cc: <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] Matrix-free method in PETSc<br>
    >  <br>
    > Thank you, that is very helpful information indeed! I will try it and send you my code when it works.<br>
    > <br>
    > Best regards,<br>
    > Yuyun<br>
    > From: Smith, Barry F. <<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>><br>
    > Sent: Saturday, February 15, 2020 10:02 PM<br>
    > To: Yuyun Yang <<a href="mailto:yyang85@stanford.edu" target="_blank">yyang85@stanford.edu</a>><br>
    > Cc: <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] Matrix-free method in PETSc<br>
    >  <br>
    >   Yuyun,<br>
    > <br>
    >     If you are speaking about using a finite difference stencil on a structured grid where you provide the Jacobian vector products yourself by looping over the grid doing the stencil operation we unfortunately do not have exactly that kind of example. <br>
    > <br>
    >     But it is actually not difficult. I suggest starting with src/ts/examples/tests/ex22.c It computes the sparse matrix explicitly with FormIJacobian() <br>
    > <br>
    >     What you need to do is instead in main() use MatCreateShell() and MatShellSetOperation(,MATOP_MULT,(void (*)(void))MyMatMult) then provide the routine MyMatMult() to do your stencil based matrix free product; note that you can create this new routine by taking the structure of IFunction() and reorganizing it to do the Jacobian product instead. You will need to get the information about the shell matrix size on each process by calling DMDAGetCorners(). <br>
    > <br>
    >     You will then remove the explicit computation of the Jacobian, and also remove the Event stuff since you don't need it.<br>
    > <br>
    >      Extending to 2 and 3d is straight forward. <br>
    > <br>
    >      Any questions let us know.<br>
    > <br>
    >    Barry<br>
    > <br>
    >    If you like this would make a great merge request with your code to improve our examples.<br>
    > <br>
    > <br>
    > > On Feb 15, 2020, at 9:42 PM, Yuyun Yang <<a href="mailto:yyang85@stanford.edu" target="_blank">yyang85@stanford.edu</a>> wrote:<br>
    > > <br>
    > > Hello team,<br>
    > > <br>
    > > I wanted to apply the Krylov subspace method to a matrix-free implementation of a stencil, such that the iterative method acts on the operation without ever constructing the matrix explicitly (for example, when doing backward Euler).<br>
    > > <br>
    > > I'm not sure whether there is already an example for that somewhere. If so, could you point me to a relevant example?<br>
    > > <br>
    > > Thank you!<br>
    > > <br>
    > > Best regards,<br>
    > > Yuyun<br>
<br>
<br>
<br>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>