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 = {5.6}

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

double lambda_tf = {5.0}

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

double mu_minus = {-0.32}

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

uint8_t base = {3}

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 – SiDB 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.

inline std::vector<charge_distribution_surface<Lyt>> groundstates() const noexcept

This function computes the ground state of the charge distributions.

Note

If degenerate states exist in the simulation result, this function will return multiple ground states that all possess the same system energy.

Returns:

A vector of charge distributions with the minimal energy.

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 in seconds.

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.

uint64_t timeout = std::numeric_limits<uint64_t>::max()

Timeout limit (in ms).

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

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.

Note

QuickSim currently does not support atomic defect simulation.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – The layout to simulate.

  • ps – QuickSim parameters.

Returns:

sidb_simulation_result is returned if the simulation was successful, otherwise std::nullopt.

Exhaustive Ground State Simulation

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

template<typename CellType = offset::ucoord_t>
struct quickexact_params

This struct stores the parameters for the QuickExact algorithm.

Template Parameters:

CellType – Cell type.

Public Types

enum class automatic_base_number_detection : uint8_t

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/clustercomplete.hpp

template<typename CellType = offset::ucoord_t>
struct clustercomplete_params

The struct containing the parameters both passed on to pre-simulator Ground State Space, and used during simulation.

Public Types

enum class ground_state_space_reporting : uint8_t

This enum class provides meaningful options for configuring the reporting of the Ground State Space statistics. These statistic may be used especially to configure the validity witness partitioning options for Ground State Space, that may impair runtimes when set too high, but could provide a large benefit to the complexity of the unfolding process of large simulation problems by performing more involved pruning procedures in the construction stage.

Values:

enumerator ON

Enabling this option will output Ground State Space statistics to the standard output.

enumerator OFF

Disabling this option will suppress the output of Ground State Space statistics.

Public Members

sidb_simulation_parameters simulation_parameters = {}

Physical simulation parameters.

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.

uint64_t validity_witness_partitioning_max_cluster_size_gss = 6

This specifies the maximum cluster size for which Ground State Space will solve an NP-complete sub-problem exhaustively. The sets of SiDBs that witness local population stability for each respective charge state may be partitioned into disjoint sets such that the number of required witnesses for each respective charge state is satisfied. If no such partition exists, the multiset charge configuration associated with the requirements may be rejected.

uint64_t num_overlapping_witnesses_limit_gss = 6

The complexity is of validity witness partitioning bounded by a factorial in the number of overlapping witnesses. This parameter thus allows the validity witness partitioning procedure to perform the reduction to overlapping witnesses for larger cluster sizes that could be runtime-impairing, then limiting specifically the length of the input to the factorial call.

uint64_t available_threads = std::thread::hardware_concurrency()

Number of threads to make available to ClusterComplete for the unfolding stage.

ground_state_space_reporting report_gss_stats = ground_state_space_reporting::OFF

Option to decide if the Ground State Space statistics are reported to the standard output. By default, this option is disabled.

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

ClusterComplete is an instantiation of a general solution to exhaustive state assignment searching for which all local predicates hold, given respective local evaluations that may be aggregated from individual inter-variable interactions. Applied to the problem of exact physical simulation of SiDBs, it is able to efficiently consider positive charges that are rare to occur, but drastically blow up exact simulation runtimes when hierarchical pruning methods are not applied. In fact, the exponential growth in problem complexity for added SiDBs is tamed by ClusterComplete, as SiDB layouts to simulate in practise amount to a high pruning efficacy, resulting in a layout-dependent reduction of the simulation base. This amounts to an effective simulation base in the real number range \([1,b]\), where \(b\in\{2,3\}\) is the given simulation base.

The part of the ClusterComplete algorithm that is implemented in this file is the destructive phase of the procedure that employs the duality of construction and destruction, folding and unfolding. The phase preceding it is the key ingredient to the achieved efficiency: the Ground State Space algorithm, which constructs a minimized hierarchical search space of charge configurations that adhere to the critical population stability criterion. In particular, it generalizes physically informed space pruning that contributes to the capabilities of the QuickExact simulator, now applying to all charge states equally, and, most importantly, it lifts the associated potential equations to higher order, allowing us to reason over potential bounds in a cluster hierarchy.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – Layout to simulate.

  • params – Parameter required for both the invocation of Ground State Space, and the simulation following.

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 and ClusterComplete due to the much better runtimes and more functionality.

Template Parameters:

Lyt – SiDB 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 : uint8_t

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.

enumerator CLUSTERCOMPLETE

ClusterComplete is a novel exact simulation engine that requires exponential runtime, though, depending on the simulation problem, it effectively reduces the base number by a real number, thus allowing problem sizes that were previously considered astronomical in size. Inherent to the simulation methodology that does not depend on the simulation base, it simulates very effectively for either base number (2 or 3).

enum class fiction::exact_sidb_simulation_engine : uint8_t

Selector exclusively for exact 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.

enumerator CLUSTERCOMPLETE

ClusterComplete is a novel exact simulation engine that requires exponential runtime, though, depending on the simulation problem, it effectively reduces the base number by a real number, thus allowing problem sizes that were previously considered astronomical. Inherent to the simulation methodology that does not depend on the simulation base, it simulates very effectively for either base number (2 or 3).

enum class fiction::heuristic_sidb_simulation_engine : uint8_t

Selector exclusively for heuristic SiDB simulation engines.

Values:

enumerator QUICKSIM

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

template<typename EngineType>
std::string_view fiction::sidb_simulation_engine_name(const EngineType &engine) noexcept

Returns the name of the given simulation engine.

Template Parameters:

EngineType – The type of the SiDB simulation engine (exhaustive/heuristic/generic).

Parameters:

engine – An SiDB simulation engine.

Returns:

The name of the simulation engine.

Energy Calculation

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

struct energy_state

This struct stores the energy state of an SiDB layout. The energy state consists of the electrostatic potential energy and the degeneracy of the state.

Public Functions

inline energy_state(const double electrostatic_potential_energy, const uint64_t degeneracy)

Default constructor.

Parameters:
  • electrostatic_potential_energy – The electrostatic potential energy of the charge distribution (eV).

  • degeneracy – The degeneracy of the state.

Public Members

double electrostatic_potential_energy

The electrostatic potential energy of the charge distribution (eV).

uint64_t degeneracy

The degeneracy of the state.

class energy_distribution

This class is used to store the energy distribution of an SiDB layout. The energy distribution is a map that contains the electrostatic potential as a key and its degeneracy as a value. To be more precise, if two different charge distributions occur with the same energy, the degeneracy value of the energy state is 2.

Public Functions

energy_distribution() = default

Default constructor.

inline std::optional<energy_state> get_nth_state(const uint64_t state_index) const noexcept

Returns the nth state (energy + degeneracy) in the energy distribution.

Parameters:

state_index – The index of the state to be retrieved.

Returns:

Energy state. If the index is out of range, std::nullopt is returned instead.

inline std::optional<uint64_t> degeneracy(const double energy) const noexcept

Returns the degeneracy value (number of states) with the given energy value.

Parameters:

energy – The energy value for which the excited state number is to be determined.

Returns:

The degeneracy of the given energy. If the energy value is not found, std::nullopt is returned instead.

inline void add_energy_state(const energy_state &state) noexcept

Adds a state to the energy distribution.

Parameters:

state – The energy state to be added.

inline std::size_t size() const noexcept

Returns the number of energy states in the energy distribution.

Returns:

The number of energy states in the energy distribution.

inline bool empty() const noexcept

Checks if the energy distribution is empty.

Returns:

true if the energy distribution is empty, false otherwise.

inline double max_energy() const noexcept

Returns the maximum energy value in the energy distribution.

Returns:

The maximum energy value in the energy distribution.

inline double min_energy() const noexcept

Returns the minimum energy value in the energy distribution.

Returns:

The minimum energy value in the energy distribution.

template<typename Fn>
inline void for_each(Fn &&fn) const

Applies a function to all states in the energy distribution.

Template Parameters:

Fn – Functor type.

Parameters:

fn – Functor to apply to each key-value pair.

template<typename Lyt>
energy_distribution fiction::calculate_energy_distribution(const std::vector<charge_distribution_surface<Lyt>> &charge_distributions)

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. To compare two energy values for equality, the comparison uses a tolerance specified by constants::ERROR_MARGIN.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:

charge_distributions – A vector of charge_distribution_surface objects for which the energy distribution is computed.

Returns:

Energy distribution.

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> &exact_results) noexcept

This function checks if the elstrostatic ground state of an SiDB layout is found by a heuristic for the physical simulation (e.g., QuickSim or SimAnneal).

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • heuristic_results – Simulation results obtained from a heuristic physical simulation.

  • exact_results – Simulation results obtained from an exact physical simulation.

Returns:

Returns true if the ground state is contained in the simulation result provided by the heuristic physical simulation. 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 Members

is_operational_params operational_params = {}

The parameters used to determine if a layout is operational or non-operational.

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).

uint64_t iteration_steps = {80}

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

double alpha = {0.7}

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

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 – Type of the truth table.

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::calculate_boltzmann_factor(const double energy, const double min_energy, const double temperature) noexcept

This function computes the Boltzmann factor for a given energy, the minimal energy, and the temperature.

Parameters:
  • energy – Boltzmann factor for the given energy is calculated.

  • min_energy – Minimum energy of the system.

  • temperature – Temperature of the system (unit: K).

Returns:

Boltzmann factor.

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 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

enum class fiction::state_type

Label to categorize ground and excited states of an SiDB layout.

Values:

enumerator ACCEPTED

A state is accepted if the charge distribution encodes the desired logic.

enumerator REJECTED

A state is rejected if the charge distributiion does not encode the desired logic. Moreover, if kinks are rejected, a charge distribution that encodes the logic, but does show kinks, is rejected.

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

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_with_kinks_accepted(const energy_distribution &energy_distribution, const std::vector<charge_distribution_surface<Lyt>> &valid_charge_distributions, const std::vector<bdl_pair<cell<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) while kinks are accepted, meaning a state with kinks is considered 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.

template<typename Lyt, typename TT>
sidb_energy_and_state_type fiction::calculate_energy_and_state_type_with_kinks_rejected(const energy_distribution &energy_distribution, const std::vector<charge_distribution_surface<Lyt>> &valid_charge_distributions, const std::vector<TT> &spec, const uint64_t input_index, const std::vector<bdl_wire<Lyt>> &input_bdl_wires, std::vector<bdl_wire<Lyt>> &output_bdl_wires) noexcept

This function takes in an SiDB energy distribution. For each charge distribution, the state type is determined (i.e. erroneous, transparent) while kinks are rejected, meaning a state with kinks is considered erroneous.

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.

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

  • input_index – The index of the current input configuration.

  • input_bdl_wires – Input BDL wires.

  • output_bdl_wires – Output BDL wires.

Returns:

Electrostatic potential energy of all charge distributions with state type.

Ground State Space Construction

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

struct ground_state_space_results

This struct is used to store the results of the Ground State Space construction.

Public Functions

inline void report(std::ostream &os = std::cout) const noexcept

Report Ground State Space statistics. A quick heuristic to assess the quality of the pruning is captured by the size of the charge space of the top cluster, which depends on the charge spaces of all clusters below it.

Parameters:

os – The output stream to write to (default: standard output).

Returns:

Prints the runtime and the number of pruned top level multisets versus the total amount possible.

Public Members

const sidb_cluster_ptr top_cluster = {nullptr}

The top cluster is the root of the cluster hierarchy. It therefore allows access to the entire cluster hierarchy, including the charge spaces of each cluster.

const std::chrono::duration<double> runtime = {}

The runtime of the construction is stored.

const uint64_t maximum_top_level_multisets = {}

The maximum size of the charge space of the top cluster, given the simulation base, can be inferred by the “stars and bars” combinatorial idea: the solution to this analogous problem determines the maximum amount of multisets of size \(N\) (where \(N\) is the number of SiDBs in the layout, and therefore in the top cluster) for the given base \(b\). In particular, the analogy is as follows: any such multiset can be seen as \(N\) stars and \(b - 1\) bars separating those stars. Then, the \(b\) partitions forming from these \(b - 1\) separators each have a respective size, adding up to \(N\). Therefore each partition is associated with an amount of one of the charge states of the multiset. Now we may compute total number of possible multisets for the top cluster as the number of combinations of \(N\) stars and \(b - 1\) bars. Hence this is computed with the following combinatorial formula: \(\binom{N + b - 1}{b - 1}\).

const uint64_t projector_state_count = {}

The total number of distinct projector states is counted. At each merge, the projector states in charge space compositions in the charge spaces of the clusters to merge are locked in the final construct, and can therefore be counted. This may be used to estimate the time it would take ClusterComplete to unfold the hierarchy.

template<typename Lyt>
ground_state_space_results fiction::ground_state_space(const Lyt &lyt, const ground_state_space_params &params = {}) noexcept

The purely constructive Ground State Space algorithm is the key ingredient of the ClusterComplete

exact SiDB simulator that lifts exact SiDB simulation to permit multiple gates in connection. It uses iterative “loop until

fixpoint” concepts to prune the simulation search space for not only a flat layout of SiDBs, but rather generalizes, and lifts the physically informed space pruning technique introduced with

QuickExact to higher order, allowing Ground State Space to prune multiset charge state configurations at any level in a cluster hierarchy.

The role of the cluster hierarchy is to rank interactions between groups, or clusters of SiDBs that together make up the whole layout, such that the variation in electrostatic potential under different charge state assignments is highest between the children clusters of clusters that low in the hierarchy. Thereby, the structure allows us to consider the most charge state assignment-dependent interaction in a more detailed physically informed space pruning analysis, enabling high pruning efficacy for the few pruning tests (with respect to the exponential search space).

Starting at a clustering of all singleton clusters, the charge spaces, ie. a set of multiset charge configurations (initially { {{-}}, {{0}}, {{+}} } or omitting the singleton multiset {{+}} in the case of base 2 pre-simulation), are pruned iteratively through potential bound analysis. Through merges, ie., replacing a set of children in the clustering with their parent, we may inspect the most crucially dependant interactions in the layout separately. The procedure finishes when the charge spaces have been folded all the way up to the top cluster, parent of all, which then contains all information resulting from the construction. ClusterComplete, without much trickery, now simply unfolds this result, allowing simulation of problems that were previously seen as astronomical, due to the (base 2 or 3) exponential growth in the number of SiDBs.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – Layout to construct the Ground State Space of.

  • params – The parameters that Ground State Space will use throughout the construction. The physical parameters that Ground State Space will use to prune simulation search space are stored in there. In particular, the user may configure parameters that decide limits on the problem sizes of pruning by validity witness partitioning. By default, these are set to avoid runtimes from being affected, as these sub-problems may scale factorially. Thereby, these parameters are especially useful for large simulation problems that could benefit from extra intensive pruning before ClusterComplete unfolds the constructed hierarchical charge space.

Returns:

The results of the construction, which include the top cluster which parents all other clusters, and thereby contains the charge spaces of each cluster.

Time-to-Solution (TTS) Statistics

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

struct time_to_solution_params

Public Members

exact_sidb_simulation_engine engine = exact_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

The confidence level represents the probability that the confidence interval calculated from the simulation contains the true value. For example, a 99.7 % (0.997) confidence level means that if the simulation were repeated many times, approximately 997 out of 1000 of the calculated confidence intervals would contain the true value.

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) const

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 in %.

double mean_single_runtime = {}

Average single simulation runtime in seconds.

double single_runtime_exact = {}

Single simulation runtime of the exact ground state simulation algorithm.

std::string algorithm

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

template<typename Lyt>
void fiction::time_to_solution(const 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 – SiDB 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.

template<typename Lyt>
void fiction::time_to_solution_for_given_simulation_results(const sidb_simulation_result<Lyt> &results_exact, const std::vector<sidb_simulation_result<Lyt>> &results_heuristic, const double confidence_level = 0.997, time_to_solution_stats *ps = nullptr) noexcept

This function calculates the Time-to-Solution (TTS) by analyzing the simulation results of a heuristic algorithm in comparison to those of an exact algorithm. It provides further statistical metrics, including the accuracy of the heuristic algorithm, and individual runtimes.

Template Parameters:

Lyt – SiDB ell-level layout type.

Parameters:
  • results_exact – Simulation results of the exact algorithm.

  • results_heuristic – Simulation of the heuristic for which the TTS is determined.

  • confidence_level – Confidence level for the TTS computation. The confidence level represents the probability that the confidence interval calculated from the simulation contains the true value. For example, a 95 % (0.95) confidence level means that if the simulation were repeated many times, approximately 95 out of 100 of the calculated confidence intervals would contain the true value.

  • ps – Pointer to a struct where the statistics of this function call (time_to_solution, acc, single runtime) are to be 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 : uint8_t

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).

enumerator MAY_OCCUR

Positive charges can occur, which means that the can_positive_charges_occur function returns true.

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.

sidb_simulation_parameters simulation_parameters = {}

Simulation parameters.

uint64_t maximal_attempts = static_cast<uint64_t>(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 = 1'000'000

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>
std::optional<Lyt> fiction::generate_random_sidb_layout(const generate_random_sidb_layout_params<coordinate<Lyt>> &params, const std::optional<Lyt> &skeleton = std::nullopt) noexcept

Generates a layout featuring a random arrangement of SiDBs. These randomly placed dots can be incorporated into an existing layout skeleton that may be optionally provided.

Template Parameters:

Lyt – SiDB cell-level SiDB layout type.

Parameters:
  • params – The parameters for generating the random layout.

  • skeleton – Optional layout to which random dots are added.

Returns:

A randomly generated SiDB layout, or std::nullopt if the process failed due to conflicting parameters.

template<typename Lyt>
std::optional<std::vector<Lyt>> fiction::generate_multiple_random_sidb_layouts(const generate_random_sidb_layout_params<coordinate<Lyt>> &params, const std::optional<Lyt> &skeleton = std::nullopt) noexcept

Generates multiple random layouts featuring a random arrangement of SiDBs. These randomly placed dots can be incorporated into an existing layout skeleton that may be optionally provided.

Template Parameters:

Lyt – SiDB cell-level SiDB layout type.

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

  • skeleton – Optional layout to which random dots are added.

Returns:

A vector containing the unique randomly generated SiDB layouts. If the design is impossible, std::nullopt

Operational Domain Computation

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

enum class fiction::operational_status : uint8_t

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 Types

enum class operational_condition : uint8_t

Condition to decide whether a layout is operational or non-operational.

Values:

enumerator TOLERATE_KINKS

Even if the I/O pins show kinks, the layout is still considered as operational.

enumerator REJECT_KINKS

The I/O pins are not allowed to show kinks. If kinks exist, the layout is considered non-operational.

enum class operational_analysis_strategy : uint8_t

Simulation method to determine if the layout is operational or non-operational. There are three possible strategies:

  • SIMULATION_ONLY: This setting does not apply any filtering strategies to determine if the layout is operational. Instead, it relies solely on physical simulation to make this determination.

  • FILTER_ONLY: This setting does only apply filtering strategies to determine if the layout is non-operational. If the layout passes all filtering strategies, it is considered operational. This is only an approximation. It may be possible that the layout is non-operational, but the filtering strategies do not detect it.

  • FILTER_THEN_SIMULATION: Before a physical simulation is conducted, the algorithm checks if filtering strategies have detected whether the layout is non-operational. This only provides any runtime benefits if kinks are rejected.

Values:

enumerator SIMULATION_ONLY

Do not apply filter strategies to determine whether the layout is operational. Instead, rely solely on physical simulation.

enumerator FILTER_ONLY

Apply filtering exclusively to determine whether the layout is non-operational. If the layout passes all filter steps, it is considered operational.

Note

This is an extremely fast approximation that may sometimes lead to false positives.

enumerator FILTER_THEN_SIMULATION

Before a physical simulation is conducted, the algorithm checks if filter strategies can determine that the layout is non-operational. This only provides any runtime benefits if kinks are rejected.

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.

bdl_input_iterator_params input_bdl_iterator_params = {}

Parameters for the BDL input iterator.

operational_condition op_condition = operational_condition::TOLERATE_KINKS

Condition to decide whether a layout is operational or non-operational.

operational_analysis_strategy strategy_to_analyze_operational_status = operational_analysis_strategy::SIMULATION_ONLY

Strategy to determine whether a layout is operational or non-operational.

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 SiDB layout using the is_operational algorithm. It determines whether the SiDB layout is operational and returns the correct result for all \(2^n\) input combinations.

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

  • TT – Type of the truth table.

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 SiDB layout (either OPERATIONAL or NON_OPERATIONAL) and the number of input combinations tested.

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, const std::vector<bdl_wire<Lyt>> &input_bdl_wire, const std::vector<bdl_wire<Lyt>> &output_bdl_wire, const std::optional<Lyt> &canvas_lyt = std::nullopt) noexcept

Determine the operational status of an SiDB layout.

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

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

  • TT – Type of the truth table.

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.

  • input_bdl_wire – Optional BDL input wires of lyt.

  • output_bdl_wire – Optional BDL output wires of lyt.

  • canvas_lyt – Optional canvas layout.

Returns:

A pair containing the operational status of the SiDB 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 layout is operational.

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

  • TT – Type of the truth table.

Parameters:
  • lyt – The SiDB layout.

  • spec – Vector of truth table specifications.

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

Returns:

The operational input combinations.

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, const std::vector<bdl_wire<Lyt>> &input_bdl_wire, const std::vector<bdl_wire<Lyt>> &output_bdl_wire, const std::optional<Lyt> &canvas_lyt = std::nullopt) noexcept

This function determines the input combinations for which the layout is operational.

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

  • TT – Type of the truth table.

Parameters:
  • lyt – The SiDB layout.

  • spec – Vector of truth table specifications.

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

  • input_bdl_wire – Optional BDL input wires of lyt.

  • output_bdl_wire – Optional BDL output wires of lyt.

  • canvas_lyt – Optional canvas layout.

Returns:

The count of operational input combinations.

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

This function determines if the layout is only considered non-operational because of kinks. This means that the layout would be considered as operational, if kinks were accepted.

Note

“Kink induced non-operational” refers to the non-operational status being exclusively caused by kinks with an otherwise correct logic match.

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

  • TT – Type of the truth table.

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:

Bool that indicates whether kinks induce the layout to become non-operational. true if the layout is non-operational due to kinks, false otherwise.

template<typename Lyt, typename TT>
bool fiction::is_kink_induced_non_operational(const Lyt &lyt, const std::vector<TT> &spec, const is_operational_params &params, const std::vector<bdl_wire<Lyt>> &input_bdl_wire, const std::vector<bdl_wire<Lyt>> &output_bdl_wire, const std::optional<Lyt> &canvas_lyt = std::nullopt) noexcept

This function determines if the layout is only considered non-operational because of kinks. This means that the layout would be considered as operational, if kinks were accepted.

Note

“Kink induced non-operational” refers to the non-operational status being exclusively caused by kinks with an otherwise correct logic match.

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

  • TT – Type of the truth table.

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.

  • input_bdl_wire – Optional BDL input wires of lyt.

  • output_bdl_wire – Optional BDL output wires of lyt.

  • canvas_lyt – Optional canvas layout.

Returns:

Bool that indicates whether kinks induce the layout to become non-operational. true if the layout is non-operational due to kinks, false otherwise.

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

This function determines all input combinations for which kinks induce the SiDB layout to become non-operational. This means that the layout is operational if kinks would be accepted.

Note

“Kink induced non-operational” refers to the non-operational status being exclusively caused by kinks with an otherwise correct logic match.

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

  • TT – Type of the truth table.

Parameters:
  • lyt – The SiDB layout.

  • spec – Vector of truth table specifications.

  • params – Parameters for the is_operational algorithm.

Returns:

The input combinations where kinks induce the SiDB layout to become non-operational.

template<typename Lyt, typename TT>
std::set<uint64_t> fiction::kink_induced_non_operational_input_patterns(const Lyt &lyt, const std::vector<TT> &spec, const is_operational_params &params, const std::vector<bdl_wire<Lyt>> &input_bdl_wire, const std::vector<bdl_wire<Lyt>> &output_bdl_wire, const std::optional<Lyt> &canvas_lyt = std::nullopt) noexcept

This function determines all input combinations for which kinks induce the SiDB layout to become non-operational. This means that the layout is operational if kinks would be accepted.

Note

“Kink induced non-operational” refers to the non-operational status being exclusively caused by kinks with an otherwise correct logic match.

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

  • TT – Type of the truth table.

Parameters:
  • lyt – The SiDB layout.

  • spec – Vector of truth table specifications.

  • params – Parameters for the is_operational algorithm.

  • input_bdl_wire – Optional BDL input wires of lyt.

  • output_bdl_wire – Optional BDL output wires of lyt.

  • canvas_lyt – Optional canvas layout.

Returns:

The input combinations where kinks induce the SiDB layout to become non-operational.

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

struct parameter_point

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

Public Functions

parameter_point() = default

Default constructor.

inline explicit parameter_point(const std::vector<double> &values)

Standard constructor.

Parameters:

values – Parameter values for each dimension.

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

Equality operator. Checks if this parameter point is equal to another point within a specified tolerance. The tolerance is defined by constants::ERROR_MARGIN.

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 if the parameter points are not equal.

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

Support for structured bindings.

Template Parameters:

I – Index of the parameter value to be returned.

Throws:

std::out_of_range – if the index is out of bounds.

Returns:

The parameter value at the specified index.

inline const std::vector<double> &get_parameters() const noexcept

Returns the parameter values for each dimension.

Returns:

The parameter values for each dimension.

enum class fiction::sweep_parameter : uint8_t

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.

class operational_domain : public fiction::sidb_simulation_domain<parameter_point, operational_status>

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 Functions

operational_domain() = default

Default constructor.

inline explicit operational_domain(const std::vector<sweep_parameter> &dims)

Standard constructor.

Parameters:

dims – Dimensions.

inline void add_dimension(const sweep_parameter &dim)

Adds a dimension to sweep over. The first dimension is the x dimension, the second dimension is the y dimension, etc.

Parameters:

dim – The dimension to add.

inline const sweep_parameter &get_dimension(const std::size_t index) const

Returns a specific dimension by index.

Parameters:

index – The index of the dimension to return.

Throws:

std::out_of_range – if the index is out of range.

Returns:

The dimension at the specified index.

inline std::size_t get_number_of_dimensions() const noexcept

Returns the number of dimensions to sweep over.

Returns:

The number of dimensions to sweep over.

class critical_temperature_domain : public fiction::sidb_simulation_domain<parameter_point, operational_status, double>

The critical_temperature_domain class collects the critical temperature and the operational status for a range of different physical parameters of a given SiDB layout. It allows for the evaluation of how the critical temperature depends on variations in the underlying parameter points. This enables simulations to explore the critical temperature’s behavior across different conditions and configurations.

Public Functions

critical_temperature_domain() = default

Default constructor.

inline explicit critical_temperature_domain(const std::vector<sweep_parameter> &dims)

Standard constructor.

Parameters:

dims – Dimensions.

inline void add_dimension(const sweep_parameter &param)

Adds a dimension to sweep over. The first dimension is the x dimension, the second dimension is the y dimension, etc.

Parameters:

param – The dimension to add.

inline const sweep_parameter &get_dimension(const std::size_t index) const

Returns a specific dimension by index.

Parameters:

index – The index of the dimension to return.

Throws:

std::out_of_range – if the index is out of range.

Returns:

The dimension at the specified index.

inline std::size_t get_number_of_dimensions() const noexcept

Returns the number of dimensions to sweep over.

Returns:

The number of dimensions to sweep over.

inline double minimum_ct() const noexcept

Finds the minimum critical temperature in the domain.

Returns:

The minimum critical temperature.

inline double maximum_ct() const noexcept

Finds the maximum critical temperature in the domain.

Returns:

The maximum critical temperature.

struct operational_domain_value_range

A range of values for a dimension sweep. The range is defined by a minimum value, a maximum value and a step size.

Public Members

sweep_parameter dimension

The sweep parameter of the dimension.

double min = {1.0}

The minimum value of the dimension sweep.

double max = {10.0}

The maximum value of the dimension sweep.

double step = {0.1}

The step size of the dimension sweep.

struct operational_domain_params

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

Public Members

is_operational_params operational_params = {}

The parameters used to determine if a layout is operational or non-operational.

std::vector<operational_domain_value_range> sweep_dimensions = {operational_domain_value_range{sweep_parameter::EPSILON_R, 1.0, 10.0, 0.1}, operational_domain_value_range{sweep_parameter::LAMBDA_TF, 1.0, 10.0, 0.1}}

The dimensions to sweep over together with their value ranges, ordered by priority. The first dimension is the x dimension, the second dimension is the y dimension, etc.

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.

std::size_t num_total_parameter_points = {0}

Total number of parameter points in the parameter space.

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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

Returns:

The 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 determine all operational “islands” 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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

Returns:

The 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 set of 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. This is repeated for all initially sampled points that do not lie within a contour. The algorithm is guaranteed to determine the contours of all operational “islands” 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. For each thusly discovered operational island, it performs another number of samples equal to the distance to an edge of each operational area. Finally, it performs up to 8 samples for each contour point (however, the actual number is usually 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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

Returns:

The operational domain of the layout.

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

Computes the critical temperature domain of the given SiDB cell-level layout. The critical temperature domain consists of all parameter combinations for which the layout is logically operational, along with the critical temperature for each specific parameter point.

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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

Returns:

The critical temperature domain of the layout.

template<typename Lyt, typename TT>
critical_temperature_domain fiction::critical_temperature_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 critical temperature domain of the given SiDB cell-level layout. The critical temperature domain consists of all parameter combinations for which the layout is logically operational, along with the critical temperature for each specific parameter point.

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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

Returns:

The critical temperature domain of the layout.

template<typename Lyt, typename TT>
critical_temperature_domain fiction::critical_temperature_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 critical temperature domain of the given SiDB cell-level layout. The critical temperature domain consists of all parameter combinations for which the layout is logically operational, along with the critical temperature for each specific parameter point.

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 determine all operational “islands” 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.

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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

Returns:

The critical temperature domain of the layout.

template<typename Lyt, typename TT>
critical_temperature_domain fiction::critical_temperature_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 critical temperature domain of the given SiDB cell-level layout. The critical temperature domain consists of all parameter combinations for which the layout is logically operational, along with the critical temperature for each specific parameter point.nt.

This algorithm first uses random sampling to find a set of 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. This is repeated for all initially sampled points that do not lie within a contour. The algorithm is guaranteed to determine the contours of all operational “islands” if the initial random sampling found at least one operational point within them. Thereby, this algorithm works for disconnected operational domains. The critical temperature is computed for each operational point.

It performs samples uniformly-distributed random samples within the parameter range. For each thusly discovered operational island, it performs another number of samples equal to the distance to an edge of each operational area. Finally, it performs up to 8 samples for each contour point (however, the actual number is usually 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.

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.

Throws:

std::invalid_argument – if the given sweep parameters are invalid.

Returns:

The critical temperature domain of the layout.

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

struct operational_domain_ratio_params

Parameters for computing the ratio of operational parameter points around a specified parameter point to the total number of parameter points in the given parameter space.

Public Members

operational_domain_params op_domain_params = {}

Parameters for the operational domain computation.

template<typename Lyt, typename TT>
double fiction::operational_domain_ratio(const Lyt &lyt, const std::vector<TT> &spec, const parameter_point &pp, const operational_domain_ratio_params &params = {}) noexcept

Calculates the ratio of operational parameter points surrounding a specified parameter point to the total number of parameter points in the given parameter space. This function is useful for assessing how robust a gate design is to variations in its parameters.

A ratio close to 1 indicates that the gate is robust, meaning it functions correctly across a broad range of parameter values. A ratio close to 0 indicates that the gate is highly sensitive to parameter variations and may fail to operate correctly.

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

  • TT – Truth table type.

Parameters:
  • lyt – The SiDB layout for which to compute the ratio of operational parameter points surrounding a specified parameter point to the total number of parameter points.

  • spec – The expected Boolean function of the layout, provided as a multi-output truth table.

  • params – Parameters.

  • pp – The specific parameter point around which the operational ratio is computed.

Returns:

The ratio of operational parameter points to the total number of parameter points in the parameter space.

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

template<typename Lyt, typename TT>
operational_status fiction::verify_logic_match(const charge_distribution_surface<Lyt> &cds, const is_operational_params &params, const std::vector<TT> &spec, const uint64_t input_pattern, const std::vector<bdl_wire<Lyt>> &input_wires, const std::vector<bdl_wire<Lyt>> &output_wires) noexcept

Checks if a given charge distribution correctly encodes the expected logic for a specified input pattern, based on a provided truth table.

Example: In the ground state charge distribution of an AND gate, kinks are rejected for the gate to be considered operational. Given an input pattern of 01, this function will:

  • Verify that the left input wire encodes 0.

  • Verify that the right input wire encodes 1.

  • Verify that the output wire encodes 0.

Note

Kinks are rejected.

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

  • TT – Truth table type.

Parameters:
  • cds – Charge distribution surface, containing charge state information for each SiDB.

  • params – The parameters used to determine if a layout is operational or non-operational.

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

  • input_pattern – The specific input pattern of the given charge distribution surface.

  • input_wires – Input BDL wires.

  • output_wires – Output BDL wires.

Returns:

The operational status indicating if the charge distribution matches the logic for the given input pattern.

Physically Valid Parameters

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

template<typename Lyt>
sidb_simulation_domain<parameter_point, uint64_t> fiction::physically_valid_parameters(Lyt &cds, const operational_domain_params &params = {}) noexcept

This function computes the physical parameters necessary for ensuring the physical validity of a given charge distribution and determines the corresponding excited state number. The ground state is denoted by zero, with each subsequent excited state incrementally numbered.

This function is designed to derive the physical parameters from charge distribution measurements of SiDB layouts, often acquired through Atomic Force Microscopy (AFM). Given a specific charge distribution, the function typically yields several physically valid parameters.

As more SiDB layouts with corresponding charge distributions are recorded, the number of physically valid parameters for all layouts decreases. Consequently, this enables a more precise determination of the physical parameters present on the surface.

Template Parameters:

Lyt – The charge distribution surface type.

Parameters:
  • cds – The charge distribution surface for which physical parameters are to be determined.

  • params – Operational domain parameters.

Returns:

Physically valid parameters with the corresponding excited state number of the given charge distribution surface for each parameter point.

Displacement Robustness Domain

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

template<typename Lyt>
struct displacement_robustness_domain

During fabrication, SiDBs may not align precisely with their intended atomic positions, resulting in displacement. This means that an SiDB is fabricated close to the desired one, typically one or a few H-Si positions away. Consequently, depending on the fabrication speed, a certain number of SiDBs may experience displacement. To address and analyze this occurrence, we introduce the Displacement Robustness Domain. This domain consists of SiDB layouts derived from an original layout, each showing displaced SiDBs, together with the operational or non-operational status, based on the specified logic.

Public Members

std::vector<std::pair<Lyt, operational_status>> operational_values = {}

Represents a domain of displacement robustness for layouts resulting from applying a displacement to a given SiDB layout.

template<typename CellType>
struct displacement_robustness_domain_params

Parameters for the determine_displacement_robustness_domain and determine_probability_of_fabricating_operational_gate algorithms.

Param CellType:

SiDB layout cell type.

Public Types

enum class displacement_analysis_mode : uint8_t

Possible modes to determine the displacement robustness domain.

Values:

enumerator EXHAUSTIVE

All possible displacements are analyzed.

enumerator RANDOM

A certain amount of all possible displacements is analyzed randomly. Defined by percentage_of_analyzed_displaced_layouts.

enum class dimer_displacement_policy : uint8_t

Specifies the allowed displacement range options for SiDB fabrication simulation.

Values:

enumerator STAY_ON_ORIGINAL_DIMER

In this mode, any displacement of SiDBs must remain within the boundaries of the initial dimer they are placed on.

enumerator ALLOW_OTHER_DIMER

In this mode, SiDBs are allowed to be displaced from the original dimer to any other dimer within the layout.

Public Members

displacement_analysis_mode analysis_mode = {displacement_analysis_mode::EXHAUSTIVE}

This parameter defines the mode of the displacement. If EXHAUSTIVE, all possible displacements are analyzed. Otherwise, a certain amount of all possible displacements is analyzed randomly.

double percentage_of_analyzed_displaced_layouts = {1.0}

This parameter defines the percentage of all possible displaced SiDB layouts that are analyzed. The default value is 1.0 (100 %), which means that all possible displacements are covered.

std::pair<uint64_t, uint64_t> displacement_variations = {1, 0}

Possible displacement range of H-Si positions in the x- and y-directions. The default value is (1, 0), which means that displacements of ±1 position in the x-direction are analyzed, with no displacement in the y-direction.

is_operational_params operational_params = {}

Parameters to check the operational status of the SiDB layout.

std::set<CellType> fixed_sidbs = {}

SiDBs in the given layout which shall not be affected by variations.

dimer_displacement_policy dimer_policy = {dimer_displacement_policy::STAY_ON_ORIGINAL_DIMER}

This flag controls whether the displacement in the y-direction can lead to changes in the Si dimer.

struct displacement_robustness_domain_stats

Statistics for the displacement robustness domain computation.

Public Members

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

Total runtime in seconds to determine the robustness of the passed SiDB layout.

std::size_t num_operational_sidb_displacements = {0}

The number of operational SiDB layouts resulting from the given layout by displacements.

std::size_t num_non_operational_sidb_displacements = {0}

The number of non-operational SiDB layouts resulting from the given layout by displacements.

template<typename Lyt>
struct displacement_robustness_domain

During fabrication, SiDBs may not align precisely with their intended atomic positions, resulting in displacement. This means that an SiDB is fabricated close to the desired one, typically one or a few H-Si positions away. Consequently, depending on the fabrication speed, a certain number of SiDBs may experience displacement. To address and analyze this occurrence, we introduce the Displacement Robustness Domain. This domain consists of SiDB layouts derived from an original layout, each showing displaced SiDBs, together with the operational or non-operational status, based on the specified logic.

Public Members

std::vector<std::pair<Lyt, operational_status>> operational_values = {}

Represents a domain of displacement robustness for layouts resulting from applying a displacement to a given SiDB layout.

template<typename Lyt, typename TT>
displacement_robustness_domain<Lyt> fiction::determine_displacement_robustness_domain(const Lyt &layout, const std::vector<TT> &spec, const displacement_robustness_domain_params<cell<Lyt>> &params = {}, displacement_robustness_domain_stats *stats = nullptr)

During fabrication, SiDBs may not align precisely with their intended atomic positions, resulting in displacement. This means that an SiDB is fabricated close to the desired one, typically one or a few H-Si positions away. Consequently, depending on the fabrication speed, a certain number of SiDBs may experience displacement.

This function determines the operational status of all possible displacements of the SiDBs of the given SiDB layout, based on the provided truth table specification and displacement robustness computation parameters. The number of displacements grows exponentially with the number of SiDBs. For small layouts, all displacements can be analyzed. For larger layouts, random sampling can be applied, controllable by the analysis_mode and percentage_of_analyzed_displaced_layouts in `params.

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

  • TT – Truth table type.

Parameters:
  • spec – Vector of truth table specifications.

  • params – Parameters for the displacement robustness computation.

  • stats – Statistics related to the displacement robustness computation.

Returns:

The displacement robustness domain of the SiDB layout.

template<typename Lyt, typename TT>
double fiction::determine_probability_of_fabricating_operational_gate(const Lyt &layout, const std::vector<TT> &spec, const displacement_robustness_domain_params<cell<Lyt>> &params = {}, const double fabrication_error_rate = 1.0)

During fabrication, SiDBs may not align precisely with their intended atomic positions, resulting in displacement. This means that an SiDB is fabricated close to the desired one, typically one or a few H-Si positions away. The percentage of displaced SiDBs depends on the fabrication speed. Therefore, SiDB layouts with high displacement tolerance are preferred to speed up the fabrication process.

This function calculates the probability of fabricating an operational SiDB layout for an originally given SiDB layout and a given fabrication error rate. A fabrication error rate of 0.0 or negative indicates that the SiDB layout is designed without displacement.

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

  • TT – The type of the truth table.

Parameters:
  • layout – The SiDB cell-level layout which is analyzed.

  • spec – Vector of truth table specifications.

  • params – Parameters for the displacement robustness computation.

  • fabrication_error_rate – The fabrication error rate. For example, 0.1 describes that 10% of all manufactured SiDBs have a slight displacement.

Returns:

The probability of fabricating an operational SiDB 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.

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 CellType>
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:

CellType – Cell type.

Public Functions

bdl_pair() = default

Standard constructor for empty BDL pairs.

inline bdl_pair(const sidb_technology::cell_type t, const CellType &u, const CellType &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.

inline constexpr bool operator==(const bdl_pair<CellType> &other) const noexcept

Equality operator.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if this BDL pair is equal to the other, false otherwise.

inline constexpr bool operator!=(const bdl_pair<CellType> &other) const noexcept

Inequality operator.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if this BDL pair is not equal to the other, false otherwise.

inline constexpr bool operator<(const bdl_pair<CellType> &other) const noexcept

Less than operator.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if this BDL pair is less than the other, false otherwise.

inline constexpr bool operator<=(const bdl_pair<CellType> &other) const noexcept

Less than or equal to operator.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if this BDL pair is less than or equal to the other, false otherwise.

inline constexpr bool operator>(const bdl_pair<CellType> &other) const

Greater than operator.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if this BDL pair is greater than the other, false otherwise.

inline constexpr bool operator>=(const bdl_pair<CellType> &other) const noexcept

Greater-than-or-equal-to operator.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if this BDL pair is greater than or equal to the other, otherwise false.

inline constexpr bool equal_ignore_type(const bdl_pair<CellType> &other) const noexcept

Equality check ignoring type.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if the upper and lower attributes are equal, otherwise false.

inline constexpr bool not_equal_ignore_type(const bdl_pair<CellType> &other) const noexcept

Inequality check ignoring type.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if the upper and lower attributes are not equal, otherwise false.

inline constexpr bool has_same_y_coordinate(const bdl_pair<CellType> &other) const noexcept

Checks if the upper and lower SiDBs in this BDL pair have the same y-coordinate as the corresponding SiDBs in the other BDL pair.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if both the upper and lower SiDBs in this pair have the same y-coordinate as the corresponding SiDBs in the other pair, otherwise false.

inline constexpr bool has_same_x_coordinate(const bdl_pair<CellType> &other) const noexcept

Checks if the upper and lower SiDBs in this BDL pair have the same x-coordinate as the corresponding SiDBs in the other BDL pair.

Parameters:

other – The other BDL pair to compare with.

Returns:

true if both the upper and lower SiDBs in this pair have the same x-coordinate as the corresponding SiDBs in the other pair, otherwise false.

Public Members

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.

CellType upper = {}

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

CellType 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 SiDBs 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<cell<Lyt>>> fiction::detect_bdl_pairs(const Lyt &lyt, const std::optional<typename technology<Lyt>::cell_type> &type = std::nullopt, 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 – Optional parameter to specify the SiDB type for which BDL pairs should be detected. If omitted, the function will detect BDL pairs for all types. Valid types include INPUT, OUTPUT, NORMAL, among others.

  • params – Parameters for the BDL pair detection algorithm.

Returns:

A vector of BDL pairs.

Binary-dot Logic (BDL) Wire Detection

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

enum class fiction::bdl_wire_selection : std::uint8_t

An enumeration of the selection of different types of wires.

Values:

enumerator ALL

Select all BDL wires.

enumerator INPUT

Select only BDL wires that start with input cells.

enumerator OUTPUT

Select only BDL wires that end with output cells.

struct detect_bdl_wires_params

This struct encapsulates parameters used for detecting BDL wires.

Public Members

double threshold_bdl_interdistance = 2.0

A distance threshold, which is used to determine if two pairs of BDLs are part of the same wire. (unit: nm).

detect_bdl_pairs_params bdl_pairs_params = {}

Parameters for the detect_bdl_pairs algorithm.

template<typename Lyt>
struct bdl_wire

This struct encapsulates a vector of bdl_pair objects, representing the pairs of SiDBs in the BDL wire.

Template Parameters:

Lyt – SiDB cell-level layout type.

Public Functions

inline bdl_wire() noexcept

Default constructor for an empty BDL wire.

inline explicit bdl_wire(const std::vector<bdl_pair<cell<Lyt>>> &p) noexcept

Constructor to initialize the BDL wire with a given vector of BDL pairs.

Also updates the start and end BDL pairs based on the given vector.

Parameters:

p – The vector of BDL pairs to initialize the wire with.

inline bdl_wire(const bdl_wire &other) noexcept

Copy constructor.

Creates a new bdl_wire object as a copy of another bdl_wire object.

Parameters:

other – The bdl_wire object to copy from.

inline bdl_wire(bdl_wire &&other) noexcept

Move constructor.

Transfers ownership of the BDL pairs, port, and start/end pairs from another bdl_wire object.

Parameters:

other – The bdl_wire object to move from.

inline bdl_wire &operator=(bdl_wire &&other) noexcept

Move assignment operator.

Transfers ownership of the BDL pairs, port, and start/end pairs from another bdl_wire object.

Parameters:

other – The bdl_wire object to move from.

Returns:

A reference to the updated object.

inline bdl_wire &operator=(const bdl_wire &other) noexcept

Copy assignment operator.

Copies the content of another bdl_wire object, including start and end pairs.

Parameters:

other – The bdl_wire object to copy from.

Returns:

A reference to the updated object.

~bdl_wire() noexcept = default

Destructor for bdl_wire.

inline void add_bdl_pair(const bdl_pair<cell<Lyt>> &pair) noexcept

Add a BDL pair to the wire.

Parameters:

pair – The BDL pair to add.

inline void erase_bdl_pair(const bdl_pair<cell<Lyt>> &pair) noexcept

Erase a specific BDL pair from the wire.

Parameters:

pair – The BDL pair to remove. The pair is compared using the equality operator (operator==).

inline std::optional<bdl_pair<cell<Lyt>>> find_bdl_pair_by_type(const sidb_technology::cell_type &t) const noexcept

Find the first Binary-dot Logic (BDL) pair of a specified type in the wire.

Parameters:

t – Type of BDL pair to search for (sidb_technology::cell_type::INPUT, sidb_technology::cell_type::OUTPUT, etc.).

Returns:

Optional containing the first BDL pair with the specified type t, or std::nullopt if no such BDL pair is found.

Public Members

std::vector<bdl_pair<cell<Lyt>>> pairs = {}

Vector of BDL pairs representing the wire.

port_direction port = {port_direction::NONE}

Port of the BDL wire.

std::optional<bdl_pair<cell<Lyt>>> first_bdl_pair = {}

First BDL pair of the wire.

std::optional<bdl_pair<cell<Lyt>>> last_bdl_pair = {}

Last BDL pair of the wire.

template<typename Lyt>
std::vector<bdl_wire<Lyt>> fiction::detect_bdl_wires(const Lyt &lyt, const detect_bdl_wires_params &params = {}, const bdl_wire_selection wire_selection = bdl_wire_selection::ALL) noexcept

This function identifies BDL wires in a given SiDB cell-level layout by detecting BDL pairs and linking them based on their spatial relationships. The function supports selection of different types of wires, such as input wires, output wires, or all wires.

Template Parameters:

Lyt – SiDB cell-level layout type.

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

  • params – Parameters used for detecting BDL wires.

  • wire_selection – The type of wires to detect, specified by the bdl_wire_selection enum. Default is bdl_wire_selection::ALL.

Returns:

A vector of BDL wires, where each wire is represented as a vector of BDL pairs.

Assess Population Stability

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

enum class fiction::transition_type : uint8_t

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.

std::unordered_map<transition_type, std::pair<typename Lyt::cell, double>> transition_potentials = {}

This map collects all charge transition types, the corresponding critical cells and the required electrostatic potential (unit: V) required to conduct the transition.

std::unordered_map<transition_type, double> distance_corresponding_to_potential

This map collects for all charge transition types, the electrostatic potential difference which is required to conduct a charge change as a distance in nanometer. This is possible since the electrostatic potential is connected to the distance.

double system_energy = {}

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

struct physical_population_stability_params

This struct stores the parameters required to simulate 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::physical_population_stability(const Lyt &lyt, const physical_population_stability_params &params) noexcept

This function simulates 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 simulated.

  • params – Parameters used to simulate the population stability.

Returns:

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

Band-Bending Resilience

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

struct band_bending_resilience_params

This struct stores the parameters required to simulate the band bending resilience of an SiDB layout

Public Members

physical_population_stability_params assess_population_stability_params = {}

Parameters for the assessing physical population stability simulation

bdl_input_iterator_params bdl_iterator_params = {}

Parameters for the input BDL iterator.

template<typename Lyt, typename TT>
double fiction::band_bending_resilience(const Lyt &lyt, const std::vector<TT> &spec, const band_bending_resilience_params &params = {}, const std::optional<transition_type> transition_type = std::nullopt) noexcept

Calculates the band bending resilience. This is the minimum electrostatic potential required to induce a charge change in an SiDB layout among all possible input combinations which was proposed in “Unifying Figures of Merit: A Versatile Cost Function for Silicon Dangling Bond Logic” by J. Drewniok, M. Walter, S. S. H. Ng, K. Walus, and R. Wille in IEEE NANO 2024 (https://ieeexplore.ieee.org/abstract/document/10628671).

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

  • TT – Truth table type.

Parameters:
  • lyt – Layout for which the band bending resilience is calculated.

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

  • params – Parameters for assessing physical population stability.

  • transition_type – The optional type of charge transition to consider. This can be used if one is only interested in a specific type of charge transition.

Returns:

The minimum potential (in V) required for charge change across all input combinations.

Convert Potential to Distance

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

inline double fiction::potential_to_distance_conversion(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.

Fabrication Defects

A collection of tools to simulate defects that can occur during the fabrication process of FCN technologies.

SiDB Defect Types

Header: fiction/technology/sidb_defects.hpp

enum class fiction::sidb_defect_type : uint8_t

Specifies the types of fabrication defects that can occur on the H-Si(100) 2x1 surface according to “Atomic defect classification of the H–Si(100) surface through multi-mode scanning probe microscopy” by Jeremiah Croshaw, Thomas Dienel, Taleana Huff, and Robert Wolkow in Journal of Nanotechnology in 2020.

Values:

enumerator NONE

Defect-free H-Si.

enumerator DB

A stray dangling bond.

enumerator SI_VACANCY

A missing silicon atom.

enumerator SINGLE_DIHYDRIDE

Double hydrogen passivation.

enumerator DIHYDRIDE_PAIR

A missing bond between dimers that leads to two double hydrogen passivations.

enumerator ONE_BY_ONE

A collection of dihydride pairs.

enumerator THREE_BY_ONE

A collection of 1 by 1’s.

enumerator SILOXANE

An oxidized dimer.

enumerator RAISED_SI

A raised silicon dimer.

enumerator MISSING_DIMER

The dimer is missing altogether.

enumerator ETCH_PIT

A collection of missing dimers.

enumerator STEP_EDGE

A step edge, which is a break in the surface reconstruction.

enumerator GUNK

Residual material.

enumerator UNKNOWN

Unknown defect.

struct sidb_defect

In accordance with the paper mentioned above, the sidb_defect struct is used to represent a specific defect on the H-Si(100) 2x1 surface that has a charge as well as relative permittivity (epsilon_r) and Thomas-Fermi screening distance (lambda_tf) values associated to it.

See “SiQAD: A Design and Simulation Tool for Atomic Silicon Quantum Dot Circuits” by S. S. H. Ng, J. Retallick, H. N. Chiu, R. Lupoiu, L. Livadaru, T. Huff, M. Rashidi, W. Vine, T. Dienel, R. A. Wolkow, and K. Walus in IEEE Transactions on Nanotechnology for more details on these values.

Public Functions

inline explicit constexpr sidb_defect(const sidb_defect_type defect_type = sidb_defect_type::UNKNOWN, const int64_t electric_charge = 0.0, const double relative_permittivity = 0.0, const double screening_distance = 0.0) noexcept

Standard constructor.

inline constexpr bool operator==(const sidb_defect &rhs) const noexcept

This operator compares two sidb_defect instances for equality. It checks if the type, charge, epsilon_r, and lambda_tf members of the two instances are equal.

Parameters:

rhssidb_defect instance to compare against.

inline constexpr bool operator!=(const sidb_defect &rhs) const noexcept

This operator compares two sidb_defect instances for inequality. It uses the operator== to check if the two instances are equal and returns the negation of the result.

Parameters:

rhssidb_defect instance to compare against.

Public Members

sidb_defect_type type

Type of defect.

int64_t charge

Electrical charge in units of the elementary charge e (e.g., 1 ^= 1*e, -2 ^= -2*e).

double epsilon_r

Electric permittivity (unitless).

double lambda_tf

Thomas-Fermi screening distance in nm.

static constexpr bool fiction::is_charged_defect_type(const sidb_defect &defect) noexcept

Checks whether the given defect type is a charged one. DB and SI_VACANCY types are charged. Those charged defects are to be avoided by a larger distance.

Parameters:

defect – Defect to check.

Returns:

true iff defect is of a charged type.

static constexpr bool fiction::is_neutral_defect_type(const sidb_defect &defect) noexcept

Checks whether the given defect type is not a charged one. Neutral defects are to be avoided as well, but not by such a large distance. Even though the NONE defect type is technically neutral, it is not a defect per se which is why this function returns false on the NONE defect input.

Parameters:

defect – Defect to check.

Returns:

true iff defect is not of a charged type.

static constexpr bool fiction::is_positively_charged_defect(const sidb_defect &defect) noexcept

Checks whether the given defect has a positive charge value assigned to it. This function is irrespective of the associated defect type.

Parameters:

defect – Defect to check.

Returns:

true iff defect has a positive charge value.

static constexpr bool fiction::is_negatively_charged_defect(const sidb_defect &defect) noexcept

Checks whether the given defect has a negative charge value assigned to it. This function is irrespective of the associated defect type.

Parameters:

defect – Defect to check.

Returns:

true iff defect has a negative charge value.

static constexpr bool fiction::is_neutrally_charged_defect(const sidb_defect &defect) noexcept

Checks whether the given defect has a neutral charge value, i.e., 0, assigned to it. This function is irrespective of the associated defect type.

Parameters:

defect – Defect to check.

Returns:

true iff defect has a neutral charge value.

constexpr const uint16_t fiction::SIDB_CHARGED_DEFECT_HORIZONTAL_SPACING = 26u

Horizontal distance to keep from charged SiDB defects. The value is to be understood as the number of DB positions rather than the number of dimers. This is true even though each defect always affects the entire dimer.

constexpr const uint16_t fiction::SIDB_CHARGED_DEFECT_VERTICAL_SPACING = 13u

Vertical distance to keep from charged SiDB defects. The value is to be understood as the number of DB positions rather than the number of dimers. This is true even though each defect always affects the entire dimer.

constexpr const uint16_t fiction::SIDB_NEUTRAL_DEFECT_HORIZONTAL_SPACING = 1u

Horizontal distance to keep from neutral SiDB defects. The value is to be understood as the number of DB positions rather than the number of dimers. This is true even though each defect always affects the entire dimer.

constexpr const uint16_t fiction::SIDB_NEUTRAL_DEFECT_VERTICAL_SPACING = 0u

Vertical distance to keep from neutral SiDB defects. The value is to be understood as the number of DB positions rather than the number of dimers. This is true even though each defect always affects the entire dimer.

static constexpr std::pair<uint16_t, uint16_t> fiction::defect_extent(const sidb_defect &defect, const std::optional<std::pair<uint16_t, uint16_t>> &charged_defect_spacing_overwrite = std::nullopt, const std::optional<std::pair<uint16_t, uint16_t>> &neutral_defect_spacing_overwrite = std::nullopt) noexcept

Returns the extent of a defect as a pair of SiDB distances in the horizontal and vertical directions. If the defect type is NONE, {0, 0} is returned.

Parameters:
  • defect – Defect type to evaluate.

  • charged_defect_spacing_overwrite – Override the default influence distance of charged atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

  • neutral_defect_spacing_overwrite – Override the default influence distance of neutral atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

Returns:

A pair of uint16_t values representing the number of horizontal and vertical SiDBs affected by the given defect type.

SiDB Defect Surface

Header: fiction/technology/sidb_defect_surface.hpp

A layout type to layer on top of any SiDB cell-level layout. It implements an interface to store and access fabrication defects on the H-Si(100) 2x1 surface.

template<typename Lyt, bool has_sidb_defect_surface = std::conjunction_v<has_assign_sidb_defect<Lyt>, has_get_sidb_defect<Lyt>>>
class sidb_defect_surface : public Lyt

A layout type to layer on top of any SiDB cell-level layout. It implements an interface to store and access fabrication defects on the H-Si(100) 2x1 surface.

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

  • has_sidb_defect_surface – Automatically determines whether a defect interface is already present.

template<typename Lyt>
class sidb_defect_surface<Lyt, true> : public Lyt
template<typename Lyt>
class sidb_defect_surface<Lyt, false> : public Lyt

Public Functions

inline explicit sidb_defect_surface(const sidb_defect_surface_params &ps = {})

Standard constructor for empty layouts.

Parameters:

ps – SiDB defect surface parameters.

inline explicit sidb_defect_surface(const typename Lyt::aspect_ratio &ar, const sidb_defect_surface_params &ps = {})

Standard constructor that layers the SiDB defect interface onto a layout with aspect ratio ar as input.

Parameters:
  • ar – aspect ratio of the layout.

  • ps – SiDB defect surface parameters.

inline explicit sidb_defect_surface(const Lyt &lyt, const sidb_defect_surface_params &ps = {})

Standard constructor that layers the SiDB defect interface onto an existing layout.

Parameters:
  • lyt – Existing layout that is to be extended by an SiDB defect interface.

  • ps – SiDB defect surface parameters.

inline sidb_defect_surface clone() const noexcept

Clones the layout returning a deep copy.

Returns:

Deep copy of the layout.

inline void assign_sidb_defect(const typename Lyt::coordinate &c, const sidb_defect &d) noexcept

Assigns a given defect type to the given coordinate.

Parameters:
  • c – Coordinate to assign defect d to.

  • d – Defect to assign to coordinate c.

inline void move_sidb_defect(const typename Lyt::coordinate &source, const typename Lyt::coordinate &target) noexcept

Moves an SiDB defect from one cell to another.

Parameters:
  • source – Source coordinate of the defect.

  • target – Target coordinate to move the defect to.

inline sidb_defect get_sidb_defect(const typename Lyt::coordinate &c) const noexcept

Returns the given coordinate’s assigned defect type. If no defect type has been assigned, NONE is returned.

Parameters:

c – Coordinate to check.

Returns:

Defect type previously assigned to c or NONE if no defect was yet assigned.

inline uint64_t num_defects() const noexcept

Returns the number of defective coordinates on the surface.

Returns:

Number of defective coordinates.

inline std::size_t num_positively_charged_defects() const noexcept

Number of positively charged defects on the surface.

Returns:

Number of positively charged defects.

inline std::size_t num_negatively_charged_defects() const noexcept

Returns the number of negatively charged defects.

Returns:

Number of negatively charged defects.

inline std::size_t num_charged_defects() const noexcept

Returns the number of charged defects.

Returns:

Number of charged defects.

inline std::size_t num_neutral_defects() const noexcept

Returns the number of neutral defects.

Returns:

Number of neutral defects.

template<typename Fn>
inline void foreach_sidb_defect(Fn &&fn) const

Applies a function to all defects on the surface. Since the defects are fetched directly from the storage map, the given function has to receive a pair of a coordinate and a defect type as its parameter.

Template Parameters:

Fn – Functor type that has to comply with the restrictions imposed by mockturtle::foreach_element.

Parameters:

fn – Functor to apply to each defect.

inline std::unordered_set<typename Lyt::coordinate> affected_sidbs(const typename Lyt::coordinate &c, const std::optional<std::pair<uint16_t, uint16_t>> &charged_defect_spacing_overwrite = std::nullopt, const std::optional<std::pair<uint16_t, uint16_t>> &neutral_defect_spacing_overwrite = std::nullopt) const noexcept

Returns all SiDB positions affected by the defect at the given coordinate. This function relies on the defect_extent function defined in sidb_defects.hpp that computes the extent of charged and neutral defect types.

If the given coordinate is defect-free, the empty set is returned.

Parameters:
  • c – Coordinate whose defect extent is to be determined.

  • charged_defect_spacing_overwrite – Override the default influence distance of charged atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

  • neutral_defect_spacing_overwrite – Override the default influence distance of neutral atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

Returns:

All SiDB positions affected by the defect at coordinate c.

inline std::unordered_set<typename Lyt::coordinate> all_affected_sidbs(const std::optional<std::pair<uint64_t, uint64_t>> &charged_defect_spacing_overwrite = std::nullopt, const std::optional<std::pair<uint64_t, uint64_t>> &neutral_defect_spacing_overwrite = std::nullopt) const noexcept

Returns all SiDB positions affected by any defect on the surface. This function relies on the defect_extent function defined in sidb_defects.hpp that computes the extent of charged and neutral defect types.

If the given surface is defect-free, the empty set is returned.

Parameters:
  • charged_defect_spacing_overwrite – Override the default influence distance of charged atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

  • neutral_defect_spacing_overwrite – Override the default influence distance of neutral atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

Returns:

All SiDB positions affected by any defect on the surface.

template<>
struct sidb_surface_storage

SiDB Defect Analysis

Header: fiction/technology/sidb_surface_analysis.hpp

template<typename GateLibrary, typename GateLyt, typename CellLyt>
auto fiction::sidb_surface_analysis(const GateLyt &gate_lyt, const CellLyt &surface, const std::optional<std::pair<uint64_t, uint64_t>> &charged_defect_spacing_overwrite = std::nullopt, const std::optional<std::pair<uint64_t, uint64_t>> &neutral_defect_spacing_overwrite = std::nullopt) noexcept

Analyzes a given defective SiDB surface and matches it against gate tiles provided by a library. Any gate type that cannot be realized on a certain tile due to disturbances caused by defects gets blacklisted on said tile. The black list is then returned by this function.

Note

The given gate library must implement both the get_functional_implementations() and get_gate_ports() functions.

Template Parameters:
  • GateLibrary – FCN gate library type to fetch the gate descriptions from.

  • GateLyt – Gate-level layout type that specifies the tiling of the SiDB surface.

  • CellLyt – SiDB cell-level layout type that is underlying to the SiDB defect surface.

Parameters:
  • gate_lyt – Gate-level layout instance that specifies the aspect ratio.

  • surface – SiDB surface that instantiates the defects.

  • charged_defect_spacing_overwrite – Override the default influence distance of charged atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

  • neutral_defect_spacing_overwrite – Override the default influence distance of neutral atomic defects on SiDBs with an optional pair of horizontal and vertical distances.

Returns:

A black list of gate functions associated with tiles.

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

template<typename CellType>
struct defect_influence_params

Parameters to determine the defect influence.

Template Parameters:

CellType – Type of the cell.

Public Types

enum class influence_definition : uint8_t

Definition of defect influence.

Values:

enumerator OPERATIONALITY_CHANGE

Influence is considered as the ability to change the operational status of the layout.

enumerator GROUND_STATE_CHANGE

Influence is considered as the ability to change the ground state of the layout.

Public Members

sidb_defect defect = {}

The defect to calculate the defect influence for.

is_operational_params operational_params = {}

Parameters for the is_operational algorithm.

CellType additional_scanning_area = {50, 6}

Area around the layout for additional defect scanning. This describes the additional space around the bounding box of the layout.

influence_definition influence_def = {influence_definition::OPERATIONALITY_CHANGE}

Definition of defect influence.

enum class fiction::defect_influence_status : uint8_t

Defines whether the influence of a defect is present at a particular position in the layout. It can be used to classify positions as having an influence or not.

Values:

enumerator INFLUENTIAL

This indicates that the defect is actively influencing the layout at this position. It implies that some form of impact, such as a change in operational status or ground state, is being caused by the defect at this position.

enumerator NON_INFLUENTIAL

This indicates that the defect does not influence the layout at this position. It implies that the layout remains unaffected by the defect at this location, meaning there is no change in the operational status or the ground state.

template<typename Lyt>
class defect_influence_domain : public fiction::sidb_simulation_domain<Lyt::cell, defect_influence_status>

A defect_influence_domain defines for each defect position the influence of the defect on the layout. Depending on the chosen definition of influence, this can either mean that the operational status or the ground state of the layout is changed due to the presence of the defect.

struct defect_influence_stats

Statistics.

template<typename Lyt, typename TT>
defect_influence_domain<Lyt> fiction::defect_influence_grid_search(const Lyt &lyt, const std::vector<TT> &spec, const defect_influence_params<cell<Lyt>> &params = {}, const std::size_t step_size = 1, defect_influence_stats *stats = nullptr)

This algorithm uses a grid search to determine the defect influence domain. The grid search is performed by exhaustively sweeping all possible atomic defect positions in x and y dimensions.

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

  • TT – Truth table type.

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

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

  • step_size – The parameter specifying the interval between consecutive defect positions to be evaluated.

  • params – Defect influence domain computation parameters.

  • stats – Statistics.

Returns:

The defect influence domain of the layout.

template<typename Lyt>
defect_influence_domain<Lyt> fiction::defect_influence_grid_search(const Lyt &lyt, const defect_influence_params<cell<Lyt>> &params = {}, const std::size_t step_size = 1, defect_influence_stats *stats = nullptr)

This algorithm uses a grid search to determine the defect influence domain. The grid search is performed by exhaustively sweeping all possible atomic defect positions in x and y dimensions.

Template Parameters:

Lyt – SiDB cell-level layout type.

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

  • step_size – The parameter specifying the interval between consecutive defect positions to be evaluated.

  • params – Defect influence domain computation parameters.

  • stats – Statistics.

Returns:

The defect influence domain of the layout.

template<typename Lyt, typename TT>
defect_influence_domain<Lyt> fiction::defect_influence_random_sampling(const Lyt &lyt, const std::vector<TT> &spec, std::size_t samples, const defect_influence_params<cell<Lyt>> &params = {}, defect_influence_stats *stats = nullptr)

This algorithm uses random sampling to find a part of the defect influence domain that might not be complete. It performs a total of `samples uniformly-distributed random samples within the specified area.

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

  • TT – Truth table type.

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

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

  • samples – Number of random samples to perform.

  • params – Defect influence domain computation parameters.

  • stats – Statistics.

Returns:

The (partial) defect influence domain of the layout.

template<typename Lyt>
defect_influence_domain<Lyt> fiction::defect_influence_random_sampling(const Lyt &lyt, std::size_t samples, const defect_influence_params<cell<Lyt>> &params = {}, defect_influence_stats *stats = nullptr)

This algorithm uses random sampling to find a part of the defect influence domain that might not be complete. It performs a total of `samples uniformly-distributed random samples within the specified area.

Template Parameters:

Lyt – SiDB cell-level layout type.

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

  • samples – Number of random samples to perform.

  • params – Defect influence domain computation parameters.

  • stats – Statistics.

Returns:

The (partial) defect influence domain of the layout.

template<typename Lyt, typename TT>
defect_influence_domain<Lyt> fiction::defect_influence_quicktrace(const Lyt &lyt, const std::vector<TT> &spec, const std::size_t samples, const defect_influence_params<cell<Lyt>> &params = {}, defect_influence_stats *stats = nullptr)

Applies contour tracing to identify the boundary (contour) between influencing and non-influencing defect positions for a given SiDB layout.

The algorithm leverages the concept of a screened Coulomb potential, where the electrostatic interaction weakens as distance increases. If a defect at position p causes the SiDB layout to be non-influential, then defects further away from the layout are also likely to have no influence on the layout’s functionality or performance. Conversely, defects closer to the layout may cause it to fail. This behavior allows for efficient contour tracing of the transition between influential and non-influential states.

The process is as follows:

  1. Initialization: Randomly select samples initial defect positions several nanometers away from the layout where they are unlikely to influence the layout.

  2. Contour Tracing: For each position, perform a defect-aware physical simulation to identify adjacent positions along the x-axis that influence the layout.

  3. Contour Following: Trace the contour of non-influential positions until the starting point is reached again, thereby closing the contour.

  4. Repetition: Repeat steps 1-3 for multiple initial heights to identify additional contours, since multiple influential-to-non-influential contours may exist. This process helps to detect all relevant transitions in the layout. This algorithm uses contour tracing to identify the transition between influencing and non-influencing defect positions of the SiDB layout. It starts by searching for defect locations on the left side (bounding_box + additional scanning area). The y-coordinate for these positions is chosen randomly. The number of samples is determined by the samples parameter. Then, the algorithm moves each defect position to the right, searching for the first last non-influencing defect position.

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

  • TT – Truth table type.

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

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

  • samples – Number of samples to perform.

  • params – Defect influence domain computation parameters.

  • stats – Defect influence computation statistics.

Returns:

The (partial) defect influence domain of the layout.

template<typename Lyt>
defect_influence_domain<Lyt> fiction::defect_influence_quicktrace(const Lyt &lyt, const std::size_t samples, const defect_influence_params<cell<Lyt>> &params = {}, defect_influence_stats *stats = nullptr)

Applies contour tracing to identify the boundary (contour) between influencing and non-influencing defect positions for a given SiDB layout.

The algorithm leverages the concept of a screened Coulomb potential, where the electrostatic interaction weakens as distance increases. If a defect at position p causes the SiDB layout to be non-influential, then defects further away from the layout are also likely to have no influence on the layout’s functionality or performance. Conversely, defects closer to the layout may cause it to fail. This behavior allows for efficient contour tracing of the transition between influential and non-influential states.

The process is as follows:

  1. Initialization: Randomly select samples initial defect positions several nanometers away from the layout where they are unlikely to influence the layout.

  2. Contour Tracing: For each position, perform a defect-aware physical simulation to identify adjacent positions along the x-axis that influence the layout.

  3. Contour Following: Trace the contour of non-influential positions until the starting point is reached again, thereby closing the contour.

  4. Repetition: Repeat steps 1-3 for multiple initial heights to identify additional contours, since multiple influential-to-non-influential contours may exist. This process helps to detect all relevant transitions in the layout.

Template Parameters:

Lyt – SiDB cell-level layout type.

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

  • samples – Number of samples to perform.

  • params – Defect influence domain computation parameters.

  • stats – Defect influence computation statistics.

Returns:

The (partial) defect influence domain of the layout.

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

template<typename CellType>
struct defect_clearance

Results of the defect clearance calculation.

Template Parameters:

CellType – Cell type of the layout.

Public Members

CellType defect_position = {}

Position with maximum distance to the SiDB layout at which the placement of an SiDB defect still causes the gate to fail.

double defect_clearance_distance = {}

The maximum of the minimum distances between any SiDB of the layout and the defect responsible for gate failure (unit: nm).

template<typename Lyt>
defect_clearance<cell<Lyt>> fiction::calculate_defect_clearance(const Lyt &lyt, const defect_influence_domain<Lyt> &defect_inf_domain) noexcept

Computes the defect clearance for a given SiDB layout based on a defect influence domain. The defect clearance is the maximum distance at which a defect can influence the layout. It calculates the minimum distance from each SiDB to any influential defect position.

Template Parameters:

Lyt – SiDB cell-level layout type.

Parameters:
  • lyt – SiDB layout for which the defect clearance is computed.

  • defect_inf_domain – Defect influence domain of the given SiDB layout.

Returns:

Defect clearance.