Input Parameters

The input parameters is a JSON file, which holds a single dictionary. Due to the modular nature of StructOpt, the input file is a dictionary of embedded dictionary, where dictionary key and values often determine the function name and kwargs, respectively. An example Au nanoparticle input file is given below. Each part will be discussed in the following sections.

Example:

{
  "seed": 0,
  "structure_type": "cluster",
  "generators": {
      "fcc": {"number_of_individuals": 5,
              "kwargs": {"atomlist": [["Au", 55]],
                         "orientation": "100",
                         "cell": [20, 20, 20],
                         "a": 4.08}}
  },
  "fitnesses": {
      "LAMMPS": {"weight": 1.0,
                 "kwargs": {"use_mpi4py": false,
                            "MPMD": 0,
                            "keep_files": true,
                            "min_style": "cg",
                            "min_modify": "line quadratic",
                            "minimize": "1e-8 1e-8 5000 10000",
                            "pair_style": "eam",
                            "potential_file": "$STRUCTOPT_HOME/potentials/Au_u3.eam",
                            "thermo_steps": 0,
                            "reference": {"Au": -3.930}}}
  },
  "relaxations": {
      "LAMMPS": {"order": 0,
                 "kwargs": {"use_mpi4py": false,
                            "MPMD": 0,
                            "keep_files": true,
                            "min_style": "cg",
                            "min_modify": "line quadratic",
                            "minimize": "1e-8 1e-8 5000 10000",
                            "pair_style": "eam",
                            "potential_file": "$STRUCTOPT_HOME/potentials/Au_u3.eam",
                            "thermo_steps": 0}}
  },
  "convergence": {
      "max_generations": 10
  },
  "mutations": {
      "move_surface_atoms": {"probability": 0.0},
      "move_atoms": {"probability": 0.0},
      "move_atoms_group": {"probability": 0.0},
      "rotate_atoms": {"probability": 0.0},
      "rotate_cluster": {"probability": 0.0},
      "rotate_all": {"probability": 0.0},
      "move_surface_defects": {"probability": 1.0}
  },
  "crossovers": {
      "rotate": {"probability": 0.7,
                 "kwargs": {"repair_composition": true}}
  },
  "predators": {
      "fuss": {"probability": 0.9},
      "diversify_module": {"probability": 0.1,
                           "kwargs": {"module": "LAMMPS",
                                      "min_diff": 0.01}}
  },
  "selections": {
      "rank": {"probability": 1.0}
  }
}

Global Parameters

Global parameters are those that determine how the optimizer should run.

structure_type

structure_type (str) is a key parameter for determining operations will be run. Currently, only cluster is supported, but StructOpt is written in a way that how the operations are carried out is seperated from the operations themselves. Hence, one can easily write new crossover, mutation, selection, and predator operations that are unique to their structure, test them, and incorporate them inside StructOpt seamlessly.

seed

seed (int): seed for the psuedo-random number generator. Two runs with the exact same input and seed should run exactly the same way. Almost all operations use random number generators. See should be an int.

convergence

convergence (dict): Convergence is a dictionary that determines when to stop the calculation. Currently, the only convergence criteria is max_generations, which is set to an int. For example, the setting below runs the optimizer for 200 generations.

Example:

"convergence": {
    "max_generations": 200
}

In the future, more convergence options will be added.

post_processing

post_processing (dict): Determines what is output as the optimizer is run. Currently, the only option is XYZs, which determines how frequentely the xyz files of each generation should be printed. The rules for this are as follows.

  • XYZs = 0: all generations are kept
  • XYZs > 0: every XYZs generation is kept
  • XYZs < 0: the last XYZs generations are kept

The below example is a run where only the last generation is kept. This behavior is by default and encouraged for saving space.

Example:

"post_processing": {
    "XYZs": -1
}

Generators

Generators are functions for initializing the population. These are pseudo-random generators that depend on the seed global parameter.

Generators are given as a dictionary entry defined by the generators key in the input file. The structure of the generators dictionary with N desired generators is given below.

Example:

"generators": {
    generator_1: {"number_of_individuals": n_1,
                  "kwargs": kwargs_1}
    generator_2: {"number_of_individuals": n_2,
                  "kwargs": kwargs_2}
    generator_3: {"number_of_individuals": n_3,
                  "kwargs": kwargs_3}
    ...
    generator_N: {"number_of_individuals": n_N,
                  "kwargs": kwargs_N}
}

The string for generator_i, is the name of the generator one wants to use. The number of individuals that generator should generate is determined by the integer n_i. The sum of all n_i values determines the total size of the population, which is fixed throughout the run. kwargs_i are dictionaries that input the kwargs to the generator function one is using. These will be specific to the function and can be found in their help function, show below.

structopt.cluster.individual.generators.ellipsoid(atomlist, fill_factor=0.74, radii=None, ratio=[1, 1, 1], cell=None)

Generates a random ellipsoid by rejection sampling.

Parameters:
  • atomlist (list) – A list of [sym, n] pairs where sym is the chemical symbol and n is the number of of sym’s to include in the individual
  • fill_factor (float) – Determines how “close” the atoms are packed. See structopt.tools.get_particle_radius for description
  • radii (list) – The size, in angstroms, of the ellipsoid in the x, y and z direction. If None, with ratio parameters and the average atomic radii, radii is automatically calculated.
  • ratio (list) – The ratio of the dimensions of the ellipsoid in the x, y and z direction.
  • cell (list) – The size, in angstroms, of the dimensions that holds the atoms object. Must be an orthogonal box.
structopt.cluster.individual.generators.sphere(atomlist, cell, fill_factor=0.74, radius=None)

Generates a random sphere of particles given an atomlist and radius.

Parameters:
  • atomlist (list) – A list of [sym, n] pairs where sym is the chemical symbol and n is the number of of sym’s to include in the individual
  • cell (list) – The x, y, and z dimensions of the cell which holds the particle.
  • fill_factor (float) – How densely packed the sphere should be. Ranges from 0 to 1.
  • radius (float) – The radius of the sphere. If None, estimated from the atomic radii
structopt.cluster.individual.generators.fcc(atomlist, cell, a, shape=[1, 1, 1], orientation=None, size=21, roundness=0.5, alpha=10, v=None, angle=None)

Generates a fcc nanoparticle of varying shape and orientation. For multi-component particles, elements are randomly distributed.

Parameters:
  • atomlist (list) – A list of [sym, n] pairs where sym is the chemical symbol and n is the number of of sym’s to include in the individual
  • cell (list) – A 3 element list that defines the x, y, and z dimensions of the simulation box
  • a (float) – fcc lattice constant
  • shape (list) – The ratio of the x, y, and z dimensions of the particle.
  • orientation (str) – The facet that is parallel to the xy plane. This is useful for LAMMPS+STEM calculations where one already knows the orientation. If None, a random orientation is chosen.
  • size (int) – Size of the fcc grid for building the nanoparticle. For example, a size of 21 means 21 x 21 x 21 supercell of fcc primitive cells will be used. Note, if the size is too small, the function will automatically expand the cell to fit the particle.
  • roundness (float) – Determines the “roundness” of the particle. A more round particle will have a smaller surface area to volume ratio, but more undercoordinated surface sites. A less round particle will take more and more the form of an octahedron.
  • alpha (float) – Parameter for determining how defective the particle will be. Higher alpha means less defective particle.
  • v (list) – Used for a custom orientation. V is the vector in which to rotate the particle. Requires angle parameters to be entered. All rotations are done with respect to the 100 plane.
  • angle (float) – Angle, in radians, to rotate atoms around vector v. All rotations are done with respect to the 100 plane.

Crossovers

Crossovers are operations for mating two individuals chosen by a selection algorithm. The purpose of the crossover is to intelligently combine (mate) different individuals (parents) in a way to create new individuals (children) that have the features of the current best individuals in the population.

Crossovers are given as a dictionary entry defined by the crossovers key in the input file. The structure of the crossovers dictionary with N desired selections is given below.

Example:

"crossovers": {
    crossover_1: {"probability": p_1,
                  "kwargs": kwargs_1}
    crossover_2: {"probability": p_2,
                  "kwargs": kwargs_2}
    crossover_3: {"probability": p_3,
                  "kwargs": kwargs_3}
    ...
    crossover_N: {"probability": p_N,
                  "kwargs": kwargs_N}
}

The string for crossover_i, is the name of the crossover one wants to use. The probability p_i is the probability of the crossover occuring if a mate is determined to happen in the population. p_i values should sum to 1. kwargs_i are dictionaries that input the kwargs to the crossover function one is using. These will be specific to the function and can be found in their help function.

Currently the only crossover in use in the algorithm is the cut-and-splice operator introduced by Deaven and Ho. The description is shown below.

structopt.cluster.population.crossovers.rotate(individual1, individual2, center_at_atom=True, repair_composition=True)

Rotates the two individuals around their centers of mass, splits them in half at the xy-plane, then splices them together. Maintains number of atoms. Note, both individuals are rotated in the same way.

Parameters:
  • individual1 (Individual) – The first parent
  • individual2 (Individual) – The second parent
  • center_at_atom (bool) – This centers the cut at an atom. This is particularly useful when one desires a crystalline solution. If both parents are crystalline, the children will likely not have grain boundaries.
  • repair_composition (bool) – If True, conserves composition. For crossovers that create children with more (less) atoms, atoms are taken from (added to) the surface of the particle. For incorrect atomic ratios, atomic symbols are randomly interchanged throughout the particle
Returns:

  • Individual (The first child)
  • Individual (The second child)

Selections

Selections are operations for choosing which individuals to “mate” in producing new individuals. Individuals are chosen based on their fitness, and different selection functions determine how diverse the children should be.

Selections are given as a dictionary entry defined by the selections key in the input file. The structure of the selections dictionary with N desired selections is given below.

Example:

"selections": {
    selection_1: {"probability": p_1,
                  "kwargs": kwargs_1}
    selection_2: {"probability": p_2,
                  "kwargs": kwargs_2}
    selection_3: {"probability": p_3,
                  "kwargs": kwargs_3}
    ...
    selection_N: {"probability": p_N,
                  "kwargs": kwargs_N}
}

The string for selection_i, is the name of the selection one wants to use. The probability p_i is the probability of the selection occuring if a mate is determined to happen in the population. p_i values should sum to 1. kwargs_i are dictionaries that input the kwargs to the selection function one is using. These will be specific to the function and can be found in their help function.

structopt.common.population.selections.tournament(population, fits, tournament_size=5, unique_pairs=False, unique_parents=False, keep_best=False)

Selects pairs in seperate “tournaments”, where a subset of the population are randomly selected and the highest fitness allowed to pass. In addition to a population, their fits, and end population size, takes in a tournament size parameter.

Parameters:
  • population (Population) – The population of individuals needed to be trimmed
  • fits (list) – List of fitnesses that correspond to the population.
  • tournament_size (int) – The number of individuals in each tournament. If 1, tournament is the same as random selection. If len(population), corresponds to the “best” selection process
  • unique_pairs (bool) – If True, all combinations of parents are unique, though parents can show up in different pairs. True increases the diversity of the population.
  • unique_parents (bool) – If True, all parents can only mate with on other individual. True increases the diversity of the population.
structopt.common.population.selections.random_selection(population, fits)

Randomly selects parents

Parameters:
  • population (Population) – An population of individuals
  • fits (list) – Fitnesses that corresponds to population
structopt.common.population.selections.best(population, fits)

Deterministic selection function that chooses adjacently ranked individuals as pairs.

Parameters:
  • population (Population) – An population of individuals
  • fits (list) – Fitnesses that corresponds to population
structopt.common.population.selections.rank(population, fits, p_min=None, unique_pairs=False, unique_parents=False)

Selection function that chooses pairs of structures based on linear ranking.

See “Grefenstette and Baker 1989 Whitley 1989”.

Parameters:
  • population (Population) – An object inherited from list that contains StructOpt individual objects.
  • fits (list) – A list of fitnesses of the population
  • p_min (float) – The probability of choosing the lowest ranked individual. Given population of size N, this should be below 1/nindiv. The probability of selecting rank N (worst) to rank 1 (best) increases from p_min to (2/N - p_min) in even, (1/N - p_min) increments. Defaults to (1/N)^2.
  • unique_pairs (bool) – If True, all combinations of parents are unique. True increases the diveristy of the population.
  • unique_parents (bool) – If True, all parents can only mate with on other individual. True increases the diversity of the population.
structopt.common.population.selections.roulette(population, fits, unique_pairs=False, unique_parents=False)

Selection function that chooses pairs of structures based on their fitness. Fitnesses are normalized from 0 to 1.

See “Grefenstette and Baker 1989 Whitley 1989”.

Parameters:
  • population (StructOpt population object) – An object inherited from list that contains StructOpt individual objects.
  • fits (list) – A list of fitnesses of the population
  • unique_pairs (bool) – If True, all combinations of parents are unique. True increases the diveristy of the population.
  • unique_parents (bool) – If True, all parents can only mate with on other individual. True increases the diversity of the population.

Predators

Similar to selections, predators are selection processes that selects individuals based on their fitness. The distinction is that while selections select individuals with positive features to duplicate in children, predators select which individuals to keep in the next generation. Note, this must be done because crossovers and (sometimes) mutations increase the population every generation, and hence each generation requires a predator step.

Predators are given as a dictionary entry defined by the predators key in the input file. The structure of the predators dictionary with N desired predators is given below

Example:

"predators": {
    predator_1: {"probability": p_1,
                 "kwargs": kwargs_1}
    predator_2: {"probability": p_2,
                 "kwargs": kwargs_2}
    predator_3: {"probability": p_3,
                 "kwargs": kwargs_3}
    ...
    predator_N: {"probability": p_N,
                 "kwargs": kwargs_N}
}

The string for predator_i, is the name of the predator one wants to use. The probability p_i is the probability of the predator occuring on every individual in the population. p_i values should sum to 1. kwargs_i are dictionaries that input the kwargs to the predator function one is using. These will be specific to the function and can be found in their help function.

structopt.common.population.predators.best(fits, nkeep)

Sorts individuals by fitness and keeps the top nkeep fitnesses.

Parameters:
  • fits (dict<int, float>) – Dictionary of <individual.id, fitness> pairs.
  • nkeep (int) – The number of individuals to keep. In a GA run, corresponds to the sum of each generators number_of_individuals
structopt.common.population.predators.roulette(fits, nkeep, T=None)

Select individuals with a probability proportional to their fitness. Fitnesses are renormalized from 0 - 1, which means minimum fitness individual is never included in in the new population.

Parameters:
  • fits (dict<int, float>) – Dictionary of <individual.id, fitness> pairs.
  • nkeep (int) – The number of individuals to keep. In a GA run, corresponds to the sum of each generators number_of_individuals
  • T (float) – If T is not None, a boltzman-like transformation is applied to all fitness values with T.
structopt.common.population.predators.rank(fits, nkeep, p_min=None)

Selection function that chooses pairs of structures based on linear ranking.

See “Grefenstette and Baker 1989 Whitley 1989”.

Parameters:
  • fits (dict<int, float>) – Dictionary of <individual.id, fitness> pairs.
  • nkeep (int) – The number of individuals to keep. In a GA run, corresponds to the sum of each generators number_of_individuals
  • p_min (float) – The probability of choosing the lowest ranked individual. Given population of size N, this should be below 1/nindiv. The probability of selecting rank N (worst) to rank 1 (best) increases from p_min to (2/N - p_min) in even, (1/N - p_min) increments. Defaults to (1/N)^2.
structopt.common.population.predators.fuss(fits, nkeep, nbest=0, fusslimit=10)

Fixed uniform selection scheme. Aimed at maintaining diversity in the population. In the case where low fit is the highest fitness, selects a fitness between min(fits) and min(fits) + fusslimit, if the difference between the min(fit) and max(fit) is larger than fusslimit.

Parameters:
  • fits (dict<int, float>) – Dictionary of <individual.id, fitness> pairs.
  • nkeep (int) – The number of individuals to keep. In a GA run, corresponds to the sum of each generators number_of_individuals
  • nbest (int) – The top n individuals to always keep (default 0)
  • fusslimit (float) – Individuals that have fitness fusslimit worse than the max fitness will not be considered
structopt.common.population.predators.tournament(fits, nkeep, tournament_size=5)

Selects individuals in seperate “tournaments”, where a subset of the population are randomly selected and the highest fitness allowed to pass. In addition to a population, their fits, and end population size, takes in a tournament size parameter.

Parameters:
  • fits (dict<int, float>) – Dictionary of <individual.id, fitness> pairs.
  • nkeep (int) – The number of individuals to keep. In a GA run, corresponds to the sum of each generators number_of_individuals
  • tournament_size (int) – The number of individuals in each tournament. If 1, tournament is the same as random selection. If len(population), corresponds to the “best” selection process

Mutations

Mutations are operations applied to individuals that change its structure and composition. It is a local search operation, though the mutation itself can be written to perform small or larger changes.

Mutations are given as a dictionary entry defined by the mutations key in the input file. The structure of the mutations dictionary with N desired mutations is given below

Example:

"mutations": {
    "preserve_best": "true" or "false",
    "keep_original": "true" or "false",
    "keep_original_best": "true" or "false,
    mutation_1: {"probability": p_1,
                 "kwargs": kwargs_1}
    mutation_2: {"probability": p_2,
                 "kwargs": kwargs_2}
    mutation_3: {"probability": p_3,
                 "kwargs": kwargs_3}
    ...
    mutation_N: {"probability": p_N,
                 "kwargs": kwargs_N}
}

The string for mutation_i, is the name of the mutation one wants to use. The probability p_i is the probability of the mutation occuring on every individual in the population. p_i values should sum to any value between 0 and 1. kwargs_i are dictionaries that input the kwargs to the mutation function one is using. These will be specific to the function and can be found in their help function.

In addition to specifying the mutations you want to use, the mutations dictionary takes three special kwargs: preserve_best, keep_original, and keep_original_best. Setting preserve_best to true, means the highest fitness individual will never be mutated. Setting keep_original to true means mutations will be applied to copies of individuals, not the individual itself. This means, the original individual is not changed through a mutation. keep_original_best applies keep_original to only the best individual.

The currently implemented mutations are shown below. Note in all functions, the first argument is the atomic structure, which inserted by the optimizer. The user defines all of the other kwargs after the first input.

structopt.common.individual.mutations.swap_positions(individual, max_natoms=0.2)

Randomly swaps the positions atoms within the individual (in place).

Parameters:
  • individual (Individual) – an individual
  • max_natoms (float or int) – if float, the maximum number of atoms whose positions will be swapped is max_natoms*len(individual) if int, the maximum number of atoms whose positions will be swapped is max_natoms if the number of atoms to be swapped is (or evaluates to) an odd integer, it is rounded down to an even integer max_natoms corresponds to the maximum number of atoms whose positions will change default: 0.20
structopt.common.individual.mutations.swap_species(individual, max_natoms=0.2)

Randomly swaps the species of atoms within the individual (in place).

Parameters:
  • individual (Individual) – an individual
  • max_natoms (float or int) – if float, the maximum number of atoms that will be swapped is max_natoms*len(individual) if int, the maximum number of atoms that will be swapped is max_natoms if the number of atoms to be swapped is (or evaluates to) an odd integer, it is rounded down to an even integer max_natoms corresponds to the maximum number of atoms whose species will change default: 0.20
structopt.common.individual.mutations.rotate_atoms(individual, max_natoms=0.2)

Randomly rotates a number of random atoms within the individual (in place).

Parameters:
  • individual (Individual) – an individual
  • max_natoms (float or int) – if float, the maximum number of atoms that will be rotated is max_natoms*len(individual) if int, the maximum number of atoms that will be rotated is max_natoms default: 0.20
structopt.common.individual.mutations.rotate_all(atoms, vector=None, angle=None, center=None)

Rotate all atoms around a single point. Most suitable for cluster calculations.

Parameters:
  • individual (Individual) – An individual.
  • vector (string or list) – The list of axes in which to rotate the atoms around. If None, is a randomly chosen direction. If ‘random’ in list, a random vector can be chosen.
  • angle (string or list) – A list of angles that will be chosen to rotate. If None, is randomly generated. Angle must be given in radians. If ‘random’ in list, a random angle is included.
  • center (string or xyz iterable) – The center in which to rotate the atoms around. If None, defaults to center of mass. Acceptable strings are COM = center of mass COP = center of positions COU = center of cell
structopt.common.individual.mutations.permutation(individual)

Swaps the chemical symbol between two elements

Parameters:individual (Individual) – An individual or atoms object.
structopt.common.individual.mutations.rattle(individual, stdev=0.5, x_avg_bond=True)

Randomly displace all atoms in a random direction with a magnitude drawn from a gaussian distribution.

Parameters:
  • individual (Individual) – An individual
  • stdev (float) – The standard deviation of the gaussian distribution to rattle all the atoms. If x_avg_bond is set to True, given as the fraction of the average bond length of the material.
  • x_avg_bond (bool) – If True, the gaussian distributions standard deviation is stdev * avg_bond_length. Note, this only applies to fcc, hcp, or bcc materials.
structopt.cluster.individual.mutations.move_atoms(individual, max_natoms=0.2)

Randomly moves atoms within a cluster.

Parameters:
  • individual (Individual) – An individual
  • max_natoms (float or int) – if float, the maximum number of atoms that will be moved is max_natoms*len(individual). if int, the maximum number of atoms that will be moved is max_natoms default: 0.20
structopt.cluster.individual.mutations.move_surface_atoms(individual, max_natoms=0.2, move_CN=11, surf_CN=11)

Randomly moves atoms at the surface to other surface sites

Parameters:
  • individual (Individual) – The individual object to be modified in place
  • max_natoms (float or int) – if float, the maximum number of atoms that will be moved is max_natoms*len(individual) if int, the maximum number of atoms that will be moved is max_natoms default: 0.20
  • move_CN (int) – The coordination number to determine which atoms can move moved. Any atom with coordination number above move_CN will not be moved
  • surf_CN (int) – The coordination number to determine which atoms are considered surface atoms. Surface atoms are used to estimating new surface sites
structopt.cluster.individual.mutations.rotate_cluster(individual, max_natoms=0.2)

Chooses a random number of atoms nearest to a random point in the cluster. These atoms are then rotated randomly around this point

Parameters:
  • individual (Individual) – An individual object
  • max_natoms (float) – The fraction of the total atoms to rotate
structopt.cluster.individual.mutations.twist(individual, max_radius=0.9)

Splits the particle randomly in half and rotates one half.

Parameters:
  • individual (structopt.Individual object) – Individual to be mutated
  • max_natoms (float) – That maximum relative distance from the center of the particle the twist is initiated
structopt.cluster.individual.mutations.swap_core_shell(individual, surf_CN=11)

Swaps atoms on the surface with an atom in the core. Only does it for different element types.

Parameters:
  • individual (Individual) – An individual
  • surf_CN (int) – The maximum coordination number of an atom to be considered surface
structopt.cluster.individual.mutations.rich2poor(individual)

Used for multi-component systems. Swaps atoms A and B so that atom A moves from a region with a high number of A-A bonds to a low number of A-A bonds.

Parameters:individual (Individual) – An individual
structopt.cluster.individual.mutations.poor2rich(individual)

Used for multi-component systems. Swaps atoms A and B so that atom A moves from a region with a low number of A-A bonds to a high number of A-A bonds.

Parameters:individual (Individual) – An individual
structopt.cluster.individual.mutations.move_surface_defects(individual, surf_CN=11)

Moves atoms around on the surface based on coordination number Moves a surface atom with a low CN to an atom with a high CN

Parameters:
  • individual (Individual) – The individual object to be modified in place
  • surf_CN (int) – The maximum coordination number to considered a surface atom
structopt.cluster.individual.mutations.enrich_surface(individual, surf_CN=11, species=None)

Mutation that selectively enriches the surface with a species.

Parameters:
  • individual (Individual) – An individual
  • surf_CN (int) – The maximum coordination number of an atom to be considered surface
  • species (str) – The surface to enrich with. If None, takes the lowest concentration species
structopt.cluster.individual.mutations.enrich_bulk(individual, surf_CN=11, species=None)

Mutation that selectively enriches the bulk with a species

Parameters:
  • individual (Individual) – An individual
  • surf_CN (int) – The maximum coordination number of an atom to be considered surface
  • species (str) – The surface to enrich with. If None, takes the lowest concentration species
structopt.cluster.individual.mutations.enrich_surface_defects(individual, surf_CN=11, species=None)

Mutation that selectively enriches defects with a species. Defects are defined as atoms atoms with lower coordination numbers

Parameters:
  • individual (Individual) – An individual
  • surf_CN (int) – The maximum coordination number of an atom to be considered surface
  • species (str) – The surface to enrich with. If None, takes the lowest concentration
structopt.cluster.individual.mutations.enrich_surface_facets(individual, surf_CN=11, species=None)

Mutation that selectively enriches facets with a species. Facets are defined as atoms atoms with higher coordination numbers.

Parameters:
  • surf_CN (int) – The maximum coordination number of an atom to be considered surface
  • species (str) – The surface to enrich with. If None, takes the lowest concentration

Relaxations

Relaxations performs a local relaxation to the atomic structure before evaluating their fitness. This is typically done after crossover and mutation operators are applied.

Relaxations differ than the previous operations in that they require varying amounts of resources. Hence, a subsequent section, Parallelization, will introduce ways to run your job with varying levels of parallel performance.

Relaxations are given as a dictionary entry defined by the relaxations key in the input file. The structure of these dictionaries is shown below.

Example:

"relaxations": {
    relaxation_1: {"order": o_1,
                   "kwargs": kwargs_1}
    relaxation_2: {"order": o_2,
                   "kwargs": kwargs_2}
    relaxation_3: {"order": o_3,
                   "kwargs": kwargs_3}
    ...
    relaxation_N: {"order": o_N,
                   "kwargs": kwargs_N}
}

The string for relaxation_i, is the name of the relaxation one wants to use. The order o_i is the order of the relaxation occuring on every individual in the population. kwargs_i are dictionaries that input the kwargs to the relaxation function one is using. These will be specific to the function. More details of each relaxation module will be given in the following subsections

LAMMPS

The LAMMPS relaxation module calls LAMMPS to relax according to some potential. Most of the kwargs can be found from the LAMMPS documentation.

class structopt.common.individual.relaxations.LAMMPS(parameters)

LAMMPS class for running LAMMPS on a single individual. Takes a dictionary, where the key: value are the parameters for running LAMMPs.

Parameters:
  • min_style (str) – The minimization scheme for running LAMMPS. See LAMMPS doc.
  • min_modify (str) – Parameters for min_style energy minimization algorithm. See LAMMPS doc.
  • minimize (str) – Convergence criteria for minimization algorithm. See LAMMPS doc.
  • pair_style (str) – Type of potential used. See LAMMPS doc.
  • potential_file (str) – The path to the potential_file. Should be absolute.
  • thermo_steps (int) – How much output to print of thermodynamic information. If set to 0, only the last step is printed.See LAMMPS doc.
  • keep_file (bool) – Will keep all of the LAMMPS input and output files for each individual. Use with caution.
  • repair (bool) – Determines whether to run an algorithm to make sure no atoms are in “space”. Atoms can be in space due to a mutation or crossover that results in a large force that shoots the atom outside of the particle.

The potential files available to use are listed below and are from the default potentials included from LAMMPS. Given a potential, enter in the potential_file kwarg as “$STRUCTOPT_HOME/potentials/<name>”. Note also that different potentials will have different lines of the pair_style kwarg. If the user would like to use an unavailable potential file, please submit an email to zxu39@wisc.edu, and the potential will be added.

AlCu.eam.alloy: Aluminum and copper alloy EAM (Cai and Ye, Phys Rev B, 54, 8398-8410 (1996))

Au_u3.eam: Gold EAM (SM Foiles et al, PRB, 33, 7983 (1986))

ZrCuAl2011.eam.alloy: Zirconium, copper, and aluminum glass (Howard Sheng at GMU. (hsheng@gmu.edu))

Fitnesses

Fitness evaluates how closely the individual satisfies the minimization criteria. One typical minimization criteria is the stability of a structure, and the formation energy is the fitness. Note, all fitness modules operate so that the lower the fitness value the more fit it is.

Fitnesses differ than the previous operations in that they require varying amounts of resources. Hence, a subsequent section, Parallelization, will introduce ways to run your job with varying levels of parallel performance.

Fitnesses are given as a dictionary entry defined by fitnesses key in the input file. The structure of these dictionaries is shown below.

Example:

"fitnesses": {
    fitness_1: {"weight": w_1,
                 "kwargs": kwargs_1}
    fitness_2: {"weight": w_2,
                 "kwargs": kwargs_2}
    fitness_3: {"weight": w_3,
                 "kwargs": kwargs_3}
    ...
    fitness_N: {"weight": w_N,
                 "kwargs": kwargs_N}
}

The string for fitness_i, is the name of the fitness one wants to use. The weight w_i is the constant to multiply the fitness value returned by the fitness_i module. Not that all selections and predators operate on the total fitness, which is a sum of each fitness and their weight. kwargs_i are dictionaries that input the kwargs to the fitness function one is using. These will be specific to the function. More details of each fitness module will be given in the following subsections

LAMMPS

The LAMMPS fitness module calls LAMMPS to calculate the potential energy of the structure. Most of the kwargs can be found from the LAMMPS documentation. In addition, most of the kwargs are the same as relaxations, except the fitness module of LAMMPS has a number of normalization options for returning the potential energy. These are described below.

class structopt.common.individual.fitnesses.LAMMPS(parameters)

LAMMPS class for running LAMMPS on a single individual. Takes a dictionary, where the key: value are the parameters for running LAMMPs.

Parameters:
  • min_style (str) – The minimization scheme for running LAMMPS. See LAMMPS doc.
  • min_modify (str) – Parameters for min_style energy minimization algorithm. See LAMMPS doc.
  • minimize (str) – Convergence criteria for minimization algorithm. Note for fitness values, the last two values are set to 0, so no relaxation is done. See LAMMPS doc.
  • pair_style (str) – Type of potential used. See LAMMPS doc.
  • potential_file (str) – The path to the potential_file. Should be absolute.
  • thermo_steps (int) – How much output to print of thermodynamic information. If set to 0, only the last step is printed.See LAMMPS doc.
  • keep_file (bool) – Will keep all of the LAMMPS input and output files for each individual. Use with caution.
  • reference (dict) – Reference energies of the particle. These are values to subtract from the values returned by LAMMPS. Given as a dictionary of {sym : E} pairs, where sym is a str denoating the the element, while E is the value to be subtracted per sym. This is typically the pure component formation energy calculated with LAMMPS. Note since this is merely a fixed subtraction, should not change the performance in constant composition runs.

The potential files available to use are listed below and are from the default potentials included from LAMMPS. Given a potential, enter in the potential_file kwarg as “$STRUCTOPT_HOME/potentials/<name>”. Note also that different potentials will have different lines of the pair_style kwarg. If the user would like to use an unavailable potential file, please submit an email to zxu39@wisc.edu, and the potential will be added.

AlCu.eam.alloy: Aluminum and copper alloy EAM (Cai and Ye, Phys Rev B, 54, 8398-8410 (1996))

Au_u3.eam: Gold EAM (SM Foiles et al, PRB, 33, 7983 (1986))

ZrCuAl2011.eam.alloy: Zirconium, copper, and aluminum glass (Howard Sheng at GMU. (hsheng@gmu.edu))

Parallelization

In addition to the module-specific parameters, each module requires two parallelization entries: use_mpi4py and MPMD_cores_per_structure. These two entries are mutually exclusive, meaning that only one can be turned on at a time. use_mpi4py can take two values, true or false depending on whether the module should use the `one-structure-per-core <>`_ parallelization.

MPMD_cores_per_structure can be disabled (if use_mpi4py is true) by setting it to 0, but otherwise specifies the number of cores that each process/structure should be allocated within the MPI_Comm_spawn_multiple command. There are two types of valid values for this parameter: 1) an integer specifying the number of cores per structure, or 2) a string of two integers separated by a dash specifying the minimum and maximum number of cores allowed (e.g. "4-16"). MPMD_cores_per_structure can also take the value of "any", and StructOpt will use as many cores as it can to run each individual.

Example:

"relaxations": {
    "LAMMPS": {"order": 0,
               "kwargs": {"use_mpi4py": true,
                          "MPMD_cores_per_structure": 0,
                          "keep_files": true,
                          "min_style": "cg\nmin_modify line quadratic",
                          "minimize": "1e-8 1e-8 5000 10000",
                          "pair_style": "eam/alloy",
                          "potential_file": "$STRUCTOPT_HOME/potentials/ZrCuAl2011.eam.alloy",
                          "thermo_steps": 1000}}
}