Skip to content

Solving Optimal Control Problems with Symbolic Universal Differential Equations

This tutorial is based on SciMLSensitivity.jl tutorial. Instead of using a classical NN architecture, here we will combine the NN with a symbolic expression from DynamicExpressions.jl (the symbolic engine behind SymbolicRegression.jl and PySR).

Here we will solve a classic optimal control problem with a universal differential equation. Let

x=u3(t)

where we want to optimize our controller u(t) such that the following is minimized:

L(θ)=i(4x(ti)2+2x(ti)2+u(ti)2)

where i is measured on (0,8) at 0.01 intervals. To do this, we rewrite the ODE in first order form:

x=vv=u3(t)

and thus

L(θ)=i(4x(ti)2+2v(ti)2+u(ti)2)

is our loss function on the first order system. We thus choose a neural network form for u and optimize the equation with respect to this loss. Note that we will first reduce control cost (the last term) by 10x in order to bump the network out of a local minimum. This looks like:

Package Imports

julia
using Lux, ComponentArrays, OrdinaryDiffEq, Optimization, OptimizationOptimJL,
      OptimizationOptimisers, SciMLSensitivity, Statistics, Printf, Random
using DynamicExpressions, SymbolicRegression, MLJ, SymbolicUtils, Latexify
using CairoMakie

Helper Functions

julia
function plot_dynamics(sol, us, ts)
    fig = Figure()
    ax = CairoMakie.Axis(fig[1, 1]; xlabel=L"t")
    ylims!(ax, (-6, 6))

    lines!(ax, ts, sol[1, :]; label=L"u_1(t)", linewidth=3)
    lines!(ax, ts, sol[2, :]; label=L"u_2(t)", linewidth=3)

    lines!(ax, ts, vec(us); label=L"u(t)", linewidth=3)

    axislegend(ax; position=:rb)

    return fig
end
plot_dynamics (generic function with 1 method)

Training a Neural Network based UDE

Let's setup the neural network. For the first part, we won't do any symbolic regression. We will plain and simple train a neural network to solve the optimal control problem.

julia
rng = Xoshiro(0)
tspan = (0.0, 8.0)

mlp = Chain(Dense(1 => 4, gelu), Dense(4 => 4, gelu), Dense(4 => 1))

function construct_ude(mlp, solver; kwargs...)
    return @compact(; mlp, solver, kwargs...) do x_in, ps
        x, ts, ret_sol = x_in

        function dudt(du, u, p, t)
            u₁, u₂ = u
            du[1] = u₂
            du[2] = mlp([t], p)[1]^3
            return
        end

        prob = ODEProblem{true}(dudt, x, extrema(ts), ps.mlp)

        sol = solve(prob, solver; saveat=ts,
            sensealg=QuadratureAdjoint(; autojacvec=ReverseDiffVJP(true)), kwargs...)

        us = mlp(reshape(ts, 1, :), ps.mlp)
        ret_sol === Val(true) && return sol, us
        return Array(sol), us
    end
end

ude = construct_ude(mlp, Vern9(); abstol=1e-10, reltol=1e-10);
┌ Warning: No @return macro found in the function body. This will lead to the generation of inefficient code.
└ @ Lux /var/lib/buildkite-agent/builds/gpuci-1/julialang/lux-dot-jl/src/helpers/compact.jl:348

Here we are going to tuse the same configuration for testing, but this is to show that we can setup them up with different ode solve configurations

julia
ude_test = construct_ude(mlp, Vern9(); abstol=1e-10, reltol=1e-10);

function train_model_1(ude, rng, ts_)
    ps, st = Lux.setup(rng, ude)
    ps = ComponentArray{Float64}(ps)
    stateful_ude = StatefulLuxLayer(ude, st)

    ts = collect(ts_)

    function loss_adjoint(θ)
        x, us = stateful_ude(([-4.0, 0.0], ts, Val(false)), θ)
        return mean(abs2, 4 .- x[1, :]) + 2 * mean(abs2, x[2, :]) + 0.1 * mean(abs2, us)
    end

    callback = function (state, l)
        state.iter % 50 == 1 && @printf "Iteration: %5d\tLoss: %10g\n" state.iter l
        return false
    end

    optf = OptimizationFunction((x, p) -> loss_adjoint(x), AutoZygote())
    optprob = OptimizationProblem(optf, ps)
    res1 = solve(optprob, Optimisers.Adam(0.001); callback, maxiters=500)

    optprob = OptimizationProblem(optf, res1.u)
    res2 = solve(optprob, LBFGS(); callback, maxiters=100)

    return StatefulLuxLayer{true}(ude, res2.u, st)
end

trained_ude = train_model_1(ude, rng, 0.0:0.01:8.0)
┌ Warning: Lux.apply(m::Lux.AbstractExplicitLayer, x::AbstractArray{<:ReverseDiff.TrackedReal}, ps, st) input was corrected to Lux.apply(m::Lux.AbstractExplicitLayer, x::ReverseDiff.TrackedArray}, ps, st).

│ 1. If this was not the desired behavior overload the dispatch on `m`.

│ 2. This might have performance implications. Check which layer was causing this problem using `Lux.Experimental.@debug_mode`.
└ @ LuxReverseDiffExt /var/lib/buildkite-agent/builds/gpuci-1/julialang/lux-dot-jl/ext/LuxReverseDiffExt.jl:25
Iteration:     1	Loss:    63.9746
Iteration:    51	Loss:    62.6926
Iteration:   101	Loss:    52.1835
Iteration:   151	Loss:    32.4412
Iteration:   201	Loss:    30.8709
Iteration:   251	Loss:    29.8951
Iteration:   301	Loss:    28.9481
Iteration:   351	Loss:    27.9728
Iteration:   401	Loss:     27.021
Iteration:   451	Loss:     26.129
Iteration:     1	Loss:    25.1138
Iteration:    51	Loss:    12.2671
julia
sol, us = ude_test(([-4.0, 0.0], 0.0:0.01:8.0, Val(true)), trained_ude.ps, trained_ude.st)[1];
plot_dynamics(sol, us, 0.0:0.01:8.0)

Now that the system is in a better behaved part of parameter space, we return to the original loss function to finish the optimization:

julia
function train_model_2(stateful_ude::StatefulLuxLayer, ts_)
    ts = collect(ts_)

    function loss_adjoint(θ)
        x, us = stateful_ude(([-4.0, 0.0], ts, Val(false)), θ)
        return mean(abs2, 4 .- x[1, :]) .+ 2 * mean(abs2, x[2, :]) .+ mean(abs2, us)
    end

    callback = function (state, l)
        state.iter % 10 == 1 && @printf "Iteration: %5d\tLoss: %10g\n" state.iter l
        return false
    end

    optf = OptimizationFunction((x, p) -> loss_adjoint(x), AutoZygote())
    optprob = OptimizationProblem(optf, stateful_ude.ps)
    res2 = solve(optprob, LBFGS(); callback, maxiters=100)

    return StatefulLuxLayer(stateful_ude.model, res2.u, stateful_ude.st)
end

trained_ude = train_model_2(trained_ude, 0.0:0.01:8.0)
┌ Warning: Lux.apply(m::Lux.AbstractExplicitLayer, x::AbstractArray{<:ReverseDiff.TrackedReal}, ps, st) input was corrected to Lux.apply(m::Lux.AbstractExplicitLayer, x::ReverseDiff.TrackedArray}, ps, st).

│ 1. If this was not the desired behavior overload the dispatch on `m`.

│ 2. This might have performance implications. Check which layer was causing this problem using `Lux.Experimental.@debug_mode`.
└ @ LuxReverseDiffExt /var/lib/buildkite-agent/builds/gpuci-1/julialang/lux-dot-jl/ext/LuxReverseDiffExt.jl:25
Iteration:     1	Loss:    12.7237
Iteration:    11	Loss:    12.7014
Iteration:    21	Loss:     12.696
Iteration:    31	Loss:    12.6788
Iteration:    41	Loss:    12.6518
Iteration:    51	Loss:    12.6165
Iteration:    61	Loss:    12.6026
Iteration:    71	Loss:    12.5638
Iteration:    81	Loss:    12.5562
Iteration:    91	Loss:    12.5425
julia
sol, us = ude_test(([-4.0, 0.0], 0.0:0.01:8.0, Val(true)), trained_ude.ps, trained_ude.st)[1];
plot_dynamics(sol, us, 0.0:0.01:8.0)

Symbolic Regression

Ok so now we have a trained neural network that solves the optimal control problem. But can we replace Dense(4 => 4, gelu) with a symbolic expression? Let's try!

Data Generation for Symbolic Regression

First, we need to generate data for the symbolic regression.

julia
ts = reshape(collect(0.0:0.1:8.0), 1, :)

X_train = mlp[1](ts, trained_ude.ps.mlp.layer_1, trained_ude.st.mlp.layer_1)[1]
4×81 Matrix{Float64}:
 -0.0301877  -0.0164291  -0.00818525  -0.00371752   -0.00153096   -0.000568211  -0.000188803  -5.57746e-5  -1.45436e-5   -3.323e-6     -6.60373e-7   -1.13291e-7   -1.6652e-8   -2.08116e-9  -2.19483e-10  -1.93832e-11  -1.42247e-12  -8.60818e-14  -4.2626e-15   -1.71386e-16  -5.55197e-18  -1.43788e-19  -2.95414e-21  -4.7774e-23   -6.03431e-25  -5.90684e-27  -4.44621e-29  -2.55356e-31  -1.11028e-33  -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0         -0.0          -0.0          -0.0          -0.0
  0.332219    0.229961    0.138639     0.0591181    -0.00807984   -0.0628106    -0.105319     -0.136222    -0.156465     -0.167257     -0.169993     -0.166168     -0.157291    -0.144809    -0.130044     -0.114146     -0.0980639    -0.0825397    -0.0681093    -0.0551212    -0.0437618    -0.034085     -0.0260429    -0.019516     -0.0143398    -0.0103273    -0.00728653   -0.0050342    -0.00340383   -0.00225098   -0.00145498   -0.00091862   -0.000566108  -0.000340279  -0.000199353  -0.000113745  -6.31591e-5  -3.41033e-5   -1.78927e-5   -9.11446e-6   -4.50423e-6   -2.15774e-6  -1.0012e-6    -4.49607e-7   -1.9525e-7    -8.19294e-8   -3.31918e-8   -1.29721e-8   -4.8868e-9    -1.77304e-9  -6.1907e-10   -2.07841e-10  -6.70401e-11  -2.07585e-11  -6.16536e-12  -1.75494e-12  -4.78355e-13  -1.24757e-13  -3.11061e-14  -7.40859e-15  -1.68412e-15  -3.6509e-16  -7.54147e-17  -1.48314e-17  -2.77468e-18  -4.93394e-19  -8.33224e-20  -1.33523e-20  -2.02868e-21  -2.91994e-22  -3.97811e-23  -5.12577e-24  -6.24108e-25  -7.17492e-26  -7.7816e-27  -7.95523e-28  -7.6596e-29  -6.94012e-30  -5.91252e-31  -4.73217e-32  -3.55524e-33
 -0.0049004  -0.0027596  -0.0014828   -0.000758704  -0.000368884  -0.000170051  -7.41582e-5   -3.05227e-5  -1.18292e-5   -4.30649e-6   -1.46924e-6   -4.68618e-7   -1.39396e-7  -3.85773e-8  -9.90839e-9   -2.35614e-9   -5.17439e-10  -1.04691e-10  -1.94659e-11  -3.31807e-12  -5.17209e-13  -7.35422e-14  -9.51521e-15  -1.11746e-15  -1.18821e-16  -1.14109e-17  -9.87249e-19  -7.67594e-20  -5.34991e-21  -3.33416e-22  -1.85337e-23  -9.16617e-25  -4.02322e-26  -1.56326e-27  -5.36376e-29  -1.62107e-30  -4.3046e-32  -1.00179e-33  -0.0          -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0         -0.0          -0.0          -0.0          -0.0
 -0.0157456  -0.0115007  -0.00823904  -0.00578652   -0.00398224   -0.0026839    -0.00177045   -0.00114237  -0.000720545  -0.000443974  -0.000267055  -0.000156707  -8.96437e-5  -4.99552e-5  -2.70995e-5   -1.43002e-5   -7.3352e-6    -3.65465e-6   -1.76736e-6   -8.28948e-7   -3.76814e-7   -1.65883e-7   -7.06681e-8   -2.91116e-8   -1.15878e-8   -4.45352e-9   -1.65136e-9   -5.90316e-10  -2.03284e-10  -6.73859e-11  -2.14856e-11  -6.58431e-12  -1.93786e-12  -5.47334e-13  -1.48241e-13  -3.84712e-14  -9.5592e-15  -2.27244e-15  -5.16434e-16  -1.12113e-16  -2.32314e-17  -4.5914e-18  -8.64823e-19  -1.55127e-19  -2.64783e-20  -4.29736e-21  -6.62649e-22  -9.70069e-23  -1.34717e-23  -1.7734e-24  -2.21115e-25  -2.6093e-26   -2.91196e-27  -3.07091e-28  -3.05796e-29  -2.87306e-30  -2.54487e-31  -2.12355e-32  -1.668e-33    -0.0          -0.0          -0.0         -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0          -0.0         -0.0          -0.0         -0.0          -0.0          -0.0          -0.0

This is the training input data. Now we generate the targets

julia
Y_train = mlp[2](X_train, trained_ude.ps.mlp.layer_2, trained_ude.st.mlp.layer_2)[1]
4×81 Matrix{Float64}:
  0.38156    0.227224   0.11173    0.0286932  -0.0289608  -0.0677555  -0.0930538  -0.10891    -0.11819    -0.122804    -0.123971    -0.122457    -0.11877    -0.113302   -0.106419   -0.098503   -0.0899583  -0.0811913  -0.0725808  -0.0644478  -0.0570352  -0.0504989  -0.0449111  -0.0402719  -0.0365263  -0.0335822  -0.0313275  -0.0296443  -0.0284189  -0.0275489  -0.0269465  -0.0265397  -0.0262721  -0.0261005  -0.0259933  -0.0259282  -0.0258897  -0.0258676  -0.0258553  -0.0258486  -0.0258451  -0.0258433  -0.0258425  -0.025842   -0.0258418  -0.0258418  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417  -0.0258417
 -0.164269  -0.169948  -0.161416  -0.141234   -0.113743   -0.0837596  -0.0555671  -0.0323158  -0.0158109  -0.00659323  -0.00420112  -0.00750861  -0.0150603  -0.0253513  -0.0370271  -0.0489965  -0.0604707  -0.0709469  -0.0801606  -0.088025   -0.0945736  -0.0999118  -0.104181   -0.107536   -0.110126   -0.112092   -0.113557   -0.114628   -0.115397   -0.115937   -0.116308   -0.116557   -0.11672    -0.116825   -0.11689    -0.116929   -0.116953   -0.116966   -0.116974   -0.116978   -0.11698    -0.116981   -0.116981   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982   -0.116982
  0.744146   0.438785   0.213599   0.0605459  -0.0363051  -0.0939244  -0.126408   -0.143826   -0.152618   -0.15651     -0.15742     -0.15617     -0.152997   -0.147911   -0.140925   -0.132186   -0.122016   -0.11089    -0.0993637  -0.0879931  -0.0772614  -0.0675329  -0.0590335  -0.0518574  -0.0459882  -0.0413297  -0.0377361  -0.035039   -0.033068   -0.0316647  -0.0306912  -0.0300331  -0.0295996  -0.0293216  -0.0291479  -0.0290423  -0.0289799  -0.0289441  -0.0289241  -0.0289132  -0.0289076  -0.0289047  -0.0289032  -0.0289026  -0.0289022  -0.0289021  -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902   -0.028902
  0.107454   0.120282   0.128643   0.133777    0.136815    0.138592    0.139644    0.140279    0.140662    0.140882     0.140989     0.141016     0.140986    0.140918    0.140826    0.140721    0.140612    0.140506    0.140406    0.140316    0.140237    0.14017     0.140114    0.140069    0.140033    0.140005    0.139984    0.139968    0.139957    0.139949    0.139943    0.139939    0.139937    0.139935    0.139934    0.139934    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933    0.139933

Fitting the Symbolic Expression

We will follow the example from SymbolicRegression.jl docs to fit the symbolic expression.

julia
srmodel = MultitargetSRRegressor(;
    binary_operators=[+, -, *, /], niterations=100, save_to_file=false);

One important note here is to transpose the data because that is how MLJ expects the data to be structured (this is in contrast to how Lux or SymbolicRegression expects the data)

julia
mach = machine(srmodel, X_train', Y_train')
fit!(mach; verbosity=0)
r = report(mach)
best_eq = [r.equations[1][r.best_idx[1]], r.equations[2][r.best_idx[2]],
    r.equations[3][r.best_idx[3]], r.equations[4][r.best_idx[4]]]
4-element Vector{DynamicExpressions.EquationModule.Node{Float64}}:
 (((((x₃ + 0.7475997630866819) + (x₂ + x₁)) * x₂) + -0.025904544947234737) - x₁) - x₁
 (((x₂ + x₁) * ((x₁ - 0.47333932991838423) + (x₂ / 0.9114683907120399))) + -0.1170202364047835) + x₃
 (((x₂ * (((x₂ + x₁) * 2.955724917341977) - -1.2513889484908958)) - x₁) - 0.028824686531356093) - x₁
 (((0.13994958435874832 - (((x₄ + (x₄ - -0.005212641736023803)) * x₂) + (x₃ / -2.1890378056697166))) + x₃) + x₃) + x₁

Let's see the expressions that SymbolicRegression.jl found. In case you were wondering, these expressions are not hardcoded, it is live updated from the output of the code above using Latexify.jl and the integration of SymbolicUtils.jl with DynamicExpressions.jl.

(x3+0.7476+x2+x1)x20.025905x1x1(x2+x1)(x10.47334+x20.91147)0.11702+x3x2((x2+x1)2.9557+1.2514)x10.028825x10.13995((x4+x4+0.0052126)x2+x32.189)+x3+x3+x1

Combining the Neural Network with the Symbolic Expression

Now that we have the symbolic expression, we can combine it with the neural network to solve the optimal control problem. but we do need to perform some finetuning.

julia
hybrid_mlp = Chain(Dense(1 => 4, gelu),
    DynamicExpressionsLayer(OperatorEnum(; binary_operators=[+, -, *, /]), best_eq),
    Dense(4 => 1); disable_optimizations=true)
Chain(
    layer_1 = Dense(1 => 4, gelu),      # 8 parameters
    layer_2 = Chain(
        layer_1 = Parallel(
            layer_1 = DynamicExpressionNode((((((x₃ + 0.7475997630866819) + (x₂ + x₁)) * x₂) + -0.025904544947234737) - x₁) - x₁),  # 2 parameters
            layer_2 = DynamicExpressionNode((((x₂ + x₁) * ((x₁ - 0.47333932991838423) + (x₂ / 0.9114683907120399))) + -0.1170202364047835) + x₃),  # 3 parameters
            layer_3 = DynamicExpressionNode((((x₂ * (((x₂ + x₁) * 2.955724917341977) - -1.2513889484908958)) - x₁) - 0.028824686531356093) - x₁),  # 3 parameters
            layer_4 = DynamicExpressionNode((((0.13994958435874832 - (((x₄ + (x₄ - -0.005212641736023803)) * x₂) + (x₃ / -2.1890378056697166))) + x₃) + x₃) + x₁),  # 3 parameters
        ),
        layer_2 = WrappedFunction(__stack1),
    ),
    layer_3 = Dense(4 => 1),            # 5 parameters
)         # Total: 24 parameters,
          #        plus 0 states.

There you have it! It is that easy to take the fitted Symbolic Expression and combine it with a neural network. Let's see how it performs before fintetuning.

julia
hybrid_ude = construct_ude(hybrid_mlp, Vern9(); abstol=1e-10, reltol=1e-10);

We want to reuse the trained neural network parameters, so we will copy them over to the new model

julia
st = Lux.initialstates(rng, hybrid_ude)
ps = (;
    mlp=(; layer_1=trained_ude.ps.mlp.layer_1,
        layer_2=Lux.initialparameters(rng, hybrid_mlp[2]),
        layer_3=trained_ude.ps.mlp.layer_3))
ps = ComponentArray(ps)

sol, us = hybrid_ude(([-4.0, 0.0], 0.0:0.01:8.0, Val(true)), ps, st)[1];
plot_dynamics(sol, us, 0.0:0.01:8.0)

Now that does perform well! But we could finetune this model very easily. We will skip that part on CI, but you can do it by using the same training code as above.

Appendix

julia
using InteractiveUtils
InteractiveUtils.versioninfo()

if @isdefined(LuxCUDA) && CUDA.functional()
    println()
    CUDA.versioninfo()
end

if @isdefined(LuxAMDGPU) && LuxAMDGPU.functional()
    println()
    AMDGPU.versioninfo()
end
Julia Version 1.10.3
Commit 0b4590a5507 (2024-04-30 10:59 UTC)
Build Info:
  Official https://julialang.org/ release
Platform Info:
  OS: Linux (x86_64-linux-gnu)
  CPU: 48 × AMD EPYC 7402 24-Core Processor
  WORD_SIZE: 64
  LIBM: libopenlibm
  LLVM: libLLVM-15.0.7 (ORCJIT, znver2)
Threads: 8 default, 0 interactive, 4 GC (on 2 virtual cores)
Environment:
  JULIA_CPU_THREADS = 2
  JULIA_AMDGPU_LOGGING_ENABLED = true
  JULIA_DEPOT_PATH = /root/.cache/julia-buildkite-plugin/depots/01872db4-8c79-43af-ab7d-12abac4f24f6
  LD_LIBRARY_PATH = /usr/local/nvidia/lib:/usr/local/nvidia/lib64
  JULIA_PKG_SERVER = 
  JULIA_NUM_THREADS = 8
  JULIA_CUDA_HARD_MEMORY_LIMIT = 100%
  JULIA_PKG_PRECOMPILE_AUTO = 0
  JULIA_DEBUG = Literate

This page was generated using Literate.jl.