What is Julia?#

Julia is a new-ish programming language targeting scientific and high-performance computing applications. Like Python or MATLAB, it is a dynamic language—you don’t need to declare the types of variables like you would in C, C++, Fortran, etc.. Unlike Python, it is a compiled language, as every function is compiled to machine code the first time it is called. The aim of Julia is to provide users with a programming language that combines the convience of Python or MATLAB with the speed of C or Fortran. The official Julia documentation has a nice page comparing Julia to other popular languages if you’d like to know more.

One of the most attractive features of Julia for OpenMDAO users is its good support for automatic differentiation (AD). There are many (perhaps too many!) AD libraries available for the language, all with particular strengths and weaknesses. Happily, an excellent library called DifferentiationInterface.jl exists which allows us to use any of these numerious AD library via a uniform API, which makes switching from one to another a matter of changing just one or two lines of code.

What is OpenMDAO.jl?#

NOTE This page is a brief, Python-centric introduction to using Julia within OpenMDAO and OpenMDAO.jl in particular. Please check out the OpenMDAO.jl docs for all the details.

OpenMDAO.jl is a Python+Julia software package that allows users to create OpenMDAO Components out of Julia code, and optionally use a Julia AD package to differentiate the Component. There are currently two possible approaches to using AD with Julia:

  1. Manual Approach: Users can write the equivalent of a compute_partials or linearize method in Julia that calls the desired AD library manually. This is quite flexible and conceptually straightforward, but requires more work on the part of the user to manage the complexities of translating the inputs and outputs vectors OpenMDAO provides to something that the AD library can work with.

  2. Automatic Approach: Users can write a relatively simple “callback function” that implements the math behind the Component they wish to create, and allow OpenMDAO.jl to call the AD for them. This is much less work for the user, but does require learning a bit more about the Julia language. This approach is also currently limited to creating ExplicitComponents (though the callback function can contain implicitness via the ImplicitAD.jl package).

This documentation will focus on approach number 2, the Automatic Approach, which consists of essentially three steps:

  1. Create a user-defined callback function

  2. Create an OpenMDAO.jl AbstractADExplicitComp struct that will keep track the user-defined function, the desired AD library, and other optional features like units, tags, etc..

  3. Pass the AbstractADExplicitComp to omjlcomps (“OpenMDAO Julia Components”), a Python library that knows how to create OpenMDAO Components from an OpenMDAO.jl AbstractADExplicitComp.

After completing those three steps, you’ll be left with an OpenMDAO Component that you can include anywhere in an OpenMDAO model. We’ll go through each of the three steps below. For more details, have a look at the OpenMDAO.jl docs.

Guided Tutorial#

Step 1: The User-Defined Function#

The user-defined function must follow one of two forms: either it can be an “in-place” function that writes its outputs to an output vector, or it can be an “out-of-place” function that returns a single output vector. The function must be written in Julia, of course, and anything that happens inside the function is up to you, of course. Both the in-place and out-of-place form must also have a params argument that will contain inputs that are needed for the calculation, but won’t be differentiated (the equivalent of OpenMDAO options). So, an example of an in-place function would be

function f_in_place!(Y, X, params)
   # calculate stuff with X and params, storing result in Y
   return nothing
end

where X is the input vector and Y is the output vector. (The function doesn’t have to return nothing, but any returned value will be ignored, so I like to include return nothing to make it clear that the return value doesn’t matter.) An out-of-place function would look like

function f_out_of_place(X, params)
   # calculate stuff with X and params, returning Y
   return Y
end

where again X is the input vector and Y is the output vector.

Now, the X and Y arguments of those functions must not be plain Julia Vectors, but ComponentVectors from the ComponentArrays.jl package. What are those? They are objects provided by the ComponentArrays.jl package that act like Vectors, but allow the user to define names for each part (“component”) of the vector. For example, we can do this from the Julia REPL:

julia> using ComponentArrays: ComponentVector

julia> x1 = ComponentVector(foo=-1.0, bar=-2.0, baz=-3.0)
ComponentVector{Float64}(foo = -1.0, bar = -2.0, baz = -3.0)

julia> @show x1 x1[3] x1.foo x1[:foo]
x1 = (foo = -1.0, bar = -2.0, baz = -3.0)
x1[3] = -3.0
x1.foo = -1.0
x1[:foo] = -1.0
-1.0

julia> 

Notice that we can get, say, the third value of x1 the usual way (x1[3]), but also by referring to the foo field value via x1.foo and by indexing the ComponentVector with the symbol :foo (x1[:foo]).

Each of the components in x1 are scalars, but they don’t have to be:

julia> x2 = ComponentVector(foo=-1.0, bar=1:4, baz=reshape(5:10, 2, 3))
ComponentVector{Float64}(foo = -1.0, bar = [1.0, 2.0, 3.0, 4.0], baz = [5.0 7.0 9.0; 6.0 8.0 10.0])

julia> @show x2 x2[:foo] x2[:bar] x2[:baz]
x2 = (foo = -1.0, bar = [1.0, 2.0, 3.0, 4.0], baz = [5.0 7.0 9.0; 6.0 8.0 10.0])
x2[:foo] = -1.0
x2[:bar] = [1.0, 2.0, 3.0, 4.0]
x2[:baz] = [5.0 7.0 9.0; 6.0 8.0 10.0]
2×3 Matrix{Float64}:
 5.0  7.0   9.0
 6.0  8.0  10.0

julia> 

In x2, the foo component is a scalar, bar refers to a Vector (aka a 1D Array) and baz refers to a Matrix (aka a 2D Array). But x2 still “looks like” a Vector:

julia> @show x2[3]  # will give the third value of `x2`, which happens to be the second value of x2[:bar]
x2[3] = 2.0
2.0

julia> @show ndims(x2)  # Should be 1, since a Vector is 1-dimensional
ndims(x2) = 1
1

julia> @show length(x2)  # length(x2) gives the total number of entries in `x2`, aka 1 + 4 + 2*3 = 11
length(x2) = 11
11

julia> @show size(x2)  # size is a length-1 tuple since a Vector has just one dimension
size(x2) = (11,)
(11,)

julia> 

Now, how will we use ComponentVectors here? We’ll use them to define the names and sizes of all the inputs and outputs to our component. For example, if we wanted to implement the Paraboloid Example, which minimizes the function

\[ f(x,y) = (x-3.0)^2 + x \cdot y + (y+4.0)^2 - 3.0 \]

using OpenMDAO.jl, we could create an input vector called X_ca that looks like this:

julia> X_ca = ComponentVector(x=1.0, y=1.0)
ComponentVector{Float64}(x = 1.0, y = 1.0)

julia> 

That has two scalar components, x and y, that correspond to the two inputs to the paraboloid function. We could similarly use an output vector Y_ca that looks like this:

julia> Y_ca = ComponentVector(f_xy=0.0)
ComponentVector{Float64}(f_xy = 0.0)

julia> 

with the single output f_xy.

We could type in a callback function that implements the paraboloid equation above into the REPL if we want:

julia> function f_paraboloid!(Y, X, params)
           # Get the inputs:
           x = @view(X[:x])
           y = @view(X[:y])
           # Could also do this:
           # x = X.x
           # y = X.y
           # or even this
           # (; x, y) = X

           # Get the output:
           f_xy = @view(Y[:f_xy])
           # Again, could also do this:
           # f_xy = Y.f_xy
           # or
           # (; f_xy) = Y

           # Do the calculation:
           @. f_xy = (x - 3.0)^2 + x*y + (y + 4.0)^2 - 3.0

           # No return value for in-place callback function.
           return nothing
       end
f_paraboloid! (generic function with 1 method)

julia> 

And then test it out with the X_ca and Y_ca ComponentVectors we created earlier:

julia> f_paraboloid!(Y_ca, X_ca, nothing)

julia> X_ca
ComponentVector{Float64}(x = 1.0, y = 1.0)

julia> Y_ca
ComponentVector{Float64}(f_xy = 27.0)

julia> 

Some remarks about the callback function:

  • The @view macro is used when extracting the inputs and outputs from the X_ca and Y_ca ComponentVectors. This creates a view into the original ComponentVector, instead of a new array with a copy of the original data, which avoids unnecessary allocations and (for the outputs) allows modifications to the view to be reflected in the Y_ca array. In this example everything is a scalar, so no allocations would have happened anyway. But it doesn’t hurt to use @view: it’s a good habit to get into, and it allows us to use the @. broadcasting macro with the scalar f_xy output.

  • The params argument is not used in this example, but it is still required, since the code in OpenMDAO.jl will expect it. We provided a nothing value for it, but we could have given it anything and not changed the result of the function, of course.

An out-of-place version of f_paraboloid! would look like this:

julia> function f_paraboloid(X, params)
           x = @view(X[:x])
           y = @view(X[:y])
           f_xy = @. (x - 3.0)^2 + x*y + (y + 4.0)^2 - 3.0
           return ComponentVector(f_xy=f_xy)
       end
f_paraboloid (generic function with 1 method)

julia> f_paraboloid(X_ca, nothing)
ComponentVector{Float64}(f_xy = 27.0)

julia> 

But we only need one or the other.

Step 2: The OpenMDAO.jl AbstractADExplicitComp#

So far we’ve done a lot of typing in the Julia REPL, which is great for testing, but not great for writing code that you’ll use more than once. Where should we store this Julia code that we’re creating, and how will we call it from Python? This section will explain my favorite way to do that. Along the way we’ll learn how to choose a Julia AD library, and create an AbstractADExplicitComp that will eventually be used to instantiate an ExplicitComponent that can be added to a OpenMDAO model.

A Julia Package for our Paraboloid#

We’ll store our Julia code for the paraboloid component in a brand new Julia package called MyParaboloidPackage. Creating a new Julia package is very easy from the REPL: just type ] to enter the Pkg REPL (look for the pkg> in the prompt to ensure you’re in “Pkg mode”), then type generate MyParaboloidPackage to create the package:

(openmdao_jl_dev) pkg> generate MyParaboloidPackage
  Generating  project MyParaboloidPackage:
    MyParaboloidPackage/Project.toml
    MyParaboloidPackage/src/MyParaboloidPackage.jl

(openmdao_jl_dev) pkg> 

(Keep this REPL open, as we’ll be using it again later.) As the output suggests, running that command creates a directory called MyParaboloidPackage in the current working directory, along with two files. The Project.toml file contains metadata about the Julia project (name, version number, a UUID to distinguish it from other Julia packages, dependencies, etc.). The src/MyParaboloidPackage.jl contains the source code of the package, and at the moment looks like this:

module MyParaboloidPackage

greet() = print("Hello World!")

end # module MyParaboloidPackage

We’ll put our paraboloid code in src/MyParaboloidPackage.jl by copy-pasting the definition of f_paraboloid! into the file. Now it should look like this (with the greet function removed):

module MyParaboloidPackage

function f_paraboloid!(Y, X, params)
    x = @view(X[:x])
    y = @view(X[:y])
    f_xy = @view(Y[:f_xy])
    @. f_xy = (x - 3.0)^2 + x*y + (y + 4.0)^2 - 3.0
    return nothing
end

end # module MyParaboloidPackage

Now, let’s see if we can create that AbstractADExplicitComp that we’ve been talking about. There are a few different flavors of AbstractADExplicitComp structs that we could use, but for this example, we’ll use DenseADExplicitComp, one that’s designed for explicit components with a dense Jacobian. If we check out the documentation for DenseADExplicitComp constructor, we’ll see that we need four things to get this done:

  • An ad_backend, to be explained shortly

  • A callback function f!

  • Y_ca, a ComponentVector of outputs

  • X_ca, a ComponentVector of inputs

In the previous section we learned how to create X_ca and Y_ca, and our callback function f! will be the f_paraboloid!. So the only missing ingrediant is the ad_backend.

The ad_backend is a Julia object that will tell the DenseADExplicitComp which AD library we want to use to differentiate our component. The object must be a subtype of AbstractADType, an abstract type provided by the ADTypes.jl package. That package is not part of the Julia Standard Library, so we’ll need to install it. How do we do that? The easiest way is from the Julia Pkg REPL, aka the thing we used to generate the MyParaboloidPackage earlier. First, type ] to ensure we’re in the Pkg REPL (again, look for pkg> in the prompt), then do this:

(openmdao_jl_dev) pkg> activate ./MyParaboloidPackage
  Activating project at `~/projects/openmdao_jl_dev/ved/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage`

(MyParaboloidPackage) pkg> 

That “activates” the environment associated with MyParaboloidPackage, which tells the Julia Pkg package manager that we want to mess around with the dependencies of it and not some other environment or package. (The Pkg REPL reminds us which environment is “active” by changing the prompt to (MyParaboloidPackage) pkg>.) We can check the status of MyParaboloidPackage via the status command:

(MyParaboloidPackage) pkg> status
Project MyParaboloidPackage v0.1.0
Status `~/projects/openmdao_jl_dev/ved/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage/Project.toml` (empty project)

(MyParaboloidPackage) pkg> 

Unsuprisingly there isn’t much going on with it yet.

The next step is to add the dependencies we need. We’ve already decided that we want ADTypes.jl, but we’ll also need:

  • ComponentArrays.jl: for the ComponentVectors we’ll provide to the DenseADExplicitComp constructor

  • OpenMDAOCore.jl: for the DenseADExplicitComp itself

  • ForwardDiff.jl: a forward-mode AD library that we’ll use for this example

So, let’s install those using the add command from the Pkg REPL:

(MyParaboloidPackage) pkg> add ADTypes ComponentArrays OpenMDAOCore ForwardDiff
   Resolving package versions...
   Installed DifferentiationInterface  v0.7.14
      Compat entries added for ADTypes, ComponentArrays, OpenMDAOCore, ForwardDiff
    Updating `~/projects/openmdao_jl_dev/ved/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage/Project.toml`
  [47edcb42] + ADTypes v1.21.0
  [b0b7db55] + ComponentArrays v0.15.31
  [f6369f11] + ForwardDiff v1.3.1
  [24d19c10] + OpenMDAOCore v0.3.3
    Updating `~/projects/openmdao_jl_dev/ved/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage/Manifest.toml`

# Lots of output removed for clarity
(MyParaboloidPackage) pkg> 

Now if we try the status command again, we’ll see this:

(MyParaboloidPackage) pkg> status
Project MyParaboloidPackage v0.1.0
Status `~/projects/openmdao_jl_dev/ved/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage/Project.toml`
  [47edcb42] ADTypes v1.21.0
  [b0b7db55] ComponentArrays v0.15.31
  [f6369f11] ForwardDiff v1.3.1
  [24d19c10] OpenMDAOCore v0.3.3

(MyParaboloidPackage) pkg> 

That shows the dependencies of MyParaboloidPackage, and their version numbers. Excellent.

Now we should have everything we need to create the DenseADExplicitComp. I think it’s easiest to create a small “helper function” that we can call later that returns our DenseADExplicitComp that looks like this:

function get_paraboloid_comp()
    ad_backend = ADTypes.AutoForwardDiff()
    X_ca = ComponentVector(x=1.0, y=1.0)
    Y_ca = ComponentVector(f_xy=0.0)
    comp = OpenMDAOCore.DenseADExplicitComp(ad_backend, f_paraboloid!, Y_ca, X_ca)
    return comp
end

The ADTypes.AutoForwardDiff() means that we’ll use ForwardDiff.jl, a forward-mode AD package, to differentiate our paraboloid.

We’ll add that function to src/MyParaboloidPackage.jl, along with some import statements that give us access to the packages we just installed. So now it should look like this:

module MyParaboloidPackage

using ADTypes: ADTypes
using ComponentArrays: ComponentVector
using ForwardDiff: ForwardDiff
using OpenMDAOCore: OpenMDAOCore

function f_paraboloid!(Y, X, params)
    x = @view(X[:x])
    y = @view(X[:y])
    f_xy = @view(Y[:f_xy])
    @. f_xy = (x - 3.0)^2 + x*y + (y + 4.0)^2 - 3.0
    return nothing
end

function get_paraboloid_comp()
    ad_backend = ADTypes.AutoForwardDiff()
    X_ca = ComponentVector(x=1.0, y=1.0)
    Y_ca = ComponentVector(f_xy=0.0)
    comp = OpenMDAOCore.DenseADExplicitComp(ad_backend, f_paraboloid!, Y_ca, X_ca)
    return comp
end

end # module MyParaboloidPackage

We can try out the get_paraboloid_comp function from the normal Julia REPL (not the Pkg REPL—just backspace until you get the julia> prompt back):

julia> using MyParaboloidPackage

julia> comp = MyParaboloidPackage.get_paraboloid_comp()

You’ll probibly see a bunch of output describing in excruciating detail all the type information associated with the DenseADExplicitComp we just created. But as long as you see something that starts with OpenMDAOCore.DenseADExplicitComp, we should be good.

Step 3: The omjlcomps.JuliaExplicitComp#

So far we have created a Julia callback function and a DenseADExplicitComp for our paraboloid. The last step is to turn that DenseADExplicitComp into something we can actually incorporate into an OpenMDAO model. We’ll do that using omjlcomps, a Python package that’s part of OpenMDAO.jl.

But first, we need to figure out a nice way of telling Python about our Julia code. The best way of doing that is using juliapkg, a Python package that allows us to add Julia dependencies to a Python package. We could install it manually ourselves using pip install juliapkg, but if instead we install OpenMDAO with the “julia” extra, we’ll get it automatically:

$ pip install openmdao[julia]

We’ll do that from here in a Jupyter notebook using this command:

!pip install openmdao[julia]
Requirement already satisfied: openmdao[julia] in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (3.42.1.dev0)
Requirement already satisfied: networkx>=3.3 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from openmdao[julia]) (3.6.1)
Requirement already satisfied: numpy!=2.4.0,>=2.0 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from openmdao[julia]) (2.3.5)
Requirement already satisfied: packaging in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from openmdao[julia]) (26.0)
Requirement already satisfied: requests in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from openmdao[julia]) (2.32.5)
Requirement already satisfied: scipy>=1.13 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from openmdao[julia]) (1.17.0)
Collecting omjlcomps>=0.2.6 (from openmdao[julia])
  Downloading omjlcomps-0.2.6-py3-none-any.whl.metadata (1.5 kB)
Collecting juliapkg~=0.1.10 (from omjlcomps>=0.2.6->openmdao[julia])
  Downloading juliapkg-0.1.22-py3-none-any.whl.metadata (6.8 kB)
Collecting juliacall~=0.9.13 (from omjlcomps>=0.2.6->openmdao[julia])
  Downloading juliacall-0.9.31-py3-none-any.whl.metadata (4.5 kB)
Requirement already satisfied: filelock<4.0,>=3.16 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from juliapkg~=0.1.10->omjlcomps>=0.2.6->openmdao[julia]) (3.20.3)
Collecting semver<4.0,>=3.0 (from juliapkg~=0.1.10->omjlcomps>=0.2.6->openmdao[julia])
  Downloading semver-3.0.4-py3-none-any.whl.metadata (6.8 kB)
Requirement already satisfied: tomli<3.0,>=2.0 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from juliapkg~=0.1.10->omjlcomps>=0.2.6->openmdao[julia]) (2.4.0)
Collecting tomlkit<0.14,>=0.13.3 (from juliapkg~=0.1.10->omjlcomps>=0.2.6->openmdao[julia])
  Downloading tomlkit-0.13.3-py3-none-any.whl.metadata (2.8 kB)
Requirement already satisfied: charset_normalizer<4,>=2 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from requests->openmdao[julia]) (3.4.4)
Requirement already satisfied: idna<4,>=2.5 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from requests->openmdao[julia]) (3.11)
Requirement already satisfied: urllib3<3,>=1.21.1 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from requests->openmdao[julia]) (2.6.3)
Requirement already satisfied: certifi>=2017.4.17 in /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages (from requests->openmdao[julia]) (2026.1.4)
Downloading omjlcomps-0.2.6-py3-none-any.whl (17 kB)
Downloading juliacall-0.9.31-py3-none-any.whl (12 kB)
Downloading juliapkg-0.1.22-py3-none-any.whl (21 kB)
Downloading semver-3.0.4-py3-none-any.whl (17 kB)
Downloading tomlkit-0.13.3-py3-none-any.whl (38 kB)
Installing collected packages: tomlkit, semver, juliapkg, juliacall, omjlcomps
?25l
   ━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━ 5/5 [omjlcomps]
?25h

Successfully installed juliacall-0.9.31 juliapkg-0.1.22 omjlcomps-0.2.6 semver-3.0.4 tomlkit-0.13.3

Now, we need to use juliapkg to tell our Python environment about MyParaboloidPackage. We do that by creating a file called juliapkg.json in our working directory that points to MyParaboloidPackage:

{"packages": {
    "MyParaboloidPackage": {"uuid": "a013c6e3-dc76-40b3-894d-62618f778e16", "dev": true, "path": "./MyParaboloidPackage"}
    }
}

(To get the correct value for the "uuid" field, copy what’s in the MyParaboloidPackage/Project.toml file. It will be different from what is listed above if you’re trying out these commands yourself.)

That tells juliapkg about MyParaboloidPackage. The "dev": true means that the package will be installed in “develop” mode, meaning any changes in the package will be reflected in code that calls it—it’s the equivalent of pip install -e or pip install --editable in Python.

juliapkg searches for juliapkg.json files in, among other places, each entry in sys.path. This is great for adding Julia dependencies to Python packages: you can add a juliapkg.json to your Python source, and juliapkg will find it automatically if your Python package is installed in the current environment. It also works well for Python scripts run from the command line: since Python automatically adds the directory of the script to sys.path, we could place the juliapkg.json file in the same directory as the run script and it should also be found by juliapkg. But it’s not-so-great for Jupyter notebooks like this one. So we’ll have to manually add what we think is the current working directory to sys.path and cross our fingers.

import os
import sys
d = os.getcwd()
sys.path.append(d)
print(sys.path[-1])
/home/runner/work/OpenMDAO/OpenMDAO/openmdao/docs/openmdao_book/features/experimental

Next, we’ll create a small julia file in the current working directory that will import the get_paraboloid_comp function that we just added to MyParaboloidPackage called paraboloid.jl:

from IPython.display import display, Code
display(Code(os.path.join(d, "paraboloid.jl")))
using MyParaboloidPackage: get_paraboloid_comp

Now we’ll create another small file, this time a Python one, called paraboloid.py that will include that paraboloid.jl julia file using JuliaCall, a Python package for calling Julia code that was automatically installed when we did pip install openmdao[julia]:

display(Code(os.path.join(d, "paraboloid.py")))
import os

# Create a new Julia module that will hold all the Julia code imported into this Python module.
import juliacall
jl = juliacall.newmodule("ParaboloidComponentsStub")

# Get the directory this file is in, then include the `paraboloid.jl` Julia source code.
d = os.path.dirname(os.path.abspath(__file__))
jl.include(os.path.join(d, "paraboloid.jl"))
# Now we have access to everything in `paraboloid.jl` in the `jl` object.

get_paraboloid_comp = jl.get_paraboloid_comp

Now we should be able to use omjlcomps (another Python package that was installed as a side effect of pip install openmdao[julia]) to create a JuliaExplicitComp. We do this by first calling get_paraboloid_comp to get the DenseADExplicitComp struct we created, then pass it as the jlcomp argument to JuliaExplicitComp:

import omjlcomps
from paraboloid import get_paraboloid_comp
jlcomp = get_paraboloid_comp()
comp = omjlcomps.JuliaExplicitComp(jlcomp=jlcomp)
[juliapkg] Found dependencies: /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages/omjlcomps/juliapkg.json
[juliapkg] Found dependencies: /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages/juliacall/juliapkg.json
[juliapkg] Found dependencies: /home/runner/work/OpenMDAO/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/juliapkg.json
[juliapkg] Found dependencies: /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/lib/python3.13/site-packages/juliapkg/juliapkg.json
[juliapkg] Locating Julia ^1.10.3
[juliapkg] Using Julia 1.12.4 at /usr/bin/julia
[juliapkg] Using Julia project at /home/runner/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/julia_env
[juliapkg] Writing Project.toml:
           | [deps]
           | OpenMDAOCore = "24d19c10-6eee-420f-95df-4537264b2753"
           | LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
           | PythonCall = "6099a3de-0909-46bc-b1f4-468b9a2dfc0d"
           | OpenSSL_jll = "458c3c95-2e84-50aa-8efc-19380b2a3a95"
           | MyParaboloidPackage = "a013c6e3-dc76-40b3-894d-62618f778e16"
           | 
           | [compat]
           | OpenMDAOCore = "^0.3.2"
           | PythonCall = "=0.9.31"
           | OpenSSL_jll = "3.0.0 - 3.6"
[juliapkg] Installing packages:
           | import Pkg
           | Pkg.Registry.update()
           | Pkg.develop([
           |   Pkg.PackageSpec(name="MyParaboloidPackage", uuid="a013c6e3-dc76-40b3-894d-62618f778e16", path=raw"/home/runner/work/OpenMDAO/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage"),
           | ])
           | Pkg.add([
           |   Pkg.PackageSpec(name="OpenMDAOCore", uuid="24d19c10-6eee-420f-95df-4537264b2753"),
           |   Pkg.PackageSpec(name="LinearAlgebra", uuid="37e2e46d-f89d-539d-b4ee-838fcccc9c8e"),
           |   Pkg.PackageSpec(name="PythonCall", uuid="6099a3de-0909-46bc-b1f4-468b9a2dfc0d"),
           |   Pkg.PackageSpec(name="OpenSSL_jll", uuid="458c3c95-2e84-50aa-8efc-19380b2a3a95"),
           | ])
           | Pkg.resolve()
           | Pkg.precompile()
  Installing known registries into `~/.julia`
       Added `General` registry to ~/.julia/registries
    Updating registry at `~/.julia/registries/General.toml`
   Resolving package versions...
   Installed IrrationalConstants ───────── v0.2.6
   Installed MicroMamba ────────────────── v0.1.14
   Installed SciMLPublic ───────────────── v1.0.1
   Installed Adapt ─────────────────────── v4.4.0
   Installed DiffRules ─────────────────── v1.15.1
   Installed Scratch ───────────────────── v1.3.0
   Installed Functors ──────────────────── v0.5.2
   Installed PythonCall ────────────────── v0.9.31
   Installed JSON3 ─────────────────────── v1.14.3
   Installed TableTraits ───────────────── v1.0.1
   Installed ADTypes ───────────────────── v1.21.0
   Installed Preferences ───────────────── v1.5.1
   Installed DiffResults ───────────────── v1.1.0
   Installed Parsers ───────────────────── v2.8.3
   Installed Tables ────────────────────── v1.12.1
   Installed SpecialFunctions ──────────── v2.7.1
   Installed Unitful ───────────────────── v1.28.0
   Installed IfElse ────────────────────── v0.1.1
   Installed DataAPI ───────────────────── v1.16.0
   Installed Pidfile ───────────────────── v1.3.0
   Installed micromamba_jll ────────────── v1.5.12+0
   Installed StaticArraysCore ──────────── v1.4.4
   Installed JLLWrappers ───────────────── v1.7.1
   Installed NaNMath ───────────────────── v1.1.3
   Installed IteratorInterfaceExtensions ─ v1.0.0
   Installed UnitfulAngles ─────────────── v0.7.2
   Installed PrecompileTools ───────────── v1.3.3
   Installed ConstructionBase ──────────── v1.6.0
   Installed DataValueInterfaces ───────── v1.0.0
   Installed OrderedCollections ────────── v1.8.1
   Installed ArrayInterface ────────────── v7.22.0
   Installed SparseMatrixColorings ─────── v0.4.23
   Installed CommonSubexpressions ──────── v0.3.1
   Installed DataStructures ────────────── v0.18.22
   Installed ForwardDiff ───────────────── v1.3.2
   Installed ChainRulesCore ────────────── v1.26.0
   Installed LogExpFunctions ───────────── v0.3.29
   Installed Requires ──────────────────── v1.3.1
   Installed UnsafePointers ────────────── v1.0.0
   Installed Compat ────────────────────── v4.18.1
   Installed Static ────────────────────── v1.3.1
   Installed MacroTools ────────────────── v0.5.16
   Installed OpenSpecFun_jll ───────────── v0.5.6+0
   Installed CommonWorldInvalidations ──── v1.0.0
   Installed StaticArrayInterface ──────── v1.9.0
   Installed pixi_jll ──────────────────── v0.41.3+0
   Installed StructTypes ───────────────── v1.11.0
   Installed DocStringExtensions ───────── v0.9.5
   Installed CondaPkg ──────────────────── v0.2.33
   Installed DifferentiationInterface ──── v0.7.16
   Installed ComponentArrays ───────────── v0.15.32
   Installed OpenMDAOCore ──────────────── v0.3.3
  Installing 1 artifacts
   Installed artifact OpenSpecFun 194.9 KiB
    Updating `~/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/julia_env/Project.toml`
  [a013c6e3] + MyParaboloidPackage v0.1.0 `~/work/OpenMDAO/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage`
  [24d19c10] + OpenMDAOCore v0.3.3
  [6099a3de] + PythonCall v0.9.31
  [37e2e46d] ~ LinearAlgebra ⇒ v1.12.0
  [458c3c95] ~ OpenSSL_jll ⇒ v3.5.4+0
    Updating `~/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/julia_env/Manifest.toml`
  [47edcb42] + ADTypes v1.21.0
  [79e6a3ab] + Adapt v4.4.0
  [4fba245c] + ArrayInterface v7.22.0
  [d360d2e6] + ChainRulesCore v1.26.0
  [bbf7d656] + CommonSubexpressions v0.3.1
  [f70d9fcc] + CommonWorldInvalidations v1.0.0
  [34da2185] + Compat v4.18.1
  [b0b7db55] + ComponentArrays v0.15.32
  [992eb4ea] + CondaPkg v0.2.33
  [187b0558] + ConstructionBase v1.6.0
  [9a962f9c] + DataAPI v1.16.0
 [864edb3b] + DataStructures v0.18.22
  [e2d170a0] + DataValueInterfaces v1.0.0
  [163ba53b] + DiffResults v1.1.0
  [b552c78f] + DiffRules v1.15.1
  [a0c0ee7d] + DifferentiationInterface v0.7.16
  [ffbed154] + DocStringExtensions v0.9.5
  [f6369f11] + ForwardDiff v1.3.2
  [d9f16b24] + Functors v0.5.2
  [615f187c] + IfElse v0.1.1
  [92d709cd] + IrrationalConstants v0.2.6
  [82899510] + IteratorInterfaceExtensions v1.0.0
  [692b3bcd] + JLLWrappers v1.7.1
  [0f8b85d8] + JSON3 v1.14.3
  [2ab3a3ac] + LogExpFunctions v0.3.29
  [1914dd2f] + MacroTools v0.5.16
  [0b3b1443] + MicroMamba v0.1.14
  [a013c6e3] + MyParaboloidPackage v0.1.0 `~/work/OpenMDAO/OpenMDAO/openmdao/docs/openmdao_book/features/experimental/MyParaboloidPackage`
  [77ba4419] + NaNMath v1.1.3
  [24d19c10] + OpenMDAOCore v0.3.3
  [bac558e1] + OrderedCollections v1.8.1
  [69de0a69] + Parsers v2.8.3
  [fa939f87] + Pidfile v1.3.0
  [aea7be01] + PrecompileTools v1.3.3
  [21216c6a] + Preferences v1.5.1
  [6099a3de] + PythonCall v0.9.31
  [ae029012] + Requires v1.3.1
  [431bcebd] + SciMLPublic v1.0.1
  [6c6a2e73] + Scratch v1.3.0
  [0a514795] + SparseMatrixColorings v0.4.23
  [276daf66] + SpecialFunctions v2.7.1
  [aedffcd0] + Static v1.3.1
  [0d7ed370] + StaticArrayInterface v1.9.0
  [1e83bf80] + StaticArraysCore v1.4.4
  [856f2bd8] + StructTypes v1.11.0
  [3783bdb8] + TableTraits v1.0.1
  [bd369af6] + Tables v1.12.1
  [1986cc42] + Unitful v1.28.0
  [6fb2a4bd] + UnitfulAngles v0.7.2
  [e17b2a0c] + UnsafePointers v1.0.0
  [efe28fd5] + OpenSpecFun_jll v0.5.6+0
 [f8abcde7] + micromamba_jll v1.5.12+0
  [4d7b5844] + pixi_jll v0.41.3+0
  [0dad84c5] + ArgTools v1.1.2
  [56f22d72] + Artifacts v1.11.0
  [2a0f44e3] + Base64 v1.11.0
  [ade2ca70] + Dates v1.11.0
  [f43a241f] + Downloads v1.7.0
  [7b1f6079] + FileWatching v1.11.0
  [b77e0a4c] + InteractiveUtils v1.11.0
  [ac6e5ff7] + JuliaSyntaxHighlighting v1.12.0
  [4af54fe1] + LazyArtifacts v1.11.0
  [b27032c2] + LibCURL v0.6.4
  [76f85450] + LibGit2 v1.11.0
  [8f399da3] + Libdl v1.11.0
  [37e2e46d] ~ LinearAlgebra ⇒ v1.12.0
  [56ddb016] + Logging v1.11.0
  [d6f4376e] + Markdown v1.11.0
  [a63ad114] + Mmap v1.11.0
  [ca575930] + NetworkOptions v1.3.0
  [44cfe95a] + Pkg v1.12.1
  [de0858da] + Printf v1.11.0
  [9a3f8284] + Random v1.11.0
  [ea8e919c] + SHA v0.7.0
  [9e88b42a] + Serialization v1.11.0
  [2f01184e] + SparseArrays v1.12.0
  [f489334b] + StyledStrings v1.11.0
  [fa267f1f] + TOML v1.0.3
  [a4e569a6] + Tar v1.10.0
  [8dfed614] + Test v1.11.0
  [cf7118a7] + UUIDs v1.11.0
  [4ec0a83e] + Unicode v1.11.0
  [e66e0078] + CompilerSupportLibraries_jll v1.3.0+1
  [deac9b47] + LibCURL_jll v8.15.0+0
  [e37daf67] + LibGit2_jll v1.9.0+0
  [29816b5a] + LibSSH2_jll v1.11.3+1
  [14a3606d] + MozillaCACerts_jll v2025.11.4
  [4536629a] + OpenBLAS_jll v0.3.29+0
  [05823500] + OpenLibm_jll v0.8.7+0
  [458c3c95] ~ OpenSSL_jll ⇒ v3.5.4+0
  [bea87d4a] + SuiteSparse_jll v7.8.3+2
  [83775a58] + Zlib_jll v1.3.1+2
  [8e850b90] + libblastrampoline_jll v5.15.0+0
  [8e850ede] + nghttp2_jll v1.64.0+1
  [3f19e933] + p7zip_jll v17.7.0+0
        Info Packages marked with  have new versions available but compatibility constraints restrict them from upgrading. To see why use `status --outdated -m`
   Resolving package versions...
     Project No packages added to or removed from `~/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/julia_env/Project.toml`
    Manifest No packages added to or removed from `~/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/julia_env/Manifest.toml`
Precompiling
 packages...
   1529.4 msIfElse
   1532.1 msDataValueInterfaces
   1530.8 msDocStringExtensions
    725.5 msUnsafePointers
    718.4 msIteratorInterfaceExtensions
   2684.3 msOrderedCollections
    739.8 msDataAPI
    794.2 msCommonWorldInvalidations
    784.0 msSciMLPublic
   4096.1 msIrrationalConstants
   1209.1 msADTypes
   3048.5 msConstructionBase
   1133.2 msStaticArraysCore
   1729.3 msRequires
   2721.7 msNaNMath
    898.1 msScratch
   1243.2 msPreferences
   3984.0 msStructTypes
   1348.7 msPidfile
    931.6 msTableTraits
   6961.0 msMacroTools
   3233.3 msCompat
    810.1 msConstructionBase → ConstructionBaseLinearAlgebraExt
   1108.0 msLogExpFunctions
    838.7 msADTypes → ADTypesConstructionBaseExt
    857.0 msDiffResults
   1170.3 msAdapt
    871.0 msPrecompileTools
    964.5 msJLLWrappers
   1045.1 msCommonSubexpressions
    973.1 msCompat → CompatLinearAlgebraExt
   3976.4 msDifferentiationInterface
   1471.0 msAdapt → AdaptSparseArraysExt
   3578.3 msTables
   2669.5 msArrayInterface
   2005.8 msmicromamba_jll
   2146.8 mspixi_jll
   1450.4 msOpenSpecFun_jll
   9424.3 msStatic
   3276.5 msDataStructures
   3990.2 msFunctors
   5224.3 msChainRulesCore
   1406.9 msDifferentiationInterface → DifferentiationInterfaceSparseArraysExt
   1573.7 msArrayInterface → ArrayInterfaceSparseArraysExt
    910.3 msArrayInterface → ArrayInterfaceStaticArraysCoreExt
  15894.8 msSparseMatrixColorings
   4087.3 msMicroMamba
   5630.8 msSpecialFunctions
   1002.2 msADTypes → ADTypesChainRulesCoreExt
   3678.8 msChainRulesCore → ChainRulesCoreSparseArraysExt
   1165.5 msDifferentiationInterface → DifferentiationInterfaceChainRulesCoreExt
   2924.5 msLogExpFunctions → LogExpFunctionsChainRulesCoreExt
   1275.0 msArrayInterface → ArrayInterfaceChainRulesCoreExt
  27906.0 msParsers
   1684.1 msDifferentiationInterface → DifferentiationInterfaceSparseMatrixColoringsExt
   1567.8 msDiffRules
   3643.9 msSpecialFunctions → SpecialFunctionsChainRulesCoreExt
   8434.0 msForwardDiff
  46952.7 msUnitful
   1200.0 msUnitful → ConstructionBaseUnitfulExt
   1269.0 msUnitful → PrintfExt
   1042.3 msUnitful → NaNMathExt
   1232.3 msUnitful → ForwardDiffExt
   3393.9 msDifferentiationInterface → DifferentiationInterfaceForwardDiffExt
   2221.3 msUnitfulAngles
  15073.7 msJSON3
  35035.8 msStaticArrayInterface
   1011.2 msComponentArrays
   1829.6 msOpenMDAOCore
   1373.6 msMyParaboloidPackage
  21877.7 msCondaPkg
  20622.9 msPythonCall
  72 dependencies successfully precompiled in 101 seconds. 34 already precompiled.
  2 dependencies had output during precompilation:
MicroMamba
 Downloading artifact: micromamba

CondaPkg
 Downloading artifact: pixi

     Project No packages added to or removed from `~/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/julia_env/Project.toml`
    Manifest No packages added to or removed from `~/work/OpenMDAO/OpenMDAO/.pixi/envs/dev/julia_env/Manifest.toml`
Detected IPython. Loading juliacall extension. See https://juliapy.github.io/PythonCall.jl/stable/compat/#IPython

Now we have a OpenMDAO Component that we can use. Let’s try it out, following the Paraboloid example from the OpenMDAO docs:

import openmdao.api as om

model = om.Group()
model.add_subsystem('parab_comp', comp)

prob = om.Problem(model)

prob.model.add_design_var("parab_comp.x")
prob.model.add_design_var("parab_comp.y")
prob.model.add_objective("parab_comp.f_xy")

prob.driver = om.ScipyOptimizeDriver()
prob.driver.options["optimizer"] = "SLSQP"

prob.setup(force_alloc_complex=True)

prob.set_val("parab_comp.x", 3.0)
prob.set_val("parab_comp.y", -4.0)

prob.run_model()
print(prob["parab_comp.f_xy"])  # Should print `[-15.]`

prob.set_val("parab_comp.x", 5.0)
prob.set_val("parab_comp.y", -2.0)

prob.run_model()
print(prob.get_val("parab_comp.f_xy"))  # Should print `[-5.]`
-15.0
-5.0

Looks good so far. A few things to note:

  • We added the JuliaExplicitComp to the Group in the normal way.

  • The variable names exposed to OpenMDAO match what was used in the input and output ComponentVectors passed to the DenseADExplicitComp constructor.

  • prob.set_val and prob.get_val worked just like usual.

We can also check our derivatives using the usual OpenMDAO methods. We can do that with the finite difference method:

print(prob.check_partials(method="fd"))
------------------------------------------------------------------------------
Component: JuliaExplicitComp 'parab_comp'
------------------------------------------------------------------------------

  parab_comp: 'f_xy' wrt 'x'

    Max Tolerance Violation (Jfwd - Jfd) - (atol + rtol * Jfd) : (-9.996325e-07)
      abs error: 1.000368e-06
      rel error: 5.001840e-07
      fwd value @ max viol: 2.000000e+00
      fd value @ max viol: 2.000001e+00 (fd:forward)

    Raw Forward Derivative (Jfwd)
    [[ 2.000000000000e+00]]

    Raw FD Derivative (Jfd)
    [[ 2.000001000368e+00]]

 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

  parab_comp: 'f_xy' wrt 'y'

    Max Tolerance Violation (Jfwd - Jfd) - (atol + rtol * Jfd) : (-7.999542e-06)
      abs error: 1.000459e-06
      rel error: 1.111621e-07
      fwd value @ max viol: 9.000000e+00
      fd value @ max viol: 9.000001e+00 (fd:forward)

    Raw Forward Derivative (Jfwd)
    [[ 9.000000000000e+00]]

    Raw FD Derivative (Jfd)
    [[ 9.000001000459e+00]]

 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

{'parab_comp': {('f_xy', 'x'): {'J_fwd': array([[2.]]), 'J_fd': array([[2.000001]]), 'rows': None, 'cols': None, 'tol violation': _ErrorData(forward=-9.99632543852158e-07, reverse=None, fwd_rev=None), 'magnitude': _MagnitudeData(forward=2.0, reverse=0.0, fd=2.0000010003684565), 'vals_at_max_error': _ErrorData(forward=(np.float64(2.0), np.float64(2.0000010003684565)), reverse=None, fwd_rev=None), 'abs error': _ErrorData(forward=1.0003684565162985e-06, reverse=None, fwd_rev=None), 'rel error': _ErrorData(forward=5.001839780740122e-07, reverse=None, fwd_rev=None)}, ('f_xy', 'y'): {'J_fwd': array([[9.]]), 'J_fd': array([[9.000001]]), 'rows': None, 'cols': None, 'tol violation': _ErrorData(forward=-7.999542276593274e-06, reverse=None, fwd_rev=None), 'magnitude': _MagnitudeData(forward=9.0, reverse=0.0, fd=9.000001000458724), 'vals_at_max_error': _ErrorData(forward=(np.float64(9.0), np.float64(9.000001000458724)), reverse=None, fwd_rev=None), 'abs error': _ErrorData(forward=1.0004587238654494e-06, reverse=None, fwd_rev=None), 'rel error': _ErrorData(forward=1.1116206807248763e-07, reverse=None, fwd_rev=None)}}}

Or we could use the complex-step method:

print(prob.check_partials(method="cs"))
------------------------------------------------------------------------------
Component: JuliaExplicitComp 'parab_comp'
------------------------------------------------------------------------------

  parab_comp: 'f_xy' wrt 'x'

    Max Tolerance Violation (Jfwd - Jfd) - (atol + rtol * Jfd) : (-2.000000e-06)
      abs error: 0.000000e+00
      rel error: 0.000000e+00
      fwd value @ max viol: 2.000000e+00
      fd value @ max viol: 2.000000e+00 (cs:None)

    Raw Forward Derivative (Jfwd)
    [[ 2.000000000000e+00]]

    Raw CS Derivative (Jfd)
    [[ 2.000000000000e+00]]

 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

  parab_comp: 'f_xy' wrt 'y'

    Max Tolerance Violation (Jfwd - Jfd) - (atol + rtol * Jfd) : (-9.000000e-06)
      abs error: 0.000000e+00
      rel error: 0.000000e+00
      fwd value @ max viol: 9.000000e+00
      fd value @ max viol: 9.000000e+00 (cs:None)

    Raw Forward Derivative (Jfwd)
    [[ 9.000000000000e+00]]

    Raw CS Derivative (Jfd)
    [[ 9.000000000000e+00]]

 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

{'parab_comp': {('f_xy', 'x'): {'J_fwd': array([[2.]]), 'J_fd': array([[2.]]), 'rows': None, 'cols': None, 'tol violation': _ErrorData(forward=-2e-06, reverse=None, fwd_rev=None), 'magnitude': _MagnitudeData(forward=2.0, reverse=0.0, fd=2.0), 'vals_at_max_error': _ErrorData(forward=(np.float64(2.0), np.float64(2.0)), reverse=None, fwd_rev=None), 'abs error': _ErrorData(forward=0.0, reverse=None, fwd_rev=None), 'rel error': _ErrorData(forward=0.0, reverse=None, fwd_rev=None)}, ('f_xy', 'y'): {'J_fwd': array([[9.]]), 'J_fd': array([[9.]]), 'rows': None, 'cols': None, 'tol violation': _ErrorData(forward=-9e-06, reverse=None, fwd_rev=None), 'magnitude': _MagnitudeData(forward=9.0, reverse=0.0, fd=9.0), 'vals_at_max_error': _ErrorData(forward=(np.float64(9.0), np.float64(9.0)), reverse=None, fwd_rev=None), 'abs error': _ErrorData(forward=0.0, reverse=None, fwd_rev=None), 'rel error': _ErrorData(forward=0.0, reverse=None, fwd_rev=None)}}}

The derivatives look great, so let’s try an optimization:

prob.run_driver()
print(f"f_xy = {prob.get_val('parab_comp.f_xy')}")  # Should print `[-27.333333]`
print(f"x = {prob.get_val('parab_comp.x')}")  # Should print `[6.666666]`
print(f"y = {prob.get_val('parab_comp.y')}")  # Should print `[-7.333333]`
Optimization terminated successfully    (Exit mode 0)
            Current function value: -27.333333333333336
            Iterations: 4
            Function evaluations: 5
            Gradient evaluations: 4
Optimization Complete
-----------------------------------
f_xy = -27.333333333333336
x = 6.666666666666666
y = -7.333333333333333

Victory!

Next Steps#

We’ve shown how to put together a simple OpenMDAO Component using Julia and OpenMDAO.jl. OpenMDAO.jl has more features that you can explore in docs. Some things we haven’t covered: