Code Monkey home page Code Monkey logo

dpfehm.jl's Introduction

DPFEHM: A Differentiable Subsurface Physics Simulator

codecov

Description

DPFEHM is a Julia module that includes differentiable numerical models with a focus on the Earth's subsurface, especially fluid flow. Currently it supports the groundwater flow equations (single phase flow), Richards equation (air/water), the advection-dispersion equation, and the 2d wave equation. Since it is differentiable, it can easily be combined with machine learning models in a workflow such as this:

Pressure management with DPFEHM and machine learning

This workflow shows how to train a machine learning model to mitigate problems with injecting fluid into the earth's subsurface (such as induced seismicity or leakage of carbon dioxide). More details on this workflow are available here.

Installation

Within Julia, you can install DPFEHM and test that it works by running

import Pkg
Pkg.add("DPFEHM")
Pkg.test("DPFEHM")

Basic Usage

The examples are a good place to get started to see how to use DPFEHM. Two examples will be described in detail here that illustrate the basic usage patterns via an examples of steady-state single-phase flow and transient Richards equation.

Here, we solve a steady-state single phase flow problem . Let's start by importing several libraries that we will use.

import DPFEHM
import GaussianRandomFields
import Optim
import PyPlot
import Random
import Zygote
Random.seed!(0)#set the seed so we can reproduce the same results with each run

Next, we'll set up the grid. Here, we use a regular grid with 100,000 nodes that covers a domain that is 50 meters by 50 meters by 5 meters.

mins = [0, 0, 0]; maxs = [50, 50, 5]#size of the domain, in meters
ns = [100, 100, 10]#number of nodes on the grid
coords, neighbors, areasoverlengths, _ = DPFEHM.regulargrid3d(mins, maxs, ns)#build the grid

The result of this grid-building is three variables that we will use. The, coords is a matrix describing the coordinates of the cell centers on the grid. The second, neighbors, is an array describing which cells neighbor other cells. The third, areasoverlengths, is another array whose length is equal to the length of neighbors and describes the area of the interface between two neighboring cells dividing by the length between the cell centers. The last variable is dumped to _ and gives the volumes of the cells. The volumes of the cells are not needed for steady state problems, so they are not used in this example.

Now we set up the boundary conditions.

Qs = zeros(size(coords, 2))
injectionnode = 1#inject in the lower left corner
Qs[injectionnode] = 1e-4#m^3/s
dirichletnodes = Int[size(coords, 2)]#fix the pressure in the upper right corner
dirichleths = zeros(size(coords, 2))
dirichleths[size(coords, 2)] = 0.0

The variable Qs describes fluid sources/sinks -- the amount of fluid injected at cell i on the grid is given by Qs[i]. In this example, the only place were we inject fluid is at node 1. Another variable, dirichletnodes is a list of cells at which the pressure will be fixed. In this example, the pressure is fixed at the last cell, which is cell number size(coords, 2). The variable dirichleths describes the pressures (or heads in hydrology jargon) that the cells are fixed at. Note that the length of dirichleths is size(coords, 2), but these values are ignored except at the indices that appear in dirichletnodes.

The final set-up step before moving on to solving the equations is to construct a heterogeneous conductivity field. Here, we use the package GaussianRandomFields to construct a conductivity field with a correlation length of 50 meters. The mean of the log-conductivity is 1e-4meters/second (note that we use a natural logarithm when defining this), and the standard deviation of the log-conductivity is 1. GaussianRandomFields is used to construct a field in 2 dimensions and then it is copied through each of the layers, so that the heterogeneity only exists in the x and y coordinate directions, but not in the z direction.

lambda = 50.0#meters -- correlation length of log-conductivity
sigma = 1.0#standard deviation of log-conductivity
mu = -9.0#mean of log conductivity -- ~1e-4 m/s, like clean sand here https://en.wikipedia.org/wiki/Hydraulic_conductivity#/media/File:Groundwater_Freeze_and_Cherry_1979_Table_2-2.png
cov = GaussianRandomFields.CovarianceFunction(2, GaussianRandomFields.Matern(lambda, 1; σ=sigma))
x_pts = range(mins[1], maxs[1]; length=ns[1])
y_pts = range(mins[2], maxs[2]; length=ns[2])
num_eigenvectors = 200
grf = GaussianRandomFields.GaussianRandomField(cov, GaussianRandomFields.KarhunenLoeve(num_eigenvectors), x_pts, y_pts)
logKs = zeros(reverse(ns)...)
logKs2d = mu .+ GaussianRandomFields.sample(grf)'#generate a random realization of the log-conductivity field
for i = 1:ns[3]#copy the 2d field to each of the 3d layers
	v = view(logKs, i, :, :)
	v .= logKs2d
end

The conductivity field is shown:

Conductivity field

Now, we look to solve the flow problem. First, we define a helper function, logKs2Ks_neighbors. This function is needed because the flow solver wants to know the conductivity on the interface between two cells, but our previous construction defined the conductivities at the cells themselves. It also converts from log-conductivity to conductivity and uses the geometric mean to move from the cells to the interfaces. The heart of this code is the call to DPFEHM.groundwater_steadystate, which solves the single phase steady-state flow problem that we pose. The solveforh function calls this function and returns the result after reshaping.

logKs2Ks_neighbors(Ks) = exp.(0.5 * (Ks[map(p->p[1], neighbors)] .+ Ks[map(p->p[2], neighbors)]))#convert from permeabilities at the nodes to permeabilities connecting the nodes
function solveforh(logKs, dirichleths)
	@assert length(logKs) == length(Qs)
	Ks_neighbors = logKs2Ks_neighbors(logKs)
	return reshape(DPFEHM.groundwater_steadystate(Ks_neighbors, neighbors, areasoverlengths, dirichletnodes, dirichleths, Qs), reverse(ns)...)
	end
end

With this function in hand, we can solve the problem using the solveforh wrapper function we previously defined. This function requires us to explicitly pass in logKs (the hydraulic conductivity) and dirichleths (the fixed-head boundary condition), but the other inputs to DPFEHM.groundwater_steadystate are fixed to global values.

h = solveforh(logKs, dirichleths)#solve for the head

The head at the bottom layer of the domain is shown (note the pressure is higher in the lower corner, where there is fluid injection, than in the rest of the domain):

Head field

DPFEHM also allows us to compute the gradient of functions involving DPFEHM.groundwater_steadystate using Zygote.gradient or Zygote.pullback.

isfreenode, nodei2freenodei, freenodei2nodei = DPFEHM.getfreenodes(length(dirichleths), dirichletnodes)
gradient_node = nodei2freenodei[div(size(coords, 2), 2) + 500]
gradient_node_x = coords[1, gradient_node]
gradient_node_y = coords[2, gradient_node]
grad = Zygote.gradient((x, y)->solveforh(x, y)[gradient_node], logKs, dirichleths)#calculate the gradient (which involves a redundant calculation of the forward pass)
function_evaluation, back = Zygote.pullback((x, y)->solveforh(x, y)[gradient_node], logKs, dirichleths)#this pullback thing lets us not redo the forward pass
print("gradient time")
grad2 = back(1.0)#compute the gradient of a function involving solveforh using Zygote.pullback

Note that the function DPFEHM.getfreenodes allows one to map indices between the free nodes (i.e., the ones that do not have fixed-pressure boundary conditions) and all nodes. The gradient of logK at the bottom layer of the domain is shown:

Gradient field

Now, we consider an example using DPFEHM's solver for Richards equation, which can be used to model flow in a porous medium where both air and water fill the pores (i.e., unsaturated flow). This example is similar to the previous example and we again start by importing several libraries, setting up the grid (lower resolution this time), the boundary conditions, and the conductivity.

import DifferentiableBackwardEuler
import DPFEHM
import GaussianRandomFields
import PyPlot
import Random
import Zygote
Random.seed!(0)#set the seed so we get the same permeability over and over

#set up the grid
mins = [0, 0, 0]; maxs = [50, 50, 5]#size of the domain, in meters
ns = [30, 30, 3]#number of nodes on the grid
coords, neighbors, areasoverlengths, volumes = DPFEHM.regulargrid3d(mins, maxs, ns)#build the grid

#set up the boundary conditions
Qs = zeros(size(coords, 2))
injectionnode = 1#inject in the lower left corner
Qs[injectionnode] = 1e-4#m^3/s
dirichletnodes = Int[size(coords, 2)]#fix the pressure in the upper right corner
dirichleths = zeros(size(coords, 2))
dirichleths[size(coords, 2)] = 0.0

#set up the conductivity field
lambda = 50.0#meters -- correlation length of log-conductivity
sigma = 1.0#standard deviation of log-conductivity
mu = -9.0#mean of log conductivity -- ~1e-4 m/s, like clean sand here https://en.wikipedia.org/wiki/Hydraulic_conductivity#/media/File:Groundwater_Freeze_and_Cherry_1979_Table_2-2.png
cov = GaussianRandomFields.CovarianceFunction(2, GaussianRandomFields.Matern(lambda, 1; σ=sigma))
x_pts = range(mins[1], maxs[1]; length=ns[1])
y_pts = range(mins[2], maxs[2]; length=ns[2])
num_eigenvectors = 200
grf = GaussianRandomFields.GaussianRandomField(cov, GaussianRandomFields.KarhunenLoeve(num_eigenvectors), x_pts, y_pts)
logKs = zeros(reverse(ns)...)
logKs2d = mu .+ GaussianRandomFields.sample(grf)'#generate a random realization of the log-conductivity field
for i = 1:ns[3]#copy the 2d field to each of the 3d layers
	v = view(logKs, i, :, :)
	v .= logKs2d
end

Since we'll be solving a time-dependent problem this time, we must set the initial condition and define a storage parameter. Note that in unsaturated flows the storage parameter is often neglected (and this can be achieved by setting specificstorage to an array of ones), but this problem involves a saturated flow so we include it here. Since this is a multi-phase flow problem, we also need to define a couple parameters that control the hydraulic conductivity's dependence on the saturation. The conductivity is given by the conductivity when saturated multiplied by a relative permeability (which depends on the saturation and varies between 0 and 1).

#set up the initial condition, the storage, and the van genuchten parameters for relative permeability
h0 = zeros(size(coords, 2))#initial condition
specificstorage = fill(0.1, size(coords, 2))#storage
alphas = fill(0.5, length(neighbors))#van genuchten relative permeability parameters
Ns = fill(1.25, length(neighbors))

With the basic description of the problem complete, now we can start to write the code for solving the equations. Note that the solveforh function does not call DPFEHM.richards_steadystate to solve the equations, and instead calls DifferentiableBackwardEuler.steps. The first argument to DifferentiableBackwardEuler.steps is the initial condition, but only at the nodes that are not controlled by the Dirichlet boundary conditions. The most important parts of this call are f_richards, f_richards_u, and f_richards_p, which we will describe in the next paragraph. The argument 0.0 is the initial time, and 60 * 60 * 24 * 365 * 1 gives the simulation time (in seconds, so 1 year). The keyword arguments will eventually be passed to DifferentialEquations.solve. The last step adds the boundary conditions back into the solution, which is needed since DifferentiableBackwardEuler.steps only solves the equations on the free nodes (i.e, the nodes where the pressure is not fixed by the boundary conditions).

logKs2Ks_neighbors(Ks) = exp.(0.5 * (Ks[map(p->p[1], neighbors)] .+ Ks[map(p->p[2], neighbors)]))#convert from permeabilities at the nodes to permeabilities at the interface between nodes using the geometric mean
isfreenode, nodei2freenodei, freenodei2nodei = DPFEHM.getfreenodes(length(Qs), dirichletnodes)
function solveforh(logKs, dirichleths)
	@assert length(logKs) == length(Qs)
	Ks_neighbors = logKs2Ks_neighbors(logKs)
	p = [Ks_neighbors; dirichleths]
	h_richards = DifferentiableBackwardEuler.steps(h0[isfreenode], f_richards, f_richards_u, f_richards_p, f_richards_t, p, 0.0, 60 * 60 * 24 * 365 * 1; abstol=1e-1, reltol=1e-1)
	h_with_bcs = hcat(map(i->DPFEHM.addboundaryconditions(h_richards[:, i], dirichletnodes, dirichleths, isfreenode, nodei2freenodei), 1:size(h_richards, 2))...)#add the dirichlet boundary conditions back
	return h_with_bcs
end

Now, we define the key functions f_richards, f_richards_u, and f_richards_p from the previous paragraph. The function f_richards basically tells DifferentiableBackwardEuler to solve du/dt=f_richards and this is the function richards_residuals. The function f_richards_u is the Jacobian of f_richards with respect to the variable that is being integrated by DifferentiableBackwardEuler. We can compute the Jacobian of richards_residuals with respect to any of its inputs using the function DPFEHM.richards_XYZ where XZY is the name of the argument (as defined within DPFEHM). In the jargon of DPFEHM's Richards equation solve, the variable we are solving for is named psi, so f_richards_u just unpacks the parameters and calls DPFEHM.richards_psi. The function f_richards_p is the Jacobian of f_richards with respect to the parameter p. Since p consists of Ks and dirichleths (or dirichletpsis in the jargon of DPFEHM's Richards solver), we concatenate the two Jacobians DPFEHM.richards_Ks and DPFEHM.richards_dirichletpsis. The last function f_richards_t is currently unused, but in principle should give the Jacobian of f_richards with respect to t.

#set up some functions needed by DifferentiableBackwardEuler
function f_richards(u, p, t)#tells DifferentiableBackwardEuler to solve du/dt=f_richards
	Ks, dirichleths = unpack(p)
	return DPFEHM.richards_residuals(u, Ks, neighbors, areasoverlengths, dirichletnodes, dirichleths, coords, alphas, Ns, Qs, specificstorage, volumes)
end
function f_richards_u(u, p, t)#give DifferentiableBackwardEuler the derivative of f_richards with respect to u -- needed for the backward euler method that we use
	Ks, dirichleths = unpack(p)
	return DPFEHM.richards_psi(u, Ks, neighbors, areasoverlengths, dirichletnodes, dirichleths, coords, alphas, Ns, Qs, specificstorage, volumes)
end
function f_richards_p(u, p, t)#give DifferentiableBackwardEuler the derivative of f_richards with respect to p -- needed for computing gradients with respect to p of functions involving the richards equation solution
	Ks, dirichleths = unpack(p)
	J1 = DPFEHM.richards_Ks(u, Ks, neighbors, areasoverlengths, dirichletnodes, dirichleths, coords, alphas, Ns, Qs, specificstorage, volumes)
	J2 = DPFEHM.richards_dirichletpsis(u, Ks, neighbors, areasoverlengths, dirichletnodes, dirichleths, coords, alphas, Ns, Qs, specificstorage, volumes)
	return hcat(J1, J2)
end
f_richards_t(u, p, t) = zeros(length(u))#the DifferentiableBackwardEuler API requires this but it currently isn't used

These functions use a helper function, unpack, which unpacks the "parameters" p from a big array into two smaller arrays. Here, unpack converts p into the conductivities, Ks and boundary conditions, dirichleths. We can also think of packing the parameters by doing `p = [Ks; dirichleths]'. This packing/unpacking is needed because DifferentiableBackwardEuler needs the parameters to be in a single array.

function unpack(p)
	@assert length(p) == length(neighbors) + size(coords, 2)
	Ks = p[1:length(neighbors)]
	dirichleths = p[length(neighbors) + 1:length(neighbors) + size(coords, 2)]
	return Ks, dirichleths
end

Now, we can solve the equations using solveforh.

h = solveforh(logKs, dirichleths)#solve for the head

Pressure

Finally, we can compute the gradient of functions involving the solution of these equations using Zygote.gradient or Zygote.pullback.

hflat2h3d(h) = reshape(h, reverse(ns)...)
gradient_node = div(size(coords, 2) + ns[3] * ns[2], 2)
gradient_node_x = coords[1, gradient_node]
gradient_node_y = coords[2, gradient_node]
grad = Zygote.gradient((x, y)->hflat2h3d(solveforh(x, y)[:, end])[gradient_node], logKs, dirichleths)#calculate the gradient (which involves a redundant calculation of the forward pass)

Gradient

Advanced usage

The examples illustrate more advanced usage including inverse problems, combining DPFEHM with a neural network, flow on discrete fracture networks, as well as solving the advection-dispersion and wave equations.

Note on the numerical methods

DPFEHM generally uses a two-point flux approximation to discretize the equations. This means, for example, that when the fluid flux between cell $i$ and cell $j$ that is proportional to the pressure gradient, the pressure gradient is approximated as being proportional to $p_i-p_j$ where $p_i$ and $p_j$ are the pressures in cells $i$ and $j$, respectively.

License

DPFEHM is provided under a BSD style license. See LICENSE.md file for the full text.

This package is part of the Orchard suite, known internally as C20086 Orchard.

Development and questions

If you would like to contribute to DPFEHM, please for the repo and make a pull request. If you have any questions, please contact Daniel O'Malley, [email protected].

dpfehm.jl's People

Contributors

apachalieva avatar kevin-mattheus-moerman avatar montyvesselinov avatar omalled avatar sygreer avatar

Stargazers

 avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar  avatar

Watchers

 avatar  avatar  avatar  avatar  avatar  avatar

dpfehm.jl's Issues

document examples

Apart from 2 examples, non of them seem to be documented. Just a few lines of in-code documentation to explain what the example does would be really useful to readers.

Suggested changes to JOSS paper

DPFEHM authors: Thank you for both making this code publicly available and for submitting a description of it to JOSS. This is very nice work! I sincerely apologize for being months behind on reviewing it.

It appears that many concerns that I might have had with the original submission have already been addressed by the other referees and the DPFEHM authors. The only real concerns I have now have to do with the JOSS paper itself (the PDF "paper", I mean). I think that it is find, overall, but it makes a few statements that I think need to be softened or qualified, or that need some additional context or information:

Item 1: The paper begins by naming some subsurface flow and transport codes and then makes the statement that "...these models cannot be efficiently integrated with machine learning frameworks. They cannot because they do not support automatic differentiation...". While it is certainly easiest and ideal to have automatic differentiation for such integration, I disagree that it is not possible to use such models in a machine learning framework; I point to a recent paper by Satish Karra et al. (https://arxiv.org/pdf/2109.03956.pdf) where PFLOTRAN is integrated -- with at least enough efficiency to be practicable -- into a machine learning workflow. The efficiency could be further improved by using a discrete adjoint approach. Although this is not done in the Karra et al. paper, I am fairly confident that such an approach could be implemented without any heroic coding effort in PFLOTRAN by using the TSAdjoint component of PETSc (https://arxiv.org/abs/1912.07696), since PFLOTRAN already has some support for the PETSc TS timesteppers and because Jacobians for some relevant adjoint systems are already coded. Given the need to code up the adjoint Jacobians, though, a nicely working automatic differentiation system (like you have with DPFEHM) seems like the ideal solution. But I suggest revising the language a bit to state that there are difficulties with using machine learning workflows with the aforementioned codes, rather than making the statement that it is impossible.

Item 2: The paper mentions other subsurface flow/transport models, but does not make it clear (for those who are less familiar with subsurface science models) that DPFEHM serves a different purpose and has different target usage than a code like PFLOTRAN or FEHM (or CrunchFlow, or TOUGH2, or several other codes). These codes support a very diverse set of features and process models that are not found in DPFEHM, which aims to support a much leaner set of features, but also opens up several interesting use cases by virtue of its differentiability and its full integration with the Julia ecosystem. It would be helpful to have just a sentence or two to make this clear.

Item 3: The paper says "DPFEHM uses a two-point flux approximation finite volume scheme. This means that an orthogonal grid is required to ensure convergence, similar to other codes such as FEHM and PFLOTRAN Alternative codes such as Amanzi-ATS (Mercer-Smith, 2020), which use a more advanced mimetic finite difference discretization, do not require orthogonal meshes". This is not an incorrect statement, but is somewhat misleading insofar as it only mentions MFD discretizations as ways to deal with non-orthogonal (I would suggest saying "non-K-orthogonal", actually) grids. There are, however, many approaches that have been developed to deal with such meshes, including multi-point flux approximations, least-squares gradient reconstruction, many different mixed finite element approaches, etc. It should be made clear that MFD approaches are not the only way to deal with such grids. I also suggest noting that virtually all commercial reservoir simulator codes use the two-point flux approximation, and that there are some good reasons that they do so, choosing simplicity and low computational expense over more robust discretizations such as MFD that involve higher complexity of implementation and greater computational expense. Indeed, PFLOTRAN used to support MFD, but this support was dropped since the cost of maintaining the MFD code was hard to justify when users found that other approaches (generally just working to get reasonably close to K-orthogonal grids) were satisfactory for accomplishing their science or engineering goals.

Item 4: This really is linked to item 1. In the second to last sentence, the paper states that "For non-differentiable models, the cost is equal to performing a number of simulations that is proportional to the number of parameters". This is not true for codes that can use the discrete adjoint solvers. (Though it is true that implementing support for this requires some developer work.)

I hope that the authors of the paper do not interpret my verbose comments as being overly critical of this work. I think this work is very nice, actually; I just want to ensure that the accompanying paper does not over-generalize or mislead readers.

I note one minor typographical issue:

  • There is a period missing between these sentences: "This means that an orthogonal grid is required to ensure convergence, similar to other codes such as FEHM and PFLOTRAN Alternative codes such as Amanzi-ATS (Mercer-Smith, 2020), which use a more advanced mimetic finite difference discretization, do not require orthogonal meshes"

Refs openjournals/joss-reviews#4560

ensure examples all run

It would be nice to have a script/action/whatever to run all the examples. Preferably, this could be run on every PullRequest, or at least, manually run to check for regressions. Currently, i'm going through all the examples, and many of them don't run without warnings/errors

fracture_network_solver_scaling/... often complains about "Package FEHM does not have DelimitedFiles in its dependencies..."
gw_multigrid_inverse/inputdeck.jl:30 file not found
piml/theislike_piml_setup.jl:73 params not defined

I haven't checked others

Refs openjournals/joss-reviews#4560

Documentation

In the current documentation (the README file), it is stated that "DPFEHM supports the groundwater flow equations (single phase flow), Richards equation (air/water), the advection-dispersion equation, and the 2d wave equation".

  • I think it would be helpful to add those equations to the documentation. Only Richards equation is linked to the wikipedia which doesn't match the description. For example, in README you talk about relative permeability and storage, but in that link they are not defined.
  • It is mentioned in the paper "DPFEHM uses a two-point flux approximation finite volume scheme". It would be also helpful to add this to the documentation and perhaps its discretized formulation.
  • It would be good to add purpose, input, output and variables that are needed inside each function in the DPFEHM.jl/src files. And adding some inline comments would help the user to contribute in the future.

Refs openjournals/joss-reviews#4560

Cannot run examples/groundwater_time_dependent.jl

Hi,

This is likely due to my inexperience, but i get the following size error:

> julia groundwater_time_dependent.jl 
PyPlot.Figure(PyObject <Figure size 640x480 with 2 Axes>)
forward solve timeERROR: LoadError: DimensionMismatch("arrays could not be broadcast to a common size; got a dimension with lengths 2699 and 2700")
Stacktrace:
  [1] _bcs1
    @ ./broadcast.jl:516 [inlined]
  [2] _bcs
    @ ./broadcast.jl:510 [inlined]
  [3] broadcast_shape
    @ ./broadcast.jl:504 [inlined]
  [4] combine_axes
    @ ./broadcast.jl:499 [inlined]
  [5] _axes
    @ ./broadcast.jl:224 [inlined]
  [6] axes
    @ ./broadcast.jl:222 [inlined]
  [7] combine_axes
    @ ./broadcast.jl:499 [inlined]
  [8] instantiate
    @ ./broadcast.jl:281 [inlined]
  [9] materialize(bc::Base.Broadcast.Broadcasted{Base.Broadcast.DefaultArrayStyle{1}, Nothing, typeof(*), Tuple{Base.Broadcast.Broadcasted{Base.Broadcast.DefaultArrayStyle{1}, Nothing, typeof(/), Tuple{Vector{Float64}, Vector{Float64}}}, Float64}})
    @ Base.Broadcast ./broadcast.jl:860
 [10] ode_determine_initdt(u0::Vector{Float64}, t::Float64, tdir::Float64, dtmax::Float64, abstol::Float64, reltol::Float64, internalnorm::typeof(DiffEqBase.ODE_DEFAULT_NORM), prob::SciMLBase.ODEProblem{Vector{Float64}, Tuple{Float64, Float64}, false, Vector{Float64}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, SciMLBase.StandardODEProblem}, integrator::OrdinaryDiffEq.ODEIntegrator{OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}, false, Vector{Float64}, Nothing, Float64, Vector{Float64}, Float64, Float64, Float64, Float64, Vector{Vector{Float64}}, SciMLBase.ODESolution{Float64, 2, Vector{Vector{Float64}}, Nothing, Nothing, Vector{Float64}, Vector{Vector{Vector{Float64}}}, SciMLBase.ODEProblem{Vector{Float64}, Tuple{Float64, Float64}, false, Vector{Float64}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, SciMLBase.StandardODEProblem}, OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}, OrdinaryDiffEq.InterpolationData{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Vector{Vector{Float64}}, Vector{Float64}, Vector{Vector{Vector{Float64}}}, OrdinaryDiffEq.ImplicitEulerConstantCache{OrdinaryDiffEq.NLSolver{DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, false, Vector{Float64}, Int64, Nothing, Float64, OrdinaryDiffEq.NLNewtonConstantCache{Float64, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, OrdinaryDiffEq.WOperator{false, Float64, LinearAlgebra.UniformScaling{Bool}, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, Vector{Float64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing}, SciMLBase.UDerivativeWrapper{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Float64, Vector{Float64}}}}}}, DiffEqBase.DEStats}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, OrdinaryDiffEq.ImplicitEulerConstantCache{OrdinaryDiffEq.NLSolver{DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, false, Vector{Float64}, Int64, Nothing, Float64, OrdinaryDiffEq.NLNewtonConstantCache{Float64, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, OrdinaryDiffEq.WOperator{false, Float64, LinearAlgebra.UniformScaling{Bool}, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, Vector{Float64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing}, SciMLBase.UDerivativeWrapper{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Float64, Vector{Float64}}}}}, OrdinaryDiffEq.DEOptions{Float64, Float64, Float64, Float64, OrdinaryDiffEq.PIController{Rational{Int64}}, typeof(DiffEqBase.ODE_DEFAULT_NORM), typeof(LinearAlgebra.opnorm), Nothing, SciMLBase.CallbackSet{Tuple{}, Tuple{}}, typeof(DiffEqBase.ODE_DEFAULT_ISOUTOFDOMAIN), typeof(DiffEqBase.ODE_DEFAULT_PROG_MESSAGE), typeof(DiffEqBase.ODE_DEFAULT_UNSTABLE_CHECK), DataStructures.BinaryHeap{Float64, DataStructures.FasterForward}, DataStructures.BinaryHeap{Float64, DataStructures.FasterForward}, Nothing, Nothing, Int64, Tuple{}, Tuple{}, Tuple{}}, Vector{Float64}, Float64, Nothing, OrdinaryDiffEq.DefaultInit})
    @ OrdinaryDiffEq ~/.julia/packages/OrdinaryDiffEq/UG9Mz/src/initdt.jl:237
 [11] auto_dt_reset!(integrator::OrdinaryDiffEq.ODEIntegrator{OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}, false, Vector{Float64}, Nothing, Float64, Vector{Float64}, Float64, Float64, Float64, Float64, Vector{Vector{Float64}}, SciMLBase.ODESolution{Float64, 2, Vector{Vector{Float64}}, Nothing, Nothing, Vector{Float64}, Vector{Vector{Vector{Float64}}}, SciMLBase.ODEProblem{Vector{Float64}, Tuple{Float64, Float64}, false, Vector{Float64}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, SciMLBase.StandardODEProblem}, OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}, OrdinaryDiffEq.InterpolationData{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Vector{Vector{Float64}}, Vector{Float64}, Vector{Vector{Vector{Float64}}}, OrdinaryDiffEq.ImplicitEulerConstantCache{OrdinaryDiffEq.NLSolver{DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, false, Vector{Float64}, Int64, Nothing, Float64, OrdinaryDiffEq.NLNewtonConstantCache{Float64, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, OrdinaryDiffEq.WOperator{false, Float64, LinearAlgebra.UniformScaling{Bool}, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, Vector{Float64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing}, SciMLBase.UDerivativeWrapper{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Float64, Vector{Float64}}}}}}, DiffEqBase.DEStats}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, OrdinaryDiffEq.ImplicitEulerConstantCache{OrdinaryDiffEq.NLSolver{DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, false, Vector{Float64}, Int64, Nothing, Float64, OrdinaryDiffEq.NLNewtonConstantCache{Float64, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, OrdinaryDiffEq.WOperator{false, Float64, LinearAlgebra.UniformScaling{Bool}, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, Vector{Float64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing}, SciMLBase.UDerivativeWrapper{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Float64, Vector{Float64}}}}}, OrdinaryDiffEq.DEOptions{Float64, Float64, Float64, Float64, OrdinaryDiffEq.PIController{Rational{Int64}}, typeof(DiffEqBase.ODE_DEFAULT_NORM), typeof(LinearAlgebra.opnorm), Nothing, SciMLBase.CallbackSet{Tuple{}, Tuple{}}, typeof(DiffEqBase.ODE_DEFAULT_ISOUTOFDOMAIN), typeof(DiffEqBase.ODE_DEFAULT_PROG_MESSAGE), typeof(DiffEqBase.ODE_DEFAULT_UNSTABLE_CHECK), DataStructures.BinaryHeap{Float64, DataStructures.FasterForward}, DataStructures.BinaryHeap{Float64, DataStructures.FasterForward}, Nothing, Nothing, Int64, Tuple{}, Tuple{}, Tuple{}}, Vector{Float64}, Float64, Nothing, OrdinaryDiffEq.DefaultInit})
    @ OrdinaryDiffEq ~/.julia/packages/OrdinaryDiffEq/UG9Mz/src/integrators/integrator_interface.jl:356
 [12] handle_dt!(integrator::OrdinaryDiffEq.ODEIntegrator{OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}, false, Vector{Float64}, Nothing, Float64, Vector{Float64}, Float64, Float64, Float64, Float64, Vector{Vector{Float64}}, SciMLBase.ODESolution{Float64, 2, Vector{Vector{Float64}}, Nothing, Nothing, Vector{Float64}, Vector{Vector{Vector{Float64}}}, SciMLBase.ODEProblem{Vector{Float64}, Tuple{Float64, Float64}, false, Vector{Float64}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, SciMLBase.StandardODEProblem}, OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}, OrdinaryDiffEq.InterpolationData{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Vector{Vector{Float64}}, Vector{Float64}, Vector{Vector{Vector{Float64}}}, OrdinaryDiffEq.ImplicitEulerConstantCache{OrdinaryDiffEq.NLSolver{DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, false, Vector{Float64}, Int64, Nothing, Float64, OrdinaryDiffEq.NLNewtonConstantCache{Float64, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, OrdinaryDiffEq.WOperator{false, Float64, LinearAlgebra.UniformScaling{Bool}, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, Vector{Float64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing}, SciMLBase.UDerivativeWrapper{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Float64, Vector{Float64}}}}}}, DiffEqBase.DEStats}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, OrdinaryDiffEq.ImplicitEulerConstantCache{OrdinaryDiffEq.NLSolver{DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, false, Vector{Float64}, Int64, Nothing, Float64, OrdinaryDiffEq.NLNewtonConstantCache{Float64, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, OrdinaryDiffEq.WOperator{false, Float64, LinearAlgebra.UniformScaling{Bool}, Float64, SciMLBase.DiffEqArrayOperator{Float64, SparseArrays.SparseMatrixCSC{Float64, Int64}, typeof(SciMLBase.DEFAULT_UPDATE_FUNC)}, Vector{Float64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing}, SciMLBase.UDerivativeWrapper{SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Float64, Vector{Float64}}}}}, OrdinaryDiffEq.DEOptions{Float64, Float64, Float64, Float64, OrdinaryDiffEq.PIController{Rational{Int64}}, typeof(DiffEqBase.ODE_DEFAULT_NORM), typeof(LinearAlgebra.opnorm), Nothing, SciMLBase.CallbackSet{Tuple{}, Tuple{}}, typeof(DiffEqBase.ODE_DEFAULT_ISOUTOFDOMAIN), typeof(DiffEqBase.ODE_DEFAULT_PROG_MESSAGE), typeof(DiffEqBase.ODE_DEFAULT_UNSTABLE_CHECK), DataStructures.BinaryHeap{Float64, DataStructures.FasterForward}, DataStructures.BinaryHeap{Float64, DataStructures.FasterForward}, Nothing, Nothing, Int64, Tuple{}, Tuple{}, Tuple{}}, Vector{Float64}, Float64, Nothing, OrdinaryDiffEq.DefaultInit})
    @ OrdinaryDiffEq ~/.julia/packages/OrdinaryDiffEq/UG9Mz/src/solve.jl:503
 [13] __init(prob::SciMLBase.ODEProblem{Vector{Float64}, Tuple{Float64, Float64}, false, Vector{Float64}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, SciMLBase.StandardODEProblem}, alg::OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}, timeseries_init::Tuple{}, ts_init::Tuple{}, ks_init::Tuple{}, recompile::Type{Val{true}}; saveat::Tuple{}, tstops::Tuple{}, d_discontinuities::Tuple{}, save_idxs::Nothing, save_everystep::Bool, save_on::Bool, save_start::Bool, save_end::Nothing, callback::Nothing, dense::Bool, calck::Bool, dt::Float64, dtmin::Nothing, dtmax::Float64, force_dtmin::Bool, adaptive::Bool, gamma::Rational{Int64}, abstol::Float64, reltol::Float64, qmin::Rational{Int64}, qmax::Int64, qsteady_min::Int64, qsteady_max::Rational{Int64}, beta1::Nothing, beta2::Nothing, qoldinit::Rational{Int64}, controller::Nothing, fullnormalize::Bool, failfactor::Int64, maxiters::Int64, internalnorm::typeof(DiffEqBase.ODE_DEFAULT_NORM), internalopnorm::typeof(LinearAlgebra.opnorm), isoutofdomain::typeof(DiffEqBase.ODE_DEFAULT_ISOUTOFDOMAIN), unstable_check::typeof(DiffEqBase.ODE_DEFAULT_UNSTABLE_CHECK), verbose::Bool, timeseries_errors::Bool, dense_errors::Bool, advance_to_tstop::Bool, stop_at_next_tstop::Bool, initialize_save::Bool, progress::Bool, progress_steps::Int64, progress_name::String, progress_message::typeof(DiffEqBase.ODE_DEFAULT_PROG_MESSAGE), userdata::Nothing, allow_extrapolation::Bool, initialize_integrator::Bool, alias_u0::Bool, alias_du0::Bool, initializealg::OrdinaryDiffEq.DefaultInit, kwargs::Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}})
    @ OrdinaryDiffEq ~/.julia/packages/OrdinaryDiffEq/UG9Mz/src/solve.jl:466
 [14] __solve(::SciMLBase.ODEProblem{Vector{Float64}, Tuple{Float64, Float64}, false, Vector{Float64}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, SciMLBase.StandardODEProblem}, ::OrdinaryDiffEq.ImplicitEuler{12, true, LinearSolve.UMFPACKFactorization, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}; kwargs::Base.Pairs{Symbol, Float64, Tuple{Symbol, Symbol}, NamedTuple{(:abstol, :reltol), Tuple{Float64, Float64}}})
    @ OrdinaryDiffEq ~/.julia/packages/OrdinaryDiffEq/UG9Mz/src/solve.jl:4
 [15] #solve_call#28
    @ ~/.julia/packages/DiffEqBase/aAyno/src/solve.jl:433 [inlined]
 [16] solve_up(prob::SciMLBase.ODEProblem{Vector{Float64}, Tuple{Float64, Float64}, false, Vector{Float64}, SciMLBase.ODEFunction{false, typeof(f_gw), LinearAlgebra.UniformScaling{Bool}, Nothing, Nothing, typeof(f_gw_u), Nothing, Nothing, SparseArrays.SparseMatrixCSC{Float64, Int64}, SparseArrays.SparseMatrixCSC{Float64, Int64}, Nothing, Nothing, Nothing, Nothing, Nothing, typeof(SciMLBase.DEFAULT_OBSERVED), Nothing}, Base.Pairs{Symbol, Union{}, Tuple{}, NamedTuple{(), Tuple{}}}, SciMLBase.StandardODEProblem}, sensealg::Nothing, u0::Vector{Float64}, p::Vector{Float64}, args::OrdinaryDiffEq.ImplicitEuler{0, true, Nothing, DiffEqBase.NLNewton{Rational{Int64}, Rational{Int64}, Rational{Int64}}, typeof(OrdinaryDiffEq.DEFAULT_PRECS), Val{:forward}, true, nothing}; kwargs::Base.Pairs{Symbol, Float64, Tuple{Symbol, Symbol}, NamedTuple{(:abstol, :reltol), Tuple{Float64, Float64}}})
    @ DiffEqBase ~/.julia/packages/DiffEqBase/aAyno/src/solve.jl:772
 [17] #solve#33
    @ ~/.julia/packages/DiffEqBase/aAyno/src/solve.jl:756 [inlined]
 [18] dae_steps_diffeq(y0::Vector{Float64}, f::Function, f_y::typeof(f_gw_u), f_p::Function, f_t::Function, p::Vector{Float64}, t0::Float64, tfinal::Float64, M::LinearAlgebra.UniformScaling{Bool}; kwargs::Base.Pairs{Symbol, Float64, Tuple{Symbol, Symbol}, NamedTuple{(:abstol, :reltol), Tuple{Float64, Float64}}})
    @ DifferentiableBackwardEuler ~/.julia/packages/DifferentiableBackwardEuler/PB0Mz/src/DifferentiableBackwardEuler.jl:108
 [19] dae_steps(y0::Vector{Float64}, f::Function, f_y::Function, f_p::Function, f_t::Function, p::Vector{Float64}, t0::Float64, tfinal::Float64, M::LinearAlgebra.UniformScaling{Bool}; kwargs::Base.Pairs{Symbol, Float64, Tuple{Symbol, Symbol}, NamedTuple{(:abstol, :reltol), Tuple{Float64, Float64}}})
    @ DifferentiableBackwardEuler ~/.julia/packages/DifferentiableBackwardEuler/PB0Mz/src/DifferentiableBackwardEuler.jl:115
 [20] #steps#23
    @ ~/.julia/packages/DifferentiableBackwardEuler/PB0Mz/src/DifferentiableBackwardEuler.jl:112 [inlined]
 [21] solveforh(logKs::Array{Float64, 3}, dirichleths::Vector{Float64})
    @ Main .../DPFEHM.jl/examples/groundwater_time_dependent.jl:61
 [22] top-level scope
    @ ./timing.jl:220
in expression starting at .../DPFEHM.jl/examples/groundwater_time_dependent.jl:92

(the "..." is a path you don't need to know)

TagBot trigger issue

This issue is used to trigger TagBot; feel free to unsubscribe.

If you haven't already, you should update your TagBot.yml to include issue comment triggers.
Please see this post on Discourse for instructions and more details.

If you'd like for me to do this for you, comment TagBot fix on this issue.
I'll open a PR within a few hours, please be patient!

code coverage

It'd be great if the readme.md had a link to a page that showed which lines of code were covered by the tests.

Recommend Projects

  • React photo React

    A declarative, efficient, and flexible JavaScript library for building user interfaces.

  • Vue.js photo Vue.js

    🖖 Vue.js is a progressive, incrementally-adoptable JavaScript framework for building UI on the web.

  • Typescript photo Typescript

    TypeScript is a superset of JavaScript that compiles to clean JavaScript output.

  • TensorFlow photo TensorFlow

    An Open Source Machine Learning Framework for Everyone

  • Django photo Django

    The Web framework for perfectionists with deadlines.

  • D3 photo D3

    Bring data to life with SVG, Canvas and HTML. 📊📈🎉

Recommend Topics

  • javascript

    JavaScript (JS) is a lightweight interpreted programming language with first-class functions.

  • web

    Some thing interesting about web. New door for the world.

  • server

    A server is a program made to process requests and deliver data to clients.

  • Machine learning

    Machine learning is a way of modeling and interpreting data that allows a piece of software to respond intelligently.

  • Game

    Some thing interesting about game, make everyone happy.

Recommend Org

  • Facebook photo Facebook

    We are working to build community through open source technology. NB: members must have two-factor auth.

  • Microsoft photo Microsoft

    Open source projects and samples from Microsoft.

  • Google photo Google

    Google ❤️ Open Source for everyone.

  • D3 photo D3

    Data-Driven Documents codes.