Optimize ACE hyper-parameters: minimize force time and fitting error.

Setup experiment

Load packages.

using AtomsBase, InteratomicPotentials, PotentialLearning
using Unitful, UnitfulAtomic
using LinearAlgebra, Random, DisplayAs
using DataFrames, Hyperopt

Define paths.

base_path = haskey(ENV, "BASE_PATH") ? ENV["BASE_PATH"] : "../../"
ds_path   = "$base_path/examples/data/a-HfO2/a-HfO2-300K-NVT-6000.extxyz"
res_path  = "$base_path/examples/Opt-ACE-aHfO2/results/";

Load utility functions.

include("$base_path/examples/utils/utils.jl");

Create experiment folder.

run(`mkdir -p $res_path`);

Load datasets

Load atomistic dataset: atomistic configurations (atom positions, geometry, etc.) + DFT data (energies, forces, etc.)

ds = load_data(ds_path, uparse("eV"), uparse("Å"))[1:1000]; # Load first 1K samples.

Split atomistic dataset into training and test

n_train, n_test = 50, 50 # Only 50 samples per dataset are used in this example.
conf_train, conf_test = split(ds, n_train, n_test)
(DataSet{num_configs = 50} 
	 Configuration{S, AtomsBase.FlexibleSystem{3, AtomsBase.Atom, Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}, Forces, Energy}
	 Configuration{S, AtomsBase.FlexibleSystem{3, AtomsBase.Atom, Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}, Forces, Energy}
	 ⋮
	 Configuration{S, AtomsBase.FlexibleSystem{3, AtomsBase.Atom, Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}, Forces, Energy}, DataSet{num_configs = 50} 
	 Configuration{S, AtomsBase.FlexibleSystem{3, AtomsBase.Atom, Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}, Forces, Energy}
	 Configuration{S, AtomsBase.FlexibleSystem{3, AtomsBase.Atom, Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}, Forces, Energy}
	 ⋮
	 Configuration{S, AtomsBase.FlexibleSystem{3, AtomsBase.Atom, Unitful.Quantity{Float64, 𝐋, Unitful.FreeUnits{(Å,), 𝐋, nothing}}}, Forces, Energy})

Optimize hyper-parameters

Define a custom loss function. Here, we minimize fitting error and force calculation time. Possible metrics are e_mae, e_rmse, e_rsq, f_mae, f_rmse, f_rsq, and time_us.

function custom_loss(
    metrics::OrderedDict
)
    e_mae     = metrics[:e_mae]
    f_mae     = metrics[:f_mae]
    time_us   = metrics[:time_us]
    e_mae_max = 0.05 # eV/atom
    f_mae_max = 0.05 # eV/Å
    w_e       = e_mae/e_mae_max
    w_f       = f_mae/f_mae_max
    w_t       = 1.0E-3
    loss = w_e * e_mae + w_f * e_mae + w_t * time_us
    return loss
end;

Define model and hyper-parameter value ranges to be optimized.

model = ACE
pars = OrderedDict( :body_order        => [2, 3, 4],
                    :polynomial_degree => [3, 4, 5],
                    :rcutoff           => LinRange(4, 6, 10),
                    :wL                => LinRange(0.5, 1.5, 10),
                    :csp               => LinRange(0.5, 1.5, 10),
                    :r0                => LinRange(0.5, 1.5, 10));

Use latin hypercube sampling to find the optimal hyper-parameters. Alternatively, use random sampling (sampler = RandomSampler()).

sampler = CLHSampler(dims=[Categorical(3), Categorical(3), Continuous(),
                           Continuous(), Continuous(), Continuous()])
iap, res = hyperlearn!(model, pars, conf_train;
                       n_samples = 10, sampler = sampler,
                       loss = custom_loss, ws = [1.0, 1.0], int = true);
E_MAE:0.184 eV/atom, F_MAE:0.251 eV/Å, Time per force per atom:83.751 µs
E_MAE:0.137 eV/atom, F_MAE:0.217 eV/Å, Time per force per atom:111.008 µs
E_MAE:0.183 eV/atom, F_MAE:0.313 eV/Å, Time per force per atom:47.608 µs
E_MAE:0.173 eV/atom, F_MAE:0.302 eV/Å, Time per force per atom:56.41 µs
E_MAE:0.135 eV/atom, F_MAE:0.111 eV/Å, Time per force per atom:328.504 µs
E_MAE:0.196 eV/atom, F_MAE:0.175 eV/Å, Time per force per atom:93.402 µs
E_MAE:0.097 eV/atom, F_MAE:0.11 eV/Å, Time per force per atom:349.337 µs
E_MAE:0.134 eV/atom, F_MAE:0.22 eV/Å, Time per force per atom:86.827 µs
E_MAE:0.22 eV/atom, F_MAE:0.301 eV/Å, Time per force per atom:115.86 µs
E_MAE:0.127 eV/atom, F_MAE:0.207 eV/Å, Time per force per atom:149.277 µs

Post-process results

Save and show results.

@save_var res_path iap.β
@save_var res_path iap.β0
@save_var res_path iap.basis
@save_dataframe res_path res
res
10×13 DataFrame
Rowe_maee_rmsee_rsqf_maef_rmsef_rsqtime_usbody_orderpolynomial_degreercutoffwLcspr0
AnyAnyAnyAnyAnyAnyAnyAnyAnyAnyAnyAnyAny
10.09706510.1137970.8378460.1098060.1411260.945789349.3372.04.04.444440.8333330.50.722222
20.1271960.159330.6821230.2074070.2682140.804191149.2773.05.06.01.166670.7222221.05556
30.1337090.1694420.6404930.2196470.2866430.77635886.82722.04.04.444440.8333330.50.722222
40.1346790.1513320.7132360.1106390.1413790.945595328.5042.04.04.00.7222221.277781.16667
50.1374390.172790.6261470.216640.2825340.782724111.0084.04.04.888891.388891.388891.38889
60.1725970.2194370.3970430.3021790.3961940.57274656.40984.03.05.111110.6111110.9444440.611111
70.1826950.2163320.4139850.3125260.4095850.54337647.60844.03.05.111110.6111110.9444440.611111
80.1838230.236560.2992750.2512160.3263660.71007983.75074.04.04.888891.388891.388891.38889
90.1962390.2244150.3693780.1752720.2251640.86200493.4022.04.04.00.7222221.277781.16667
100.2204470.2778090.03359770.3009630.3929720.579666115.863.05.06.01.166670.7222221.05556

Plot error vs time.

err_time = plot_err_time(res)
@save_fig res_path err_time
DisplayAs.PNG(err_time)
Example block output

This page was generated using Literate.jl.