JAXVacua: overview#

What’s in this notebook? This notebook provides an overview of JAXVacua’s core capabilities.

  • Quick Start — a self-contained end-to-end example: load a model, find a SUSY flux vacuum with Newton’s method, extract \(W_0\), \(g_s\), and the mass spectrum.

  • periods — period vector, gauge kinetic matrix, ISD matrix; LCS prepotential; worldsheet instanton corrections and GV invariants.

  • css — Kähler potential and metric in affine coordinates; gauge choices; visualisation of the Kähler potential landscape.

  • FluxEFT — GVW superpotential, F-terms, scalar potential, Hessian; visualisation of the scalar potential landscape.

  • sampling — uniform and ISD sampling of flux/moduli configurations.

(Created: Andreas Schachner, June 2024)

Outline#

  1. Setup

  2. Quick Start

  3. Introduction to JAXVacua

  4. Sub-modules and basic capabilities

  5. Take-aways

  6. Further reading

Setup#

# General imports
import warnings
import numpy as np
from tqdm.auto import tqdm

# JAX imports
import jax
import jax.numpy as jnp
jax.config.update("jax_enable_x64", True)

# Plotting
import seaborn as sn
import matplotlib.pyplot as plt
cmap = sn.color_palette("viridis", as_cmap=True)

# JAXVacua
import jaxvacua as jvc

warnings.filterwarnings('ignore')

How to read this overview#

JAXVacua is not only a Newton solver. It is a stack of composable tools for flux compactifications:

Layer

What JAXVacua provides

Where to go next

Geometry and periods

Calabi-Yau topological data, period vectors, Kähler potentials, gauge-kinetic and ISD matrices

NB03, NB09, API pages for lcs, periods, css

Flux effective theory

Superpotential, F-terms, scalar potential, Hessians, masses and physical observables

NB05, NB14, jaxvacua.flux_eft

Sampling and search

Initial-condition samplers, ISD-inspired flux completion, scipy and internal Newton refinement

NB04, NB06, NB07, jaxvacua.sampling, jaxvacua.flux_vacua_finder

Bound and enumerate

Finite-region flux bounding, enumeration, stochastic bounded sampling and cluster helpers

NB08, jaxvacua.flux_bounding

Special limits and reductions

coni-LCS prepotentials, PFV workflows, conifold freezing and reduced EFTs

NB09, NB10, NB12

The rest of this notebook gives a compact tour through these layers. It is intended as orientation; the neighbouring tutorials provide the detailed workflows.

Quick Start#

A complete end-to-end workflow: load a model → find a SUSY flux vacuum → extract physical observables. The detailed sections below explain each step.

Step

What happens

Load FluxVacuaFinder

Reads topological data, builds the LCS period vector, and exposes the Newton solver

Newton refinement

Solves F-term conditions \(D_I W = 0\) via Newton’s method starting from a nearby point

Physical observables

Superpotential \(W_0\), string coupling \(g_s = 1/\mathrm{Im}(\tau)\), D3 tadpole \(N_\mathrm{flux}\), mass spectrum

print(f"Default: precision={jvc.precision}, jnp.array(1.0).dtype={jnp.array(1.0).dtype}")

jvc.set_precision("float32")
print(f"float32: precision={jvc.precision}, jnp.array(1.0).dtype={jnp.array(1.0).dtype}")

jvc.set_precision("float64")
print(f"float64: precision={jvc.precision}, jnp.array(1.0).dtype={jnp.array(1.0).dtype}")
# 2. Load model and define a starting point near a known vacuum
# h12=2 Kreuzer-Skarke model at Large Complex Structure, with worldsheet instanton corrections
# FluxVacuaFinder extends FluxEFT with Newton-refinement helpers.
model = jvc.FluxVacuaFinder(h12=2, model_ID=1, maximum_degree=2)


# Flux configuration from arXiv:1912.10047 (nearby starting point, not exact vacuum)
fluxes  = jnp.array([7, 3, -24, 0, -16, 50, 0, 3, -4, 0, 0, 0])
z_start = jnp.array([2.7j, 2.0j])
tau_start = 6.9j
# 3. Refine to exact SUSY vacuum (D_I W = 0) using Newton's method
z_sol, tau_sol, res = model.newton_method_flux_vacua(
    z_start, tau_start, fluxes, mode="SUSY", step_size_Newton=0.1, tol=1e-12)

print(f"Converged moduli:   z = {z_sol}")
print(f"Axio-dilaton:       \u03c4 = {tau_sol:.8f}")
print(f"F-term residual:    max|D_I W| = {res:.2e}")
# 4. Extract physical observables at the SUSY vacuum
W0     = model.superpotential(z_sol, tau_sol, fluxes)
g_s    = float(1.0 / jnp.imag(tau_sol))
N_flux = int(model.tadpole(fluxes))
H      = model.hessian(z_sol, jnp.conj(z_sol), tau_sol, jnp.conj(tau_sol), fluxes, mode="SUSY")
masses = jnp.sort(jnp.real(jnp.linalg.eigvalsh(H)))

print(f"Superpotential:              |W\u2080|= {np.abs(W0):.8f}")
print(f"String coupling:             g_s = 1/Im(\u03c4) = {g_s:.4f}")
print(f"D3-tadpole:                  N_flux = {N_flux}")
print(f"Hessian eigenvalues (\u221d m\u00b2): {masses}")

Introduction to JAXVacua#

JAXVacua is a software package designed for the efficient exploration of string theory vacua, leveraging the computational power of JAX. Developed as a tool for high-dimensional optimisation and machine learning applications in theoretical physics, JAXVacua provides a scalable and flexible framework for analysing complex energy landscapes in string theory and related fields.

Key features of JAXVacua include:

  • Automatic Differentiation: Utilising JAX’s autograd capabilities to compute gradients and Hessians efficiently, aiding in the search for vacua.

  • Just-In-Time (JIT) Compilation: Enhancing performance through XLA-based compilation for CPU, GPU, and TPU acceleration.

  • Vectorisation (vmap): Enabling simultaneous evaluation of multiple candidate vacua, significantly improving computational efficiency.

  • Parallelism (pmap): Allowing large-scale parallel computations across multiple devices, making it suitable for extensive searches in high-dimensional parameter spaces.

JAXVacua is particularly well-suited for studies in flux compactifications, moduli stabilisation, and other problems requiring rapid and scalable numerical analysis. By integrating modern machine learning techniques with traditional methods in string theory, it provides a powerful computational framework for advancing research in fundamental physics.

To reiterate, moduli stabilisation requires finding vacua where the scalar potential, derived from fluxes, branes, and other ingredients in string theory, has local minima. As we argue below, automatic differentiation can play a crucial role in this process by efficiently computing gradients and Hessians of the scalar potential. The first derivatives (gradients) help identify critical points where the potential is stationary, while second derivatives (Hessians) determine their stability by checking for positive eigenvalues (indicating local minima). Given the high-dimensional and complex nature of the moduli space, manual differentiation is infeasible, and numerical methods like finite differences introduce errors. AD provides exact derivatives up to machine precision without symbolic computation overhead, making it highly efficient.

Frameworks like JAX further enhance this process by enabling rapid gradient computation (jax.grad), second-order differentiation (jax.hessian), and Just-In-Time (JIT) compilation for acceleration. These features allow for efficient optimisation and exploration of the moduli space, making AD a powerful tool for string compactifications and moduli stabilisation. We will introduce this process further below.

Sub-modules and basic capabilities#

Let us describe one example in some more detail. First, we provide different options to load examples:

  • From file: we provide files for selected models from the Kreuzer-Skarke (KS) list for varying \(h^{1,2}\) as well as for all Complete Intersection CY threefolds (CICYs). These models can be accessed by specifying \(h^{1,2}\) and a corresponding model_ID.

  • From input: if either periods or other relevant information specifying the periods are known, they can be provided directly. E.g., at LCS, the intersection numbers of the mirror CY as well as some other topological information suffices to completely specify the periods. See NB09: Moduli-space limits.

  • From CYTools: for KS models at LCS, we can simply provide a cytools.CalabiYau object to compute the relevant topological data of the mirror CY. See NB03: CYTools Interface.

  • From the database: models can be loaded directly from StringForge’s TDFDatabase or CICYDatabase using a HuggingFace-hosted parquet database. See the StringForge documentation.

The package consists of the following modules:

Module

Class / entry point

Description

jaxvacua.periods

periods

Period vector, prepotential, gauge kinetic matrix, ISD matrix

jaxvacua.css

css

Kähler potential, metric, F-terms in affine coordinates

jaxvacua.flux_eft

FluxEFT

GVW superpotential, scalar potential, Hessian, and masses

jaxvacua.flux_vacua_finder

FluxVacuaFinder

Newton refinement and higher-level vacuum-search workflows

jaxvacua.sampling

data_sampler

Moduli/flux sampling, ISD sampling

jaxvacua.flux_bounding

bounded_fluxes

Bounding boxes and enumeration for integer flux search

jaxvacua.flux_eft

flux_eft

EFT validity checks, instanton convergence

jaxvacua.flux_utils

utility functions

Standalone PFV algebra (split, convert fluxes)

jaxvacua.freezer

ConifoldFreezer

Freeze out heavy (conifold) moduli in a reduced EFT

jaxvacua.lcs

lcs_tree

Topological data container for LCS models

jaxvacua.hypergeometric_models

HypergeometricModels

Closed-form one-modulus prepotentials at LCS / K-point / C-point

jaxvacua.cytools_interface

utility functions

Extract topological data from CYTools objects

Module composition: periodscssFluxEFTsampling / flux_bounding / flux_eft / freezer

These classes are useful to compute periods for a large fraction of models; see in particular NB03: CYTools Interface, and the CICY tutorial in the StringForge documentation. Beyond periods at LCS, our implementation allows users to easily include their own custom period class for different regions in (complex structure) moduli space — see NB09: Moduli-space limits.

periods sub-module#

This sub-module contains periods class implementing standard formulas in terms of the periods.

The holomorphic \(3\)-form \(\Omega\) on \(X\) can be parametrised by the real, symplectic basis of \(3\)-forms \((\alpha_I, \beta^J)\in H_{-}^{3}(X_{3},\mathbb{Z})\), \(I,J\in \{0, \ldots, h^{1,2}_-(X)\}\), together with a dual basis of \(3\)-cycles \((A_{I},B^{J})\in H^{-}_{3}(X_{3},\mathbb{Z})\) so that

\[ \Omega_{3}= X^I \, \alpha_I - \, \mathcal{F}_{J} \, \beta^J\; ,\quad X^{I}= \int_{A_{I}}\, \Omega_{3}\; ,\quad F_{J}= \int_{B^{J}}\, \Omega_{3} \, . \]

Here, \((X^{I},\mathcal{F}_{J})\) are the so-called \emph{periods} of \(\Omega_{3}\) which are usually arranged in the \emph{period vector}\index{Period vector}

\[\begin{split} \Pi=\left (\begin{array}{c}\mathcal{F}_{J}\\[0.3em] X^{I}\end{array} \right )\, . \end{split}\]

Due to the special Kähler structure of \(\mathcal{M}_{\text{cs}}(X)\), half of the periods can be obtained from a prepotential \(F(X^I)\), that is,

\[\begin{split} \Pi=\left (\begin{array}{c} \partial_{X^I}F\\ X^{I} \end{array} \right )\, . \end{split}\]

Another readily available function that can be computed from knowing the periods is the gauge kinetic matrix which is defined as

\[ \mathcal{N}_{IJ} = (\mathcal{F}_{I} ,\, D_{\bar{\imath}}\overline{\mathcal{F}}_I )\, \cdot\, (X^{J},\, D_{\bar{\jmath}}\overline{X}^J)^{-1}= P_{IK} \, (Q^{-1})^{K}{}_{J}\quad , \; P_{IJ} = (\mathcal{F}_{I} ,\, D_{\bar{\imath}}\overline{\mathcal{F}}_I )_{J}\quad , \; Q^{I}{}_{J} = (X^{I},\, D_{\bar{\jmath}}\overline{X}^I)_{J}\, , \]
see e.g. Sect. 2 above Eq. (2.12) in 2110.05511 or Eq. (3.5) in 2310.06040. If we are working with a prepotential \(F\), the gauge kinetic matrix can also be computed using
\[ \mathcal{N}_{I J}=\overline{F}_{IJ}+2\text{i}\, \dfrac{\text{Im}(F_{I L})X^{L} \, \text{Im}(F_{J K})X^{K}}{X^{M}\text{Im}(F_{MN})X^{N}}\; ,\quad F_{IJ}= \partial_{X^{I}} \partial_{X^{J}}F\, . \]
Both options are implemented in our code, see periods.gauge_kinetic_matrix_periods and periods.gauge_kinetic_matrix_prepotential.

By writing \(\mathcal{N}=\mathcal{R}+\text{i} \mathcal{I}\), the Hodge-star matrix \(\mathcal{M}\) is given by

\[\begin{split} \mathcal{M}=\left (\begin{array}{cc} -\mathcal{I}^{-1} \hspace{0.3cm}&\mathcal{I}^{-1}\mathcal{R} \\ \mathcal{R}\mathcal{I}^{-1} & -\mathcal{I}-\mathcal{R}\mathcal{I}^{-1}\mathcal{R} \end{array} \right )\, . \end{split}\]
This function can be called via periods.ISD_matrix. Note that our convention differs slightly from that of Eq. (3.7) in 2310.06040 due to the ordering of periods in the period vector. Specifically, we have
\[ \mathcal{M}_{\mathrm{here}} = \mathcal{M}_{\mathrm{there}}^{-1}\, . \]

Large Complex Structure limits#

At Large Complex Structure (LCS), the prepotential can be expressed as

\[ F_{\mathrm{LCS}}(X) = F_{\mathrm{poly}}(X) + F_{\mathrm{inst}}(X)\, . \]

Here, the polynomial part \(F_{\mathrm{poly}}\) of the LCS prepotential \(F_{\mathrm{LCS}}\) can be expressed in terms of the periods \(X^I=(X^0,X^i)\) as

\[ F_{\mathrm{poly}}(X) = -\frac{1}{6X^0}\widetilde{\kappa}_{ijk}X^iX^jX^k+\frac{1}{2}a_{ij}X^iX^j +b_{i}X^i\, X^0 + \dfrac{\text{i}}{2}\tilde{\xi}\, (X^0)^2\, , \]

Here, \(\widetilde{\kappa}_{ijk}\) are the triple intersection numbers of the mirror dual Calabi-Yau threefold \(\widetilde{X}\). Here, we defined

\[\begin{split} a_{ij} = \dfrac{1}{2}\begin{cases} \widetilde{\kappa}_{iij} & i\geq j\\[0.3em] \widetilde{\kappa}_{ijj} & i<j \end{cases} \, , \quad b_i = \dfrac{1}{24} \int_{\tilde{D}^i}\, c_2(\widetilde{X})\, , \quad \tilde{\xi}=\frac{\zeta(3)\, \chi(\widetilde{X})}{(2\pi)^3}\, . \end{split}\]

The instanton part \(F_{\mathrm{inst}}\) of the LCS prepotential \(F_{\mathrm{LCS}}\) can be expressed in terms of the periods \(X^I=(X^0,X^i)\) as

\[ F_{\mathrm{inst}}(X) = -\frac{(X^0)^2}{(2\pi\mathrm{i})^3}\, \sum_{q\in\mathcal{M}(\widetilde{X})}\, n_q^{0}\, \text{Li}_3\left (\text{e}^{2\pi \text{i} q_i X^i / X^0}\right )\; , \quad \text{Li}_3\left (x\right )=\sum_{m=1}^{\infty}\, \dfrac{x^{m}}{m^{3}}\, . \]

Here the sum is performed over all effective curve classes \(q\in\mathcal{M}(\widetilde{X})\) in the Mori cone \(\mathcal{M}(\widetilde{X})\) of the mirror dual manifold \(\widetilde{X}\). Here, the \(n_q^{0}\) are the genus-0 Gopakumar-Vafa (GV) invariants which can be computed systematically using methods described in hep-th/9308122.

The infinite sum appearing in the poly-logarithm \(\text{Li}_3\) can be rewritten to arrive at

\[ \sum_{q\in\mathcal{M}(\widetilde{X})}\, n_q^{0}\, \text{Li}_3\left (\text{e}^{2\pi \text{i} q_i X^i / X^0}\right ) = \sum_{q\in\mathcal{M}(\widetilde{X})}\, N_q\, \text{e}^{2\pi \text{i} q_i X^i / X^0} \]

in terms of genus-0 Gromov-Witten (GW) invariants \(N_q\). We typically work with the latter to simplify the calculation.

The above functions are implemented in the periods class. When initialised with limit="LCS", the period computation uses the LCS prepotential constructed from the topological data stored in the lcs_tree.

Choosing limit and maximum_degree#

The periods and FluxEFT constructors accept a limit argument that selects the region of moduli space:

limit

When to use

What changes

"LCS"

Near the Large Complex Structure point; topological data (intersection numbers, GV invariants) known

Uses the polynomial + instanton LCS prepotential \(F_\mathrm{LCS} = F_\mathrm{poly} + F_\mathrm{inst}\)

"coniLCS"

Near a conifold locus within the LCS region

Applies a basis change to align the conifold charge direction; see NB09 and NB10

None / custom

User-supplied prepotential or period vector

Accepts prepotential_input or period_input callables; see NB09

The maximum_degree parameter controls how many Gromov-Witten terms are included in the instanton sum \(F_\mathrm{inst}\). Setting maximum_degree=0 uses only the polynomial part \(F_\mathrm{poly}\) (faster); increasing it adds higher-degree worldsheet corrections.

Example for periods class#

As an example, let us collect a simple model at LCS with \(h^{1,2}=2\), model_ID=1 and model_type="KS":

periods = jvc.periods(h12=2,model_ID=1,model_type="KS",limit="LCS")
periods

This example has mirror triple intersection numbers

periods.lcs_tree.intnums

It turns out that this is precisely the model discussed e.g. in https://inspirehep.net/literature/1772253.

Let us perform a small computation here by verifying that the ISD matrix \(\mathcal{M}\) representing the Hodge star \(\star\) as a matrix for a given integral basis is indeed symplectic. First, we evaluate the ISD matrix \(\mathcal{M}\) for the projective coordinates \(X = (1, i, i)\):

X = 1j*jnp.ones(periods.h12+1)
X = X.at[0].set(1.)
M = periods.ISD_matrix(X,jnp.conj(X))
M

As expected, the output is real as it should be. Now, we verify the following three conditions:

\[ \mathcal{M} = \mathcal{M}^T\; ,\quad \mathcal{M}^T\Sigma\mathcal{M} = \Sigma\; ,\quad \mathcal{M}^{-1}=\Sigma^T\mathcal{M}\Sigma\, . \]

Note that the second and third condition are not independent, but we test both regardless. These three conditions can be written in jax.numpy as:

cond1 = M - M.T
cond2 = jnp.matmul(M.T,jnp.matmul(periods.sigma,M))- periods.sigma
cond3 = jnp.linalg.inv(M)-jnp.matmul(periods.sigma.T,jnp.matmul(M,periods.sigma))

The tests then proceed as follows:

tol = 1e-10
test1 = bool(jnp.max(jnp.abs(cond1)) < tol)
test2 = bool(jnp.max(jnp.abs(cond2)) < tol)
test3 = bool(jnp.max(jnp.abs(cond3)) < tol)
test1, test2, test3

We see that all three conditions are indeed satisfied. Below, we will use the ISD matrix to test that the ISD condition for Type IIB GKP-type vacua is indeed satisfied.

A full list of methods and attributes can be found in the accompanying documentation.

GV invariants and worldsheet instanton corrections#

At LCS, the prepotential splits as \(F_\mathrm{LCS} = F_\mathrm{poly} + F_\mathrm{inst}\). The instanton part \(F_\mathrm{inst}\) encodes genus-0 Gopakumar-Vafa (GV) invariants through the poly-logarithm sum

\[ F_\mathrm{inst}(z) = -\frac{1}{(2\pi\mathrm{i})^3} \sum_{q \in \mathcal{M}(\widetilde{X})} n_q^0\, \mathrm{Li}_3\!\left(e^{2\pi\mathrm{i}\, q_i z^i}\right). \]
Deep in the LCS region (\(\mathrm{Im}(z^i) \gg 1\)), \(F_\mathrm{inst}\) is exponentially suppressed relative to \(F_\mathrm{poly}\). The maximum_degree parameter truncates the sum over curve classes \(q\) at a maximum degree \(\sum_i q_i \leq \) maximum_degree.

# Compare polynomial and instanton parts of the prepotential at the known vacuum
# css inherits F_LCS_poly and F_inst from the periods LCS implementation
css_inst = jvc.css(h12=2, model_ID=1, model_type="KS", maximum_degree=3)

z_vac = jnp.array([2.74215j, 2.05662j])   # from arXiv:1912.10047
F_poly = css_inst.F_LCS_poly(z_vac)
F_inst = css_inst.F_inst(z_vac)

print(f"|F_poly| = {abs(F_poly):.6f}")
print(f"|F_inst| = {abs(F_inst):.6f}")
print(f"|F_inst / F_poly| = {abs(F_inst / F_poly):.2e}   (instanton corrections are suppressed at LCS)")

css sub-module (complex structure sector)#

This sub-module contains css class containing functions to handle the Kähler geometry of the complex structure moduli space.

A parametrisation of the complex structure moduli space \(\mathcal{M}_{\text{cs}}(X)\) is conveniently obtained by choosing half of the periods, say \(X^{I}\) from above, as projective coordinates. Away from the locus \(X^0=0\), one then defines local affine coordinates

\[ z^{i}=\dfrac{X^{i}}{X^{0}}\; ,\quad i=1,\ldots ,h^{1,2}_{-}\, . \]

Typically, we normalise the 3-form \(\Omega\) such that \(X^0=1\), but our implementation allows for arbitrary choices for the normalisation. From the period class above, we obtain functions depending only on the complex structure moduli \(z^i\) by introducing suitable coordinates. The period vector in terms of these coordinates is again given in terms of the prepotential \(F(z)\) through (setting \(X^{0}=1\))

\[\begin{split} \Pi=\left (\begin{array}{c} 2F-z^{i} \partial_{i}F\\ \partial_{i}F\\ 1 \\ z^{i} \end{array} \right )\, . \end{split}\]

From this, we can define e.g. the Kähler potential

\[ K = -\log\biggl (-\mathrm{i}\Pi^{\dagger}\Sigma\Pi\biggl )-\log(-\mathrm{i}(\tau-\bar{\tau}) \]

as well as all of its first and second derivatives \(\partial_{z^i}K,\partial_{\tau}K,\partial_{z^i}\partial_{z^j}K,\partial_{\tau}\partial_{\tau}K\) etc. such as the Kähler metric

\[ K_{i\bar{\jmath}} = \partial_{z^i}\partial_{\overline{z}^j}K\; ,\quad K_{\tau\bar{\tau}} = \partial_{\tau}\partial_{\overline{\tau}}K\, . \]

The css class inherits some functions from the parent class jaxvacua.periods. E.g., the function jaxvacua.periods.ISD_matrix (which takes as input the periods and their complex conjugates) decends to jaxvacua.css.ISD_matrix which now requires as input \(z^{i}\) and \(\bar{z}^{i}\). Other functions inherited in this way are the derivatives of the ISD matrix, the gauge kinetic matrix \(\mathcal{N}\) (see above) and its derivatives.

Choice of gauge#

Let us look at the example from above:

css = jvc.css(h12=2,model_ID = 1,model_type="KS")

Per default, the gauge fixing \(X^0 = 1\)

z = 1j*jnp.ones(css.h12)
M = css.ISD_matrix(z,jnp.conj(z))
# Alias
#M = css.M(z,jnp.conj(z))
M

This should match

X = 1j*jnp.ones(periods.h12+1)
X = X.at[0].set(1.)
M_periods = periods.ISD_matrix(X,jnp.conj(X))
M-M_periods

We can also try a different gauge

# Define different gauge
X0 = 2+1j

css_X0 = jvc.css(h12=2,model_ID = 1,model_type="KS",gauge_choice=X0)

z = 1j*jnp.ones(css.h12)
M = css_X0.ISD_matrix(z,jnp.conj(z))

X = 1j*jnp.ones(periods.h12+1)*X0
X = X.at[0].set(X0)
M_periods = periods.ISD_matrix(X,jnp.conj(X))
M-M_periods

We find agreement again!

Computing the Kähler metric#

To compute the Käher metric, we use

# Choose point
z = 1j*jnp.ones(css.h12)
tau = 1j

# Evaluate metric
km = css.kahler_metric(z,jnp.conj(z),tau,jnp.conj(tau))
km

The inverse metric is computed as follows:

ikm = css.inverse_kahler_metric(z,jnp.conj(z),tau,jnp.conj(tau))
ikm

We can check that it indeed matches the expected form by computing the inverse of the Kähler metric km from above:

ikm-jnp.linalg.inv(km)

Kähler potential landscape#

To build geometric intuition, we visualise how the Kähler potential \(K\) varies as a function of the modulus \(z_1\) along the purely imaginary direction \(z_1 = \mathrm{i}\,\mathrm{Im}(z_1)\) (the axion \(\mathrm{Re}(z_1) = 0\)), with \(z_2 = 2\mathrm{i}\) and \(\tau = 3\mathrm{i}\) fixed. The Kähler potential grows at small \(\mathrm{Im}(z_1)\) (approaching the boundary of moduli space) and decreases logarithmically for large values, as expected from \(K \approx -\log(\kappa_{ijk}\mathrm{Im}(z)^i\mathrm{Im}(z)^j\mathrm{Im}(z)^k)\) at LCS.

# Visualise the Kähler potential K(z1) along the imaginary axis
# Fix z2 = 2i and tau = 3i, sweep Im(z1) from 1.2 to 5.5
im_z1_vals = jnp.linspace(1.2, 5.5, 120)
z2_ref     = 2.0j
tau_ref    = 3.0j

def K_along_z1(im_z1):
    z = jnp.array([im_z1 * 1j, z2_ref])
    return jnp.real(css.kahler_potential(z, jnp.conj(z), tau_ref, jnp.conj(tau_ref)))

K_vals = jax.vmap(K_along_z1)(im_z1_vals)

fig, ax = plt.subplots(figsize=(6, 3.5))
ax.plot(im_z1_vals, K_vals, lw=2, color='steelblue')
ax.axvline(2.742, color='crimson', ls='--', lw=1.5, label=r'vacuum $\mathrm{Im}(z_1^*)$')
ax.set_xlabel(r'$\mathrm{Im}(z_1)$', fontsize=13)
ax.set_ylabel(r'$K(\mathrm{i}\,\mathrm{Im}(z_1),\,2\mathrm{i},\,3\mathrm{i})$', fontsize=12)
ax.set_title(r'Kähler potential along the purely imaginary axis', fontsize=13)
ax.legend(fontsize=11)
sn.despine()
plt.tight_layout()
plt.show()

FluxEFT sub-module#

Flux compactifications of Type IIB string theory#

This sub-module contains FluxEFT class for computations of the flux scalar potential induced by 3-form flux backgrounds. Flux compactifications are a central framework for connecting string theory to lower-dimensional physics. They involve the compactification of the extra spatial dimensions of string theory on a compact manifold, with additional background fluxes stabilising the moduli fields associated with the geometry. In the context of Type IIB string theory, compactification on Calabi-Yau manifolds with background fluxes provides a mechanism to stabilise moduli, a crucial step towards constructing realistic four-dimensional (4D) effective theories.

Type IIB string theory in ten dimensions contains two types of bosonic fields relevant for flux compactifications: the Neveu-Schwarz–Neveu-Schwarz (NSNS) 2-form field \(B_2\) and the Ramond-Ramond (RR) fields \(C_0\), \(C_2\), and \(C_4\). The dynamics of these fields are encoded in the 10D Type IIB supergravity action. Compactifying the theory on a six-dimensional Calabi-Yau threefold \(X\) leads to an effective 4D theory.

A Calabi-Yau threefold is a complex three-dimensional Kähler manifold with a vanishing first Chern class. Such manifolds admit a unique holomorphic \((3,0)\)-form \(\Omega\) and a Kähler form \(J\). The metric on the internal space can be written as \begin{equation} \mathrm{d}s^2_{10} = \mathrm{e}^{2A(y)} g_{\mu\nu} \mathrm{d}x^\mu \mathrm{d}x^\nu + \mathrm{e}^{-2A(y)} g_{mn} \mathrm{d}y^m \mathrm{d}y^n, , \end{equation} where \(g_{\mu\nu}\) is the 4D metric, \(g_{mn}\) is the metric on \(X\), and \(A(y)\) is the warp factor. Luckily, for what we are about to do, information about the CY metric \( g_{mn} \) and the warp factor \( A \) are mostly irrelevant.

It should be noted that there has been significant progress in constructing CY metrics \( g_{mn} \) numerically using machine learning tools. Many open-source software packages have been developed over time, such as cymetric, cyjax, and cymyc. This is another interesting use case for applying neural networks (NNs) to the string landscape that, unfortunately, we will not have time to cover in this lecture.

In flux compactifications, the background fluxes are introduced as integrals of field strengths over the cycles \(\Sigma_\alpha\) of the Calabi-Yau manifold. The relevant field strength is the combined 3-form flux \begin{equation} G_3 = F_3 - \tau H_3, \end{equation} where \(F_3 = \mathrm{d}C_2\) is the RR 3-form flux, \(H_3 = \mathrm{d}B_2\) is the NSNS 3-form flux, and \(\tau = C_0 + \mathrm{i} \, \mathrm{e}^{-\phi}\) is the axio-dilaton. The associated NSNS- and RR-fluxes \((h_\alpha,f_\alpha)\) must satisfy the quantisation condition \begin{equation} h_{\alpha}=\frac{1}{(2\pi)^2 \alpha’} \int_{\Sigma_\alpha} H_3 \in \mathbb{Z}, , \quad f_{\alpha}=\frac{1}{(2\pi)^2 \alpha’} \int_{\Sigma_\alpha} F_3 \in \mathbb{Z}, , \end{equation} where \(\Sigma_\alpha\) are 3-cycles in \(X\). For general Calabi-Yau compactifications, the number of such fluxes is given by \(4(h^{1,2}+1)\) where \(h^{1,2}\) is a topological invariant, a so-called Hodge number, which counts the number of complex structure deformations of the Calabi-Yau threefold. Each such deformation mode is associated with a separate scalar field, i.e., a modulus \(z^i\) which appears in the EFT.

The dynamics of moduli for given fluxes are governed by the Gukov-Vafa-Witten (GVW) superpotential
\begin{equation} W_{\text{GVW}}(\tau,Z)= \int_{X}, G_{3}\wedge\Omega_{3} = \left (f-\tau h\right )^{T}\cdot \Sigma\cdot \Pi(Z), . \end{equation} The associated Kähler potential is given by
\begin{equation} K = -\log\left(-\mathrm{i} \int_X \Omega \wedge \bar{\Omega} \right) - \log\left(-\mathrm{i}(\tau - \bar{\tau})\right). \end{equation} The scalar potential is then given by \begin{equation} V_{\text{Flux}}= \mathrm{e}^{K} \left( K^{\tau\overline{\tau}}D_{\tau} W, D_{\overline{\tau}}\overline{ W}+ K^{i\bar\jmath}D_{i} W, D_{\bar\jmath}\overline{ W} \right), , \end{equation} where \begin{equation} D_\tau W = \partial_\tau W + (\partial_\tau K)W , , \quad D_i W = \partial_i W + (\partial_i K)W , , \end{equation} are the Kähler covariant derivatives. A vacuum solution is a minimum (or ground state) of the potential \(V\) with \(\partial_\tau V=\partial_{i}V=0\). Supersymmetric vacua correspond to solutions to the \(F\)-term equations \begin{equation} D_\tau W = 0 , ,\quad D_i W = 0 , . \end{equation} These conditions can equivalently written as ISD condition, namely

\[ m_{J}-\tau n_{J}-\overline{\mathcal{N}}_{JI}\left (m^{I}-\tau n^{I}\right ) =0\quad\Leftrightarrow\quad f=\left (\mathcal{M}\Sigma s+1c_{0}\right )h \]
For solutions of this type, the vacuum energy \(V_0=\langle V\rangle=0\) vanishes, i.e., we find a flat Minkowski minimum.

Example for FluxEFT class#

Let us again look at our example from above:

model = jvc.FluxEFT(h12=2,model_ID = 1,model_type="KS",maximum_degree=0)
model.lcs_tree.a_matrix = jnp.array([[4.5,1.5],[1.5,0.]])
model

The flux sector is a subclass of css (the complex-structure sector) and thus has access to objects like the Kähler metric etc. In addition, as the name suggests, it implements functions that depend on the 3-form fluxes like the GVW superpotential, the flux scalar potential and more. The flux sector still has the period sector as an attribute which can be called via model.periods which e.g. in the LCS limit contains also information about the topology of the mirror CY (like triple intersection numbers, GVs, etc.).

Now, to showcase the basic functionality, let us look at a specific flux vacuum obtained in 1912.10047

fluxes=jnp.array([7, 3, -24, 0, -16, 50,0, 3, -4, 0, 0, 0])
u1sol=2.74215479602462524879172086700112955631003945168828832743217138983767*1j
u2sol=2.05661613496943436323419976712599580262262253939859294519039244649420*1j
tau0 = 6.85540179778358427172610564536555609784128313762349971439377181031816*1j

z0 = jnp.array([u1sol,u2sol])
ctau0 = jnp.conj(tau0)
cz0 = jnp.conj(z0)

We can check that this solution corresponds to an actual vacuum satisfying the \(F\)-term conditions

model.DW(z0,cz0,tau0,ctau0,fluxes)

We see that the \(F\)-term conditions are approximately satisfied. This is because above we neglected corrections from worldsheet instanton corrections. These can be included by increasing the value for maximum_degree

model_inst = jvc.FluxEFT(h12=2,model_ID = 1,model_type="KS",maximum_degree=2)
model_inst.lcs_tree.a_matrix = jnp.array([[4.5,1.5],[1.5,0.]])
model_inst.DW(z0,cz0,tau0,ctau0,fluxes)

We see that the \(F\)-term conditions including instanton corrections are better approximated by the above solutions.

Other functions pre-implemented include \(\partial_i D_j W\)

model.dDW(z0,cz0,tau0,ctau0,fluxes)

or \(D_i D_j W\)

model.DDW(z0,cz0,tau0,ctau0,fluxes)

Moreover, we implemented the standard \(F\)-term scalar potential induced by 3-form fluxes through the GVW superpotential

model.scalar_potential(z0,cz0,tau0,ctau0,fluxes)

Moreover, we can compute the Hessian and its eigenvalues

H = model.hessian(z0,cz0,tau0,ctau0,fluxes)
H_evs = jnp.linalg.eigvalsh(H)
H_evs

Two eigenvalues are highly suppressed because at the classical level one complex direction remains flat.

Scalar potential landscape#

The flux scalar potential \(V\) is generally a complicated function of the moduli. As a simple visualisation, we plot \(|V(z_1)|\) as a function of \(\mathrm{Im}(z_1)\) with all other fields fixed at the vacuum values \(z_2 = z_2^*\), \(\tau = \tau^*\). The minimum at \(z_1 = z_1^*\) is clearly visible — the potential vanishes exactly there since this is a SUSY (\(D_IW = 0\)) Minkowski vacuum with \(\langle V \rangle = 0\).

# Scalar potential slice: fix z2, tau at the vacuum and sweep Im(z1)
fluxes_plot = jnp.array([7, 3, -24, 0, -16, 50, 0, 3, -4, 0, 0, 0])
z2_vac  = 2.05662j
tau_vac = 6.85540j
im_z1_scan = jnp.linspace(1.5, 5.0, 200)

def V_along_z1(im_z1):
    z = jnp.array([im_z1 * 1j, z2_vac])
    return jnp.real(model.scalar_potential(z, jnp.conj(z), tau_vac, jnp.conj(tau_vac), fluxes_plot))

V_vals = jax.vmap(V_along_z1)(im_z1_scan)

fig, ax = plt.subplots(figsize=(6, 3.5))
ax.semilogy(im_z1_scan, jnp.abs(V_vals), lw=2, color='darkorange')
ax.axvline(2.742, color='crimson', ls='--', lw=1.5, label=r'vacuum $\mathrm{Im}(z_1^*)$')
ax.set_xlabel(r'$\mathrm{Im}(z_1)$', fontsize=13)
ax.set_ylabel(r'$|V(z_1)|$', fontsize=13)
ax.set_title(r'Scalar potential slice with $z_2 = z_2^*$, $\tau = \tau^*$', fontsize=13)
ax.legend(fontsize=11)
sn.despine()
plt.tight_layout()
plt.show()

jaxvacua.sampling sub-module#

The data_sampler class provides utilities for generating random starting points in moduli space and for sampling flux configurations. It wraps a FluxEFT model and respects the Kähler cone constraints on the moduli.

Sampling modes for moduli (moduli_sampling_mode):

Mode

Description

"box"

Uniform sampling in \([\mathrm{Im}(z)_\mathrm{min}, \mathrm{Im}(z)_\mathrm{max}]^{h^{1,2}}\)

"cone"

Samples inside the stretched Kähler cone

"tip_ray"

Samples near the tip of the cone along random rays

"random_ray"

Uniform along a random ray in the cone

For a detailed comparison of these modes, see Notebook 4: Sampling Module.

sampler = jvc.data_sampler(model,flux_bounds=[-5,5])

We can then sample initial guesses by calling

moduli,tau,fluxes = sampler.initial_guesses(10)
moduli,tau,fluxes

ISD sampling#

ISD sampling is the most efficient strategy for finding SUSY flux vacua. The key insight is that any SUSY vacuum must satisfy the imaginary self-dual (ISD) condition \(\star G_3 = \mathrm{i}G_3\), which is equivalent to \(D_IW = 0\). ISD sampling exploits this by:

  1. Drawing a random point \((z^i, \tau)\) in moduli space.

  2. Choosing half of the flux quanta randomly (e.g. the NSNS fluxes \(h\)).

  3. Solving analytically for the remaining half (RR fluxes \(f\)) such that the ISD condition \(f = (s\,\mathcal{M}\Sigma + c_0\mathbf{1})h\) is satisfied, where \(\tau = c_0 + \mathrm{i}s\).

The resulting (rounded) flux is close to ISD, giving a starting point near a true vacuum. Subsequent Newton refinement then converges rapidly. For the algorithmic walkthrough, see NB06: ISD sampling principle; for batched campaigns, see NB07: ISD sampling workflows.

# ISD sampling: generate flux guesses that approximately satisfy the ISD condition
# initial_guesses_ISD samples (z, tau) randomly and solves for the ISD-completing fluxes
moduli_isd, tau_isd, fluxes_isd = sampler.initial_guesses_ISD(
    20, mode="H",
    minval_moduli=1.5, maxval_moduli=4.0,
    minval_dilaton=1.5, maxval_dilaton=5.0,
    minval_fluxes=-5, maxval_fluxes=5
)

# Evaluate F-term residuals — without any optimisation, ISD sampling already
# produces starting points close to a SUSY vacuum
DW_isd = jax.vmap(model.DW)(
    moduli_isd, jnp.conj(moduli_isd), tau_isd, jnp.conj(tau_isd), fluxes_isd)
residuals = jnp.max(jnp.abs(DW_isd), axis=-1)
print(f"Residuals max|D_I W| for 20 ISD-sampled starting points:")
print(residuals)
print(f"\nMean residual: {jnp.mean(residuals):.4f}")

Take-aways#

This notebook introduced the four main classes of JAXVacua — periods, css, FluxEFT, and data_sampler — and demonstrated the basic workflow from loading a model to finding a SUSY flux vacuum.

Further reading#

Next steps:

Notebook

Topic

3: CYTools Interface

Load models directly from a CYTools CalabiYau object

CICY tutorial (in stringforge docs)

Work with Complete Intersection Calabi-Yau threefolds

5: Finding Flux Vacua

Newton’s method and SciPy optimisers in detail

4: Sampling Module

Moduli sampling modes compared

6: ISD Sampling Principle

Algorithmic ISD completion for a single vacuum

7: ISD Sampling Workflows

Batched moduli-first and flux-first campaigns

8: Flux Bounding

Systematic flux enumeration within a D3 tadpole bound

9: Moduli Space Limits

LCS and coni-LCS limits

10: coni-LCS Pipeline

Conifold PFV workflow and freezer setup

11: Hypergeometric Models

Closed-form one-modulus prepotentials