pyMBE.pyMBE

   1#
   2# Copyright (C) 2023-2026 pyMBE-dev team
   3#
   4# This file is part of pyMBE.
   5#
   6# pyMBE is free software: you can redistribute it and/or modify
   7# it under the terms of the GNU General Public License as published by
   8# the Free Software Foundation, either version 3 of the License, or
   9# (at your option) any later version.
  10#
  11# pyMBE is distributed in the hope that it will be useful,
  12# but WITHOUT ANY WARRANTY; without even the implied warranty of
  13# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  14# GNU General Public License for more details.
  15#
  16# You should have received a copy of the GNU General Public License
  17# along with this program.  If not, see <http://www.gnu.org/licenses/>.
  18
  19import re
  20import json
  21import pint
  22import numpy as np
  23import pandas as pd
  24import scipy.constants
  25import scipy.optimize
  26import logging
  27import importlib.resources
  28
  29# Database
  30from pyMBE.storage.manager import Manager
  31from pyMBE.storage.pint_quantity import PintQuantity
  32## Templates
  33from pyMBE.storage.templates.particle import ParticleTemplate, ParticleStateTemplate
  34from pyMBE.storage.templates.residue import ResidueTemplate
  35from pyMBE.storage.templates.molecule import MoleculeTemplate
  36from pyMBE.storage.templates.peptide import PeptideTemplate
  37from pyMBE.storage.templates.protein import ProteinTemplate
  38from pyMBE.storage.templates.hydrogel import HydrogelTemplate, HydrogelNode, HydrogelChain
  39from pyMBE.storage.templates.bond import BondTemplate
  40from pyMBE.storage.templates.angle import AngleTemplate
  41from pyMBE.storage.templates.lj import LJInteractionTemplate
  42## Instances
  43from pyMBE.storage.instances.particle import ParticleInstance
  44from pyMBE.storage.instances.residue import ResidueInstance
  45from pyMBE.storage.instances.molecule import MoleculeInstance
  46from pyMBE.storage.instances.peptide import PeptideInstance
  47from pyMBE.storage.instances.protein import ProteinInstance
  48from pyMBE.storage.instances.bond import BondInstance
  49from pyMBE.storage.instances.angle import AngleInstance
  50from pyMBE.storage.instances.hydrogel import HydrogelInstance
  51## Reactions
  52from pyMBE.storage.reactions.reaction import Reaction, ReactionParticipant
  53# Utilities
  54import pyMBE.lib.handy_functions as hf
  55import pyMBE.storage.io as io
  56
  57class pymbe_library():
  58    """
  59    Core library of the Molecular Builder for ESPResSo (pyMBE).
  60
  61    Attributes:
  62        N_A ('pint.Quantity'):
  63            Avogadro number.
  64
  65        kB ('pint.Quantity'):
  66            Boltzmann constant.
  67
  68        e ('pint.Quantity'):
  69            Elementary charge.
  70
  71        kT ('pint.Quantity'):
  72            Thermal energy corresponding to the set temperature.
  73
  74        Kw ('pint.Quantity'):
  75            Ionic product of water, used in G-RxMC and Donnan-related calculations.
  76
  77        db ('Manager'):
  78            Database manager holding all pyMBE templates, instances and reactions.
  79
  80        rng ('numpy.random.Generator'):
  81            Random number generator initialized with the provided seed.
  82
  83        units ('pint.UnitRegistry'):
  84            Pint unit registry used for unit-aware calculations.
  85
  86        lattice_builder ('pyMBE.lib.lattice.LatticeBuilder'):
  87            Optional lattice builder object (initialized as ''None'').
  88            
  89        root ('importlib.resources.abc.Traversable'):
  90            Root path to the pyMBE package resources.
  91    """
  92
  93    def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None):
  94        """
  95        Initializes the pyMBE library.
  96
  97        Args:
  98            seed ('int'):
  99                Seed for the random number generator.
 100
 101            temperature ('pint.Quantity', optional):
 102                Simulation temperature. If ''None'', defaults to 298.15 K.
 103
 104            unit_length ('pint.Quantity', optional):
 105                Reference length for reduced units. If ''None'', defaults to
 106                0.355 nm.
 107
 108            unit_charge ('pint.Quantity', optional):
 109                Reference charge for reduced units. If ''None'', defaults to
 110                one elementary charge.
 111
 112            Kw ('pint.Quantity', optional):
 113                Ionic product of water (typically in mol²/L²). If ''None'',
 114                defaults to 1e-14 mol²/L².
 115        """
 116        # Seed and RNG
 117        self.seed=seed
 118        self.rng = np.random.default_rng(seed)
 119        self.units=pint.UnitRegistry()
 120        self.N_A=scipy.constants.N_A / self.units.mol
 121        self.kB=scipy.constants.k * self.units.J / self.units.K
 122        self.e=scipy.constants.e * self.units.C
 123        self.set_reduced_units(unit_length=unit_length, 
 124                               unit_charge=unit_charge,
 125                               temperature=temperature, 
 126                               Kw=Kw)
 127        
 128        self.db = Manager(units=self.units)
 129        self.lattice_builder = None
 130        self.root = importlib.resources.files(__package__)
 131
 132    def _check_bond_inputs(self, bond_type, bond_parameters):
 133        """
 134        Checks that the input bond parameters are valid within the current pyMBE implementation.
 135
 136        Args:
 137            bond_type ('str'): 
 138                label to identify the potential to model the bond.
 139            
 140            bond_parameters ('dict'): 
 141                parameters of the potential of the bond.
 142        """
 143        valid_bond_types   = ["harmonic", "FENE"] 
 144        if bond_type not in valid_bond_types:
 145            raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}")
 146        required_parameters = {"harmonic": ["r_0","k"],
 147                                "FENE": ["r_0","k","d_r_max"]}
 148        for required_parameter in required_parameters[bond_type]:
 149            if required_parameter not in bond_parameters.keys():
 150                raise ValueError(f"Missing required parameter {required_parameter} for {bond_type} bond")
 151            
 152    def _check_dimensionality(self, variable, expected_dimensionality):
 153        """
 154        Checks if the dimensionality of 'variable' matches 'expected_dimensionality'.
 155
 156        Args:
 157            variable ('pint.Quantity'): 
 158                Quantity to be checked.
 159
 160            expected_dimensionality ('str'): 
 161                Expected dimension of the variable.
 162
 163        Returns:
 164            ('bool'): 
 165                'True' if the variable if of the expected dimensionality, 'False' otherwise.
 166
 167        Notes:
 168            - 'expected_dimensionality' takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality).
 169            - For example, to check for a variable corresponding to a velocity 'expected_dimensionality = "[length]/[time]"'    
 170        """
 171        correct_dimensionality=variable.check(f"{expected_dimensionality}")      
 172        if not correct_dimensionality:
 173            raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}")
 174        return correct_dimensionality   
 175            
 176    def _check_pka_set(self, pka_set):
 177        """
 178        Checks that 'pka_set' has the formatting expected by pyMBE.
 179       
 180        Args:
 181            pka_set ('dict'): 
 182                {"name" : {"pka_value": pka, "acidity": acidity}}
 183        """
 184        required_keys=['pka_value','acidity']
 185        for required_key in required_keys:
 186            for pka_name, pka_entry in pka_set.items():
 187                if required_key not in pka_entry:
 188                    raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")')
 189        return
 190
 191    def _create_espresso_bond_instance(self, bond_type, bond_parameters):
 192        """
 193        Creates an ESPResSo bond instance.
 194
 195        Args:
 196            bond_type ('str'): 
 197                label to identify the potential to model the bond.
 198
 199            bond_parameters ('dict'): 
 200                parameters of the potential of the bond.
 201
 202        Notes:
 203            Currently, only HARMONIC and FENE bonds are supported.
 204
 205            For a HARMONIC bond the dictionary must contain:
 206                - k ('Pint.Quantity')      : Magnitude of the bond. It should have units of energy/length**2 
 207                using the 'pmb.units' UnitRegistry.
 208                - r_0 ('Pint.Quantity')    : Equilibrium bond length. It should have units of length using 
 209                the 'pmb.units' UnitRegistry.
 210           
 211            For a FENE bond the dictionary must additionally contain:
 212                - d_r_max ('Pint.Quantity'): Maximal stretching length for FENE. It should have 
 213                units of length using the 'pmb.units' UnitRegistry. Default 'None'.
 214
 215        Returns:
 216            ('espressomd.interactions'): instance of an ESPResSo bond object
 217        """
 218        from espressomd import interactions
 219        self._check_bond_inputs(bond_parameters=bond_parameters,
 220                                bond_type=bond_type)
 221        if bond_type == 'harmonic':
 222            bond_instance = interactions.HarmonicBond(k = bond_parameters["k"].m_as("reduced_energy/reduced_length**2"),
 223                                                      r_0 = bond_parameters["r_0"].m_as("reduced_length"))
 224        elif bond_type == 'FENE':
 225            bond_instance    = interactions.FeneBond(k = bond_parameters["k"].m_as("reduced_energy/reduced_length**2"),
 226                                                      r_0 = bond_parameters["r_0"].m_as("reduced_length"),
 227                                                      d_r_max = bond_parameters["d_r_max"].m_as("reduced_length"))    
 228        return bond_instance
 229
 230    def _create_hydrogel_chain(self, hydrogel_chain, nodes, espresso_system, use_default_bond=False, gen_angle=False):
 231        """
 232        Creates a chain between two nodes of a hydrogel.
 233
 234        Args:
 235            hydrogel_chain ('HydrogelChain'): 
 236                template of a hydrogel chain
 237            nodes ('dict'): 
 238                {node_index: {"name": node_particle_name, "pos": node_position, "id": node_particle_instance_id}}
 239
 240            espresso_system ('espressomd.system.System'): 
 241                ESPResSo system object where the hydrogel chain will be created.
 242
 243            use_default_bond ('bool', optional): 
 244                If True, use a default bond template if no specific template exists. Defaults to False.
 245
 246            gen_angle ('bool', optional):
 247                If True, generate the angle potentials internal to the created
 248                chain molecule. Junction angles near the hydrogel crosslinkers
 249                are handled separately at the hydrogel level.
 250
 251        Return:
 252            ('int'): 
 253                molecule_id of the created hydrogel chian.
 254
 255        Notes:
 256            - If the chain is defined between node_start = ''[0 0 0]'' and node_end = ''[1 1 1]'', the chain will be placed between these two nodes.
 257            - The chain will be placed in the direction of the vector between 'node_start' and 'node_end'. 
 258        """
 259        if self.lattice_builder is None:
 260            raise ValueError("LatticeBuilder is not initialized. Use 'initialize_lattice_builder' first.")
 261        molecule_tpl = self.db.get_template(pmb_type="molecule",
 262                                            name=hydrogel_chain.molecule_name)
 263        residue_list = molecule_tpl.residue_list
 264        molecule_name = molecule_tpl.name
 265        node_start = hydrogel_chain.node_start
 266        node_end = hydrogel_chain.node_end
 267        node_start_label = self.lattice_builder._create_node_label(node_start)
 268        node_end_label = self.lattice_builder._create_node_label(node_end)
 269        _, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end)
 270        if node_start == node_end and residue_list != residue_list[::-1]:
 271            raise ValueError(f"Aborted creation of hydrogel chain between '{node_start}' and '{node_end}' because pyMBE could not resolve a unique topology for that chain")
 272        if reverse:
 273            reverse_residue_order=True
 274        else:
 275            reverse_residue_order=False
 276        start_node_id = nodes[node_start_label]["id"]
 277        end_node_id = nodes[node_end_label]["id"]
 278        # Finding a backbone vector between node_start and node_end
 279        vec_between_nodes = np.array(nodes[node_end_label]["pos"]) - np.array(nodes[node_start_label]["pos"])
 280        vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l)
 281        backbone_vector = vec_between_nodes / np.linalg.norm(vec_between_nodes)
 282        if reverse_residue_order:
 283            vec_between_nodes *= -1.0
 284        # Calculate the start position of the chain
 285        chain_residues = self.db.get_template(pmb_type="molecule",
 286                                              name=molecule_name).residue_list
 287        part_start_chain_name = self.db.get_template(pmb_type="residue",
 288                                                     name=chain_residues[0]).central_bead
 289        lj_parameters = self.get_lj_parameters(particle_name1=nodes[node_start_label]["name"],
 290                                               particle_name2=part_start_chain_name)
 291        bond_tpl = self.get_bond_template(particle_name1=nodes[node_start_label]["name"],
 292                                          particle_name2=part_start_chain_name,
 293                                          use_default_bond=use_default_bond)
 294        l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
 295                                              bond_type=bond_tpl.bond_type,
 296                                              bond_parameters=bond_tpl.get_parameters(ureg=self.units))
 297        first_bead_pos = np.array((nodes[node_start_label]["pos"])) + np.array(backbone_vector)*l0
 298        mol_id = self.create_molecule(name=molecule_name,  # Use the name defined earlier
 299                                      number_of_molecules=1,  # Creating one chain
 300                                      espresso_system=espresso_system,
 301                                      list_of_first_residue_positions=[first_bead_pos.tolist()], #Start at the first node
 302                                      backbone_vector=np.array(backbone_vector)/l0,
 303                                      use_default_bond=use_default_bond,
 304                                      reverse_residue_order=reverse_residue_order,
 305                                      gen_angle=gen_angle)[0]
 306        # Bond chain to the hydrogel nodes
 307        chain_pids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
 308                                                             attribute="molecule_id",
 309                                                             value=mol_id)
 310        self.create_bond(particle_id1=start_node_id,
 311                         particle_id2=chain_pids[0],
 312                         espresso_system=espresso_system,
 313                         use_default_bond=use_default_bond)
 314        self.create_bond(particle_id1=chain_pids[-1],
 315                         particle_id2=end_node_id,
 316                         espresso_system=espresso_system,
 317                         use_default_bond=use_default_bond)
 318        return mol_id
 319
 320    def _generate_hydrogel_crosslinker_angles(self, espresso_system, central_particle_ids):
 321        """
 322        Generate hydrogel angles centered on crosslinkers and adjacent terminal beads.
 323
 324        If the user defines any explicit angle template for such junction
 325        triplets, then all required junction triplets must be defined. If none
 326        are defined, hydrogel construction proceeds without crosslinker-adjacent
 327        angles.
 328        """
 329        particle_instances = self.db.get_instances(pmb_type="particle")
 330        bonded_neighbors = {}
 331        for bond in self.db.get_instances(pmb_type="bond").values():
 332            bonded_neighbors.setdefault(bond.particle_id1, set()).add(bond.particle_id2)
 333            bonded_neighbors.setdefault(bond.particle_id2, set()).add(bond.particle_id1)
 334
 335        triplets = []
 336        for central_particle_id in sorted(set(central_particle_ids)):
 337            neighbors = sorted(bonded_neighbors.get(central_particle_id, set()))
 338            central_name = particle_instances[central_particle_id].name
 339            for idx_i in range(len(neighbors)):
 340                for idx_k in range(idx_i + 1, len(neighbors)):
 341                    side_particle_id1 = neighbors[idx_i]
 342                    side_particle_id3 = neighbors[idx_k]
 343                    side_name1 = particle_instances[side_particle_id1].name
 344                    side_name3 = particle_instances[side_particle_id3].name
 345                    angle_key = AngleTemplate.make_angle_key(side1=side_name1,
 346                                                             central=central_name,
 347                                                             side2=side_name3)
 348                    triplets.append((side_particle_id1,
 349                                     central_particle_id,
 350                                     side_particle_id3,
 351                                     angle_key))
 352
 353        defined_angle_templates = self.db.get_templates(pmb_type="angle")
 354        defined_angle_keys = {
 355            angle_key
 356            for _, _, _, angle_key in triplets
 357            if angle_key in defined_angle_templates
 358        }
 359        if not defined_angle_keys:
 360            logging.warning("No angle templates defined for hydrogel crosslinkers")
 361            return
 362
 363        missing_angle_keys = sorted({
 364            angle_key
 365            for _, _, _, angle_key in triplets
 366            if angle_key not in defined_angle_keys
 367        })
 368        if missing_angle_keys:
 369            raise ValueError(
 370                "Hydrogel crosslinker-adjacent angle templates must be defined for all required triplets. "
 371                f"Missing definitions for: {missing_angle_keys}"
 372            )
 373
 374        for side_particle_id1, central_particle_id, side_particle_id3, _ in triplets:
 375            self.create_angular_potential(particle_id1=side_particle_id1,
 376                              particle_id2=central_particle_id,
 377                              particle_id3=side_particle_id3,
 378                              espresso_system=espresso_system,
 379                              use_default_angle=False)
 380
 381    def _create_hydrogel_node(self, node_index, node_name, espresso_system):
 382        """
 383        Set a node residue type.
 384        
 385        Args:
 386            node_index ('str'): 
 387                Lattice node index in the form of a string, e.g. "[0 0 0]".
 388
 389            node_name ('str'): 
 390                name of the node particle defined in pyMBE.
 391
 392            espresso_system (espressomd.system.System): 
 393                ESPResSo system object where the hydrogel node will be created.
 394
 395        Returns:
 396            ('tuple(list,int)'):
 397                ('list'): Position of the node in the lattice.
 398                ('int'): Particle ID of the node.
 399        """
 400        if self.lattice_builder is None:
 401            raise ValueError("LatticeBuilder is not initialized. Use 'initialize_lattice_builder' first.")
 402        node_position = np.array(node_index)*0.25*self.lattice_builder.box_l
 403        p_id = self.create_particle(name = node_name,
 404                                    espresso_system=espresso_system,
 405                                    number_of_particles=1,
 406                                    position = [node_position])
 407        key = self.lattice_builder._get_node_by_label(f"[{node_index[0]} {node_index[1]} {node_index[2]}]")
 408        self.lattice_builder.nodes[key] = node_name
 409        return node_position.tolist(), p_id[0]
 410
 411    def _get_espresso_bond_instance(self, bond_template, espresso_system):
 412        """
 413        Retrieve or create a bond instance in an ESPResSo system for a given pair of particle names.
 414
 415        Args:
 416            bond_template ('BondTemplate'): 
 417                BondTemplate object from the pyMBE database.
 418            espresso_system ('espressomd.system.System'): 
 419                An ESPResSo system object where the bond will be added or retrieved.
 420
 421        Returns:
 422            ('espressomd.interactions.BondedInteraction'): 
 423                The ESPResSo bond instance object.
 424
 425        Notes:
 426            When a new bond instance is created, it is not added to the ESPResSo system.
 427        """
 428        if bond_template.name in self.db.espresso_bond_instances.keys():
 429            bond_inst = self.db.espresso_bond_instances[bond_template.name]
 430        else:   
 431            # Create an instance of the bond 
 432            bond_inst = self._create_espresso_bond_instance(bond_type=bond_template.bond_type,
 433                                                            bond_parameters=bond_template.get_parameters(self.units))
 434            self.db.espresso_bond_instances[bond_template.name]= bond_inst
 435            espresso_system.bonded_inter.add(bond_inst)
 436        return bond_inst
 437
 438    def _get_label_id_map(self, pmb_type):
 439        """
 440        Returns the key used to access the particle ID map for a given pyMBE object type.
 441
 442        Args:
 443            pmb_type ('str'):
 444                pyMBE object type for which the particle ID map label is requested.
 445
 446        Returns:
 447            'str':
 448                Label identifying the appropriate particle ID map. 
 449        """
 450        if pmb_type in self.db._assembly_like_types:
 451            label="assembly_map"
 452        elif pmb_type in self.db._molecule_like_types:
 453            label="molecule_map"
 454        else:
 455            label=f"{pmb_type}_map"
 456        return label
 457
 458    def _get_residue_list_from_sequence(self, sequence):
 459        """
 460        Convenience function to get a 'residue_list' from a protein or peptide 'sequence'.
 461
 462        Args:
 463            sequence ('lst'): 
 464                 Sequence of the peptide or protein.
 465
 466        Returns:
 467            residue_list ('list' of 'str'): 
 468                List of the 'name's of the 'residue's  in the sequence of the 'molecule'.             
 469        """
 470        residue_list = []
 471        for item in sequence:
 472            residue_name='AA-'+item
 473            residue_list.append(residue_name)
 474        return residue_list
 475    
 476    def _get_template_type(self, name, allowed_types):
 477        """
 478        Validate that a template name resolves unambiguously to exactly one
 479        allowed pmb_type in the pyMBE database and return it.
 480
 481        Args:
 482            name ('str'): 
 483                Name of the template to validate.
 484
 485            allowed_types ('set[str]'):  
 486                Set of allowed pmb_type values (e.g. {"molecule", "peptide"}).
 487
 488        Returns:
 489            ('str'): 
 490                Resolved pmb_type.
 491
 492        Notess:
 493            - This method does *not* return the template itself, only the validated pmb_type. 
 494        """
 495        registered_pmb_types_with_name = self.db._find_template_types(name=name)
 496        filtered_types = allowed_types.intersection(registered_pmb_types_with_name)
 497        if len(filtered_types) > 1:
 498            raise ValueError(f"Ambiguous template name '{name}': found {len(filtered_types)} templates in the pyMBE database. Molecule creation aborted.")  
 499        if len(filtered_types) == 0:
 500            raise ValueError(f"No {allowed_types} template found with name '{name}'. Found templates of types: {filtered_types}.")
 501        return next(iter(filtered_types))
 502
 503    def _delete_particles_from_espresso(self, particle_ids, espresso_system):
 504        """
 505        Remove a list of particles from an ESPResSo simulation system.
 506
 507        Args:
 508            particle_ids  ('Iterable[int]'):
 509                A list (or other iterable) of ESPResSo particle IDs to remove.
 510
 511            espresso_system ('espressomd.system.System'):
 512                The ESPResSo simulation system from which the particles
 513                will be removed.
 514
 515        Notess:
 516            - This method removes particles only from the ESPResSo simulation,
 517            **not** from the pyMBE database. Database cleanup must be handled
 518            separately by the caller.
 519            - Attempting to remove a non-existent particle ID will raise
 520            an ESPResSo error.
 521        """
 522        for pid in particle_ids:
 523            espresso_system.part.by_id(pid).remove()
 524
 525    def calculate_center_of_mass(self, instance_id, pmb_type, espresso_system):
 526        """
 527        Calculates the center of mass of a pyMBE object instance in an ESPResSo system.
 528
 529        Args:
 530            instance_id ('int'):
 531                pyMBE instance ID of the object whose center of mass is calculated.
 532
 533            pmb_type ('str'):
 534                Type of the pyMBE object. Must correspond to a particle-aggregating
 535                template type (e.g. '"molecule"', '"residue"', '"peptide"', '"protein"').
 536
 537            espresso_system ('espressomd.system.System'):
 538                ESPResSo system containing the particle instances.
 539
 540        Returns:
 541            ('numpy.ndarray'):
 542                Array of shape '(3,)' containing the Cartesian coordinates of the
 543                center of mass.
 544
 545        Notes:
 546            - This method assumes equal mass for all particles.
 547            - Periodic boundary conditions are *not* unfolded; positions are taken
 548            directly from ESPResSo particle coordinates.
 549        """
 550        center_of_mass = np.zeros(3)
 551        axis_list = [0,1,2]
 552        inst = self.db.get_instance(pmb_type=pmb_type,
 553                                    instance_id=instance_id)
 554        particle_id_list = self.get_particle_id_map(object_name=inst.name)["all"]
 555        for pid in particle_id_list:
 556            for axis in axis_list:
 557                center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis]
 558        center_of_mass = center_of_mass / len(particle_id_list)
 559        return center_of_mass
 560
 561    def calculate_HH(self, template_name, pH_list=None, pka_set=None):
 562        """
 563        Calculates the charge in the template object according to the ideal  Henderson–Hasselbalch titration curve.
 564
 565        Args:
 566            template_name ('str'):
 567                Name of the template.
 568
 569            pH_list ('list[float]', optional):
 570                pH values at which the charge is evaluated.
 571                Defaults to 50 values between 2 and 12.
 572
 573            pka_set ('dict', optional):
 574                Mapping: {particle_name: {"pka_value": 'float', "acidity": "acidic"|"basic"}}
 575
 576        Returns:
 577            'list[float]':
 578                Net molecular charge at each pH value.
 579        """
 580        if pH_list is None:
 581            pH_list = np.linspace(2, 12, 50)
 582        if pka_set is None:
 583            pka_set = self.get_pka_set()
 584        self._check_pka_set(pka_set=pka_set)
 585        particle_counts = self.db.get_particle_templates_under(template_name=template_name,
 586                                                               return_counts=True)
 587        if not particle_counts:
 588            return [None] * len(pH_list)
 589        charge_number_map = self.get_charge_number_map()
 590        def formal_charge(particle_name):
 591            tpl = self.db.get_template(name=particle_name, 
 592                                       pmb_type="particle")
 593            state = self.db.get_template(name=tpl.initial_state,
 594                                         pmb_type="particle_state")
 595            return charge_number_map[state.es_type]
 596        Z_HH = []
 597        for pH in pH_list:
 598            Z = 0.0
 599            for particle, multiplicity in particle_counts.items():
 600                if particle in pka_set:
 601                    pka = pka_set[particle]["pka_value"]
 602                    acidity = pka_set[particle]["acidity"]
 603                    if acidity == "acidic":
 604                        psi = -1
 605                    elif acidity == "basic":
 606                        psi = +1
 607                    else:
 608                        raise ValueError(f"Unknown acidity '{acidity}' for particle '{particle}'")
 609                    charge = psi / (1.0 + 10.0 ** (psi * (pH - pka)))
 610                    Z += multiplicity * charge
 611                else:
 612                    Z += multiplicity * formal_charge(particle)
 613            Z_HH.append(Z)
 614        return Z_HH   
 615
 616    def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
 617        """
 618        Computes macromolecular charges using the Henderson–Hasselbalch equation
 619        coupled to ideal Donnan partitioning.
 620
 621        Args:
 622            c_macro ('dict'):
 623                Mapping of macromolecular species names to their concentrations
 624                in the system:
 625                '{molecule_name: concentration}'.
 626                
 627            c_salt ('float' or 'pint.Quantity'):
 628                Salt concentration in the reservoir.
 629
 630            pH_list ('list[float]', optional):
 631                List of pH values in the reservoir at which the calculation is
 632                performed. If 'None', 50 equally spaced values between 2 and 12
 633                are used.
 634
 635            pka_set ('dict', optional):
 636                Dictionary defining the acid–base properties of titratable particle
 637                types:
 638                '{particle_name: {"pka_value": float, "acidity": "acidic" | "basic"}}'.
 639                If 'None', the pKa set is taken from the pyMBE database.
 640
 641        Returns:
 642            'dict':
 643                Dictionary containing:
 644                - '"charges_dict"' ('dict'):
 645                    Mapping '{molecule_name: list}' of Henderson–Hasselbalch–Donnan
 646                    charges evaluated at each pH value.
 647                - '"pH_system_list"' ('list[float]'):
 648                    Effective pH values inside the system phase after Donnan
 649                    partitioning.
 650                - '"partition_coefficients"' ('list[float]'):
 651                    Partition coefficients of monovalent cations at each pH value.
 652
 653        Notes:
 654            - This method assumes **ideal Donnan equilibrium** and **monovalent salt**.
 655            - The ionic strength of the reservoir includes both salt and
 656            pH-dependent H⁺/OH⁻ contributions.
 657            - All charged macromolecular species present in the system must be
 658            included in 'c_macro'; missing species will lead to incorrect results.
 659            - The nonlinear Donnan equilibrium equation is solved using a scalar
 660            root finder ('brentq') in logarithmic form for numerical stability.
 661            - This method is intended for **two-phase systems**; for single-phase
 662            systems use 'calculate_HH' instead.
 663        """
 664        if pH_list is None:
 665            pH_list=np.linspace(2,12,50)
 666        if pka_set is None:
 667            pka_set=self.get_pka_set() 
 668        self._check_pka_set(pka_set=pka_set)
 669        partition_coefficients_list = []
 670        pH_system_list = []
 671        Z_HH_Donnan={}
 672        for key in c_macro:
 673            Z_HH_Donnan[key] = []
 674        def calc_charges(c_macro, pH):
 675            """
 676            Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation.
 677
 678            Args:
 679                c_macro ('dict'): 
 680                    {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 681
 682                pH ('float'): 
 683                    pH-value that is used in the HH equation.
 684
 685            Returns:
 686                ('dict'): 
 687                    {"molecule_name": charge}
 688            """
 689            charge = {}
 690            for name in c_macro:
 691                charge[name] = self.calculate_HH(name, [pH], pka_set)[0]
 692            return charge
 693
 694        def calc_partition_coefficient(charge, c_macro):
 695            """
 696            Calculates the partition coefficients of positive ions according to the ideal Donnan theory.
 697
 698            Args:
 699                charge ('dict'): 
 700                    {"molecule_name": charge}
 701
 702                c_macro ('dict'): 
 703                    {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 704            """
 705            nonlocal ionic_strength_res
 706            charge_density = 0.0
 707            for key in charge:
 708                charge_density += charge[key] * c_macro[key]
 709            return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude
 710        for pH_value in pH_list:    
 711            # calculate the ionic strength of the reservoir
 712            if pH_value <= 7.0:
 713                ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 
 714            elif pH_value > 7.0:
 715                ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt
 716            #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations
 717            #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation.
 718            #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects.
 719            equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro))
 720            logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq")
 721            partition_coefficient = 10**logxi.root
 722            charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient))
 723            for key in c_macro:
 724                Z_HH_Donnan[key].append(charges_temp[key])
 725            pH_system_list.append(pH_value - np.log10(partition_coefficient))
 726            partition_coefficients_list.append(partition_coefficient)
 727        return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
 728
 729    def calculate_net_charge(self,espresso_system,object_name,pmb_type,dimensionless=False):
 730        """
 731        Calculates the net charge per instance of a given pmb object type.
 732
 733        Args:
 734            espresso_system (espressomd.system.System):
 735                ESPResSo system containing the particles.
 736            object_name (str):
 737                Name of the object (e.g. molecule, residue, peptide, protein).
 738            pmb_type (str):
 739                Type of object to analyze. Must be molecule-like.
 740            dimensionless (bool, optional):
 741                If True, return charge as a pure number.
 742                If False, return a quantity with reduced_charge units.
 743
 744        Returns:
 745            dict:
 746                {"mean": mean_net_charge, "instances": {instance_id: net_charge}}
 747        """
 748        id_map = self.get_particle_id_map(object_name=object_name)
 749        label = self._get_label_id_map(pmb_type=pmb_type)
 750        instance_map = id_map[label]
 751        charges = {}
 752        for instance_id, particle_ids in instance_map.items():
 753            if dimensionless:
 754                net_charge = 0.0
 755            else:
 756                net_charge = 0 * self.units.Quantity(1, "reduced_charge")
 757            for pid in particle_ids:
 758                q = espresso_system.part.by_id(pid).q
 759                if not dimensionless:
 760                    q *= self.units.Quantity(1, "reduced_charge")
 761                net_charge += q
 762            charges[instance_id] = net_charge
 763        # Mean charge
 764        if dimensionless:
 765            mean_charge = float(np.mean(list(charges.values())))
 766        else:
 767            mean_charge = (np.mean([q.magnitude for q in charges.values()])* self.units.Quantity(1, "reduced_charge"))
 768        return {"mean": mean_charge, "instances": charges}
 769
 770    def center_object_in_simulation_box(self, instance_id, espresso_system, pmb_type):
 771        """
 772        Centers a pyMBE object instance in the simulation box of an ESPResSo system.
 773        The object is translated such that its center of mass coincides with the
 774        geometric center of the ESPResSo simulation box.
 775
 776        Args:
 777            instance_id ('int'):
 778                ID of the pyMBE object instance to be centered.
 779
 780            pmb_type ('str'):
 781                Type of the pyMBE object.
 782
 783            espresso_system ('espressomd.system.System'):
 784                ESPResSo system object in which the particles are defined.
 785
 786        Notes:
 787            - Works for both cubic and non-cubic simulation boxes.
 788        """
 789        inst = self.db.get_instance(instance_id=instance_id,
 790                                    pmb_type=pmb_type)
 791        center_of_mass = self.calculate_center_of_mass(instance_id=instance_id,
 792                                                       espresso_system=espresso_system,
 793                                                       pmb_type=pmb_type)
 794        box_center = [espresso_system.box_l[0]/2.0,
 795                      espresso_system.box_l[1]/2.0,
 796                      espresso_system.box_l[2]/2.0]
 797        particle_id_list = self.get_particle_id_map(object_name=inst.name)["all"]
 798        for pid in particle_id_list:
 799            es_pos = espresso_system.part.by_id(pid).pos
 800            espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center
 801
 802    def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt):    
 803        """
 804        Creates a 'c_salt' concentration of 'cation_name' and 'anion_name' ions into the 'espresso_system'.
 805
 806        Args:
 807            espresso_system('espressomd.system.System'): instance of an espresso system object.
 808            cation_name('str'): 'name' of a particle with a positive charge.
 809            anion_name('str'): 'name' of a particle with a negative charge.
 810            c_salt('float'): Salt concentration.
 811            
 812        Returns:
 813            c_salt_calculated('float'): Calculated salt concentration added to 'espresso_system'.
 814        """ 
 815        cation_tpl = self.db.get_template(pmb_type="particle",
 816                                          name=cation_name)
 817        cation_state = self.db.get_template(pmb_type="particle_state",
 818                                            name=cation_tpl.initial_state)
 819        cation_charge = cation_state.z
 820        anion_tpl = self.db.get_template(pmb_type="particle",
 821                                          name=anion_name)
 822        anion_state = self.db.get_template(pmb_type="particle_state",
 823                                            name=anion_tpl.initial_state)
 824        anion_charge = anion_state.z
 825        if cation_charge <= 0:
 826            raise ValueError(f'ERROR cation charge must be positive, charge {cation_charge}')
 827        if anion_charge >= 0:
 828            raise ValueError(f'ERROR anion charge must be negative, charge {anion_charge}')
 829        # Calculate the number of ions in the simulation box
 830        volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3')
 831        if c_salt.check('[substance] [length]**-3'):
 832            N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude)
 833            c_salt_calculated=N_ions/(volume*self.N_A)
 834        elif c_salt.check('[length]**-3'):
 835            N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude)
 836            c_salt_calculated=N_ions/volume
 837        else:
 838            raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt)
 839        N_cation = N_ions*abs(anion_charge)
 840        N_anion = N_ions*abs(cation_charge)
 841        self.create_particle(espresso_system=espresso_system, 
 842                             name=cation_name, 
 843                             number_of_particles=N_cation)
 844        self.create_particle(espresso_system=espresso_system, 
 845                             name=anion_name, 
 846                             number_of_particles=N_anion)
 847        if c_salt_calculated.check('[substance] [length]**-3'):
 848            logging.info(f"added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions")
 849        elif c_salt_calculated.check('[length]**-3'):
 850            logging.info(f"added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions")
 851        return c_salt_calculated
 852
 853    def create_bond(self, particle_id1, particle_id2, espresso_system, use_default_bond=False):
 854        """
 855        Creates a bond between two particle instances in an ESPResSo system and registers it in the pyMBE database.
 856
 857        This method performs the following steps:
 858            1. Retrieves the particle instances corresponding to 'particle_id1' and 'particle_id2' from the database.
 859            2. Retrieves or creates the corresponding ESPResSo bond instance using the bond template.
 860            3. Adds the ESPResSo bond instance to the ESPResSo system if it was newly created.
 861            4. Adds the bond to the first particle's bond list in ESPResSo.
 862            5. Creates a 'BondInstance' in the database and registers it.
 863
 864        Args:
 865            particle_id1 ('int'): 
 866                pyMBE and ESPResSo ID of the first particle.
 867
 868            particle_id2 ('int'): 
 869                pyMBE and ESPResSo ID of the second particle.
 870
 871            espresso_system ('espressomd.system.System'): 
 872                ESPResSo system object where the bond will be created.
 873
 874            use_default_bond ('bool', optional): 
 875                If True, use a default bond template if no specific template exists. Defaults to False.
 876
 877        Returns:
 878            ('int'): 
 879                bond_id of the bond instance created in the pyMBE database.
 880        """
 881        particle_inst_1 = self.db.get_instance(pmb_type="particle",
 882                                               instance_id=particle_id1)
 883        particle_inst_2 = self.db.get_instance(pmb_type="particle",
 884                                               instance_id=particle_id2)
 885        bond_tpl = self.get_bond_template(particle_name1=particle_inst_1.name,
 886                                          particle_name2=particle_inst_2.name,
 887                                          use_default_bond=use_default_bond)
 888        bond_inst = self._get_espresso_bond_instance(bond_template=bond_tpl,
 889                                                    espresso_system=espresso_system)
 890        espresso_system.part.by_id(particle_id1).add_bond((bond_inst, particle_id2))
 891        bond_id = self.db._propose_instance_id(pmb_type="bond")
 892        pmb_bond_instance = BondInstance(bond_id=bond_id,
 893                                         name=bond_tpl.name,
 894                                         particle_id1=particle_id1,
 895                                         particle_id2=particle_id2)
 896        self.db._register_instance(instance=pmb_bond_instance)
 897
 898    def create_counterions(self, object_name, cation_name, anion_name, espresso_system):
 899        """
 900        Creates particles of 'cation_name' and 'anion_name' in 'espresso_system' to counter the net charge of 'object_name'.
 901        
 902        Args:
 903            object_name ('str'): 
 904                'name' of a pyMBE object.
 905
 906            espresso_system ('espressomd.system.System'): 
 907                Instance of a system object from the espressomd library.
 908
 909            cation_name ('str'): 
 910                'name' of a particle with a positive charge.
 911
 912            anion_name ('str'): 
 913                'name' of a particle with a negative charge.
 914
 915        Returns: 
 916            ('dict'): 
 917                {"name": number}
 918
 919        Notes:
 920            This function currently does not support the creation of counterions for hydrogels.
 921        """ 
 922        cation_tpl = self.db.get_template(pmb_type="particle",
 923                                          name=cation_name)
 924        cation_state = self.db.get_template(pmb_type="particle_state",
 925                                            name=cation_tpl.initial_state)
 926        cation_charge = cation_state.z
 927        anion_tpl = self.db.get_template(pmb_type="particle",
 928                                          name=anion_name)
 929        anion_state = self.db.get_template(pmb_type="particle_state",
 930                                            name=anion_tpl.initial_state)
 931        anion_charge = anion_state.z
 932        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
 933        counterion_number={}
 934        object_charge={}
 935        for name in ['positive', 'negative']:
 936            object_charge[name]=0
 937        for id in object_ids:
 938            if espresso_system.part.by_id(id).q > 0:
 939                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 940            elif espresso_system.part.by_id(id).q < 0:
 941                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 942        if object_charge['positive'] % abs(anion_charge) == 0:
 943            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
 944        else:
 945            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
 946        if object_charge['negative'] % abs(cation_charge) == 0:
 947            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
 948        else:
 949            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
 950        if counterion_number[cation_name] > 0: 
 951            self.create_particle(espresso_system=espresso_system, 
 952                                 name=cation_name, 
 953                                 number_of_particles=counterion_number[cation_name])
 954        else:
 955            counterion_number[cation_name]=0
 956        if counterion_number[anion_name] > 0:
 957            self.create_particle(espresso_system=espresso_system, 
 958                                 name=anion_name, 
 959                                 number_of_particles=counterion_number[anion_name])
 960        else:
 961            counterion_number[anion_name] = 0
 962        logging.info('the following counter-ions have been created: ')
 963        for name in counterion_number.keys():
 964            logging.info(f'Ion type: {name} created number: {counterion_number[name]}')
 965        return counterion_number
 966
 967    def create_hydrogel(self, name, espresso_system, use_default_bond=False, gen_angle=False):
 968        """ 
 969        Creates a hydrogel in espresso_system using a pyMBE hydrogel template given by 'name'
 970
 971        Args:
 972            name ('str'): 
 973                name of the hydrogel template in the pyMBE database.
 974
 975            espresso_system ('espressomd.system.System'): 
 976                ESPResSo system object where the hydrogel will be created.
 977
 978            use_default_bond ('bool', optional): 
 979                If True, use a default bond template if no specific template exists. Defaults to False.
 980
 981            gen_angle ('bool', optional):
 982                If True, generate angle potentials for the internal hydrogel
 983                chains and, when explicitly defined, for all crosslinker-adjacent
 984                triplets. Defaults to False.
 985
 986        Returns:
 987            ('int'): id of the hydrogel instance created.
 988        """
 989        if not self.db._has_template(name=name, pmb_type="hydrogel"):
 990            raise ValueError(f"Hydrogel template with name '{name}' is not defined in the pyMBE database.")
 991        hydrogel_tpl = self.db.get_template(pmb_type="hydrogel",
 992                                            name=name)
 993        assembly_id = self.db._propose_instance_id(pmb_type="hydrogel")
 994        # Create the nodes
 995        nodes = {}
 996        hydrogel_angle_centers = set()
 997        node_topology = hydrogel_tpl.node_map
 998        for node in node_topology:
 999            node_index = node.lattice_index
1000            node_name = node.particle_name
1001            node_pos, node_id = self._create_hydrogel_node(node_index=node_index,
1002                                                          node_name=node_name,
1003                                                          espresso_system=espresso_system)
1004            node_label = self.lattice_builder._create_node_label(node_index=node_index)
1005            nodes[node_label] = {"name": node_name, "id": node_id, "pos": node_pos} 
1006            self.db._update_instance(instance_id=node_id,
1007                                     pmb_type="particle",
1008                                     attribute="assembly_id",
1009                                     value=assembly_id)
1010        for hydrogel_chain in hydrogel_tpl.chain_map:
1011            molecule_id = self._create_hydrogel_chain(hydrogel_chain=hydrogel_chain,
1012                                                      nodes=nodes, 
1013                                                      espresso_system=espresso_system,
1014                                                      use_default_bond=use_default_bond,
1015                                                      gen_angle=gen_angle)
1016            self.db._update_instance(instance_id=molecule_id,
1017                                     pmb_type="molecule",
1018                                     attribute="assembly_id",
1019                                     value=assembly_id)
1020            if gen_angle:
1021                residue_ids = self.db._find_instance_ids_by_attribute(pmb_type="residue",
1022                                                                      attribute="molecule_id",
1023                                                                      value=molecule_id)
1024                first_residue_id = min(residue_ids)
1025                last_residue_id = max(residue_ids)
1026                first_residue = self.db.get_instance(pmb_type="residue",
1027                                                     instance_id=first_residue_id)
1028                last_residue = self.db.get_instance(pmb_type="residue",
1029                                                    instance_id=last_residue_id)
1030                first_central_bead_name = self.db.get_template(pmb_type="residue",
1031                                                               name=first_residue.name).central_bead
1032                last_central_bead_name = self.db.get_template(pmb_type="residue",
1033                                                              name=last_residue.name).central_bead
1034                particle_instances = self.db.get_instances(pmb_type="particle")
1035                first_residue_particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1036                                                                                      attribute="residue_id",
1037                                                                                      value=first_residue_id)
1038                last_residue_particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1039                                                                                     attribute="residue_id",
1040                                                                                     value=last_residue_id)
1041                first_bead_id = None
1042                for particle_id in first_residue_particle_ids:
1043                    if particle_instances[particle_id].name == first_central_bead_name:
1044                        first_bead_id = particle_id
1045                        break
1046
1047                last_bead_id = None
1048                for particle_id in last_residue_particle_ids:
1049                    if particle_instances[particle_id].name == last_central_bead_name:
1050                        last_bead_id = particle_id
1051                        break
1052                node_start_label = self.lattice_builder._create_node_label(hydrogel_chain.node_start)
1053                node_end_label = self.lattice_builder._create_node_label(hydrogel_chain.node_end)
1054                hydrogel_angle_centers.update({
1055                    nodes[node_start_label]["id"],
1056                    nodes[node_end_label]["id"],
1057                    first_bead_id,
1058                    last_bead_id,
1059                })
1060        self.db._propagate_id(root_type="hydrogel", 
1061                                root_id=assembly_id, 
1062                                attribute="assembly_id", 
1063                                value=assembly_id)
1064        if gen_angle:
1065            self._generate_hydrogel_crosslinker_angles(espresso_system=espresso_system,
1066                                                       central_particle_ids=hydrogel_angle_centers)
1067        # Register an hydrogel instance in the pyMBE databasegit 
1068        self.db._register_instance(HydrogelInstance(name=name,
1069                                                    assembly_id=assembly_id))
1070        return assembly_id
1071
1072    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False, reverse_residue_order = False, gen_angle=False):
1073        """
1074        Creates instances of a given molecule template name into ESPResSo.
1075
1076        Args:
1077            name ('str'): 
1078                Label of the molecule type to be created. 'name'.
1079
1080            espresso_system ('espressomd.system.System'): 
1081                Instance of a system object from espressomd library.
1082
1083            number_of_molecules ('int'): 
1084                Number of molecules or peptides of type 'name' to be created.
1085
1086            list_of_first_residue_positions ('list', optional): 
1087                List of coordinates where the central bead of the first_residue_position will be created, random by default.
1088
1089            backbone_vector ('list' of 'float'): 
1090                Backbone vector of the molecule, random by default. Central beads of the residues in the 'residue_list' are placed along this vector. 
1091
1092            use_default_bond('bool', optional): 
1093                Controls if a bond of type 'default' is used to bond particles with undefined bonds in the pyMBE database.
1094
1095            reverse_residue_order('bool', optional): 
1096                Creates residues in reverse sequential order than the one defined in the molecule template. Defaults to False.
1097
1098        Returns:
1099            ('list' of 'int'): 
1100                List with the 'molecule_id' of the pyMBE molecule instances created into 'espresso_system'.
1101
1102        Notes:
1103            - This function can be used to create both molecules and peptides.    
1104        """
1105        pmb_type = self._get_template_type(name=name,
1106                                           allowed_types={"molecule", "peptide"})
1107        if number_of_molecules <= 0:
1108            return {}
1109        if list_of_first_residue_positions is not None:
1110            for item in list_of_first_residue_positions:
1111                if not isinstance(item, list):
1112                    raise ValueError("The provided input position is not a nested list. Should be a nested list with elements of 3D lists, corresponding to xyz coord.")
1113                elif len(item) != 3:
1114                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
1115
1116            if len(list_of_first_residue_positions) != number_of_molecules:
1117                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
1118        # Generate an arbitrary random unit vector
1119        if backbone_vector is None:
1120            backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
1121                                                                      radius=1, 
1122                                                                      n_samples=1,
1123                                                                      on_surface=True)[0]
1124        else:
1125            backbone_vector = np.array(backbone_vector)
1126        first_residue = True
1127        molecule_tpl = self.db.get_template(pmb_type=pmb_type,
1128                                            name=name)
1129        if reverse_residue_order:
1130            residue_list = molecule_tpl.residue_list[::-1]
1131        else:
1132            residue_list = molecule_tpl.residue_list
1133        pos_index = 0 
1134        molecule_ids = []
1135        for _ in range(number_of_molecules):        
1136            molecule_id = self.db._propose_instance_id(pmb_type=pmb_type)
1137            for residue in residue_list:
1138                if first_residue:
1139                    if list_of_first_residue_positions is None:
1140                        central_bead_pos = None
1141                    else:
1142                        for item in list_of_first_residue_positions:
1143                            central_bead_pos = [np.array(list_of_first_residue_positions[pos_index])]
1144                            
1145                    residue_id = self.create_residue(name=residue,
1146                                                     espresso_system=espresso_system, 
1147                                                     central_bead_position=central_bead_pos,  
1148                                                     use_default_bond= use_default_bond, 
1149                                                     backbone_vector=backbone_vector)
1150                    
1151                    # Add molecule_id to the residue instance and all particles associated
1152                    self.db._propagate_id(root_type="residue", 
1153                                          root_id=residue_id,
1154                                          attribute="molecule_id", 
1155                                          value=molecule_id)
1156                    particle_ids_in_residue = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1157                                                                                      attribute="residue_id",
1158                                                                                      value=residue_id)
1159                    prev_central_bead_id = particle_ids_in_residue[0]
1160                    prev_central_bead_name = self.db.get_instance(pmb_type="particle", 
1161                                                                  instance_id=prev_central_bead_id).name
1162                    prev_central_bead_pos = espresso_system.part.by_id(prev_central_bead_id).pos
1163                    first_residue = False          
1164                else:
1165                    
1166                    # Calculate the starting position of the new residue
1167                    residue_tpl = self.db.get_template(pmb_type="residue",
1168                                                       name=residue)
1169                    lj_parameters = self.get_lj_parameters(particle_name1=prev_central_bead_name,
1170                                                           particle_name2=residue_tpl.central_bead)
1171                    bond_tpl = self.get_bond_template(particle_name1=prev_central_bead_name,
1172                                                      particle_name2=residue_tpl.central_bead,
1173                                                      use_default_bond=use_default_bond)
1174                    l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1175                                                          bond_type=bond_tpl.bond_type,
1176                                                          bond_parameters=bond_tpl.get_parameters(ureg=self.units))
1177                    central_bead_pos = prev_central_bead_pos+backbone_vector*l0
1178                    # Create the residue
1179                    residue_id = self.create_residue(name=residue, 
1180                                                     espresso_system=espresso_system, 
1181                                                     central_bead_position=[central_bead_pos],
1182                                                     use_default_bond= use_default_bond, 
1183                                                     backbone_vector=backbone_vector)
1184                    # Add molecule_id to the residue instance and all particles associated
1185                    self.db._propagate_id(root_type="residue", 
1186                                          root_id=residue_id, 
1187                                          attribute="molecule_id", 
1188                                          value=molecule_id)
1189                    particle_ids_in_residue = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1190                                                                                      attribute="residue_id",
1191                                                                                      value=residue_id)
1192                    central_bead_id = particle_ids_in_residue[0]
1193
1194                    # Bond the central beads of the new and previous residues
1195                    self.create_bond(particle_id1=prev_central_bead_id,
1196                                     particle_id2=central_bead_id,
1197                                     espresso_system=espresso_system,
1198                                     use_default_bond=use_default_bond)
1199                    
1200                    prev_central_bead_id = central_bead_id                    
1201                    prev_central_bead_name = self.db.get_instance(pmb_type="particle", instance_id=central_bead_id).name
1202                    prev_central_bead_pos =central_bead_pos
1203            # Create a Peptide or Molecule instance and register it on the pyMBE database
1204            if pmb_type == "molecule":
1205                inst = MoleculeInstance(molecule_id=molecule_id,
1206                                        name=name)
1207            elif pmb_type == "peptide":
1208                inst = PeptideInstance(name=name,
1209                                       molecule_id=molecule_id)
1210            self.db._register_instance(inst)
1211            if gen_angle:
1212                self._generate_angles_for_entity(
1213                    espresso_system=espresso_system,
1214                    entity_id=molecule_id,
1215                    entity_id_col='molecule_id')
1216            first_residue = True
1217            pos_index+=1
1218            molecule_ids.append(molecule_id)
1219        return molecule_ids
1220    
1221    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
1222        """
1223        Creates one or more particles in an ESPResSo system based on the particle template in the pyMBE database.
1224        
1225        Args:
1226            name ('str'): 
1227                Label of the particle template in the pyMBE database. 
1228
1229            espresso_system ('espressomd.system.System'): 
1230                Instance of a system object from the espressomd library.
1231
1232            number_of_particles ('int'): 
1233                Number of particles to be created.
1234
1235            position (list of ['float','float','float'], optional): 
1236                Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
1237
1238            fix ('bool', optional): 
1239                Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
1240
1241        Returns:
1242            ('list' of 'int'): 
1243                List with the ids of the particles created into 'espresso_system'.
1244        """       
1245        if number_of_particles <=0:
1246            return []
1247        if not self.db._has_template(name=name, pmb_type="particle"):
1248            raise ValueError(f"Particle template with name '{name}' is not defined in the pyMBE database.")
1249        
1250        part_tpl = self.db.get_template(pmb_type="particle",
1251                                        name=name)
1252        part_state = self.db.get_template(pmb_type="particle_state",
1253                                         name=part_tpl.initial_state)
1254        z = part_state.z
1255        es_type = part_state.es_type
1256        # Create the new particles into  ESPResSo 
1257        created_pid_list=[]
1258        for index in range(number_of_particles):
1259            if position is None:
1260                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1261            else:
1262                particle_position = position[index]
1263            
1264            particle_id = self.db._propose_instance_id(pmb_type="particle")
1265            created_pid_list.append(particle_id)
1266            kwargs = dict(id=particle_id, pos=particle_position, type=es_type, q=z)
1267            if fix:
1268                kwargs["fix"] = 3 * [fix]
1269            espresso_system.part.add(**kwargs)
1270            part_inst = ParticleInstance(name=name,
1271                                         particle_id=particle_id,
1272                                         initial_state=part_state.name)
1273            self.db._register_instance(part_inst)
1274                              
1275        return created_pid_list
1276
1277    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1278        """
1279        Creates one or more protein molecules in an ESPResSo system based on the 
1280        protein template in the pyMBE database and a provided topology.
1281
1282        Args:
1283            name (str):
1284                Name of the protein template stored in the pyMBE database.
1285            
1286            number_of_proteins (int):
1287                Number of protein molecules to generate.  
1288            
1289            espresso_system (espressomd.system.System):
1290                The ESPResSo simulation system where the protein molecules will be created.
1291            
1292            topology_dict (dict):
1293                Dictionary defining the internal structure of the protein. Expected format:
1294                    {"ResidueName1": {"initial_pos": np.ndarray,
1295                                      "chain_id": int,
1296                                      "radius": float},
1297                     "ResidueName2": { ... },
1298                        ...
1299                    }
1300                The '"initial_pos"' entry is required and represents the residue’s
1301                reference coordinates before shifting to the protein's center-of-mass.
1302
1303        Returns:
1304            ('list' of 'int'): 
1305                List of the molecule_id of the Protein instances created into ESPResSo.
1306
1307        Notes:
1308            - Particles are created using 'create_particle()' with 'fix=True',
1309            meaning they are initially immobilized.
1310            - The function assumes all residues in 'topology_dict' correspond to
1311            particle templates already defined in the pyMBE database.
1312            - Bonds between residues are not created here; it assumes a rigid body representation of the protein.
1313        """
1314        if number_of_proteins <= 0:
1315            return
1316        if not self.db._has_template(name=name, pmb_type="protein"):
1317            raise ValueError(f"Protein template with name '{name}' is not defined in the pyMBE database.")
1318        protein_tpl = self.db.get_template(pmb_type="protein", name=name)
1319        box_half = espresso_system.box_l[0] / 2.0
1320        # Create protein
1321        mol_ids = []
1322        for _ in range(number_of_proteins):
1323            # create a molecule identifier in pyMBE
1324            molecule_id = self.db._propose_instance_id(pmb_type="protein")
1325            # place protein COM randomly
1326            protein_center = self.generate_coordinates_outside_sphere(radius=1,
1327                                                                      max_dist=box_half,
1328                                                                      n_samples=1,
1329                                                                      center=[box_half]*3)[0]
1330            residues = hf.get_residues_from_topology_dict(topology_dict=topology_dict,
1331                                                         model=protein_tpl.model)
1332            # CREATE RESIDUES + PARTICLES
1333            for _, rdata in residues.items():
1334                base_resname = rdata["resname"]  
1335                residue_name = f"AA-{base_resname}"
1336                # residue instance ID
1337                residue_id = self.db._propose_instance_id("residue")
1338                # register ResidueInstance
1339                self.db._register_instance(ResidueInstance(name=residue_name,
1340                                                           residue_id=residue_id,
1341                                                           molecule_id=molecule_id))
1342                # PARTICLE CREATION
1343                for bead_id in rdata["beads"]:
1344                    bead_type = re.split(r'\d+', bead_id)[0]
1345                    relative_pos = topology_dict[bead_id]["initial_pos"]
1346                    absolute_pos = relative_pos + protein_center
1347                    particle_id = self.create_particle(name=bead_type,
1348                                                       espresso_system=espresso_system,
1349                                                       number_of_particles=1,
1350                                                       position=[absolute_pos],
1351                                                       fix=True)[0]
1352                    # update metadata
1353                    self.db._update_instance(instance_id=particle_id,
1354                                             pmb_type="particle",
1355                                             attribute="molecule_id",
1356                                             value=molecule_id)
1357                    self.db._update_instance(instance_id=particle_id,
1358                                             pmb_type="particle",
1359                                             attribute="residue_id",
1360                                             value=residue_id)
1361            protein_inst = ProteinInstance(name=name,
1362                                           molecule_id=molecule_id)
1363            self.db._register_instance(protein_inst)
1364            mol_ids.append(molecule_id)
1365        return mol_ids
1366
1367    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None, gen_angle=False):
1368        """
1369        Creates a residue  into ESPResSo.
1370
1371        Args:
1372            name ('str'): 
1373                Label of the residue type to be created. 
1374
1375            espresso_system ('espressomd.system.System'): 
1376                Instance of a system object from espressomd library.
1377
1378            central_bead_position ('list' of 'float'): 
1379                Position of the central bead.
1380
1381            use_default_bond ('bool'): 
1382                Switch to control if a bond of type 'default' is used to bond a particle whose bonds types are not defined in the pyMBE database.
1383
1384            backbone_vector ('list' of 'float'): 
1385                Backbone vector of the molecule. All side chains are created perpendicularly to 'backbone_vector'.
1386
1387        Returns:
1388            (int): 
1389                residue_id of the residue created.
1390        """
1391        if not self.db._has_template(name=name, pmb_type="residue"):
1392            raise ValueError(f"Residue template with name '{name}' is not defined in the pyMBE database.")
1393        res_tpl = self.db.get_template(pmb_type="residue",
1394                                       name=name)
1395        # Assign a residue_id
1396        residue_id = self.db._propose_instance_id(pmb_type="residue")
1397        res_inst = ResidueInstance(name=name,
1398                                   residue_id=residue_id)
1399        self.db._register_instance(res_inst)
1400        # create the principal bead   
1401        central_bead_name = res_tpl.central_bead 
1402        central_bead_id = self.create_particle(name=central_bead_name,
1403                                               espresso_system=espresso_system,
1404                                               position=central_bead_position,
1405                                               number_of_particles = 1)[0]
1406        
1407        central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1408        # Assigns residue_id to the central_bead particle created.
1409        self.db._update_instance(pmb_type="particle",
1410                                 instance_id=central_bead_id,
1411                                 attribute="residue_id",
1412                                 value=residue_id)
1413        
1414        # create the lateral beads  
1415        side_chain_list = res_tpl.side_chains
1416        side_chain_beads_ids = []
1417        for side_chain_name in side_chain_list:
1418            pmb_type = self._get_template_type(name=side_chain_name,
1419                                               allowed_types={"particle", "residue"})
1420            if pmb_type == 'particle':
1421                lj_parameters = self.get_lj_parameters(particle_name1=central_bead_name,
1422                                                       particle_name2=side_chain_name)
1423                bond_tpl = self.get_bond_template(particle_name1=central_bead_name,
1424                                                  particle_name2=side_chain_name,
1425                                                  use_default_bond=use_default_bond)
1426                l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1427                                                      bond_type=bond_tpl.bond_type,
1428                                                      bond_parameters=bond_tpl.get_parameters(ureg=self.units))               
1429                if backbone_vector is None:
1430                    bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1431                                                                radius=l0, 
1432                                                                n_samples=1,
1433                                                                on_surface=True)[0]
1434                else:
1435                    bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector),
1436                                                                                                magnitude=l0)
1437                    
1438                side_bead_id = self.create_particle(name=side_chain_name, 
1439                                                    espresso_system=espresso_system,
1440                                                    position=[bead_position], 
1441                                                    number_of_particles=1)[0]
1442                side_chain_beads_ids.append(side_bead_id)
1443                self.db._update_instance(pmb_type="particle",
1444                                         instance_id=side_bead_id,
1445                                         attribute="residue_id",
1446                                         value=residue_id)
1447                self.create_bond(particle_id1=central_bead_id,
1448                                 particle_id2=side_bead_id,
1449                                 espresso_system=espresso_system,
1450                                 use_default_bond=use_default_bond)
1451            elif pmb_type == 'residue':
1452                side_residue_tpl = self.db.get_template(name=side_chain_name,
1453                                                        pmb_type=pmb_type)
1454                central_bead_side_chain = side_residue_tpl.central_bead
1455                lj_parameters = self.get_lj_parameters(particle_name1=central_bead_name,
1456                                                       particle_name2=central_bead_side_chain)
1457                bond_tpl = self.get_bond_template(particle_name1=central_bead_name,
1458                                                  particle_name2=central_bead_side_chain,
1459                                                  use_default_bond=use_default_bond)
1460                l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1461                                                      bond_type=bond_tpl.bond_type,
1462                                                      bond_parameters=bond_tpl.get_parameters(ureg=self.units))
1463                if backbone_vector is None:
1464                    residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1465                                                                radius=l0, 
1466                                                                n_samples=1,
1467                                                                on_surface=True)[0]
1468                else:
1469                    residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1470                                                                                                    magnitude=l0)
1471                side_residue_id = self.create_residue(name=side_chain_name, 
1472                                                      espresso_system=espresso_system,
1473                                                      central_bead_position=[residue_position],
1474                                                      use_default_bond=use_default_bond)
1475                # Find particle ids of the inner residue
1476                side_chain_beads_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1477                                                                               attribute="residue_id",
1478                                                                               value=side_residue_id)
1479                # Change the residue_id of the residue in the side chain to the one of the outer residue
1480                for particle_id in side_chain_beads_ids:
1481                    self.db._update_instance(instance_id=particle_id,
1482                                             pmb_type="particle",
1483                                             attribute="residue_id",
1484                                             value=residue_id)
1485                # Remove the instance of the inner residue
1486                self.db.delete_instance(pmb_type="residue",
1487                                        instance_id=side_residue_id)
1488                self.create_bond(particle_id1=central_bead_id,
1489                                 particle_id2=side_chain_beads_ids[0],
1490                                 espresso_system=espresso_system,
1491                                 use_default_bond=use_default_bond)
1492        if gen_angle:
1493            self._generate_angles_for_entity(espresso_system=espresso_system,
1494                                             entity_id=residue_id,
1495                                             entity_id_col="residue_id")
1496        return  residue_id
1497
1498    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1499        """
1500        Defines bond templates for each particle pair in 'particle_pairs' in the pyMBE database.
1501
1502        Args:
1503            bond_type ('str'): 
1504                label to identify the potential to model the bond.
1505
1506            bond_parameters ('dict'): 
1507                parameters of the potential of the bond.
1508
1509            particle_pairs ('lst'): 
1510                list of the 'names' of the 'particles' to be bonded.
1511
1512        Notes:
1513            -Currently, only HARMONIC and FENE bonds are supported.
1514            - For a HARMONIC bond the dictionary must contain the following parameters:
1515                - k ('pint.Quantity')      : Magnitude of the bond. It should have units of energy/length**2 
1516                using the 'pmb.units' UnitRegistry.
1517                - r_0 ('pint.Quantity')    : Equilibrium bond length. It should have units of length using 
1518                the 'pmb.units' UnitRegistry.
1519           - For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:              
1520                - d_r_max ('pint.Quantity'): Maximal stretching length for FENE. It should have 
1521                units of length using the 'pmb.units' UnitRegistry. Default 'None'.
1522        """
1523        self._check_bond_inputs(bond_parameters=bond_parameters,
1524                                bond_type=bond_type)
1525        parameters_expected_dimensions={"r_0": "length",
1526                                        "k": "energy/length**2",
1527                                        "d_r_max": "length"}
1528
1529        parameters_tpl = {}
1530        for key in bond_parameters.keys():
1531            parameters_tpl[key]= PintQuantity.from_quantity(q=bond_parameters[key],
1532                                                            expected_dimension=parameters_expected_dimensions[key],
1533                                                            ureg=self.units)
1534
1535        bond_names=[]
1536        for particle_name1, particle_name2 in particle_pairs:
1537            
1538            tpl = BondTemplate(particle_name1=particle_name1,
1539                               particle_name2=particle_name2,
1540                               parameters=parameters_tpl,
1541                               bond_type=bond_type)
1542            tpl._make_name()
1543            if tpl.name in bond_names:
1544                raise RuntimeError(f"Bond {tpl.name} has already been defined, please check the list of particle pairs")
1545            bond_names.append(tpl.name)
1546            self.db._register_template(tpl)
1547
1548    
1549    def define_default_bond(self, bond_type, bond_parameters):
1550        """
1551        Defines a bond template as a "default" template in the pyMBE database.
1552        
1553        Args:
1554            bond_type ('str'): 
1555                label to identify the potential to model the bond.
1556
1557            bond_parameters ('dict'): 
1558                parameters of the potential of the bond.
1559            
1560        Notes:
1561            - Currently, only harmonic and FENE bonds are supported. 
1562        """
1563        self._check_bond_inputs(bond_parameters=bond_parameters,
1564                                bond_type=bond_type)
1565        parameters_expected_dimensions={"r_0": "length",
1566                                        "k": "energy/length**2",
1567                                        "d_r_max": "length"}
1568        parameters_tpl = {}
1569        for key in bond_parameters.keys():
1570            parameters_tpl[key]= PintQuantity.from_quantity(q=bond_parameters[key],
1571                                                            expected_dimension=parameters_expected_dimensions[key],
1572                                                            ureg=self.units)
1573        tpl = BondTemplate(parameters=parameters_tpl,
1574                               bond_type=bond_type)
1575        tpl.name = "default"
1576        self.db._register_template(tpl)
1577
1578    def define_angular_potential(self, angle_type, angle_parameters, particle_triplets):
1579        """
1580        Defines angle potential templates for each particle triplet in `particle_triplets`.
1581
1582        Args:
1583            angle_type ('str'):
1584                Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
1585
1586            angle_parameters ('dict'):
1587                Parameters of the angle potential. Must contain:
1588                    - "k" ('pint.Quantity'): Bending stiffness with dimensions of energy.
1589                    - "phi_0" ('float'): Equilibrium angle in radians.
1590
1591            particle_triplets ('list[tuple[str,str,str]]'):
1592                List of (side_particle1, central_particle, side_particle2) triplets.
1593        """
1594        valid_angle_types = ["harmonic", "cosine", "harmonic_cosine"]
1595        if angle_type not in valid_angle_types:
1596            raise NotImplementedError(f"Angle potential type '{angle_type}' currently not implemented in pyMBE, accepted types are {valid_angle_types}")
1597
1598        if "k" not in angle_parameters:
1599            raise ValueError("Magnitude of the angle potential (k) is missing")
1600        if "phi_0" not in angle_parameters:
1601            raise ValueError("Equilibrium angle (phi_0) is missing")
1602
1603        parameters_tpl = {
1604            "k": PintQuantity.from_quantity(q=angle_parameters["k"],
1605                                            expected_dimension="energy",
1606                                            ureg=self.units),
1607            "phi_0": PintQuantity.from_quantity(q=angle_parameters["phi_0"],
1608                                                expected_dimension="dimensionless",
1609                                                ureg=self.units),
1610        }
1611
1612        angle_names = []
1613        for side1, central, side2 in particle_triplets:
1614            tpl = AngleTemplate(side_particle1=side1,
1615                                central_particle=central,
1616                                side_particle2=side2,
1617                                parameters=parameters_tpl,
1618                                angle_type=angle_type)
1619            tpl._make_name()
1620            if tpl.name in angle_names:
1621                raise RuntimeError(f"Angle {tpl.name} has already been defined, please check the list of particle triplets")
1622            angle_names.append(tpl.name)
1623            self.db._register_template(tpl)
1624
1625    def define_default_angular_potential(self, angle_type, angle_parameters):
1626        """
1627        Defines an angle template as a "default" template in the pyMBE database.
1628
1629        Args:
1630            angle_type ('str'):
1631                Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
1632
1633            angle_parameters ('dict'):
1634                Parameters of the angle potential (k, phi_0).
1635        """
1636        valid_angle_types = ["harmonic", "cosine", "harmonic_cosine"]
1637        if angle_type not in valid_angle_types:
1638            raise NotImplementedError(f"Angle potential type '{angle_type}' currently not implemented in pyMBE, accepted types are {valid_angle_types}")
1639
1640        if "k" not in angle_parameters:
1641            raise ValueError("Magnitude of the angle potential (k) is missing")
1642        if "phi_0" not in angle_parameters:
1643            raise ValueError("Equilibrium angle (phi_0) is missing")
1644
1645        parameters_tpl = {
1646            "k": PintQuantity.from_quantity(q=angle_parameters["k"],
1647                                            expected_dimension="energy",
1648                                            ureg=self.units),
1649            "phi_0": PintQuantity.from_quantity(q=angle_parameters["phi_0"],
1650                                                expected_dimension="dimensionless",
1651                                                ureg=self.units),
1652        }
1653        tpl = AngleTemplate(parameters=parameters_tpl,
1654                            angle_type=angle_type)
1655        tpl.name = "default"
1656        self.db._register_template(tpl)
1657
1658    def create_angular_potential(self, particle_id1, particle_id2, particle_id3, espresso_system, use_default_angle=False):
1659        """
1660        Creates an angle between three particle instances in an ESPResSo system
1661        and registers it in the pyMBE database.
1662
1663        Args:
1664            particle_id1 ('int'): ID of the first side particle.
1665            particle_id2 ('int'): ID of the central particle.
1666            particle_id3 ('int'): ID of the second side particle.
1667            espresso_system ('espressomd.system.System'): ESPResSo system.
1668            use_default_angle ('bool', optional): If True, use the default angle if no specific one is found.
1669        """
1670        particle_inst_1 = self.db.get_instance(pmb_type="particle", instance_id=particle_id1)
1671        particle_inst_2 = self.db.get_instance(pmb_type="particle", instance_id=particle_id2)
1672        particle_inst_3 = self.db.get_instance(pmb_type="particle", instance_id=particle_id3)
1673
1674        # Verify that bonds exist between side particles and central particle
1675        bond_instances = self.db.get_instances(pmb_type="bond")
1676        bonded_pairs = set()
1677        for bond in bond_instances.values():
1678            pair = frozenset([bond.particle_id1, bond.particle_id2])
1679            bonded_pairs.add(pair)
1680        if frozenset([particle_id1, particle_id2]) not in bonded_pairs:
1681            raise ValueError(f"Cannot create angle: no bond exists between particle {particle_id1} and central particle {particle_id2}.")
1682        if frozenset([particle_id3, particle_id2]) not in bonded_pairs:
1683            raise ValueError(f"Cannot create angle: no bond exists between particle {particle_id3} and central particle {particle_id2}.")
1684
1685        angle_tpl = self.get_angle_template(side_name1=particle_inst_1.name,
1686                                            central_name=particle_inst_2.name,
1687                                            side_name2=particle_inst_3.name,
1688                                            use_default_angle=use_default_angle)
1689        angle_inst = self._get_espresso_angle_instance(angle_template=angle_tpl, espresso_system=espresso_system)
1690
1691        # ESPResSo angle bonds are added to the central particle
1692        espresso_system.part.by_id(particle_id2).add_bond((angle_inst, particle_id1, particle_id3))
1693
1694        angle_id = self.db._propose_instance_id(pmb_type="angle")
1695        pmb_angle_instance = AngleInstance(angle_id=angle_id,
1696                                           name=angle_tpl.name,
1697                                           particle_id1=particle_id1,
1698                                           particle_id2=particle_id2,
1699                                           particle_id3=particle_id3)
1700        self.db._register_instance(instance=pmb_angle_instance)
1701
1702    def get_angle_template(self, side_name1, central_name, side_name2, use_default_angle=False):
1703        """
1704        Retrieves an angle template connecting three particle templates.
1705
1706        Args:
1707            side_name1 ('str'): Name of the first side particle.
1708            central_name ('str'): Name of the central particle.
1709            side_name2 ('str'): Name of the second side particle.
1710            use_default_angle ('bool', optional): If True, fall back to the default angle template.
1711
1712        Returns:
1713            ('AngleTemplate'): The matching angle template.
1714        """
1715        angle_key = AngleTemplate.make_angle_key(side1=side_name1, central=central_name, side2=side_name2)
1716        try:
1717            return self.db.get_template(name=angle_key, pmb_type="angle")
1718        except ValueError:
1719            pass
1720
1721        if use_default_angle:
1722            return self.db.get_template(name="default", pmb_type="angle")
1723
1724        raise ValueError(f"No angle template found for '{side_name1}-{central_name}-{side_name2}', and default angles are deactivated.")
1725
1726    def _get_espresso_angle_instance(self, angle_template, espresso_system):
1727        """
1728        Retrieve or create an angle interaction in an ESPResSo system for a given angle template.
1729
1730        Args:
1731            angle_template ('AngleTemplate'): The angle template to use.
1732            espresso_system ('espressomd.system.System'): ESPResSo system.
1733
1734        Returns:
1735            ('espressomd.interactions.BondedInteraction'): The ESPResSo angle interaction object.
1736        """
1737        if angle_template.name in self.db.espresso_angle_instances:
1738            return self.db.espresso_angle_instances[angle_template.name]
1739        angle_inst = self._create_espresso_angle_instance(angle_type=angle_template.angle_type,
1740                                                          angle_parameters=angle_template.get_parameters(self.units))
1741        self.db.espresso_angle_instances[angle_template.name] = angle_inst
1742        espresso_system.bonded_inter.add(angle_inst)
1743        return angle_inst
1744
1745    def _create_espresso_angle_instance(self, angle_type, angle_parameters):
1746        """
1747        Creates an ESPResSo angle interaction object.
1748
1749        Args:
1750            angle_type ('str'): Type of angle potential ("harmonic", "cosine", "harmonic_cosine").
1751            angle_parameters ('dict'): Parameters of the angle potential (k, phi_0).
1752
1753        Returns:
1754            ('espressomd.interactions.BondedInteraction'): The ESPResSo angle interaction object.
1755        """
1756        from espressomd import interactions
1757
1758        k = angle_parameters["k"].m_as("reduced_energy")
1759        phi_0 = float(angle_parameters["phi_0"].magnitude)
1760
1761        if angle_type == "harmonic":
1762            return interactions.AngleHarmonic(bend=k, phi0=phi_0)
1763        elif angle_type == "cosine":
1764            return interactions.AngleCosine(bend=k, phi0=phi_0)
1765        elif angle_type == "harmonic_cosine":
1766            return interactions.AngleCossquare(bend=k, phi0=phi_0)
1767
1768    def _generate_angles_for_entity(self, espresso_system, entity_id, entity_id_col):
1769        """
1770        Auto-generates angles from bond topology for an entity (molecule or residue).
1771
1772        For each particle in the entity that has two or more bonded neighbors,
1773        this method finds all neighbor pairs and applies any matching angle potential.
1774
1775        Args:
1776            espresso_system ('espressomd.system.System'): ESPResSo system.
1777            entity_id ('int'): The molecule_id or residue_id to generate angles for.
1778            entity_id_col ('str'): Either "molecule_id" or "residue_id".
1779        """
1780        # Get all particle IDs for this entity
1781        particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1782                                                               attribute=entity_id_col,
1783                                                               value=entity_id)
1784        if not particle_ids:
1785            return
1786
1787        # Build neighbor map from bond instances
1788        neighbors = {pid: set() for pid in particle_ids}
1789        pid_set = set(particle_ids)
1790        bond_instances = self.db.get_instances(pmb_type="bond")
1791        for bond in bond_instances.values():
1792            i, j = bond.particle_id1, bond.particle_id2
1793            if i in pid_set and j in pid_set:
1794                neighbors[i].add(j)
1795                neighbors[j].add(i)
1796
1797        # For each particle with 2+ neighbors, generate angles
1798        for j in particle_ids:
1799            nbs = sorted(neighbors[j])
1800            if len(nbs) < 2:
1801                continue
1802
1803            for idx_i in range(len(nbs)):
1804                for idx_k in range(idx_i + 1, len(nbs)):
1805                    i = nbs[idx_i]
1806                    k = nbs[idx_k]
1807                    try:
1808                        self.create_angular_potential(particle_id1=i,
1809                                          particle_id2=j,
1810                                          particle_id3=k,
1811                                          espresso_system=espresso_system,
1812                                          use_default_angle=True)
1813                    except ValueError:
1814                        # No angle template defined for this triplet — skip
1815                        continue
1816
1817    def define_hydrogel(self, name, node_map, chain_map):
1818        """
1819        Defines a hydrogel template in the pyMBE database.
1820
1821        Args:
1822            name ('str'): 
1823                Unique label that identifies the 'hydrogel'.
1824
1825            node_map ('list of dict'): 
1826                [{"particle_name": , "lattice_index": }, ... ]
1827
1828            chain_map ('list of dict'): 
1829                [{"node_start": , "node_end": , "residue_list": , ... ]
1830        """
1831        # Sanity tests
1832        node_indices = {tuple(entry['lattice_index']) for entry in node_map}                
1833        chain_map_connectivity = set()
1834        for entry in chain_map:
1835            start = self.lattice_builder.node_labels[entry['node_start']]
1836            end = self.lattice_builder.node_labels[entry['node_end']]
1837            chain_map_connectivity.add((start,end))
1838        if self.lattice_builder.lattice.connectivity != chain_map_connectivity:
1839            raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs")
1840        diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices}
1841        if node_indices != diamond_indices:
1842            raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ")
1843        # Register information in the pyMBE database
1844        nodes=[]
1845        for entry in node_map:
1846            nodes.append(HydrogelNode(particle_name=entry["particle_name"],
1847                                      lattice_index=entry["lattice_index"]))
1848        chains=[]
1849        for chain in chain_map:
1850            chains.append(HydrogelChain(node_start=chain["node_start"],
1851                                        node_end=chain["node_end"],
1852                                        molecule_name=chain["molecule_name"]))
1853        tpl = HydrogelTemplate(name=name,
1854                               node_map=nodes,
1855                               chain_map=chains)
1856        self.db._register_template(tpl)
1857
1858    def define_molecule(self, name, residue_list):
1859        """
1860        Defines a molecule template in the pyMBE database.
1861
1862        Args:
1863            name('str'): 
1864                Unique label that identifies the 'molecule'.
1865
1866            residue_list ('list' of 'str'): 
1867                List of the 'name's of the 'residue's  in the sequence of the 'molecule'.  
1868        """
1869        tpl = MoleculeTemplate(name=name,
1870                               residue_list=residue_list)
1871        self.db._register_template(tpl)
1872
1873    def define_monoprototic_acidbase_reaction(self, particle_name, pka, acidity, metadata=None):
1874        """
1875        Defines an acid-base reaction for a monoprototic particle in the pyMBE database.
1876
1877        Args:
1878            particle_name ('str'): 
1879                Unique label that identifies the particle template. 
1880
1881            pka ('float'): 
1882                pka-value of the acid or base.
1883
1884            acidity ('str'): 
1885                Identifies whether if the particle is 'acidic' or 'basic'.
1886
1887            metadata ('dict', optional): 
1888                Additional information to be stored in the reaction. Defaults to None.
1889        """
1890        supported_acidities = ["acidic", "basic"]
1891        if acidity not in supported_acidities:
1892            raise ValueError(f"Unsupported acidity '{acidity}' for particle '{particle_name}'. Supported acidities are {supported_acidities}.")
1893        reaction_type = "monoprotic"
1894        if acidity == "basic":
1895            reaction_type += "_base"
1896        else:
1897            reaction_type += "_acid"
1898        reaction = Reaction(participants=[ReactionParticipant(particle_name=particle_name,
1899                                                              state_name=f"{particle_name}H", 
1900                                                              coefficient=-1),
1901                                          ReactionParticipant(particle_name=particle_name,
1902                                                              state_name=f"{particle_name}",
1903                                                              coefficient=1)],
1904                            reaction_type=reaction_type,
1905                            pK=pka,
1906                            metadata=metadata)
1907        self.db._register_reaction(reaction)
1908
1909    def define_monoprototic_particle_states(self, particle_name, acidity):
1910        """
1911        Defines particle states for a monoprotonic particle template including the charges in each of its possible states. 
1912
1913        Args:
1914            particle_name ('str'): 
1915                Unique label that identifies the particle template. 
1916
1917            acidity ('str'): 
1918                Identifies whether the particle is 'acidic' or 'basic'.
1919        """
1920        acidity_valid_keys = ['acidic', 'basic']
1921        if not pd.isna(acidity):
1922            if acidity not in acidity_valid_keys:
1923                raise ValueError(f"Acidity {acidity} provided for particle name  {particle_name} is not supported. Valid keys are: {acidity_valid_keys}")
1924        if acidity == "acidic":
1925            states = [{"name": f"{particle_name}H", "z": 0}, 
1926                      {"name": f"{particle_name}",  "z": -1}]
1927            
1928        elif acidity == "basic":
1929            states = [{"name": f"{particle_name}H", "z": 1}, 
1930                      {"name": f"{particle_name}",  "z": 0}]
1931        self.define_particle_states(particle_name=particle_name, 
1932                                    states=states)
1933
1934    def define_particle(self, name,  sigma, epsilon, z=0, acidity=pd.NA, pka=pd.NA, cutoff=pd.NA, offset=pd.NA):
1935        """
1936        Defines a particle template in the pyMBE database.
1937
1938        Args:
1939            name('str'):
1940                 Unique label that identifies this particle type.  
1941
1942            sigma('pint.Quantity'): 
1943                Sigma parameter used to set up Lennard-Jones interactions for this particle type. 
1944
1945            epsilon('pint.Quantity'): 
1946                Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe.
1947
1948            z('int', optional): 
1949                Permanent charge number of this particle type. Defaults to 0.
1950
1951            acidity('str', optional): 
1952                Identifies whether if the particle is 'acidic' or 'basic', used to setup constant pH simulations. Defaults to pd.NA.
1953
1954            pka('float', optional):
1955                If 'particle' is an acid or a base, it defines its  pka-value. Defaults to pd.NA.
1956
1957            cutoff('pint.Quantity', optional): 
1958                Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1959
1960            offset('pint.Quantity', optional): 
1961                Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1962            
1963        Notes:
1964            - 'sigma', 'cutoff' and 'offset' must have a dimensitonality of '[length]' and should be defined using pmb.units.
1965            - 'epsilon' must have a dimensitonality of '[energy]' and should be defined using pmb.units.
1966            - 'cutoff' defaults to '2**(1./6.) reduced_length'. 
1967            - 'offset' defaults to 0.
1968            - For more information on 'sigma', 'epsilon', 'cutoff' and 'offset' check 'pmb.setup_lj_interactions()'.
1969        """ 
1970        # If 'cutoff' and 'offset' are not defined, default them to the following values
1971        if pd.isna(cutoff):
1972            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1973        if pd.isna(offset):
1974            offset=self.units.Quantity(0, "reduced_length")
1975        # Define particle states
1976        if acidity is pd.NA:
1977            states = [{"name": f"{name}",  "z": z}]
1978            self.define_particle_states(particle_name=name, 
1979                                        states=states)
1980            initial_state = name
1981        else:
1982            self.define_monoprototic_particle_states(particle_name=name,
1983                                                  acidity=acidity)
1984            initial_state = f"{name}H"
1985            if pka is not pd.NA:
1986                self.define_monoprototic_acidbase_reaction(particle_name=name,
1987                                                           acidity=acidity,
1988                                                           pka=pka)
1989        tpl = ParticleTemplate(name=name, 
1990                               sigma=PintQuantity.from_quantity(q=sigma, expected_dimension="length", ureg=self.units), 
1991                               epsilon=PintQuantity.from_quantity(q=epsilon, expected_dimension="energy", ureg=self.units),
1992                               cutoff=PintQuantity.from_quantity(q=cutoff, expected_dimension="length", ureg=self.units), 
1993                               offset=PintQuantity.from_quantity(q=offset, expected_dimension="length", ureg=self.units),
1994                               initial_state=initial_state)
1995        self.db._register_template(tpl)
1996    
1997    def define_particle_states(self, particle_name, states):
1998        """
1999        Define the chemical states of an existing particle template.
2000
2001        Args:
2002            particle_name ('str'):
2003                Name of a particle template. 
2004
2005            states ('list' of 'dict'):
2006                List of dictionaries defining the particle states. Each dictionary
2007                must contain:
2008                - 'name' ('str'): Name of the particle state (e.g. '"H"', '"-"',
2009                '"neutral"').
2010                - 'z' ('int'): Charge number of the particle in this state.
2011                Example:
2012                states = [{"name": "AH", "z": 0},     # protonated
2013                         {"name": "A-", "z": -1}]    # deprotonated
2014        Notes:
2015            - Each state is assigned a unique Espresso 'es_type' automatically.
2016            - Chemical reactions (e.g. acid–base equilibria) are **not** created by
2017            this method and must be defined separately (e.g. via
2018            'set_particle_acidity()' or custom reaction definitions).
2019            - Particles without explicitly defined states are assumed to have a
2020            single, implicit state with their default charge.
2021        """
2022        for s in states:
2023            state = ParticleStateTemplate(particle_name=particle_name,
2024                                          name=s["name"],
2025                                          z=s["z"],
2026                                          es_type=self.propose_unused_type())
2027            self.db._register_template(state)
2028
2029    def define_peptide(self, name, sequence, model):
2030        """
2031        Defines a peptide template in the pyMBE database.
2032
2033        Args:
2034            name ('str'): 
2035                Unique label that identifies the peptide.
2036
2037            sequence ('str'): 
2038                Sequence of the peptide.
2039
2040            model ('str'): 
2041                Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
2042        """
2043        valid_keys = ['1beadAA','2beadAA']
2044        if model not in valid_keys:
2045            raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
2046        clean_sequence = hf.protein_sequence_parser(sequence=sequence)    
2047        residue_list = self._get_residue_list_from_sequence(sequence=clean_sequence)
2048        tpl = PeptideTemplate(name=name,
2049                            residue_list=residue_list,
2050                            model=model,
2051                            sequence=sequence)
2052        self.db._register_template(tpl)        
2053    
2054    def define_protein(self, name, sequence, model):
2055        """
2056        Defines a protein template in the pyMBE database.
2057
2058        Args:
2059            name ('str'): 
2060                Unique label that identifies the protein.
2061
2062            sequence ('str'): 
2063                Sequence of the protein.
2064
2065            model ('string'): 
2066                Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
2067
2068        Notes:
2069            - Currently, only 'lj_setup_mode="wca"' is supported. This corresponds to setting up the WCA potential.
2070        """
2071        valid_model_keys = ['1beadAA','2beadAA']
2072        if model not in valid_model_keys:
2073            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
2074        
2075        residue_list = self._get_residue_list_from_sequence(sequence=sequence)
2076        tpl = ProteinTemplate(name=name,
2077                              model=model,
2078                              residue_list=residue_list,
2079                              sequence=sequence)
2080        self.db._register_template(tpl)
2081    
2082    def define_residue(self, name, central_bead, side_chains):
2083        """
2084        Defines a residue template in the pyMBE database.
2085
2086        Args:
2087            name ('str'): 
2088                Unique label that identifies the residue.
2089
2090            central_bead ('str'): 
2091                'name' of the 'particle' to be placed as central_bead of the residue.
2092
2093            side_chains('list' of 'str'): 
2094                List of 'name's of the pmb_objects to be placed as side_chains of the residue. Currently, only pyMBE objects of type 'particle' or 'residue' are supported.
2095        """
2096        tpl = ResidueTemplate(name=name,
2097                              central_bead=central_bead,
2098                              side_chains=side_chains)
2099        self.db._register_template(tpl)
2100
2101    def delete_instances_in_system(self, instance_id, pmb_type, espresso_system):
2102        """
2103        Deletes the instance with instance_id from the ESPResSo system. 
2104        Related assembly, molecule, residue, particles and bond instances will also be deleted from the pyMBE dataframe.
2105
2106        Args:
2107            instance_id ('int'): 
2108                id of the assembly to be deleted. 
2109
2110            pmb_type ('str'): 
2111                the instance type to be deleted. 
2112
2113            espresso_system ('espressomd.system.System'): 
2114                Instance of a system class from espressomd library.
2115        """
2116        if pmb_type == "particle":
2117            instance_identifier = "particle_id"
2118        elif pmb_type == "residue":
2119            instance_identifier = "residue_id"
2120        elif pmb_type in self.db._molecule_like_types:
2121            instance_identifier = "molecule_id"
2122        elif pmb_type in self.db._assembly_like_types:
2123            instance_identifier = "assembly_id"
2124        particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
2125                                                               attribute=instance_identifier,
2126                                                               value=instance_id)
2127        self._delete_particles_from_espresso(particle_ids=particle_ids,
2128                                             espresso_system=espresso_system)
2129        self.db.delete_instance(pmb_type=pmb_type,
2130                                instance_id=instance_id)
2131
2132    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
2133        """
2134        Determines ionic concentrations in the reservoir at fixed pH and salt concentration.
2135
2136        Args:
2137            pH_res ('float'):
2138                Target pH value in the reservoir.
2139
2140            c_salt_res ('pint.Quantity'):
2141                Concentration of monovalent salt (e.g., NaCl) in the reservoir.
2142
2143            activity_coefficient_monovalent_pair ('callable'):
2144                Function returning the activity coefficient of a monovalent ion pair
2145                as a function of ionic strength:
2146                'gamma = activity_coefficient_monovalent_pair(I)'.
2147
2148            max_number_sc_runs ('int', optional):
2149                Maximum number of self-consistent iterations allowed before
2150                convergence is enforced. Defaults to 200.
2151
2152        Returns:
2153            tuple:
2154                (cH_res, cOH_res, cNa_res, cCl_res)
2155                - cH_res ('pint.Quantity'): Concentration of H⁺ ions.
2156                - cOH_res ('pint.Quantity'): Concentration of OH⁻ ions.
2157                - cNa_res ('pint.Quantity'): Concentration of Na⁺ ions.
2158                - cCl_res ('pint.Quantity'): Concentration of Cl⁻ ions.
2159
2160        Notess:
2161            - The algorithm enforces electroneutrality in the reservoir.
2162            - Water autodissociation is included via the equilibrium constant 'Kw'.
2163            - Non-ideal effects enter through activity coefficients depending on
2164            ionic strength.
2165            - The implementation follows the self-consistent scheme described in
2166            Landsgesell (PhD thesis, Sec. 5.3, doi:10.18419/opus-10831), adapted
2167            from the original code (doi:10.18419/darus-2237).
2168        """
2169        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
2170            """
2171            Iteratively determines reservoir ion concentrations self-consistently.
2172
2173            Args:
2174                cH_res ('pint.Quantity'):
2175                    Current estimate of the H⁺ concentration.
2176                c_salt_res ('pint.Quantity'):
2177                    Concentration of monovalent salt in the reservoir.
2178
2179            Returns:
2180                'tuple':
2181                    (cH_res, cOH_res, cNa_res, cCl_res)
2182            """
2183            # Initial ideal estimate
2184            cOH_res = self.Kw / cH_res
2185            if cOH_res >= cH_res:
2186                cNa_res = c_salt_res + (cOH_res - cH_res)
2187                cCl_res = c_salt_res
2188            else:
2189                cCl_res = c_salt_res + (cH_res - cOH_res)
2190                cNa_res = c_salt_res
2191            # Self-consistent iteration
2192            for _ in range(max_number_sc_runs):
2193                ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2194                cOH_new = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
2195                if cOH_new >= cH_res:
2196                    cNa_new = c_salt_res + (cOH_new - cH_res)
2197                    cCl_new = c_salt_res
2198                else:
2199                    cCl_new = c_salt_res + (cH_res - cOH_new)
2200                    cNa_new = c_salt_res
2201                # Update values
2202                cOH_res = cOH_new
2203                cNa_res = cNa_new
2204                cCl_res = cCl_new
2205            return cH_res, cOH_res, cNa_res, cCl_res
2206        # Initial guess for H+ concentration from target pH
2207        cH_res = 10 ** (-pH_res) * self.units.mol / self.units.l
2208        # First self-consistent solve
2209        cH_res, cOH_res, cNa_res, cCl_res = (determine_reservoir_concentrations_selfconsistently(cH_res, 
2210                                                                                                 c_salt_res))
2211        ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2212        determined_pH = -np.log10(cH_res.to("mol/L").magnitude* np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
2213        # Outer loop to enforce target pH
2214        while abs(determined_pH - pH_res) > 1e-6:
2215            if determined_pH > pH_res:
2216                cH_res *= 1.005
2217            else:
2218                cH_res /= 1.003
2219            cH_res, cOH_res, cNa_res, cCl_res = (determine_reservoir_concentrations_selfconsistently(cH_res, 
2220                                                                                                     c_salt_res))
2221            ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2222            determined_pH = -np.log10(cH_res.to("mol/L").magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
2223        return cH_res, cOH_res, cNa_res, cCl_res
2224
2225    def enable_motion_of_rigid_object(self, instance_id, pmb_type, espresso_system):
2226        """
2227        Enables translational and rotational motion of a rigid pyMBE object instance
2228        in an ESPResSo system.This method creates a rigid-body center particle at the center of mass of
2229        the specified pyMBE object and attaches all constituent particles to it
2230        using ESPResSo virtual sites. The resulting rigid object can translate and
2231        rotate as a single body.
2232
2233        Args:
2234            instance_id ('int'):
2235                Instance ID of the pyMBE object whose rigid-body motion is enabled.
2236
2237            pmb_type ('str'):
2238                pyMBE object type of the instance (e.g. '"molecule"', '"peptide"',
2239                '"protein"', or any assembly-like type).
2240
2241            espresso_system ('espressomd.system.System'):
2242                ESPResSo system in which the rigid object is defined.
2243
2244        Notess:
2245            - This method requires ESPResSo to be compiled with the following
2246            features enabled:
2247                - '"VIRTUAL_SITES_RELATIVE"'
2248                - '"MASS"'
2249            - A new ESPResSo particle is created to represent the rigid-body center.
2250            - The mass of the rigid-body center is set to the number of particles
2251            belonging to the object.
2252            - The rotational inertia tensor is approximated from the squared
2253            distances of the particles to the center of mass.
2254        """
2255        logging.info('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
2256        inst = self.db.get_instance(pmb_type=pmb_type,
2257                                    instance_id=instance_id)
2258        label = self._get_label_id_map(pmb_type=pmb_type)
2259        particle_ids_list = self.get_particle_id_map(object_name=inst.name)[label][instance_id]
2260        center_of_mass = self.calculate_center_of_mass (instance_id=instance_id,
2261                                                        espresso_system=espresso_system,
2262                                                        pmb_type=pmb_type)
2263        rigid_object_center = espresso_system.part.add(pos=center_of_mass,
2264                                                        rotation=[True,True,True], 
2265                                                        type=self.propose_unused_type())
2266        rigid_object_center.mass = len(particle_ids_list)
2267        momI = 0
2268        for pid in particle_ids_list:
2269            momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
2270        rigid_object_center.rinertia = np.ones(3) * momI        
2271        for particle_id in particle_ids_list:
2272            pid = espresso_system.part.by_id(particle_id)
2273            pid.vs_auto_relate_to(rigid_object_center.id)
2274
2275    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2276        """
2277        Generates random coordinates outside a sphere and inside a larger bounding sphere.
2278
2279        Args:
2280            center ('array-like'):
2281                Coordinates of the center of the spheres.
2282
2283            radius ('float'):
2284                Radius of the inner exclusion sphere. Must be positive.
2285
2286            max_dist ('float'):
2287                Radius of the outer sampling sphere. Must be larger than 'radius'.
2288
2289            n_samples ('int'):
2290                Number of coordinates to generate.
2291
2292        Returns:
2293            'list' of 'numpy.ndarray':
2294                List of coordinates lying outside the inner sphere and inside the
2295                outer sphere.
2296
2297        Notess:
2298            - Points are uniformly sampled inside a sphere of radius 'max_dist' centered at 'center' 
2299            and only those with a distance greater than or equal to 'radius' from the center are retained.
2300        """
2301        if not radius > 0: 
2302            raise ValueError (f'The value of {radius} must be a positive value')
2303        if not radius < max_dist:
2304            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
2305        coord_list = []
2306        counter = 0
2307        while counter<n_samples:
2308            coord = self.generate_random_points_in_a_sphere(center=center, 
2309                                            radius=max_dist,
2310                                            n_samples=1)[0]
2311            if np.linalg.norm(coord-np.asarray(center))>=radius:
2312                coord_list.append (coord)
2313                counter += 1
2314        return coord_list
2315    
2316    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2317        """
2318        Generates uniformly distributed random points inside or on the surface of a sphere.
2319
2320        Args:
2321            center ('array-like'):
2322                Coordinates of the center of the sphere.
2323
2324            radius ('float'):
2325                Radius of the sphere.
2326
2327            n_samples ('int'):
2328                Number of sample points to generate.
2329
2330            on_surface ('bool', optional):
2331                If True, points are uniformly sampled on the surface of the sphere.
2332                If False, points are uniformly sampled within the sphere volume.
2333                Defaults to False.
2334
2335        Returns:
2336            'numpy.ndarray':
2337                Array of shape '(n_samples, d)' containing the generated coordinates,
2338                where 'd' is the dimensionality of 'center'.
2339        Notes:
2340            - Points are sampled in a space whose dimensionality is inferred 
2341            from the length of 'center'.
2342        """
2343        # initial values
2344        center=np.array(center)
2345        d = center.shape[0]
2346        # sample n_samples points in d dimensions from a standard normal distribution
2347        samples = self.rng.normal(size=(n_samples, d))
2348        # make the samples lie on the surface of the unit hypersphere
2349        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
2350        samples /= normalize_radii
2351        if not on_surface:
2352            # make the samples lie inside the hypersphere with the correct density
2353            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
2354            new_radii = np.power(uniform_points, 1/d)
2355            samples *= new_radii
2356        # scale the points to have the correct radius and center
2357        samples = samples * radius + center
2358        return samples 
2359
2360    def generate_trial_perpendicular_vector(self,vector,magnitude):
2361        """
2362        Generates a random vector perpendicular to a given vector.
2363
2364        Args:
2365            vector ('array-like'):
2366                Reference vector to which the generated vector will be perpendicular.
2367
2368            magnitude ('float'):
2369                Desired magnitude of the perpendicular vector.
2370
2371        Returns:
2372            'numpy.ndarray':
2373                Vector orthogonal to 'vector' with norm equal to 'magnitude'.
2374        """ 
2375        np_vec = np.array(vector) 
2376        if np.all(np_vec == 0):
2377            raise ValueError('Zero vector')
2378        np_vec /= np.linalg.norm(np_vec) 
2379        # Generate a random vector 
2380        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
2381                                                                center=[0,0,0],
2382                                                                n_samples=1, 
2383                                                                on_surface=True)[0]
2384        # Project the random vector onto the input vector and subtract the projection
2385        projection = np.dot(random_vector, np_vec) * np_vec
2386        perpendicular_vector = random_vector - projection
2387        # Normalize the perpendicular vector to have the same magnitude as the input vector
2388        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
2389        return perpendicular_vector*magnitude
2390            
2391    def get_bond_template(self, particle_name1, particle_name2, use_default_bond=False) :
2392        """
2393        Retrieves a bond template connecting two particle templates.
2394
2395        Args:
2396            particle_name1 ('str'):
2397                Name of the first particle template.
2398
2399            particle_name2 ('str'):
2400                Name of the second particle template.
2401
2402            use_default_bond ('bool', optional):
2403                If True, returns the default bond template when no specific bond
2404                template is found. Defaults to False.
2405
2406        Returns:
2407            'BondTemplate':
2408                Bond template object retrieved from the pyMBE database.
2409            
2410        Notes:
2411            - This method searches the pyMBE database for a bond template defined between particle templates with names 'particle_name1' and 'particle_name2'. 
2412            - If no specific bond template is found and 'use_default_bond' is enabled, a default bond template is returned instead.
2413        """
2414        # Try to find a specific bond template
2415        bond_key = BondTemplate.make_bond_key(pn1=particle_name1,
2416                                              pn2=particle_name2)
2417        try:
2418            return self.db.get_template(name=bond_key, 
2419                                        pmb_type="bond")
2420        except ValueError:
2421            pass
2422
2423        #  Fallback to default bond if allowed
2424        if use_default_bond:
2425            return self.db.get_template(name="default", 
2426                                        pmb_type="bond")
2427
2428        # No bond template found
2429        raise ValueError(f"No bond template found between '{particle_name1}' and '{particle_name2}', and default bonds are deactivated.")
2430    
2431    def get_charge_number_map(self):
2432        """
2433        Construct a mapping from ESPResSo particle types to their charge numbers.
2434
2435        Returns:
2436            'dict[int, float]':
2437                Dictionary mapping ESPResSo particle types to charge numbers,
2438                ''{es_type: z}''.
2439
2440        Notess:
2441            - The mapping is built from particle *states*, not instances.
2442            - If multiple templates define states with the same ''es_type'',
2443            the last encountered definition will overwrite previous ones.
2444            This behavior is intentional and assumes database consistency.
2445            - Neutral particles (''z = 0'') are included in the map.
2446        """
2447        charge_number_map = {}
2448        particle_templates = self.db.get_templates("particle")
2449        for tpl in particle_templates.values():
2450            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
2451                charge_number_map[state.es_type] = state.z
2452        return charge_number_map
2453
2454    def get_instances_df(self, pmb_type):
2455        """
2456        Returns a dataframe with all instances of type 'pmb_type' in the pyMBE database.
2457
2458        Args:
2459            pmb_type ('str'): 
2460                pmb type to search instances in the pyMBE database.
2461        
2462        Returns:
2463            ('Pandas.Dataframe'): 
2464                Dataframe with all instances of type 'pmb_type'.
2465        """
2466        return self.db._get_instances_df(pmb_type=pmb_type)
2467
2468    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2469        """
2470        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2471        'particle_name1' and 'particle_name2' in the pyMBE database, calculated according to the provided combining rule.
2472
2473        Args:
2474            particle_name1 ('str'): 
2475                label of the type of the first particle type
2476
2477            particle_name2 ('str'): 
2478                label of the type of the second particle type
2479
2480            combining_rule ('string', optional): 
2481                combining rule used to calculate 'sigma' and 'epsilon' for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2482
2483        Returns:
2484            ('dict'):
2485                {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2486
2487        Notes:
2488            - Currently, the only 'combining_rule' supported is Lorentz-Berthelot.
2489            - If the sigma value of 'particle_name1' or 'particle_name2' is 0, the function will return an empty dictionary. No LJ interactions are set up for particles with sigma = 0.
2490        """
2491        supported_combining_rules=["Lorentz-Berthelot"]
2492        if combining_rule not in supported_combining_rules:
2493            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2494        part_tpl1 = self.db.get_template(name=particle_name1,
2495                                         pmb_type="particle")
2496        part_tpl2 = self.db.get_template(name=particle_name2,
2497                                         pmb_type="particle")
2498        lj_parameters1 = part_tpl1.get_lj_parameters(ureg=self.units)
2499        lj_parameters2 = part_tpl2.get_lj_parameters(ureg=self.units)
2500
2501        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2502        if part_tpl1.sigma.magnitude == 0 or part_tpl2.sigma.magnitude == 0:
2503            return {}
2504        # Apply combining rule
2505        if combining_rule == 'Lorentz-Berthelot':
2506            sigma=(lj_parameters1["sigma"]+lj_parameters2["sigma"])/2
2507            cutoff=(lj_parameters1["cutoff"]+lj_parameters2["cutoff"])/2
2508            offset=(lj_parameters1["offset"]+lj_parameters2["offset"])/2
2509            epsilon=np.sqrt(lj_parameters1["epsilon"]*lj_parameters2["epsilon"])
2510        return {"sigma": sigma, "cutoff": cutoff, "offset": offset, "epsilon": epsilon}    
2511
2512    def get_particle_id_map(self, object_name):
2513        """
2514        Collect all particle IDs associated with an object of given name in the
2515        pyMBE database. 
2516
2517        Args:
2518            object_name ('str'): 
2519                Name of the object.
2520
2521        Returns:
2522            ('dict'): 
2523                {"all": [particle_ids],
2524                 "residue_map": {residue_id: [particle_ids]},
2525                 "molecule_map": {molecule_id: [particle_ids]},
2526                 "assembly_map": {assembly_id: [particle_ids]},}
2527
2528        Notess:
2529            - Works for all supported pyMBE templates.
2530            - Relies in the internal method Manager.get_particle_id_map, see method for the detailed code.
2531        """
2532        return self.db.get_particle_id_map(object_name=object_name)
2533
2534    def get_pka_set(self):
2535        """
2536        Retrieve the pKa set for all titratable particles in the pyMBE database.
2537
2538        Returns:
2539            ('dict'): 
2540                Dictionary of the form:
2541                {"particle_name": {"pka_value": float,
2542                                   "acidity": "acidic" | "basic"}}
2543        Notes:
2544            - If a particle participates in multiple acid/base reactions, an error is raised.
2545        """
2546        pka_set = {}
2547        supported_reactions = ["monoprotic_acid",
2548                               "monoprotic_base"]
2549        for reaction in self.db._reactions.values():
2550            if reaction.reaction_type not in supported_reactions:
2551                continue
2552            # Identify involved particle(s)
2553            particle_names = {participant.particle_name for participant in reaction.participants}
2554            particle_name = particle_names.pop()
2555            if particle_name in pka_set:
2556                raise ValueError(f"Multiple acid/base reactions found for particle '{particle_name}'.")
2557            pka_set[particle_name] = {"pka_value": reaction.pK}
2558            if reaction.reaction_type == "monoprotic_acid":
2559                acidity = "acidic"
2560            elif reaction.reaction_type == "monoprotic_base":
2561                acidity = "basic"
2562            pka_set[particle_name]["acidity"] = acidity
2563        return pka_set
2564    
2565    def get_radius_map(self, dimensionless=True):
2566        """
2567        Gets the effective radius of each particle defined in the pyMBE database. 
2568
2569        Args:
2570            dimensionless ('bool'):
2571                If ``True``, return magnitudes expressed in ``reduced_length``.
2572                If ``False``, return Pint quantities with units.
2573        
2574        Returns:
2575            ('dict'): 
2576                {espresso_type: radius}.
2577
2578        Notes:
2579            - The radius corresponds to (sigma+offset)/2
2580        """
2581        if "particle" not in self.db._templates:
2582            return {}          
2583        result = {}
2584        for _, tpl in self.db._templates["particle"].items():
2585            radius = (tpl.sigma.to_quantity(self.units) + tpl.offset.to_quantity(self.units))/2.0
2586            if dimensionless:
2587                magnitude_reduced_length = radius.m_as("reduced_length")
2588                radius = magnitude_reduced_length
2589            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
2590                result[state.es_type] = radius
2591        return result
2592
2593    def get_reactions_df(self):
2594        """
2595        Returns a dataframe with all reaction templates in the pyMBE database.
2596
2597        Returns:
2598            (Pandas.Dataframe): 
2599                Dataframe with all  reaction templates.
2600        """
2601        return self.db._get_reactions_df()
2602
2603    def get_reduced_units(self):
2604        """
2605        Returns the  current set of reduced units defined in pyMBE.
2606
2607        Returns:
2608            reduced_units_text ('str'): 
2609                text with information about the current set of reduced units.
2610
2611        """
2612        unit_length=self.units.Quantity(1,'reduced_length')
2613        unit_energy=self.units.Quantity(1,'reduced_energy')
2614        unit_charge=self.units.Quantity(1,'reduced_charge')
2615        reduced_units_text = "\n".join(["Current set of reduced units:",
2616                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2617                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2618                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2619                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2620                                        ])   
2621        return reduced_units_text
2622
2623    def get_templates_df(self, pmb_type):
2624        """
2625        Returns a dataframe with all templates of type 'pmb_type' in the pyMBE database.
2626
2627        Args:
2628            pmb_type ('str'): 
2629                pmb type to search templates in the pyMBE database.
2630        
2631        Returns:
2632            ('Pandas.Dataframe'): 
2633                Dataframe with all templates of type given by 'pmb_type'.
2634        """
2635        return self.db._get_templates_df(pmb_type=pmb_type)
2636
2637    def get_type_map(self):
2638        """
2639        Return the mapping of ESPResSo types for all particle states defined in the pyMBE database.
2640        
2641        Returns:
2642            'dict[str, int]':
2643                A dictionary mapping each particle state to its corresponding ESPResSo type:
2644                {state_name: es_type, ...}
2645        """
2646        
2647        return self.db.get_es_types_map()
2648
2649    def initialize_lattice_builder(self, diamond_lattice):
2650        """
2651        Initialize the lattice builder with the DiamondLattice object.
2652
2653        Args:
2654            diamond_lattice ('DiamondLattice'): 
2655                DiamondLattice object from the 'lib/lattice' module to be used in the LatticeBuilder.
2656        """
2657        from .lib.lattice import LatticeBuilder, DiamondLattice
2658        if not isinstance(diamond_lattice, DiamondLattice):
2659            raise TypeError("Currently only DiamondLattice objects are supported.")
2660        self.lattice_builder = LatticeBuilder(lattice=diamond_lattice)
2661        logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}")
2662        return self.lattice_builder
2663
2664    def load_database(self, folder, format='csv'):
2665        """
2666        Loads a pyMBE database stored in 'folder'.
2667
2668        Args:
2669            folder ('str' or 'Path'): 
2670                Path to the folder where the pyMBE database was stored.
2671
2672            format ('str', optional): 
2673                Format of the database to be loaded. Defaults to 'csv'.
2674
2675        Return:
2676            ('dict'): 
2677                metadata with additional information about the source of the information in the database.
2678
2679        Notes:
2680            - The folder must contain the files generated by 'pmb.save_database()'.
2681            - Currently, only 'csv' format is supported.
2682        """
2683        supported_formats = ['csv']
2684        if format not in supported_formats:
2685            raise ValueError(f"Format {format} not supported. Supported formats are {supported_formats}")
2686        if format == 'csv':
2687            metadata =io._load_database_csv(self.db, 
2688                                            folder=folder)
2689        return metadata
2690        
2691    
2692    def load_pka_set(self, filename):
2693        """
2694        Load a pKa set and attach chemical states and acid–base reactions
2695        to existing particle templates.
2696
2697        Args:
2698            filename ('str'): 
2699                Path to a JSON file containing the pKa set. Expected format:
2700                {"metadata": {...},
2701                  "data": {"A": {"acidity": "acidic", "pka_value": 4.5},
2702                           "B": {"acidity": "basic",  "pka_value": 9.8}}}
2703
2704        Returns:
2705            ('dict'): 
2706                Dictionary with bibliographic metadata about the original work were the pKa set was determined.
2707
2708        Notes:
2709            - This method is designed for monoprotic acids and bases only.
2710        """
2711        with open(filename, "r") as f:
2712            pka_data = json.load(f)
2713        pka_set = pka_data["data"]
2714        metadata = pka_data.get("metadata", {})
2715        self._check_pka_set(pka_set)
2716        for particle_name, entry in pka_set.items():
2717            acidity = entry["acidity"]
2718            pka = entry["pka_value"]
2719            self.define_monoprototic_acidbase_reaction(particle_name=particle_name,
2720                                                       pka=pka,
2721                                                       acidity=acidity,
2722                                                       metadata=metadata)
2723        return metadata
2724            
2725    def propose_unused_type(self):
2726        """
2727        Propose an unused ESPResSo particle type.
2728
2729        Returns:
2730            ('int'): 
2731                The next available integer ESPResSo type. Returns ''0'' if no integer types are currently defined.
2732        """
2733        type_map = self.get_type_map()
2734        # Flatten all es_type values across all particles and states
2735        all_types = []
2736        for es_type in type_map.values():
2737            all_types.append(es_type)
2738        # If no es_types exist, start at 0
2739        if not all_types:
2740            return 0
2741        return max(all_types) + 1
2742       
2743    def read_protein_vtf(self, filename, unit_length=None):
2744        """
2745        Loads a coarse-grained protein model from a VTF file.
2746
2747        Args:
2748            filename ('str'): 
2749                Path to the VTF file.
2750                
2751            unit_length ('Pint.Quantity'): 
2752                Unit of length for coordinates (pyMBE UnitRegistry). Defaults to Angstrom.
2753
2754        Returns:
2755            ('tuple'):
2756                ('dict'): Particle topology.    
2757                ('str'):  One-letter amino-acid sequence (including n/c ends).
2758        """
2759        logging.info(f"Loading protein coarse-grain model file: {filename}")
2760        if unit_length is None:
2761            unit_length = 1 * self.units.angstrom
2762        atoms = {}        # atom_id -> atom info
2763        coords = []       # ordered coordinates
2764        residues = {}     # resid -> resname (first occurrence)
2765        has_n_term = False
2766        has_c_term = False
2767        aa_3to1 = {"ALA": "A", "ARG": "R", "ASN": "N", "ASP": "D",
2768                   "CYS": "C", "GLU": "E", "GLN": "Q", "GLY": "G",
2769                   "HIS": "H", "ILE": "I", "LEU": "L", "LYS": "K",
2770                   "MET": "M", "PHE": "F", "PRO": "P", "SER": "S",
2771                   "THR": "T", "TRP": "W", "TYR": "Y", "VAL": "V",
2772                   "n": "n", "c": "c"}
2773        # --- parse VTF ---
2774        with open(filename, "r") as f:
2775            for line in f:
2776                fields = line.split()
2777                if not fields:
2778                    continue
2779                if fields[0] == "atom":
2780                    atom_id = int(fields[1])
2781                    atom_name = fields[3]
2782                    resname = fields[5]
2783                    resid = int(fields[7])
2784                    chain_id = fields[9]
2785                    radius = float(fields[11]) * unit_length
2786                    atoms[atom_id] = {"name": atom_name,
2787                                     "resname": resname,
2788                                     "resid": resid,
2789                                     "chain_id": chain_id,
2790                                     "radius": radius}
2791                    if resname == "n":
2792                        has_n_term = True
2793                    elif resname == "c":
2794                        has_c_term = True
2795                    # register residue 
2796                    if resid not in residues:
2797                        residues[resid] = resname
2798                elif fields[0].isnumeric():
2799                    xyz = [(float(x) * unit_length).to("reduced_length").magnitude
2800                        for x in fields[1:4]]
2801                    coords.append(xyz)
2802        sequence = ""
2803        # N-terminus
2804        if has_n_term:
2805            sequence += "n"
2806        # protein residues only
2807        protein_resids = sorted(resid for resid, resname in residues.items()  if resname not in ("n", "c", "Ca"))
2808        for resid in protein_resids:
2809            resname = residues[resid]
2810            try:
2811                sequence += aa_3to1[resname]
2812            except KeyError:
2813                raise ValueError(f"Unknown residue name '{resname}' in VTF file")
2814        # C-terminus
2815        if has_c_term:
2816            sequence += "c"
2817        last_resid = max(protein_resids)
2818        # --- build topology ---
2819        topology_dict = {}
2820        for atom_id in sorted(atoms.keys()):
2821            atom = atoms[atom_id]
2822            resname = atom["resname"]
2823            resid = atom["resid"]
2824            # apply labeling rules
2825            if resname == "n":
2826                label_resid = 0
2827            elif resname == "c":
2828                label_resid = last_resid + 1
2829            elif resname == "Ca":
2830                label_resid = last_resid + 2
2831            else:
2832                label_resid = resid  # preserve original resid 
2833            label = f"{atom['name']}{label_resid}"
2834            if label in topology_dict:
2835                raise ValueError(f"Duplicate particle label '{label}'. Check VTF residue definitions.")
2836            topology_dict[label] = {"initial_pos": coords[atom_id - 1], "chain_id": atom["chain_id"], "radius": atom["radius"],}
2837        return topology_dict, sequence
2838
2839    
2840    def save_database(self, folder, format='csv'):
2841        """
2842        Saves the current pyMBE database into a file 'filename'.
2843
2844        Args:
2845            folder ('str' or 'Path'): 
2846                Path to the folder where the database files will be saved.
2847
2848        """
2849        supported_formats = ['csv']
2850        if format not in supported_formats:
2851            raise ValueError(f"Format {format} not supported. Supported formats are: {supported_formats}")
2852        if format == 'csv':
2853            io._save_database_csv(self.db, 
2854                                folder=folder)
2855
2856    def set_particle_initial_state(self, particle_name, state_name):
2857        """
2858        Sets the default initial state of a particle template defined in the pyMBE database.
2859
2860        Args:
2861            particle_name ('str'): 
2862                Unique label that identifies the particle template. 
2863
2864            state_name ('str'): 
2865                Name of the state to be set as default initial state.
2866        """
2867        part_tpl = self.db.get_template(name=particle_name,
2868        
2869                                        pmb_type="particle")
2870        part_tpl.initial_state = state_name
2871        logging.info(f"Default initial state of particle {particle_name} set to {state_name}.")
2872
2873    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None):
2874        """
2875        Sets the set of reduced units used by pyMBE.units and it prints it.
2876
2877        Args:
2878            unit_length ('pint.Quantity', optional): 
2879                Reduced unit of length defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2880
2881            unit_charge ('pint.Quantity', optional): 
2882                Reduced unit of charge defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2883
2884            temperature ('pint.Quantity', optional): 
2885                Temperature of the system, defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2886
2887            Kw ('pint.Quantity', optional): 
2888                Ionic product of water in mol^2/l^2. Defaults to None. 
2889
2890        Notes:
2891            - If no 'temperature' is given, a value of 298.15 K is assumed by default.
2892            - If no 'unit_length' is given, a value of 0.355 nm is assumed by default.
2893            - If no 'unit_charge' is given, a value of 1 elementary charge is assumed by default. 
2894            - If no 'Kw' is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2895        """
2896        if unit_length is None:
2897            unit_length= 0.355*self.units.nm
2898        if temperature is None:
2899            temperature = 298.15 * self.units.K
2900        if unit_charge is None:
2901            unit_charge = scipy.constants.e * self.units.C
2902        if Kw is None:
2903            Kw = 1e-14
2904        # Sanity check
2905        variables=[unit_length,temperature,unit_charge]
2906        dimensionalities=["[length]","[temperature]","[charge]"]
2907        for variable,dimensionality in zip(variables,dimensionalities):
2908            self._check_dimensionality(variable,dimensionality)
2909        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2910        self.kT=temperature*self.kB
2911        self.units._build_cache()
2912        self.units.define(f'reduced_energy = {self.kT} ')
2913        self.units.define(f'reduced_length = {unit_length}')
2914        self.units.define(f'reduced_charge = {unit_charge}')
2915        logging.info(self.get_reduced_units())
2916
2917    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, use_exclusion_radius_per_type = False):
2918        """
2919        Sets up the Acid/Base reactions for acidic/basic particles defined in the pyMBE database
2920        to be sampled in the constant pH ensemble. 
2921
2922        Args:
2923            counter_ion ('str'): 
2924                'name' of the counter_ion 'particle'.
2925
2926            constant_pH ('float'): 
2927                pH-value.
2928
2929            exclusion_range ('pint.Quantity', optional): 
2930                Below this value, no particles will be inserted.
2931
2932            use_exclusion_radius_per_type ('bool', optional): 
2933                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
2934
2935        Returns:
2936            ('reaction_methods.ConstantpHEnsemble'): 
2937                Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2938        """
2939        from espressomd import reaction_methods
2940        if exclusion_range is None:
2941            exclusion_range = max(self.get_radius_map().values())*2.0
2942        if use_exclusion_radius_per_type:
2943            exclusion_radius_per_type = self.get_radius_map()
2944        else:
2945            exclusion_radius_per_type = {}
2946        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2947                                                exclusion_range=exclusion_range, 
2948                                                seed=self.seed, 
2949                                                constant_pH=constant_pH,
2950                                                exclusion_radius_per_type = exclusion_radius_per_type)
2951        conterion_tpl = self.db.get_template(name=counter_ion,
2952                                             pmb_type="particle")
2953        conterion_state = self.db.get_template(name=conterion_tpl.initial_state,
2954                                               pmb_type="particle_state")
2955        for reaction in self.db.get_reactions():
2956            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
2957                continue
2958            default_charges = {}
2959            reactant_types  = []
2960            product_types   = []
2961            for participant in reaction.participants:
2962                state_tpl = self.db.get_template(name=participant.state_name,
2963                                                 pmb_type="particle_state")
2964                default_charges[state_tpl.es_type] = state_tpl.z
2965                if participant.coefficient < 0:
2966                    reactant_types.append(state_tpl.es_type)
2967                elif participant.coefficient > 0:
2968                    product_types.append(state_tpl.es_type)
2969            # Add counterion to the products
2970            if conterion_state.es_type not in product_types:
2971                product_types.append(conterion_state.es_type)
2972                default_charges[conterion_state.es_type] = conterion_state.z
2973                reaction.add_participant(particle_name=counter_ion,
2974                                         state_name=conterion_tpl.initial_state,
2975                                         coefficient=1)
2976            gamma=10**-reaction.pK
2977            RE.add_reaction(gamma=gamma,
2978                            reactant_types=reactant_types,
2979                            product_types=product_types,
2980                            default_charges=default_charges)
2981            reaction.add_simulation_method(simulation_method="cpH")
2982        return RE
2983
2984    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2985        """
2986        Sets up grand-canonical coupling to a reservoir of salt.
2987        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2988
2989        Args:
2990            c_salt_res ('pint.Quantity'): 
2991                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2992
2993            salt_cation_name ('str'): 
2994                Name of the salt cation (e.g. Na+) particle.
2995
2996            salt_anion_name ('str'): 
2997                Name of the salt anion (e.g. Cl-) particle.
2998
2999            activity_coefficient ('callable'): 
3000                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3001
3002            exclusion_range('pint.Quantity', optional): 
3003                For distances shorter than this value, no particles will be inserted.
3004
3005            use_exclusion_radius_per_type('bool',optional): 
3006                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
3007
3008        Returns:
3009            ('reaction_methods.ReactionEnsemble'): 
3010                Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
3011        """
3012        from espressomd import reaction_methods
3013        if exclusion_range is None:
3014            exclusion_range = max(self.get_radius_map().values())*2.0
3015        if use_exclusion_radius_per_type:
3016            exclusion_radius_per_type = self.get_radius_map()
3017        else:
3018            exclusion_radius_per_type = {}
3019        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3020                                               exclusion_range=exclusion_range, 
3021                                               seed=self.seed, 
3022                                               exclusion_radius_per_type = exclusion_radius_per_type)
3023        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3024        determined_activity_coefficient = activity_coefficient(c_salt_res)
3025        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
3026        cation_tpl = self.db.get_template(pmb_type="particle",
3027                                          name=salt_cation_name)
3028        cation_state = self.db.get_template(pmb_type="particle_state",
3029                                            name=cation_tpl.initial_state)
3030        anion_tpl = self.db.get_template(pmb_type="particle",
3031                                          name=salt_anion_name)
3032        anion_state = self.db.get_template(pmb_type="particle_state",
3033                                            name=anion_tpl.initial_state)
3034        salt_cation_es_type = cation_state.es_type
3035        salt_anion_es_type = anion_state.es_type     
3036        salt_cation_charge = cation_state.z
3037        salt_anion_charge = anion_state.z
3038        if salt_cation_charge <= 0:
3039            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3040        if salt_anion_charge >= 0:
3041            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3042        # Grand-canonical coupling to the reservoir
3043        RE.add_reaction(gamma = K_salt.magnitude,
3044                        reactant_types = [],
3045                        reactant_coefficients = [],
3046                        product_types = [ salt_cation_es_type, salt_anion_es_type ],
3047                        product_coefficients = [ 1, 1 ],
3048                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3049                                           salt_anion_es_type: salt_anion_charge})
3050        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3051                                                            state_name=cation_state.name,
3052                                                            coefficient=1),
3053                                        ReactionParticipant(particle_name=salt_anion_name,
3054                                                            state_name=anion_state.name,
3055                                                            coefficient=1)],
3056                           pK=-np.log10(K_salt.magnitude),
3057                           reaction_type="ion_insertion",
3058                           simulation_method="GCMC")
3059        self.db._register_reaction(rx_tpl)
3060        return RE
3061
3062    def setup_grxmc_reactions(self, pH_res, c_salt_res, proton_name, hydroxide_name, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
3063        """
3064        Sets up acid/base reactions for acidic/basic monoprotic particles defined in the pyMBE database, 
3065        as well as a grand-canonical coupling to a reservoir of small ions. 
3066        
3067        Args:
3068            pH_res ('float'): 
3069                pH-value in the reservoir.
3070
3071            c_salt_res ('pint.Quantity'): 
3072                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3073
3074            proton_name ('str'): 
3075                Name of the proton (H+) particle.
3076
3077            hydroxide_name ('str'): 
3078                Name of the hydroxide (OH-) particle.
3079
3080            salt_cation_name ('str'): 
3081                Name of the salt cation (e.g. Na+) particle.
3082
3083            salt_anion_name ('str'): 
3084                Name of the salt anion (e.g. Cl-) particle.
3085
3086            activity_coefficient ('callable'): 
3087                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3088
3089            exclusion_range('pint.Quantity', optional): 
3090                For distances shorter than this value, no particles will be inserted.
3091
3092            use_exclusion_radius_per_type('bool', optional): 
3093                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
3094
3095        Returns:
3096            'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':
3097
3098                'reaction_methods.ReactionEnsemble':  
3099                    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.
3100                
3101                'pint.Quantity': 
3102                    Ionic strength of the reservoir (useful for calculating partition coefficients).
3103
3104        Notess:
3105            - This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
3106
3107        [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3108        """
3109        from espressomd import reaction_methods
3110        if exclusion_range is None:
3111            exclusion_range = max(self.get_radius_map().values())*2.0
3112        if use_exclusion_radius_per_type:
3113            exclusion_radius_per_type = self.get_radius_map()
3114        else:
3115            exclusion_radius_per_type = {}
3116        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3117                                               exclusion_range=exclusion_range, 
3118                                               seed=self.seed, 
3119                                               exclusion_radius_per_type = exclusion_radius_per_type)
3120        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3121        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3122        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3123        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3124        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3125        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3126        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3127        cation_tpl = self.db.get_template(pmb_type="particle",
3128                                          name=salt_cation_name)
3129        cation_state = self.db.get_template(pmb_type="particle_state",
3130                                            name=cation_tpl.initial_state)
3131        anion_tpl = self.db.get_template(pmb_type="particle",
3132                                          name=salt_anion_name)
3133        anion_state = self.db.get_template(pmb_type="particle_state",
3134                                            name=anion_tpl.initial_state)
3135        proton_tpl = self.db.get_template(pmb_type="particle",
3136                                          name=proton_name)
3137        proton_state = self.db.get_template(pmb_type="particle_state",
3138                                            name=proton_tpl.initial_state)
3139        hydroxide_tpl = self.db.get_template(pmb_type="particle",
3140                                             name=hydroxide_name)
3141        hydroxide_state = self.db.get_template(pmb_type="particle_state",
3142                                               name=hydroxide_tpl.initial_state)
3143        proton_es_type = proton_state.es_type
3144        hydroxide_es_type = hydroxide_state.es_type
3145        salt_cation_es_type = cation_state.es_type
3146        salt_anion_es_type = anion_state.es_type
3147        proton_charge = proton_state.z
3148        hydroxide_charge = hydroxide_state.z          
3149        salt_cation_charge = cation_state.z
3150        salt_anion_charge = anion_state.z      
3151        if proton_charge <= 0:
3152            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
3153        if salt_cation_charge <= 0:
3154            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3155        if hydroxide_charge >= 0:
3156            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
3157        if salt_anion_charge >= 0:
3158            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3159        # Grand-canonical coupling to the reservoir
3160        # 0 = H+ + OH-
3161        RE.add_reaction(gamma = K_W.magnitude,
3162                        reactant_types = [],
3163                        reactant_coefficients = [],
3164                        product_types = [ proton_es_type, hydroxide_es_type ],
3165                        product_coefficients = [ 1, 1 ],
3166                        default_charges = {proton_es_type: proton_charge, 
3167                                           hydroxide_es_type: hydroxide_charge})
3168        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3169                                                            state_name=proton_state.name,
3170                                                            coefficient=1),
3171                                        ReactionParticipant(particle_name=hydroxide_name,
3172                                                            state_name=hydroxide_state.name,
3173                                                            coefficient=1)],
3174                           pK=-np.log10(K_W.magnitude),
3175                           reaction_type="ion_insertion",
3176                           simulation_method="GRxMC")
3177        self.db._register_reaction(rx_tpl)
3178        # 0 = Na+ + Cl-
3179        RE.add_reaction(gamma = K_NACL.magnitude,
3180                        reactant_types = [],
3181                        reactant_coefficients = [],
3182                        product_types = [ salt_cation_es_type, salt_anion_es_type ],
3183                        product_coefficients = [ 1, 1 ],
3184                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3185                                        salt_anion_es_type: salt_anion_charge})
3186        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3187                                                            state_name=cation_state.name,
3188                                                            coefficient=1),
3189                                        ReactionParticipant(particle_name=salt_anion_name,
3190                                                            state_name=anion_state.name,
3191                                                            coefficient=1)],
3192                           pK=-np.log10(K_NACL.magnitude),
3193                           reaction_type="ion_insertion",
3194                           simulation_method="GRxMC")
3195        self.db._register_reaction(rx_tpl)
3196        # 0 = Na+ + OH-
3197        RE.add_reaction(gamma = (K_NACL * K_W / K_HCL).magnitude,
3198                        reactant_types = [],
3199                        reactant_coefficients = [],
3200                        product_types = [ salt_cation_es_type, hydroxide_es_type ],
3201                        product_coefficients = [ 1, 1 ],
3202                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3203                                           hydroxide_es_type: hydroxide_charge})
3204        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3205                                                            state_name=cation_state.name,
3206                                                            coefficient=1),
3207                                        ReactionParticipant(particle_name=hydroxide_name,
3208                                                            state_name=hydroxide_state.name,
3209                                                            coefficient=1)],
3210                           pK=-np.log10((K_NACL * K_W / K_HCL).magnitude),
3211                           reaction_type="ion_insertion",
3212                           simulation_method="GRxMC")
3213        self.db._register_reaction(rx_tpl)
3214        # 0 = H+ + Cl-
3215        RE.add_reaction(gamma = K_HCL.magnitude,
3216                        reactant_types = [],
3217                        reactant_coefficients = [],
3218                        product_types = [ proton_es_type, salt_anion_es_type ],
3219                        product_coefficients = [ 1, 1 ],
3220                        default_charges = {proton_es_type: proton_charge, 
3221                                           salt_anion_es_type: salt_anion_charge})
3222        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3223                                                            state_name=proton_state.name,
3224                                                            coefficient=1),
3225                                        ReactionParticipant(particle_name=salt_anion_name,
3226                                                            state_name=anion_state.name,
3227                                                            coefficient=1)],
3228                           pK=-np.log10(K_HCL.magnitude),
3229                           reaction_type="ion_insertion",
3230                           simulation_method="GRxMC")
3231        self.db._register_reaction(rx_tpl)
3232        # Annealing moves to ensure sufficient sampling
3233        # Cation annealing H+ = Na+
3234        RE.add_reaction(gamma = (K_NACL / K_HCL).magnitude,
3235                        reactant_types = [proton_es_type],
3236                        reactant_coefficients = [ 1 ],
3237                        product_types = [ salt_cation_es_type ],
3238                        product_coefficients = [ 1 ],
3239                        default_charges = {proton_es_type: proton_charge, 
3240                                           salt_cation_es_type: salt_cation_charge})
3241        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3242                                                            state_name=proton_state.name,
3243                                                            coefficient=-1),
3244                                        ReactionParticipant(particle_name=salt_cation_name,
3245                                                            state_name=cation_state.name,
3246                                                            coefficient=1)],
3247                           pK=-np.log10((K_NACL / K_HCL).magnitude),
3248                           reaction_type="particle replacement",
3249                           simulation_method="GRxMC")
3250        self.db._register_reaction(rx_tpl)
3251        # Anion annealing OH- = Cl- 
3252        RE.add_reaction(gamma = (K_HCL / K_W).magnitude,
3253                        reactant_types = [hydroxide_es_type],
3254                        reactant_coefficients = [ 1 ],
3255                        product_types = [ salt_anion_es_type ],
3256                        product_coefficients = [ 1 ],
3257            default_charges = {hydroxide_es_type: hydroxide_charge, 
3258                               salt_anion_es_type: salt_anion_charge})
3259        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=hydroxide_name,
3260                                                            state_name=hydroxide_state.name,
3261                                                            coefficient=-1),
3262                                        ReactionParticipant(particle_name=salt_anion_name,
3263                                                            state_name=anion_state.name,
3264                                                            coefficient=1)],
3265                           pK=-np.log10((K_HCL / K_W).magnitude),
3266                           reaction_type="particle replacement",
3267                           simulation_method="GRxMC")
3268        self.db._register_reaction(rx_tpl)
3269        for reaction in self.db.get_reactions():
3270            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
3271                continue
3272            default_charges = {}
3273            reactant_types  = []
3274            product_types   = []
3275            for participant in reaction.participants:
3276                state_tpl = self.db.get_template(name=participant.state_name,
3277                                                 pmb_type="particle_state")
3278                default_charges[state_tpl.es_type] = state_tpl.z
3279                if participant.coefficient < 0:
3280                    reactant_types.append(state_tpl.es_type)
3281                    reactant_name=state_tpl.particle_name
3282                    reactant_state_name=state_tpl.name
3283                elif participant.coefficient > 0:
3284                    product_types.append(state_tpl.es_type)
3285                    product_name=state_tpl.particle_name
3286                    product_state_name=state_tpl.name
3287
3288            Ka = (10**-reaction.pK * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3289            # Reaction in terms of proton: HA = A + H+
3290            RE.add_reaction(gamma=Ka.magnitude,
3291                            reactant_types=reactant_types,
3292                            reactant_coefficients=[1],
3293                            product_types=product_types+[proton_es_type],
3294                            product_coefficients=[1, 1],
3295                            default_charges= default_charges | {proton_es_type: proton_charge})
3296            reaction.add_participant(particle_name=proton_name,
3297                                     state_name=proton_state.name,
3298                                     coefficient=1)
3299            reaction.add_simulation_method("GRxMC")
3300            # Reaction in terms of salt cation: HA = A + Na+
3301            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
3302                            reactant_types=reactant_types,
3303                            reactant_coefficients=[1],
3304                            product_types=product_types+[salt_cation_es_type],
3305                            product_coefficients=[1, 1],
3306                            default_charges=default_charges | {salt_cation_es_type: salt_cation_charge})
3307            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3308                                                                state_name=reactant_state_name,
3309                                                                coefficient=-1),
3310                                            ReactionParticipant(particle_name=product_name,
3311                                                                state_name=product_state_name,
3312                                                                coefficient=1),
3313                                            ReactionParticipant(particle_name=salt_cation_name,
3314                                                                state_name=cation_state.name,
3315                                                                coefficient=1),],
3316                              pK=-np.log10((Ka * K_NACL / K_HCL).magnitude),
3317                              reaction_type=reaction.reaction_type+"_salt",
3318                              simulation_method="GRxMC")
3319            self.db._register_reaction(rx_tpl)
3320            # Reaction in terms of hydroxide: OH- + HA = A
3321            RE.add_reaction(gamma=(Ka / K_W).magnitude,
3322                            reactant_types=reactant_types+[hydroxide_es_type],
3323                            reactant_coefficients=[1, 1],
3324                            product_types=product_types,
3325                            product_coefficients=[1],
3326                            default_charges=default_charges | {hydroxide_es_type: hydroxide_charge})
3327            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3328                                                                state_name=reactant_state_name,
3329                                                                coefficient=-1),
3330                                            ReactionParticipant(particle_name=product_name,
3331                                                                state_name=product_state_name,
3332                                                                coefficient=1),
3333                                            ReactionParticipant(particle_name=hydroxide_name,
3334                                                                state_name=hydroxide_state.name,
3335                                                                coefficient=-1),],
3336                              pK=-np.log10((Ka / K_W).magnitude),
3337                              reaction_type=reaction.reaction_type+"_conjugate",
3338                              simulation_method="GRxMC")
3339            self.db._register_reaction(rx_tpl)
3340            # Reaction in terms of salt anion: Cl- + HA = A
3341            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
3342                            reactant_types=reactant_types+[salt_anion_es_type],
3343                            reactant_coefficients=[1, 1],
3344                            product_types=product_types,
3345                            product_coefficients=[1],
3346                            default_charges=default_charges | {salt_anion_es_type: salt_anion_charge})
3347            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3348                                                                state_name=reactant_state_name,
3349                                                                coefficient=-1),
3350                                            ReactionParticipant(particle_name=product_name,
3351                                                                state_name=product_state_name,
3352                                                                coefficient=1),
3353                                            ReactionParticipant(particle_name=salt_anion_name,
3354                                                                state_name=anion_state.name,
3355                                                                coefficient=-1),],
3356                              pK=-np.log10((Ka / K_HCL).magnitude),
3357                              reaction_type=reaction.reaction_type+"_salt",
3358                              simulation_method="GRxMC")
3359            self.db._register_reaction(rx_tpl)
3360        return RE, ionic_strength_res
3361
3362    def setup_grxmc_unified(self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
3363        """
3364        Sets up acid/base reactions for acidic/basic 'particles' defined in the pyMBE database, as well as a grand-canonical coupling to a 
3365        reservoir of small ions using a unified formulation for small ions.
3366
3367        Args:
3368            pH_res ('float'): 
3369                pH-value in the reservoir.
3370
3371            c_salt_res ('pint.Quantity'): 
3372                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3373
3374            cation_name ('str'): 
3375                Name of the cationic particle.
3376
3377            anion_name ('str'): 
3378                Name of the anionic particle.
3379
3380            activity_coefficient ('callable'): 
3381                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3382
3383            exclusion_range('pint.Quantity', optional): 
3384                Below this value, no particles will be inserted.
3385            
3386            use_exclusion_radius_per_type('bool', optional): 
3387                Controls if one exclusion_radius per each espresso_type. Defaults to 'False'.
3388
3389        Returns:
3390            'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':
3391
3392                'reaction_methods.ReactionEnsemble':  
3393                    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.
3394                
3395                'pint.Quantity': 
3396                    Ionic strength of the reservoir (useful for calculating partition coefficients).
3397
3398        Notes:
3399            - This implementation uses the formulation of the grand-reaction method by Curk et al. [1], which relies on "unified" ion types X+ = {H+, Na+} and X- = {OH-, Cl-}. 
3400            - A function that implements the original version of the grand-reaction method by Landsgesell et al. [2] is also available under the name 'setup_grxmc_reactions'.
3401
3402        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3403        [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3404        """
3405        from espressomd import reaction_methods
3406        if exclusion_range is None:
3407            exclusion_range = max(self.get_radius_map().values())*2.0
3408        if use_exclusion_radius_per_type:
3409            exclusion_radius_per_type = self.get_radius_map()
3410        else:
3411            exclusion_radius_per_type = {}
3412        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3413                                               exclusion_range=exclusion_range, 
3414                                               seed=self.seed, 
3415                                               exclusion_radius_per_type = exclusion_radius_per_type)
3416        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3417        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3418        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3419        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3420        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3421        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3422        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3423        K_XX = a_cation * a_anion
3424        cation_tpl = self.db.get_template(pmb_type="particle",
3425                                          name=cation_name)
3426        cation_state = self.db.get_template(pmb_type="particle_state",
3427                                            name=cation_tpl.initial_state)
3428        anion_tpl = self.db.get_template(pmb_type="particle",
3429                                          name=anion_name)
3430        anion_state = self.db.get_template(pmb_type="particle_state",
3431                                            name=anion_tpl.initial_state)
3432        cation_es_type = cation_state.es_type
3433        anion_es_type = anion_state.es_type     
3434        cation_charge = cation_state.z
3435        anion_charge = anion_state.z
3436        if cation_charge <= 0:
3437            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3438        if anion_charge >= 0:
3439            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3440        # Coupling to the reservoir: 0 = X+ + X-
3441        RE.add_reaction(gamma = K_XX.magnitude,
3442                        reactant_types = [],
3443                        reactant_coefficients = [],
3444                        product_types = [ cation_es_type, anion_es_type ],
3445                        product_coefficients = [ 1, 1 ],
3446                        default_charges = {cation_es_type: cation_charge, 
3447                                           anion_es_type: anion_charge})
3448        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=cation_name,
3449                                                            state_name=cation_state.name,
3450                                                            coefficient=1),
3451                                        ReactionParticipant(particle_name=anion_name,
3452                                                            state_name=anion_state.name,
3453                                                            coefficient=1)],
3454                           pK=-np.log10(K_XX.magnitude),
3455                           reaction_type="ion_insertion",
3456                           simulation_method="GCMC")
3457        self.db._register_reaction(rx_tpl)
3458        for reaction in self.db.get_reactions():
3459            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
3460                continue
3461            default_charges = {}
3462            reactant_types  = []
3463            product_types   = []
3464            for participant in reaction.participants:
3465                state_tpl = self.db.get_template(name=participant.state_name,
3466                                                 pmb_type="particle_state")
3467                default_charges[state_tpl.es_type] = state_tpl.z
3468                if participant.coefficient < 0:
3469                    reactant_types.append(state_tpl.es_type)
3470                    reactant_name=state_tpl.particle_name
3471                    reactant_state_name=state_tpl.name
3472                elif participant.coefficient > 0:
3473                    product_types.append(state_tpl.es_type)
3474                    product_name=state_tpl.particle_name
3475                    product_state_name=state_tpl.name
3476
3477            Ka = (10**-reaction.pK * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3478            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3479            # Reaction in terms of small cation: HA = A + X+
3480            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3481                            reactant_types=reactant_types,
3482                            reactant_coefficients=[1],
3483                            product_types=product_types+[cation_es_type],
3484                            product_coefficients=[1, 1],
3485                            default_charges=default_charges|{cation_es_type: cation_charge})
3486            reaction.add_participant(particle_name=cation_name,
3487                                     state_name=cation_state.name,
3488                                     coefficient=1)
3489            reaction.add_simulation_method("GRxMC")
3490            # Reaction in terms of small anion: X- + HA = A
3491            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3492                            reactant_types=reactant_types+[anion_es_type],
3493                            reactant_coefficients=[1, 1],
3494                            product_types=product_types,
3495                            product_coefficients=[1],
3496                            default_charges=default_charges|{anion_es_type: anion_charge})
3497            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3498                                                                state_name=reactant_state_name,
3499                                                                coefficient=-1),
3500                                            ReactionParticipant(particle_name=product_name,
3501                                                                state_name=product_state_name,
3502                                                                coefficient=1),
3503                                            ReactionParticipant(particle_name=anion_name,
3504                                                                state_name=anion_state.name,
3505                                                                coefficient=-1),],
3506                              pK=-np.log10(gamma_K_AX.magnitude / K_XX.magnitude),
3507                              reaction_type=reaction.reaction_type+"_conjugate",
3508                              simulation_method="GRxMC")
3509            self.db._register_reaction(rx_tpl)
3510        return RE, ionic_strength_res
3511
3512    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot'):
3513        """
3514        Sets up the Lennard-Jones (LJ) potential between all pairs of particle states defined in the pyMBE database.
3515
3516        Args:
3517            espresso_system('espressomd.system.System'): 
3518                Instance of a system object from the espressomd library.
3519
3520            shift_potential('bool', optional): 
3521                If True, a shift will be automatically computed such that the potential is continuous at the cutoff radius. Otherwise, no shift will be applied. Defaults to True.
3522
3523            combining_rule('string', optional): 
3524                combining rule used to calculate 'sigma' and 'epsilon' for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3525
3526            warning('bool', optional): 
3527                switch to activate/deactivate warning messages. Defaults to True.
3528
3529        Notes:
3530            - Currently, the only 'combining_rule' supported is Lorentz-Berthelot.
3531            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3532
3533        """
3534        from itertools import combinations_with_replacement
3535        particle_templates = self.db.get_templates("particle")
3536        shift = "auto" if shift_potential else 0
3537        if shift == "auto":
3538            shift_tpl = shift
3539        else:
3540            shift_tpl = PintQuantity.from_quantity(q=shift*self.units.reduced_length,
3541                                                   expected_dimension="length",
3542                                                   ureg=self.units)
3543        # Get all particle states registered in pyMBE
3544        state_entries = []
3545        for tpl in particle_templates.values():
3546            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
3547                state_entries.append((tpl, state))
3548
3549        # Iterate over all unique state pairs
3550        for (tpl1, state1), (tpl2, state2) in combinations_with_replacement(state_entries, 2):
3551
3552            lj_parameters = self.get_lj_parameters(particle_name1=tpl1.name,
3553                                                   particle_name2=tpl2.name,
3554                                                   combining_rule=combining_rule)
3555            if not lj_parameters:
3556                continue
3557
3558            espresso_system.non_bonded_inter[state1.es_type, state2.es_type].lennard_jones.set_params(
3559                epsilon=lj_parameters["epsilon"].to("reduced_energy").magnitude,
3560                sigma=lj_parameters["sigma"].to("reduced_length").magnitude,
3561                cutoff=lj_parameters["cutoff"].to("reduced_length").magnitude,
3562                offset=lj_parameters["offset"].to("reduced_length").magnitude,
3563                shift=shift)
3564                
3565            lj_template = LJInteractionTemplate(state1=state1.name,
3566                                                state2=state2.name,
3567                                                sigma=PintQuantity.from_quantity(q=lj_parameters["sigma"],
3568                                                                                 expected_dimension="length",
3569                                                                                 ureg=self.units),
3570                                                epsilon=PintQuantity.from_quantity(q=lj_parameters["epsilon"],
3571                                                                                   expected_dimension="energy",
3572                                                                                   ureg=self.units),
3573                                                cutoff=PintQuantity.from_quantity(q=lj_parameters["cutoff"],
3574                                                                                  expected_dimension="length",
3575                                                                                  ureg=self.units),
3576                                                offset=PintQuantity.from_quantity(q=lj_parameters["offset"],
3577                                                                                  expected_dimension="length",
3578                                                                                  ureg=self.units),
3579                                                shift=shift_tpl)
3580            self.db._register_template(lj_template)
class pymbe_library:
  58class pymbe_library():
  59    """
  60    Core library of the Molecular Builder for ESPResSo (pyMBE).
  61
  62    Attributes:
  63        N_A ('pint.Quantity'):
  64            Avogadro number.
  65
  66        kB ('pint.Quantity'):
  67            Boltzmann constant.
  68
  69        e ('pint.Quantity'):
  70            Elementary charge.
  71
  72        kT ('pint.Quantity'):
  73            Thermal energy corresponding to the set temperature.
  74
  75        Kw ('pint.Quantity'):
  76            Ionic product of water, used in G-RxMC and Donnan-related calculations.
  77
  78        db ('Manager'):
  79            Database manager holding all pyMBE templates, instances and reactions.
  80
  81        rng ('numpy.random.Generator'):
  82            Random number generator initialized with the provided seed.
  83
  84        units ('pint.UnitRegistry'):
  85            Pint unit registry used for unit-aware calculations.
  86
  87        lattice_builder ('pyMBE.lib.lattice.LatticeBuilder'):
  88            Optional lattice builder object (initialized as ''None'').
  89            
  90        root ('importlib.resources.abc.Traversable'):
  91            Root path to the pyMBE package resources.
  92    """
  93
  94    def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None):
  95        """
  96        Initializes the pyMBE library.
  97
  98        Args:
  99            seed ('int'):
 100                Seed for the random number generator.
 101
 102            temperature ('pint.Quantity', optional):
 103                Simulation temperature. If ''None'', defaults to 298.15 K.
 104
 105            unit_length ('pint.Quantity', optional):
 106                Reference length for reduced units. If ''None'', defaults to
 107                0.355 nm.
 108
 109            unit_charge ('pint.Quantity', optional):
 110                Reference charge for reduced units. If ''None'', defaults to
 111                one elementary charge.
 112
 113            Kw ('pint.Quantity', optional):
 114                Ionic product of water (typically in mol²/L²). If ''None'',
 115                defaults to 1e-14 mol²/L².
 116        """
 117        # Seed and RNG
 118        self.seed=seed
 119        self.rng = np.random.default_rng(seed)
 120        self.units=pint.UnitRegistry()
 121        self.N_A=scipy.constants.N_A / self.units.mol
 122        self.kB=scipy.constants.k * self.units.J / self.units.K
 123        self.e=scipy.constants.e * self.units.C
 124        self.set_reduced_units(unit_length=unit_length, 
 125                               unit_charge=unit_charge,
 126                               temperature=temperature, 
 127                               Kw=Kw)
 128        
 129        self.db = Manager(units=self.units)
 130        self.lattice_builder = None
 131        self.root = importlib.resources.files(__package__)
 132
 133    def _check_bond_inputs(self, bond_type, bond_parameters):
 134        """
 135        Checks that the input bond parameters are valid within the current pyMBE implementation.
 136
 137        Args:
 138            bond_type ('str'): 
 139                label to identify the potential to model the bond.
 140            
 141            bond_parameters ('dict'): 
 142                parameters of the potential of the bond.
 143        """
 144        valid_bond_types   = ["harmonic", "FENE"] 
 145        if bond_type not in valid_bond_types:
 146            raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}")
 147        required_parameters = {"harmonic": ["r_0","k"],
 148                                "FENE": ["r_0","k","d_r_max"]}
 149        for required_parameter in required_parameters[bond_type]:
 150            if required_parameter not in bond_parameters.keys():
 151                raise ValueError(f"Missing required parameter {required_parameter} for {bond_type} bond")
 152            
 153    def _check_dimensionality(self, variable, expected_dimensionality):
 154        """
 155        Checks if the dimensionality of 'variable' matches 'expected_dimensionality'.
 156
 157        Args:
 158            variable ('pint.Quantity'): 
 159                Quantity to be checked.
 160
 161            expected_dimensionality ('str'): 
 162                Expected dimension of the variable.
 163
 164        Returns:
 165            ('bool'): 
 166                'True' if the variable if of the expected dimensionality, 'False' otherwise.
 167
 168        Notes:
 169            - 'expected_dimensionality' takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality).
 170            - For example, to check for a variable corresponding to a velocity 'expected_dimensionality = "[length]/[time]"'    
 171        """
 172        correct_dimensionality=variable.check(f"{expected_dimensionality}")      
 173        if not correct_dimensionality:
 174            raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}")
 175        return correct_dimensionality   
 176            
 177    def _check_pka_set(self, pka_set):
 178        """
 179        Checks that 'pka_set' has the formatting expected by pyMBE.
 180       
 181        Args:
 182            pka_set ('dict'): 
 183                {"name" : {"pka_value": pka, "acidity": acidity}}
 184        """
 185        required_keys=['pka_value','acidity']
 186        for required_key in required_keys:
 187            for pka_name, pka_entry in pka_set.items():
 188                if required_key not in pka_entry:
 189                    raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")')
 190        return
 191
 192    def _create_espresso_bond_instance(self, bond_type, bond_parameters):
 193        """
 194        Creates an ESPResSo bond instance.
 195
 196        Args:
 197            bond_type ('str'): 
 198                label to identify the potential to model the bond.
 199
 200            bond_parameters ('dict'): 
 201                parameters of the potential of the bond.
 202
 203        Notes:
 204            Currently, only HARMONIC and FENE bonds are supported.
 205
 206            For a HARMONIC bond the dictionary must contain:
 207                - k ('Pint.Quantity')      : Magnitude of the bond. It should have units of energy/length**2 
 208                using the 'pmb.units' UnitRegistry.
 209                - r_0 ('Pint.Quantity')    : Equilibrium bond length. It should have units of length using 
 210                the 'pmb.units' UnitRegistry.
 211           
 212            For a FENE bond the dictionary must additionally contain:
 213                - d_r_max ('Pint.Quantity'): Maximal stretching length for FENE. It should have 
 214                units of length using the 'pmb.units' UnitRegistry. Default 'None'.
 215
 216        Returns:
 217            ('espressomd.interactions'): instance of an ESPResSo bond object
 218        """
 219        from espressomd import interactions
 220        self._check_bond_inputs(bond_parameters=bond_parameters,
 221                                bond_type=bond_type)
 222        if bond_type == 'harmonic':
 223            bond_instance = interactions.HarmonicBond(k = bond_parameters["k"].m_as("reduced_energy/reduced_length**2"),
 224                                                      r_0 = bond_parameters["r_0"].m_as("reduced_length"))
 225        elif bond_type == 'FENE':
 226            bond_instance    = interactions.FeneBond(k = bond_parameters["k"].m_as("reduced_energy/reduced_length**2"),
 227                                                      r_0 = bond_parameters["r_0"].m_as("reduced_length"),
 228                                                      d_r_max = bond_parameters["d_r_max"].m_as("reduced_length"))    
 229        return bond_instance
 230
 231    def _create_hydrogel_chain(self, hydrogel_chain, nodes, espresso_system, use_default_bond=False, gen_angle=False):
 232        """
 233        Creates a chain between two nodes of a hydrogel.
 234
 235        Args:
 236            hydrogel_chain ('HydrogelChain'): 
 237                template of a hydrogel chain
 238            nodes ('dict'): 
 239                {node_index: {"name": node_particle_name, "pos": node_position, "id": node_particle_instance_id}}
 240
 241            espresso_system ('espressomd.system.System'): 
 242                ESPResSo system object where the hydrogel chain will be created.
 243
 244            use_default_bond ('bool', optional): 
 245                If True, use a default bond template if no specific template exists. Defaults to False.
 246
 247            gen_angle ('bool', optional):
 248                If True, generate the angle potentials internal to the created
 249                chain molecule. Junction angles near the hydrogel crosslinkers
 250                are handled separately at the hydrogel level.
 251
 252        Return:
 253            ('int'): 
 254                molecule_id of the created hydrogel chian.
 255
 256        Notes:
 257            - If the chain is defined between node_start = ''[0 0 0]'' and node_end = ''[1 1 1]'', the chain will be placed between these two nodes.
 258            - The chain will be placed in the direction of the vector between 'node_start' and 'node_end'. 
 259        """
 260        if self.lattice_builder is None:
 261            raise ValueError("LatticeBuilder is not initialized. Use 'initialize_lattice_builder' first.")
 262        molecule_tpl = self.db.get_template(pmb_type="molecule",
 263                                            name=hydrogel_chain.molecule_name)
 264        residue_list = molecule_tpl.residue_list
 265        molecule_name = molecule_tpl.name
 266        node_start = hydrogel_chain.node_start
 267        node_end = hydrogel_chain.node_end
 268        node_start_label = self.lattice_builder._create_node_label(node_start)
 269        node_end_label = self.lattice_builder._create_node_label(node_end)
 270        _, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end)
 271        if node_start == node_end and residue_list != residue_list[::-1]:
 272            raise ValueError(f"Aborted creation of hydrogel chain between '{node_start}' and '{node_end}' because pyMBE could not resolve a unique topology for that chain")
 273        if reverse:
 274            reverse_residue_order=True
 275        else:
 276            reverse_residue_order=False
 277        start_node_id = nodes[node_start_label]["id"]
 278        end_node_id = nodes[node_end_label]["id"]
 279        # Finding a backbone vector between node_start and node_end
 280        vec_between_nodes = np.array(nodes[node_end_label]["pos"]) - np.array(nodes[node_start_label]["pos"])
 281        vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l)
 282        backbone_vector = vec_between_nodes / np.linalg.norm(vec_between_nodes)
 283        if reverse_residue_order:
 284            vec_between_nodes *= -1.0
 285        # Calculate the start position of the chain
 286        chain_residues = self.db.get_template(pmb_type="molecule",
 287                                              name=molecule_name).residue_list
 288        part_start_chain_name = self.db.get_template(pmb_type="residue",
 289                                                     name=chain_residues[0]).central_bead
 290        lj_parameters = self.get_lj_parameters(particle_name1=nodes[node_start_label]["name"],
 291                                               particle_name2=part_start_chain_name)
 292        bond_tpl = self.get_bond_template(particle_name1=nodes[node_start_label]["name"],
 293                                          particle_name2=part_start_chain_name,
 294                                          use_default_bond=use_default_bond)
 295        l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
 296                                              bond_type=bond_tpl.bond_type,
 297                                              bond_parameters=bond_tpl.get_parameters(ureg=self.units))
 298        first_bead_pos = np.array((nodes[node_start_label]["pos"])) + np.array(backbone_vector)*l0
 299        mol_id = self.create_molecule(name=molecule_name,  # Use the name defined earlier
 300                                      number_of_molecules=1,  # Creating one chain
 301                                      espresso_system=espresso_system,
 302                                      list_of_first_residue_positions=[first_bead_pos.tolist()], #Start at the first node
 303                                      backbone_vector=np.array(backbone_vector)/l0,
 304                                      use_default_bond=use_default_bond,
 305                                      reverse_residue_order=reverse_residue_order,
 306                                      gen_angle=gen_angle)[0]
 307        # Bond chain to the hydrogel nodes
 308        chain_pids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
 309                                                             attribute="molecule_id",
 310                                                             value=mol_id)
 311        self.create_bond(particle_id1=start_node_id,
 312                         particle_id2=chain_pids[0],
 313                         espresso_system=espresso_system,
 314                         use_default_bond=use_default_bond)
 315        self.create_bond(particle_id1=chain_pids[-1],
 316                         particle_id2=end_node_id,
 317                         espresso_system=espresso_system,
 318                         use_default_bond=use_default_bond)
 319        return mol_id
 320
 321    def _generate_hydrogel_crosslinker_angles(self, espresso_system, central_particle_ids):
 322        """
 323        Generate hydrogel angles centered on crosslinkers and adjacent terminal beads.
 324
 325        If the user defines any explicit angle template for such junction
 326        triplets, then all required junction triplets must be defined. If none
 327        are defined, hydrogel construction proceeds without crosslinker-adjacent
 328        angles.
 329        """
 330        particle_instances = self.db.get_instances(pmb_type="particle")
 331        bonded_neighbors = {}
 332        for bond in self.db.get_instances(pmb_type="bond").values():
 333            bonded_neighbors.setdefault(bond.particle_id1, set()).add(bond.particle_id2)
 334            bonded_neighbors.setdefault(bond.particle_id2, set()).add(bond.particle_id1)
 335
 336        triplets = []
 337        for central_particle_id in sorted(set(central_particle_ids)):
 338            neighbors = sorted(bonded_neighbors.get(central_particle_id, set()))
 339            central_name = particle_instances[central_particle_id].name
 340            for idx_i in range(len(neighbors)):
 341                for idx_k in range(idx_i + 1, len(neighbors)):
 342                    side_particle_id1 = neighbors[idx_i]
 343                    side_particle_id3 = neighbors[idx_k]
 344                    side_name1 = particle_instances[side_particle_id1].name
 345                    side_name3 = particle_instances[side_particle_id3].name
 346                    angle_key = AngleTemplate.make_angle_key(side1=side_name1,
 347                                                             central=central_name,
 348                                                             side2=side_name3)
 349                    triplets.append((side_particle_id1,
 350                                     central_particle_id,
 351                                     side_particle_id3,
 352                                     angle_key))
 353
 354        defined_angle_templates = self.db.get_templates(pmb_type="angle")
 355        defined_angle_keys = {
 356            angle_key
 357            for _, _, _, angle_key in triplets
 358            if angle_key in defined_angle_templates
 359        }
 360        if not defined_angle_keys:
 361            logging.warning("No angle templates defined for hydrogel crosslinkers")
 362            return
 363
 364        missing_angle_keys = sorted({
 365            angle_key
 366            for _, _, _, angle_key in triplets
 367            if angle_key not in defined_angle_keys
 368        })
 369        if missing_angle_keys:
 370            raise ValueError(
 371                "Hydrogel crosslinker-adjacent angle templates must be defined for all required triplets. "
 372                f"Missing definitions for: {missing_angle_keys}"
 373            )
 374
 375        for side_particle_id1, central_particle_id, side_particle_id3, _ in triplets:
 376            self.create_angular_potential(particle_id1=side_particle_id1,
 377                              particle_id2=central_particle_id,
 378                              particle_id3=side_particle_id3,
 379                              espresso_system=espresso_system,
 380                              use_default_angle=False)
 381
 382    def _create_hydrogel_node(self, node_index, node_name, espresso_system):
 383        """
 384        Set a node residue type.
 385        
 386        Args:
 387            node_index ('str'): 
 388                Lattice node index in the form of a string, e.g. "[0 0 0]".
 389
 390            node_name ('str'): 
 391                name of the node particle defined in pyMBE.
 392
 393            espresso_system (espressomd.system.System): 
 394                ESPResSo system object where the hydrogel node will be created.
 395
 396        Returns:
 397            ('tuple(list,int)'):
 398                ('list'): Position of the node in the lattice.
 399                ('int'): Particle ID of the node.
 400        """
 401        if self.lattice_builder is None:
 402            raise ValueError("LatticeBuilder is not initialized. Use 'initialize_lattice_builder' first.")
 403        node_position = np.array(node_index)*0.25*self.lattice_builder.box_l
 404        p_id = self.create_particle(name = node_name,
 405                                    espresso_system=espresso_system,
 406                                    number_of_particles=1,
 407                                    position = [node_position])
 408        key = self.lattice_builder._get_node_by_label(f"[{node_index[0]} {node_index[1]} {node_index[2]}]")
 409        self.lattice_builder.nodes[key] = node_name
 410        return node_position.tolist(), p_id[0]
 411
 412    def _get_espresso_bond_instance(self, bond_template, espresso_system):
 413        """
 414        Retrieve or create a bond instance in an ESPResSo system for a given pair of particle names.
 415
 416        Args:
 417            bond_template ('BondTemplate'): 
 418                BondTemplate object from the pyMBE database.
 419            espresso_system ('espressomd.system.System'): 
 420                An ESPResSo system object where the bond will be added or retrieved.
 421
 422        Returns:
 423            ('espressomd.interactions.BondedInteraction'): 
 424                The ESPResSo bond instance object.
 425
 426        Notes:
 427            When a new bond instance is created, it is not added to the ESPResSo system.
 428        """
 429        if bond_template.name in self.db.espresso_bond_instances.keys():
 430            bond_inst = self.db.espresso_bond_instances[bond_template.name]
 431        else:   
 432            # Create an instance of the bond 
 433            bond_inst = self._create_espresso_bond_instance(bond_type=bond_template.bond_type,
 434                                                            bond_parameters=bond_template.get_parameters(self.units))
 435            self.db.espresso_bond_instances[bond_template.name]= bond_inst
 436            espresso_system.bonded_inter.add(bond_inst)
 437        return bond_inst
 438
 439    def _get_label_id_map(self, pmb_type):
 440        """
 441        Returns the key used to access the particle ID map for a given pyMBE object type.
 442
 443        Args:
 444            pmb_type ('str'):
 445                pyMBE object type for which the particle ID map label is requested.
 446
 447        Returns:
 448            'str':
 449                Label identifying the appropriate particle ID map. 
 450        """
 451        if pmb_type in self.db._assembly_like_types:
 452            label="assembly_map"
 453        elif pmb_type in self.db._molecule_like_types:
 454            label="molecule_map"
 455        else:
 456            label=f"{pmb_type}_map"
 457        return label
 458
 459    def _get_residue_list_from_sequence(self, sequence):
 460        """
 461        Convenience function to get a 'residue_list' from a protein or peptide 'sequence'.
 462
 463        Args:
 464            sequence ('lst'): 
 465                 Sequence of the peptide or protein.
 466
 467        Returns:
 468            residue_list ('list' of 'str'): 
 469                List of the 'name's of the 'residue's  in the sequence of the 'molecule'.             
 470        """
 471        residue_list = []
 472        for item in sequence:
 473            residue_name='AA-'+item
 474            residue_list.append(residue_name)
 475        return residue_list
 476    
 477    def _get_template_type(self, name, allowed_types):
 478        """
 479        Validate that a template name resolves unambiguously to exactly one
 480        allowed pmb_type in the pyMBE database and return it.
 481
 482        Args:
 483            name ('str'): 
 484                Name of the template to validate.
 485
 486            allowed_types ('set[str]'):  
 487                Set of allowed pmb_type values (e.g. {"molecule", "peptide"}).
 488
 489        Returns:
 490            ('str'): 
 491                Resolved pmb_type.
 492
 493        Notess:
 494            - This method does *not* return the template itself, only the validated pmb_type. 
 495        """
 496        registered_pmb_types_with_name = self.db._find_template_types(name=name)
 497        filtered_types = allowed_types.intersection(registered_pmb_types_with_name)
 498        if len(filtered_types) > 1:
 499            raise ValueError(f"Ambiguous template name '{name}': found {len(filtered_types)} templates in the pyMBE database. Molecule creation aborted.")  
 500        if len(filtered_types) == 0:
 501            raise ValueError(f"No {allowed_types} template found with name '{name}'. Found templates of types: {filtered_types}.")
 502        return next(iter(filtered_types))
 503
 504    def _delete_particles_from_espresso(self, particle_ids, espresso_system):
 505        """
 506        Remove a list of particles from an ESPResSo simulation system.
 507
 508        Args:
 509            particle_ids  ('Iterable[int]'):
 510                A list (or other iterable) of ESPResSo particle IDs to remove.
 511
 512            espresso_system ('espressomd.system.System'):
 513                The ESPResSo simulation system from which the particles
 514                will be removed.
 515
 516        Notess:
 517            - This method removes particles only from the ESPResSo simulation,
 518            **not** from the pyMBE database. Database cleanup must be handled
 519            separately by the caller.
 520            - Attempting to remove a non-existent particle ID will raise
 521            an ESPResSo error.
 522        """
 523        for pid in particle_ids:
 524            espresso_system.part.by_id(pid).remove()
 525
 526    def calculate_center_of_mass(self, instance_id, pmb_type, espresso_system):
 527        """
 528        Calculates the center of mass of a pyMBE object instance in an ESPResSo system.
 529
 530        Args:
 531            instance_id ('int'):
 532                pyMBE instance ID of the object whose center of mass is calculated.
 533
 534            pmb_type ('str'):
 535                Type of the pyMBE object. Must correspond to a particle-aggregating
 536                template type (e.g. '"molecule"', '"residue"', '"peptide"', '"protein"').
 537
 538            espresso_system ('espressomd.system.System'):
 539                ESPResSo system containing the particle instances.
 540
 541        Returns:
 542            ('numpy.ndarray'):
 543                Array of shape '(3,)' containing the Cartesian coordinates of the
 544                center of mass.
 545
 546        Notes:
 547            - This method assumes equal mass for all particles.
 548            - Periodic boundary conditions are *not* unfolded; positions are taken
 549            directly from ESPResSo particle coordinates.
 550        """
 551        center_of_mass = np.zeros(3)
 552        axis_list = [0,1,2]
 553        inst = self.db.get_instance(pmb_type=pmb_type,
 554                                    instance_id=instance_id)
 555        particle_id_list = self.get_particle_id_map(object_name=inst.name)["all"]
 556        for pid in particle_id_list:
 557            for axis in axis_list:
 558                center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis]
 559        center_of_mass = center_of_mass / len(particle_id_list)
 560        return center_of_mass
 561
 562    def calculate_HH(self, template_name, pH_list=None, pka_set=None):
 563        """
 564        Calculates the charge in the template object according to the ideal  Henderson–Hasselbalch titration curve.
 565
 566        Args:
 567            template_name ('str'):
 568                Name of the template.
 569
 570            pH_list ('list[float]', optional):
 571                pH values at which the charge is evaluated.
 572                Defaults to 50 values between 2 and 12.
 573
 574            pka_set ('dict', optional):
 575                Mapping: {particle_name: {"pka_value": 'float', "acidity": "acidic"|"basic"}}
 576
 577        Returns:
 578            'list[float]':
 579                Net molecular charge at each pH value.
 580        """
 581        if pH_list is None:
 582            pH_list = np.linspace(2, 12, 50)
 583        if pka_set is None:
 584            pka_set = self.get_pka_set()
 585        self._check_pka_set(pka_set=pka_set)
 586        particle_counts = self.db.get_particle_templates_under(template_name=template_name,
 587                                                               return_counts=True)
 588        if not particle_counts:
 589            return [None] * len(pH_list)
 590        charge_number_map = self.get_charge_number_map()
 591        def formal_charge(particle_name):
 592            tpl = self.db.get_template(name=particle_name, 
 593                                       pmb_type="particle")
 594            state = self.db.get_template(name=tpl.initial_state,
 595                                         pmb_type="particle_state")
 596            return charge_number_map[state.es_type]
 597        Z_HH = []
 598        for pH in pH_list:
 599            Z = 0.0
 600            for particle, multiplicity in particle_counts.items():
 601                if particle in pka_set:
 602                    pka = pka_set[particle]["pka_value"]
 603                    acidity = pka_set[particle]["acidity"]
 604                    if acidity == "acidic":
 605                        psi = -1
 606                    elif acidity == "basic":
 607                        psi = +1
 608                    else:
 609                        raise ValueError(f"Unknown acidity '{acidity}' for particle '{particle}'")
 610                    charge = psi / (1.0 + 10.0 ** (psi * (pH - pka)))
 611                    Z += multiplicity * charge
 612                else:
 613                    Z += multiplicity * formal_charge(particle)
 614            Z_HH.append(Z)
 615        return Z_HH   
 616
 617    def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
 618        """
 619        Computes macromolecular charges using the Henderson–Hasselbalch equation
 620        coupled to ideal Donnan partitioning.
 621
 622        Args:
 623            c_macro ('dict'):
 624                Mapping of macromolecular species names to their concentrations
 625                in the system:
 626                '{molecule_name: concentration}'.
 627                
 628            c_salt ('float' or 'pint.Quantity'):
 629                Salt concentration in the reservoir.
 630
 631            pH_list ('list[float]', optional):
 632                List of pH values in the reservoir at which the calculation is
 633                performed. If 'None', 50 equally spaced values between 2 and 12
 634                are used.
 635
 636            pka_set ('dict', optional):
 637                Dictionary defining the acid–base properties of titratable particle
 638                types:
 639                '{particle_name: {"pka_value": float, "acidity": "acidic" | "basic"}}'.
 640                If 'None', the pKa set is taken from the pyMBE database.
 641
 642        Returns:
 643            'dict':
 644                Dictionary containing:
 645                - '"charges_dict"' ('dict'):
 646                    Mapping '{molecule_name: list}' of Henderson–Hasselbalch–Donnan
 647                    charges evaluated at each pH value.
 648                - '"pH_system_list"' ('list[float]'):
 649                    Effective pH values inside the system phase after Donnan
 650                    partitioning.
 651                - '"partition_coefficients"' ('list[float]'):
 652                    Partition coefficients of monovalent cations at each pH value.
 653
 654        Notes:
 655            - This method assumes **ideal Donnan equilibrium** and **monovalent salt**.
 656            - The ionic strength of the reservoir includes both salt and
 657            pH-dependent H⁺/OH⁻ contributions.
 658            - All charged macromolecular species present in the system must be
 659            included in 'c_macro'; missing species will lead to incorrect results.
 660            - The nonlinear Donnan equilibrium equation is solved using a scalar
 661            root finder ('brentq') in logarithmic form for numerical stability.
 662            - This method is intended for **two-phase systems**; for single-phase
 663            systems use 'calculate_HH' instead.
 664        """
 665        if pH_list is None:
 666            pH_list=np.linspace(2,12,50)
 667        if pka_set is None:
 668            pka_set=self.get_pka_set() 
 669        self._check_pka_set(pka_set=pka_set)
 670        partition_coefficients_list = []
 671        pH_system_list = []
 672        Z_HH_Donnan={}
 673        for key in c_macro:
 674            Z_HH_Donnan[key] = []
 675        def calc_charges(c_macro, pH):
 676            """
 677            Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation.
 678
 679            Args:
 680                c_macro ('dict'): 
 681                    {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 682
 683                pH ('float'): 
 684                    pH-value that is used in the HH equation.
 685
 686            Returns:
 687                ('dict'): 
 688                    {"molecule_name": charge}
 689            """
 690            charge = {}
 691            for name in c_macro:
 692                charge[name] = self.calculate_HH(name, [pH], pka_set)[0]
 693            return charge
 694
 695        def calc_partition_coefficient(charge, c_macro):
 696            """
 697            Calculates the partition coefficients of positive ions according to the ideal Donnan theory.
 698
 699            Args:
 700                charge ('dict'): 
 701                    {"molecule_name": charge}
 702
 703                c_macro ('dict'): 
 704                    {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 705            """
 706            nonlocal ionic_strength_res
 707            charge_density = 0.0
 708            for key in charge:
 709                charge_density += charge[key] * c_macro[key]
 710            return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude
 711        for pH_value in pH_list:    
 712            # calculate the ionic strength of the reservoir
 713            if pH_value <= 7.0:
 714                ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 
 715            elif pH_value > 7.0:
 716                ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt
 717            #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations
 718            #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation.
 719            #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects.
 720            equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro))
 721            logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq")
 722            partition_coefficient = 10**logxi.root
 723            charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient))
 724            for key in c_macro:
 725                Z_HH_Donnan[key].append(charges_temp[key])
 726            pH_system_list.append(pH_value - np.log10(partition_coefficient))
 727            partition_coefficients_list.append(partition_coefficient)
 728        return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
 729
 730    def calculate_net_charge(self,espresso_system,object_name,pmb_type,dimensionless=False):
 731        """
 732        Calculates the net charge per instance of a given pmb object type.
 733
 734        Args:
 735            espresso_system (espressomd.system.System):
 736                ESPResSo system containing the particles.
 737            object_name (str):
 738                Name of the object (e.g. molecule, residue, peptide, protein).
 739            pmb_type (str):
 740                Type of object to analyze. Must be molecule-like.
 741            dimensionless (bool, optional):
 742                If True, return charge as a pure number.
 743                If False, return a quantity with reduced_charge units.
 744
 745        Returns:
 746            dict:
 747                {"mean": mean_net_charge, "instances": {instance_id: net_charge}}
 748        """
 749        id_map = self.get_particle_id_map(object_name=object_name)
 750        label = self._get_label_id_map(pmb_type=pmb_type)
 751        instance_map = id_map[label]
 752        charges = {}
 753        for instance_id, particle_ids in instance_map.items():
 754            if dimensionless:
 755                net_charge = 0.0
 756            else:
 757                net_charge = 0 * self.units.Quantity(1, "reduced_charge")
 758            for pid in particle_ids:
 759                q = espresso_system.part.by_id(pid).q
 760                if not dimensionless:
 761                    q *= self.units.Quantity(1, "reduced_charge")
 762                net_charge += q
 763            charges[instance_id] = net_charge
 764        # Mean charge
 765        if dimensionless:
 766            mean_charge = float(np.mean(list(charges.values())))
 767        else:
 768            mean_charge = (np.mean([q.magnitude for q in charges.values()])* self.units.Quantity(1, "reduced_charge"))
 769        return {"mean": mean_charge, "instances": charges}
 770
 771    def center_object_in_simulation_box(self, instance_id, espresso_system, pmb_type):
 772        """
 773        Centers a pyMBE object instance in the simulation box of an ESPResSo system.
 774        The object is translated such that its center of mass coincides with the
 775        geometric center of the ESPResSo simulation box.
 776
 777        Args:
 778            instance_id ('int'):
 779                ID of the pyMBE object instance to be centered.
 780
 781            pmb_type ('str'):
 782                Type of the pyMBE object.
 783
 784            espresso_system ('espressomd.system.System'):
 785                ESPResSo system object in which the particles are defined.
 786
 787        Notes:
 788            - Works for both cubic and non-cubic simulation boxes.
 789        """
 790        inst = self.db.get_instance(instance_id=instance_id,
 791                                    pmb_type=pmb_type)
 792        center_of_mass = self.calculate_center_of_mass(instance_id=instance_id,
 793                                                       espresso_system=espresso_system,
 794                                                       pmb_type=pmb_type)
 795        box_center = [espresso_system.box_l[0]/2.0,
 796                      espresso_system.box_l[1]/2.0,
 797                      espresso_system.box_l[2]/2.0]
 798        particle_id_list = self.get_particle_id_map(object_name=inst.name)["all"]
 799        for pid in particle_id_list:
 800            es_pos = espresso_system.part.by_id(pid).pos
 801            espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center
 802
 803    def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt):    
 804        """
 805        Creates a 'c_salt' concentration of 'cation_name' and 'anion_name' ions into the 'espresso_system'.
 806
 807        Args:
 808            espresso_system('espressomd.system.System'): instance of an espresso system object.
 809            cation_name('str'): 'name' of a particle with a positive charge.
 810            anion_name('str'): 'name' of a particle with a negative charge.
 811            c_salt('float'): Salt concentration.
 812            
 813        Returns:
 814            c_salt_calculated('float'): Calculated salt concentration added to 'espresso_system'.
 815        """ 
 816        cation_tpl = self.db.get_template(pmb_type="particle",
 817                                          name=cation_name)
 818        cation_state = self.db.get_template(pmb_type="particle_state",
 819                                            name=cation_tpl.initial_state)
 820        cation_charge = cation_state.z
 821        anion_tpl = self.db.get_template(pmb_type="particle",
 822                                          name=anion_name)
 823        anion_state = self.db.get_template(pmb_type="particle_state",
 824                                            name=anion_tpl.initial_state)
 825        anion_charge = anion_state.z
 826        if cation_charge <= 0:
 827            raise ValueError(f'ERROR cation charge must be positive, charge {cation_charge}')
 828        if anion_charge >= 0:
 829            raise ValueError(f'ERROR anion charge must be negative, charge {anion_charge}')
 830        # Calculate the number of ions in the simulation box
 831        volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3')
 832        if c_salt.check('[substance] [length]**-3'):
 833            N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude)
 834            c_salt_calculated=N_ions/(volume*self.N_A)
 835        elif c_salt.check('[length]**-3'):
 836            N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude)
 837            c_salt_calculated=N_ions/volume
 838        else:
 839            raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt)
 840        N_cation = N_ions*abs(anion_charge)
 841        N_anion = N_ions*abs(cation_charge)
 842        self.create_particle(espresso_system=espresso_system, 
 843                             name=cation_name, 
 844                             number_of_particles=N_cation)
 845        self.create_particle(espresso_system=espresso_system, 
 846                             name=anion_name, 
 847                             number_of_particles=N_anion)
 848        if c_salt_calculated.check('[substance] [length]**-3'):
 849            logging.info(f"added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions")
 850        elif c_salt_calculated.check('[length]**-3'):
 851            logging.info(f"added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions")
 852        return c_salt_calculated
 853
 854    def create_bond(self, particle_id1, particle_id2, espresso_system, use_default_bond=False):
 855        """
 856        Creates a bond between two particle instances in an ESPResSo system and registers it in the pyMBE database.
 857
 858        This method performs the following steps:
 859            1. Retrieves the particle instances corresponding to 'particle_id1' and 'particle_id2' from the database.
 860            2. Retrieves or creates the corresponding ESPResSo bond instance using the bond template.
 861            3. Adds the ESPResSo bond instance to the ESPResSo system if it was newly created.
 862            4. Adds the bond to the first particle's bond list in ESPResSo.
 863            5. Creates a 'BondInstance' in the database and registers it.
 864
 865        Args:
 866            particle_id1 ('int'): 
 867                pyMBE and ESPResSo ID of the first particle.
 868
 869            particle_id2 ('int'): 
 870                pyMBE and ESPResSo ID of the second particle.
 871
 872            espresso_system ('espressomd.system.System'): 
 873                ESPResSo system object where the bond will be created.
 874
 875            use_default_bond ('bool', optional): 
 876                If True, use a default bond template if no specific template exists. Defaults to False.
 877
 878        Returns:
 879            ('int'): 
 880                bond_id of the bond instance created in the pyMBE database.
 881        """
 882        particle_inst_1 = self.db.get_instance(pmb_type="particle",
 883                                               instance_id=particle_id1)
 884        particle_inst_2 = self.db.get_instance(pmb_type="particle",
 885                                               instance_id=particle_id2)
 886        bond_tpl = self.get_bond_template(particle_name1=particle_inst_1.name,
 887                                          particle_name2=particle_inst_2.name,
 888                                          use_default_bond=use_default_bond)
 889        bond_inst = self._get_espresso_bond_instance(bond_template=bond_tpl,
 890                                                    espresso_system=espresso_system)
 891        espresso_system.part.by_id(particle_id1).add_bond((bond_inst, particle_id2))
 892        bond_id = self.db._propose_instance_id(pmb_type="bond")
 893        pmb_bond_instance = BondInstance(bond_id=bond_id,
 894                                         name=bond_tpl.name,
 895                                         particle_id1=particle_id1,
 896                                         particle_id2=particle_id2)
 897        self.db._register_instance(instance=pmb_bond_instance)
 898
 899    def create_counterions(self, object_name, cation_name, anion_name, espresso_system):
 900        """
 901        Creates particles of 'cation_name' and 'anion_name' in 'espresso_system' to counter the net charge of 'object_name'.
 902        
 903        Args:
 904            object_name ('str'): 
 905                'name' of a pyMBE object.
 906
 907            espresso_system ('espressomd.system.System'): 
 908                Instance of a system object from the espressomd library.
 909
 910            cation_name ('str'): 
 911                'name' of a particle with a positive charge.
 912
 913            anion_name ('str'): 
 914                'name' of a particle with a negative charge.
 915
 916        Returns: 
 917            ('dict'): 
 918                {"name": number}
 919
 920        Notes:
 921            This function currently does not support the creation of counterions for hydrogels.
 922        """ 
 923        cation_tpl = self.db.get_template(pmb_type="particle",
 924                                          name=cation_name)
 925        cation_state = self.db.get_template(pmb_type="particle_state",
 926                                            name=cation_tpl.initial_state)
 927        cation_charge = cation_state.z
 928        anion_tpl = self.db.get_template(pmb_type="particle",
 929                                          name=anion_name)
 930        anion_state = self.db.get_template(pmb_type="particle_state",
 931                                            name=anion_tpl.initial_state)
 932        anion_charge = anion_state.z
 933        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
 934        counterion_number={}
 935        object_charge={}
 936        for name in ['positive', 'negative']:
 937            object_charge[name]=0
 938        for id in object_ids:
 939            if espresso_system.part.by_id(id).q > 0:
 940                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 941            elif espresso_system.part.by_id(id).q < 0:
 942                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 943        if object_charge['positive'] % abs(anion_charge) == 0:
 944            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
 945        else:
 946            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
 947        if object_charge['negative'] % abs(cation_charge) == 0:
 948            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
 949        else:
 950            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
 951        if counterion_number[cation_name] > 0: 
 952            self.create_particle(espresso_system=espresso_system, 
 953                                 name=cation_name, 
 954                                 number_of_particles=counterion_number[cation_name])
 955        else:
 956            counterion_number[cation_name]=0
 957        if counterion_number[anion_name] > 0:
 958            self.create_particle(espresso_system=espresso_system, 
 959                                 name=anion_name, 
 960                                 number_of_particles=counterion_number[anion_name])
 961        else:
 962            counterion_number[anion_name] = 0
 963        logging.info('the following counter-ions have been created: ')
 964        for name in counterion_number.keys():
 965            logging.info(f'Ion type: {name} created number: {counterion_number[name]}')
 966        return counterion_number
 967
 968    def create_hydrogel(self, name, espresso_system, use_default_bond=False, gen_angle=False):
 969        """ 
 970        Creates a hydrogel in espresso_system using a pyMBE hydrogel template given by 'name'
 971
 972        Args:
 973            name ('str'): 
 974                name of the hydrogel template in the pyMBE database.
 975
 976            espresso_system ('espressomd.system.System'): 
 977                ESPResSo system object where the hydrogel will be created.
 978
 979            use_default_bond ('bool', optional): 
 980                If True, use a default bond template if no specific template exists. Defaults to False.
 981
 982            gen_angle ('bool', optional):
 983                If True, generate angle potentials for the internal hydrogel
 984                chains and, when explicitly defined, for all crosslinker-adjacent
 985                triplets. Defaults to False.
 986
 987        Returns:
 988            ('int'): id of the hydrogel instance created.
 989        """
 990        if not self.db._has_template(name=name, pmb_type="hydrogel"):
 991            raise ValueError(f"Hydrogel template with name '{name}' is not defined in the pyMBE database.")
 992        hydrogel_tpl = self.db.get_template(pmb_type="hydrogel",
 993                                            name=name)
 994        assembly_id = self.db._propose_instance_id(pmb_type="hydrogel")
 995        # Create the nodes
 996        nodes = {}
 997        hydrogel_angle_centers = set()
 998        node_topology = hydrogel_tpl.node_map
 999        for node in node_topology:
1000            node_index = node.lattice_index
1001            node_name = node.particle_name
1002            node_pos, node_id = self._create_hydrogel_node(node_index=node_index,
1003                                                          node_name=node_name,
1004                                                          espresso_system=espresso_system)
1005            node_label = self.lattice_builder._create_node_label(node_index=node_index)
1006            nodes[node_label] = {"name": node_name, "id": node_id, "pos": node_pos} 
1007            self.db._update_instance(instance_id=node_id,
1008                                     pmb_type="particle",
1009                                     attribute="assembly_id",
1010                                     value=assembly_id)
1011        for hydrogel_chain in hydrogel_tpl.chain_map:
1012            molecule_id = self._create_hydrogel_chain(hydrogel_chain=hydrogel_chain,
1013                                                      nodes=nodes, 
1014                                                      espresso_system=espresso_system,
1015                                                      use_default_bond=use_default_bond,
1016                                                      gen_angle=gen_angle)
1017            self.db._update_instance(instance_id=molecule_id,
1018                                     pmb_type="molecule",
1019                                     attribute="assembly_id",
1020                                     value=assembly_id)
1021            if gen_angle:
1022                residue_ids = self.db._find_instance_ids_by_attribute(pmb_type="residue",
1023                                                                      attribute="molecule_id",
1024                                                                      value=molecule_id)
1025                first_residue_id = min(residue_ids)
1026                last_residue_id = max(residue_ids)
1027                first_residue = self.db.get_instance(pmb_type="residue",
1028                                                     instance_id=first_residue_id)
1029                last_residue = self.db.get_instance(pmb_type="residue",
1030                                                    instance_id=last_residue_id)
1031                first_central_bead_name = self.db.get_template(pmb_type="residue",
1032                                                               name=first_residue.name).central_bead
1033                last_central_bead_name = self.db.get_template(pmb_type="residue",
1034                                                              name=last_residue.name).central_bead
1035                particle_instances = self.db.get_instances(pmb_type="particle")
1036                first_residue_particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1037                                                                                      attribute="residue_id",
1038                                                                                      value=first_residue_id)
1039                last_residue_particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1040                                                                                     attribute="residue_id",
1041                                                                                     value=last_residue_id)
1042                first_bead_id = None
1043                for particle_id in first_residue_particle_ids:
1044                    if particle_instances[particle_id].name == first_central_bead_name:
1045                        first_bead_id = particle_id
1046                        break
1047
1048                last_bead_id = None
1049                for particle_id in last_residue_particle_ids:
1050                    if particle_instances[particle_id].name == last_central_bead_name:
1051                        last_bead_id = particle_id
1052                        break
1053                node_start_label = self.lattice_builder._create_node_label(hydrogel_chain.node_start)
1054                node_end_label = self.lattice_builder._create_node_label(hydrogel_chain.node_end)
1055                hydrogel_angle_centers.update({
1056                    nodes[node_start_label]["id"],
1057                    nodes[node_end_label]["id"],
1058                    first_bead_id,
1059                    last_bead_id,
1060                })
1061        self.db._propagate_id(root_type="hydrogel", 
1062                                root_id=assembly_id, 
1063                                attribute="assembly_id", 
1064                                value=assembly_id)
1065        if gen_angle:
1066            self._generate_hydrogel_crosslinker_angles(espresso_system=espresso_system,
1067                                                       central_particle_ids=hydrogel_angle_centers)
1068        # Register an hydrogel instance in the pyMBE databasegit 
1069        self.db._register_instance(HydrogelInstance(name=name,
1070                                                    assembly_id=assembly_id))
1071        return assembly_id
1072
1073    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False, reverse_residue_order = False, gen_angle=False):
1074        """
1075        Creates instances of a given molecule template name into ESPResSo.
1076
1077        Args:
1078            name ('str'): 
1079                Label of the molecule type to be created. 'name'.
1080
1081            espresso_system ('espressomd.system.System'): 
1082                Instance of a system object from espressomd library.
1083
1084            number_of_molecules ('int'): 
1085                Number of molecules or peptides of type 'name' to be created.
1086
1087            list_of_first_residue_positions ('list', optional): 
1088                List of coordinates where the central bead of the first_residue_position will be created, random by default.
1089
1090            backbone_vector ('list' of 'float'): 
1091                Backbone vector of the molecule, random by default. Central beads of the residues in the 'residue_list' are placed along this vector. 
1092
1093            use_default_bond('bool', optional): 
1094                Controls if a bond of type 'default' is used to bond particles with undefined bonds in the pyMBE database.
1095
1096            reverse_residue_order('bool', optional): 
1097                Creates residues in reverse sequential order than the one defined in the molecule template. Defaults to False.
1098
1099        Returns:
1100            ('list' of 'int'): 
1101                List with the 'molecule_id' of the pyMBE molecule instances created into 'espresso_system'.
1102
1103        Notes:
1104            - This function can be used to create both molecules and peptides.    
1105        """
1106        pmb_type = self._get_template_type(name=name,
1107                                           allowed_types={"molecule", "peptide"})
1108        if number_of_molecules <= 0:
1109            return {}
1110        if list_of_first_residue_positions is not None:
1111            for item in list_of_first_residue_positions:
1112                if not isinstance(item, list):
1113                    raise ValueError("The provided input position is not a nested list. Should be a nested list with elements of 3D lists, corresponding to xyz coord.")
1114                elif len(item) != 3:
1115                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
1116
1117            if len(list_of_first_residue_positions) != number_of_molecules:
1118                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
1119        # Generate an arbitrary random unit vector
1120        if backbone_vector is None:
1121            backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
1122                                                                      radius=1, 
1123                                                                      n_samples=1,
1124                                                                      on_surface=True)[0]
1125        else:
1126            backbone_vector = np.array(backbone_vector)
1127        first_residue = True
1128        molecule_tpl = self.db.get_template(pmb_type=pmb_type,
1129                                            name=name)
1130        if reverse_residue_order:
1131            residue_list = molecule_tpl.residue_list[::-1]
1132        else:
1133            residue_list = molecule_tpl.residue_list
1134        pos_index = 0 
1135        molecule_ids = []
1136        for _ in range(number_of_molecules):        
1137            molecule_id = self.db._propose_instance_id(pmb_type=pmb_type)
1138            for residue in residue_list:
1139                if first_residue:
1140                    if list_of_first_residue_positions is None:
1141                        central_bead_pos = None
1142                    else:
1143                        for item in list_of_first_residue_positions:
1144                            central_bead_pos = [np.array(list_of_first_residue_positions[pos_index])]
1145                            
1146                    residue_id = self.create_residue(name=residue,
1147                                                     espresso_system=espresso_system, 
1148                                                     central_bead_position=central_bead_pos,  
1149                                                     use_default_bond= use_default_bond, 
1150                                                     backbone_vector=backbone_vector)
1151                    
1152                    # Add molecule_id to the residue instance and all particles associated
1153                    self.db._propagate_id(root_type="residue", 
1154                                          root_id=residue_id,
1155                                          attribute="molecule_id", 
1156                                          value=molecule_id)
1157                    particle_ids_in_residue = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1158                                                                                      attribute="residue_id",
1159                                                                                      value=residue_id)
1160                    prev_central_bead_id = particle_ids_in_residue[0]
1161                    prev_central_bead_name = self.db.get_instance(pmb_type="particle", 
1162                                                                  instance_id=prev_central_bead_id).name
1163                    prev_central_bead_pos = espresso_system.part.by_id(prev_central_bead_id).pos
1164                    first_residue = False          
1165                else:
1166                    
1167                    # Calculate the starting position of the new residue
1168                    residue_tpl = self.db.get_template(pmb_type="residue",
1169                                                       name=residue)
1170                    lj_parameters = self.get_lj_parameters(particle_name1=prev_central_bead_name,
1171                                                           particle_name2=residue_tpl.central_bead)
1172                    bond_tpl = self.get_bond_template(particle_name1=prev_central_bead_name,
1173                                                      particle_name2=residue_tpl.central_bead,
1174                                                      use_default_bond=use_default_bond)
1175                    l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1176                                                          bond_type=bond_tpl.bond_type,
1177                                                          bond_parameters=bond_tpl.get_parameters(ureg=self.units))
1178                    central_bead_pos = prev_central_bead_pos+backbone_vector*l0
1179                    # Create the residue
1180                    residue_id = self.create_residue(name=residue, 
1181                                                     espresso_system=espresso_system, 
1182                                                     central_bead_position=[central_bead_pos],
1183                                                     use_default_bond= use_default_bond, 
1184                                                     backbone_vector=backbone_vector)
1185                    # Add molecule_id to the residue instance and all particles associated
1186                    self.db._propagate_id(root_type="residue", 
1187                                          root_id=residue_id, 
1188                                          attribute="molecule_id", 
1189                                          value=molecule_id)
1190                    particle_ids_in_residue = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1191                                                                                      attribute="residue_id",
1192                                                                                      value=residue_id)
1193                    central_bead_id = particle_ids_in_residue[0]
1194
1195                    # Bond the central beads of the new and previous residues
1196                    self.create_bond(particle_id1=prev_central_bead_id,
1197                                     particle_id2=central_bead_id,
1198                                     espresso_system=espresso_system,
1199                                     use_default_bond=use_default_bond)
1200                    
1201                    prev_central_bead_id = central_bead_id                    
1202                    prev_central_bead_name = self.db.get_instance(pmb_type="particle", instance_id=central_bead_id).name
1203                    prev_central_bead_pos =central_bead_pos
1204            # Create a Peptide or Molecule instance and register it on the pyMBE database
1205            if pmb_type == "molecule":
1206                inst = MoleculeInstance(molecule_id=molecule_id,
1207                                        name=name)
1208            elif pmb_type == "peptide":
1209                inst = PeptideInstance(name=name,
1210                                       molecule_id=molecule_id)
1211            self.db._register_instance(inst)
1212            if gen_angle:
1213                self._generate_angles_for_entity(
1214                    espresso_system=espresso_system,
1215                    entity_id=molecule_id,
1216                    entity_id_col='molecule_id')
1217            first_residue = True
1218            pos_index+=1
1219            molecule_ids.append(molecule_id)
1220        return molecule_ids
1221    
1222    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
1223        """
1224        Creates one or more particles in an ESPResSo system based on the particle template in the pyMBE database.
1225        
1226        Args:
1227            name ('str'): 
1228                Label of the particle template in the pyMBE database. 
1229
1230            espresso_system ('espressomd.system.System'): 
1231                Instance of a system object from the espressomd library.
1232
1233            number_of_particles ('int'): 
1234                Number of particles to be created.
1235
1236            position (list of ['float','float','float'], optional): 
1237                Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
1238
1239            fix ('bool', optional): 
1240                Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
1241
1242        Returns:
1243            ('list' of 'int'): 
1244                List with the ids of the particles created into 'espresso_system'.
1245        """       
1246        if number_of_particles <=0:
1247            return []
1248        if not self.db._has_template(name=name, pmb_type="particle"):
1249            raise ValueError(f"Particle template with name '{name}' is not defined in the pyMBE database.")
1250        
1251        part_tpl = self.db.get_template(pmb_type="particle",
1252                                        name=name)
1253        part_state = self.db.get_template(pmb_type="particle_state",
1254                                         name=part_tpl.initial_state)
1255        z = part_state.z
1256        es_type = part_state.es_type
1257        # Create the new particles into  ESPResSo 
1258        created_pid_list=[]
1259        for index in range(number_of_particles):
1260            if position is None:
1261                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1262            else:
1263                particle_position = position[index]
1264            
1265            particle_id = self.db._propose_instance_id(pmb_type="particle")
1266            created_pid_list.append(particle_id)
1267            kwargs = dict(id=particle_id, pos=particle_position, type=es_type, q=z)
1268            if fix:
1269                kwargs["fix"] = 3 * [fix]
1270            espresso_system.part.add(**kwargs)
1271            part_inst = ParticleInstance(name=name,
1272                                         particle_id=particle_id,
1273                                         initial_state=part_state.name)
1274            self.db._register_instance(part_inst)
1275                              
1276        return created_pid_list
1277
1278    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1279        """
1280        Creates one or more protein molecules in an ESPResSo system based on the 
1281        protein template in the pyMBE database and a provided topology.
1282
1283        Args:
1284            name (str):
1285                Name of the protein template stored in the pyMBE database.
1286            
1287            number_of_proteins (int):
1288                Number of protein molecules to generate.  
1289            
1290            espresso_system (espressomd.system.System):
1291                The ESPResSo simulation system where the protein molecules will be created.
1292            
1293            topology_dict (dict):
1294                Dictionary defining the internal structure of the protein. Expected format:
1295                    {"ResidueName1": {"initial_pos": np.ndarray,
1296                                      "chain_id": int,
1297                                      "radius": float},
1298                     "ResidueName2": { ... },
1299                        ...
1300                    }
1301                The '"initial_pos"' entry is required and represents the residue’s
1302                reference coordinates before shifting to the protein's center-of-mass.
1303
1304        Returns:
1305            ('list' of 'int'): 
1306                List of the molecule_id of the Protein instances created into ESPResSo.
1307
1308        Notes:
1309            - Particles are created using 'create_particle()' with 'fix=True',
1310            meaning they are initially immobilized.
1311            - The function assumes all residues in 'topology_dict' correspond to
1312            particle templates already defined in the pyMBE database.
1313            - Bonds between residues are not created here; it assumes a rigid body representation of the protein.
1314        """
1315        if number_of_proteins <= 0:
1316            return
1317        if not self.db._has_template(name=name, pmb_type="protein"):
1318            raise ValueError(f"Protein template with name '{name}' is not defined in the pyMBE database.")
1319        protein_tpl = self.db.get_template(pmb_type="protein", name=name)
1320        box_half = espresso_system.box_l[0] / 2.0
1321        # Create protein
1322        mol_ids = []
1323        for _ in range(number_of_proteins):
1324            # create a molecule identifier in pyMBE
1325            molecule_id = self.db._propose_instance_id(pmb_type="protein")
1326            # place protein COM randomly
1327            protein_center = self.generate_coordinates_outside_sphere(radius=1,
1328                                                                      max_dist=box_half,
1329                                                                      n_samples=1,
1330                                                                      center=[box_half]*3)[0]
1331            residues = hf.get_residues_from_topology_dict(topology_dict=topology_dict,
1332                                                         model=protein_tpl.model)
1333            # CREATE RESIDUES + PARTICLES
1334            for _, rdata in residues.items():
1335                base_resname = rdata["resname"]  
1336                residue_name = f"AA-{base_resname}"
1337                # residue instance ID
1338                residue_id = self.db._propose_instance_id("residue")
1339                # register ResidueInstance
1340                self.db._register_instance(ResidueInstance(name=residue_name,
1341                                                           residue_id=residue_id,
1342                                                           molecule_id=molecule_id))
1343                # PARTICLE CREATION
1344                for bead_id in rdata["beads"]:
1345                    bead_type = re.split(r'\d+', bead_id)[0]
1346                    relative_pos = topology_dict[bead_id]["initial_pos"]
1347                    absolute_pos = relative_pos + protein_center
1348                    particle_id = self.create_particle(name=bead_type,
1349                                                       espresso_system=espresso_system,
1350                                                       number_of_particles=1,
1351                                                       position=[absolute_pos],
1352                                                       fix=True)[0]
1353                    # update metadata
1354                    self.db._update_instance(instance_id=particle_id,
1355                                             pmb_type="particle",
1356                                             attribute="molecule_id",
1357                                             value=molecule_id)
1358                    self.db._update_instance(instance_id=particle_id,
1359                                             pmb_type="particle",
1360                                             attribute="residue_id",
1361                                             value=residue_id)
1362            protein_inst = ProteinInstance(name=name,
1363                                           molecule_id=molecule_id)
1364            self.db._register_instance(protein_inst)
1365            mol_ids.append(molecule_id)
1366        return mol_ids
1367
1368    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None, gen_angle=False):
1369        """
1370        Creates a residue  into ESPResSo.
1371
1372        Args:
1373            name ('str'): 
1374                Label of the residue type to be created. 
1375
1376            espresso_system ('espressomd.system.System'): 
1377                Instance of a system object from espressomd library.
1378
1379            central_bead_position ('list' of 'float'): 
1380                Position of the central bead.
1381
1382            use_default_bond ('bool'): 
1383                Switch to control if a bond of type 'default' is used to bond a particle whose bonds types are not defined in the pyMBE database.
1384
1385            backbone_vector ('list' of 'float'): 
1386                Backbone vector of the molecule. All side chains are created perpendicularly to 'backbone_vector'.
1387
1388        Returns:
1389            (int): 
1390                residue_id of the residue created.
1391        """
1392        if not self.db._has_template(name=name, pmb_type="residue"):
1393            raise ValueError(f"Residue template with name '{name}' is not defined in the pyMBE database.")
1394        res_tpl = self.db.get_template(pmb_type="residue",
1395                                       name=name)
1396        # Assign a residue_id
1397        residue_id = self.db._propose_instance_id(pmb_type="residue")
1398        res_inst = ResidueInstance(name=name,
1399                                   residue_id=residue_id)
1400        self.db._register_instance(res_inst)
1401        # create the principal bead   
1402        central_bead_name = res_tpl.central_bead 
1403        central_bead_id = self.create_particle(name=central_bead_name,
1404                                               espresso_system=espresso_system,
1405                                               position=central_bead_position,
1406                                               number_of_particles = 1)[0]
1407        
1408        central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1409        # Assigns residue_id to the central_bead particle created.
1410        self.db._update_instance(pmb_type="particle",
1411                                 instance_id=central_bead_id,
1412                                 attribute="residue_id",
1413                                 value=residue_id)
1414        
1415        # create the lateral beads  
1416        side_chain_list = res_tpl.side_chains
1417        side_chain_beads_ids = []
1418        for side_chain_name in side_chain_list:
1419            pmb_type = self._get_template_type(name=side_chain_name,
1420                                               allowed_types={"particle", "residue"})
1421            if pmb_type == 'particle':
1422                lj_parameters = self.get_lj_parameters(particle_name1=central_bead_name,
1423                                                       particle_name2=side_chain_name)
1424                bond_tpl = self.get_bond_template(particle_name1=central_bead_name,
1425                                                  particle_name2=side_chain_name,
1426                                                  use_default_bond=use_default_bond)
1427                l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1428                                                      bond_type=bond_tpl.bond_type,
1429                                                      bond_parameters=bond_tpl.get_parameters(ureg=self.units))               
1430                if backbone_vector is None:
1431                    bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1432                                                                radius=l0, 
1433                                                                n_samples=1,
1434                                                                on_surface=True)[0]
1435                else:
1436                    bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector),
1437                                                                                                magnitude=l0)
1438                    
1439                side_bead_id = self.create_particle(name=side_chain_name, 
1440                                                    espresso_system=espresso_system,
1441                                                    position=[bead_position], 
1442                                                    number_of_particles=1)[0]
1443                side_chain_beads_ids.append(side_bead_id)
1444                self.db._update_instance(pmb_type="particle",
1445                                         instance_id=side_bead_id,
1446                                         attribute="residue_id",
1447                                         value=residue_id)
1448                self.create_bond(particle_id1=central_bead_id,
1449                                 particle_id2=side_bead_id,
1450                                 espresso_system=espresso_system,
1451                                 use_default_bond=use_default_bond)
1452            elif pmb_type == 'residue':
1453                side_residue_tpl = self.db.get_template(name=side_chain_name,
1454                                                        pmb_type=pmb_type)
1455                central_bead_side_chain = side_residue_tpl.central_bead
1456                lj_parameters = self.get_lj_parameters(particle_name1=central_bead_name,
1457                                                       particle_name2=central_bead_side_chain)
1458                bond_tpl = self.get_bond_template(particle_name1=central_bead_name,
1459                                                  particle_name2=central_bead_side_chain,
1460                                                  use_default_bond=use_default_bond)
1461                l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1462                                                      bond_type=bond_tpl.bond_type,
1463                                                      bond_parameters=bond_tpl.get_parameters(ureg=self.units))
1464                if backbone_vector is None:
1465                    residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1466                                                                radius=l0, 
1467                                                                n_samples=1,
1468                                                                on_surface=True)[0]
1469                else:
1470                    residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1471                                                                                                    magnitude=l0)
1472                side_residue_id = self.create_residue(name=side_chain_name, 
1473                                                      espresso_system=espresso_system,
1474                                                      central_bead_position=[residue_position],
1475                                                      use_default_bond=use_default_bond)
1476                # Find particle ids of the inner residue
1477                side_chain_beads_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1478                                                                               attribute="residue_id",
1479                                                                               value=side_residue_id)
1480                # Change the residue_id of the residue in the side chain to the one of the outer residue
1481                for particle_id in side_chain_beads_ids:
1482                    self.db._update_instance(instance_id=particle_id,
1483                                             pmb_type="particle",
1484                                             attribute="residue_id",
1485                                             value=residue_id)
1486                # Remove the instance of the inner residue
1487                self.db.delete_instance(pmb_type="residue",
1488                                        instance_id=side_residue_id)
1489                self.create_bond(particle_id1=central_bead_id,
1490                                 particle_id2=side_chain_beads_ids[0],
1491                                 espresso_system=espresso_system,
1492                                 use_default_bond=use_default_bond)
1493        if gen_angle:
1494            self._generate_angles_for_entity(espresso_system=espresso_system,
1495                                             entity_id=residue_id,
1496                                             entity_id_col="residue_id")
1497        return  residue_id
1498
1499    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1500        """
1501        Defines bond templates for each particle pair in 'particle_pairs' in the pyMBE database.
1502
1503        Args:
1504            bond_type ('str'): 
1505                label to identify the potential to model the bond.
1506
1507            bond_parameters ('dict'): 
1508                parameters of the potential of the bond.
1509
1510            particle_pairs ('lst'): 
1511                list of the 'names' of the 'particles' to be bonded.
1512
1513        Notes:
1514            -Currently, only HARMONIC and FENE bonds are supported.
1515            - For a HARMONIC bond the dictionary must contain the following parameters:
1516                - k ('pint.Quantity')      : Magnitude of the bond. It should have units of energy/length**2 
1517                using the 'pmb.units' UnitRegistry.
1518                - r_0 ('pint.Quantity')    : Equilibrium bond length. It should have units of length using 
1519                the 'pmb.units' UnitRegistry.
1520           - For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:              
1521                - d_r_max ('pint.Quantity'): Maximal stretching length for FENE. It should have 
1522                units of length using the 'pmb.units' UnitRegistry. Default 'None'.
1523        """
1524        self._check_bond_inputs(bond_parameters=bond_parameters,
1525                                bond_type=bond_type)
1526        parameters_expected_dimensions={"r_0": "length",
1527                                        "k": "energy/length**2",
1528                                        "d_r_max": "length"}
1529
1530        parameters_tpl = {}
1531        for key in bond_parameters.keys():
1532            parameters_tpl[key]= PintQuantity.from_quantity(q=bond_parameters[key],
1533                                                            expected_dimension=parameters_expected_dimensions[key],
1534                                                            ureg=self.units)
1535
1536        bond_names=[]
1537        for particle_name1, particle_name2 in particle_pairs:
1538            
1539            tpl = BondTemplate(particle_name1=particle_name1,
1540                               particle_name2=particle_name2,
1541                               parameters=parameters_tpl,
1542                               bond_type=bond_type)
1543            tpl._make_name()
1544            if tpl.name in bond_names:
1545                raise RuntimeError(f"Bond {tpl.name} has already been defined, please check the list of particle pairs")
1546            bond_names.append(tpl.name)
1547            self.db._register_template(tpl)
1548
1549    
1550    def define_default_bond(self, bond_type, bond_parameters):
1551        """
1552        Defines a bond template as a "default" template in the pyMBE database.
1553        
1554        Args:
1555            bond_type ('str'): 
1556                label to identify the potential to model the bond.
1557
1558            bond_parameters ('dict'): 
1559                parameters of the potential of the bond.
1560            
1561        Notes:
1562            - Currently, only harmonic and FENE bonds are supported. 
1563        """
1564        self._check_bond_inputs(bond_parameters=bond_parameters,
1565                                bond_type=bond_type)
1566        parameters_expected_dimensions={"r_0": "length",
1567                                        "k": "energy/length**2",
1568                                        "d_r_max": "length"}
1569        parameters_tpl = {}
1570        for key in bond_parameters.keys():
1571            parameters_tpl[key]= PintQuantity.from_quantity(q=bond_parameters[key],
1572                                                            expected_dimension=parameters_expected_dimensions[key],
1573                                                            ureg=self.units)
1574        tpl = BondTemplate(parameters=parameters_tpl,
1575                               bond_type=bond_type)
1576        tpl.name = "default"
1577        self.db._register_template(tpl)
1578
1579    def define_angular_potential(self, angle_type, angle_parameters, particle_triplets):
1580        """
1581        Defines angle potential templates for each particle triplet in `particle_triplets`.
1582
1583        Args:
1584            angle_type ('str'):
1585                Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
1586
1587            angle_parameters ('dict'):
1588                Parameters of the angle potential. Must contain:
1589                    - "k" ('pint.Quantity'): Bending stiffness with dimensions of energy.
1590                    - "phi_0" ('float'): Equilibrium angle in radians.
1591
1592            particle_triplets ('list[tuple[str,str,str]]'):
1593                List of (side_particle1, central_particle, side_particle2) triplets.
1594        """
1595        valid_angle_types = ["harmonic", "cosine", "harmonic_cosine"]
1596        if angle_type not in valid_angle_types:
1597            raise NotImplementedError(f"Angle potential type '{angle_type}' currently not implemented in pyMBE, accepted types are {valid_angle_types}")
1598
1599        if "k" not in angle_parameters:
1600            raise ValueError("Magnitude of the angle potential (k) is missing")
1601        if "phi_0" not in angle_parameters:
1602            raise ValueError("Equilibrium angle (phi_0) is missing")
1603
1604        parameters_tpl = {
1605            "k": PintQuantity.from_quantity(q=angle_parameters["k"],
1606                                            expected_dimension="energy",
1607                                            ureg=self.units),
1608            "phi_0": PintQuantity.from_quantity(q=angle_parameters["phi_0"],
1609                                                expected_dimension="dimensionless",
1610                                                ureg=self.units),
1611        }
1612
1613        angle_names = []
1614        for side1, central, side2 in particle_triplets:
1615            tpl = AngleTemplate(side_particle1=side1,
1616                                central_particle=central,
1617                                side_particle2=side2,
1618                                parameters=parameters_tpl,
1619                                angle_type=angle_type)
1620            tpl._make_name()
1621            if tpl.name in angle_names:
1622                raise RuntimeError(f"Angle {tpl.name} has already been defined, please check the list of particle triplets")
1623            angle_names.append(tpl.name)
1624            self.db._register_template(tpl)
1625
1626    def define_default_angular_potential(self, angle_type, angle_parameters):
1627        """
1628        Defines an angle template as a "default" template in the pyMBE database.
1629
1630        Args:
1631            angle_type ('str'):
1632                Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
1633
1634            angle_parameters ('dict'):
1635                Parameters of the angle potential (k, phi_0).
1636        """
1637        valid_angle_types = ["harmonic", "cosine", "harmonic_cosine"]
1638        if angle_type not in valid_angle_types:
1639            raise NotImplementedError(f"Angle potential type '{angle_type}' currently not implemented in pyMBE, accepted types are {valid_angle_types}")
1640
1641        if "k" not in angle_parameters:
1642            raise ValueError("Magnitude of the angle potential (k) is missing")
1643        if "phi_0" not in angle_parameters:
1644            raise ValueError("Equilibrium angle (phi_0) is missing")
1645
1646        parameters_tpl = {
1647            "k": PintQuantity.from_quantity(q=angle_parameters["k"],
1648                                            expected_dimension="energy",
1649                                            ureg=self.units),
1650            "phi_0": PintQuantity.from_quantity(q=angle_parameters["phi_0"],
1651                                                expected_dimension="dimensionless",
1652                                                ureg=self.units),
1653        }
1654        tpl = AngleTemplate(parameters=parameters_tpl,
1655                            angle_type=angle_type)
1656        tpl.name = "default"
1657        self.db._register_template(tpl)
1658
1659    def create_angular_potential(self, particle_id1, particle_id2, particle_id3, espresso_system, use_default_angle=False):
1660        """
1661        Creates an angle between three particle instances in an ESPResSo system
1662        and registers it in the pyMBE database.
1663
1664        Args:
1665            particle_id1 ('int'): ID of the first side particle.
1666            particle_id2 ('int'): ID of the central particle.
1667            particle_id3 ('int'): ID of the second side particle.
1668            espresso_system ('espressomd.system.System'): ESPResSo system.
1669            use_default_angle ('bool', optional): If True, use the default angle if no specific one is found.
1670        """
1671        particle_inst_1 = self.db.get_instance(pmb_type="particle", instance_id=particle_id1)
1672        particle_inst_2 = self.db.get_instance(pmb_type="particle", instance_id=particle_id2)
1673        particle_inst_3 = self.db.get_instance(pmb_type="particle", instance_id=particle_id3)
1674
1675        # Verify that bonds exist between side particles and central particle
1676        bond_instances = self.db.get_instances(pmb_type="bond")
1677        bonded_pairs = set()
1678        for bond in bond_instances.values():
1679            pair = frozenset([bond.particle_id1, bond.particle_id2])
1680            bonded_pairs.add(pair)
1681        if frozenset([particle_id1, particle_id2]) not in bonded_pairs:
1682            raise ValueError(f"Cannot create angle: no bond exists between particle {particle_id1} and central particle {particle_id2}.")
1683        if frozenset([particle_id3, particle_id2]) not in bonded_pairs:
1684            raise ValueError(f"Cannot create angle: no bond exists between particle {particle_id3} and central particle {particle_id2}.")
1685
1686        angle_tpl = self.get_angle_template(side_name1=particle_inst_1.name,
1687                                            central_name=particle_inst_2.name,
1688                                            side_name2=particle_inst_3.name,
1689                                            use_default_angle=use_default_angle)
1690        angle_inst = self._get_espresso_angle_instance(angle_template=angle_tpl, espresso_system=espresso_system)
1691
1692        # ESPResSo angle bonds are added to the central particle
1693        espresso_system.part.by_id(particle_id2).add_bond((angle_inst, particle_id1, particle_id3))
1694
1695        angle_id = self.db._propose_instance_id(pmb_type="angle")
1696        pmb_angle_instance = AngleInstance(angle_id=angle_id,
1697                                           name=angle_tpl.name,
1698                                           particle_id1=particle_id1,
1699                                           particle_id2=particle_id2,
1700                                           particle_id3=particle_id3)
1701        self.db._register_instance(instance=pmb_angle_instance)
1702
1703    def get_angle_template(self, side_name1, central_name, side_name2, use_default_angle=False):
1704        """
1705        Retrieves an angle template connecting three particle templates.
1706
1707        Args:
1708            side_name1 ('str'): Name of the first side particle.
1709            central_name ('str'): Name of the central particle.
1710            side_name2 ('str'): Name of the second side particle.
1711            use_default_angle ('bool', optional): If True, fall back to the default angle template.
1712
1713        Returns:
1714            ('AngleTemplate'): The matching angle template.
1715        """
1716        angle_key = AngleTemplate.make_angle_key(side1=side_name1, central=central_name, side2=side_name2)
1717        try:
1718            return self.db.get_template(name=angle_key, pmb_type="angle")
1719        except ValueError:
1720            pass
1721
1722        if use_default_angle:
1723            return self.db.get_template(name="default", pmb_type="angle")
1724
1725        raise ValueError(f"No angle template found for '{side_name1}-{central_name}-{side_name2}', and default angles are deactivated.")
1726
1727    def _get_espresso_angle_instance(self, angle_template, espresso_system):
1728        """
1729        Retrieve or create an angle interaction in an ESPResSo system for a given angle template.
1730
1731        Args:
1732            angle_template ('AngleTemplate'): The angle template to use.
1733            espresso_system ('espressomd.system.System'): ESPResSo system.
1734
1735        Returns:
1736            ('espressomd.interactions.BondedInteraction'): The ESPResSo angle interaction object.
1737        """
1738        if angle_template.name in self.db.espresso_angle_instances:
1739            return self.db.espresso_angle_instances[angle_template.name]
1740        angle_inst = self._create_espresso_angle_instance(angle_type=angle_template.angle_type,
1741                                                          angle_parameters=angle_template.get_parameters(self.units))
1742        self.db.espresso_angle_instances[angle_template.name] = angle_inst
1743        espresso_system.bonded_inter.add(angle_inst)
1744        return angle_inst
1745
1746    def _create_espresso_angle_instance(self, angle_type, angle_parameters):
1747        """
1748        Creates an ESPResSo angle interaction object.
1749
1750        Args:
1751            angle_type ('str'): Type of angle potential ("harmonic", "cosine", "harmonic_cosine").
1752            angle_parameters ('dict'): Parameters of the angle potential (k, phi_0).
1753
1754        Returns:
1755            ('espressomd.interactions.BondedInteraction'): The ESPResSo angle interaction object.
1756        """
1757        from espressomd import interactions
1758
1759        k = angle_parameters["k"].m_as("reduced_energy")
1760        phi_0 = float(angle_parameters["phi_0"].magnitude)
1761
1762        if angle_type == "harmonic":
1763            return interactions.AngleHarmonic(bend=k, phi0=phi_0)
1764        elif angle_type == "cosine":
1765            return interactions.AngleCosine(bend=k, phi0=phi_0)
1766        elif angle_type == "harmonic_cosine":
1767            return interactions.AngleCossquare(bend=k, phi0=phi_0)
1768
1769    def _generate_angles_for_entity(self, espresso_system, entity_id, entity_id_col):
1770        """
1771        Auto-generates angles from bond topology for an entity (molecule or residue).
1772
1773        For each particle in the entity that has two or more bonded neighbors,
1774        this method finds all neighbor pairs and applies any matching angle potential.
1775
1776        Args:
1777            espresso_system ('espressomd.system.System'): ESPResSo system.
1778            entity_id ('int'): The molecule_id or residue_id to generate angles for.
1779            entity_id_col ('str'): Either "molecule_id" or "residue_id".
1780        """
1781        # Get all particle IDs for this entity
1782        particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1783                                                               attribute=entity_id_col,
1784                                                               value=entity_id)
1785        if not particle_ids:
1786            return
1787
1788        # Build neighbor map from bond instances
1789        neighbors = {pid: set() for pid in particle_ids}
1790        pid_set = set(particle_ids)
1791        bond_instances = self.db.get_instances(pmb_type="bond")
1792        for bond in bond_instances.values():
1793            i, j = bond.particle_id1, bond.particle_id2
1794            if i in pid_set and j in pid_set:
1795                neighbors[i].add(j)
1796                neighbors[j].add(i)
1797
1798        # For each particle with 2+ neighbors, generate angles
1799        for j in particle_ids:
1800            nbs = sorted(neighbors[j])
1801            if len(nbs) < 2:
1802                continue
1803
1804            for idx_i in range(len(nbs)):
1805                for idx_k in range(idx_i + 1, len(nbs)):
1806                    i = nbs[idx_i]
1807                    k = nbs[idx_k]
1808                    try:
1809                        self.create_angular_potential(particle_id1=i,
1810                                          particle_id2=j,
1811                                          particle_id3=k,
1812                                          espresso_system=espresso_system,
1813                                          use_default_angle=True)
1814                    except ValueError:
1815                        # No angle template defined for this triplet — skip
1816                        continue
1817
1818    def define_hydrogel(self, name, node_map, chain_map):
1819        """
1820        Defines a hydrogel template in the pyMBE database.
1821
1822        Args:
1823            name ('str'): 
1824                Unique label that identifies the 'hydrogel'.
1825
1826            node_map ('list of dict'): 
1827                [{"particle_name": , "lattice_index": }, ... ]
1828
1829            chain_map ('list of dict'): 
1830                [{"node_start": , "node_end": , "residue_list": , ... ]
1831        """
1832        # Sanity tests
1833        node_indices = {tuple(entry['lattice_index']) for entry in node_map}                
1834        chain_map_connectivity = set()
1835        for entry in chain_map:
1836            start = self.lattice_builder.node_labels[entry['node_start']]
1837            end = self.lattice_builder.node_labels[entry['node_end']]
1838            chain_map_connectivity.add((start,end))
1839        if self.lattice_builder.lattice.connectivity != chain_map_connectivity:
1840            raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs")
1841        diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices}
1842        if node_indices != diamond_indices:
1843            raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ")
1844        # Register information in the pyMBE database
1845        nodes=[]
1846        for entry in node_map:
1847            nodes.append(HydrogelNode(particle_name=entry["particle_name"],
1848                                      lattice_index=entry["lattice_index"]))
1849        chains=[]
1850        for chain in chain_map:
1851            chains.append(HydrogelChain(node_start=chain["node_start"],
1852                                        node_end=chain["node_end"],
1853                                        molecule_name=chain["molecule_name"]))
1854        tpl = HydrogelTemplate(name=name,
1855                               node_map=nodes,
1856                               chain_map=chains)
1857        self.db._register_template(tpl)
1858
1859    def define_molecule(self, name, residue_list):
1860        """
1861        Defines a molecule template in the pyMBE database.
1862
1863        Args:
1864            name('str'): 
1865                Unique label that identifies the 'molecule'.
1866
1867            residue_list ('list' of 'str'): 
1868                List of the 'name's of the 'residue's  in the sequence of the 'molecule'.  
1869        """
1870        tpl = MoleculeTemplate(name=name,
1871                               residue_list=residue_list)
1872        self.db._register_template(tpl)
1873
1874    def define_monoprototic_acidbase_reaction(self, particle_name, pka, acidity, metadata=None):
1875        """
1876        Defines an acid-base reaction for a monoprototic particle in the pyMBE database.
1877
1878        Args:
1879            particle_name ('str'): 
1880                Unique label that identifies the particle template. 
1881
1882            pka ('float'): 
1883                pka-value of the acid or base.
1884
1885            acidity ('str'): 
1886                Identifies whether if the particle is 'acidic' or 'basic'.
1887
1888            metadata ('dict', optional): 
1889                Additional information to be stored in the reaction. Defaults to None.
1890        """
1891        supported_acidities = ["acidic", "basic"]
1892        if acidity not in supported_acidities:
1893            raise ValueError(f"Unsupported acidity '{acidity}' for particle '{particle_name}'. Supported acidities are {supported_acidities}.")
1894        reaction_type = "monoprotic"
1895        if acidity == "basic":
1896            reaction_type += "_base"
1897        else:
1898            reaction_type += "_acid"
1899        reaction = Reaction(participants=[ReactionParticipant(particle_name=particle_name,
1900                                                              state_name=f"{particle_name}H", 
1901                                                              coefficient=-1),
1902                                          ReactionParticipant(particle_name=particle_name,
1903                                                              state_name=f"{particle_name}",
1904                                                              coefficient=1)],
1905                            reaction_type=reaction_type,
1906                            pK=pka,
1907                            metadata=metadata)
1908        self.db._register_reaction(reaction)
1909
1910    def define_monoprototic_particle_states(self, particle_name, acidity):
1911        """
1912        Defines particle states for a monoprotonic particle template including the charges in each of its possible states. 
1913
1914        Args:
1915            particle_name ('str'): 
1916                Unique label that identifies the particle template. 
1917
1918            acidity ('str'): 
1919                Identifies whether the particle is 'acidic' or 'basic'.
1920        """
1921        acidity_valid_keys = ['acidic', 'basic']
1922        if not pd.isna(acidity):
1923            if acidity not in acidity_valid_keys:
1924                raise ValueError(f"Acidity {acidity} provided for particle name  {particle_name} is not supported. Valid keys are: {acidity_valid_keys}")
1925        if acidity == "acidic":
1926            states = [{"name": f"{particle_name}H", "z": 0}, 
1927                      {"name": f"{particle_name}",  "z": -1}]
1928            
1929        elif acidity == "basic":
1930            states = [{"name": f"{particle_name}H", "z": 1}, 
1931                      {"name": f"{particle_name}",  "z": 0}]
1932        self.define_particle_states(particle_name=particle_name, 
1933                                    states=states)
1934
1935    def define_particle(self, name,  sigma, epsilon, z=0, acidity=pd.NA, pka=pd.NA, cutoff=pd.NA, offset=pd.NA):
1936        """
1937        Defines a particle template in the pyMBE database.
1938
1939        Args:
1940            name('str'):
1941                 Unique label that identifies this particle type.  
1942
1943            sigma('pint.Quantity'): 
1944                Sigma parameter used to set up Lennard-Jones interactions for this particle type. 
1945
1946            epsilon('pint.Quantity'): 
1947                Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe.
1948
1949            z('int', optional): 
1950                Permanent charge number of this particle type. Defaults to 0.
1951
1952            acidity('str', optional): 
1953                Identifies whether if the particle is 'acidic' or 'basic', used to setup constant pH simulations. Defaults to pd.NA.
1954
1955            pka('float', optional):
1956                If 'particle' is an acid or a base, it defines its  pka-value. Defaults to pd.NA.
1957
1958            cutoff('pint.Quantity', optional): 
1959                Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1960
1961            offset('pint.Quantity', optional): 
1962                Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1963            
1964        Notes:
1965            - 'sigma', 'cutoff' and 'offset' must have a dimensitonality of '[length]' and should be defined using pmb.units.
1966            - 'epsilon' must have a dimensitonality of '[energy]' and should be defined using pmb.units.
1967            - 'cutoff' defaults to '2**(1./6.) reduced_length'. 
1968            - 'offset' defaults to 0.
1969            - For more information on 'sigma', 'epsilon', 'cutoff' and 'offset' check 'pmb.setup_lj_interactions()'.
1970        """ 
1971        # If 'cutoff' and 'offset' are not defined, default them to the following values
1972        if pd.isna(cutoff):
1973            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1974        if pd.isna(offset):
1975            offset=self.units.Quantity(0, "reduced_length")
1976        # Define particle states
1977        if acidity is pd.NA:
1978            states = [{"name": f"{name}",  "z": z}]
1979            self.define_particle_states(particle_name=name, 
1980                                        states=states)
1981            initial_state = name
1982        else:
1983            self.define_monoprototic_particle_states(particle_name=name,
1984                                                  acidity=acidity)
1985            initial_state = f"{name}H"
1986            if pka is not pd.NA:
1987                self.define_monoprototic_acidbase_reaction(particle_name=name,
1988                                                           acidity=acidity,
1989                                                           pka=pka)
1990        tpl = ParticleTemplate(name=name, 
1991                               sigma=PintQuantity.from_quantity(q=sigma, expected_dimension="length", ureg=self.units), 
1992                               epsilon=PintQuantity.from_quantity(q=epsilon, expected_dimension="energy", ureg=self.units),
1993                               cutoff=PintQuantity.from_quantity(q=cutoff, expected_dimension="length", ureg=self.units), 
1994                               offset=PintQuantity.from_quantity(q=offset, expected_dimension="length", ureg=self.units),
1995                               initial_state=initial_state)
1996        self.db._register_template(tpl)
1997    
1998    def define_particle_states(self, particle_name, states):
1999        """
2000        Define the chemical states of an existing particle template.
2001
2002        Args:
2003            particle_name ('str'):
2004                Name of a particle template. 
2005
2006            states ('list' of 'dict'):
2007                List of dictionaries defining the particle states. Each dictionary
2008                must contain:
2009                - 'name' ('str'): Name of the particle state (e.g. '"H"', '"-"',
2010                '"neutral"').
2011                - 'z' ('int'): Charge number of the particle in this state.
2012                Example:
2013                states = [{"name": "AH", "z": 0},     # protonated
2014                         {"name": "A-", "z": -1}]    # deprotonated
2015        Notes:
2016            - Each state is assigned a unique Espresso 'es_type' automatically.
2017            - Chemical reactions (e.g. acid–base equilibria) are **not** created by
2018            this method and must be defined separately (e.g. via
2019            'set_particle_acidity()' or custom reaction definitions).
2020            - Particles without explicitly defined states are assumed to have a
2021            single, implicit state with their default charge.
2022        """
2023        for s in states:
2024            state = ParticleStateTemplate(particle_name=particle_name,
2025                                          name=s["name"],
2026                                          z=s["z"],
2027                                          es_type=self.propose_unused_type())
2028            self.db._register_template(state)
2029
2030    def define_peptide(self, name, sequence, model):
2031        """
2032        Defines a peptide template in the pyMBE database.
2033
2034        Args:
2035            name ('str'): 
2036                Unique label that identifies the peptide.
2037
2038            sequence ('str'): 
2039                Sequence of the peptide.
2040
2041            model ('str'): 
2042                Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
2043        """
2044        valid_keys = ['1beadAA','2beadAA']
2045        if model not in valid_keys:
2046            raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
2047        clean_sequence = hf.protein_sequence_parser(sequence=sequence)    
2048        residue_list = self._get_residue_list_from_sequence(sequence=clean_sequence)
2049        tpl = PeptideTemplate(name=name,
2050                            residue_list=residue_list,
2051                            model=model,
2052                            sequence=sequence)
2053        self.db._register_template(tpl)        
2054    
2055    def define_protein(self, name, sequence, model):
2056        """
2057        Defines a protein template in the pyMBE database.
2058
2059        Args:
2060            name ('str'): 
2061                Unique label that identifies the protein.
2062
2063            sequence ('str'): 
2064                Sequence of the protein.
2065
2066            model ('string'): 
2067                Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
2068
2069        Notes:
2070            - Currently, only 'lj_setup_mode="wca"' is supported. This corresponds to setting up the WCA potential.
2071        """
2072        valid_model_keys = ['1beadAA','2beadAA']
2073        if model not in valid_model_keys:
2074            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
2075        
2076        residue_list = self._get_residue_list_from_sequence(sequence=sequence)
2077        tpl = ProteinTemplate(name=name,
2078                              model=model,
2079                              residue_list=residue_list,
2080                              sequence=sequence)
2081        self.db._register_template(tpl)
2082    
2083    def define_residue(self, name, central_bead, side_chains):
2084        """
2085        Defines a residue template in the pyMBE database.
2086
2087        Args:
2088            name ('str'): 
2089                Unique label that identifies the residue.
2090
2091            central_bead ('str'): 
2092                'name' of the 'particle' to be placed as central_bead of the residue.
2093
2094            side_chains('list' of 'str'): 
2095                List of 'name's of the pmb_objects to be placed as side_chains of the residue. Currently, only pyMBE objects of type 'particle' or 'residue' are supported.
2096        """
2097        tpl = ResidueTemplate(name=name,
2098                              central_bead=central_bead,
2099                              side_chains=side_chains)
2100        self.db._register_template(tpl)
2101
2102    def delete_instances_in_system(self, instance_id, pmb_type, espresso_system):
2103        """
2104        Deletes the instance with instance_id from the ESPResSo system. 
2105        Related assembly, molecule, residue, particles and bond instances will also be deleted from the pyMBE dataframe.
2106
2107        Args:
2108            instance_id ('int'): 
2109                id of the assembly to be deleted. 
2110
2111            pmb_type ('str'): 
2112                the instance type to be deleted. 
2113
2114            espresso_system ('espressomd.system.System'): 
2115                Instance of a system class from espressomd library.
2116        """
2117        if pmb_type == "particle":
2118            instance_identifier = "particle_id"
2119        elif pmb_type == "residue":
2120            instance_identifier = "residue_id"
2121        elif pmb_type in self.db._molecule_like_types:
2122            instance_identifier = "molecule_id"
2123        elif pmb_type in self.db._assembly_like_types:
2124            instance_identifier = "assembly_id"
2125        particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
2126                                                               attribute=instance_identifier,
2127                                                               value=instance_id)
2128        self._delete_particles_from_espresso(particle_ids=particle_ids,
2129                                             espresso_system=espresso_system)
2130        self.db.delete_instance(pmb_type=pmb_type,
2131                                instance_id=instance_id)
2132
2133    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
2134        """
2135        Determines ionic concentrations in the reservoir at fixed pH and salt concentration.
2136
2137        Args:
2138            pH_res ('float'):
2139                Target pH value in the reservoir.
2140
2141            c_salt_res ('pint.Quantity'):
2142                Concentration of monovalent salt (e.g., NaCl) in the reservoir.
2143
2144            activity_coefficient_monovalent_pair ('callable'):
2145                Function returning the activity coefficient of a monovalent ion pair
2146                as a function of ionic strength:
2147                'gamma = activity_coefficient_monovalent_pair(I)'.
2148
2149            max_number_sc_runs ('int', optional):
2150                Maximum number of self-consistent iterations allowed before
2151                convergence is enforced. Defaults to 200.
2152
2153        Returns:
2154            tuple:
2155                (cH_res, cOH_res, cNa_res, cCl_res)
2156                - cH_res ('pint.Quantity'): Concentration of H⁺ ions.
2157                - cOH_res ('pint.Quantity'): Concentration of OH⁻ ions.
2158                - cNa_res ('pint.Quantity'): Concentration of Na⁺ ions.
2159                - cCl_res ('pint.Quantity'): Concentration of Cl⁻ ions.
2160
2161        Notess:
2162            - The algorithm enforces electroneutrality in the reservoir.
2163            - Water autodissociation is included via the equilibrium constant 'Kw'.
2164            - Non-ideal effects enter through activity coefficients depending on
2165            ionic strength.
2166            - The implementation follows the self-consistent scheme described in
2167            Landsgesell (PhD thesis, Sec. 5.3, doi:10.18419/opus-10831), adapted
2168            from the original code (doi:10.18419/darus-2237).
2169        """
2170        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
2171            """
2172            Iteratively determines reservoir ion concentrations self-consistently.
2173
2174            Args:
2175                cH_res ('pint.Quantity'):
2176                    Current estimate of the H⁺ concentration.
2177                c_salt_res ('pint.Quantity'):
2178                    Concentration of monovalent salt in the reservoir.
2179
2180            Returns:
2181                'tuple':
2182                    (cH_res, cOH_res, cNa_res, cCl_res)
2183            """
2184            # Initial ideal estimate
2185            cOH_res = self.Kw / cH_res
2186            if cOH_res >= cH_res:
2187                cNa_res = c_salt_res + (cOH_res - cH_res)
2188                cCl_res = c_salt_res
2189            else:
2190                cCl_res = c_salt_res + (cH_res - cOH_res)
2191                cNa_res = c_salt_res
2192            # Self-consistent iteration
2193            for _ in range(max_number_sc_runs):
2194                ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2195                cOH_new = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
2196                if cOH_new >= cH_res:
2197                    cNa_new = c_salt_res + (cOH_new - cH_res)
2198                    cCl_new = c_salt_res
2199                else:
2200                    cCl_new = c_salt_res + (cH_res - cOH_new)
2201                    cNa_new = c_salt_res
2202                # Update values
2203                cOH_res = cOH_new
2204                cNa_res = cNa_new
2205                cCl_res = cCl_new
2206            return cH_res, cOH_res, cNa_res, cCl_res
2207        # Initial guess for H+ concentration from target pH
2208        cH_res = 10 ** (-pH_res) * self.units.mol / self.units.l
2209        # First self-consistent solve
2210        cH_res, cOH_res, cNa_res, cCl_res = (determine_reservoir_concentrations_selfconsistently(cH_res, 
2211                                                                                                 c_salt_res))
2212        ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2213        determined_pH = -np.log10(cH_res.to("mol/L").magnitude* np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
2214        # Outer loop to enforce target pH
2215        while abs(determined_pH - pH_res) > 1e-6:
2216            if determined_pH > pH_res:
2217                cH_res *= 1.005
2218            else:
2219                cH_res /= 1.003
2220            cH_res, cOH_res, cNa_res, cCl_res = (determine_reservoir_concentrations_selfconsistently(cH_res, 
2221                                                                                                     c_salt_res))
2222            ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2223            determined_pH = -np.log10(cH_res.to("mol/L").magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
2224        return cH_res, cOH_res, cNa_res, cCl_res
2225
2226    def enable_motion_of_rigid_object(self, instance_id, pmb_type, espresso_system):
2227        """
2228        Enables translational and rotational motion of a rigid pyMBE object instance
2229        in an ESPResSo system.This method creates a rigid-body center particle at the center of mass of
2230        the specified pyMBE object and attaches all constituent particles to it
2231        using ESPResSo virtual sites. The resulting rigid object can translate and
2232        rotate as a single body.
2233
2234        Args:
2235            instance_id ('int'):
2236                Instance ID of the pyMBE object whose rigid-body motion is enabled.
2237
2238            pmb_type ('str'):
2239                pyMBE object type of the instance (e.g. '"molecule"', '"peptide"',
2240                '"protein"', or any assembly-like type).
2241
2242            espresso_system ('espressomd.system.System'):
2243                ESPResSo system in which the rigid object is defined.
2244
2245        Notess:
2246            - This method requires ESPResSo to be compiled with the following
2247            features enabled:
2248                - '"VIRTUAL_SITES_RELATIVE"'
2249                - '"MASS"'
2250            - A new ESPResSo particle is created to represent the rigid-body center.
2251            - The mass of the rigid-body center is set to the number of particles
2252            belonging to the object.
2253            - The rotational inertia tensor is approximated from the squared
2254            distances of the particles to the center of mass.
2255        """
2256        logging.info('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
2257        inst = self.db.get_instance(pmb_type=pmb_type,
2258                                    instance_id=instance_id)
2259        label = self._get_label_id_map(pmb_type=pmb_type)
2260        particle_ids_list = self.get_particle_id_map(object_name=inst.name)[label][instance_id]
2261        center_of_mass = self.calculate_center_of_mass (instance_id=instance_id,
2262                                                        espresso_system=espresso_system,
2263                                                        pmb_type=pmb_type)
2264        rigid_object_center = espresso_system.part.add(pos=center_of_mass,
2265                                                        rotation=[True,True,True], 
2266                                                        type=self.propose_unused_type())
2267        rigid_object_center.mass = len(particle_ids_list)
2268        momI = 0
2269        for pid in particle_ids_list:
2270            momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
2271        rigid_object_center.rinertia = np.ones(3) * momI        
2272        for particle_id in particle_ids_list:
2273            pid = espresso_system.part.by_id(particle_id)
2274            pid.vs_auto_relate_to(rigid_object_center.id)
2275
2276    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2277        """
2278        Generates random coordinates outside a sphere and inside a larger bounding sphere.
2279
2280        Args:
2281            center ('array-like'):
2282                Coordinates of the center of the spheres.
2283
2284            radius ('float'):
2285                Radius of the inner exclusion sphere. Must be positive.
2286
2287            max_dist ('float'):
2288                Radius of the outer sampling sphere. Must be larger than 'radius'.
2289
2290            n_samples ('int'):
2291                Number of coordinates to generate.
2292
2293        Returns:
2294            'list' of 'numpy.ndarray':
2295                List of coordinates lying outside the inner sphere and inside the
2296                outer sphere.
2297
2298        Notess:
2299            - Points are uniformly sampled inside a sphere of radius 'max_dist' centered at 'center' 
2300            and only those with a distance greater than or equal to 'radius' from the center are retained.
2301        """
2302        if not radius > 0: 
2303            raise ValueError (f'The value of {radius} must be a positive value')
2304        if not radius < max_dist:
2305            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
2306        coord_list = []
2307        counter = 0
2308        while counter<n_samples:
2309            coord = self.generate_random_points_in_a_sphere(center=center, 
2310                                            radius=max_dist,
2311                                            n_samples=1)[0]
2312            if np.linalg.norm(coord-np.asarray(center))>=radius:
2313                coord_list.append (coord)
2314                counter += 1
2315        return coord_list
2316    
2317    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2318        """
2319        Generates uniformly distributed random points inside or on the surface of a sphere.
2320
2321        Args:
2322            center ('array-like'):
2323                Coordinates of the center of the sphere.
2324
2325            radius ('float'):
2326                Radius of the sphere.
2327
2328            n_samples ('int'):
2329                Number of sample points to generate.
2330
2331            on_surface ('bool', optional):
2332                If True, points are uniformly sampled on the surface of the sphere.
2333                If False, points are uniformly sampled within the sphere volume.
2334                Defaults to False.
2335
2336        Returns:
2337            'numpy.ndarray':
2338                Array of shape '(n_samples, d)' containing the generated coordinates,
2339                where 'd' is the dimensionality of 'center'.
2340        Notes:
2341            - Points are sampled in a space whose dimensionality is inferred 
2342            from the length of 'center'.
2343        """
2344        # initial values
2345        center=np.array(center)
2346        d = center.shape[0]
2347        # sample n_samples points in d dimensions from a standard normal distribution
2348        samples = self.rng.normal(size=(n_samples, d))
2349        # make the samples lie on the surface of the unit hypersphere
2350        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
2351        samples /= normalize_radii
2352        if not on_surface:
2353            # make the samples lie inside the hypersphere with the correct density
2354            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
2355            new_radii = np.power(uniform_points, 1/d)
2356            samples *= new_radii
2357        # scale the points to have the correct radius and center
2358        samples = samples * radius + center
2359        return samples 
2360
2361    def generate_trial_perpendicular_vector(self,vector,magnitude):
2362        """
2363        Generates a random vector perpendicular to a given vector.
2364
2365        Args:
2366            vector ('array-like'):
2367                Reference vector to which the generated vector will be perpendicular.
2368
2369            magnitude ('float'):
2370                Desired magnitude of the perpendicular vector.
2371
2372        Returns:
2373            'numpy.ndarray':
2374                Vector orthogonal to 'vector' with norm equal to 'magnitude'.
2375        """ 
2376        np_vec = np.array(vector) 
2377        if np.all(np_vec == 0):
2378            raise ValueError('Zero vector')
2379        np_vec /= np.linalg.norm(np_vec) 
2380        # Generate a random vector 
2381        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
2382                                                                center=[0,0,0],
2383                                                                n_samples=1, 
2384                                                                on_surface=True)[0]
2385        # Project the random vector onto the input vector and subtract the projection
2386        projection = np.dot(random_vector, np_vec) * np_vec
2387        perpendicular_vector = random_vector - projection
2388        # Normalize the perpendicular vector to have the same magnitude as the input vector
2389        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
2390        return perpendicular_vector*magnitude
2391            
2392    def get_bond_template(self, particle_name1, particle_name2, use_default_bond=False) :
2393        """
2394        Retrieves a bond template connecting two particle templates.
2395
2396        Args:
2397            particle_name1 ('str'):
2398                Name of the first particle template.
2399
2400            particle_name2 ('str'):
2401                Name of the second particle template.
2402
2403            use_default_bond ('bool', optional):
2404                If True, returns the default bond template when no specific bond
2405                template is found. Defaults to False.
2406
2407        Returns:
2408            'BondTemplate':
2409                Bond template object retrieved from the pyMBE database.
2410            
2411        Notes:
2412            - This method searches the pyMBE database for a bond template defined between particle templates with names 'particle_name1' and 'particle_name2'. 
2413            - If no specific bond template is found and 'use_default_bond' is enabled, a default bond template is returned instead.
2414        """
2415        # Try to find a specific bond template
2416        bond_key = BondTemplate.make_bond_key(pn1=particle_name1,
2417                                              pn2=particle_name2)
2418        try:
2419            return self.db.get_template(name=bond_key, 
2420                                        pmb_type="bond")
2421        except ValueError:
2422            pass
2423
2424        #  Fallback to default bond if allowed
2425        if use_default_bond:
2426            return self.db.get_template(name="default", 
2427                                        pmb_type="bond")
2428
2429        # No bond template found
2430        raise ValueError(f"No bond template found between '{particle_name1}' and '{particle_name2}', and default bonds are deactivated.")
2431    
2432    def get_charge_number_map(self):
2433        """
2434        Construct a mapping from ESPResSo particle types to their charge numbers.
2435
2436        Returns:
2437            'dict[int, float]':
2438                Dictionary mapping ESPResSo particle types to charge numbers,
2439                ''{es_type: z}''.
2440
2441        Notess:
2442            - The mapping is built from particle *states*, not instances.
2443            - If multiple templates define states with the same ''es_type'',
2444            the last encountered definition will overwrite previous ones.
2445            This behavior is intentional and assumes database consistency.
2446            - Neutral particles (''z = 0'') are included in the map.
2447        """
2448        charge_number_map = {}
2449        particle_templates = self.db.get_templates("particle")
2450        for tpl in particle_templates.values():
2451            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
2452                charge_number_map[state.es_type] = state.z
2453        return charge_number_map
2454
2455    def get_instances_df(self, pmb_type):
2456        """
2457        Returns a dataframe with all instances of type 'pmb_type' in the pyMBE database.
2458
2459        Args:
2460            pmb_type ('str'): 
2461                pmb type to search instances in the pyMBE database.
2462        
2463        Returns:
2464            ('Pandas.Dataframe'): 
2465                Dataframe with all instances of type 'pmb_type'.
2466        """
2467        return self.db._get_instances_df(pmb_type=pmb_type)
2468
2469    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2470        """
2471        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2472        'particle_name1' and 'particle_name2' in the pyMBE database, calculated according to the provided combining rule.
2473
2474        Args:
2475            particle_name1 ('str'): 
2476                label of the type of the first particle type
2477
2478            particle_name2 ('str'): 
2479                label of the type of the second particle type
2480
2481            combining_rule ('string', optional): 
2482                combining rule used to calculate 'sigma' and 'epsilon' for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2483
2484        Returns:
2485            ('dict'):
2486                {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2487
2488        Notes:
2489            - Currently, the only 'combining_rule' supported is Lorentz-Berthelot.
2490            - If the sigma value of 'particle_name1' or 'particle_name2' is 0, the function will return an empty dictionary. No LJ interactions are set up for particles with sigma = 0.
2491        """
2492        supported_combining_rules=["Lorentz-Berthelot"]
2493        if combining_rule not in supported_combining_rules:
2494            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2495        part_tpl1 = self.db.get_template(name=particle_name1,
2496                                         pmb_type="particle")
2497        part_tpl2 = self.db.get_template(name=particle_name2,
2498                                         pmb_type="particle")
2499        lj_parameters1 = part_tpl1.get_lj_parameters(ureg=self.units)
2500        lj_parameters2 = part_tpl2.get_lj_parameters(ureg=self.units)
2501
2502        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2503        if part_tpl1.sigma.magnitude == 0 or part_tpl2.sigma.magnitude == 0:
2504            return {}
2505        # Apply combining rule
2506        if combining_rule == 'Lorentz-Berthelot':
2507            sigma=(lj_parameters1["sigma"]+lj_parameters2["sigma"])/2
2508            cutoff=(lj_parameters1["cutoff"]+lj_parameters2["cutoff"])/2
2509            offset=(lj_parameters1["offset"]+lj_parameters2["offset"])/2
2510            epsilon=np.sqrt(lj_parameters1["epsilon"]*lj_parameters2["epsilon"])
2511        return {"sigma": sigma, "cutoff": cutoff, "offset": offset, "epsilon": epsilon}    
2512
2513    def get_particle_id_map(self, object_name):
2514        """
2515        Collect all particle IDs associated with an object of given name in the
2516        pyMBE database. 
2517
2518        Args:
2519            object_name ('str'): 
2520                Name of the object.
2521
2522        Returns:
2523            ('dict'): 
2524                {"all": [particle_ids],
2525                 "residue_map": {residue_id: [particle_ids]},
2526                 "molecule_map": {molecule_id: [particle_ids]},
2527                 "assembly_map": {assembly_id: [particle_ids]},}
2528
2529        Notess:
2530            - Works for all supported pyMBE templates.
2531            - Relies in the internal method Manager.get_particle_id_map, see method for the detailed code.
2532        """
2533        return self.db.get_particle_id_map(object_name=object_name)
2534
2535    def get_pka_set(self):
2536        """
2537        Retrieve the pKa set for all titratable particles in the pyMBE database.
2538
2539        Returns:
2540            ('dict'): 
2541                Dictionary of the form:
2542                {"particle_name": {"pka_value": float,
2543                                   "acidity": "acidic" | "basic"}}
2544        Notes:
2545            - If a particle participates in multiple acid/base reactions, an error is raised.
2546        """
2547        pka_set = {}
2548        supported_reactions = ["monoprotic_acid",
2549                               "monoprotic_base"]
2550        for reaction in self.db._reactions.values():
2551            if reaction.reaction_type not in supported_reactions:
2552                continue
2553            # Identify involved particle(s)
2554            particle_names = {participant.particle_name for participant in reaction.participants}
2555            particle_name = particle_names.pop()
2556            if particle_name in pka_set:
2557                raise ValueError(f"Multiple acid/base reactions found for particle '{particle_name}'.")
2558            pka_set[particle_name] = {"pka_value": reaction.pK}
2559            if reaction.reaction_type == "monoprotic_acid":
2560                acidity = "acidic"
2561            elif reaction.reaction_type == "monoprotic_base":
2562                acidity = "basic"
2563            pka_set[particle_name]["acidity"] = acidity
2564        return pka_set
2565    
2566    def get_radius_map(self, dimensionless=True):
2567        """
2568        Gets the effective radius of each particle defined in the pyMBE database. 
2569
2570        Args:
2571            dimensionless ('bool'):
2572                If ``True``, return magnitudes expressed in ``reduced_length``.
2573                If ``False``, return Pint quantities with units.
2574        
2575        Returns:
2576            ('dict'): 
2577                {espresso_type: radius}.
2578
2579        Notes:
2580            - The radius corresponds to (sigma+offset)/2
2581        """
2582        if "particle" not in self.db._templates:
2583            return {}          
2584        result = {}
2585        for _, tpl in self.db._templates["particle"].items():
2586            radius = (tpl.sigma.to_quantity(self.units) + tpl.offset.to_quantity(self.units))/2.0
2587            if dimensionless:
2588                magnitude_reduced_length = radius.m_as("reduced_length")
2589                radius = magnitude_reduced_length
2590            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
2591                result[state.es_type] = radius
2592        return result
2593
2594    def get_reactions_df(self):
2595        """
2596        Returns a dataframe with all reaction templates in the pyMBE database.
2597
2598        Returns:
2599            (Pandas.Dataframe): 
2600                Dataframe with all  reaction templates.
2601        """
2602        return self.db._get_reactions_df()
2603
2604    def get_reduced_units(self):
2605        """
2606        Returns the  current set of reduced units defined in pyMBE.
2607
2608        Returns:
2609            reduced_units_text ('str'): 
2610                text with information about the current set of reduced units.
2611
2612        """
2613        unit_length=self.units.Quantity(1,'reduced_length')
2614        unit_energy=self.units.Quantity(1,'reduced_energy')
2615        unit_charge=self.units.Quantity(1,'reduced_charge')
2616        reduced_units_text = "\n".join(["Current set of reduced units:",
2617                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2618                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2619                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2620                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2621                                        ])   
2622        return reduced_units_text
2623
2624    def get_templates_df(self, pmb_type):
2625        """
2626        Returns a dataframe with all templates of type 'pmb_type' in the pyMBE database.
2627
2628        Args:
2629            pmb_type ('str'): 
2630                pmb type to search templates in the pyMBE database.
2631        
2632        Returns:
2633            ('Pandas.Dataframe'): 
2634                Dataframe with all templates of type given by 'pmb_type'.
2635        """
2636        return self.db._get_templates_df(pmb_type=pmb_type)
2637
2638    def get_type_map(self):
2639        """
2640        Return the mapping of ESPResSo types for all particle states defined in the pyMBE database.
2641        
2642        Returns:
2643            'dict[str, int]':
2644                A dictionary mapping each particle state to its corresponding ESPResSo type:
2645                {state_name: es_type, ...}
2646        """
2647        
2648        return self.db.get_es_types_map()
2649
2650    def initialize_lattice_builder(self, diamond_lattice):
2651        """
2652        Initialize the lattice builder with the DiamondLattice object.
2653
2654        Args:
2655            diamond_lattice ('DiamondLattice'): 
2656                DiamondLattice object from the 'lib/lattice' module to be used in the LatticeBuilder.
2657        """
2658        from .lib.lattice import LatticeBuilder, DiamondLattice
2659        if not isinstance(diamond_lattice, DiamondLattice):
2660            raise TypeError("Currently only DiamondLattice objects are supported.")
2661        self.lattice_builder = LatticeBuilder(lattice=diamond_lattice)
2662        logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}")
2663        return self.lattice_builder
2664
2665    def load_database(self, folder, format='csv'):
2666        """
2667        Loads a pyMBE database stored in 'folder'.
2668
2669        Args:
2670            folder ('str' or 'Path'): 
2671                Path to the folder where the pyMBE database was stored.
2672
2673            format ('str', optional): 
2674                Format of the database to be loaded. Defaults to 'csv'.
2675
2676        Return:
2677            ('dict'): 
2678                metadata with additional information about the source of the information in the database.
2679
2680        Notes:
2681            - The folder must contain the files generated by 'pmb.save_database()'.
2682            - Currently, only 'csv' format is supported.
2683        """
2684        supported_formats = ['csv']
2685        if format not in supported_formats:
2686            raise ValueError(f"Format {format} not supported. Supported formats are {supported_formats}")
2687        if format == 'csv':
2688            metadata =io._load_database_csv(self.db, 
2689                                            folder=folder)
2690        return metadata
2691        
2692    
2693    def load_pka_set(self, filename):
2694        """
2695        Load a pKa set and attach chemical states and acid–base reactions
2696        to existing particle templates.
2697
2698        Args:
2699            filename ('str'): 
2700                Path to a JSON file containing the pKa set. Expected format:
2701                {"metadata": {...},
2702                  "data": {"A": {"acidity": "acidic", "pka_value": 4.5},
2703                           "B": {"acidity": "basic",  "pka_value": 9.8}}}
2704
2705        Returns:
2706            ('dict'): 
2707                Dictionary with bibliographic metadata about the original work were the pKa set was determined.
2708
2709        Notes:
2710            - This method is designed for monoprotic acids and bases only.
2711        """
2712        with open(filename, "r") as f:
2713            pka_data = json.load(f)
2714        pka_set = pka_data["data"]
2715        metadata = pka_data.get("metadata", {})
2716        self._check_pka_set(pka_set)
2717        for particle_name, entry in pka_set.items():
2718            acidity = entry["acidity"]
2719            pka = entry["pka_value"]
2720            self.define_monoprototic_acidbase_reaction(particle_name=particle_name,
2721                                                       pka=pka,
2722                                                       acidity=acidity,
2723                                                       metadata=metadata)
2724        return metadata
2725            
2726    def propose_unused_type(self):
2727        """
2728        Propose an unused ESPResSo particle type.
2729
2730        Returns:
2731            ('int'): 
2732                The next available integer ESPResSo type. Returns ''0'' if no integer types are currently defined.
2733        """
2734        type_map = self.get_type_map()
2735        # Flatten all es_type values across all particles and states
2736        all_types = []
2737        for es_type in type_map.values():
2738            all_types.append(es_type)
2739        # If no es_types exist, start at 0
2740        if not all_types:
2741            return 0
2742        return max(all_types) + 1
2743       
2744    def read_protein_vtf(self, filename, unit_length=None):
2745        """
2746        Loads a coarse-grained protein model from a VTF file.
2747
2748        Args:
2749            filename ('str'): 
2750                Path to the VTF file.
2751                
2752            unit_length ('Pint.Quantity'): 
2753                Unit of length for coordinates (pyMBE UnitRegistry). Defaults to Angstrom.
2754
2755        Returns:
2756            ('tuple'):
2757                ('dict'): Particle topology.    
2758                ('str'):  One-letter amino-acid sequence (including n/c ends).
2759        """
2760        logging.info(f"Loading protein coarse-grain model file: {filename}")
2761        if unit_length is None:
2762            unit_length = 1 * self.units.angstrom
2763        atoms = {}        # atom_id -> atom info
2764        coords = []       # ordered coordinates
2765        residues = {}     # resid -> resname (first occurrence)
2766        has_n_term = False
2767        has_c_term = False
2768        aa_3to1 = {"ALA": "A", "ARG": "R", "ASN": "N", "ASP": "D",
2769                   "CYS": "C", "GLU": "E", "GLN": "Q", "GLY": "G",
2770                   "HIS": "H", "ILE": "I", "LEU": "L", "LYS": "K",
2771                   "MET": "M", "PHE": "F", "PRO": "P", "SER": "S",
2772                   "THR": "T", "TRP": "W", "TYR": "Y", "VAL": "V",
2773                   "n": "n", "c": "c"}
2774        # --- parse VTF ---
2775        with open(filename, "r") as f:
2776            for line in f:
2777                fields = line.split()
2778                if not fields:
2779                    continue
2780                if fields[0] == "atom":
2781                    atom_id = int(fields[1])
2782                    atom_name = fields[3]
2783                    resname = fields[5]
2784                    resid = int(fields[7])
2785                    chain_id = fields[9]
2786                    radius = float(fields[11]) * unit_length
2787                    atoms[atom_id] = {"name": atom_name,
2788                                     "resname": resname,
2789                                     "resid": resid,
2790                                     "chain_id": chain_id,
2791                                     "radius": radius}
2792                    if resname == "n":
2793                        has_n_term = True
2794                    elif resname == "c":
2795                        has_c_term = True
2796                    # register residue 
2797                    if resid not in residues:
2798                        residues[resid] = resname
2799                elif fields[0].isnumeric():
2800                    xyz = [(float(x) * unit_length).to("reduced_length").magnitude
2801                        for x in fields[1:4]]
2802                    coords.append(xyz)
2803        sequence = ""
2804        # N-terminus
2805        if has_n_term:
2806            sequence += "n"
2807        # protein residues only
2808        protein_resids = sorted(resid for resid, resname in residues.items()  if resname not in ("n", "c", "Ca"))
2809        for resid in protein_resids:
2810            resname = residues[resid]
2811            try:
2812                sequence += aa_3to1[resname]
2813            except KeyError:
2814                raise ValueError(f"Unknown residue name '{resname}' in VTF file")
2815        # C-terminus
2816        if has_c_term:
2817            sequence += "c"
2818        last_resid = max(protein_resids)
2819        # --- build topology ---
2820        topology_dict = {}
2821        for atom_id in sorted(atoms.keys()):
2822            atom = atoms[atom_id]
2823            resname = atom["resname"]
2824            resid = atom["resid"]
2825            # apply labeling rules
2826            if resname == "n":
2827                label_resid = 0
2828            elif resname == "c":
2829                label_resid = last_resid + 1
2830            elif resname == "Ca":
2831                label_resid = last_resid + 2
2832            else:
2833                label_resid = resid  # preserve original resid 
2834            label = f"{atom['name']}{label_resid}"
2835            if label in topology_dict:
2836                raise ValueError(f"Duplicate particle label '{label}'. Check VTF residue definitions.")
2837            topology_dict[label] = {"initial_pos": coords[atom_id - 1], "chain_id": atom["chain_id"], "radius": atom["radius"],}
2838        return topology_dict, sequence
2839
2840    
2841    def save_database(self, folder, format='csv'):
2842        """
2843        Saves the current pyMBE database into a file 'filename'.
2844
2845        Args:
2846            folder ('str' or 'Path'): 
2847                Path to the folder where the database files will be saved.
2848
2849        """
2850        supported_formats = ['csv']
2851        if format not in supported_formats:
2852            raise ValueError(f"Format {format} not supported. Supported formats are: {supported_formats}")
2853        if format == 'csv':
2854            io._save_database_csv(self.db, 
2855                                folder=folder)
2856
2857    def set_particle_initial_state(self, particle_name, state_name):
2858        """
2859        Sets the default initial state of a particle template defined in the pyMBE database.
2860
2861        Args:
2862            particle_name ('str'): 
2863                Unique label that identifies the particle template. 
2864
2865            state_name ('str'): 
2866                Name of the state to be set as default initial state.
2867        """
2868        part_tpl = self.db.get_template(name=particle_name,
2869        
2870                                        pmb_type="particle")
2871        part_tpl.initial_state = state_name
2872        logging.info(f"Default initial state of particle {particle_name} set to {state_name}.")
2873
2874    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None):
2875        """
2876        Sets the set of reduced units used by pyMBE.units and it prints it.
2877
2878        Args:
2879            unit_length ('pint.Quantity', optional): 
2880                Reduced unit of length defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2881
2882            unit_charge ('pint.Quantity', optional): 
2883                Reduced unit of charge defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2884
2885            temperature ('pint.Quantity', optional): 
2886                Temperature of the system, defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2887
2888            Kw ('pint.Quantity', optional): 
2889                Ionic product of water in mol^2/l^2. Defaults to None. 
2890
2891        Notes:
2892            - If no 'temperature' is given, a value of 298.15 K is assumed by default.
2893            - If no 'unit_length' is given, a value of 0.355 nm is assumed by default.
2894            - If no 'unit_charge' is given, a value of 1 elementary charge is assumed by default. 
2895            - If no 'Kw' is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2896        """
2897        if unit_length is None:
2898            unit_length= 0.355*self.units.nm
2899        if temperature is None:
2900            temperature = 298.15 * self.units.K
2901        if unit_charge is None:
2902            unit_charge = scipy.constants.e * self.units.C
2903        if Kw is None:
2904            Kw = 1e-14
2905        # Sanity check
2906        variables=[unit_length,temperature,unit_charge]
2907        dimensionalities=["[length]","[temperature]","[charge]"]
2908        for variable,dimensionality in zip(variables,dimensionalities):
2909            self._check_dimensionality(variable,dimensionality)
2910        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2911        self.kT=temperature*self.kB
2912        self.units._build_cache()
2913        self.units.define(f'reduced_energy = {self.kT} ')
2914        self.units.define(f'reduced_length = {unit_length}')
2915        self.units.define(f'reduced_charge = {unit_charge}')
2916        logging.info(self.get_reduced_units())
2917
2918    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, use_exclusion_radius_per_type = False):
2919        """
2920        Sets up the Acid/Base reactions for acidic/basic particles defined in the pyMBE database
2921        to be sampled in the constant pH ensemble. 
2922
2923        Args:
2924            counter_ion ('str'): 
2925                'name' of the counter_ion 'particle'.
2926
2927            constant_pH ('float'): 
2928                pH-value.
2929
2930            exclusion_range ('pint.Quantity', optional): 
2931                Below this value, no particles will be inserted.
2932
2933            use_exclusion_radius_per_type ('bool', optional): 
2934                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
2935
2936        Returns:
2937            ('reaction_methods.ConstantpHEnsemble'): 
2938                Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2939        """
2940        from espressomd import reaction_methods
2941        if exclusion_range is None:
2942            exclusion_range = max(self.get_radius_map().values())*2.0
2943        if use_exclusion_radius_per_type:
2944            exclusion_radius_per_type = self.get_radius_map()
2945        else:
2946            exclusion_radius_per_type = {}
2947        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2948                                                exclusion_range=exclusion_range, 
2949                                                seed=self.seed, 
2950                                                constant_pH=constant_pH,
2951                                                exclusion_radius_per_type = exclusion_radius_per_type)
2952        conterion_tpl = self.db.get_template(name=counter_ion,
2953                                             pmb_type="particle")
2954        conterion_state = self.db.get_template(name=conterion_tpl.initial_state,
2955                                               pmb_type="particle_state")
2956        for reaction in self.db.get_reactions():
2957            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
2958                continue
2959            default_charges = {}
2960            reactant_types  = []
2961            product_types   = []
2962            for participant in reaction.participants:
2963                state_tpl = self.db.get_template(name=participant.state_name,
2964                                                 pmb_type="particle_state")
2965                default_charges[state_tpl.es_type] = state_tpl.z
2966                if participant.coefficient < 0:
2967                    reactant_types.append(state_tpl.es_type)
2968                elif participant.coefficient > 0:
2969                    product_types.append(state_tpl.es_type)
2970            # Add counterion to the products
2971            if conterion_state.es_type not in product_types:
2972                product_types.append(conterion_state.es_type)
2973                default_charges[conterion_state.es_type] = conterion_state.z
2974                reaction.add_participant(particle_name=counter_ion,
2975                                         state_name=conterion_tpl.initial_state,
2976                                         coefficient=1)
2977            gamma=10**-reaction.pK
2978            RE.add_reaction(gamma=gamma,
2979                            reactant_types=reactant_types,
2980                            product_types=product_types,
2981                            default_charges=default_charges)
2982            reaction.add_simulation_method(simulation_method="cpH")
2983        return RE
2984
2985    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2986        """
2987        Sets up grand-canonical coupling to a reservoir of salt.
2988        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2989
2990        Args:
2991            c_salt_res ('pint.Quantity'): 
2992                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2993
2994            salt_cation_name ('str'): 
2995                Name of the salt cation (e.g. Na+) particle.
2996
2997            salt_anion_name ('str'): 
2998                Name of the salt anion (e.g. Cl-) particle.
2999
3000            activity_coefficient ('callable'): 
3001                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3002
3003            exclusion_range('pint.Quantity', optional): 
3004                For distances shorter than this value, no particles will be inserted.
3005
3006            use_exclusion_radius_per_type('bool',optional): 
3007                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
3008
3009        Returns:
3010            ('reaction_methods.ReactionEnsemble'): 
3011                Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
3012        """
3013        from espressomd import reaction_methods
3014        if exclusion_range is None:
3015            exclusion_range = max(self.get_radius_map().values())*2.0
3016        if use_exclusion_radius_per_type:
3017            exclusion_radius_per_type = self.get_radius_map()
3018        else:
3019            exclusion_radius_per_type = {}
3020        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3021                                               exclusion_range=exclusion_range, 
3022                                               seed=self.seed, 
3023                                               exclusion_radius_per_type = exclusion_radius_per_type)
3024        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3025        determined_activity_coefficient = activity_coefficient(c_salt_res)
3026        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
3027        cation_tpl = self.db.get_template(pmb_type="particle",
3028                                          name=salt_cation_name)
3029        cation_state = self.db.get_template(pmb_type="particle_state",
3030                                            name=cation_tpl.initial_state)
3031        anion_tpl = self.db.get_template(pmb_type="particle",
3032                                          name=salt_anion_name)
3033        anion_state = self.db.get_template(pmb_type="particle_state",
3034                                            name=anion_tpl.initial_state)
3035        salt_cation_es_type = cation_state.es_type
3036        salt_anion_es_type = anion_state.es_type     
3037        salt_cation_charge = cation_state.z
3038        salt_anion_charge = anion_state.z
3039        if salt_cation_charge <= 0:
3040            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3041        if salt_anion_charge >= 0:
3042            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3043        # Grand-canonical coupling to the reservoir
3044        RE.add_reaction(gamma = K_salt.magnitude,
3045                        reactant_types = [],
3046                        reactant_coefficients = [],
3047                        product_types = [ salt_cation_es_type, salt_anion_es_type ],
3048                        product_coefficients = [ 1, 1 ],
3049                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3050                                           salt_anion_es_type: salt_anion_charge})
3051        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3052                                                            state_name=cation_state.name,
3053                                                            coefficient=1),
3054                                        ReactionParticipant(particle_name=salt_anion_name,
3055                                                            state_name=anion_state.name,
3056                                                            coefficient=1)],
3057                           pK=-np.log10(K_salt.magnitude),
3058                           reaction_type="ion_insertion",
3059                           simulation_method="GCMC")
3060        self.db._register_reaction(rx_tpl)
3061        return RE
3062
3063    def setup_grxmc_reactions(self, pH_res, c_salt_res, proton_name, hydroxide_name, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
3064        """
3065        Sets up acid/base reactions for acidic/basic monoprotic particles defined in the pyMBE database, 
3066        as well as a grand-canonical coupling to a reservoir of small ions. 
3067        
3068        Args:
3069            pH_res ('float'): 
3070                pH-value in the reservoir.
3071
3072            c_salt_res ('pint.Quantity'): 
3073                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3074
3075            proton_name ('str'): 
3076                Name of the proton (H+) particle.
3077
3078            hydroxide_name ('str'): 
3079                Name of the hydroxide (OH-) particle.
3080
3081            salt_cation_name ('str'): 
3082                Name of the salt cation (e.g. Na+) particle.
3083
3084            salt_anion_name ('str'): 
3085                Name of the salt anion (e.g. Cl-) particle.
3086
3087            activity_coefficient ('callable'): 
3088                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3089
3090            exclusion_range('pint.Quantity', optional): 
3091                For distances shorter than this value, no particles will be inserted.
3092
3093            use_exclusion_radius_per_type('bool', optional): 
3094                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
3095
3096        Returns:
3097            'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':
3098
3099                'reaction_methods.ReactionEnsemble':  
3100                    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.
3101                
3102                'pint.Quantity': 
3103                    Ionic strength of the reservoir (useful for calculating partition coefficients).
3104
3105        Notess:
3106            - This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
3107
3108        [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3109        """
3110        from espressomd import reaction_methods
3111        if exclusion_range is None:
3112            exclusion_range = max(self.get_radius_map().values())*2.0
3113        if use_exclusion_radius_per_type:
3114            exclusion_radius_per_type = self.get_radius_map()
3115        else:
3116            exclusion_radius_per_type = {}
3117        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3118                                               exclusion_range=exclusion_range, 
3119                                               seed=self.seed, 
3120                                               exclusion_radius_per_type = exclusion_radius_per_type)
3121        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3122        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3123        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3124        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3125        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3126        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3127        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3128        cation_tpl = self.db.get_template(pmb_type="particle",
3129                                          name=salt_cation_name)
3130        cation_state = self.db.get_template(pmb_type="particle_state",
3131                                            name=cation_tpl.initial_state)
3132        anion_tpl = self.db.get_template(pmb_type="particle",
3133                                          name=salt_anion_name)
3134        anion_state = self.db.get_template(pmb_type="particle_state",
3135                                            name=anion_tpl.initial_state)
3136        proton_tpl = self.db.get_template(pmb_type="particle",
3137                                          name=proton_name)
3138        proton_state = self.db.get_template(pmb_type="particle_state",
3139                                            name=proton_tpl.initial_state)
3140        hydroxide_tpl = self.db.get_template(pmb_type="particle",
3141                                             name=hydroxide_name)
3142        hydroxide_state = self.db.get_template(pmb_type="particle_state",
3143                                               name=hydroxide_tpl.initial_state)
3144        proton_es_type = proton_state.es_type
3145        hydroxide_es_type = hydroxide_state.es_type
3146        salt_cation_es_type = cation_state.es_type
3147        salt_anion_es_type = anion_state.es_type
3148        proton_charge = proton_state.z
3149        hydroxide_charge = hydroxide_state.z          
3150        salt_cation_charge = cation_state.z
3151        salt_anion_charge = anion_state.z      
3152        if proton_charge <= 0:
3153            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
3154        if salt_cation_charge <= 0:
3155            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3156        if hydroxide_charge >= 0:
3157            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
3158        if salt_anion_charge >= 0:
3159            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3160        # Grand-canonical coupling to the reservoir
3161        # 0 = H+ + OH-
3162        RE.add_reaction(gamma = K_W.magnitude,
3163                        reactant_types = [],
3164                        reactant_coefficients = [],
3165                        product_types = [ proton_es_type, hydroxide_es_type ],
3166                        product_coefficients = [ 1, 1 ],
3167                        default_charges = {proton_es_type: proton_charge, 
3168                                           hydroxide_es_type: hydroxide_charge})
3169        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3170                                                            state_name=proton_state.name,
3171                                                            coefficient=1),
3172                                        ReactionParticipant(particle_name=hydroxide_name,
3173                                                            state_name=hydroxide_state.name,
3174                                                            coefficient=1)],
3175                           pK=-np.log10(K_W.magnitude),
3176                           reaction_type="ion_insertion",
3177                           simulation_method="GRxMC")
3178        self.db._register_reaction(rx_tpl)
3179        # 0 = Na+ + Cl-
3180        RE.add_reaction(gamma = K_NACL.magnitude,
3181                        reactant_types = [],
3182                        reactant_coefficients = [],
3183                        product_types = [ salt_cation_es_type, salt_anion_es_type ],
3184                        product_coefficients = [ 1, 1 ],
3185                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3186                                        salt_anion_es_type: salt_anion_charge})
3187        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3188                                                            state_name=cation_state.name,
3189                                                            coefficient=1),
3190                                        ReactionParticipant(particle_name=salt_anion_name,
3191                                                            state_name=anion_state.name,
3192                                                            coefficient=1)],
3193                           pK=-np.log10(K_NACL.magnitude),
3194                           reaction_type="ion_insertion",
3195                           simulation_method="GRxMC")
3196        self.db._register_reaction(rx_tpl)
3197        # 0 = Na+ + OH-
3198        RE.add_reaction(gamma = (K_NACL * K_W / K_HCL).magnitude,
3199                        reactant_types = [],
3200                        reactant_coefficients = [],
3201                        product_types = [ salt_cation_es_type, hydroxide_es_type ],
3202                        product_coefficients = [ 1, 1 ],
3203                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3204                                           hydroxide_es_type: hydroxide_charge})
3205        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3206                                                            state_name=cation_state.name,
3207                                                            coefficient=1),
3208                                        ReactionParticipant(particle_name=hydroxide_name,
3209                                                            state_name=hydroxide_state.name,
3210                                                            coefficient=1)],
3211                           pK=-np.log10((K_NACL * K_W / K_HCL).magnitude),
3212                           reaction_type="ion_insertion",
3213                           simulation_method="GRxMC")
3214        self.db._register_reaction(rx_tpl)
3215        # 0 = H+ + Cl-
3216        RE.add_reaction(gamma = K_HCL.magnitude,
3217                        reactant_types = [],
3218                        reactant_coefficients = [],
3219                        product_types = [ proton_es_type, salt_anion_es_type ],
3220                        product_coefficients = [ 1, 1 ],
3221                        default_charges = {proton_es_type: proton_charge, 
3222                                           salt_anion_es_type: salt_anion_charge})
3223        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3224                                                            state_name=proton_state.name,
3225                                                            coefficient=1),
3226                                        ReactionParticipant(particle_name=salt_anion_name,
3227                                                            state_name=anion_state.name,
3228                                                            coefficient=1)],
3229                           pK=-np.log10(K_HCL.magnitude),
3230                           reaction_type="ion_insertion",
3231                           simulation_method="GRxMC")
3232        self.db._register_reaction(rx_tpl)
3233        # Annealing moves to ensure sufficient sampling
3234        # Cation annealing H+ = Na+
3235        RE.add_reaction(gamma = (K_NACL / K_HCL).magnitude,
3236                        reactant_types = [proton_es_type],
3237                        reactant_coefficients = [ 1 ],
3238                        product_types = [ salt_cation_es_type ],
3239                        product_coefficients = [ 1 ],
3240                        default_charges = {proton_es_type: proton_charge, 
3241                                           salt_cation_es_type: salt_cation_charge})
3242        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3243                                                            state_name=proton_state.name,
3244                                                            coefficient=-1),
3245                                        ReactionParticipant(particle_name=salt_cation_name,
3246                                                            state_name=cation_state.name,
3247                                                            coefficient=1)],
3248                           pK=-np.log10((K_NACL / K_HCL).magnitude),
3249                           reaction_type="particle replacement",
3250                           simulation_method="GRxMC")
3251        self.db._register_reaction(rx_tpl)
3252        # Anion annealing OH- = Cl- 
3253        RE.add_reaction(gamma = (K_HCL / K_W).magnitude,
3254                        reactant_types = [hydroxide_es_type],
3255                        reactant_coefficients = [ 1 ],
3256                        product_types = [ salt_anion_es_type ],
3257                        product_coefficients = [ 1 ],
3258            default_charges = {hydroxide_es_type: hydroxide_charge, 
3259                               salt_anion_es_type: salt_anion_charge})
3260        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=hydroxide_name,
3261                                                            state_name=hydroxide_state.name,
3262                                                            coefficient=-1),
3263                                        ReactionParticipant(particle_name=salt_anion_name,
3264                                                            state_name=anion_state.name,
3265                                                            coefficient=1)],
3266                           pK=-np.log10((K_HCL / K_W).magnitude),
3267                           reaction_type="particle replacement",
3268                           simulation_method="GRxMC")
3269        self.db._register_reaction(rx_tpl)
3270        for reaction in self.db.get_reactions():
3271            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
3272                continue
3273            default_charges = {}
3274            reactant_types  = []
3275            product_types   = []
3276            for participant in reaction.participants:
3277                state_tpl = self.db.get_template(name=participant.state_name,
3278                                                 pmb_type="particle_state")
3279                default_charges[state_tpl.es_type] = state_tpl.z
3280                if participant.coefficient < 0:
3281                    reactant_types.append(state_tpl.es_type)
3282                    reactant_name=state_tpl.particle_name
3283                    reactant_state_name=state_tpl.name
3284                elif participant.coefficient > 0:
3285                    product_types.append(state_tpl.es_type)
3286                    product_name=state_tpl.particle_name
3287                    product_state_name=state_tpl.name
3288
3289            Ka = (10**-reaction.pK * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3290            # Reaction in terms of proton: HA = A + H+
3291            RE.add_reaction(gamma=Ka.magnitude,
3292                            reactant_types=reactant_types,
3293                            reactant_coefficients=[1],
3294                            product_types=product_types+[proton_es_type],
3295                            product_coefficients=[1, 1],
3296                            default_charges= default_charges | {proton_es_type: proton_charge})
3297            reaction.add_participant(particle_name=proton_name,
3298                                     state_name=proton_state.name,
3299                                     coefficient=1)
3300            reaction.add_simulation_method("GRxMC")
3301            # Reaction in terms of salt cation: HA = A + Na+
3302            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
3303                            reactant_types=reactant_types,
3304                            reactant_coefficients=[1],
3305                            product_types=product_types+[salt_cation_es_type],
3306                            product_coefficients=[1, 1],
3307                            default_charges=default_charges | {salt_cation_es_type: salt_cation_charge})
3308            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3309                                                                state_name=reactant_state_name,
3310                                                                coefficient=-1),
3311                                            ReactionParticipant(particle_name=product_name,
3312                                                                state_name=product_state_name,
3313                                                                coefficient=1),
3314                                            ReactionParticipant(particle_name=salt_cation_name,
3315                                                                state_name=cation_state.name,
3316                                                                coefficient=1),],
3317                              pK=-np.log10((Ka * K_NACL / K_HCL).magnitude),
3318                              reaction_type=reaction.reaction_type+"_salt",
3319                              simulation_method="GRxMC")
3320            self.db._register_reaction(rx_tpl)
3321            # Reaction in terms of hydroxide: OH- + HA = A
3322            RE.add_reaction(gamma=(Ka / K_W).magnitude,
3323                            reactant_types=reactant_types+[hydroxide_es_type],
3324                            reactant_coefficients=[1, 1],
3325                            product_types=product_types,
3326                            product_coefficients=[1],
3327                            default_charges=default_charges | {hydroxide_es_type: hydroxide_charge})
3328            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3329                                                                state_name=reactant_state_name,
3330                                                                coefficient=-1),
3331                                            ReactionParticipant(particle_name=product_name,
3332                                                                state_name=product_state_name,
3333                                                                coefficient=1),
3334                                            ReactionParticipant(particle_name=hydroxide_name,
3335                                                                state_name=hydroxide_state.name,
3336                                                                coefficient=-1),],
3337                              pK=-np.log10((Ka / K_W).magnitude),
3338                              reaction_type=reaction.reaction_type+"_conjugate",
3339                              simulation_method="GRxMC")
3340            self.db._register_reaction(rx_tpl)
3341            # Reaction in terms of salt anion: Cl- + HA = A
3342            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
3343                            reactant_types=reactant_types+[salt_anion_es_type],
3344                            reactant_coefficients=[1, 1],
3345                            product_types=product_types,
3346                            product_coefficients=[1],
3347                            default_charges=default_charges | {salt_anion_es_type: salt_anion_charge})
3348            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3349                                                                state_name=reactant_state_name,
3350                                                                coefficient=-1),
3351                                            ReactionParticipant(particle_name=product_name,
3352                                                                state_name=product_state_name,
3353                                                                coefficient=1),
3354                                            ReactionParticipant(particle_name=salt_anion_name,
3355                                                                state_name=anion_state.name,
3356                                                                coefficient=-1),],
3357                              pK=-np.log10((Ka / K_HCL).magnitude),
3358                              reaction_type=reaction.reaction_type+"_salt",
3359                              simulation_method="GRxMC")
3360            self.db._register_reaction(rx_tpl)
3361        return RE, ionic_strength_res
3362
3363    def setup_grxmc_unified(self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
3364        """
3365        Sets up acid/base reactions for acidic/basic 'particles' defined in the pyMBE database, as well as a grand-canonical coupling to a 
3366        reservoir of small ions using a unified formulation for small ions.
3367
3368        Args:
3369            pH_res ('float'): 
3370                pH-value in the reservoir.
3371
3372            c_salt_res ('pint.Quantity'): 
3373                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3374
3375            cation_name ('str'): 
3376                Name of the cationic particle.
3377
3378            anion_name ('str'): 
3379                Name of the anionic particle.
3380
3381            activity_coefficient ('callable'): 
3382                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3383
3384            exclusion_range('pint.Quantity', optional): 
3385                Below this value, no particles will be inserted.
3386            
3387            use_exclusion_radius_per_type('bool', optional): 
3388                Controls if one exclusion_radius per each espresso_type. Defaults to 'False'.
3389
3390        Returns:
3391            'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':
3392
3393                'reaction_methods.ReactionEnsemble':  
3394                    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.
3395                
3396                'pint.Quantity': 
3397                    Ionic strength of the reservoir (useful for calculating partition coefficients).
3398
3399        Notes:
3400            - This implementation uses the formulation of the grand-reaction method by Curk et al. [1], which relies on "unified" ion types X+ = {H+, Na+} and X- = {OH-, Cl-}. 
3401            - A function that implements the original version of the grand-reaction method by Landsgesell et al. [2] is also available under the name 'setup_grxmc_reactions'.
3402
3403        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3404        [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3405        """
3406        from espressomd import reaction_methods
3407        if exclusion_range is None:
3408            exclusion_range = max(self.get_radius_map().values())*2.0
3409        if use_exclusion_radius_per_type:
3410            exclusion_radius_per_type = self.get_radius_map()
3411        else:
3412            exclusion_radius_per_type = {}
3413        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3414                                               exclusion_range=exclusion_range, 
3415                                               seed=self.seed, 
3416                                               exclusion_radius_per_type = exclusion_radius_per_type)
3417        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3418        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3419        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3420        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3421        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3422        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3423        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3424        K_XX = a_cation * a_anion
3425        cation_tpl = self.db.get_template(pmb_type="particle",
3426                                          name=cation_name)
3427        cation_state = self.db.get_template(pmb_type="particle_state",
3428                                            name=cation_tpl.initial_state)
3429        anion_tpl = self.db.get_template(pmb_type="particle",
3430                                          name=anion_name)
3431        anion_state = self.db.get_template(pmb_type="particle_state",
3432                                            name=anion_tpl.initial_state)
3433        cation_es_type = cation_state.es_type
3434        anion_es_type = anion_state.es_type     
3435        cation_charge = cation_state.z
3436        anion_charge = anion_state.z
3437        if cation_charge <= 0:
3438            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3439        if anion_charge >= 0:
3440            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3441        # Coupling to the reservoir: 0 = X+ + X-
3442        RE.add_reaction(gamma = K_XX.magnitude,
3443                        reactant_types = [],
3444                        reactant_coefficients = [],
3445                        product_types = [ cation_es_type, anion_es_type ],
3446                        product_coefficients = [ 1, 1 ],
3447                        default_charges = {cation_es_type: cation_charge, 
3448                                           anion_es_type: anion_charge})
3449        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=cation_name,
3450                                                            state_name=cation_state.name,
3451                                                            coefficient=1),
3452                                        ReactionParticipant(particle_name=anion_name,
3453                                                            state_name=anion_state.name,
3454                                                            coefficient=1)],
3455                           pK=-np.log10(K_XX.magnitude),
3456                           reaction_type="ion_insertion",
3457                           simulation_method="GCMC")
3458        self.db._register_reaction(rx_tpl)
3459        for reaction in self.db.get_reactions():
3460            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
3461                continue
3462            default_charges = {}
3463            reactant_types  = []
3464            product_types   = []
3465            for participant in reaction.participants:
3466                state_tpl = self.db.get_template(name=participant.state_name,
3467                                                 pmb_type="particle_state")
3468                default_charges[state_tpl.es_type] = state_tpl.z
3469                if participant.coefficient < 0:
3470                    reactant_types.append(state_tpl.es_type)
3471                    reactant_name=state_tpl.particle_name
3472                    reactant_state_name=state_tpl.name
3473                elif participant.coefficient > 0:
3474                    product_types.append(state_tpl.es_type)
3475                    product_name=state_tpl.particle_name
3476                    product_state_name=state_tpl.name
3477
3478            Ka = (10**-reaction.pK * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3479            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3480            # Reaction in terms of small cation: HA = A + X+
3481            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3482                            reactant_types=reactant_types,
3483                            reactant_coefficients=[1],
3484                            product_types=product_types+[cation_es_type],
3485                            product_coefficients=[1, 1],
3486                            default_charges=default_charges|{cation_es_type: cation_charge})
3487            reaction.add_participant(particle_name=cation_name,
3488                                     state_name=cation_state.name,
3489                                     coefficient=1)
3490            reaction.add_simulation_method("GRxMC")
3491            # Reaction in terms of small anion: X- + HA = A
3492            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3493                            reactant_types=reactant_types+[anion_es_type],
3494                            reactant_coefficients=[1, 1],
3495                            product_types=product_types,
3496                            product_coefficients=[1],
3497                            default_charges=default_charges|{anion_es_type: anion_charge})
3498            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3499                                                                state_name=reactant_state_name,
3500                                                                coefficient=-1),
3501                                            ReactionParticipant(particle_name=product_name,
3502                                                                state_name=product_state_name,
3503                                                                coefficient=1),
3504                                            ReactionParticipant(particle_name=anion_name,
3505                                                                state_name=anion_state.name,
3506                                                                coefficient=-1),],
3507                              pK=-np.log10(gamma_K_AX.magnitude / K_XX.magnitude),
3508                              reaction_type=reaction.reaction_type+"_conjugate",
3509                              simulation_method="GRxMC")
3510            self.db._register_reaction(rx_tpl)
3511        return RE, ionic_strength_res
3512
3513    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot'):
3514        """
3515        Sets up the Lennard-Jones (LJ) potential between all pairs of particle states defined in the pyMBE database.
3516
3517        Args:
3518            espresso_system('espressomd.system.System'): 
3519                Instance of a system object from the espressomd library.
3520
3521            shift_potential('bool', optional): 
3522                If True, a shift will be automatically computed such that the potential is continuous at the cutoff radius. Otherwise, no shift will be applied. Defaults to True.
3523
3524            combining_rule('string', optional): 
3525                combining rule used to calculate 'sigma' and 'epsilon' for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3526
3527            warning('bool', optional): 
3528                switch to activate/deactivate warning messages. Defaults to True.
3529
3530        Notes:
3531            - Currently, the only 'combining_rule' supported is Lorentz-Berthelot.
3532            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3533
3534        """
3535        from itertools import combinations_with_replacement
3536        particle_templates = self.db.get_templates("particle")
3537        shift = "auto" if shift_potential else 0
3538        if shift == "auto":
3539            shift_tpl = shift
3540        else:
3541            shift_tpl = PintQuantity.from_quantity(q=shift*self.units.reduced_length,
3542                                                   expected_dimension="length",
3543                                                   ureg=self.units)
3544        # Get all particle states registered in pyMBE
3545        state_entries = []
3546        for tpl in particle_templates.values():
3547            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
3548                state_entries.append((tpl, state))
3549
3550        # Iterate over all unique state pairs
3551        for (tpl1, state1), (tpl2, state2) in combinations_with_replacement(state_entries, 2):
3552
3553            lj_parameters = self.get_lj_parameters(particle_name1=tpl1.name,
3554                                                   particle_name2=tpl2.name,
3555                                                   combining_rule=combining_rule)
3556            if not lj_parameters:
3557                continue
3558
3559            espresso_system.non_bonded_inter[state1.es_type, state2.es_type].lennard_jones.set_params(
3560                epsilon=lj_parameters["epsilon"].to("reduced_energy").magnitude,
3561                sigma=lj_parameters["sigma"].to("reduced_length").magnitude,
3562                cutoff=lj_parameters["cutoff"].to("reduced_length").magnitude,
3563                offset=lj_parameters["offset"].to("reduced_length").magnitude,
3564                shift=shift)
3565                
3566            lj_template = LJInteractionTemplate(state1=state1.name,
3567                                                state2=state2.name,
3568                                                sigma=PintQuantity.from_quantity(q=lj_parameters["sigma"],
3569                                                                                 expected_dimension="length",
3570                                                                                 ureg=self.units),
3571                                                epsilon=PintQuantity.from_quantity(q=lj_parameters["epsilon"],
3572                                                                                   expected_dimension="energy",
3573                                                                                   ureg=self.units),
3574                                                cutoff=PintQuantity.from_quantity(q=lj_parameters["cutoff"],
3575                                                                                  expected_dimension="length",
3576                                                                                  ureg=self.units),
3577                                                offset=PintQuantity.from_quantity(q=lj_parameters["offset"],
3578                                                                                  expected_dimension="length",
3579                                                                                  ureg=self.units),
3580                                                shift=shift_tpl)
3581            self.db._register_template(lj_template)

Core library of the Molecular Builder for ESPResSo (pyMBE).

Attributes:
  • N_A ('pint.Quantity'): Avogadro number.
  • kB ('pint.Quantity'): Boltzmann constant.
  • e ('pint.Quantity'): Elementary charge.
  • kT ('pint.Quantity'): Thermal energy corresponding to the set temperature.
  • Kw ('pint.Quantity'): Ionic product of water, used in G-RxMC and Donnan-related calculations.
  • db ('Manager'): Database manager holding all pyMBE templates, instances and reactions.
  • rng ('numpy.random.Generator'): Random number generator initialized with the provided seed.
  • units ('pint.UnitRegistry'): Pint unit registry used for unit-aware calculations.
  • lattice_builder ('pyMBE.lib.lattice.LatticeBuilder'): Optional lattice builder object (initialized as ''None'').
  • root ('importlib.resources.abc.Traversable'): Root path to the pyMBE package resources.
pymbe_library(seed, temperature=None, unit_length=None, unit_charge=None, Kw=None)
 94    def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None):
 95        """
 96        Initializes the pyMBE library.
 97
 98        Args:
 99            seed ('int'):
100                Seed for the random number generator.
101
102            temperature ('pint.Quantity', optional):
103                Simulation temperature. If ''None'', defaults to 298.15 K.
104
105            unit_length ('pint.Quantity', optional):
106                Reference length for reduced units. If ''None'', defaults to
107                0.355 nm.
108
109            unit_charge ('pint.Quantity', optional):
110                Reference charge for reduced units. If ''None'', defaults to
111                one elementary charge.
112
113            Kw ('pint.Quantity', optional):
114                Ionic product of water (typically in mol²/L²). If ''None'',
115                defaults to 1e-14 mol²/L².
116        """
117        # Seed and RNG
118        self.seed=seed
119        self.rng = np.random.default_rng(seed)
120        self.units=pint.UnitRegistry()
121        self.N_A=scipy.constants.N_A / self.units.mol
122        self.kB=scipy.constants.k * self.units.J / self.units.K
123        self.e=scipy.constants.e * self.units.C
124        self.set_reduced_units(unit_length=unit_length, 
125                               unit_charge=unit_charge,
126                               temperature=temperature, 
127                               Kw=Kw)
128        
129        self.db = Manager(units=self.units)
130        self.lattice_builder = None
131        self.root = importlib.resources.files(__package__)

Initializes the pyMBE library.

Arguments:
  • seed ('int'): Seed for the random number generator.
  • temperature ('pint.Quantity', optional): Simulation temperature. If ''None'', defaults to 298.15 K.
  • unit_length ('pint.Quantity', optional): Reference length for reduced units. If ''None'', defaults to 0.355 nm.
  • unit_charge ('pint.Quantity', optional): Reference charge for reduced units. If ''None'', defaults to one elementary charge.
  • Kw ('pint.Quantity', optional): Ionic product of water (typically in mol²/L²). If ''None'', defaults to 1e-14 mol²/L².
seed
rng
units
N_A
kB
e
db
lattice_builder
root
def calculate_center_of_mass(self, instance_id, pmb_type, espresso_system):
526    def calculate_center_of_mass(self, instance_id, pmb_type, espresso_system):
527        """
528        Calculates the center of mass of a pyMBE object instance in an ESPResSo system.
529
530        Args:
531            instance_id ('int'):
532                pyMBE instance ID of the object whose center of mass is calculated.
533
534            pmb_type ('str'):
535                Type of the pyMBE object. Must correspond to a particle-aggregating
536                template type (e.g. '"molecule"', '"residue"', '"peptide"', '"protein"').
537
538            espresso_system ('espressomd.system.System'):
539                ESPResSo system containing the particle instances.
540
541        Returns:
542            ('numpy.ndarray'):
543                Array of shape '(3,)' containing the Cartesian coordinates of the
544                center of mass.
545
546        Notes:
547            - This method assumes equal mass for all particles.
548            - Periodic boundary conditions are *not* unfolded; positions are taken
549            directly from ESPResSo particle coordinates.
550        """
551        center_of_mass = np.zeros(3)
552        axis_list = [0,1,2]
553        inst = self.db.get_instance(pmb_type=pmb_type,
554                                    instance_id=instance_id)
555        particle_id_list = self.get_particle_id_map(object_name=inst.name)["all"]
556        for pid in particle_id_list:
557            for axis in axis_list:
558                center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis]
559        center_of_mass = center_of_mass / len(particle_id_list)
560        return center_of_mass

Calculates the center of mass of a pyMBE object instance in an ESPResSo system.

Arguments:
  • instance_id ('int'): pyMBE instance ID of the object whose center of mass is calculated.
  • pmb_type ('str'): Type of the pyMBE object. Must correspond to a particle-aggregating template type (e.g. '"molecule"', '"residue"', '"peptide"', '"protein"').
  • espresso_system ('espressomd.system.System'): ESPResSo system containing the particle instances.
Returns:

('numpy.ndarray'): Array of shape '(3,)' containing the Cartesian coordinates of the center of mass.

Notes:
  • This method assumes equal mass for all particles.
  • Periodic boundary conditions are not unfolded; positions are taken directly from ESPResSo particle coordinates.
def calculate_HH(self, template_name, pH_list=None, pka_set=None):
562    def calculate_HH(self, template_name, pH_list=None, pka_set=None):
563        """
564        Calculates the charge in the template object according to the ideal  Henderson–Hasselbalch titration curve.
565
566        Args:
567            template_name ('str'):
568                Name of the template.
569
570            pH_list ('list[float]', optional):
571                pH values at which the charge is evaluated.
572                Defaults to 50 values between 2 and 12.
573
574            pka_set ('dict', optional):
575                Mapping: {particle_name: {"pka_value": 'float', "acidity": "acidic"|"basic"}}
576
577        Returns:
578            'list[float]':
579                Net molecular charge at each pH value.
580        """
581        if pH_list is None:
582            pH_list = np.linspace(2, 12, 50)
583        if pka_set is None:
584            pka_set = self.get_pka_set()
585        self._check_pka_set(pka_set=pka_set)
586        particle_counts = self.db.get_particle_templates_under(template_name=template_name,
587                                                               return_counts=True)
588        if not particle_counts:
589            return [None] * len(pH_list)
590        charge_number_map = self.get_charge_number_map()
591        def formal_charge(particle_name):
592            tpl = self.db.get_template(name=particle_name, 
593                                       pmb_type="particle")
594            state = self.db.get_template(name=tpl.initial_state,
595                                         pmb_type="particle_state")
596            return charge_number_map[state.es_type]
597        Z_HH = []
598        for pH in pH_list:
599            Z = 0.0
600            for particle, multiplicity in particle_counts.items():
601                if particle in pka_set:
602                    pka = pka_set[particle]["pka_value"]
603                    acidity = pka_set[particle]["acidity"]
604                    if acidity == "acidic":
605                        psi = -1
606                    elif acidity == "basic":
607                        psi = +1
608                    else:
609                        raise ValueError(f"Unknown acidity '{acidity}' for particle '{particle}'")
610                    charge = psi / (1.0 + 10.0 ** (psi * (pH - pka)))
611                    Z += multiplicity * charge
612                else:
613                    Z += multiplicity * formal_charge(particle)
614            Z_HH.append(Z)
615        return Z_HH   

Calculates the charge in the template object according to the ideal Henderson–Hasselbalch titration curve.

Arguments:
  • template_name ('str'): Name of the template.
  • pH_list ('list[float]', optional): pH values at which the charge is evaluated. Defaults to 50 values between 2 and 12.
  • pka_set ('dict', optional): Mapping: {particle_name: {"pka_value": 'float', "acidity": "acidic"|"basic"}}
Returns:

'list[float]': Net molecular charge at each pH value.

def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
617    def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
618        """
619        Computes macromolecular charges using the Henderson–Hasselbalch equation
620        coupled to ideal Donnan partitioning.
621
622        Args:
623            c_macro ('dict'):
624                Mapping of macromolecular species names to their concentrations
625                in the system:
626                '{molecule_name: concentration}'.
627                
628            c_salt ('float' or 'pint.Quantity'):
629                Salt concentration in the reservoir.
630
631            pH_list ('list[float]', optional):
632                List of pH values in the reservoir at which the calculation is
633                performed. If 'None', 50 equally spaced values between 2 and 12
634                are used.
635
636            pka_set ('dict', optional):
637                Dictionary defining the acid–base properties of titratable particle
638                types:
639                '{particle_name: {"pka_value": float, "acidity": "acidic" | "basic"}}'.
640                If 'None', the pKa set is taken from the pyMBE database.
641
642        Returns:
643            'dict':
644                Dictionary containing:
645                - '"charges_dict"' ('dict'):
646                    Mapping '{molecule_name: list}' of Henderson–Hasselbalch–Donnan
647                    charges evaluated at each pH value.
648                - '"pH_system_list"' ('list[float]'):
649                    Effective pH values inside the system phase after Donnan
650                    partitioning.
651                - '"partition_coefficients"' ('list[float]'):
652                    Partition coefficients of monovalent cations at each pH value.
653
654        Notes:
655            - This method assumes **ideal Donnan equilibrium** and **monovalent salt**.
656            - The ionic strength of the reservoir includes both salt and
657            pH-dependent H⁺/OH⁻ contributions.
658            - All charged macromolecular species present in the system must be
659            included in 'c_macro'; missing species will lead to incorrect results.
660            - The nonlinear Donnan equilibrium equation is solved using a scalar
661            root finder ('brentq') in logarithmic form for numerical stability.
662            - This method is intended for **two-phase systems**; for single-phase
663            systems use 'calculate_HH' instead.
664        """
665        if pH_list is None:
666            pH_list=np.linspace(2,12,50)
667        if pka_set is None:
668            pka_set=self.get_pka_set() 
669        self._check_pka_set(pka_set=pka_set)
670        partition_coefficients_list = []
671        pH_system_list = []
672        Z_HH_Donnan={}
673        for key in c_macro:
674            Z_HH_Donnan[key] = []
675        def calc_charges(c_macro, pH):
676            """
677            Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation.
678
679            Args:
680                c_macro ('dict'): 
681                    {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
682
683                pH ('float'): 
684                    pH-value that is used in the HH equation.
685
686            Returns:
687                ('dict'): 
688                    {"molecule_name": charge}
689            """
690            charge = {}
691            for name in c_macro:
692                charge[name] = self.calculate_HH(name, [pH], pka_set)[0]
693            return charge
694
695        def calc_partition_coefficient(charge, c_macro):
696            """
697            Calculates the partition coefficients of positive ions according to the ideal Donnan theory.
698
699            Args:
700                charge ('dict'): 
701                    {"molecule_name": charge}
702
703                c_macro ('dict'): 
704                    {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
705            """
706            nonlocal ionic_strength_res
707            charge_density = 0.0
708            for key in charge:
709                charge_density += charge[key] * c_macro[key]
710            return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude
711        for pH_value in pH_list:    
712            # calculate the ionic strength of the reservoir
713            if pH_value <= 7.0:
714                ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 
715            elif pH_value > 7.0:
716                ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt
717            #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations
718            #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation.
719            #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects.
720            equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro))
721            logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq")
722            partition_coefficient = 10**logxi.root
723            charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient))
724            for key in c_macro:
725                Z_HH_Donnan[key].append(charges_temp[key])
726            pH_system_list.append(pH_value - np.log10(partition_coefficient))
727            partition_coefficients_list.append(partition_coefficient)
728        return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}

Computes macromolecular charges using the Henderson–Hasselbalch equation coupled to ideal Donnan partitioning.

Arguments:
  • c_macro ('dict'): Mapping of macromolecular species names to their concentrations in the system: '{molecule_name: concentration}'.
  • c_salt ('float' or 'pint.Quantity'): Salt concentration in the reservoir.
  • pH_list ('list[float]', optional): List of pH values in the reservoir at which the calculation is performed. If 'None', 50 equally spaced values between 2 and 12 are used.
  • pka_set ('dict', optional): Dictionary defining the acid–base properties of titratable particle types: '{particle_name: {"pka_value": float, "acidity": "acidic" | "basic"}}'. If 'None', the pKa set is taken from the pyMBE database.
Returns:

'dict': Dictionary containing: - '"charges_dict"' ('dict'): Mapping '{molecule_name: list}' of Henderson–Hasselbalch–Donnan charges evaluated at each pH value. - '"pH_system_list"' ('list[float]'): Effective pH values inside the system phase after Donnan partitioning. - '"partition_coefficients"' ('list[float]'): Partition coefficients of monovalent cations at each pH value.

Notes:
  • This method assumes ideal Donnan equilibrium and monovalent salt.
  • The ionic strength of the reservoir includes both salt and pH-dependent H⁺/OH⁻ contributions.
  • All charged macromolecular species present in the system must be included in 'c_macro'; missing species will lead to incorrect results.
  • The nonlinear Donnan equilibrium equation is solved using a scalar root finder ('brentq') in logarithmic form for numerical stability.
  • This method is intended for two-phase systems; for single-phase systems use 'calculate_HH' instead.
def calculate_net_charge(self, espresso_system, object_name, pmb_type, dimensionless=False):
730    def calculate_net_charge(self,espresso_system,object_name,pmb_type,dimensionless=False):
731        """
732        Calculates the net charge per instance of a given pmb object type.
733
734        Args:
735            espresso_system (espressomd.system.System):
736                ESPResSo system containing the particles.
737            object_name (str):
738                Name of the object (e.g. molecule, residue, peptide, protein).
739            pmb_type (str):
740                Type of object to analyze. Must be molecule-like.
741            dimensionless (bool, optional):
742                If True, return charge as a pure number.
743                If False, return a quantity with reduced_charge units.
744
745        Returns:
746            dict:
747                {"mean": mean_net_charge, "instances": {instance_id: net_charge}}
748        """
749        id_map = self.get_particle_id_map(object_name=object_name)
750        label = self._get_label_id_map(pmb_type=pmb_type)
751        instance_map = id_map[label]
752        charges = {}
753        for instance_id, particle_ids in instance_map.items():
754            if dimensionless:
755                net_charge = 0.0
756            else:
757                net_charge = 0 * self.units.Quantity(1, "reduced_charge")
758            for pid in particle_ids:
759                q = espresso_system.part.by_id(pid).q
760                if not dimensionless:
761                    q *= self.units.Quantity(1, "reduced_charge")
762                net_charge += q
763            charges[instance_id] = net_charge
764        # Mean charge
765        if dimensionless:
766            mean_charge = float(np.mean(list(charges.values())))
767        else:
768            mean_charge = (np.mean([q.magnitude for q in charges.values()])* self.units.Quantity(1, "reduced_charge"))
769        return {"mean": mean_charge, "instances": charges}

Calculates the net charge per instance of a given pmb object type.

Arguments:
  • espresso_system (espressomd.system.System): ESPResSo system containing the particles.
  • object_name (str): Name of the object (e.g. molecule, residue, peptide, protein).
  • pmb_type (str): Type of object to analyze. Must be molecule-like.
  • dimensionless (bool, optional): If True, return charge as a pure number. If False, return a quantity with reduced_charge units.
Returns:

dict: {"mean": mean_net_charge, "instances": {instance_id: net_charge}}

def center_object_in_simulation_box(self, instance_id, espresso_system, pmb_type):
771    def center_object_in_simulation_box(self, instance_id, espresso_system, pmb_type):
772        """
773        Centers a pyMBE object instance in the simulation box of an ESPResSo system.
774        The object is translated such that its center of mass coincides with the
775        geometric center of the ESPResSo simulation box.
776
777        Args:
778            instance_id ('int'):
779                ID of the pyMBE object instance to be centered.
780
781            pmb_type ('str'):
782                Type of the pyMBE object.
783
784            espresso_system ('espressomd.system.System'):
785                ESPResSo system object in which the particles are defined.
786
787        Notes:
788            - Works for both cubic and non-cubic simulation boxes.
789        """
790        inst = self.db.get_instance(instance_id=instance_id,
791                                    pmb_type=pmb_type)
792        center_of_mass = self.calculate_center_of_mass(instance_id=instance_id,
793                                                       espresso_system=espresso_system,
794                                                       pmb_type=pmb_type)
795        box_center = [espresso_system.box_l[0]/2.0,
796                      espresso_system.box_l[1]/2.0,
797                      espresso_system.box_l[2]/2.0]
798        particle_id_list = self.get_particle_id_map(object_name=inst.name)["all"]
799        for pid in particle_id_list:
800            es_pos = espresso_system.part.by_id(pid).pos
801            espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center

Centers a pyMBE object instance in the simulation box of an ESPResSo system. The object is translated such that its center of mass coincides with the geometric center of the ESPResSo simulation box.

Arguments:
  • instance_id ('int'): ID of the pyMBE object instance to be centered.
  • pmb_type ('str'): Type of the pyMBE object.
  • espresso_system ('espressomd.system.System'): ESPResSo system object in which the particles are defined.
Notes:
  • Works for both cubic and non-cubic simulation boxes.
def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt):
803    def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt):    
804        """
805        Creates a 'c_salt' concentration of 'cation_name' and 'anion_name' ions into the 'espresso_system'.
806
807        Args:
808            espresso_system('espressomd.system.System'): instance of an espresso system object.
809            cation_name('str'): 'name' of a particle with a positive charge.
810            anion_name('str'): 'name' of a particle with a negative charge.
811            c_salt('float'): Salt concentration.
812            
813        Returns:
814            c_salt_calculated('float'): Calculated salt concentration added to 'espresso_system'.
815        """ 
816        cation_tpl = self.db.get_template(pmb_type="particle",
817                                          name=cation_name)
818        cation_state = self.db.get_template(pmb_type="particle_state",
819                                            name=cation_tpl.initial_state)
820        cation_charge = cation_state.z
821        anion_tpl = self.db.get_template(pmb_type="particle",
822                                          name=anion_name)
823        anion_state = self.db.get_template(pmb_type="particle_state",
824                                            name=anion_tpl.initial_state)
825        anion_charge = anion_state.z
826        if cation_charge <= 0:
827            raise ValueError(f'ERROR cation charge must be positive, charge {cation_charge}')
828        if anion_charge >= 0:
829            raise ValueError(f'ERROR anion charge must be negative, charge {anion_charge}')
830        # Calculate the number of ions in the simulation box
831        volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3')
832        if c_salt.check('[substance] [length]**-3'):
833            N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude)
834            c_salt_calculated=N_ions/(volume*self.N_A)
835        elif c_salt.check('[length]**-3'):
836            N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude)
837            c_salt_calculated=N_ions/volume
838        else:
839            raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt)
840        N_cation = N_ions*abs(anion_charge)
841        N_anion = N_ions*abs(cation_charge)
842        self.create_particle(espresso_system=espresso_system, 
843                             name=cation_name, 
844                             number_of_particles=N_cation)
845        self.create_particle(espresso_system=espresso_system, 
846                             name=anion_name, 
847                             number_of_particles=N_anion)
848        if c_salt_calculated.check('[substance] [length]**-3'):
849            logging.info(f"added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions")
850        elif c_salt_calculated.check('[length]**-3'):
851            logging.info(f"added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions")
852        return c_salt_calculated

Creates a 'c_salt' concentration of 'cation_name' and 'anion_name' ions into the 'espresso_system'.

Arguments:
  • espresso_system('espressomd.system.System'): instance of an espresso system object.
  • cation_name('str'): 'name' of a particle with a positive charge.
  • anion_name('str'): 'name' of a particle with a negative charge.
  • c_salt('float'): Salt concentration.
Returns:

c_salt_calculated('float'): Calculated salt concentration added to 'espresso_system'.

def create_bond( self, particle_id1, particle_id2, espresso_system, use_default_bond=False):
854    def create_bond(self, particle_id1, particle_id2, espresso_system, use_default_bond=False):
855        """
856        Creates a bond between two particle instances in an ESPResSo system and registers it in the pyMBE database.
857
858        This method performs the following steps:
859            1. Retrieves the particle instances corresponding to 'particle_id1' and 'particle_id2' from the database.
860            2. Retrieves or creates the corresponding ESPResSo bond instance using the bond template.
861            3. Adds the ESPResSo bond instance to the ESPResSo system if it was newly created.
862            4. Adds the bond to the first particle's bond list in ESPResSo.
863            5. Creates a 'BondInstance' in the database and registers it.
864
865        Args:
866            particle_id1 ('int'): 
867                pyMBE and ESPResSo ID of the first particle.
868
869            particle_id2 ('int'): 
870                pyMBE and ESPResSo ID of the second particle.
871
872            espresso_system ('espressomd.system.System'): 
873                ESPResSo system object where the bond will be created.
874
875            use_default_bond ('bool', optional): 
876                If True, use a default bond template if no specific template exists. Defaults to False.
877
878        Returns:
879            ('int'): 
880                bond_id of the bond instance created in the pyMBE database.
881        """
882        particle_inst_1 = self.db.get_instance(pmb_type="particle",
883                                               instance_id=particle_id1)
884        particle_inst_2 = self.db.get_instance(pmb_type="particle",
885                                               instance_id=particle_id2)
886        bond_tpl = self.get_bond_template(particle_name1=particle_inst_1.name,
887                                          particle_name2=particle_inst_2.name,
888                                          use_default_bond=use_default_bond)
889        bond_inst = self._get_espresso_bond_instance(bond_template=bond_tpl,
890                                                    espresso_system=espresso_system)
891        espresso_system.part.by_id(particle_id1).add_bond((bond_inst, particle_id2))
892        bond_id = self.db._propose_instance_id(pmb_type="bond")
893        pmb_bond_instance = BondInstance(bond_id=bond_id,
894                                         name=bond_tpl.name,
895                                         particle_id1=particle_id1,
896                                         particle_id2=particle_id2)
897        self.db._register_instance(instance=pmb_bond_instance)

Creates a bond between two particle instances in an ESPResSo system and registers it in the pyMBE database.

This method performs the following steps:
  1. Retrieves the particle instances corresponding to 'particle_id1' and 'particle_id2' from the database.
  2. Retrieves or creates the corresponding ESPResSo bond instance using the bond template.
  3. Adds the ESPResSo bond instance to the ESPResSo system if it was newly created.
  4. Adds the bond to the first particle's bond list in ESPResSo.
  5. Creates a 'BondInstance' in the database and registers it.
Arguments:
  • particle_id1 ('int'): pyMBE and ESPResSo ID of the first particle.
  • particle_id2 ('int'): pyMBE and ESPResSo ID of the second particle.
  • espresso_system ('espressomd.system.System'): ESPResSo system object where the bond will be created.
  • use_default_bond ('bool', optional): If True, use a default bond template if no specific template exists. Defaults to False.
Returns:

('int'): bond_id of the bond instance created in the pyMBE database.

def create_counterions(self, object_name, cation_name, anion_name, espresso_system):
899    def create_counterions(self, object_name, cation_name, anion_name, espresso_system):
900        """
901        Creates particles of 'cation_name' and 'anion_name' in 'espresso_system' to counter the net charge of 'object_name'.
902        
903        Args:
904            object_name ('str'): 
905                'name' of a pyMBE object.
906
907            espresso_system ('espressomd.system.System'): 
908                Instance of a system object from the espressomd library.
909
910            cation_name ('str'): 
911                'name' of a particle with a positive charge.
912
913            anion_name ('str'): 
914                'name' of a particle with a negative charge.
915
916        Returns: 
917            ('dict'): 
918                {"name": number}
919
920        Notes:
921            This function currently does not support the creation of counterions for hydrogels.
922        """ 
923        cation_tpl = self.db.get_template(pmb_type="particle",
924                                          name=cation_name)
925        cation_state = self.db.get_template(pmb_type="particle_state",
926                                            name=cation_tpl.initial_state)
927        cation_charge = cation_state.z
928        anion_tpl = self.db.get_template(pmb_type="particle",
929                                          name=anion_name)
930        anion_state = self.db.get_template(pmb_type="particle_state",
931                                            name=anion_tpl.initial_state)
932        anion_charge = anion_state.z
933        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
934        counterion_number={}
935        object_charge={}
936        for name in ['positive', 'negative']:
937            object_charge[name]=0
938        for id in object_ids:
939            if espresso_system.part.by_id(id).q > 0:
940                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
941            elif espresso_system.part.by_id(id).q < 0:
942                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
943        if object_charge['positive'] % abs(anion_charge) == 0:
944            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
945        else:
946            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
947        if object_charge['negative'] % abs(cation_charge) == 0:
948            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
949        else:
950            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
951        if counterion_number[cation_name] > 0: 
952            self.create_particle(espresso_system=espresso_system, 
953                                 name=cation_name, 
954                                 number_of_particles=counterion_number[cation_name])
955        else:
956            counterion_number[cation_name]=0
957        if counterion_number[anion_name] > 0:
958            self.create_particle(espresso_system=espresso_system, 
959                                 name=anion_name, 
960                                 number_of_particles=counterion_number[anion_name])
961        else:
962            counterion_number[anion_name] = 0
963        logging.info('the following counter-ions have been created: ')
964        for name in counterion_number.keys():
965            logging.info(f'Ion type: {name} created number: {counterion_number[name]}')
966        return counterion_number

Creates particles of 'cation_name' and 'anion_name' in 'espresso_system' to counter the net charge of 'object_name'.

Arguments:
  • object_name ('str'): 'name' of a pyMBE object.
  • espresso_system ('espressomd.system.System'): Instance of a system object from the espressomd library.
  • cation_name ('str'): 'name' of a particle with a positive charge.
  • anion_name ('str'): 'name' of a particle with a negative charge.

Returns: ('dict'): {"name": number}

Notes:

This function currently does not support the creation of counterions for hydrogels.

def create_hydrogel(self, name, espresso_system, use_default_bond=False, gen_angle=False):
 968    def create_hydrogel(self, name, espresso_system, use_default_bond=False, gen_angle=False):
 969        """ 
 970        Creates a hydrogel in espresso_system using a pyMBE hydrogel template given by 'name'
 971
 972        Args:
 973            name ('str'): 
 974                name of the hydrogel template in the pyMBE database.
 975
 976            espresso_system ('espressomd.system.System'): 
 977                ESPResSo system object where the hydrogel will be created.
 978
 979            use_default_bond ('bool', optional): 
 980                If True, use a default bond template if no specific template exists. Defaults to False.
 981
 982            gen_angle ('bool', optional):
 983                If True, generate angle potentials for the internal hydrogel
 984                chains and, when explicitly defined, for all crosslinker-adjacent
 985                triplets. Defaults to False.
 986
 987        Returns:
 988            ('int'): id of the hydrogel instance created.
 989        """
 990        if not self.db._has_template(name=name, pmb_type="hydrogel"):
 991            raise ValueError(f"Hydrogel template with name '{name}' is not defined in the pyMBE database.")
 992        hydrogel_tpl = self.db.get_template(pmb_type="hydrogel",
 993                                            name=name)
 994        assembly_id = self.db._propose_instance_id(pmb_type="hydrogel")
 995        # Create the nodes
 996        nodes = {}
 997        hydrogel_angle_centers = set()
 998        node_topology = hydrogel_tpl.node_map
 999        for node in node_topology:
1000            node_index = node.lattice_index
1001            node_name = node.particle_name
1002            node_pos, node_id = self._create_hydrogel_node(node_index=node_index,
1003                                                          node_name=node_name,
1004                                                          espresso_system=espresso_system)
1005            node_label = self.lattice_builder._create_node_label(node_index=node_index)
1006            nodes[node_label] = {"name": node_name, "id": node_id, "pos": node_pos} 
1007            self.db._update_instance(instance_id=node_id,
1008                                     pmb_type="particle",
1009                                     attribute="assembly_id",
1010                                     value=assembly_id)
1011        for hydrogel_chain in hydrogel_tpl.chain_map:
1012            molecule_id = self._create_hydrogel_chain(hydrogel_chain=hydrogel_chain,
1013                                                      nodes=nodes, 
1014                                                      espresso_system=espresso_system,
1015                                                      use_default_bond=use_default_bond,
1016                                                      gen_angle=gen_angle)
1017            self.db._update_instance(instance_id=molecule_id,
1018                                     pmb_type="molecule",
1019                                     attribute="assembly_id",
1020                                     value=assembly_id)
1021            if gen_angle:
1022                residue_ids = self.db._find_instance_ids_by_attribute(pmb_type="residue",
1023                                                                      attribute="molecule_id",
1024                                                                      value=molecule_id)
1025                first_residue_id = min(residue_ids)
1026                last_residue_id = max(residue_ids)
1027                first_residue = self.db.get_instance(pmb_type="residue",
1028                                                     instance_id=first_residue_id)
1029                last_residue = self.db.get_instance(pmb_type="residue",
1030                                                    instance_id=last_residue_id)
1031                first_central_bead_name = self.db.get_template(pmb_type="residue",
1032                                                               name=first_residue.name).central_bead
1033                last_central_bead_name = self.db.get_template(pmb_type="residue",
1034                                                              name=last_residue.name).central_bead
1035                particle_instances = self.db.get_instances(pmb_type="particle")
1036                first_residue_particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1037                                                                                      attribute="residue_id",
1038                                                                                      value=first_residue_id)
1039                last_residue_particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1040                                                                                     attribute="residue_id",
1041                                                                                     value=last_residue_id)
1042                first_bead_id = None
1043                for particle_id in first_residue_particle_ids:
1044                    if particle_instances[particle_id].name == first_central_bead_name:
1045                        first_bead_id = particle_id
1046                        break
1047
1048                last_bead_id = None
1049                for particle_id in last_residue_particle_ids:
1050                    if particle_instances[particle_id].name == last_central_bead_name:
1051                        last_bead_id = particle_id
1052                        break
1053                node_start_label = self.lattice_builder._create_node_label(hydrogel_chain.node_start)
1054                node_end_label = self.lattice_builder._create_node_label(hydrogel_chain.node_end)
1055                hydrogel_angle_centers.update({
1056                    nodes[node_start_label]["id"],
1057                    nodes[node_end_label]["id"],
1058                    first_bead_id,
1059                    last_bead_id,
1060                })
1061        self.db._propagate_id(root_type="hydrogel", 
1062                                root_id=assembly_id, 
1063                                attribute="assembly_id", 
1064                                value=assembly_id)
1065        if gen_angle:
1066            self._generate_hydrogel_crosslinker_angles(espresso_system=espresso_system,
1067                                                       central_particle_ids=hydrogel_angle_centers)
1068        # Register an hydrogel instance in the pyMBE databasegit 
1069        self.db._register_instance(HydrogelInstance(name=name,
1070                                                    assembly_id=assembly_id))
1071        return assembly_id

Creates a hydrogel in espresso_system using a pyMBE hydrogel template given by 'name'

Arguments:
  • name ('str'): name of the hydrogel template in the pyMBE database.
  • espresso_system ('espressomd.system.System'): ESPResSo system object where the hydrogel will be created.
  • use_default_bond ('bool', optional): If True, use a default bond template if no specific template exists. Defaults to False.
  • gen_angle ('bool', optional): If True, generate angle potentials for the internal hydrogel chains and, when explicitly defined, for all crosslinker-adjacent triplets. Defaults to False.
Returns:

('int'): id of the hydrogel instance created.

def create_molecule( self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False, reverse_residue_order=False, gen_angle=False):
1073    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False, reverse_residue_order = False, gen_angle=False):
1074        """
1075        Creates instances of a given molecule template name into ESPResSo.
1076
1077        Args:
1078            name ('str'): 
1079                Label of the molecule type to be created. 'name'.
1080
1081            espresso_system ('espressomd.system.System'): 
1082                Instance of a system object from espressomd library.
1083
1084            number_of_molecules ('int'): 
1085                Number of molecules or peptides of type 'name' to be created.
1086
1087            list_of_first_residue_positions ('list', optional): 
1088                List of coordinates where the central bead of the first_residue_position will be created, random by default.
1089
1090            backbone_vector ('list' of 'float'): 
1091                Backbone vector of the molecule, random by default. Central beads of the residues in the 'residue_list' are placed along this vector. 
1092
1093            use_default_bond('bool', optional): 
1094                Controls if a bond of type 'default' is used to bond particles with undefined bonds in the pyMBE database.
1095
1096            reverse_residue_order('bool', optional): 
1097                Creates residues in reverse sequential order than the one defined in the molecule template. Defaults to False.
1098
1099        Returns:
1100            ('list' of 'int'): 
1101                List with the 'molecule_id' of the pyMBE molecule instances created into 'espresso_system'.
1102
1103        Notes:
1104            - This function can be used to create both molecules and peptides.    
1105        """
1106        pmb_type = self._get_template_type(name=name,
1107                                           allowed_types={"molecule", "peptide"})
1108        if number_of_molecules <= 0:
1109            return {}
1110        if list_of_first_residue_positions is not None:
1111            for item in list_of_first_residue_positions:
1112                if not isinstance(item, list):
1113                    raise ValueError("The provided input position is not a nested list. Should be a nested list with elements of 3D lists, corresponding to xyz coord.")
1114                elif len(item) != 3:
1115                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
1116
1117            if len(list_of_first_residue_positions) != number_of_molecules:
1118                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
1119        # Generate an arbitrary random unit vector
1120        if backbone_vector is None:
1121            backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
1122                                                                      radius=1, 
1123                                                                      n_samples=1,
1124                                                                      on_surface=True)[0]
1125        else:
1126            backbone_vector = np.array(backbone_vector)
1127        first_residue = True
1128        molecule_tpl = self.db.get_template(pmb_type=pmb_type,
1129                                            name=name)
1130        if reverse_residue_order:
1131            residue_list = molecule_tpl.residue_list[::-1]
1132        else:
1133            residue_list = molecule_tpl.residue_list
1134        pos_index = 0 
1135        molecule_ids = []
1136        for _ in range(number_of_molecules):        
1137            molecule_id = self.db._propose_instance_id(pmb_type=pmb_type)
1138            for residue in residue_list:
1139                if first_residue:
1140                    if list_of_first_residue_positions is None:
1141                        central_bead_pos = None
1142                    else:
1143                        for item in list_of_first_residue_positions:
1144                            central_bead_pos = [np.array(list_of_first_residue_positions[pos_index])]
1145                            
1146                    residue_id = self.create_residue(name=residue,
1147                                                     espresso_system=espresso_system, 
1148                                                     central_bead_position=central_bead_pos,  
1149                                                     use_default_bond= use_default_bond, 
1150                                                     backbone_vector=backbone_vector)
1151                    
1152                    # Add molecule_id to the residue instance and all particles associated
1153                    self.db._propagate_id(root_type="residue", 
1154                                          root_id=residue_id,
1155                                          attribute="molecule_id", 
1156                                          value=molecule_id)
1157                    particle_ids_in_residue = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1158                                                                                      attribute="residue_id",
1159                                                                                      value=residue_id)
1160                    prev_central_bead_id = particle_ids_in_residue[0]
1161                    prev_central_bead_name = self.db.get_instance(pmb_type="particle", 
1162                                                                  instance_id=prev_central_bead_id).name
1163                    prev_central_bead_pos = espresso_system.part.by_id(prev_central_bead_id).pos
1164                    first_residue = False          
1165                else:
1166                    
1167                    # Calculate the starting position of the new residue
1168                    residue_tpl = self.db.get_template(pmb_type="residue",
1169                                                       name=residue)
1170                    lj_parameters = self.get_lj_parameters(particle_name1=prev_central_bead_name,
1171                                                           particle_name2=residue_tpl.central_bead)
1172                    bond_tpl = self.get_bond_template(particle_name1=prev_central_bead_name,
1173                                                      particle_name2=residue_tpl.central_bead,
1174                                                      use_default_bond=use_default_bond)
1175                    l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1176                                                          bond_type=bond_tpl.bond_type,
1177                                                          bond_parameters=bond_tpl.get_parameters(ureg=self.units))
1178                    central_bead_pos = prev_central_bead_pos+backbone_vector*l0
1179                    # Create the residue
1180                    residue_id = self.create_residue(name=residue, 
1181                                                     espresso_system=espresso_system, 
1182                                                     central_bead_position=[central_bead_pos],
1183                                                     use_default_bond= use_default_bond, 
1184                                                     backbone_vector=backbone_vector)
1185                    # Add molecule_id to the residue instance and all particles associated
1186                    self.db._propagate_id(root_type="residue", 
1187                                          root_id=residue_id, 
1188                                          attribute="molecule_id", 
1189                                          value=molecule_id)
1190                    particle_ids_in_residue = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1191                                                                                      attribute="residue_id",
1192                                                                                      value=residue_id)
1193                    central_bead_id = particle_ids_in_residue[0]
1194
1195                    # Bond the central beads of the new and previous residues
1196                    self.create_bond(particle_id1=prev_central_bead_id,
1197                                     particle_id2=central_bead_id,
1198                                     espresso_system=espresso_system,
1199                                     use_default_bond=use_default_bond)
1200                    
1201                    prev_central_bead_id = central_bead_id                    
1202                    prev_central_bead_name = self.db.get_instance(pmb_type="particle", instance_id=central_bead_id).name
1203                    prev_central_bead_pos =central_bead_pos
1204            # Create a Peptide or Molecule instance and register it on the pyMBE database
1205            if pmb_type == "molecule":
1206                inst = MoleculeInstance(molecule_id=molecule_id,
1207                                        name=name)
1208            elif pmb_type == "peptide":
1209                inst = PeptideInstance(name=name,
1210                                       molecule_id=molecule_id)
1211            self.db._register_instance(inst)
1212            if gen_angle:
1213                self._generate_angles_for_entity(
1214                    espresso_system=espresso_system,
1215                    entity_id=molecule_id,
1216                    entity_id_col='molecule_id')
1217            first_residue = True
1218            pos_index+=1
1219            molecule_ids.append(molecule_id)
1220        return molecule_ids

Creates instances of a given molecule template name into ESPResSo.

Arguments:
  • name ('str'): Label of the molecule type to be created. 'name'.
  • espresso_system ('espressomd.system.System'): Instance of a system object from espressomd library.
  • number_of_molecules ('int'): Number of molecules or peptides of type 'name' to be created.
  • list_of_first_residue_positions ('list', optional): List of coordinates where the central bead of the first_residue_position will be created, random by default.
  • backbone_vector ('list' of 'float'): Backbone vector of the molecule, random by default. Central beads of the residues in the 'residue_list' are placed along this vector.
  • use_default_bond('bool', optional): Controls if a bond of type 'default' is used to bond particles with undefined bonds in the pyMBE database.
  • reverse_residue_order('bool', optional): Creates residues in reverse sequential order than the one defined in the molecule template. Defaults to False.
Returns:

('list' of 'int'): List with the 'molecule_id' of the pyMBE molecule instances created into 'espresso_system'.

Notes:
  • This function can be used to create both molecules and peptides.
def create_particle( self, name, espresso_system, number_of_particles, position=None, fix=False):
1222    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
1223        """
1224        Creates one or more particles in an ESPResSo system based on the particle template in the pyMBE database.
1225        
1226        Args:
1227            name ('str'): 
1228                Label of the particle template in the pyMBE database. 
1229
1230            espresso_system ('espressomd.system.System'): 
1231                Instance of a system object from the espressomd library.
1232
1233            number_of_particles ('int'): 
1234                Number of particles to be created.
1235
1236            position (list of ['float','float','float'], optional): 
1237                Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
1238
1239            fix ('bool', optional): 
1240                Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
1241
1242        Returns:
1243            ('list' of 'int'): 
1244                List with the ids of the particles created into 'espresso_system'.
1245        """       
1246        if number_of_particles <=0:
1247            return []
1248        if not self.db._has_template(name=name, pmb_type="particle"):
1249            raise ValueError(f"Particle template with name '{name}' is not defined in the pyMBE database.")
1250        
1251        part_tpl = self.db.get_template(pmb_type="particle",
1252                                        name=name)
1253        part_state = self.db.get_template(pmb_type="particle_state",
1254                                         name=part_tpl.initial_state)
1255        z = part_state.z
1256        es_type = part_state.es_type
1257        # Create the new particles into  ESPResSo 
1258        created_pid_list=[]
1259        for index in range(number_of_particles):
1260            if position is None:
1261                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1262            else:
1263                particle_position = position[index]
1264            
1265            particle_id = self.db._propose_instance_id(pmb_type="particle")
1266            created_pid_list.append(particle_id)
1267            kwargs = dict(id=particle_id, pos=particle_position, type=es_type, q=z)
1268            if fix:
1269                kwargs["fix"] = 3 * [fix]
1270            espresso_system.part.add(**kwargs)
1271            part_inst = ParticleInstance(name=name,
1272                                         particle_id=particle_id,
1273                                         initial_state=part_state.name)
1274            self.db._register_instance(part_inst)
1275                              
1276        return created_pid_list

Creates one or more particles in an ESPResSo system based on the particle template in the pyMBE database.

Arguments:
  • name ('str'): Label of the particle template in the pyMBE database.
  • espresso_system ('espressomd.system.System'): Instance of a system object from the espressomd library.
  • number_of_particles ('int'): Number of particles to be created.
  • position (list of ['float','float','float'], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
  • fix ('bool', optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
Returns:

('list' of 'int'): List with the ids of the particles created into 'espresso_system'.

def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1278    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1279        """
1280        Creates one or more protein molecules in an ESPResSo system based on the 
1281        protein template in the pyMBE database and a provided topology.
1282
1283        Args:
1284            name (str):
1285                Name of the protein template stored in the pyMBE database.
1286            
1287            number_of_proteins (int):
1288                Number of protein molecules to generate.  
1289            
1290            espresso_system (espressomd.system.System):
1291                The ESPResSo simulation system where the protein molecules will be created.
1292            
1293            topology_dict (dict):
1294                Dictionary defining the internal structure of the protein. Expected format:
1295                    {"ResidueName1": {"initial_pos": np.ndarray,
1296                                      "chain_id": int,
1297                                      "radius": float},
1298                     "ResidueName2": { ... },
1299                        ...
1300                    }
1301                The '"initial_pos"' entry is required and represents the residue’s
1302                reference coordinates before shifting to the protein's center-of-mass.
1303
1304        Returns:
1305            ('list' of 'int'): 
1306                List of the molecule_id of the Protein instances created into ESPResSo.
1307
1308        Notes:
1309            - Particles are created using 'create_particle()' with 'fix=True',
1310            meaning they are initially immobilized.
1311            - The function assumes all residues in 'topology_dict' correspond to
1312            particle templates already defined in the pyMBE database.
1313            - Bonds between residues are not created here; it assumes a rigid body representation of the protein.
1314        """
1315        if number_of_proteins <= 0:
1316            return
1317        if not self.db._has_template(name=name, pmb_type="protein"):
1318            raise ValueError(f"Protein template with name '{name}' is not defined in the pyMBE database.")
1319        protein_tpl = self.db.get_template(pmb_type="protein", name=name)
1320        box_half = espresso_system.box_l[0] / 2.0
1321        # Create protein
1322        mol_ids = []
1323        for _ in range(number_of_proteins):
1324            # create a molecule identifier in pyMBE
1325            molecule_id = self.db._propose_instance_id(pmb_type="protein")
1326            # place protein COM randomly
1327            protein_center = self.generate_coordinates_outside_sphere(radius=1,
1328                                                                      max_dist=box_half,
1329                                                                      n_samples=1,
1330                                                                      center=[box_half]*3)[0]
1331            residues = hf.get_residues_from_topology_dict(topology_dict=topology_dict,
1332                                                         model=protein_tpl.model)
1333            # CREATE RESIDUES + PARTICLES
1334            for _, rdata in residues.items():
1335                base_resname = rdata["resname"]  
1336                residue_name = f"AA-{base_resname}"
1337                # residue instance ID
1338                residue_id = self.db._propose_instance_id("residue")
1339                # register ResidueInstance
1340                self.db._register_instance(ResidueInstance(name=residue_name,
1341                                                           residue_id=residue_id,
1342                                                           molecule_id=molecule_id))
1343                # PARTICLE CREATION
1344                for bead_id in rdata["beads"]:
1345                    bead_type = re.split(r'\d+', bead_id)[0]
1346                    relative_pos = topology_dict[bead_id]["initial_pos"]
1347                    absolute_pos = relative_pos + protein_center
1348                    particle_id = self.create_particle(name=bead_type,
1349                                                       espresso_system=espresso_system,
1350                                                       number_of_particles=1,
1351                                                       position=[absolute_pos],
1352                                                       fix=True)[0]
1353                    # update metadata
1354                    self.db._update_instance(instance_id=particle_id,
1355                                             pmb_type="particle",
1356                                             attribute="molecule_id",
1357                                             value=molecule_id)
1358                    self.db._update_instance(instance_id=particle_id,
1359                                             pmb_type="particle",
1360                                             attribute="residue_id",
1361                                             value=residue_id)
1362            protein_inst = ProteinInstance(name=name,
1363                                           molecule_id=molecule_id)
1364            self.db._register_instance(protein_inst)
1365            mol_ids.append(molecule_id)
1366        return mol_ids

Creates one or more protein molecules in an ESPResSo system based on the protein template in the pyMBE database and a provided topology.

Arguments:
  • name (str): Name of the protein template stored in the pyMBE database.
  • number_of_proteins (int): Number of protein molecules to generate.
  • espresso_system (espressomd.system.System): The ESPResSo simulation system where the protein molecules will be created.
  • topology_dict (dict): Dictionary defining the internal structure of the protein. Expected format: {"ResidueName1": {"initial_pos": np.ndarray, "chain_id": int, "radius": float}, "ResidueName2": { ... }, ... } The '"initial_pos"' entry is required and represents the residue’s reference coordinates before shifting to the protein's center-of-mass.
Returns:

('list' of 'int'): List of the molecule_id of the Protein instances created into ESPResSo.

Notes:
  • Particles are created using 'create_particle()' with 'fix=True', meaning they are initially immobilized.
  • The function assumes all residues in 'topology_dict' correspond to particle templates already defined in the pyMBE database.
  • Bonds between residues are not created here; it assumes a rigid body representation of the protein.
def create_residue( self, name, espresso_system, central_bead_position=None, use_default_bond=False, backbone_vector=None, gen_angle=False):
1368    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None, gen_angle=False):
1369        """
1370        Creates a residue  into ESPResSo.
1371
1372        Args:
1373            name ('str'): 
1374                Label of the residue type to be created. 
1375
1376            espresso_system ('espressomd.system.System'): 
1377                Instance of a system object from espressomd library.
1378
1379            central_bead_position ('list' of 'float'): 
1380                Position of the central bead.
1381
1382            use_default_bond ('bool'): 
1383                Switch to control if a bond of type 'default' is used to bond a particle whose bonds types are not defined in the pyMBE database.
1384
1385            backbone_vector ('list' of 'float'): 
1386                Backbone vector of the molecule. All side chains are created perpendicularly to 'backbone_vector'.
1387
1388        Returns:
1389            (int): 
1390                residue_id of the residue created.
1391        """
1392        if not self.db._has_template(name=name, pmb_type="residue"):
1393            raise ValueError(f"Residue template with name '{name}' is not defined in the pyMBE database.")
1394        res_tpl = self.db.get_template(pmb_type="residue",
1395                                       name=name)
1396        # Assign a residue_id
1397        residue_id = self.db._propose_instance_id(pmb_type="residue")
1398        res_inst = ResidueInstance(name=name,
1399                                   residue_id=residue_id)
1400        self.db._register_instance(res_inst)
1401        # create the principal bead   
1402        central_bead_name = res_tpl.central_bead 
1403        central_bead_id = self.create_particle(name=central_bead_name,
1404                                               espresso_system=espresso_system,
1405                                               position=central_bead_position,
1406                                               number_of_particles = 1)[0]
1407        
1408        central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1409        # Assigns residue_id to the central_bead particle created.
1410        self.db._update_instance(pmb_type="particle",
1411                                 instance_id=central_bead_id,
1412                                 attribute="residue_id",
1413                                 value=residue_id)
1414        
1415        # create the lateral beads  
1416        side_chain_list = res_tpl.side_chains
1417        side_chain_beads_ids = []
1418        for side_chain_name in side_chain_list:
1419            pmb_type = self._get_template_type(name=side_chain_name,
1420                                               allowed_types={"particle", "residue"})
1421            if pmb_type == 'particle':
1422                lj_parameters = self.get_lj_parameters(particle_name1=central_bead_name,
1423                                                       particle_name2=side_chain_name)
1424                bond_tpl = self.get_bond_template(particle_name1=central_bead_name,
1425                                                  particle_name2=side_chain_name,
1426                                                  use_default_bond=use_default_bond)
1427                l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1428                                                      bond_type=bond_tpl.bond_type,
1429                                                      bond_parameters=bond_tpl.get_parameters(ureg=self.units))               
1430                if backbone_vector is None:
1431                    bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1432                                                                radius=l0, 
1433                                                                n_samples=1,
1434                                                                on_surface=True)[0]
1435                else:
1436                    bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector),
1437                                                                                                magnitude=l0)
1438                    
1439                side_bead_id = self.create_particle(name=side_chain_name, 
1440                                                    espresso_system=espresso_system,
1441                                                    position=[bead_position], 
1442                                                    number_of_particles=1)[0]
1443                side_chain_beads_ids.append(side_bead_id)
1444                self.db._update_instance(pmb_type="particle",
1445                                         instance_id=side_bead_id,
1446                                         attribute="residue_id",
1447                                         value=residue_id)
1448                self.create_bond(particle_id1=central_bead_id,
1449                                 particle_id2=side_bead_id,
1450                                 espresso_system=espresso_system,
1451                                 use_default_bond=use_default_bond)
1452            elif pmb_type == 'residue':
1453                side_residue_tpl = self.db.get_template(name=side_chain_name,
1454                                                        pmb_type=pmb_type)
1455                central_bead_side_chain = side_residue_tpl.central_bead
1456                lj_parameters = self.get_lj_parameters(particle_name1=central_bead_name,
1457                                                       particle_name2=central_bead_side_chain)
1458                bond_tpl = self.get_bond_template(particle_name1=central_bead_name,
1459                                                  particle_name2=central_bead_side_chain,
1460                                                  use_default_bond=use_default_bond)
1461                l0 = hf.calculate_initial_bond_length(lj_parameters=lj_parameters,
1462                                                      bond_type=bond_tpl.bond_type,
1463                                                      bond_parameters=bond_tpl.get_parameters(ureg=self.units))
1464                if backbone_vector is None:
1465                    residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1466                                                                radius=l0, 
1467                                                                n_samples=1,
1468                                                                on_surface=True)[0]
1469                else:
1470                    residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1471                                                                                                    magnitude=l0)
1472                side_residue_id = self.create_residue(name=side_chain_name, 
1473                                                      espresso_system=espresso_system,
1474                                                      central_bead_position=[residue_position],
1475                                                      use_default_bond=use_default_bond)
1476                # Find particle ids of the inner residue
1477                side_chain_beads_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
1478                                                                               attribute="residue_id",
1479                                                                               value=side_residue_id)
1480                # Change the residue_id of the residue in the side chain to the one of the outer residue
1481                for particle_id in side_chain_beads_ids:
1482                    self.db._update_instance(instance_id=particle_id,
1483                                             pmb_type="particle",
1484                                             attribute="residue_id",
1485                                             value=residue_id)
1486                # Remove the instance of the inner residue
1487                self.db.delete_instance(pmb_type="residue",
1488                                        instance_id=side_residue_id)
1489                self.create_bond(particle_id1=central_bead_id,
1490                                 particle_id2=side_chain_beads_ids[0],
1491                                 espresso_system=espresso_system,
1492                                 use_default_bond=use_default_bond)
1493        if gen_angle:
1494            self._generate_angles_for_entity(espresso_system=espresso_system,
1495                                             entity_id=residue_id,
1496                                             entity_id_col="residue_id")
1497        return  residue_id

Creates a residue into ESPResSo.

Arguments:
  • name ('str'): Label of the residue type to be created.
  • espresso_system ('espressomd.system.System'): Instance of a system object from espressomd library.
  • central_bead_position ('list' of 'float'): Position of the central bead.
  • use_default_bond ('bool'): Switch to control if a bond of type 'default' is used to bond a particle whose bonds types are not defined in the pyMBE database.
  • backbone_vector ('list' of 'float'): Backbone vector of the molecule. All side chains are created perpendicularly to 'backbone_vector'.
Returns:

(int): residue_id of the residue created.

def define_bond(self, bond_type, bond_parameters, particle_pairs):
1499    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1500        """
1501        Defines bond templates for each particle pair in 'particle_pairs' in the pyMBE database.
1502
1503        Args:
1504            bond_type ('str'): 
1505                label to identify the potential to model the bond.
1506
1507            bond_parameters ('dict'): 
1508                parameters of the potential of the bond.
1509
1510            particle_pairs ('lst'): 
1511                list of the 'names' of the 'particles' to be bonded.
1512
1513        Notes:
1514            -Currently, only HARMONIC and FENE bonds are supported.
1515            - For a HARMONIC bond the dictionary must contain the following parameters:
1516                - k ('pint.Quantity')      : Magnitude of the bond. It should have units of energy/length**2 
1517                using the 'pmb.units' UnitRegistry.
1518                - r_0 ('pint.Quantity')    : Equilibrium bond length. It should have units of length using 
1519                the 'pmb.units' UnitRegistry.
1520           - For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:              
1521                - d_r_max ('pint.Quantity'): Maximal stretching length for FENE. It should have 
1522                units of length using the 'pmb.units' UnitRegistry. Default 'None'.
1523        """
1524        self._check_bond_inputs(bond_parameters=bond_parameters,
1525                                bond_type=bond_type)
1526        parameters_expected_dimensions={"r_0": "length",
1527                                        "k": "energy/length**2",
1528                                        "d_r_max": "length"}
1529
1530        parameters_tpl = {}
1531        for key in bond_parameters.keys():
1532            parameters_tpl[key]= PintQuantity.from_quantity(q=bond_parameters[key],
1533                                                            expected_dimension=parameters_expected_dimensions[key],
1534                                                            ureg=self.units)
1535
1536        bond_names=[]
1537        for particle_name1, particle_name2 in particle_pairs:
1538            
1539            tpl = BondTemplate(particle_name1=particle_name1,
1540                               particle_name2=particle_name2,
1541                               parameters=parameters_tpl,
1542                               bond_type=bond_type)
1543            tpl._make_name()
1544            if tpl.name in bond_names:
1545                raise RuntimeError(f"Bond {tpl.name} has already been defined, please check the list of particle pairs")
1546            bond_names.append(tpl.name)
1547            self.db._register_template(tpl)

Defines bond templates for each particle pair in 'particle_pairs' in the pyMBE database.

Arguments:
  • bond_type ('str'): label to identify the potential to model the bond.
  • bond_parameters ('dict'): parameters of the potential of the bond.
  • particle_pairs ('lst'): list of the 'names' of the 'particles' to be bonded.
Notes:

-Currently, only HARMONIC and FENE bonds are supported.

  • For a HARMONIC bond the dictionary must contain the following parameters:
    • k ('pint.Quantity') : Magnitude of the bond. It should have units of energy/length**2 using the 'pmb.units' UnitRegistry.
    • r_0 ('pint.Quantity') : Equilibrium bond length. It should have units of length using the 'pmb.units' UnitRegistry.
    • For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:
    • d_r_max ('pint.Quantity'): Maximal stretching length for FENE. It should have units of length using the 'pmb.units' UnitRegistry. Default 'None'.
def define_default_bond(self, bond_type, bond_parameters):
1550    def define_default_bond(self, bond_type, bond_parameters):
1551        """
1552        Defines a bond template as a "default" template in the pyMBE database.
1553        
1554        Args:
1555            bond_type ('str'): 
1556                label to identify the potential to model the bond.
1557
1558            bond_parameters ('dict'): 
1559                parameters of the potential of the bond.
1560            
1561        Notes:
1562            - Currently, only harmonic and FENE bonds are supported. 
1563        """
1564        self._check_bond_inputs(bond_parameters=bond_parameters,
1565                                bond_type=bond_type)
1566        parameters_expected_dimensions={"r_0": "length",
1567                                        "k": "energy/length**2",
1568                                        "d_r_max": "length"}
1569        parameters_tpl = {}
1570        for key in bond_parameters.keys():
1571            parameters_tpl[key]= PintQuantity.from_quantity(q=bond_parameters[key],
1572                                                            expected_dimension=parameters_expected_dimensions[key],
1573                                                            ureg=self.units)
1574        tpl = BondTemplate(parameters=parameters_tpl,
1575                               bond_type=bond_type)
1576        tpl.name = "default"
1577        self.db._register_template(tpl)

Defines a bond template as a "default" template in the pyMBE database.

Arguments:
  • bond_type ('str'): label to identify the potential to model the bond.
  • bond_parameters ('dict'): parameters of the potential of the bond.
Notes:
  • Currently, only harmonic and FENE bonds are supported.
def define_angular_potential(self, angle_type, angle_parameters, particle_triplets):
1579    def define_angular_potential(self, angle_type, angle_parameters, particle_triplets):
1580        """
1581        Defines angle potential templates for each particle triplet in `particle_triplets`.
1582
1583        Args:
1584            angle_type ('str'):
1585                Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
1586
1587            angle_parameters ('dict'):
1588                Parameters of the angle potential. Must contain:
1589                    - "k" ('pint.Quantity'): Bending stiffness with dimensions of energy.
1590                    - "phi_0" ('float'): Equilibrium angle in radians.
1591
1592            particle_triplets ('list[tuple[str,str,str]]'):
1593                List of (side_particle1, central_particle, side_particle2) triplets.
1594        """
1595        valid_angle_types = ["harmonic", "cosine", "harmonic_cosine"]
1596        if angle_type not in valid_angle_types:
1597            raise NotImplementedError(f"Angle potential type '{angle_type}' currently not implemented in pyMBE, accepted types are {valid_angle_types}")
1598
1599        if "k" not in angle_parameters:
1600            raise ValueError("Magnitude of the angle potential (k) is missing")
1601        if "phi_0" not in angle_parameters:
1602            raise ValueError("Equilibrium angle (phi_0) is missing")
1603
1604        parameters_tpl = {
1605            "k": PintQuantity.from_quantity(q=angle_parameters["k"],
1606                                            expected_dimension="energy",
1607                                            ureg=self.units),
1608            "phi_0": PintQuantity.from_quantity(q=angle_parameters["phi_0"],
1609                                                expected_dimension="dimensionless",
1610                                                ureg=self.units),
1611        }
1612
1613        angle_names = []
1614        for side1, central, side2 in particle_triplets:
1615            tpl = AngleTemplate(side_particle1=side1,
1616                                central_particle=central,
1617                                side_particle2=side2,
1618                                parameters=parameters_tpl,
1619                                angle_type=angle_type)
1620            tpl._make_name()
1621            if tpl.name in angle_names:
1622                raise RuntimeError(f"Angle {tpl.name} has already been defined, please check the list of particle triplets")
1623            angle_names.append(tpl.name)
1624            self.db._register_template(tpl)

Defines angle potential templates for each particle triplet in particle_triplets.

Arguments:
  • angle_type ('str'): Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
  • angle_parameters ('dict'): Parameters of the angle potential. Must contain:
    • "k" ('pint.Quantity'): Bending stiffness with dimensions of energy.
    • "phi_0" ('float'): Equilibrium angle in radians.
  • particle_triplets ('list[tuple[str,str,str]]'): List of (side_particle1, central_particle, side_particle2) triplets.
def define_default_angular_potential(self, angle_type, angle_parameters):
1626    def define_default_angular_potential(self, angle_type, angle_parameters):
1627        """
1628        Defines an angle template as a "default" template in the pyMBE database.
1629
1630        Args:
1631            angle_type ('str'):
1632                Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
1633
1634            angle_parameters ('dict'):
1635                Parameters of the angle potential (k, phi_0).
1636        """
1637        valid_angle_types = ["harmonic", "cosine", "harmonic_cosine"]
1638        if angle_type not in valid_angle_types:
1639            raise NotImplementedError(f"Angle potential type '{angle_type}' currently not implemented in pyMBE, accepted types are {valid_angle_types}")
1640
1641        if "k" not in angle_parameters:
1642            raise ValueError("Magnitude of the angle potential (k) is missing")
1643        if "phi_0" not in angle_parameters:
1644            raise ValueError("Equilibrium angle (phi_0) is missing")
1645
1646        parameters_tpl = {
1647            "k": PintQuantity.from_quantity(q=angle_parameters["k"],
1648                                            expected_dimension="energy",
1649                                            ureg=self.units),
1650            "phi_0": PintQuantity.from_quantity(q=angle_parameters["phi_0"],
1651                                                expected_dimension="dimensionless",
1652                                                ureg=self.units),
1653        }
1654        tpl = AngleTemplate(parameters=parameters_tpl,
1655                            angle_type=angle_type)
1656        tpl.name = "default"
1657        self.db._register_template(tpl)

Defines an angle template as a "default" template in the pyMBE database.

Arguments:
  • angle_type ('str'): Type of angle potential. Supported: "harmonic", "cosine", "harmonic_cosine".
  • angle_parameters ('dict'): Parameters of the angle potential (k, phi_0).
def create_angular_potential( self, particle_id1, particle_id2, particle_id3, espresso_system, use_default_angle=False):
1659    def create_angular_potential(self, particle_id1, particle_id2, particle_id3, espresso_system, use_default_angle=False):
1660        """
1661        Creates an angle between three particle instances in an ESPResSo system
1662        and registers it in the pyMBE database.
1663
1664        Args:
1665            particle_id1 ('int'): ID of the first side particle.
1666            particle_id2 ('int'): ID of the central particle.
1667            particle_id3 ('int'): ID of the second side particle.
1668            espresso_system ('espressomd.system.System'): ESPResSo system.
1669            use_default_angle ('bool', optional): If True, use the default angle if no specific one is found.
1670        """
1671        particle_inst_1 = self.db.get_instance(pmb_type="particle", instance_id=particle_id1)
1672        particle_inst_2 = self.db.get_instance(pmb_type="particle", instance_id=particle_id2)
1673        particle_inst_3 = self.db.get_instance(pmb_type="particle", instance_id=particle_id3)
1674
1675        # Verify that bonds exist between side particles and central particle
1676        bond_instances = self.db.get_instances(pmb_type="bond")
1677        bonded_pairs = set()
1678        for bond in bond_instances.values():
1679            pair = frozenset([bond.particle_id1, bond.particle_id2])
1680            bonded_pairs.add(pair)
1681        if frozenset([particle_id1, particle_id2]) not in bonded_pairs:
1682            raise ValueError(f"Cannot create angle: no bond exists between particle {particle_id1} and central particle {particle_id2}.")
1683        if frozenset([particle_id3, particle_id2]) not in bonded_pairs:
1684            raise ValueError(f"Cannot create angle: no bond exists between particle {particle_id3} and central particle {particle_id2}.")
1685
1686        angle_tpl = self.get_angle_template(side_name1=particle_inst_1.name,
1687                                            central_name=particle_inst_2.name,
1688                                            side_name2=particle_inst_3.name,
1689                                            use_default_angle=use_default_angle)
1690        angle_inst = self._get_espresso_angle_instance(angle_template=angle_tpl, espresso_system=espresso_system)
1691
1692        # ESPResSo angle bonds are added to the central particle
1693        espresso_system.part.by_id(particle_id2).add_bond((angle_inst, particle_id1, particle_id3))
1694
1695        angle_id = self.db._propose_instance_id(pmb_type="angle")
1696        pmb_angle_instance = AngleInstance(angle_id=angle_id,
1697                                           name=angle_tpl.name,
1698                                           particle_id1=particle_id1,
1699                                           particle_id2=particle_id2,
1700                                           particle_id3=particle_id3)
1701        self.db._register_instance(instance=pmb_angle_instance)

Creates an angle between three particle instances in an ESPResSo system and registers it in the pyMBE database.

Arguments:
  • particle_id1 ('int'): ID of the first side particle.
  • particle_id2 ('int'): ID of the central particle.
  • particle_id3 ('int'): ID of the second side particle.
  • espresso_system ('espressomd.system.System'): ESPResSo system.
  • use_default_angle ('bool', optional): If True, use the default angle if no specific one is found.
def get_angle_template(self, side_name1, central_name, side_name2, use_default_angle=False):
1703    def get_angle_template(self, side_name1, central_name, side_name2, use_default_angle=False):
1704        """
1705        Retrieves an angle template connecting three particle templates.
1706
1707        Args:
1708            side_name1 ('str'): Name of the first side particle.
1709            central_name ('str'): Name of the central particle.
1710            side_name2 ('str'): Name of the second side particle.
1711            use_default_angle ('bool', optional): If True, fall back to the default angle template.
1712
1713        Returns:
1714            ('AngleTemplate'): The matching angle template.
1715        """
1716        angle_key = AngleTemplate.make_angle_key(side1=side_name1, central=central_name, side2=side_name2)
1717        try:
1718            return self.db.get_template(name=angle_key, pmb_type="angle")
1719        except ValueError:
1720            pass
1721
1722        if use_default_angle:
1723            return self.db.get_template(name="default", pmb_type="angle")
1724
1725        raise ValueError(f"No angle template found for '{side_name1}-{central_name}-{side_name2}', and default angles are deactivated.")

Retrieves an angle template connecting three particle templates.

Arguments:
  • side_name1 ('str'): Name of the first side particle.
  • central_name ('str'): Name of the central particle.
  • side_name2 ('str'): Name of the second side particle.
  • use_default_angle ('bool', optional): If True, fall back to the default angle template.
Returns:

('AngleTemplate'): The matching angle template.

def define_hydrogel(self, name, node_map, chain_map):
1818    def define_hydrogel(self, name, node_map, chain_map):
1819        """
1820        Defines a hydrogel template in the pyMBE database.
1821
1822        Args:
1823            name ('str'): 
1824                Unique label that identifies the 'hydrogel'.
1825
1826            node_map ('list of dict'): 
1827                [{"particle_name": , "lattice_index": }, ... ]
1828
1829            chain_map ('list of dict'): 
1830                [{"node_start": , "node_end": , "residue_list": , ... ]
1831        """
1832        # Sanity tests
1833        node_indices = {tuple(entry['lattice_index']) for entry in node_map}                
1834        chain_map_connectivity = set()
1835        for entry in chain_map:
1836            start = self.lattice_builder.node_labels[entry['node_start']]
1837            end = self.lattice_builder.node_labels[entry['node_end']]
1838            chain_map_connectivity.add((start,end))
1839        if self.lattice_builder.lattice.connectivity != chain_map_connectivity:
1840            raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs")
1841        diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices}
1842        if node_indices != diamond_indices:
1843            raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ")
1844        # Register information in the pyMBE database
1845        nodes=[]
1846        for entry in node_map:
1847            nodes.append(HydrogelNode(particle_name=entry["particle_name"],
1848                                      lattice_index=entry["lattice_index"]))
1849        chains=[]
1850        for chain in chain_map:
1851            chains.append(HydrogelChain(node_start=chain["node_start"],
1852                                        node_end=chain["node_end"],
1853                                        molecule_name=chain["molecule_name"]))
1854        tpl = HydrogelTemplate(name=name,
1855                               node_map=nodes,
1856                               chain_map=chains)
1857        self.db._register_template(tpl)

Defines a hydrogel template in the pyMBE database.

Arguments:
  • name ('str'): Unique label that identifies the 'hydrogel'.
  • node_map ('list of dict'): [{"particle_name": , "lattice_index": }, ... ]
  • chain_map ('list of dict'): [{"node_start": , "node_end": , "residue_list": , ... ]
def define_molecule(self, name, residue_list):
1859    def define_molecule(self, name, residue_list):
1860        """
1861        Defines a molecule template in the pyMBE database.
1862
1863        Args:
1864            name('str'): 
1865                Unique label that identifies the 'molecule'.
1866
1867            residue_list ('list' of 'str'): 
1868                List of the 'name's of the 'residue's  in the sequence of the 'molecule'.  
1869        """
1870        tpl = MoleculeTemplate(name=name,
1871                               residue_list=residue_list)
1872        self.db._register_template(tpl)

Defines a molecule template in the pyMBE database.

Arguments:
  • name('str'): Unique label that identifies the 'molecule'.
  • residue_list ('list' of 'str'): List of the 'name's of the 'residue's in the sequence of the 'molecule'.
def define_monoprototic_acidbase_reaction(self, particle_name, pka, acidity, metadata=None):
1874    def define_monoprototic_acidbase_reaction(self, particle_name, pka, acidity, metadata=None):
1875        """
1876        Defines an acid-base reaction for a monoprototic particle in the pyMBE database.
1877
1878        Args:
1879            particle_name ('str'): 
1880                Unique label that identifies the particle template. 
1881
1882            pka ('float'): 
1883                pka-value of the acid or base.
1884
1885            acidity ('str'): 
1886                Identifies whether if the particle is 'acidic' or 'basic'.
1887
1888            metadata ('dict', optional): 
1889                Additional information to be stored in the reaction. Defaults to None.
1890        """
1891        supported_acidities = ["acidic", "basic"]
1892        if acidity not in supported_acidities:
1893            raise ValueError(f"Unsupported acidity '{acidity}' for particle '{particle_name}'. Supported acidities are {supported_acidities}.")
1894        reaction_type = "monoprotic"
1895        if acidity == "basic":
1896            reaction_type += "_base"
1897        else:
1898            reaction_type += "_acid"
1899        reaction = Reaction(participants=[ReactionParticipant(particle_name=particle_name,
1900                                                              state_name=f"{particle_name}H", 
1901                                                              coefficient=-1),
1902                                          ReactionParticipant(particle_name=particle_name,
1903                                                              state_name=f"{particle_name}",
1904                                                              coefficient=1)],
1905                            reaction_type=reaction_type,
1906                            pK=pka,
1907                            metadata=metadata)
1908        self.db._register_reaction(reaction)

Defines an acid-base reaction for a monoprototic particle in the pyMBE database.

Arguments:
  • particle_name ('str'): Unique label that identifies the particle template.
  • pka ('float'): pka-value of the acid or base.
  • acidity ('str'): Identifies whether if the particle is 'acidic' or 'basic'.
  • metadata ('dict', optional): Additional information to be stored in the reaction. Defaults to None.
def define_monoprototic_particle_states(self, particle_name, acidity):
1910    def define_monoprototic_particle_states(self, particle_name, acidity):
1911        """
1912        Defines particle states for a monoprotonic particle template including the charges in each of its possible states. 
1913
1914        Args:
1915            particle_name ('str'): 
1916                Unique label that identifies the particle template. 
1917
1918            acidity ('str'): 
1919                Identifies whether the particle is 'acidic' or 'basic'.
1920        """
1921        acidity_valid_keys = ['acidic', 'basic']
1922        if not pd.isna(acidity):
1923            if acidity not in acidity_valid_keys:
1924                raise ValueError(f"Acidity {acidity} provided for particle name  {particle_name} is not supported. Valid keys are: {acidity_valid_keys}")
1925        if acidity == "acidic":
1926            states = [{"name": f"{particle_name}H", "z": 0}, 
1927                      {"name": f"{particle_name}",  "z": -1}]
1928            
1929        elif acidity == "basic":
1930            states = [{"name": f"{particle_name}H", "z": 1}, 
1931                      {"name": f"{particle_name}",  "z": 0}]
1932        self.define_particle_states(particle_name=particle_name, 
1933                                    states=states)

Defines particle states for a monoprotonic particle template including the charges in each of its possible states.

Arguments:
  • particle_name ('str'): Unique label that identifies the particle template.
  • acidity ('str'): Identifies whether the particle is 'acidic' or 'basic'.
def define_particle( self, name, sigma, epsilon, z=0, acidity=<NA>, pka=<NA>, cutoff=<NA>, offset=<NA>):
1935    def define_particle(self, name,  sigma, epsilon, z=0, acidity=pd.NA, pka=pd.NA, cutoff=pd.NA, offset=pd.NA):
1936        """
1937        Defines a particle template in the pyMBE database.
1938
1939        Args:
1940            name('str'):
1941                 Unique label that identifies this particle type.  
1942
1943            sigma('pint.Quantity'): 
1944                Sigma parameter used to set up Lennard-Jones interactions for this particle type. 
1945
1946            epsilon('pint.Quantity'): 
1947                Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe.
1948
1949            z('int', optional): 
1950                Permanent charge number of this particle type. Defaults to 0.
1951
1952            acidity('str', optional): 
1953                Identifies whether if the particle is 'acidic' or 'basic', used to setup constant pH simulations. Defaults to pd.NA.
1954
1955            pka('float', optional):
1956                If 'particle' is an acid or a base, it defines its  pka-value. Defaults to pd.NA.
1957
1958            cutoff('pint.Quantity', optional): 
1959                Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1960
1961            offset('pint.Quantity', optional): 
1962                Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1963            
1964        Notes:
1965            - 'sigma', 'cutoff' and 'offset' must have a dimensitonality of '[length]' and should be defined using pmb.units.
1966            - 'epsilon' must have a dimensitonality of '[energy]' and should be defined using pmb.units.
1967            - 'cutoff' defaults to '2**(1./6.) reduced_length'. 
1968            - 'offset' defaults to 0.
1969            - For more information on 'sigma', 'epsilon', 'cutoff' and 'offset' check 'pmb.setup_lj_interactions()'.
1970        """ 
1971        # If 'cutoff' and 'offset' are not defined, default them to the following values
1972        if pd.isna(cutoff):
1973            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1974        if pd.isna(offset):
1975            offset=self.units.Quantity(0, "reduced_length")
1976        # Define particle states
1977        if acidity is pd.NA:
1978            states = [{"name": f"{name}",  "z": z}]
1979            self.define_particle_states(particle_name=name, 
1980                                        states=states)
1981            initial_state = name
1982        else:
1983            self.define_monoprototic_particle_states(particle_name=name,
1984                                                  acidity=acidity)
1985            initial_state = f"{name}H"
1986            if pka is not pd.NA:
1987                self.define_monoprototic_acidbase_reaction(particle_name=name,
1988                                                           acidity=acidity,
1989                                                           pka=pka)
1990        tpl = ParticleTemplate(name=name, 
1991                               sigma=PintQuantity.from_quantity(q=sigma, expected_dimension="length", ureg=self.units), 
1992                               epsilon=PintQuantity.from_quantity(q=epsilon, expected_dimension="energy", ureg=self.units),
1993                               cutoff=PintQuantity.from_quantity(q=cutoff, expected_dimension="length", ureg=self.units), 
1994                               offset=PintQuantity.from_quantity(q=offset, expected_dimension="length", ureg=self.units),
1995                               initial_state=initial_state)
1996        self.db._register_template(tpl)

Defines a particle template in the pyMBE database.

Arguments:
  • name('str'): Unique label that identifies this particle type.
  • sigma('pint.Quantity'): Sigma parameter used to set up Lennard-Jones interactions for this particle type.
  • epsilon('pint.Quantity'): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe.
  • z('int', optional): Permanent charge number of this particle type. Defaults to 0.
  • acidity('str', optional): Identifies whether if the particle is 'acidic' or 'basic', used to setup constant pH simulations. Defaults to pd.NA.
  • pka('float', optional): If 'particle' is an acid or a base, it defines its pka-value. Defaults to pd.NA.
  • cutoff('pint.Quantity', optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
  • offset('pint.Quantity', optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
Notes:
  • 'sigma', 'cutoff' and 'offset' must have a dimensitonality of '[length]' and should be defined using pmb.units.
  • 'epsilon' must have a dimensitonality of '[energy]' and should be defined using pmb.units.
  • 'cutoff' defaults to '2**(1./6.) reduced_length'.
  • 'offset' defaults to 0.
  • For more information on 'sigma', 'epsilon', 'cutoff' and 'offset' check 'pmb.setup_lj_interactions()'.
def define_particle_states(self, particle_name, states):
1998    def define_particle_states(self, particle_name, states):
1999        """
2000        Define the chemical states of an existing particle template.
2001
2002        Args:
2003            particle_name ('str'):
2004                Name of a particle template. 
2005
2006            states ('list' of 'dict'):
2007                List of dictionaries defining the particle states. Each dictionary
2008                must contain:
2009                - 'name' ('str'): Name of the particle state (e.g. '"H"', '"-"',
2010                '"neutral"').
2011                - 'z' ('int'): Charge number of the particle in this state.
2012                Example:
2013                states = [{"name": "AH", "z": 0},     # protonated
2014                         {"name": "A-", "z": -1}]    # deprotonated
2015        Notes:
2016            - Each state is assigned a unique Espresso 'es_type' automatically.
2017            - Chemical reactions (e.g. acid–base equilibria) are **not** created by
2018            this method and must be defined separately (e.g. via
2019            'set_particle_acidity()' or custom reaction definitions).
2020            - Particles without explicitly defined states are assumed to have a
2021            single, implicit state with their default charge.
2022        """
2023        for s in states:
2024            state = ParticleStateTemplate(particle_name=particle_name,
2025                                          name=s["name"],
2026                                          z=s["z"],
2027                                          es_type=self.propose_unused_type())
2028            self.db._register_template(state)

Define the chemical states of an existing particle template.

Arguments:
  • particle_name ('str'): Name of a particle template.
  • states ('list' of 'dict'): List of dictionaries defining the particle states. Each dictionary must contain:
    • 'name' ('str'): Name of the particle state (e.g. '"H"', '"-"', '"neutral"').
    • 'z' ('int'): Charge number of the particle in this state. Example: states = [{"name": "AH", "z": 0}, # protonated {"name": "A-", "z": -1}] # deprotonated
Notes:
  • Each state is assigned a unique Espresso 'es_type' automatically.
  • Chemical reactions (e.g. acid–base equilibria) are not created by this method and must be defined separately (e.g. via 'set_particle_acidity()' or custom reaction definitions).
  • Particles without explicitly defined states are assumed to have a single, implicit state with their default charge.
def define_peptide(self, name, sequence, model):
2030    def define_peptide(self, name, sequence, model):
2031        """
2032        Defines a peptide template in the pyMBE database.
2033
2034        Args:
2035            name ('str'): 
2036                Unique label that identifies the peptide.
2037
2038            sequence ('str'): 
2039                Sequence of the peptide.
2040
2041            model ('str'): 
2042                Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
2043        """
2044        valid_keys = ['1beadAA','2beadAA']
2045        if model not in valid_keys:
2046            raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
2047        clean_sequence = hf.protein_sequence_parser(sequence=sequence)    
2048        residue_list = self._get_residue_list_from_sequence(sequence=clean_sequence)
2049        tpl = PeptideTemplate(name=name,
2050                            residue_list=residue_list,
2051                            model=model,
2052                            sequence=sequence)
2053        self.db._register_template(tpl)        

Defines a peptide template in the pyMBE database.

Arguments:
  • name ('str'): Unique label that identifies the peptide.
  • sequence ('str'): Sequence of the peptide.
  • model ('str'): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
def define_protein(self, name, sequence, model):
2055    def define_protein(self, name, sequence, model):
2056        """
2057        Defines a protein template in the pyMBE database.
2058
2059        Args:
2060            name ('str'): 
2061                Unique label that identifies the protein.
2062
2063            sequence ('str'): 
2064                Sequence of the protein.
2065
2066            model ('string'): 
2067                Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
2068
2069        Notes:
2070            - Currently, only 'lj_setup_mode="wca"' is supported. This corresponds to setting up the WCA potential.
2071        """
2072        valid_model_keys = ['1beadAA','2beadAA']
2073        if model not in valid_model_keys:
2074            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
2075        
2076        residue_list = self._get_residue_list_from_sequence(sequence=sequence)
2077        tpl = ProteinTemplate(name=name,
2078                              model=model,
2079                              residue_list=residue_list,
2080                              sequence=sequence)
2081        self.db._register_template(tpl)

Defines a protein template in the pyMBE database.

Arguments:
  • name ('str'): Unique label that identifies the protein.
  • sequence ('str'): Sequence of the protein.
  • model ('string'): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
Notes:
  • Currently, only 'lj_setup_mode="wca"' is supported. This corresponds to setting up the WCA potential.
def define_residue(self, name, central_bead, side_chains):
2083    def define_residue(self, name, central_bead, side_chains):
2084        """
2085        Defines a residue template in the pyMBE database.
2086
2087        Args:
2088            name ('str'): 
2089                Unique label that identifies the residue.
2090
2091            central_bead ('str'): 
2092                'name' of the 'particle' to be placed as central_bead of the residue.
2093
2094            side_chains('list' of 'str'): 
2095                List of 'name's of the pmb_objects to be placed as side_chains of the residue. Currently, only pyMBE objects of type 'particle' or 'residue' are supported.
2096        """
2097        tpl = ResidueTemplate(name=name,
2098                              central_bead=central_bead,
2099                              side_chains=side_chains)
2100        self.db._register_template(tpl)

Defines a residue template in the pyMBE database.

Arguments:
  • name ('str'): Unique label that identifies the residue.
  • central_bead ('str'): 'name' of the 'particle' to be placed as central_bead of the residue.
  • side_chains('list' of 'str'): List of 'name's of the pmb_objects to be placed as side_chains of the residue. Currently, only pyMBE objects of type 'particle' or 'residue' are supported.
def delete_instances_in_system(self, instance_id, pmb_type, espresso_system):
2102    def delete_instances_in_system(self, instance_id, pmb_type, espresso_system):
2103        """
2104        Deletes the instance with instance_id from the ESPResSo system. 
2105        Related assembly, molecule, residue, particles and bond instances will also be deleted from the pyMBE dataframe.
2106
2107        Args:
2108            instance_id ('int'): 
2109                id of the assembly to be deleted. 
2110
2111            pmb_type ('str'): 
2112                the instance type to be deleted. 
2113
2114            espresso_system ('espressomd.system.System'): 
2115                Instance of a system class from espressomd library.
2116        """
2117        if pmb_type == "particle":
2118            instance_identifier = "particle_id"
2119        elif pmb_type == "residue":
2120            instance_identifier = "residue_id"
2121        elif pmb_type in self.db._molecule_like_types:
2122            instance_identifier = "molecule_id"
2123        elif pmb_type in self.db._assembly_like_types:
2124            instance_identifier = "assembly_id"
2125        particle_ids = self.db._find_instance_ids_by_attribute(pmb_type="particle",
2126                                                               attribute=instance_identifier,
2127                                                               value=instance_id)
2128        self._delete_particles_from_espresso(particle_ids=particle_ids,
2129                                             espresso_system=espresso_system)
2130        self.db.delete_instance(pmb_type=pmb_type,
2131                                instance_id=instance_id)

Deletes the instance with instance_id from the ESPResSo system. Related assembly, molecule, residue, particles and bond instances will also be deleted from the pyMBE dataframe.

Arguments:
  • instance_id ('int'): id of the assembly to be deleted.
  • pmb_type ('str'): the instance type to be deleted.
  • espresso_system ('espressomd.system.System'): Instance of a system class from espressomd library.
def determine_reservoir_concentrations( self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
2133    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
2134        """
2135        Determines ionic concentrations in the reservoir at fixed pH and salt concentration.
2136
2137        Args:
2138            pH_res ('float'):
2139                Target pH value in the reservoir.
2140
2141            c_salt_res ('pint.Quantity'):
2142                Concentration of monovalent salt (e.g., NaCl) in the reservoir.
2143
2144            activity_coefficient_monovalent_pair ('callable'):
2145                Function returning the activity coefficient of a monovalent ion pair
2146                as a function of ionic strength:
2147                'gamma = activity_coefficient_monovalent_pair(I)'.
2148
2149            max_number_sc_runs ('int', optional):
2150                Maximum number of self-consistent iterations allowed before
2151                convergence is enforced. Defaults to 200.
2152
2153        Returns:
2154            tuple:
2155                (cH_res, cOH_res, cNa_res, cCl_res)
2156                - cH_res ('pint.Quantity'): Concentration of H⁺ ions.
2157                - cOH_res ('pint.Quantity'): Concentration of OH⁻ ions.
2158                - cNa_res ('pint.Quantity'): Concentration of Na⁺ ions.
2159                - cCl_res ('pint.Quantity'): Concentration of Cl⁻ ions.
2160
2161        Notess:
2162            - The algorithm enforces electroneutrality in the reservoir.
2163            - Water autodissociation is included via the equilibrium constant 'Kw'.
2164            - Non-ideal effects enter through activity coefficients depending on
2165            ionic strength.
2166            - The implementation follows the self-consistent scheme described in
2167            Landsgesell (PhD thesis, Sec. 5.3, doi:10.18419/opus-10831), adapted
2168            from the original code (doi:10.18419/darus-2237).
2169        """
2170        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
2171            """
2172            Iteratively determines reservoir ion concentrations self-consistently.
2173
2174            Args:
2175                cH_res ('pint.Quantity'):
2176                    Current estimate of the H⁺ concentration.
2177                c_salt_res ('pint.Quantity'):
2178                    Concentration of monovalent salt in the reservoir.
2179
2180            Returns:
2181                'tuple':
2182                    (cH_res, cOH_res, cNa_res, cCl_res)
2183            """
2184            # Initial ideal estimate
2185            cOH_res = self.Kw / cH_res
2186            if cOH_res >= cH_res:
2187                cNa_res = c_salt_res + (cOH_res - cH_res)
2188                cCl_res = c_salt_res
2189            else:
2190                cCl_res = c_salt_res + (cH_res - cOH_res)
2191                cNa_res = c_salt_res
2192            # Self-consistent iteration
2193            for _ in range(max_number_sc_runs):
2194                ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2195                cOH_new = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
2196                if cOH_new >= cH_res:
2197                    cNa_new = c_salt_res + (cOH_new - cH_res)
2198                    cCl_new = c_salt_res
2199                else:
2200                    cCl_new = c_salt_res + (cH_res - cOH_new)
2201                    cNa_new = c_salt_res
2202                # Update values
2203                cOH_res = cOH_new
2204                cNa_res = cNa_new
2205                cCl_res = cCl_new
2206            return cH_res, cOH_res, cNa_res, cCl_res
2207        # Initial guess for H+ concentration from target pH
2208        cH_res = 10 ** (-pH_res) * self.units.mol / self.units.l
2209        # First self-consistent solve
2210        cH_res, cOH_res, cNa_res, cCl_res = (determine_reservoir_concentrations_selfconsistently(cH_res, 
2211                                                                                                 c_salt_res))
2212        ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2213        determined_pH = -np.log10(cH_res.to("mol/L").magnitude* np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
2214        # Outer loop to enforce target pH
2215        while abs(determined_pH - pH_res) > 1e-6:
2216            if determined_pH > pH_res:
2217                cH_res *= 1.005
2218            else:
2219                cH_res /= 1.003
2220            cH_res, cOH_res, cNa_res, cCl_res = (determine_reservoir_concentrations_selfconsistently(cH_res, 
2221                                                                                                     c_salt_res))
2222            ionic_strength_res = 0.5 * (cNa_res + cCl_res + cOH_res + cH_res)
2223            determined_pH = -np.log10(cH_res.to("mol/L").magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
2224        return cH_res, cOH_res, cNa_res, cCl_res

Determines ionic concentrations in the reservoir at fixed pH and salt concentration.

Arguments:
  • pH_res ('float'): Target pH value in the reservoir.
  • c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g., NaCl) in the reservoir.
  • activity_coefficient_monovalent_pair ('callable'): Function returning the activity coefficient of a monovalent ion pair as a function of ionic strength: 'gamma = activity_coefficient_monovalent_pair(I)'.
  • max_number_sc_runs ('int', optional): Maximum number of self-consistent iterations allowed before convergence is enforced. Defaults to 200.
Returns:

tuple: (cH_res, cOH_res, cNa_res, cCl_res) - cH_res ('pint.Quantity'): Concentration of H⁺ ions. - cOH_res ('pint.Quantity'): Concentration of OH⁻ ions. - cNa_res ('pint.Quantity'): Concentration of Na⁺ ions. - cCl_res ('pint.Quantity'): Concentration of Cl⁻ ions.

Notess:
  • The algorithm enforces electroneutrality in the reservoir.
  • Water autodissociation is included via the equilibrium constant 'Kw'.
  • Non-ideal effects enter through activity coefficients depending on ionic strength.
  • The implementation follows the self-consistent scheme described in Landsgesell (PhD thesis, Sec. 5.3, doi:10.18419/opus-10831), adapted from the original code (doi:10.18419/darus-2237).
def enable_motion_of_rigid_object(self, instance_id, pmb_type, espresso_system):
2226    def enable_motion_of_rigid_object(self, instance_id, pmb_type, espresso_system):
2227        """
2228        Enables translational and rotational motion of a rigid pyMBE object instance
2229        in an ESPResSo system.This method creates a rigid-body center particle at the center of mass of
2230        the specified pyMBE object and attaches all constituent particles to it
2231        using ESPResSo virtual sites. The resulting rigid object can translate and
2232        rotate as a single body.
2233
2234        Args:
2235            instance_id ('int'):
2236                Instance ID of the pyMBE object whose rigid-body motion is enabled.
2237
2238            pmb_type ('str'):
2239                pyMBE object type of the instance (e.g. '"molecule"', '"peptide"',
2240                '"protein"', or any assembly-like type).
2241
2242            espresso_system ('espressomd.system.System'):
2243                ESPResSo system in which the rigid object is defined.
2244
2245        Notess:
2246            - This method requires ESPResSo to be compiled with the following
2247            features enabled:
2248                - '"VIRTUAL_SITES_RELATIVE"'
2249                - '"MASS"'
2250            - A new ESPResSo particle is created to represent the rigid-body center.
2251            - The mass of the rigid-body center is set to the number of particles
2252            belonging to the object.
2253            - The rotational inertia tensor is approximated from the squared
2254            distances of the particles to the center of mass.
2255        """
2256        logging.info('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
2257        inst = self.db.get_instance(pmb_type=pmb_type,
2258                                    instance_id=instance_id)
2259        label = self._get_label_id_map(pmb_type=pmb_type)
2260        particle_ids_list = self.get_particle_id_map(object_name=inst.name)[label][instance_id]
2261        center_of_mass = self.calculate_center_of_mass (instance_id=instance_id,
2262                                                        espresso_system=espresso_system,
2263                                                        pmb_type=pmb_type)
2264        rigid_object_center = espresso_system.part.add(pos=center_of_mass,
2265                                                        rotation=[True,True,True], 
2266                                                        type=self.propose_unused_type())
2267        rigid_object_center.mass = len(particle_ids_list)
2268        momI = 0
2269        for pid in particle_ids_list:
2270            momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
2271        rigid_object_center.rinertia = np.ones(3) * momI        
2272        for particle_id in particle_ids_list:
2273            pid = espresso_system.part.by_id(particle_id)
2274            pid.vs_auto_relate_to(rigid_object_center.id)

Enables translational and rotational motion of a rigid pyMBE object instance in an ESPResSo system.This method creates a rigid-body center particle at the center of mass of the specified pyMBE object and attaches all constituent particles to it using ESPResSo virtual sites. The resulting rigid object can translate and rotate as a single body.

Arguments:
  • instance_id ('int'): Instance ID of the pyMBE object whose rigid-body motion is enabled.
  • pmb_type ('str'): pyMBE object type of the instance (e.g. '"molecule"', '"peptide"', '"protein"', or any assembly-like type).
  • espresso_system ('espressomd.system.System'): ESPResSo system in which the rigid object is defined.
Notess:
  • This method requires ESPResSo to be compiled with the following features enabled:
    • '"VIRTUAL_SITES_RELATIVE"'
    • '"MASS"'
  • A new ESPResSo particle is created to represent the rigid-body center.
  • The mass of the rigid-body center is set to the number of particles belonging to the object.
  • The rotational inertia tensor is approximated from the squared distances of the particles to the center of mass.
def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2276    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2277        """
2278        Generates random coordinates outside a sphere and inside a larger bounding sphere.
2279
2280        Args:
2281            center ('array-like'):
2282                Coordinates of the center of the spheres.
2283
2284            radius ('float'):
2285                Radius of the inner exclusion sphere. Must be positive.
2286
2287            max_dist ('float'):
2288                Radius of the outer sampling sphere. Must be larger than 'radius'.
2289
2290            n_samples ('int'):
2291                Number of coordinates to generate.
2292
2293        Returns:
2294            'list' of 'numpy.ndarray':
2295                List of coordinates lying outside the inner sphere and inside the
2296                outer sphere.
2297
2298        Notess:
2299            - Points are uniformly sampled inside a sphere of radius 'max_dist' centered at 'center' 
2300            and only those with a distance greater than or equal to 'radius' from the center are retained.
2301        """
2302        if not radius > 0: 
2303            raise ValueError (f'The value of {radius} must be a positive value')
2304        if not radius < max_dist:
2305            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
2306        coord_list = []
2307        counter = 0
2308        while counter<n_samples:
2309            coord = self.generate_random_points_in_a_sphere(center=center, 
2310                                            radius=max_dist,
2311                                            n_samples=1)[0]
2312            if np.linalg.norm(coord-np.asarray(center))>=radius:
2313                coord_list.append (coord)
2314                counter += 1
2315        return coord_list

Generates random coordinates outside a sphere and inside a larger bounding sphere.

Arguments:
  • center ('array-like'): Coordinates of the center of the spheres.
  • radius ('float'): Radius of the inner exclusion sphere. Must be positive.
  • max_dist ('float'): Radius of the outer sampling sphere. Must be larger than 'radius'.
  • n_samples ('int'): Number of coordinates to generate.
Returns:

'list' of 'numpy.ndarray': List of coordinates lying outside the inner sphere and inside the outer sphere.

Notess:
  • Points are uniformly sampled inside a sphere of radius 'max_dist' centered at 'center' and only those with a distance greater than or equal to 'radius' from the center are retained.
def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2317    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2318        """
2319        Generates uniformly distributed random points inside or on the surface of a sphere.
2320
2321        Args:
2322            center ('array-like'):
2323                Coordinates of the center of the sphere.
2324
2325            radius ('float'):
2326                Radius of the sphere.
2327
2328            n_samples ('int'):
2329                Number of sample points to generate.
2330
2331            on_surface ('bool', optional):
2332                If True, points are uniformly sampled on the surface of the sphere.
2333                If False, points are uniformly sampled within the sphere volume.
2334                Defaults to False.
2335
2336        Returns:
2337            'numpy.ndarray':
2338                Array of shape '(n_samples, d)' containing the generated coordinates,
2339                where 'd' is the dimensionality of 'center'.
2340        Notes:
2341            - Points are sampled in a space whose dimensionality is inferred 
2342            from the length of 'center'.
2343        """
2344        # initial values
2345        center=np.array(center)
2346        d = center.shape[0]
2347        # sample n_samples points in d dimensions from a standard normal distribution
2348        samples = self.rng.normal(size=(n_samples, d))
2349        # make the samples lie on the surface of the unit hypersphere
2350        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
2351        samples /= normalize_radii
2352        if not on_surface:
2353            # make the samples lie inside the hypersphere with the correct density
2354            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
2355            new_radii = np.power(uniform_points, 1/d)
2356            samples *= new_radii
2357        # scale the points to have the correct radius and center
2358        samples = samples * radius + center
2359        return samples 

Generates uniformly distributed random points inside or on the surface of a sphere.

Arguments:
  • center ('array-like'): Coordinates of the center of the sphere.
  • radius ('float'): Radius of the sphere.
  • n_samples ('int'): Number of sample points to generate.
  • on_surface ('bool', optional): If True, points are uniformly sampled on the surface of the sphere. If False, points are uniformly sampled within the sphere volume. Defaults to False.
Returns:

'numpy.ndarray': Array of shape '(n_samples, d)' containing the generated coordinates, where 'd' is the dimensionality of 'center'.

Notes:
  • Points are sampled in a space whose dimensionality is inferred from the length of 'center'.
def generate_trial_perpendicular_vector(self, vector, magnitude):
2361    def generate_trial_perpendicular_vector(self,vector,magnitude):
2362        """
2363        Generates a random vector perpendicular to a given vector.
2364
2365        Args:
2366            vector ('array-like'):
2367                Reference vector to which the generated vector will be perpendicular.
2368
2369            magnitude ('float'):
2370                Desired magnitude of the perpendicular vector.
2371
2372        Returns:
2373            'numpy.ndarray':
2374                Vector orthogonal to 'vector' with norm equal to 'magnitude'.
2375        """ 
2376        np_vec = np.array(vector) 
2377        if np.all(np_vec == 0):
2378            raise ValueError('Zero vector')
2379        np_vec /= np.linalg.norm(np_vec) 
2380        # Generate a random vector 
2381        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
2382                                                                center=[0,0,0],
2383                                                                n_samples=1, 
2384                                                                on_surface=True)[0]
2385        # Project the random vector onto the input vector and subtract the projection
2386        projection = np.dot(random_vector, np_vec) * np_vec
2387        perpendicular_vector = random_vector - projection
2388        # Normalize the perpendicular vector to have the same magnitude as the input vector
2389        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
2390        return perpendicular_vector*magnitude

Generates a random vector perpendicular to a given vector.

Arguments:
  • vector ('array-like'): Reference vector to which the generated vector will be perpendicular.
  • magnitude ('float'): Desired magnitude of the perpendicular vector.
Returns:

'numpy.ndarray': Vector orthogonal to 'vector' with norm equal to 'magnitude'.

def get_bond_template(self, particle_name1, particle_name2, use_default_bond=False):
2392    def get_bond_template(self, particle_name1, particle_name2, use_default_bond=False) :
2393        """
2394        Retrieves a bond template connecting two particle templates.
2395
2396        Args:
2397            particle_name1 ('str'):
2398                Name of the first particle template.
2399
2400            particle_name2 ('str'):
2401                Name of the second particle template.
2402
2403            use_default_bond ('bool', optional):
2404                If True, returns the default bond template when no specific bond
2405                template is found. Defaults to False.
2406
2407        Returns:
2408            'BondTemplate':
2409                Bond template object retrieved from the pyMBE database.
2410            
2411        Notes:
2412            - This method searches the pyMBE database for a bond template defined between particle templates with names 'particle_name1' and 'particle_name2'. 
2413            - If no specific bond template is found and 'use_default_bond' is enabled, a default bond template is returned instead.
2414        """
2415        # Try to find a specific bond template
2416        bond_key = BondTemplate.make_bond_key(pn1=particle_name1,
2417                                              pn2=particle_name2)
2418        try:
2419            return self.db.get_template(name=bond_key, 
2420                                        pmb_type="bond")
2421        except ValueError:
2422            pass
2423
2424        #  Fallback to default bond if allowed
2425        if use_default_bond:
2426            return self.db.get_template(name="default", 
2427                                        pmb_type="bond")
2428
2429        # No bond template found
2430        raise ValueError(f"No bond template found between '{particle_name1}' and '{particle_name2}', and default bonds are deactivated.")

Retrieves a bond template connecting two particle templates.

Arguments:
  • particle_name1 ('str'): Name of the first particle template.
  • particle_name2 ('str'): Name of the second particle template.
  • use_default_bond ('bool', optional): If True, returns the default bond template when no specific bond template is found. Defaults to False.
Returns:

'BondTemplate': Bond template object retrieved from the pyMBE database.

Notes:
  • This method searches the pyMBE database for a bond template defined between particle templates with names 'particle_name1' and 'particle_name2'.
  • If no specific bond template is found and 'use_default_bond' is enabled, a default bond template is returned instead.
def get_charge_number_map(self):
2432    def get_charge_number_map(self):
2433        """
2434        Construct a mapping from ESPResSo particle types to their charge numbers.
2435
2436        Returns:
2437            'dict[int, float]':
2438                Dictionary mapping ESPResSo particle types to charge numbers,
2439                ''{es_type: z}''.
2440
2441        Notess:
2442            - The mapping is built from particle *states*, not instances.
2443            - If multiple templates define states with the same ''es_type'',
2444            the last encountered definition will overwrite previous ones.
2445            This behavior is intentional and assumes database consistency.
2446            - Neutral particles (''z = 0'') are included in the map.
2447        """
2448        charge_number_map = {}
2449        particle_templates = self.db.get_templates("particle")
2450        for tpl in particle_templates.values():
2451            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
2452                charge_number_map[state.es_type] = state.z
2453        return charge_number_map

Construct a mapping from ESPResSo particle types to their charge numbers.

Returns:

'dict[int, float]': Dictionary mapping ESPResSo particle types to charge numbers, ''{es_type: z}''.

Notess:
  • The mapping is built from particle states, not instances.
  • If multiple templates define states with the same ''es_type'', the last encountered definition will overwrite previous ones. This behavior is intentional and assumes database consistency.
  • Neutral particles (''z = 0'') are included in the map.
def get_instances_df(self, pmb_type):
2455    def get_instances_df(self, pmb_type):
2456        """
2457        Returns a dataframe with all instances of type 'pmb_type' in the pyMBE database.
2458
2459        Args:
2460            pmb_type ('str'): 
2461                pmb type to search instances in the pyMBE database.
2462        
2463        Returns:
2464            ('Pandas.Dataframe'): 
2465                Dataframe with all instances of type 'pmb_type'.
2466        """
2467        return self.db._get_instances_df(pmb_type=pmb_type)

Returns a dataframe with all instances of type 'pmb_type' in the pyMBE database.

Arguments:
  • pmb_type ('str'): pmb type to search instances in the pyMBE database.
Returns:

('Pandas.Dataframe'): Dataframe with all instances of type 'pmb_type'.

def get_lj_parameters( self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2469    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2470        """
2471        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2472        'particle_name1' and 'particle_name2' in the pyMBE database, calculated according to the provided combining rule.
2473
2474        Args:
2475            particle_name1 ('str'): 
2476                label of the type of the first particle type
2477
2478            particle_name2 ('str'): 
2479                label of the type of the second particle type
2480
2481            combining_rule ('string', optional): 
2482                combining rule used to calculate 'sigma' and 'epsilon' for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2483
2484        Returns:
2485            ('dict'):
2486                {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2487
2488        Notes:
2489            - Currently, the only 'combining_rule' supported is Lorentz-Berthelot.
2490            - If the sigma value of 'particle_name1' or 'particle_name2' is 0, the function will return an empty dictionary. No LJ interactions are set up for particles with sigma = 0.
2491        """
2492        supported_combining_rules=["Lorentz-Berthelot"]
2493        if combining_rule not in supported_combining_rules:
2494            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2495        part_tpl1 = self.db.get_template(name=particle_name1,
2496                                         pmb_type="particle")
2497        part_tpl2 = self.db.get_template(name=particle_name2,
2498                                         pmb_type="particle")
2499        lj_parameters1 = part_tpl1.get_lj_parameters(ureg=self.units)
2500        lj_parameters2 = part_tpl2.get_lj_parameters(ureg=self.units)
2501
2502        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2503        if part_tpl1.sigma.magnitude == 0 or part_tpl2.sigma.magnitude == 0:
2504            return {}
2505        # Apply combining rule
2506        if combining_rule == 'Lorentz-Berthelot':
2507            sigma=(lj_parameters1["sigma"]+lj_parameters2["sigma"])/2
2508            cutoff=(lj_parameters1["cutoff"]+lj_parameters2["cutoff"])/2
2509            offset=(lj_parameters1["offset"]+lj_parameters2["offset"])/2
2510            epsilon=np.sqrt(lj_parameters1["epsilon"]*lj_parameters2["epsilon"])
2511        return {"sigma": sigma, "cutoff": cutoff, "offset": offset, "epsilon": epsilon}    

Returns the Lennard-Jones parameters for the interaction between the particle types given by 'particle_name1' and 'particle_name2' in the pyMBE database, calculated according to the provided combining rule.

Arguments:
  • particle_name1 ('str'): label of the type of the first particle type
  • particle_name2 ('str'): label of the type of the second particle type
  • combining_rule ('string', optional): combining rule used to calculate 'sigma' and 'epsilon' for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
Returns:

('dict'): {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}

Notes:
  • Currently, the only 'combining_rule' supported is Lorentz-Berthelot.
  • If the sigma value of 'particle_name1' or 'particle_name2' is 0, the function will return an empty dictionary. No LJ interactions are set up for particles with sigma = 0.
def get_particle_id_map(self, object_name):
2513    def get_particle_id_map(self, object_name):
2514        """
2515        Collect all particle IDs associated with an object of given name in the
2516        pyMBE database. 
2517
2518        Args:
2519            object_name ('str'): 
2520                Name of the object.
2521
2522        Returns:
2523            ('dict'): 
2524                {"all": [particle_ids],
2525                 "residue_map": {residue_id: [particle_ids]},
2526                 "molecule_map": {molecule_id: [particle_ids]},
2527                 "assembly_map": {assembly_id: [particle_ids]},}
2528
2529        Notess:
2530            - Works for all supported pyMBE templates.
2531            - Relies in the internal method Manager.get_particle_id_map, see method for the detailed code.
2532        """
2533        return self.db.get_particle_id_map(object_name=object_name)

Collect all particle IDs associated with an object of given name in the pyMBE database.

Arguments:
  • object_name ('str'): Name of the object.
Returns:

('dict'): {"all": [particle_ids], "residue_map": {residue_id: [particle_ids]}, "molecule_map": {molecule_id: [particle_ids]}, "assembly_map": {assembly_id: [particle_ids]},}

Notess:
  • Works for all supported pyMBE templates.
  • Relies in the internal method Manager.get_particle_id_map, see method for the detailed code.
def get_pka_set(self):
2535    def get_pka_set(self):
2536        """
2537        Retrieve the pKa set for all titratable particles in the pyMBE database.
2538
2539        Returns:
2540            ('dict'): 
2541                Dictionary of the form:
2542                {"particle_name": {"pka_value": float,
2543                                   "acidity": "acidic" | "basic"}}
2544        Notes:
2545            - If a particle participates in multiple acid/base reactions, an error is raised.
2546        """
2547        pka_set = {}
2548        supported_reactions = ["monoprotic_acid",
2549                               "monoprotic_base"]
2550        for reaction in self.db._reactions.values():
2551            if reaction.reaction_type not in supported_reactions:
2552                continue
2553            # Identify involved particle(s)
2554            particle_names = {participant.particle_name for participant in reaction.participants}
2555            particle_name = particle_names.pop()
2556            if particle_name in pka_set:
2557                raise ValueError(f"Multiple acid/base reactions found for particle '{particle_name}'.")
2558            pka_set[particle_name] = {"pka_value": reaction.pK}
2559            if reaction.reaction_type == "monoprotic_acid":
2560                acidity = "acidic"
2561            elif reaction.reaction_type == "monoprotic_base":
2562                acidity = "basic"
2563            pka_set[particle_name]["acidity"] = acidity
2564        return pka_set

Retrieve the pKa set for all titratable particles in the pyMBE database.

Returns:

('dict'): Dictionary of the form: {"particle_name": {"pka_value": float, "acidity": "acidic" | "basic"}}

Notes:
  • If a particle participates in multiple acid/base reactions, an error is raised.
def get_radius_map(self, dimensionless=True):
2566    def get_radius_map(self, dimensionless=True):
2567        """
2568        Gets the effective radius of each particle defined in the pyMBE database. 
2569
2570        Args:
2571            dimensionless ('bool'):
2572                If ``True``, return magnitudes expressed in ``reduced_length``.
2573                If ``False``, return Pint quantities with units.
2574        
2575        Returns:
2576            ('dict'): 
2577                {espresso_type: radius}.
2578
2579        Notes:
2580            - The radius corresponds to (sigma+offset)/2
2581        """
2582        if "particle" not in self.db._templates:
2583            return {}          
2584        result = {}
2585        for _, tpl in self.db._templates["particle"].items():
2586            radius = (tpl.sigma.to_quantity(self.units) + tpl.offset.to_quantity(self.units))/2.0
2587            if dimensionless:
2588                magnitude_reduced_length = radius.m_as("reduced_length")
2589                radius = magnitude_reduced_length
2590            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
2591                result[state.es_type] = radius
2592        return result

Gets the effective radius of each particle defined in the pyMBE database.

Arguments:
  • dimensionless ('bool'): If True, return magnitudes expressed in reduced_length. If False, return Pint quantities with units.
Returns:

('dict'): {espresso_type: radius}.

Notes:
  • The radius corresponds to (sigma+offset)/2
def get_reactions_df(self):
2594    def get_reactions_df(self):
2595        """
2596        Returns a dataframe with all reaction templates in the pyMBE database.
2597
2598        Returns:
2599            (Pandas.Dataframe): 
2600                Dataframe with all  reaction templates.
2601        """
2602        return self.db._get_reactions_df()

Returns a dataframe with all reaction templates in the pyMBE database.

Returns:

(Pandas.Dataframe): Dataframe with all reaction templates.

def get_reduced_units(self):
2604    def get_reduced_units(self):
2605        """
2606        Returns the  current set of reduced units defined in pyMBE.
2607
2608        Returns:
2609            reduced_units_text ('str'): 
2610                text with information about the current set of reduced units.
2611
2612        """
2613        unit_length=self.units.Quantity(1,'reduced_length')
2614        unit_energy=self.units.Quantity(1,'reduced_energy')
2615        unit_charge=self.units.Quantity(1,'reduced_charge')
2616        reduced_units_text = "\n".join(["Current set of reduced units:",
2617                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2618                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2619                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2620                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2621                                        ])   
2622        return reduced_units_text

Returns the current set of reduced units defined in pyMBE.

Returns:

reduced_units_text ('str'): text with information about the current set of reduced units.

def get_templates_df(self, pmb_type):
2624    def get_templates_df(self, pmb_type):
2625        """
2626        Returns a dataframe with all templates of type 'pmb_type' in the pyMBE database.
2627
2628        Args:
2629            pmb_type ('str'): 
2630                pmb type to search templates in the pyMBE database.
2631        
2632        Returns:
2633            ('Pandas.Dataframe'): 
2634                Dataframe with all templates of type given by 'pmb_type'.
2635        """
2636        return self.db._get_templates_df(pmb_type=pmb_type)

Returns a dataframe with all templates of type 'pmb_type' in the pyMBE database.

Arguments:
  • pmb_type ('str'): pmb type to search templates in the pyMBE database.
Returns:

('Pandas.Dataframe'): Dataframe with all templates of type given by 'pmb_type'.

def get_type_map(self):
2638    def get_type_map(self):
2639        """
2640        Return the mapping of ESPResSo types for all particle states defined in the pyMBE database.
2641        
2642        Returns:
2643            'dict[str, int]':
2644                A dictionary mapping each particle state to its corresponding ESPResSo type:
2645                {state_name: es_type, ...}
2646        """
2647        
2648        return self.db.get_es_types_map()

Return the mapping of ESPResSo types for all particle states defined in the pyMBE database.

Returns:

'dict[str, int]': A dictionary mapping each particle state to its corresponding ESPResSo type: {state_name: es_type, ...}

def initialize_lattice_builder(self, diamond_lattice):
2650    def initialize_lattice_builder(self, diamond_lattice):
2651        """
2652        Initialize the lattice builder with the DiamondLattice object.
2653
2654        Args:
2655            diamond_lattice ('DiamondLattice'): 
2656                DiamondLattice object from the 'lib/lattice' module to be used in the LatticeBuilder.
2657        """
2658        from .lib.lattice import LatticeBuilder, DiamondLattice
2659        if not isinstance(diamond_lattice, DiamondLattice):
2660            raise TypeError("Currently only DiamondLattice objects are supported.")
2661        self.lattice_builder = LatticeBuilder(lattice=diamond_lattice)
2662        logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}")
2663        return self.lattice_builder

Initialize the lattice builder with the DiamondLattice object.

Arguments:
  • diamond_lattice ('DiamondLattice'): DiamondLattice object from the 'lib/lattice' module to be used in the LatticeBuilder.
def load_database(self, folder, format='csv'):
2665    def load_database(self, folder, format='csv'):
2666        """
2667        Loads a pyMBE database stored in 'folder'.
2668
2669        Args:
2670            folder ('str' or 'Path'): 
2671                Path to the folder where the pyMBE database was stored.
2672
2673            format ('str', optional): 
2674                Format of the database to be loaded. Defaults to 'csv'.
2675
2676        Return:
2677            ('dict'): 
2678                metadata with additional information about the source of the information in the database.
2679
2680        Notes:
2681            - The folder must contain the files generated by 'pmb.save_database()'.
2682            - Currently, only 'csv' format is supported.
2683        """
2684        supported_formats = ['csv']
2685        if format not in supported_formats:
2686            raise ValueError(f"Format {format} not supported. Supported formats are {supported_formats}")
2687        if format == 'csv':
2688            metadata =io._load_database_csv(self.db, 
2689                                            folder=folder)
2690        return metadata

Loads a pyMBE database stored in 'folder'.

Arguments:
  • folder ('str' or 'Path'): Path to the folder where the pyMBE database was stored.
  • format ('str', optional): Format of the database to be loaded. Defaults to 'csv'.
Return:

('dict'): metadata with additional information about the source of the information in the database.

Notes:
  • The folder must contain the files generated by 'pmb.save_database()'.
  • Currently, only 'csv' format is supported.
def load_pka_set(self, filename):
2693    def load_pka_set(self, filename):
2694        """
2695        Load a pKa set and attach chemical states and acid–base reactions
2696        to existing particle templates.
2697
2698        Args:
2699            filename ('str'): 
2700                Path to a JSON file containing the pKa set. Expected format:
2701                {"metadata": {...},
2702                  "data": {"A": {"acidity": "acidic", "pka_value": 4.5},
2703                           "B": {"acidity": "basic",  "pka_value": 9.8}}}
2704
2705        Returns:
2706            ('dict'): 
2707                Dictionary with bibliographic metadata about the original work were the pKa set was determined.
2708
2709        Notes:
2710            - This method is designed for monoprotic acids and bases only.
2711        """
2712        with open(filename, "r") as f:
2713            pka_data = json.load(f)
2714        pka_set = pka_data["data"]
2715        metadata = pka_data.get("metadata", {})
2716        self._check_pka_set(pka_set)
2717        for particle_name, entry in pka_set.items():
2718            acidity = entry["acidity"]
2719            pka = entry["pka_value"]
2720            self.define_monoprototic_acidbase_reaction(particle_name=particle_name,
2721                                                       pka=pka,
2722                                                       acidity=acidity,
2723                                                       metadata=metadata)
2724        return metadata

Load a pKa set and attach chemical states and acid–base reactions to existing particle templates.

Arguments:
  • filename ('str'): Path to a JSON file containing the pKa set. Expected format: {"metadata": {...}, "data": {"A": {"acidity": "acidic", "pka_value": 4.5}, "B": {"acidity": "basic", "pka_value": 9.8}}}
Returns:

('dict'): Dictionary with bibliographic metadata about the original work were the pKa set was determined.

Notes:
  • This method is designed for monoprotic acids and bases only.
def propose_unused_type(self):
2726    def propose_unused_type(self):
2727        """
2728        Propose an unused ESPResSo particle type.
2729
2730        Returns:
2731            ('int'): 
2732                The next available integer ESPResSo type. Returns ''0'' if no integer types are currently defined.
2733        """
2734        type_map = self.get_type_map()
2735        # Flatten all es_type values across all particles and states
2736        all_types = []
2737        for es_type in type_map.values():
2738            all_types.append(es_type)
2739        # If no es_types exist, start at 0
2740        if not all_types:
2741            return 0
2742        return max(all_types) + 1

Propose an unused ESPResSo particle type.

Returns:

('int'): The next available integer ESPResSo type. Returns ''0'' if no integer types are currently defined.

def read_protein_vtf(self, filename, unit_length=None):
2744    def read_protein_vtf(self, filename, unit_length=None):
2745        """
2746        Loads a coarse-grained protein model from a VTF file.
2747
2748        Args:
2749            filename ('str'): 
2750                Path to the VTF file.
2751                
2752            unit_length ('Pint.Quantity'): 
2753                Unit of length for coordinates (pyMBE UnitRegistry). Defaults to Angstrom.
2754
2755        Returns:
2756            ('tuple'):
2757                ('dict'): Particle topology.    
2758                ('str'):  One-letter amino-acid sequence (including n/c ends).
2759        """
2760        logging.info(f"Loading protein coarse-grain model file: {filename}")
2761        if unit_length is None:
2762            unit_length = 1 * self.units.angstrom
2763        atoms = {}        # atom_id -> atom info
2764        coords = []       # ordered coordinates
2765        residues = {}     # resid -> resname (first occurrence)
2766        has_n_term = False
2767        has_c_term = False
2768        aa_3to1 = {"ALA": "A", "ARG": "R", "ASN": "N", "ASP": "D",
2769                   "CYS": "C", "GLU": "E", "GLN": "Q", "GLY": "G",
2770                   "HIS": "H", "ILE": "I", "LEU": "L", "LYS": "K",
2771                   "MET": "M", "PHE": "F", "PRO": "P", "SER": "S",
2772                   "THR": "T", "TRP": "W", "TYR": "Y", "VAL": "V",
2773                   "n": "n", "c": "c"}
2774        # --- parse VTF ---
2775        with open(filename, "r") as f:
2776            for line in f:
2777                fields = line.split()
2778                if not fields:
2779                    continue
2780                if fields[0] == "atom":
2781                    atom_id = int(fields[1])
2782                    atom_name = fields[3]
2783                    resname = fields[5]
2784                    resid = int(fields[7])
2785                    chain_id = fields[9]
2786                    radius = float(fields[11]) * unit_length
2787                    atoms[atom_id] = {"name": atom_name,
2788                                     "resname": resname,
2789                                     "resid": resid,
2790                                     "chain_id": chain_id,
2791                                     "radius": radius}
2792                    if resname == "n":
2793                        has_n_term = True
2794                    elif resname == "c":
2795                        has_c_term = True
2796                    # register residue 
2797                    if resid not in residues:
2798                        residues[resid] = resname
2799                elif fields[0].isnumeric():
2800                    xyz = [(float(x) * unit_length).to("reduced_length").magnitude
2801                        for x in fields[1:4]]
2802                    coords.append(xyz)
2803        sequence = ""
2804        # N-terminus
2805        if has_n_term:
2806            sequence += "n"
2807        # protein residues only
2808        protein_resids = sorted(resid for resid, resname in residues.items()  if resname not in ("n", "c", "Ca"))
2809        for resid in protein_resids:
2810            resname = residues[resid]
2811            try:
2812                sequence += aa_3to1[resname]
2813            except KeyError:
2814                raise ValueError(f"Unknown residue name '{resname}' in VTF file")
2815        # C-terminus
2816        if has_c_term:
2817            sequence += "c"
2818        last_resid = max(protein_resids)
2819        # --- build topology ---
2820        topology_dict = {}
2821        for atom_id in sorted(atoms.keys()):
2822            atom = atoms[atom_id]
2823            resname = atom["resname"]
2824            resid = atom["resid"]
2825            # apply labeling rules
2826            if resname == "n":
2827                label_resid = 0
2828            elif resname == "c":
2829                label_resid = last_resid + 1
2830            elif resname == "Ca":
2831                label_resid = last_resid + 2
2832            else:
2833                label_resid = resid  # preserve original resid 
2834            label = f"{atom['name']}{label_resid}"
2835            if label in topology_dict:
2836                raise ValueError(f"Duplicate particle label '{label}'. Check VTF residue definitions.")
2837            topology_dict[label] = {"initial_pos": coords[atom_id - 1], "chain_id": atom["chain_id"], "radius": atom["radius"],}
2838        return topology_dict, sequence

Loads a coarse-grained protein model from a VTF file.

Arguments:
  • filename ('str'): Path to the VTF file.
  • unit_length ('Pint.Quantity'): Unit of length for coordinates (pyMBE UnitRegistry). Defaults to Angstrom.
Returns:

('tuple'): ('dict'): Particle topology.
('str'): One-letter amino-acid sequence (including n/c ends).

def save_database(self, folder, format='csv'):
2841    def save_database(self, folder, format='csv'):
2842        """
2843        Saves the current pyMBE database into a file 'filename'.
2844
2845        Args:
2846            folder ('str' or 'Path'): 
2847                Path to the folder where the database files will be saved.
2848
2849        """
2850        supported_formats = ['csv']
2851        if format not in supported_formats:
2852            raise ValueError(f"Format {format} not supported. Supported formats are: {supported_formats}")
2853        if format == 'csv':
2854            io._save_database_csv(self.db, 
2855                                folder=folder)

Saves the current pyMBE database into a file 'filename'.

Arguments:
  • folder ('str' or 'Path'): Path to the folder where the database files will be saved.
def set_particle_initial_state(self, particle_name, state_name):
2857    def set_particle_initial_state(self, particle_name, state_name):
2858        """
2859        Sets the default initial state of a particle template defined in the pyMBE database.
2860
2861        Args:
2862            particle_name ('str'): 
2863                Unique label that identifies the particle template. 
2864
2865            state_name ('str'): 
2866                Name of the state to be set as default initial state.
2867        """
2868        part_tpl = self.db.get_template(name=particle_name,
2869        
2870                                        pmb_type="particle")
2871        part_tpl.initial_state = state_name
2872        logging.info(f"Default initial state of particle {particle_name} set to {state_name}.")

Sets the default initial state of a particle template defined in the pyMBE database.

Arguments:
  • particle_name ('str'): Unique label that identifies the particle template.
  • state_name ('str'): Name of the state to be set as default initial state.
def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None):
2874    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None):
2875        """
2876        Sets the set of reduced units used by pyMBE.units and it prints it.
2877
2878        Args:
2879            unit_length ('pint.Quantity', optional): 
2880                Reduced unit of length defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2881
2882            unit_charge ('pint.Quantity', optional): 
2883                Reduced unit of charge defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2884
2885            temperature ('pint.Quantity', optional): 
2886                Temperature of the system, defined using the 'pmb.units' UnitRegistry. Defaults to None. 
2887
2888            Kw ('pint.Quantity', optional): 
2889                Ionic product of water in mol^2/l^2. Defaults to None. 
2890
2891        Notes:
2892            - If no 'temperature' is given, a value of 298.15 K is assumed by default.
2893            - If no 'unit_length' is given, a value of 0.355 nm is assumed by default.
2894            - If no 'unit_charge' is given, a value of 1 elementary charge is assumed by default. 
2895            - If no 'Kw' is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2896        """
2897        if unit_length is None:
2898            unit_length= 0.355*self.units.nm
2899        if temperature is None:
2900            temperature = 298.15 * self.units.K
2901        if unit_charge is None:
2902            unit_charge = scipy.constants.e * self.units.C
2903        if Kw is None:
2904            Kw = 1e-14
2905        # Sanity check
2906        variables=[unit_length,temperature,unit_charge]
2907        dimensionalities=["[length]","[temperature]","[charge]"]
2908        for variable,dimensionality in zip(variables,dimensionalities):
2909            self._check_dimensionality(variable,dimensionality)
2910        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2911        self.kT=temperature*self.kB
2912        self.units._build_cache()
2913        self.units.define(f'reduced_energy = {self.kT} ')
2914        self.units.define(f'reduced_length = {unit_length}')
2915        self.units.define(f'reduced_charge = {unit_charge}')
2916        logging.info(self.get_reduced_units())

Sets the set of reduced units used by pyMBE.units and it prints it.

Arguments:
  • unit_length ('pint.Quantity', optional): Reduced unit of length defined using the 'pmb.units' UnitRegistry. Defaults to None.
  • unit_charge ('pint.Quantity', optional): Reduced unit of charge defined using the 'pmb.units' UnitRegistry. Defaults to None.
  • temperature ('pint.Quantity', optional): Temperature of the system, defined using the 'pmb.units' UnitRegistry. Defaults to None.
  • Kw ('pint.Quantity', optional): Ionic product of water in mol^2/l^2. Defaults to None.
Notes:
  • If no 'temperature' is given, a value of 298.15 K is assumed by default.
  • If no 'unit_length' is given, a value of 0.355 nm is assumed by default.
  • If no 'unit_charge' is given, a value of 1 elementary charge is assumed by default.
  • If no 'Kw' is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default.
def setup_cpH( self, counter_ion, constant_pH, exclusion_range=None, use_exclusion_radius_per_type=False):
2918    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, use_exclusion_radius_per_type = False):
2919        """
2920        Sets up the Acid/Base reactions for acidic/basic particles defined in the pyMBE database
2921        to be sampled in the constant pH ensemble. 
2922
2923        Args:
2924            counter_ion ('str'): 
2925                'name' of the counter_ion 'particle'.
2926
2927            constant_pH ('float'): 
2928                pH-value.
2929
2930            exclusion_range ('pint.Quantity', optional): 
2931                Below this value, no particles will be inserted.
2932
2933            use_exclusion_radius_per_type ('bool', optional): 
2934                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
2935
2936        Returns:
2937            ('reaction_methods.ConstantpHEnsemble'): 
2938                Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2939        """
2940        from espressomd import reaction_methods
2941        if exclusion_range is None:
2942            exclusion_range = max(self.get_radius_map().values())*2.0
2943        if use_exclusion_radius_per_type:
2944            exclusion_radius_per_type = self.get_radius_map()
2945        else:
2946            exclusion_radius_per_type = {}
2947        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2948                                                exclusion_range=exclusion_range, 
2949                                                seed=self.seed, 
2950                                                constant_pH=constant_pH,
2951                                                exclusion_radius_per_type = exclusion_radius_per_type)
2952        conterion_tpl = self.db.get_template(name=counter_ion,
2953                                             pmb_type="particle")
2954        conterion_state = self.db.get_template(name=conterion_tpl.initial_state,
2955                                               pmb_type="particle_state")
2956        for reaction in self.db.get_reactions():
2957            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
2958                continue
2959            default_charges = {}
2960            reactant_types  = []
2961            product_types   = []
2962            for participant in reaction.participants:
2963                state_tpl = self.db.get_template(name=participant.state_name,
2964                                                 pmb_type="particle_state")
2965                default_charges[state_tpl.es_type] = state_tpl.z
2966                if participant.coefficient < 0:
2967                    reactant_types.append(state_tpl.es_type)
2968                elif participant.coefficient > 0:
2969                    product_types.append(state_tpl.es_type)
2970            # Add counterion to the products
2971            if conterion_state.es_type not in product_types:
2972                product_types.append(conterion_state.es_type)
2973                default_charges[conterion_state.es_type] = conterion_state.z
2974                reaction.add_participant(particle_name=counter_ion,
2975                                         state_name=conterion_tpl.initial_state,
2976                                         coefficient=1)
2977            gamma=10**-reaction.pK
2978            RE.add_reaction(gamma=gamma,
2979                            reactant_types=reactant_types,
2980                            product_types=product_types,
2981                            default_charges=default_charges)
2982            reaction.add_simulation_method(simulation_method="cpH")
2983        return RE

Sets up the Acid/Base reactions for acidic/basic particles defined in the pyMBE database to be sampled in the constant pH ensemble.

Arguments:
  • counter_ion ('str'): 'name' of the counter_ion 'particle'.
  • constant_pH ('float'): pH-value.
  • exclusion_range ('pint.Quantity', optional): Below this value, no particles will be inserted.
  • use_exclusion_radius_per_type ('bool', optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
Returns:

('reaction_methods.ConstantpHEnsemble'): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.

def setup_gcmc( self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type=False):
2985    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2986        """
2987        Sets up grand-canonical coupling to a reservoir of salt.
2988        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2989
2990        Args:
2991            c_salt_res ('pint.Quantity'): 
2992                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2993
2994            salt_cation_name ('str'): 
2995                Name of the salt cation (e.g. Na+) particle.
2996
2997            salt_anion_name ('str'): 
2998                Name of the salt anion (e.g. Cl-) particle.
2999
3000            activity_coefficient ('callable'): 
3001                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3002
3003            exclusion_range('pint.Quantity', optional): 
3004                For distances shorter than this value, no particles will be inserted.
3005
3006            use_exclusion_radius_per_type('bool',optional): 
3007                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
3008
3009        Returns:
3010            ('reaction_methods.ReactionEnsemble'): 
3011                Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
3012        """
3013        from espressomd import reaction_methods
3014        if exclusion_range is None:
3015            exclusion_range = max(self.get_radius_map().values())*2.0
3016        if use_exclusion_radius_per_type:
3017            exclusion_radius_per_type = self.get_radius_map()
3018        else:
3019            exclusion_radius_per_type = {}
3020        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3021                                               exclusion_range=exclusion_range, 
3022                                               seed=self.seed, 
3023                                               exclusion_radius_per_type = exclusion_radius_per_type)
3024        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3025        determined_activity_coefficient = activity_coefficient(c_salt_res)
3026        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
3027        cation_tpl = self.db.get_template(pmb_type="particle",
3028                                          name=salt_cation_name)
3029        cation_state = self.db.get_template(pmb_type="particle_state",
3030                                            name=cation_tpl.initial_state)
3031        anion_tpl = self.db.get_template(pmb_type="particle",
3032                                          name=salt_anion_name)
3033        anion_state = self.db.get_template(pmb_type="particle_state",
3034                                            name=anion_tpl.initial_state)
3035        salt_cation_es_type = cation_state.es_type
3036        salt_anion_es_type = anion_state.es_type     
3037        salt_cation_charge = cation_state.z
3038        salt_anion_charge = anion_state.z
3039        if salt_cation_charge <= 0:
3040            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3041        if salt_anion_charge >= 0:
3042            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3043        # Grand-canonical coupling to the reservoir
3044        RE.add_reaction(gamma = K_salt.magnitude,
3045                        reactant_types = [],
3046                        reactant_coefficients = [],
3047                        product_types = [ salt_cation_es_type, salt_anion_es_type ],
3048                        product_coefficients = [ 1, 1 ],
3049                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3050                                           salt_anion_es_type: salt_anion_charge})
3051        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3052                                                            state_name=cation_state.name,
3053                                                            coefficient=1),
3054                                        ReactionParticipant(particle_name=salt_anion_name,
3055                                                            state_name=anion_state.name,
3056                                                            coefficient=1)],
3057                           pK=-np.log10(K_salt.magnitude),
3058                           reaction_type="ion_insertion",
3059                           simulation_method="GCMC")
3060        self.db._register_reaction(rx_tpl)
3061        return RE

Sets up grand-canonical coupling to a reservoir of salt. For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.

Arguments:
  • c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
  • salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
  • salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
  • activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
  • exclusion_range('pint.Quantity', optional): For distances shorter than this value, no particles will be inserted.
  • use_exclusion_radius_per_type('bool',optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
Returns:

('reaction_methods.ReactionEnsemble'): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.

def setup_grxmc_reactions( self, pH_res, c_salt_res, proton_name, hydroxide_name, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type=False):
3063    def setup_grxmc_reactions(self, pH_res, c_salt_res, proton_name, hydroxide_name, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
3064        """
3065        Sets up acid/base reactions for acidic/basic monoprotic particles defined in the pyMBE database, 
3066        as well as a grand-canonical coupling to a reservoir of small ions. 
3067        
3068        Args:
3069            pH_res ('float'): 
3070                pH-value in the reservoir.
3071
3072            c_salt_res ('pint.Quantity'): 
3073                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3074
3075            proton_name ('str'): 
3076                Name of the proton (H+) particle.
3077
3078            hydroxide_name ('str'): 
3079                Name of the hydroxide (OH-) particle.
3080
3081            salt_cation_name ('str'): 
3082                Name of the salt cation (e.g. Na+) particle.
3083
3084            salt_anion_name ('str'): 
3085                Name of the salt anion (e.g. Cl-) particle.
3086
3087            activity_coefficient ('callable'): 
3088                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3089
3090            exclusion_range('pint.Quantity', optional): 
3091                For distances shorter than this value, no particles will be inserted.
3092
3093            use_exclusion_radius_per_type('bool', optional): 
3094                Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
3095
3096        Returns:
3097            'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':
3098
3099                'reaction_methods.ReactionEnsemble':  
3100                    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.
3101                
3102                'pint.Quantity': 
3103                    Ionic strength of the reservoir (useful for calculating partition coefficients).
3104
3105        Notess:
3106            - This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
3107
3108        [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3109        """
3110        from espressomd import reaction_methods
3111        if exclusion_range is None:
3112            exclusion_range = max(self.get_radius_map().values())*2.0
3113        if use_exclusion_radius_per_type:
3114            exclusion_radius_per_type = self.get_radius_map()
3115        else:
3116            exclusion_radius_per_type = {}
3117        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3118                                               exclusion_range=exclusion_range, 
3119                                               seed=self.seed, 
3120                                               exclusion_radius_per_type = exclusion_radius_per_type)
3121        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3122        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3123        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3124        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3125        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3126        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3127        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3128        cation_tpl = self.db.get_template(pmb_type="particle",
3129                                          name=salt_cation_name)
3130        cation_state = self.db.get_template(pmb_type="particle_state",
3131                                            name=cation_tpl.initial_state)
3132        anion_tpl = self.db.get_template(pmb_type="particle",
3133                                          name=salt_anion_name)
3134        anion_state = self.db.get_template(pmb_type="particle_state",
3135                                            name=anion_tpl.initial_state)
3136        proton_tpl = self.db.get_template(pmb_type="particle",
3137                                          name=proton_name)
3138        proton_state = self.db.get_template(pmb_type="particle_state",
3139                                            name=proton_tpl.initial_state)
3140        hydroxide_tpl = self.db.get_template(pmb_type="particle",
3141                                             name=hydroxide_name)
3142        hydroxide_state = self.db.get_template(pmb_type="particle_state",
3143                                               name=hydroxide_tpl.initial_state)
3144        proton_es_type = proton_state.es_type
3145        hydroxide_es_type = hydroxide_state.es_type
3146        salt_cation_es_type = cation_state.es_type
3147        salt_anion_es_type = anion_state.es_type
3148        proton_charge = proton_state.z
3149        hydroxide_charge = hydroxide_state.z          
3150        salt_cation_charge = cation_state.z
3151        salt_anion_charge = anion_state.z      
3152        if proton_charge <= 0:
3153            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
3154        if salt_cation_charge <= 0:
3155            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3156        if hydroxide_charge >= 0:
3157            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
3158        if salt_anion_charge >= 0:
3159            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3160        # Grand-canonical coupling to the reservoir
3161        # 0 = H+ + OH-
3162        RE.add_reaction(gamma = K_W.magnitude,
3163                        reactant_types = [],
3164                        reactant_coefficients = [],
3165                        product_types = [ proton_es_type, hydroxide_es_type ],
3166                        product_coefficients = [ 1, 1 ],
3167                        default_charges = {proton_es_type: proton_charge, 
3168                                           hydroxide_es_type: hydroxide_charge})
3169        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3170                                                            state_name=proton_state.name,
3171                                                            coefficient=1),
3172                                        ReactionParticipant(particle_name=hydroxide_name,
3173                                                            state_name=hydroxide_state.name,
3174                                                            coefficient=1)],
3175                           pK=-np.log10(K_W.magnitude),
3176                           reaction_type="ion_insertion",
3177                           simulation_method="GRxMC")
3178        self.db._register_reaction(rx_tpl)
3179        # 0 = Na+ + Cl-
3180        RE.add_reaction(gamma = K_NACL.magnitude,
3181                        reactant_types = [],
3182                        reactant_coefficients = [],
3183                        product_types = [ salt_cation_es_type, salt_anion_es_type ],
3184                        product_coefficients = [ 1, 1 ],
3185                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3186                                        salt_anion_es_type: salt_anion_charge})
3187        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3188                                                            state_name=cation_state.name,
3189                                                            coefficient=1),
3190                                        ReactionParticipant(particle_name=salt_anion_name,
3191                                                            state_name=anion_state.name,
3192                                                            coefficient=1)],
3193                           pK=-np.log10(K_NACL.magnitude),
3194                           reaction_type="ion_insertion",
3195                           simulation_method="GRxMC")
3196        self.db._register_reaction(rx_tpl)
3197        # 0 = Na+ + OH-
3198        RE.add_reaction(gamma = (K_NACL * K_W / K_HCL).magnitude,
3199                        reactant_types = [],
3200                        reactant_coefficients = [],
3201                        product_types = [ salt_cation_es_type, hydroxide_es_type ],
3202                        product_coefficients = [ 1, 1 ],
3203                        default_charges = {salt_cation_es_type: salt_cation_charge, 
3204                                           hydroxide_es_type: hydroxide_charge})
3205        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=salt_cation_name,
3206                                                            state_name=cation_state.name,
3207                                                            coefficient=1),
3208                                        ReactionParticipant(particle_name=hydroxide_name,
3209                                                            state_name=hydroxide_state.name,
3210                                                            coefficient=1)],
3211                           pK=-np.log10((K_NACL * K_W / K_HCL).magnitude),
3212                           reaction_type="ion_insertion",
3213                           simulation_method="GRxMC")
3214        self.db._register_reaction(rx_tpl)
3215        # 0 = H+ + Cl-
3216        RE.add_reaction(gamma = K_HCL.magnitude,
3217                        reactant_types = [],
3218                        reactant_coefficients = [],
3219                        product_types = [ proton_es_type, salt_anion_es_type ],
3220                        product_coefficients = [ 1, 1 ],
3221                        default_charges = {proton_es_type: proton_charge, 
3222                                           salt_anion_es_type: salt_anion_charge})
3223        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3224                                                            state_name=proton_state.name,
3225                                                            coefficient=1),
3226                                        ReactionParticipant(particle_name=salt_anion_name,
3227                                                            state_name=anion_state.name,
3228                                                            coefficient=1)],
3229                           pK=-np.log10(K_HCL.magnitude),
3230                           reaction_type="ion_insertion",
3231                           simulation_method="GRxMC")
3232        self.db._register_reaction(rx_tpl)
3233        # Annealing moves to ensure sufficient sampling
3234        # Cation annealing H+ = Na+
3235        RE.add_reaction(gamma = (K_NACL / K_HCL).magnitude,
3236                        reactant_types = [proton_es_type],
3237                        reactant_coefficients = [ 1 ],
3238                        product_types = [ salt_cation_es_type ],
3239                        product_coefficients = [ 1 ],
3240                        default_charges = {proton_es_type: proton_charge, 
3241                                           salt_cation_es_type: salt_cation_charge})
3242        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=proton_name,
3243                                                            state_name=proton_state.name,
3244                                                            coefficient=-1),
3245                                        ReactionParticipant(particle_name=salt_cation_name,
3246                                                            state_name=cation_state.name,
3247                                                            coefficient=1)],
3248                           pK=-np.log10((K_NACL / K_HCL).magnitude),
3249                           reaction_type="particle replacement",
3250                           simulation_method="GRxMC")
3251        self.db._register_reaction(rx_tpl)
3252        # Anion annealing OH- = Cl- 
3253        RE.add_reaction(gamma = (K_HCL / K_W).magnitude,
3254                        reactant_types = [hydroxide_es_type],
3255                        reactant_coefficients = [ 1 ],
3256                        product_types = [ salt_anion_es_type ],
3257                        product_coefficients = [ 1 ],
3258            default_charges = {hydroxide_es_type: hydroxide_charge, 
3259                               salt_anion_es_type: salt_anion_charge})
3260        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=hydroxide_name,
3261                                                            state_name=hydroxide_state.name,
3262                                                            coefficient=-1),
3263                                        ReactionParticipant(particle_name=salt_anion_name,
3264                                                            state_name=anion_state.name,
3265                                                            coefficient=1)],
3266                           pK=-np.log10((K_HCL / K_W).magnitude),
3267                           reaction_type="particle replacement",
3268                           simulation_method="GRxMC")
3269        self.db._register_reaction(rx_tpl)
3270        for reaction in self.db.get_reactions():
3271            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
3272                continue
3273            default_charges = {}
3274            reactant_types  = []
3275            product_types   = []
3276            for participant in reaction.participants:
3277                state_tpl = self.db.get_template(name=participant.state_name,
3278                                                 pmb_type="particle_state")
3279                default_charges[state_tpl.es_type] = state_tpl.z
3280                if participant.coefficient < 0:
3281                    reactant_types.append(state_tpl.es_type)
3282                    reactant_name=state_tpl.particle_name
3283                    reactant_state_name=state_tpl.name
3284                elif participant.coefficient > 0:
3285                    product_types.append(state_tpl.es_type)
3286                    product_name=state_tpl.particle_name
3287                    product_state_name=state_tpl.name
3288
3289            Ka = (10**-reaction.pK * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3290            # Reaction in terms of proton: HA = A + H+
3291            RE.add_reaction(gamma=Ka.magnitude,
3292                            reactant_types=reactant_types,
3293                            reactant_coefficients=[1],
3294                            product_types=product_types+[proton_es_type],
3295                            product_coefficients=[1, 1],
3296                            default_charges= default_charges | {proton_es_type: proton_charge})
3297            reaction.add_participant(particle_name=proton_name,
3298                                     state_name=proton_state.name,
3299                                     coefficient=1)
3300            reaction.add_simulation_method("GRxMC")
3301            # Reaction in terms of salt cation: HA = A + Na+
3302            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
3303                            reactant_types=reactant_types,
3304                            reactant_coefficients=[1],
3305                            product_types=product_types+[salt_cation_es_type],
3306                            product_coefficients=[1, 1],
3307                            default_charges=default_charges | {salt_cation_es_type: salt_cation_charge})
3308            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3309                                                                state_name=reactant_state_name,
3310                                                                coefficient=-1),
3311                                            ReactionParticipant(particle_name=product_name,
3312                                                                state_name=product_state_name,
3313                                                                coefficient=1),
3314                                            ReactionParticipant(particle_name=salt_cation_name,
3315                                                                state_name=cation_state.name,
3316                                                                coefficient=1),],
3317                              pK=-np.log10((Ka * K_NACL / K_HCL).magnitude),
3318                              reaction_type=reaction.reaction_type+"_salt",
3319                              simulation_method="GRxMC")
3320            self.db._register_reaction(rx_tpl)
3321            # Reaction in terms of hydroxide: OH- + HA = A
3322            RE.add_reaction(gamma=(Ka / K_W).magnitude,
3323                            reactant_types=reactant_types+[hydroxide_es_type],
3324                            reactant_coefficients=[1, 1],
3325                            product_types=product_types,
3326                            product_coefficients=[1],
3327                            default_charges=default_charges | {hydroxide_es_type: hydroxide_charge})
3328            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3329                                                                state_name=reactant_state_name,
3330                                                                coefficient=-1),
3331                                            ReactionParticipant(particle_name=product_name,
3332                                                                state_name=product_state_name,
3333                                                                coefficient=1),
3334                                            ReactionParticipant(particle_name=hydroxide_name,
3335                                                                state_name=hydroxide_state.name,
3336                                                                coefficient=-1),],
3337                              pK=-np.log10((Ka / K_W).magnitude),
3338                              reaction_type=reaction.reaction_type+"_conjugate",
3339                              simulation_method="GRxMC")
3340            self.db._register_reaction(rx_tpl)
3341            # Reaction in terms of salt anion: Cl- + HA = A
3342            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
3343                            reactant_types=reactant_types+[salt_anion_es_type],
3344                            reactant_coefficients=[1, 1],
3345                            product_types=product_types,
3346                            product_coefficients=[1],
3347                            default_charges=default_charges | {salt_anion_es_type: salt_anion_charge})
3348            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3349                                                                state_name=reactant_state_name,
3350                                                                coefficient=-1),
3351                                            ReactionParticipant(particle_name=product_name,
3352                                                                state_name=product_state_name,
3353                                                                coefficient=1),
3354                                            ReactionParticipant(particle_name=salt_anion_name,
3355                                                                state_name=anion_state.name,
3356                                                                coefficient=-1),],
3357                              pK=-np.log10((Ka / K_HCL).magnitude),
3358                              reaction_type=reaction.reaction_type+"_salt",
3359                              simulation_method="GRxMC")
3360            self.db._register_reaction(rx_tpl)
3361        return RE, ionic_strength_res

Sets up acid/base reactions for acidic/basic monoprotic particles defined in the pyMBE database, as well as a grand-canonical coupling to a reservoir of small ions.

Arguments:
  • pH_res ('float'): pH-value in the reservoir.
  • c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
  • proton_name ('str'): Name of the proton (H+) particle.
  • hydroxide_name ('str'): Name of the hydroxide (OH-) particle.
  • salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
  • salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
  • activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
  • exclusion_range('pint.Quantity', optional): For distances shorter than this value, no particles will be inserted.
  • use_exclusion_radius_per_type('bool', optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to 'False'.
Returns:

'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':

'reaction_methods.ReactionEnsemble':  
    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.

'pint.Quantity': 
    Ionic strength of the reservoir (useful for calculating partition coefficients).
Notess:
  • This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].

[1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.

def setup_grxmc_unified( self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type=False):
3363    def setup_grxmc_unified(self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
3364        """
3365        Sets up acid/base reactions for acidic/basic 'particles' defined in the pyMBE database, as well as a grand-canonical coupling to a 
3366        reservoir of small ions using a unified formulation for small ions.
3367
3368        Args:
3369            pH_res ('float'): 
3370                pH-value in the reservoir.
3371
3372            c_salt_res ('pint.Quantity'): 
3373                Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3374
3375            cation_name ('str'): 
3376                Name of the cationic particle.
3377
3378            anion_name ('str'): 
3379                Name of the anionic particle.
3380
3381            activity_coefficient ('callable'): 
3382                A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3383
3384            exclusion_range('pint.Quantity', optional): 
3385                Below this value, no particles will be inserted.
3386            
3387            use_exclusion_radius_per_type('bool', optional): 
3388                Controls if one exclusion_radius per each espresso_type. Defaults to 'False'.
3389
3390        Returns:
3391            'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':
3392
3393                'reaction_methods.ReactionEnsemble':  
3394                    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.
3395                
3396                'pint.Quantity': 
3397                    Ionic strength of the reservoir (useful for calculating partition coefficients).
3398
3399        Notes:
3400            - This implementation uses the formulation of the grand-reaction method by Curk et al. [1], which relies on "unified" ion types X+ = {H+, Na+} and X- = {OH-, Cl-}. 
3401            - A function that implements the original version of the grand-reaction method by Landsgesell et al. [2] is also available under the name 'setup_grxmc_reactions'.
3402
3403        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3404        [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3405        """
3406        from espressomd import reaction_methods
3407        if exclusion_range is None:
3408            exclusion_range = max(self.get_radius_map().values())*2.0
3409        if use_exclusion_radius_per_type:
3410            exclusion_radius_per_type = self.get_radius_map()
3411        else:
3412            exclusion_radius_per_type = {}
3413        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3414                                               exclusion_range=exclusion_range, 
3415                                               seed=self.seed, 
3416                                               exclusion_radius_per_type = exclusion_radius_per_type)
3417        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3418        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3419        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3420        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3421        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3422        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3423        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3424        K_XX = a_cation * a_anion
3425        cation_tpl = self.db.get_template(pmb_type="particle",
3426                                          name=cation_name)
3427        cation_state = self.db.get_template(pmb_type="particle_state",
3428                                            name=cation_tpl.initial_state)
3429        anion_tpl = self.db.get_template(pmb_type="particle",
3430                                          name=anion_name)
3431        anion_state = self.db.get_template(pmb_type="particle_state",
3432                                            name=anion_tpl.initial_state)
3433        cation_es_type = cation_state.es_type
3434        anion_es_type = anion_state.es_type     
3435        cation_charge = cation_state.z
3436        anion_charge = anion_state.z
3437        if cation_charge <= 0:
3438            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3439        if anion_charge >= 0:
3440            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3441        # Coupling to the reservoir: 0 = X+ + X-
3442        RE.add_reaction(gamma = K_XX.magnitude,
3443                        reactant_types = [],
3444                        reactant_coefficients = [],
3445                        product_types = [ cation_es_type, anion_es_type ],
3446                        product_coefficients = [ 1, 1 ],
3447                        default_charges = {cation_es_type: cation_charge, 
3448                                           anion_es_type: anion_charge})
3449        rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=cation_name,
3450                                                            state_name=cation_state.name,
3451                                                            coefficient=1),
3452                                        ReactionParticipant(particle_name=anion_name,
3453                                                            state_name=anion_state.name,
3454                                                            coefficient=1)],
3455                           pK=-np.log10(K_XX.magnitude),
3456                           reaction_type="ion_insertion",
3457                           simulation_method="GCMC")
3458        self.db._register_reaction(rx_tpl)
3459        for reaction in self.db.get_reactions():
3460            if reaction.reaction_type not in ["monoprotic_acid", "monoprotic_base"]:
3461                continue
3462            default_charges = {}
3463            reactant_types  = []
3464            product_types   = []
3465            for participant in reaction.participants:
3466                state_tpl = self.db.get_template(name=participant.state_name,
3467                                                 pmb_type="particle_state")
3468                default_charges[state_tpl.es_type] = state_tpl.z
3469                if participant.coefficient < 0:
3470                    reactant_types.append(state_tpl.es_type)
3471                    reactant_name=state_tpl.particle_name
3472                    reactant_state_name=state_tpl.name
3473                elif participant.coefficient > 0:
3474                    product_types.append(state_tpl.es_type)
3475                    product_name=state_tpl.particle_name
3476                    product_state_name=state_tpl.name
3477
3478            Ka = (10**-reaction.pK * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3479            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3480            # Reaction in terms of small cation: HA = A + X+
3481            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3482                            reactant_types=reactant_types,
3483                            reactant_coefficients=[1],
3484                            product_types=product_types+[cation_es_type],
3485                            product_coefficients=[1, 1],
3486                            default_charges=default_charges|{cation_es_type: cation_charge})
3487            reaction.add_participant(particle_name=cation_name,
3488                                     state_name=cation_state.name,
3489                                     coefficient=1)
3490            reaction.add_simulation_method("GRxMC")
3491            # Reaction in terms of small anion: X- + HA = A
3492            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3493                            reactant_types=reactant_types+[anion_es_type],
3494                            reactant_coefficients=[1, 1],
3495                            product_types=product_types,
3496                            product_coefficients=[1],
3497                            default_charges=default_charges|{anion_es_type: anion_charge})
3498            rx_tpl = Reaction(participants=[ReactionParticipant(particle_name=reactant_name,
3499                                                                state_name=reactant_state_name,
3500                                                                coefficient=-1),
3501                                            ReactionParticipant(particle_name=product_name,
3502                                                                state_name=product_state_name,
3503                                                                coefficient=1),
3504                                            ReactionParticipant(particle_name=anion_name,
3505                                                                state_name=anion_state.name,
3506                                                                coefficient=-1),],
3507                              pK=-np.log10(gamma_K_AX.magnitude / K_XX.magnitude),
3508                              reaction_type=reaction.reaction_type+"_conjugate",
3509                              simulation_method="GRxMC")
3510            self.db._register_reaction(rx_tpl)
3511        return RE, ionic_strength_res

Sets up acid/base reactions for acidic/basic 'particles' defined in the pyMBE database, as well as a grand-canonical coupling to a reservoir of small ions using a unified formulation for small ions.

Arguments:
  • pH_res ('float'): pH-value in the reservoir.
  • c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
  • cation_name ('str'): Name of the cationic particle.
  • anion_name ('str'): Name of the anionic particle.
  • activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
  • exclusion_range('pint.Quantity', optional): Below this value, no particles will be inserted.
  • use_exclusion_radius_per_type('bool', optional): Controls if one exclusion_radius per each espresso_type. Defaults to 'False'.
Returns:

'tuple(reaction_methods.ReactionEnsemble,pint.Quantity)':

'reaction_methods.ReactionEnsemble':  
    espressomd reaction_methods object with all reactions necesary to run the GRxMC ensamble.

'pint.Quantity': 
    Ionic strength of the reservoir (useful for calculating partition coefficients).
Notes:
  • This implementation uses the formulation of the grand-reaction method by Curk et al. [1], which relies on "unified" ion types X+ = {H+, Na+} and X- = {OH-, Cl-}.
  • A function that implements the original version of the grand-reaction method by Landsgesell et al. [2] is also available under the name 'setup_grxmc_reactions'.

[1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4). [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.

def setup_lj_interactions( self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot'):
3513    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot'):
3514        """
3515        Sets up the Lennard-Jones (LJ) potential between all pairs of particle states defined in the pyMBE database.
3516
3517        Args:
3518            espresso_system('espressomd.system.System'): 
3519                Instance of a system object from the espressomd library.
3520
3521            shift_potential('bool', optional): 
3522                If True, a shift will be automatically computed such that the potential is continuous at the cutoff radius. Otherwise, no shift will be applied. Defaults to True.
3523
3524            combining_rule('string', optional): 
3525                combining rule used to calculate 'sigma' and 'epsilon' for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3526
3527            warning('bool', optional): 
3528                switch to activate/deactivate warning messages. Defaults to True.
3529
3530        Notes:
3531            - Currently, the only 'combining_rule' supported is Lorentz-Berthelot.
3532            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3533
3534        """
3535        from itertools import combinations_with_replacement
3536        particle_templates = self.db.get_templates("particle")
3537        shift = "auto" if shift_potential else 0
3538        if shift == "auto":
3539            shift_tpl = shift
3540        else:
3541            shift_tpl = PintQuantity.from_quantity(q=shift*self.units.reduced_length,
3542                                                   expected_dimension="length",
3543                                                   ureg=self.units)
3544        # Get all particle states registered in pyMBE
3545        state_entries = []
3546        for tpl in particle_templates.values():
3547            for state in self.db.get_particle_states_templates(particle_name=tpl.name).values():
3548                state_entries.append((tpl, state))
3549
3550        # Iterate over all unique state pairs
3551        for (tpl1, state1), (tpl2, state2) in combinations_with_replacement(state_entries, 2):
3552
3553            lj_parameters = self.get_lj_parameters(particle_name1=tpl1.name,
3554                                                   particle_name2=tpl2.name,
3555                                                   combining_rule=combining_rule)
3556            if not lj_parameters:
3557                continue
3558
3559            espresso_system.non_bonded_inter[state1.es_type, state2.es_type].lennard_jones.set_params(
3560                epsilon=lj_parameters["epsilon"].to("reduced_energy").magnitude,
3561                sigma=lj_parameters["sigma"].to("reduced_length").magnitude,
3562                cutoff=lj_parameters["cutoff"].to("reduced_length").magnitude,
3563                offset=lj_parameters["offset"].to("reduced_length").magnitude,
3564                shift=shift)
3565                
3566            lj_template = LJInteractionTemplate(state1=state1.name,
3567                                                state2=state2.name,
3568                                                sigma=PintQuantity.from_quantity(q=lj_parameters["sigma"],
3569                                                                                 expected_dimension="length",
3570                                                                                 ureg=self.units),
3571                                                epsilon=PintQuantity.from_quantity(q=lj_parameters["epsilon"],
3572                                                                                   expected_dimension="energy",
3573                                                                                   ureg=self.units),
3574                                                cutoff=PintQuantity.from_quantity(q=lj_parameters["cutoff"],
3575                                                                                  expected_dimension="length",
3576                                                                                  ureg=self.units),
3577                                                offset=PintQuantity.from_quantity(q=lj_parameters["offset"],
3578                                                                                  expected_dimension="length",
3579                                                                                  ureg=self.units),
3580                                                shift=shift_tpl)
3581            self.db._register_template(lj_template)

Sets up the Lennard-Jones (LJ) potential between all pairs of particle states defined in the pyMBE database.

Arguments:
  • espresso_system('espressomd.system.System'): Instance of a system object from the espressomd library.
  • shift_potential('bool', optional): If True, a shift will be automatically computed such that the potential is continuous at the cutoff radius. Otherwise, no shift will be applied. Defaults to True.
  • combining_rule('string', optional): combining rule used to calculate 'sigma' and 'epsilon' for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
  • warning('bool', optional): switch to activate/deactivate warning messages. Defaults to True.
Notes: