<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, May 11, 2018 at 12:25 AM, Smith, Barry F. <span dir="ltr"><<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@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"><div class="HOEnZb"><div class="h5"><br>
<br>
> On May 10, 2018, at 4:12 PM, Jed Brown <<a href="mailto:jed@jedbrown.org">jed@jedbrown.org</a>> wrote:<br>
> <br>
> "Zhang, Hong" <<a href="mailto:hongzhang@anl.gov">hongzhang@anl.gov</a>> writes:<br>
> <br>
>> Dear PETSc folks,<br>
>> <br>
>> Current TS APIs (IFunction/IJacobian+<wbr>RHSFunction/RHSJacobian) were designed for the fully implicit formulation F(t,U,Udot) = G(t,U).<br>
>> Shampine's paper (<a href="https://www.sciencedirect.com/science/article/pii/S0377042706004110?via%3Dihub" rel="noreferrer" target="_blank">https://www.sciencedirect.<wbr>com/science/article/pii/<wbr>S0377042706004110?via%3Dihub</a><<a href="https://www.sciencedirect.com/science/article/pii/S0377042706004110?via=ihub" rel="noreferrer" target="_blank">h<wbr>ttps://www.sciencedirect.com/<wbr>science/article/pii/<wbr>S0377042706004110?via=ihub</a>>) explains some reasoning behind it.<br>
>> <br>
>> Our formulation is general enough to cover all the following common cases<br>
>> <br>
>>  *   Udot = G(t,U) (classic ODE)<br>
>>  *   M Udot = G(t,U)  (ODEs/DAEs for mechanical and electronic systems)<br>
>>  *   M(t,U) Udot = G(t,U) (PDEs)<br>
>> <br>
>> Yet the TS APIs provide the capability to solve both simple problems and complicated problems. However, we are not doing well to make TS easy to use and efficient especially for simple problems. Over the years, we have clearly seen the main drawbacks including:<br>
>> 1. The shift parameter exposed in IJacobian is terribly confusing, especially to new users. Also it is not conceptually straightforward when using AD or finite differences on IFunction to approximate IJacobian.<br>
> <br>
> What isn't straightforward about AD or FD on the IFunction?  That one<br>
> bit of chain rule?<br>
> <br>
>> 2. It is difficult to switch from fully implicit to fully explicit. Users cannot use explicit methods when they provide IFunction/IJacobian.<br>
> <br>
> This is a real issue, but it's extremely common for PDE to have boundary<br>
> conditions enforced as algebraic constraints, thus yielding a DAE.<br>
> <br>
>> 3. The structure of mass matrix is completely invisible to TS. This means giving up all the opportunities to improve efficiency. For example, when M is constant or weekly dependent on U, we might not want to evaluate/update it every time IJacobian is called. If M is diagonal, the Jacobian can be shifted more efficiently than just using MatAXPY().<br>
> <br>
> I don't understand <br>
> <br>
>> 4. Reshifting the Jacobian is unnecessarily expensive and sometimes buggy.<br>
> <br>
> Why is "reshifting" needed?  After a step is rejected and when the step<br>
> size changes for a linear constant operator?<br>
> <br>
>> Consider the scenario below.<br>
>> shift = a;<br>
>>  TSComputeIJacobian()<br>
>>  shift = b;<br>
>>  TSComputeIJacobian() // with the same U and Udot as last call<br>
>> Changing the shift parameter requires the Jacobian function to be evaluated again. If users provide only RHSJacobian, the Jacobian will not be updated/reshifted in the second call because TSComputeRHSJacobian() finds out that U has not been changed. This issue is fixable by adding more logic into the already convoluted implementation of TSComputeIJacobian(), but the intention here is to illustrate the cumbersomeness of current IJacobian and the growing complications in TSComputeIJacobian() that IJacobian causes.<br>
>> <br>
>> So I propose that we have two separate matrices dF/dUdot and dF/dU, and remove the shift parameter from IJacobian. dF/dU will be calculated by IJacobian; dF/dUdot will be calculated by a new callback function and default to an identity matrix if it is not provided by users. Then the users do not need to assemble the shifted Jacobian since TS will handle the shifting whenever needed. And knowing the structure of dF/dUdot (the mass matrix), TS will become more flexible. In particular, we will have<br>
>> <br>
>>  *   easy access to the unshifted Jacobian dF/dU (this simplifies the adjoint implementation a lot),<br>
> <br>
> How does this simplify the adjoint?<br>
> <br>
>>  *   plenty of opportunities to optimize TS when the mass matrix is diagonal or constant or weekly dependent on U (which accounts for almost all use cases in practice),<br>
> <br>
> But longer critical path,<br>
<br>
</div></div>   What do you mean by longer critical path?<br>
<span class=""><br>
> more storage required, and more data motion.<br>
<br>
</span>  The extra storage needed is related to the size of the mass matrix correct? And the extra data motion is related to the size of the mass matrix correct?<br>
<br>
   Is the extra work coming from a needed call to MatAXPY (to combine the scaled mass matrix with the Jacobian) in Hong's approach? While, in theory, the user can avoid the MatAXPY in the current code if they have custom code that assembles directly the scaled mass matrix and Jacobian? But surely most users would not write such custom code and would themselves keep a copy of the mass matrix (likely constant) and use MatAXPY() to combine the copy of the mass matrix with the Jacobian they compute at each timestep/stage?  Or am I missing something?</blockquote><div><br></div><div>I assemble the combined system directly.</div><div><br></div><div>   Matt</div><div><br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><span class="">
> And if the mass matrix is simple, won't it take a very small fraction of<br>
> time, thus have little gains from "optimizing it"?<br>
<br>
</span>   Is the Function approach only theoretically much more efficient than Hong's approach when the mass matrix is nontrivial? That is the mass matrix has a nonzero structure similar to the Jacobian? <br>
<div class="HOEnZb"><div class="h5">> <br>
>>  *   easy conversion from fully implicit to fully explicit,<br>
>>  *   an IJacobian without shift parameter that is easy to understand and easy to port to other software.<br>
> <br>
> Note that even CVODE has an interface similar to PETSc; e.g., gamma<br>
> parameter in CVSpilsPrecSetupFn.<br>
> <br>
>> Regarding the changes on the user side, most IJacobian users should not have problem splitting the old Jacobian if they compute dF/dUdot and dF/dU explicitly. If dF/dUdot is too complicated to build, matrix-free is an alternative option.<br>
>> <br>
>> While this proposal is somehow related to Barry's idea of having a TSSetMassMatrix() last year, I hope it provides more details for your information. Any of your comments and feedback would be appreciated.<br>
>> <br>
>> Thanks,<br>
>> Hong<br>
<br>
</div></div></blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature" data-smartmail="gmail_signature"><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.caam.rice.edu/~mk51/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div>
</div></div>