Physical Simulation of Silicon Dangling Bond Logic

These headers provide functions for physically simulating an SiDB layout, which is a crucial step in the physical design flow of SiDB layouts, as they are used to validate their functionality.

Physical Parameters

Header: fiction/algorithms/simulation/sidb/sidb_simulation_parameters.hpp

struct sidb_simulation_parameters

This struct collects all physical parameters for physical SiDB simulations. It can be useful to adjust them, especially when experiments create new insights. However, the default values are commonly used.

Public Functions

inline explicit constexpr sidb_simulation_parameters(const uint8_t base_number = 3, const double mu = -0.32, const double relative_permittivity = 5.6, const double screening_distance = 5.0)

Default constructor.

Parameters:
  • base_number – simulation can be conducted with 2 and 3 charge states. 2 = (Negative, Neutral), 3 = (Negative, Neutral, Positive).

  • mu – (µ-) is the energy transition level (0/-) in eV.

  • relative_permittivity – it describes the electric field reduction due to polarization.

  • screening_distance – also known as “Thomas-Fermi screening” and it describes the electric field screening due to free charges in nm.

inline double k() const noexcept

k is the Coulomb constant K_E divided by epsilon_r (unit: \(N \cdot m^{2} \cdot C^{-2}\)).

inline double mu_plus() const noexcept

mu_plus (µ+) is the energy transition level (+/0) (unit: eV).

Public Members

double epsilon_r

epsilon_r is the electric permittivity. It is a material specific number (unit-less).

double lambda_tf

lambda_tf is the Thomas-Fermi screening distance (unit: nm).

double mu_minus

mu_minus (µ-) is the energy transition level (0/-) (unit: eV).

uint8_t base

base can be either 2 or 3 and describes the assumed number of charge states of one SiDB. It often makes sense to assume only negatively and neutrally charged SiDBs.

Simulation Result

Header: fiction/algorithms/simulation/sidb/sidb_simulation_result.hpp

template<typename Lyt>
struct sidb_simulation_result

This struct defines a unified return type for all SiDB simulation algorithms. It contains the name of the algorithm, the total simulation runtime, the charge distributions determined by the algorithm, the physical parameters used in the simulation, and (optional) algorithm-specific named simulation parameters.

Template Parameters:

Lyt – Cell-level layout type.

Public Functions

inline sidb_simulation_result() noexcept

Default constructor. It only exists to allow for the use of static_assert statements that restrict the type of Lyt.

Public Members

std::string algorithm_name = {}

Name of the algorithm used to determine the charge distributions.

std::chrono::duration<double> simulation_runtime = {}

Total simulation runtime.

std::vector<charge_distribution_surface<Lyt>> charge_distributions = {}

Charge distributions determined by the algorithm.

sidb_simulation_parameters simulation_parameters = {}

Physical parameters used in the simulation.

std::unordered_map<std::string, std::any> additional_simulation_parameters = {}

Additional named simulation parameters. This is used to store algorithm-dependent parameters that are not part of the sidb_simulation_parameters struct.

The key of the map is the name of the parameter, the element is the value of the parameter.

Heuristic Ground State Simulation

Header: fiction/algorithms/simulation/sidb/quicksim.hpp

struct quicksim_params

This struct stores the parameters for the QuickSim algorithm.

Public Members

sidb_simulation_parameters simulation_parameters = {}

Simulation parameters for the simulation of the physical SiDB system.

uint64_t iteration_steps = {80}

Number of iterations to run the simulation for.

double alpha = {0.7}

alpha parameter for the QuickSim algorithm (should be reduced if no result is found).

uint64_t number_threads = {std::thread::hardware_concurrency()}

Number of threads to spawn. By default the number of threads is set to the number of available hardware threads.

template<typename Lyt>
sidb_simulation_result<Lyt> fiction::quicksim(const Lyt &lyt, const quicksim_params &ps = quicksim_params{})

The QuickSim algorithm which was proposed in “QuickSim: Efficient and Accurate Physical Simulation of Silicon Dangling Bond Logic” by J. Drewniok, M. Walter, S. S. H. Ng, K. Walus, and R. Wille in IEEE NANO 2023 (https://ieeexplore.ieee.org/document/10231266) is an electrostatic ground state simulation algorithm for SiDB layouts. It determines physically valid charge configurations (with minimal energy) of a given (already initialized) charge distribution layout. Depending on the simulation parameters, the ground state is found with a certain probability after one run.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – The layout to simulate.

  • ps – Physical parameters. They are material-specific and may vary from experiment to experiment.

Returns:

sidb_simulation_result is returned with all results.

Exhaustive Ground State Simulation

Header: fiction/algorithms/simulation/sidb/quickexact.hpp

template<typename CellType>
struct quickexact_params

This struct stores the parameters for the QuickExact algorithm.

Public Types

enum class automatic_base_number_detection

Modes to use for the QuickExact algorithm.

Values:

enumerator ON

Simulation is conducted with the required base number (i.e., if positively charged SiDBs can occur, three state simulation is conducted).

enumerator OFF

The base number from the physical parameter is used for the simulation.

Public Members

sidb_simulation_parameters simulation_parameters = {}

All parameters for physical SiDB simulations.

automatic_base_number_detection base_number_detection = automatic_base_number_detection::ON

If ON, QuickExact checks which base number is required for the simulation, i.e., whether 3-state is necessary or 2-state simulation is sufficient.

std::unordered_map<CellType, double> local_external_potential = {}

Local external electrostatic potentials (e.g., locally applied electrodes).

double global_potential = 0

Global external electrostatic potential. Value is applied on each cell in the layout.

template<typename Lyt>
sidb_simulation_result<Lyt> fiction::quickexact(const Lyt &lyt, const quickexact_params<cell<Lyt>> &params = {}) noexcept

QuickExact is a quick and exact physical simulation algorithm designed specifically for SiDB layouts. It was proposed in “The Need for Speed: Efficient Exact Simulation of Silicon Dangling Bond Logic” by J. Drewniok, M. Walter, and R. Wille in ASP-DAC 2024 (https://ieeexplore.ieee.org/document/10473946). It determines all physically valid charge configurations of a given SiDB layout, providing a significant performance advantage of more than three orders of magnitude over ExGS (exhaustive_ground_state_simulation).

The performance improvement of QuickExact can be attributed to the incorporation of three key ideas:

  1. Advanced Negative SiDB Detection: QuickExact efficiently identifies SiDBs that require negative charges in a physically valid charge distribution. By pre-assigned them in advance, the search space is pruned by a factor of \(2^k\), where k is the number of found SiDBs.

  2. Dependent SiDB Selection: The algorithm selects a dependent SiDB, whose charge state is always derived from its n-1 neighbors. This dependency simplifies the computation process and contributes to the overall efficiency of QuickExact.

  3. Gray Code Representation: QuickExact employs Gray code to represent and traverse through all charge configurations. By using Gray code, only one charge state changes at a time, making the computation of the local electrostatic potential easier.

Additionally, QuickExact also considers global and local electrostatic potentials, as well as existing defects. This holistic approach ensures an accurate representation of the physical behavior of the SiDB layout.

In summary, QuickExact combines advanced SiDB charge detection, dependent SiDB selection, and the use of Gray code to achieve outstanding performance and enable efficient simulations of SiDB layouts, even in scenarios where positively-charged SiDBs occur due to small spacing.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – Layout to simulate.

  • params – Parameter required for the simulation.

Returns:

Simulation Results.

Header: fiction/algorithms/simulation/sidb/exhaustive_ground_state_simulation.hpp

template<typename Lyt>
sidb_simulation_result<Lyt> fiction::exhaustive_ground_state_simulation(const Lyt &lyt, const sidb_simulation_parameters &params = sidb_simulation_parameters{}) noexcept

Exhaustive Ground State Simulation (ExGS) which was proposed in “Computer-Aided Design of Atomic Silicon Quantum Dots and Computational Applications” by S. S. H. Ng (https://dx.doi.org/10.14288/1.0392909) computes all physically valid charge configurations of a given SiDB layout. All possible charge configurations are passed and checked for physical validity. As a consequence, its runtime grows exponentially with the number of SiDBs per layout. Therefore, only layouts with up to 30 SiDBs can be simulated in a reasonable time. However, since all charge configurations are checked for validity, 100 % simulation accuracy is guaranteed.

Note

This was the first exact simulation approach. However, it is replaced by QuickExact due to the much better runtimes and more functionality.

Template Parameters:

Lyt – Cell-level layout type.

Parameters:
  • lyt – The layout to simulate.

  • params – Simulation parameters.

  • ps – Simulation statistics.

Returns:

sidb_simulation_result is returned with all results.

Engine Selectors

Header: fiction/algorithms/simulation/sidb/sidb_simulation_engine.hpp

enum class fiction::sidb_simulation_engine

Selector for the available SiDB simulation engines.

Values:

enumerator EXGS

Exhaustive Ground State Search (EXGS) is an exact simulation engine that always has exponential runtime.

enumerator QUICKSIM

QuickSim is a heuristic simulation engine that only requires polynomial runtime.

enumerator QUICKEXACT

QuickExact is also an exact simulation engine that requires exponential runtime, but it scales a lot better than ExGS due to its effective search-space pruning.

enum class fiction::exhaustive_sidb_simulation_engine

Selector exclusively for exhaustive SiDB simulation engines.

Values:

enumerator EXGS

Exhaustive Ground State Search (EXGS) is an exact simulation engine that always has exponential runtime.

enumerator QUICKEXACT

QuickExact is also an exact simulation engine that requires exponential runtime, but it scales a lot better than ExGS due to its effective search-space pruning.

Energy Calculation

Header: fiction/algorithms/simulation/sidb/energy_distribution.hpp

using fiction::sidb_energy_distribution = std::map<double, uint64_t>

Data type to collect electrostatic potential energies (unit: eV) of charge distributions with corresponding degeneracy (i.e., how often a certain energy value occurs).

template<typename Lyt>
sidb_energy_distribution fiction::energy_distribution(const std::vector<charge_distribution_surface<Lyt>> &input_vec) noexcept

This function takes in a vector of charge_distribution_surface objects and returns a map containing the system energy and the number of occurrences of that energy in the input vector.

Template Parameters:

Lyt – Cell-level layout type.

Parameters:

input_vec – A vector of charge_distribution_surface objects for which statistics are to be computed.

Returns:

A map containing the system energy as the key and the number of occurrences of that energy in the input vector as the value.

Header: fiction/algorithms/simulation/sidb/minimum_energy.hpp

template<typename InputIt>
double fiction::minimum_energy(const InputIt first, const InputIt last) noexcept

Computes the minimum energy of a range of charge_distribution_surface objects. If the range is empty, infinity is returned.

Template Parameters:

InputIt – Must meet the requirements of LegacyInputIterator.

Parameters:
  • first – Begin of the range to examime.

  • last – End of the range to examine.

Returns:

Value of the minimum energy found in the input range (unit: eV), or infinity if the range is empty.

template<typename InputIt>
InputIt fiction::minimum_energy_distribution(const InputIt first, const InputIt last) noexcept

Returns an iterator to the charge distribution of minimum energy contained in a range of charge_distribution_surface objects. If the range is empty, last is returned.

Template Parameters:

InputIt – Must meet the requirements of LegacyInputIterator.

Parameters:
  • first – Begin of the range to examime.

  • last – End of the range to examine.

Returns:

Iterator to the minimum energy charge distribution found in the input range, or last if the range is empty.

Header: fiction/algorithms/simulation/sidb/is_ground_state.hpp

template<typename Lyt>
bool fiction::is_ground_state(const sidb_simulation_result<Lyt> &heuristic_results, const sidb_simulation_result<Lyt> &exhaustive_results) noexcept

This function checks if the ground state is found by the QuickSim algorithm.

Template Parameters:

Lyt – Cell-level layout type.

Parameters:
  • heuristic_results – All found physically valid charge distribution surfaces obtained by a heuristic algorithm.

  • exhaustive_results – All valid charge distribution surfaces determined by ExGS.

Returns:

Returns true if the relative difference between the lowest energies of the two sets is less than \(0.00001\), false otherwise.

Temperature Behavior

Header: fiction/algorithms/simulation/sidb/critical_temperature.hpp

struct critical_temperature_params

This struct stores the parameters for the Critical Temperature` algorithm.

Public Types

enum class simulation_engine : uint8_t

An enumeration of simulation modes (exact vs. approximate) to use for the Critical Temperature Simulation.

Values:

enumerator EXACT

This simulation engine computes Critical Temperature values with 100 % accuracy.

enumerator APPROXIMATE

This simulation engine quickly calculates the Critical Temperature. However, there may be deviations from the exact Critical Temperature. This mode is recommended for larger layouts (> 40 SiDBs).

Public Members

sidb_simulation_parameters simulation_parameters = {}

All parameters for physical SiDB simulations.

simulation_engine engine = simulation_engine::EXACT

Simulation mode to determine the Critical Temperature.

double confidence_level = {0.99}

Probability threshold for ground state population. The temperature at which the simulation finds the ground state to be populated with a probability of less than the given percentage, is determined to be the critical temperature. For gate-based simulation, this is the probability of erroneous calculations of the gate.

double max_temperature = {400}

Maximum simulation temperature beyond which no simulation will be conducted (~ 126 °C by default) (unit: K).

detect_bdl_pairs_params bdl_params = {}

Parameters for the BDL pair detection algorithms.

uint64_t iteration_steps = {80}

Number of iteration steps for the QuickSim algorithm (only applicable if engine == APPROXIMATE).

double alpha = {0.7}

Alpha parameter for the QuickSim algorithm (only applicable if engine == APPROXIMATE).

template<typename Lyt, typename TT>
double fiction::critical_temperature_gate_based(const Lyt &lyt, const std::vector<TT> &spec, const critical_temperature_params &params = {}, critical_temperature_stats *pst = nullptr)

This algorithm performs temperature-aware SiDB simulation as proposed in “Temperature Behavior of Silicon Dangling Bond Logic” by J. Drewniok, M. Walter, and R. Wille in IEEE NANO 2023 (https://ieeexplore.ieee.org/document/10231259). It comes in two flavors: gate-based and non-gate based.

For Gate-based Critical Temperature Simulation, the Critical Temperature is defined as follows: The temperature at which the erroneous charge distributions are populated by more than \(1 - \eta\), where \(\eta \in [0,1]\).

Template Parameters:
  • Lyt – SiDB cell-level layout type.

  • TT – The type of the truth table specifying the gate behavior.

Parameters:
  • lyt – The layout to simulate.

  • spec – Expected Boolean function of the layout given as a multi-output truth table.

  • params – Simulation and physical parameters.

  • pst – Statistics.

Returns:

The critical temperature (unit: K).

template<typename Lyt>
double fiction::critical_temperature_non_gate_based(const Lyt &lyt, const critical_temperature_params &params = {}, critical_temperature_stats *pst = nullptr)

For Non-gate-based Critical Temperature simulation, the Critical Temperature is defined as follows: The temperature at which the excited charge distributions are populated by more than \(1 - \eta\), where \(\eta \in [0,1]\) is the confidence level for the presence of a working gate.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – The layout to simulate.

  • params – Simulation and physical parameters.

  • pst – Statistics.

Returns:

The critical temperature (unit: K)

Header: fiction/algorithms/simulation/sidb/occupation_probability_of_excited_states.hpp

inline double fiction::occupation_probability_gate_based(const sidb_energy_and_state_type &energy_and_state_type, const double temperature) noexcept

This function computes the occupation probability of erroneous charge distributions (output charge does not match the expected output according the truth table) at a given temperature.

Parameters:
  • energy_and_state_type – This contains the energies of all possible charge distributions together with the information if the charge distribution (state) is transparent or erroneous.

  • temperature – System temperature to assume (unit: K).

Returns:

The occupation probability of all erroneous states is returned.

inline double fiction::occupation_probability_non_gate_based(const sidb_energy_distribution &energy_distribution, const double temperature) noexcept

This function computes the occupation probability of excited states (charge distributions with energy higher than the ground state) at a given temperature.

Parameters:
  • energy_distribution – This contains the energies in eV of all possible charge distributions with the degeneracy.

  • temperature – System temperature to assume (unit: K).

Returns:

The total occupation probability of all excited states is returned.

Header: fiction/algorithms/simulation/sidb/calculate_energy_and_state_type.hpp

using fiction::sidb_energy_and_state_type = std::vector<std::pair<double, bool>>

Data type to collect electrostatic potential energies (in eV) of charge distributions with corresponding state types (i.e., true = transparent, false = erroneous).

template<typename Lyt, typename TT>
sidb_energy_and_state_type fiction::calculate_energy_and_state_type(const sidb_energy_distribution &energy_distribution, const std::vector<charge_distribution_surface<Lyt>> &valid_charge_distributions, const std::vector<bdl_pair<Lyt>> &output_bdl_pairs, const std::vector<TT> &spec, const uint64_t input_index) noexcept

This function takes in an SiDB energy distribution. For each charge distribution, the state type is determined (i.e. erroneous, transparent).

Template Parameters:
  • Lyt – SiDB cell-level layout type.

  • TT – The type of the truth table specifying the gate behavior.

Parameters:
  • energy_distribution – Energy distribution.

  • valid_charge_distributions – Physically valid charge distributions.

  • output_bdl_pairs – Output BDL pairs.

  • spec – Expected Boolean function of the layout given as a multi-output truth table.

  • input_index – The index of the current input configuration.

Returns:

Electrostatic potential energy of all charge distributions with state type.

Maximum Defect Influence Distance

Header: fiction/algorithms/simulation/sidb/maximum_defect_influence_position_and_distance.hpp

struct maximum_defect_influence_distance_params

This struct stores the parameters for the maximum_defect_influence_position_and_distance algorithm.

Public Members

sidb_defect defect = {}

The defect to calculate the maximum defect influence distance for.

sidb_simulation_parameters simulation_parameters = {}

Physical simulation parameters.

std::pair<int32_t, int32_t> additional_scanning_area = {50, 6}

The pair describes the width and height of the area around the gate, which is also used to place defects.

Note

If SiQAD coordinates are used, the second entry describes the number of dimer rows.

template<typename Lyt>
std::pair<typename Lyt::cell, double> fiction::maximum_defect_influence_position_and_distance(const Lyt &lyt, const maximum_defect_influence_distance_params &params = {})

Calculates the maximum distance at which a given defect can influence the layout’s ground state.

This function simulates the influence of defects on a SiDB cell-level layout. It computes the maximum influence distance, defined as the minimum distance between any SiDB cell and the given defect, at which the defect can still affect the layout’s ground state, potentially altering its behavior, such as gate functionality.

Parameters:
  • lyt – The SiDB cell-level layout for which the influence distance is being determined.

  • params – Parameters used to calculate the defect’s maximum influence distance.

Returns:

Pair with the first element describing the position with maximum distance to the layout where a placed defect can still affect the ground state of the layout. The second entry describes the distance of the defect from the layout.

Time-to-Solution (TTS) Statistics

Header: fiction/algorithms/simulation/sidb/time_to_solution.hpp

struct time_to_solution_params

Public Members

exhaustive_sidb_simulation_engine engine = exhaustive_sidb_simulation_engine::QUICKEXACT

Exhaustive simulation algorithm used to simulate the ground state as reference.

uint64_t repetitions = 100

Number of iterations of the heuristic algorithm used to determine the simulation accuracy (repetitions = 100 means that accuracy is precise to 1%).

double confidence_level = 0.997

Confidence level.

struct time_to_solution_stats

This struct stores the time-to-solution, the simulation accuracy and the average single simulation runtime of QuickSim, the single runtime of the exact simulator used, and the number of valid charge configurations found by the exact algorithm.

Public Functions

inline void report(std::ostream &out = std::cout)

Print the results to the given output stream.

Parameters:

out – Output stream.

Public Members

double time_to_solution = {0}

Time-to-solution in seconds.

double acc = {}

Accuracy of the simulation.

double mean_single_runtime = {}

Average single simulation runtime in seconds.

double single_runtime_exhaustive = {}

Single simulation runtime of the exhaustive ground state searcher in seconds.

std::string algorithm

Exhaustive simulation algorithm used to simulate the ground state as reference.

template<typename Lyt>
void fiction::time_to_solution(Lyt &lyt, const quicksim_params &quicksim_params, const time_to_solution_params &tts_params = {}, time_to_solution_stats *ps = nullptr) noexcept

This function determines the time-to-solution (TTS) and the accuracy (acc) of the QuickSim algorithm.

Template Parameters:

Lyt – Cell-level layout type.

Parameters:
  • lyt – Layout that is used for the simulation.

  • quicksim_params – Parameters required for the QuickSim algorithm.

  • tts_params – Parameters used for the time-to-solution calculation.

  • ps – Pointer to a struct where the results (time_to_solution, acc, single runtime) are stored.

Random SiDB Layout Generator

Header: fiction/algorithms/simulation/sidb/random_sidb_layout_generator.hpp

template<typename CoordinateType>
struct generate_random_sidb_layout_params

This struct stores the parameters for the generate_random_sidb_layout algorithm.

Public Types

enum class positive_charges

An enumeration of modes to use for the generation of random SiDB layouts to control control the appearance of positive charges.

Values:

enumerator ALLOWED

Positive charges can occur (i.e. SiDBs can be placed right next to each other).

enumerator FORBIDDEN

Positive charges are not allowed to occur (i.e. SiDBs need to be seperated by a few lattice points).

Public Members

std::pair<CoordinateType, CoordinateType> coordinate_pair

Two coordinates that span the region where SiDBs may be placed (order is not important). The first coordinate is the upper left corner and the second coordinate is the lower right corner of the area.

uint64_t number_of_sidbs = 0

Number of SiDBs that are placed on the layout.

positive_charges positive_sidbs = positive_charges::ALLOWED

If positively charged SiDBs should be prevented, SiDBs are not placed closer than the minimal_spacing.

double minimal_spacing = 2.0

If positively charged SiDBs should be prevented, SiDBs are not placed closer than this value (Euclidean distance of two cells).

uint64_t maximal_attempts = 10E6

Maximum number of steps to place the specified number of SiDBs. Example: If the area, where SiDBs can be placed, is small and many SiDBs are to be placed, several tries are required to generate a layout with no positively charged SiDBs.

uint64_t number_of_unique_generated_layouts = 1

The desired number of unique layouts to be generated.

uint64_t maximal_attempts_for_multiple_layouts = 10E6

The maximum number of attempts allowed to generate the given number of unique layouts (default: \(10^{6}\)). Example: If the area, where SiDBs can be placed, is small and many SiDBs are to be placed, it may be difficult or even impossible to find several unique (given by number_of_unique_generated_layouts) layouts. Therefore, this parameter sets a limit for the maximum number of tries.

template<typename Lyt>
Lyt fiction::generate_random_sidb_layout(const Lyt &lyt_skeleton, const generate_random_sidb_layout_params<coordinate<Lyt>> &params)

Generates a random layout of SiDBs by adding them to the provided layout skeleton. The layout skeleton serves as the starting layout to which SiDBs are added to create the final layout.

Template Parameters:

Lyt – Cell-level SiDB layout type.

Parameters:
  • lyt_skeleton – A layout to which random cells are added to create the final layout.

  • params – The parameters for generating the random layout.

Returns:

A randomly-generated layout of SiDBs.

template<typename Lyt>
std::vector<Lyt> fiction::generate_multiple_random_sidb_layouts(const Lyt &lyt_skeleton, const generate_random_sidb_layout_params<coordinate<Lyt>> &params)

Generates multiple unique random SiDB layouts by adding them to the provided layout skeleton. The layout skeleton serves as the starting layout to which SiDBs are added to create unique SiDB layouts.

Template Parameters:

Lyt – Cell-level SiDB layout type.

Parameters:
  • lyt_skeleton – A layout to which random SiDBs are added to create unique layouts.

  • params – The parameters for generating the random SiDB layouts.

Returns:

A vector containing the unique randomly generated SiDB layouts.

Operational Domain Computation

Header: fiction/algorithms/simulation/sidb/is_operational.hpp

enum class fiction::operational_status

Possible operational status of a layout.

Values:

enumerator OPERATIONAL

The layout is operational.

enumerator NON_OPERATIONAL

The layout is non-operational.

struct is_operational_params

Parameters for the is_operational algorithm.

Public Members

sidb_simulation_parameters simulation_parameters = {}

The simulation parameters for the physical simulation of the ground state.

sidb_simulation_engine sim_engine = {sidb_simulation_engine::QUICKEXACT}

The simulation engine to be used for the operational domain computation.

detect_bdl_pairs_params bdl_params = {}

Parameters for the BDL pair detection algorithms.

template<typename Lyt, typename TT>
std::pair<operational_status, std::size_t> fiction::is_operational(const Lyt &lyt, const std::vector<TT> &spec, const is_operational_params &params = {}) noexcept

Determine the operational status of an SiDB layout.

This function checks the operational status of a given gate layout using the is_operational algorithm. It determines whether the gate layout is operational and returns the correct result for all \(2^n\) input combinations.

Template Parameters:
  • Lyt – SiDB cell-level layout type.

  • TT – The type of the truth table specifying the layout behavior.

Parameters:
  • lyt – The SiDB cell-level layout to be checked.

  • spec – Expected Boolean function of the layout given as a multi-output truth table.

  • params – Parameters for the is_operational algorithm.

Returns:

A pair containing the operational status of the gate layout (either OPERATIONAL or NON_OPERATIONAL) and the number of input combinations tested.

template<typename Lyt, typename TT>
std::set<uint64_t> fiction::operational_input_patterns(const Lyt &lyt, const std::vector<TT> &spec, const is_operational_params &params = {}) noexcept

This function determines the input combinations for which the SiDB-based logic, represented by the provided layout (lyt) and truth table specifications (spec), produces the correct output.

Template Parameters:
  • Lyt – Type of the cell-level layout.

  • TT – Type of the truth table.

Parameters:
  • lyt – The SiDB layout.

  • spec – Vector of truth table specifications.

  • params – Parameters to simualte if a input combination is operational.

Returns:

The count of operational input combinations.

Header: fiction/algorithms/simulation/sidb/operational_domain.hpp

struct operational_domain

An operational domain is a set of simulation parameter values for which a given SiDB layout is logically operational. This means that a layout is deemed operational if the layout’s ground state corresponds with a given Boolean function at the layout’s outputs for all possible input combinations. In this implementation, \(n\) BDL input wires and a single BDL output wire are assumed for a given layout. Any operational domain computation algorithm toggles through all \(2^n\) input combinations and evaluates the layout’s output behavior in accordance with the given Boolean function. The layout is only considered operational for a certain parameter combination, if the output behavior is correct for all input combinations. The operational domain can be computed by sweeping over specified simulation parameters and checking the operational status of the layout for each parameter combination. The operational domain is then defined as the set of all parameter combinations for which the layout is operational. Different techniques for performing these sweep are implemented.

Public Types

enum class sweep_parameter

Possible sweep parameters for the operational domain computation.

Values:

enumerator EPSILON_R

The relative permittivity of the dielectric material.

enumerator LAMBDA_TF

The Thomas-Fermi screening length.

enumerator MU_MINUS

The energy transition level.

Public Members

sweep_parameter x_dimension = {operational_domain::sweep_parameter::EPSILON_R}

X dimension sweep parameter.

sweep_parameter y_dimension = {operational_domain::sweep_parameter::LAMBDA_TF}

Y dimension sweep parameter.

locked_parallel_flat_hash_map<parameter_point, operational_status> operational_values = {}

The operational status of the layout for each specified parameter combination. This constitutes the operational domain. The key of the map is the parameter point, which holds the parameter values in the x and y dimension. The operational status is stored as the value of the map.

struct parameter_point

The parameter point holds parameter values in the x and y dimension.

Public Functions

parameter_point() = default

Standard default constructor.

inline parameter_point(const double x_val, const double y_val)

Standard constructor.

Parameters:
  • x_val – X dimension parameter value.

  • y_val – Y dimension parameter value.

inline bool operator==(const parameter_point &other) const noexcept

Equality operator.

Parameters:

other – Other parameter point to compare with.

Returns:

true iff the parameter points are equal.

inline bool operator!=(const parameter_point &other) const noexcept

Inequality operator.

Parameters:

other – Other parameter point to compare with.

Returns:

true iff the parameter points are not equal.

template<std::size_t I>
inline auto get() const noexcept

Support for structured bindings.

Template Parameters:

I – Index of the parameter value to be returned.

Returns:

The parameter value at the specified index.

Public Members

double x

X dimension parameter value.

double y

Y dimension parameter value.

struct operational_domain_params

Parameters for the operational domain computation. The parameters are used across the different operational domain computation algorithms.

Public Members

sidb_simulation_parameters simulation_parameters = {}

The simulation parameters for the operational domain computation. Most parameters will be kept constant across sweeps, but the sweep parameters are adjusted in each simulation step and thus overwritten in this object.

sidb_simulation_engine sim_engine = {sidb_simulation_engine::QUICKEXACT}

The simulation engine to be used for the operational domain computation.

operational_domain::sweep_parameter x_dimension = {operational_domain::sweep_parameter::EPSILON_R}

The sweep parameter for the x dimension.

double x_min = {1.0}

The minimum value of the x dimension sweep.

double x_max = {10.0}

The maximum value of the x dimension sweep.

double x_step = {0.1}

The step size of the x dimension sweep.

operational_domain::sweep_parameter y_dimension = {operational_domain::sweep_parameter::LAMBDA_TF}

The sweep parameter for the y dimension.

double y_min = {1.0}

The minimum value of the y dimension sweep.

double y_max = {10.0}

The maximum value of the y dimension sweep.

double y_step = {0.1}

The step size of the y dimension sweep.

detect_bdl_pairs_params bdl_params = {}

The parameters for the BDL pair detection, which is necessary during the operational domain computation to detect input and output BDL pairs.

struct operational_domain_stats

Statistics for the operational domain computation. The statistics are used across the different operational domain computation algorithms.

Public Members

mockturtle::stopwatch::duration time_total = {0}

The total runtime of the operational domain computation.

std::size_t num_simulator_invocations = {0}

Number of simulator invocations.

std::size_t num_evaluated_parameter_combinations = {0}

Number of evaluated parameter combinations.

std::size_t num_operational_parameter_combinations = {0}

Number of parameter combinations, for which the layout is operational.

std::size_t num_non_operational_parameter_combinations = {0}

Number of parameter combinations, for which the layout is non-operational.

template<typename Lyt, typename TT>
operational_domain fiction::operational_domain_grid_search(const Lyt &lyt, const std::vector<TT> &spec, const operational_domain_params &params = {}, operational_domain_stats *stats = nullptr)

Computes the operational domain of the given SiDB cell-level layout. The operational domain is the set of all parameter combinations for which the layout is logically operational. Logical operation is defined as the layout implementing the given truth table. The input BDL pairs of the layout are assumed to be in the same order as the inputs of the truth table.

This algorithm uses a grid search to find the operational domain. The grid search is performed by exhaustively sweeping the parameter space in the x and y dimensions. Since grid search is exhaustive, the algorithm is guaranteed to find the operational domain, if it exists within the parameter range. However, the algorithm performs a quadratic number of operational checks on the layout, where each operational check consists of up to \(2^n\) exact ground state simulations, where \(n\) is the number of inputs of the layout. Each exact ground state simulation has exponential complexity in of itself. Therefore, the algorithm is only feasible for small layouts with few inputs.

Template Parameters:
  • Lyt – SiDB cell-level layout type.

  • TT – Truth table type.

Parameters:
  • lyt – Layout to compute the operational domain for.

  • spec – Expected vector of truth tables of the layout. Each truth table represents an output of the Boolean function.

  • params – Operational domain computation parameters.

  • stats – Operational domain computation statistics.

Returns:

The operational domain of the layout.

template<typename Lyt, typename TT>
operational_domain fiction::operational_domain_random_sampling(const Lyt &lyt, const std::vector<TT> &spec, const std::size_t samples, const operational_domain_params &params = {}, operational_domain_stats *stats = nullptr)

Computes the operational domain of the given SiDB cell-level layout. The operational domain is the set of all parameter combinations for which the layout is logically operational. Logical operation is defined as the layout implementing the given truth table. The input BDL pairs of the layout are assumed to be in the same order as the inputs of the truth table.

This algorithm uses random sampling to find a part of the operational domain that might not be complete. It performs a total of samples uniformly-distributed random samples within the parameter range. For each sample, the algorithm performs one operational check on the layout, where each operational check consists of up to \(2^n\) exact ground state simulations, where \(n\) is the number of inputs of the layout. Each exact ground state simulation has exponential complexity in of itself. Therefore, the algorithm is only feasible for small layouts with few inputs.

Template Parameters:
  • Lyt – SiDB cell-level layout type.

  • TT – Truth table type.

Parameters:
  • lyt – Layout to compute the operational domain for.

  • spec – Expected Boolean function of the layout given as a multi-output truth table.

  • samples – Number of samples to perform.

  • params – Operational domain computation parameters.

  • stats – Operational domain computation statistics.

Returns:

The (partial) operational domain of the layout.

template<typename Lyt, typename TT>
operational_domain fiction::operational_domain_flood_fill(const Lyt &lyt, const std::vector<TT> &spec, const std::size_t samples, const operational_domain_params &params = {}, operational_domain_stats *stats = nullptr)

Computes the operational domain of the given SiDB cell-level layout. The operational domain is the set of all parameter combinations for which the layout is logically operational. Logical operation is defined as the layout implementing the given truth table. The input BDL pairs of the layout are assumed to be in the same order as the inputs of the truth table.

This algorithm first uses random sampling to find several operational points within the parameter range. From there, it employs the “flood fill” algorithm to explore the operational domain. The algorithm is guaranteed to find all operational areas in their entirety if the initial random sampling found at least one operational point within them. Thereby, this algorithm works for disconnected operational domains.

It performs samples uniformly-distributed random samples within the parameter range. From there, it performs another number of samples equal to the number of points within the operational domain plus the first non-operational point in each direction. For each sample, the algorithm performs one operational check on the layout, where each operational check consists of up to \(2^n\) exact ground state simulations, where \(n\) is the number of inputs of the layout. Each exact ground state simulation has exponential complexity in of itself. Therefore, the algorithm is only feasible for small layouts with few inputs.

This flavor of operational domain computation was proposed in “Reducing the Complexity of Operational Domain Computation in Silicon Dangling Bond Logic” by M. Walter, J. Drewniok, S. S. H. Ng, K. Walus, and R. Wille in NANOARCH 2023.

Template Parameters:
  • Lyt – SiDB cell-level layout type.

  • TT – Truth table type.

Parameters:
  • lyt – Layout to compute the operational domain for.

  • spec – Expected Boolean function of the layout given as a multi-output truth table.

  • samples – Number of samples to perform.

  • params – Operational domain computation parameters.

  • stats – Operational domain computation statistics.

Returns:

The (partial) operational domain of the layout.

template<typename Lyt, typename TT>
operational_domain fiction::operational_domain_contour_tracing(const Lyt &lyt, const std::vector<TT> &spec, const std::size_t samples, const operational_domain_params &params = {}, operational_domain_stats *stats = nullptr)

Computes the operational domain of the given SiDB cell-level layout. The operational domain is the set of all parameter combinations for which the layout is logically operational. Logical operation is defined as the layout implementing the given truth table. The input BDL pairs of the layout are assumed to be in the same order as the inputs of the truth table.

This algorithm first uses random sampling to find a single operational point within the parameter range. From there, it traverses outwards to find the edge of the operational area and performs Moore neighborhood contour tracing to explore the contour of the operational domain. If the operational domain is connected, the algorithm is guaranteed to find the contours of the entire operational domain within the parameter range if the initial random sampling found an operational point.

It performs up to samples uniformly-distributed random samples within the parameter range until an operational point is found. From there, it performs another number of samples equal to the distance to an edge of the operational area. Finally, it performs up to 8 samples for each contour point (however, the actual number is usually much lower). For each sample, the algorithm performs one operational check on the layout, where each operational check consists of up to \(2^n\) exact ground state simulations, where \(n\) is the number of inputs of the layout. Each exact ground state simulation has exponential complexity in of itself. Therefore, the algorithm is only feasible for small layouts with few inputs.

This flavor of operational domain computation was proposed in “Reducing the Complexity of Operational Domain Computation in Silicon Dangling Bond Logic” by M. Walter, J. Drewniok, S. S. H. Ng, K. Walus, and R. Wille in NANOARCH 2023.

Template Parameters:
  • Lyt – SiDB cell-level layout type.

  • TT – Truth table type.

Parameters:
  • lyt – Layout to compute the operational domain for.

  • spec – Expected Boolean function of the layout given as a multi-output truth table.

  • samples – Number of samples to perform.

  • params – Operational domain computation parameters.

  • stats – Operational domain computation statistics.

Returns:

The (partial) operational domain of the layout.

Utility Functions

Simulation Equivalence Checking

Header: fiction/algorithms/simulation/sidb/check_simulation_results_for_equivalence.hpp

template<typename Lyt>
bool fiction::check_simulation_results_for_equivalence(sidb_simulation_result<Lyt> result1, sidb_simulation_result<Lyt> result2)

This function compares two SiDB simulation results for equivalence. Two results are considered equivalent if they have the same number of charge distributions and if each corresponding charge distribution has the same electrostatic potential energy and charge states for all cells.

Template Parameters:

Lyt – The SiDB cell-level layout type used in the simulation results.

Parameters:
  • result1 – The first SiDB simulation result to compare.

  • result2 – The second SiDB simulation result to compare.

Returns:

true if the two simulation results are equivalent, false otherwise.

Determine the Ground State from Simulation Results

Header: fiction/algorithms/simulation/sidb/determine_groundstate_from_simulation_results.hpp

template<typename Lyt>
std::vector<charge_distribution_surface<Lyt>> fiction::determine_groundstate_from_simulation_results(const sidb_simulation_result<Lyt> &simulation_results) noexcept

This function calculates the ground state charge distributions from the provided simulation results. The ground state charge distributions are those with energy closest to the minimum energy found in the simulation results.

Note

When degenerate states exist, there are multiple ground states with the same energy.

Template Parameters:

Lyt – The layout type used in the simulation results.

Parameters:

simulation_results – The simulation results containing charge distributions.

Returns:

A vector of charge distributions with the minimal energy.

Charge Detection

Header: fiction/algorithms/simulation/sidb/can_positive_charges_occur.hpp

template<typename Lyt>
bool fiction::can_positive_charges_occur(const Lyt &lyt, const sidb_simulation_parameters &sim_params) noexcept

This algorithm determines if positively charged SiDBs can occur in a given SiDB cell-level layout due to strong electrostatic interaction.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – The layout to be analyzed.

  • sim_params – Physical parameters used to determine whether positively charged SiDBs can occur.

Binary-dot Logic (BDL) Pair Detection

Header: fiction/algorithms/simulation/sidb/detect_bdl_pairs.hpp

template<typename Lyt>
struct bdl_pair

A Binary-dot Logic (BDL) pair is a pair of SiDBs that are close to each other and, thus, most likely share a charge.

Template Parameters:

Lyt – SiDB cell-level layout type.

Public Functions

bdl_pair() = default

Standard constructor for empty BDL pairs.

inline bdl_pair(const sidb_technology::cell_type t, const cell<Lyt> &u, const cell<Lyt> &l) noexcept

Constructor for BDL pairs.

Parameters:
  • t – Type of the SiDBs in the pair.

  • u – The upper SiDB of the pair.

  • l – The lower SiDB of the pair.

Public Members

const sidb_technology::cell_type type = {}

The type of the SiDBs in the pair. BDL SiDBs must be of the same type. They can either be normal, input, or output SiDBs.

const cell<Lyt> upper = {}

The upper SiDB of the pair. Upper and lower are defined relative to each other via the operator< overload.

const cell<Lyt> lower = {}

The lower SiDB of the pair. Upper and lower are defined relative to each other via the operator< overload.

struct detect_bdl_pairs_params

Parameters for the BDL pair detection algorithms.

Public Members

double minimum_distance = {0.75}

The minimum distance between two dots to be considered a BDL pair. This is useful to prevent, e.g., SiDBs of atomic wires to be considered BDL pairs. (unit: nm).

double maximum_distance = {1.5}

The maximum distance between two dots to be considered a BDL pair. This is useful to prevent unlikely pairings of SiDBs that are far apart and to improve performance of the matching algorithm. (unit: nm).

template<typename Lyt>
std::vector<bdl_pair<Lyt>> fiction::detect_bdl_pairs(const Lyt &lyt, const typename technology<Lyt>::cell_type type, const detect_bdl_pairs_params &params = {}) noexcept

This algorithm detects BDL pairs in an SiDB layout. It does so by first collecting all dots of the given type and then uniquely pairing them up based on their distance. Lower and upper distance thresholds can be defined (defaults = 0.75 nm and 1.5 nm, respectively) to narrow down the range in which SiDBs could be considered a BDL pair. The distance between two dots is computed using the sidb_nm_distance function. The algorithm returns a vector of BDL pairs.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – The layout to detect BDL pairs in.

  • type – The type of the SiDBs to detect BDL pairs for, e.g., INPUT, OUTPUT, NORMAL.

  • params – Parameters for the BDL pair detection algorithm.

Returns:

A vector of BDL pairs.

Assess Population Stability

Header: fiction/algorithms/simulation/sidb/assess_physical_population_stability.hpp

enum class fiction::transition_type

Possible types of charge transitions that can occur in an SiDB layout. These transitions represent changes in the charge state of SiDBs, including transitions from neutral to negative, negative to neutral, neutral to positive, and positive to neutral.

Values:

enumerator NEUTRAL_TO_NEGATIVE

SiDB is neutrally charged, but is closest to being negatively charged.

enumerator NEGATIVE_TO_NEUTRAL

SiDB is negatively charged, but is closest to being neutrally charged.

enumerator NEUTRAL_TO_POSITIVE

SiDB is neutrally charged, but is closest to being positively charged.

enumerator POSITIVE_TO_NEUTRAL

SiDB is positively charged, but is closest to being neutrally charged.

template<typename Lyt>
struct population_stability_information

This struct encapsulates information related to the population stability of a charge distribution. It includes details about the SiDB closest to a charge transition (critical cell), the specific charge state transition, the electrostatic potential difference required for the transition, the corresponding distance, and the total electrostatic energy of the given charge distribution.

Template Parameters:

Lyt – SiDB cell-level layout type.

Public Members

Lyt::cell critical_cell = {}

SiDB cell which is closest to a charge transition.

transition_type transition_from_to = {}

Charge transition from the current charge state to the closest one.

double minimum_potential_difference_to_transition = {}

Absolute electrostatic potential (unit: V) required for the charge state transition.

double distance_corresponding_to_potential = {}

Distance (unit: nm) corresponding to the minimum potential difference.

double system_energy = {}

Total electrostatic energy (unit: eV) of given charge distribution.

struct assess_physical_population_stability_params

This struct stores the parameters required to assess the population stability.

Public Members

sidb_simulation_parameters simulation_parameters = {}

Parameters for the electrostatic potential.

uint64_t precision_for_distance_corresponding_to_potential = 2

The precision level for the conversion from the minimum potential difference to the corresponding distance.

template<typename Lyt>
std::vector<population_stability_information<Lyt>> fiction::assess_physical_population_stability(const Lyt &lyt, const assess_physical_population_stability_params &params) noexcept

This function assesses the population stability of each physically valid charge distributions of a given SiDB layout. It determines the minimum absolute electrostatic potential required to induce a charge distribution transition. The function also identifies the SiDB for which this is the case (critical SiDB) and the corresponding charge state transition (i.e., the change from one charge state to another).

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – The layout for which the population stability is assessed.

  • params – Parameters used to assess the population stability.

Returns:

A vector of population stability information for all physically valid charge distributions of the given SiDB layout.

Convert Potential to Distance

Header: fiction/algorithms/simulation/sidb/convert_potential_to_distance.hpp

inline double fiction::convert_potential_to_distance(const double potential, const sidb_simulation_parameters &params = sidb_simulation_parameters{}, const uint64_t precision = 2) noexcept

The electrostatic potential on hydrogen-passivated silicon is typically modeled using a screened Coulomb potential. This electrostatic potential is commonly employed to determine the electrostatic potential for a given distance (between SiDB and point under consideration) and given physical parameters. However, the function provided here serves the inverse purpose by calculating the distance for a given potential and given physical parameters.

Note

Runtime depends exponentially on the provided precision.

Parameters:
  • params – The physical parameters for a given hydrogen-passivated silicon surface.

  • potential – The electrostatic potential (unit: V) to be converted to a distance.

  • precision – The precision level for the conversion, specifying the number of decimal places.

Returns:

The distance (unit: nm) corresponding to the given electrostatic potential.