<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=us-ascii">
</head>
<body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class="">
<br class="">
<div><br class="">
<blockquote type="cite" class="">
<div class="">On May 10, 2018, at 4:12 PM, Jed Brown <<a href="mailto:jed@jedbrown.org" class="">jed@jedbrown.org</a>> wrote:</div>
<br class="Apple-interchange-newline">
<div class=""><span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">"Zhang,
 Hong" <</span><a href="mailto:hongzhang@anl.gov" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px;" class="">hongzhang@anl.gov</a><span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">>
 writes:</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
Dear PETSc folks,<br class="">
<br class="">
Current TS APIs (IFunction/IJacobian+RHSFunction/RHSJacobian) were designed for the fully implicit formulation F(t,U,Udot) = G(t,U).<br class="">
Shampine's paper (<a href="https://www.sciencedirect.com/science/article/pii/S0377042706004110?via%3Dihub%3Chttps://www.sciencedirect.com/science/article/pii/S0377042706004110?via=ihub%3E" class="">https://www.sciencedirect.com/science/article/pii/S0377042706004110?via%3Dihub<https://www.sciencedirect.com/science/article/pii/S0377042706004110?via=ihub></a>)
 explains some reasoning behind it.<br class="">
<br class="">
Our formulation is general enough to cover all the following common cases<br class="">
<br class="">
 *   Udot = G(t,U) (classic ODE)<br class="">
 *   M Udot = G(t,U)  (ODEs/DAEs for mechanical and electronic systems)<br class="">
 *   M(t,U) Udot = G(t,U) (PDEs)<br class="">
<br class="">
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 class="">
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 class="">
</blockquote>
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">What
 isn't straightforward about AD or FD on the IFunction?  That one</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">bit
 of chain rule?</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
</div>
</blockquote>
<div><br class="">
</div>
<div>Assembling shift*dF/dUdot+dF/dU could mean just a bit extra work for experienced users or experts, but could also be fairly annoying for others. Given IFunction, an AD tool generates dF/dUdot and dF/dU naturally, but it would not do that extra "bit of
 chain rule" automatically. It is our fault to force users to understand or worry about the shift business while it can be done easily by PETSc.   </div>
<br class="">
<blockquote type="cite" class="">
<div class=""><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
2. It is difficult to switch from fully implicit to fully explicit. Users cannot use explicit methods when they provide IFunction/IJacobian.<br class="">
</blockquote>
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">This
 is a real issue, but it's extremely common for PDE to have boundary</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">conditions
 enforced as algebraic constraints, thus yielding a DAE.</span></div>
</blockquote>
<div><br class="">
</div>
<div>Of course the conversion is not always possible. But it is possible for many cases if we have the mass matrix. That is one of the reasons why we want to have dF/dUdot.</div>
<br class="">
<blockquote type="cite" class="">
<div class=""><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
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 class="">
</blockquote>
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">I
 don't understand<span class="Apple-converted-space"> </span></span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
</div>
</blockquote>
<div><br class="">
</div>
<div>Whenever IJacobian is called, both dF/dUdot (the mass matrix) and dF/dU will be computed. Think about the case M Udot = G(t,U) where M is constant. For efficiency, we actually just need to update the Jacobian of G during time integration and reuse the
 mass matrix. This is just one of many cases that we can take advantage of the structure of the mass matrix, and hiding the mass matrix in IJacobian prevents us doing that.</div>
<br class="">
<blockquote type="cite" class="">
<div class=""><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
4. Reshifting the Jacobian is unnecessarily expensive and sometimes buggy.<br class="">
</blockquote>
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">Why
 is "reshifting" needed?  After a step is rejected and when the step</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">size
 changes for a linear constant operator?</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
</div>
</blockquote>
<div><br class="">
</div>
It doesn't not have to be a linear constant operator. Reusing Jacobian across stages or even time steps for nonlinear problems is not uncommon.</div>
<div><br class="">
<blockquote type="cite" class="">
<div class=""><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
Consider the scenario below.<br class="">
shift = a;<br class="">
 TSComputeIJacobian()<br class="">
 shift = b;<br class="">
 TSComputeIJacobian() // with the same U and Udot as last call<br class="">
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 class="">
<br class="">
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 class="">
<br class="">
 *   easy access to the unshifted Jacobian dF/dU (this simplifies the adjoint implementation a lot),<br class="">
</blockquote>
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">How
 does this simplify the adjoint?</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
</div>
</blockquote>
<div><br class="">
</div>
<div>For example, when computing the adjoint for M Udot = G(t,U), I often need to compute just dG/dU or multiplying M with the adjoint variables. It is doable with IJacobian, but tricky and inefficient since IJacobian always evaluates both M and dG/dU.</div>
<div><br class="">
</div>
<blockquote type="cite" class="">
<div class=""><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
 *   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 class="">
</blockquote>
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">But
 longer critical path, more storage required, and more data motion.</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">And
 if the mass matrix is simple, won't it take a very small fraction of</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">time,
 thus have little gains from "optimizing it"?</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
</div>
</blockquote>
<div><br class="">
</div>
<div>I would not worry much about the storage. In most cases, mass matrices would be diagonal, thus not eating much memory. For complicated mass matrices, most PDE solvers like libmesh already store them separately. </div>
<div><br class="">
</div>
<div>The gains may depends on not only complexity but also scale. Savings could be significant if we can avoid computing a huge constant mass matrix over and over again. And the optimization normally does not require too much effort. For example, reusing the
 mass matrix would be similar to lagging a Jacobian.</div>
<br class="">
<blockquote type="cite" class="">
<div class=""><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
 *   easy conversion from fully implicit to fully explicit,<br class="">
 *   an IJacobian without shift parameter that is easy to understand and easy to port to other software.<br class="">
</blockquote>
<br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">Note
 that even CVODE has an interface similar to PETSc; e.g., gamma</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<span style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; float: none; display: inline !important;" class="">parameter
 in CVSpilsPrecSetupFn.</span><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
</div>
</blockquote>
<div><br class="">
</div>
<div>The good thing in CVODE is that it does not have gamma parameter in the Jacobian function, see CVSpilsJacTimesSetupFn.</div>
<div><br class="">
</div>
<div>Thanks,</div>
<div>Hong</div>
<br class="">
<blockquote type="cite" class="">
<div class=""><br style="caret-color: rgb(0, 0, 0); font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
<blockquote type="cite" style="font-family: Verdana; font-size: 14px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; text-decoration: none;" class="">
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 class="">
<br class="">
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 class="">
<br class="">
Thanks,<br class="">
Hong</blockquote>
</div>
</blockquote>
</div>
<br class="">
</body>
</html>