Calculators

Calculator interfaces for molecular simulations using AIMNet2.

AIMNet2Calculator

The core calculator for running AIMNet2 inference. It handles model loading, device management, and application of long-range interactions (Coulomb and Dispersion).

Key Features

  • Format Support: Loads both legacy .jpt models and new .pt format.
  • Long-Range Interactions: Automatically attaches LRCoulomb and DFTD3 modules based on model metadata.
  • Overrides: You can force specific long-range behavior using needs_coulomb and needs_dispersion arguments.
  • Batching: Automatically batches large molecules/systems based on nb_threshold.

AIMNet2Calculator(model='aimnet2', nb_threshold=120, needs_coulomb=None, needs_dispersion=None, device=None, compile_model=False, compile_kwargs=None, train=False)

Generic AIMNet2 calculator.

A helper class to load AIMNet2 models and perform inference.

Parameters

model : str | nn.Module Model name (from registry), path to model file, or nn.Module instance. nb_threshold : int Threshold for neighbor list batching. Molecules larger than this use flattened processing. Default is 120. needs_coulomb : bool | None Whether to add external Coulomb module. If None (default), determined from model metadata. If True/False, overrides metadata. needs_dispersion : bool | None Whether to add external DFTD3 module. If None (default), determined from model metadata. If True/False, overrides metadata. device : str | None Device to run the model on ("cuda", "cpu", or specific like "cuda:0"). If None (default), auto-detects CUDA availability. compile_model : bool Whether to compile the model with torch.compile(). Default is False. compile_kwargs : dict | None Additional keyword arguments to pass to torch.compile(). Default is None. train : bool Whether to enable training mode. Default is False (inference mode). When False, all model parameters have requires_grad=False, which improves torch.compile compatibility and reduces memory usage. Set to True only when training the model.

Attributes

model : nn.Module The loaded AIMNet2 model. device : str Device the model is running on ("cuda" or "cpu"). cutoff : float Short-range cutoff distance in Angstroms. cutoff_lr : float | None Long-range cutoff distance, or None if no LR modules. external_coulomb : LRCoulomb | None External Coulomb module if attached. external_dftd3 : DFTD3 | None External DFTD3 module if attached.

Notes

External LR module behavior:

  • For file-loaded models (str): metadata is loaded from file
  • For nn.Module: metadata is read from model.metadata attribute if available
  • Explicit flags (needs_coulomb, needs_dispersion) override metadata
  • If no metadata and no explicit flags, no external LR modules are added
Source code in aimnet/calculators/calculator.py
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
def __init__(
    self,
    model: str | nn.Module = "aimnet2",
    nb_threshold: int = 120,
    needs_coulomb: bool | None = None,
    needs_dispersion: bool | None = None,
    device: str | None = None,
    compile_model: bool = False,
    compile_kwargs: dict | None = None,
    train: bool = False,
):
    # Device selection: use provided or auto-detect
    if device is not None:
        self.device = device
    else:
        self.device = "cuda" if torch.cuda.is_available() else "cpu"
    self.external_coulomb: LRCoulomb | None = None
    self.external_dftd3: DFTD3 | None = None
    # Default cutoffs for LR modules
    self._default_dsf_cutoff = 15.0
    self._default_dftd3_cutoff = 15.0
    self._default_dftd3_smoothing = 0.2

    # Load model and get metadata
    metadata: dict | None = None
    if isinstance(model, str):
        p = get_model_path(model)
        self.model, metadata = load_model(p, device=self.device)
        self.cutoff = metadata["cutoff"]
    elif isinstance(model, nn.Module):
        self.model = model.to(self.device)
        self.cutoff = getattr(self.model, "cutoff", 5.0)
        metadata = getattr(self.model, "_metadata", None)
    else:
        raise TypeError("Invalid model type/name.")

    # Compile model if requested
    if compile_model:
        kwargs = compile_kwargs or {}
        self.model = torch.compile(self.model, **kwargs)

    # Resolve final flags (explicit overrides metadata)
    final_needs_coulomb = (
        needs_coulomb
        if needs_coulomb is not None
        else (metadata.get("needs_coulomb", False) if metadata is not None else False)
    )
    final_needs_dispersion = (
        needs_dispersion
        if needs_dispersion is not None
        else (metadata.get("needs_dispersion", False) if metadata is not None else False)
    )

    # Set up external Coulomb if needed
    if final_needs_coulomb:
        sr_embedded = metadata.get("coulomb_mode") == "sr_embedded" if metadata is not None else False
        # For PBC, user can switch to DSF/Ewald via set_lrcoulomb_method()
        # When sr_embedded=True: model has SRCoulomb which subtracts SR, so external
        # should compute FULL (subtract_sr=False) to give: (NN - SR) + FULL = NN + LR
        # When sr_embedded=False: model has no SR embedded, so external should compute
        # LR only (subtract_sr=True) to avoid double-counting
        self.external_coulomb = LRCoulomb(
            key_in="charges",
            key_out="energy",
            method="simple",
            rc=metadata.get("coulomb_sr_rc", 4.6) if metadata is not None else 4.6,
            envelope=metadata.get("coulomb_sr_envelope", "exp") if metadata is not None else "exp",
            subtract_sr=not sr_embedded,
        )
        self.external_coulomb = self.external_coulomb.to(self.device)

    # Set up external DFTD3 if needed
    if final_needs_dispersion:
        d3_params = metadata.get("d3_params") if metadata else None
        if d3_params is None:
            raise ValueError(
                "needs_dispersion=True but d3_params not found in metadata. "
                "Provide d3_params in model metadata or set needs_dispersion=False."
            )
        self.external_dftd3 = DFTD3(
            s8=d3_params["s8"],
            a1=d3_params["a1"],
            a2=d3_params["a2"],
            s6=d3_params.get("s6", 1.0),
        )
        self.external_dftd3 = self.external_dftd3.to(self.device)

    # Determine if model has long-range modules (embedded or external)
    has_embedded_lr = metadata.get("has_embedded_lr", False) if metadata is not None else False
    self.lr = (
        hasattr(self.model, "cutoff_lr")
        or self.external_coulomb is not None
        or self.external_dftd3 is not None
        or has_embedded_lr
    )
    # Set cutoff_lr based on model attribute or external modules
    if hasattr(self.model, "cutoff_lr"):
        self.cutoff_lr = getattr(self.model, "cutoff_lr", float("inf"))
    elif self.external_coulomb is not None:
        # For "simple" method, use inf (all pairs). For DSF, use dsf_cutoff.
        if self.external_coulomb.method == "simple":
            self.cutoff_lr = float("inf")
        else:
            self.cutoff_lr = self._default_dsf_cutoff
    elif self.external_dftd3 is not None:
        self.cutoff_lr = self._default_dftd3_cutoff
    elif has_embedded_lr:
        # Embedded LR modules (D3TS, SRCoulomb) need nbmat_lr
        self.cutoff_lr = self._default_dftd3_cutoff
    else:
        self.cutoff_lr = None
    self.nb_threshold = nb_threshold

    # Create adaptive neighbor list instances
    self._nblist = AdaptiveNeighborList(cutoff=self.cutoff)

    # Track separate cutoffs for LR modules
    self._coulomb_cutoff: float | None = None
    self._dftd3_cutoff: float = self._default_dftd3_cutoff
    if self.external_coulomb is not None:
        if self.external_coulomb.method == "simple":
            self._coulomb_cutoff = float("inf")
        elif self.external_coulomb.method == "ewald":
            self._coulomb_cutoff = None  # Ewald manages its own cutoff
        else:
            self._coulomb_cutoff = self.external_coulomb.dsf_rc
    if self.external_dftd3 is not None:
        self._dftd3_cutoff = self.external_dftd3.smoothing_off

    # Create long-range neighbor list(s) if LR modules present
    self._nblist_lr: AdaptiveNeighborList | None = None
    self._nblist_dftd3: AdaptiveNeighborList | None = None
    self._nblist_coulomb: AdaptiveNeighborList | None = None
    self._update_lr_nblists()

    # indicator if input was flattened
    self._batch = None
    self._max_mol_size: int = 0
    # placeholder for tensors that require grad
    self._saved_for_grad = {}
    # set flag of current Coulomb method
    self._coulomb_method: str | None = None
    if self.external_coulomb is not None:
        self._coulomb_method = self.external_coulomb.method
    elif self._has_embedded_coulomb():
        # Legacy models have embedded Coulomb with "simple" method
        self._coulomb_method = "simple"

    # Set training mode (default False for inference)
    self._train = train
    self.model.train(train)
    if not train:
        # Disable gradients on all parameters for inference mode
        for param in self.model.parameters():
            param.requires_grad_(False)
        if self.external_coulomb is not None:
            for param in self.external_coulomb.parameters():
                param.requires_grad_(False)
        if self.external_dftd3 is not None:
            for param in self.external_dftd3.parameters():
                param.requires_grad_(False)

coulomb_cutoff property

Get the current Coulomb cutoff distance.

Returns

float | None The cutoff distance for Coulomb calculations, or None if not applicable. For "simple" method, this is inf. For "ewald", this is None. Use set_lrcoulomb_method() to change.

coulomb_method property

Get the current Coulomb method.

Returns

str | None One of "simple", "dsf", "ewald", or None if no external Coulomb. For legacy models with embedded Coulomb, returns None.

dftd3_cutoff property

Get the current DFTD3 cutoff distance.

Returns

float The cutoff distance for DFTD3 calculations in Angstroms.

has_external_coulomb property

Check if calculator has external Coulomb module attached.

Returns True for new-format models that were trained with Coulomb and have it externalized. For legacy models, Coulomb is embedded in the model itself, so this returns False.

has_external_dftd3 property

Check if calculator has external DFTD3 module attached.

Returns True for new-format models that were trained with DFTD3/D3BJ dispersion and have it externalized. For legacy models or D3TS models, dispersion is embedded in the model itself, so this returns False.

mol_flatten(data)

Flatten the input data for multiple molecules. Will not flatten for batched input and molecule size below threshold.

Source code in aimnet/calculators/calculator.py
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
def mol_flatten(self, data: dict[str, Tensor]) -> dict[str, Tensor]:
    """Flatten the input data for multiple molecules.
    Will not flatten for batched input and molecule size below threshold.
    """
    ndim = data["coord"].ndim
    if ndim == 2:
        self._batch = None
        if "mol_idx" not in data:
            data["mol_idx"] = torch.zeros(data["coord"].shape[0], dtype=torch.long, device=self.device)
            self._max_mol_size = data["coord"].shape[0]
        elif data["mol_idx"][-1] == 0:
            self._max_mol_size = len(data["mol_idx"])
        else:
            self._max_mol_size = data["mol_idx"].unique(return_counts=True)[1].max().item()

    elif ndim == 3:
        B, N = data["coord"].shape[:2]
        # Force flattening for PBC (cell present) to ensure make_nbmat computes proper neighbor lists with shifts
        if self.nb_threshold < N or self.device == "cpu" or data.get("cell") is not None:
            self._batch = B
            data["mol_idx"] = torch.repeat_interleave(
                torch.arange(0, B, device=self.device), torch.full((B,), N, device=self.device)
            )
            for k, v in data.items():
                if k in self.atom_feature_keys:
                    data[k] = v.flatten(0, 1)
        else:
            self._batch = None
        self._max_mol_size = N
    return data

set_dftd3_cutoff(cutoff=None, smoothing_fraction=None)

Set DFTD3 cutoff and smoothing.

Parameters

cutoff : float | None Cutoff distance in Angstroms for DFTD3 calculation. Default is _default_dftd3_cutoff (15.0). smoothing_fraction : float | None Fraction of cutoff used as smoothing width. Default is _default_dftd3_smoothing (0.2).

Notes

This method only affects external DFTD3 modules attached to new-format models. For legacy models with embedded DFTD3, the smoothing is fixed.

Updates _dftd3_cutoff and rebuilds neighbor lists.

Source code in aimnet/calculators/calculator.py
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
def set_dftd3_cutoff(self, cutoff: float | None = None, smoothing_fraction: float | None = None) -> None:
    """Set DFTD3 cutoff and smoothing.

    Parameters
    ----------
    cutoff : float | None
        Cutoff distance in Angstroms for DFTD3 calculation.
        Default is _default_dftd3_cutoff (15.0).
    smoothing_fraction : float | None
        Fraction of cutoff used as smoothing width.
        Default is _default_dftd3_smoothing (0.2).

    Notes
    -----
    This method only affects external DFTD3 modules attached to
    new-format models. For legacy models with embedded DFTD3,
    the smoothing is fixed.

    Updates _dftd3_cutoff and rebuilds neighbor lists.
    """
    if cutoff is None:
        cutoff = self._default_dftd3_cutoff
    if smoothing_fraction is None:
        smoothing_fraction = self._default_dftd3_smoothing

    self._dftd3_cutoff = cutoff
    if self.external_dftd3 is not None:
        self.external_dftd3.set_smoothing(cutoff, smoothing_fraction)
    self._update_lr_nblists()

set_lr_cutoff(cutoff)

Set the unified long-range cutoff for all LR modules.

Parameters

cutoff : float Cutoff distance in Angstroms for LR neighbor lists.

Notes

This updates both _coulomb_cutoff and _dftd3_cutoff. Ewald uses its own internal neighbor list and ignores this cutoff.

Source code in aimnet/calculators/calculator.py
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
def set_lr_cutoff(self, cutoff: float) -> None:
    """Set the unified long-range cutoff for all LR modules.

    Parameters
    ----------
    cutoff : float
        Cutoff distance in Angstroms for LR neighbor lists.

    Notes
    -----
    This updates both _coulomb_cutoff and _dftd3_cutoff.
    Ewald uses its own internal neighbor list and ignores this cutoff.
    """
    # Update both cutoffs (but not for ewald which manages its own)
    if self._coulomb_method != "ewald":
        self._coulomb_cutoff = cutoff
    self._dftd3_cutoff = cutoff
    self.cutoff_lr = cutoff
    self._update_lr_nblists()

set_lrcoulomb_method(method, cutoff=15.0, dsf_alpha=0.2, ewald_accuracy=1e-08)

Set the long-range Coulomb method.

Parameters

method : str One of "simple", "dsf", or "ewald". cutoff : float Cutoff distance for DSF neighbor list. Default is 15.0. Not used for Ewald (which computes cutoffs from accuracy). dsf_alpha : float Alpha parameter for DSF method. Default is 0.2. ewald_accuracy : float Target accuracy for Ewald summation. Controls the real-space and reciprocal-space cutoffs. Lower values give higher accuracy but require more computation. Default is 1e-8.

The Ewald cutoffs are computed as:
- eta = (V^2 / N)^(1/6) / sqrt(2*pi)
- cutoff_real = sqrt(-2 * ln(accuracy)) * eta
- cutoff_recip = sqrt(-2 * ln(accuracy)) / eta
Notes

For new-format models with external Coulomb, this updates the external module. For legacy models with embedded Coulomb, a warning is issued as those modules cannot be modified at runtime.

Source code in aimnet/calculators/calculator.py
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
def set_lrcoulomb_method(
    self,
    method: Literal["simple", "dsf", "ewald"],
    cutoff: float = 15.0,
    dsf_alpha: float = 0.2,
    ewald_accuracy: float = 1e-8,
):
    """Set the long-range Coulomb method.

    Parameters
    ----------
    method : str
        One of "simple", "dsf", or "ewald".
    cutoff : float
        Cutoff distance for DSF neighbor list. Default is 15.0.
        Not used for Ewald (which computes cutoffs from accuracy).
    dsf_alpha : float
        Alpha parameter for DSF method. Default is 0.2.
    ewald_accuracy : float
        Target accuracy for Ewald summation. Controls the real-space
        and reciprocal-space cutoffs. Lower values give higher accuracy
        but require more computation. Default is 1e-8.

        The Ewald cutoffs are computed as:
        - eta = (V^2 / N)^(1/6) / sqrt(2*pi)
        - cutoff_real = sqrt(-2 * ln(accuracy)) * eta
        - cutoff_recip = sqrt(-2 * ln(accuracy)) / eta

    Notes
    -----
    For new-format models with external Coulomb, this updates the external module.
    For legacy models with embedded Coulomb, a warning is issued as those modules
    cannot be modified at runtime.
    """
    if method not in ("simple", "dsf", "ewald"):
        raise ValueError(f"Invalid method: {method}")

    # Warn if model has embedded Coulomb (legacy models)
    if self._has_embedded_coulomb() and self.external_coulomb is None:
        warnings.warn(
            "Model has embedded Coulomb module (legacy format). "
            "set_lrcoulomb_method() only affects external Coulomb modules. "
            "For legacy models, the Coulomb method cannot be changed at runtime.",
            stacklevel=2,
        )

    # Update external LRCoulomb module if present
    if self.external_coulomb is not None:
        self.external_coulomb.method = method
        if method == "dsf":
            self.external_coulomb.dsf_alpha = dsf_alpha
            self.external_coulomb.dsf_rc = cutoff
        elif method == "ewald":
            self.external_coulomb.ewald_accuracy = ewald_accuracy

    # Update _coulomb_cutoff based on method
    if method == "simple":
        self._coulomb_cutoff = float("inf")
    elif method == "dsf":
        self._coulomb_cutoff = cutoff
    elif method == "ewald":
        self._coulomb_cutoff = None  # Ewald manages its own real-space cutoff

    # Update cutoff_lr for backward compatibility
    if self._coulomb_cutoff is not None:
        self.cutoff_lr = self._coulomb_cutoff
    else:
        # Ewald - use DFTD3 cutoff if available, else None
        self.cutoff_lr = self._dftd3_cutoff if self.external_dftd3 is not None else None

    self._coulomb_method = method
    self._update_lr_nblists()

options: show_root_heading: true show_source: true

AIMNet2ASE

ASE (Atomic Simulation Environment) calculator interface.

!!! note "Installation" Requires the ase extra: pip install aimnet[ase]

This calculator integrates with ASE's Atoms object, supporting energy, forces, stress, and dipole moment calculations. It operates in eV and Angstrom.

Usage Example

from ase.io import read
from aimnet.calculators import AIMNet2ASE

atoms = read("molecule.xyz")
atoms.calc = AIMNet2ASE("aimnet2")

print(atoms.get_potential_energy())
print(atoms.get_forces())

AIMNet2ASE(base_calc='aimnet2', charge=0, mult=1)

Bases: Calculator

Source code in aimnet/calculators/aimnet2ase.py
24
25
26
27
28
29
30
31
32
33
34
35
36
37
def __init__(self, base_calc: AIMNet2Calculator | str = "aimnet2", charge=0, mult=1):
    super().__init__()
    if isinstance(base_calc, str):
        base_calc = AIMNet2Calculator(base_calc)
    self.base_calc = base_calc
    self.reset()
    self.charge = charge
    self.mult = mult
    self.update_tensors()
    # list of implemented species
    if hasattr(base_calc, "implemented_species"):
        self.implemented_species = base_calc.implemented_species.cpu().numpy()  # type: ignore
    else:
        self.implemented_species = None

options: show_root_heading: true show_source: true

AIMNet2Pysis

PySisyphus calculator interface.

!!! note "Installation" Requires the pysis extra: pip install aimnet[pysis]

This interface adapts AIMNet2 for use with PySisyphus optimizers. It handles unit conversion automatically:

  • Input: Converts Angstrom (PySisyphus) to Angstrom (AIMNet2).
  • Output: Converts eV/Angstrom (AIMNet2) to Hartree/Bohr (PySisyphus).

AIMNet2Pysis(model='aimnet2', charge=0, mult=1, **kwargs)

Bases: Calculator

Source code in aimnet/calculators/aimnet2pysis.py
21
22
23
24
25
def __init__(self, model: AIMNet2Calculator | str = "aimnet2", charge=0, mult=1, **kwargs):
    super().__init__(charge=charge, mult=mult, **kwargs)
    if isinstance(model, str):
        model = AIMNet2Calculator(model)
    self.model = model

options: show_root_heading: true show_source: true

Model Registry

Utilities for loading pre-trained models. Models are automatically downloaded from the remote repository to a local cache (aimnet/calculators/assets/) upon first use.

CLI Command

You can clear the local model cache using the CLI:

aimnet clear_model_cache

model_registry

options: show_root_heading: true show_source: true