pyMBE

   1#
   2# Copyright (C) 2023-2025 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 sys
  21import json
  22import pint
  23import numpy as np
  24import pandas as pd
  25import scipy.constants
  26import scipy.optimize
  27from lib.lattice import LatticeBuilder, DiamondLattice
  28import logging
  29
  30class pymbe_library():
  31    """
  32    The library for the Molecular Builder for ESPResSo (pyMBE)
  33
  34    Attributes:
  35        N_A(`pint.Quantity`): Avogadro number.
  36        Kb(`pint.Quantity`): Boltzmann constant.
  37        e(`pint.Quantity`): Elementary charge.
  38        df(`Pandas.Dataframe`): Dataframe used to bookkeep all the information stored in pyMBE. Typically refered as `pmb.df`. 
  39        kT(`pint.Quantity`): Thermal energy.
  40        Kw(`pint.Quantity`): Ionic product of water. Used in the setup of the G-RxMC method.
  41    """
  42
  43    class NumpyEncoder(json.JSONEncoder):
  44        """
  45        Custom JSON encoder that converts NumPy arrays to Python lists
  46        and NumPy scalars to Python scalars.
  47        """
  48        def default(self, obj):
  49            if isinstance(obj, np.ndarray):
  50                return obj.tolist()
  51            if isinstance(obj, np.generic):
  52                return obj.item()
  53            return super().default(obj)
  54
  55    def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None):
  56        """
  57        Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 
  58        and sets up  the `pmb.df` for bookkeeping.
  59
  60        Args:
  61            temperature(`pint.Quantity`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None.
  62            unit_length(`pint.Quantity`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None.
  63            unit_charge (`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
  64            Kw (`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
  65        
  66        Note:
  67            - If no `temperature` is given, a value of 298.15 K is assumed by default.
  68            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
  69            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
  70            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
  71        """
  72        # Seed and RNG
  73        self.seed=seed
  74        self.rng = np.random.default_rng(seed)
  75        self.units=pint.UnitRegistry()
  76        self.N_A=scipy.constants.N_A / self.units.mol
  77        self.kB=scipy.constants.k * self.units.J / self.units.K
  78        self.e=scipy.constants.e * self.units.C
  79        self.set_reduced_units(unit_length=unit_length, 
  80                               unit_charge=unit_charge,
  81                               temperature=temperature, 
  82                               Kw=Kw, 
  83                               verbose=False)
  84        self.setup_df()
  85        self.lattice_builder = None
  86        return
  87
  88    def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False):
  89        """
  90        Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles.
  91
  92        Args:
  93            particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles
  94            particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles
  95            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False.
  96
  97        Returns:
  98            index(`int`): Row index where the bond information has been added in pmb.df.
  99        """
 100        particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0]
 101        particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0]
 102        
 103        bond_key = self.find_bond_key(particle_name1=particle_name1,
 104                                    particle_name2=particle_name2, 
 105                                    use_default_bond=use_default_bond)
 106        if not bond_key:
 107            return None
 108        self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1)
 109        indexs = np.where(self.df['name']==bond_key)
 110        index_list = list (indexs[0])
 111        used_bond_df = self.df.loc[self.df['particle_id2'].notnull()]
 112        #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict
 113        used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 )
 114        used_bond_index = used_bond_df.index.to_list()
 115        if not index_list:
 116            return None
 117        for index in index_list:
 118            if index not in used_bond_index:
 119                self.clean_df_row(index=int(index))
 120                self.df.at[index,'particle_id'] = particle_id1
 121                self.df.at[index,'particle_id2'] = particle_id2
 122                break
 123        return index
 124
 125    def add_bonds_to_espresso(self, espresso_system) :
 126        """
 127        Adds all bonds defined in `pmb.df` to `espresso_system`.
 128
 129        Args:
 130            espresso_system(`espressomd.system.System`): system object of espressomd library
 131        """
 132
 133        if 'bond' in self.df["pmb_type"].values:
 134            bond_df = self.df.loc[self.df ['pmb_type'] == 'bond']
 135            bond_list = bond_df.bond_object.values.tolist()
 136            for bond in bond_list:
 137                espresso_system.bonded_inter.add(bond)
 138
 139        else:
 140            print ('WARNING: There are no bonds defined in pymbe.df')
 141        
 142        return
 143
 144    def add_value_to_df(self,index,key,new_value, non_standard_value=False, overwrite=False):
 145        """
 146        Adds a value to a cell in the `pmb.df` DataFrame.
 147
 148        Args:
 149            index(`int`): index of the row to add the value to.
 150            key(`str`): the column label to add the value to.
 151            non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False.
 152            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
 153        """
 154
 155        token = "#protected:"
 156
 157        def protect(obj):
 158            if non_standard_value:
 159                return token + json.dumps(obj, cls=self.NumpyEncoder)
 160            return obj
 161
 162        def deprotect(obj):
 163            if non_standard_value and isinstance(obj, str) and obj.startswith(token):
 164                return json.loads(obj.removeprefix(token))
 165            return obj
 166
 167        # Make sure index is a scalar integer value
 168        index = int(index)
 169        assert isinstance(index, int), '`index` should be a scalar integer value.'
 170        idx = pd.IndexSlice
 171        if self.check_if_df_cell_has_a_value(index=index,key=key):
 172            old_value = self.df.loc[index,idx[key]]
 173            if not pd.Series([protect(old_value)]).equals(pd.Series([protect(new_value)])):
 174                name=self.df.loc[index,('name','')]
 175                pmb_type=self.df.loc[index,('pmb_type','')]
 176                logging.debug(f"You are attempting to redefine the properties of {name} of pmb_type {pmb_type}")    
 177                if overwrite:
 178                    logging.info(f'Overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}')
 179                if not overwrite:
 180                    logging.debug(f"pyMBE has preserved of the entry `{key}`: old_value = {old_value}. If you want to overwrite it with new_value = {new_value}, activate the switch overwrite = True ")
 181                    return
 182
 183        self.df.loc[index,idx[key]] = protect(new_value)
 184        if non_standard_value:
 185            self.df[key] = self.df[key].apply(deprotect)
 186        return
 187    
 188    def assign_molecule_id(self, molecule_index):
 189        """
 190        Assigns the `molecule_id` of the pmb object given by `pmb_type`
 191        
 192        Args:
 193            molecule_index(`int`): index of the current `pmb_object_type` to assign the `molecule_id`
 194        Returns:
 195            molecule_id(`int`): Id of the molecule
 196        """
 197
 198        self.clean_df_row(index=int(molecule_index))
 199        
 200        if self.df['molecule_id'].isnull().values.all():
 201            molecule_id = 0        
 202        else:
 203            molecule_id = self.df['molecule_id'].max() +1
 204
 205        self.add_value_to_df (key=('molecule_id',''),
 206                                index=int(molecule_index),
 207                                new_value=molecule_id)
 208
 209        return molecule_id
 210    
 211    def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system):
 212        """
 213        Calculates the center of the molecule with a given molecule_id.
 214
 215        Args:
 216            molecule_id(`int`): Id of the molecule whose center of mass is to be calculated.
 217            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 218        
 219        Returns:
 220            center_of_mass(`lst`): Coordinates of the center of mass.
 221        """
 222        center_of_mass = np.zeros(3)
 223        axis_list = [0,1,2]
 224        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
 225        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
 226        for pid in particle_id_list:
 227            for axis in axis_list:
 228                center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis]
 229        center_of_mass = center_of_mass / len(particle_id_list)
 230        return center_of_mass
 231
 232    def calculate_HH(self, molecule_name, pH_list=None, pka_set=None):
 233        """
 234        Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 
 235        for molecules with the name `molecule_name`.
 236
 237        Args:
 238            molecule_name(`str`): name of the molecule to calculate the ideal charge for
 239            pH_list(`lst`): pH-values to calculate. 
 240            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
 241
 242        Returns:
 243            Z_HH(`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list`
 244
 245        Note:
 246            - This function supports objects with pmb types: "molecule", "peptide" and "protein".
 247            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
 248            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
 249            - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan`  should be used.
 250        """
 251        if pH_list is None:
 252            pH_list=np.linspace(2,12,50)
 253        if pka_set is None:
 254            pka_set=self.get_pka_set() 
 255        self.check_pka_set(pka_set=pka_set)
 256        charge_number_map = self.get_charge_number_map()
 257        Z_HH=[]
 258        for pH_value in pH_list:    
 259            Z=0
 260            index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 
 261            residue_list = self.df.at [index,('residue_list','')]
 262            sequence = self.df.at [index,('sequence','')]
 263            if np.any(pd.isnull(sequence)):
 264                # Molecule has no sequence
 265                for residue in residue_list:
 266                    list_of_particles_in_residue = self.search_particles_in_residue(residue)
 267                    for particle in list_of_particles_in_residue:
 268                        if particle in pka_set.keys():
 269                            if pka_set[particle]['acidity'] == 'acidic':
 270                                psi=-1
 271                            elif pka_set[particle]['acidity']== 'basic':
 272                                psi=+1
 273                            else:
 274                                psi=0
 275                            Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value'])))                      
 276                Z_HH.append(Z)
 277            else:
 278                # Molecule has a sequence
 279                for name in sequence:
 280                    if name in pka_set.keys():
 281                        if pka_set[name]['acidity'] == 'acidic':
 282                            psi=-1
 283                        elif pka_set[name]['acidity']== 'basic':
 284                            psi=+1
 285                        else:
 286                            psi=0
 287                        Z+=psi/(1+10**(psi*(pH_value-pka_set[name]['pka_value'])))
 288                    else:
 289                        state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
 290                        Z+=charge_number_map[state_one_type]
 291                Z_HH.append(Z)
 292
 293        return Z_HH
 294
 295    def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
 296        """
 297        Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning.
 298
 299        Args:
 300            c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 301            c_salt('float'): Salt concentration in the reservoir.
 302            pH_list('lst'): List of pH-values in the reservoir. 
 303            pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}.
 304
 305        Returns:
 306            {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
 307            pH_system_list ('lst'): List of pH_values in the system.
 308            partition_coefficients_list ('lst'): List of partition coefficients of cations.
 309
 310        Note:
 311            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
 312            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
 313            - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result.
 314        """
 315        if pH_list is None:
 316            pH_list=np.linspace(2,12,50)
 317        if pka_set is None:
 318            pka_set=self.get_pka_set() 
 319        self.check_pka_set(pka_set=pka_set)
 320
 321        partition_coefficients_list = []
 322        pH_system_list = []
 323        Z_HH_Donnan={}
 324        for key in c_macro:
 325            Z_HH_Donnan[key] = []
 326
 327        def calc_charges(c_macro, pH):
 328            """
 329            Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation.
 330
 331            Args:
 332                c_macro('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 333                pH('float'): pH-value that is used in the HH equation.
 334
 335            Returns:
 336                charge('dict'): {"molecule_name": charge}
 337            """
 338            charge = {}
 339            for name in c_macro:
 340                charge[name] = self.calculate_HH(name, [pH], pka_set)[0]
 341            return charge
 342
 343        def calc_partition_coefficient(charge, c_macro):
 344            """
 345            Calculates the partition coefficients of positive ions according to the ideal Donnan theory.
 346
 347            Args:
 348                charge('dict'): {"molecule_name": charge}
 349                c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 350            """
 351            nonlocal ionic_strength_res
 352            charge_density = 0.0
 353            for key in charge:
 354                charge_density += charge[key] * c_macro[key]
 355            return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude
 356
 357        for pH_value in pH_list:    
 358            # calculate the ionic strength of the reservoir
 359            if pH_value <= 7.0:
 360                ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 
 361            elif pH_value > 7.0:
 362                ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt
 363
 364            #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations
 365            #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation.
 366            #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects.
 367            equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro))
 368            logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq")
 369            partition_coefficient = 10**logxi.root
 370
 371            charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient))
 372            for key in c_macro:
 373                Z_HH_Donnan[key].append(charges_temp[key])
 374
 375            pH_system_list.append(pH_value - np.log10(partition_coefficient))
 376            partition_coefficients_list.append(partition_coefficient)
 377
 378        return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
 379
 380    def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset):
 381        """
 382        Calculates the initial bond length that is used when setting up molecules,
 383        based on the minimum of the sum of bonded and short-range (LJ) interactions.
 384        
 385        Args:
 386            bond_object(`espressomd.interactions.BondedInteractions`): instance of a bond object from espressomd library
 387            bond_type(`str`): label identifying the used bonded potential
 388            epsilon(`pint.Quantity`): LJ epsilon of the interaction between the particles
 389            sigma(`pint.Quantity`): LJ sigma of the interaction between the particles
 390            cutoff(`pint.Quantity`): cutoff-radius of the LJ interaction 
 391            offset(`pint.Quantity`): offset of the LJ interaction
 392        """    
 393        def truncated_lj_potential(x, epsilon, sigma, cutoff,offset):
 394            if x>cutoff:
 395                return 0.0
 396            else:
 397                return 4*epsilon*((sigma/(x-offset))**12-(sigma/(x-offset))**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6)
 398
 399        epsilon_red=epsilon.to('reduced_energy').magnitude
 400        sigma_red=sigma.to('reduced_length').magnitude
 401        cutoff_red=cutoff.to('reduced_length').magnitude
 402        offset_red=offset.to('reduced_length').magnitude
 403
 404        if bond_type == "harmonic":
 405            r_0 = bond_object.params.get('r_0')
 406            k = bond_object.params.get('k')
 407            l0 = scipy.optimize.minimize(lambda x: 0.5*k*(x-r_0)**2 + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red, offset_red), x0=r_0).x
 408        elif bond_type == "FENE":
 409            r_0 = bond_object.params.get('r_0')
 410            k = bond_object.params.get('k')
 411            d_r_max = bond_object.params.get('d_r_max')
 412            l0 = scipy.optimize.minimize(lambda x: -0.5*k*(d_r_max**2)*np.log(1-((x-r_0)/d_r_max)**2) + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red,offset_red), x0=1.0).x
 413        return l0
 414
 415    def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False):
 416        '''
 417        Calculates the net charge per molecule of molecules with `name` = molecule_name. 
 418        Returns the net charge per molecule and a maps with the net charge per residue and molecule.
 419
 420        Args:
 421            espresso_system(`espressomd.system.System`): system information 
 422            molecule_name(`str`): name of the molecule to calculate the net charge
 423            dimensionless(`bool'): sets if the charge is returned with a dimension or not
 424
 425        Returns:
 426            (`dict`): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }}
 427
 428        Note:
 429            - The net charge of the molecule is averaged over all molecules of type `name` 
 430            - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name`
 431        '''        
 432        valid_pmb_types = ["molecule", "protein"]
 433        pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0]
 434        if pmb_type not in valid_pmb_types:
 435            raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}")      
 436
 437        id_map = self.get_particle_id_map(object_name=molecule_name)
 438        def create_charge_map(espresso_system,id_map,label):
 439            charge_number_map={}
 440            for super_id in id_map[label].keys():
 441                if dimensionless:
 442                    net_charge=0
 443                else:
 444                    net_charge=0 * self.units.Quantity(1,'reduced_charge')
 445                for pid in id_map[label][super_id]:
 446                    if dimensionless:
 447                        net_charge+=espresso_system.part.by_id(pid).q
 448                    else:
 449                        net_charge+=espresso_system.part.by_id(pid).q * self.units.Quantity(1,'reduced_charge')
 450                charge_number_map[super_id]=net_charge
 451            return charge_number_map
 452        net_charge_molecules=create_charge_map(label="molecule_map",
 453                                                espresso_system=espresso_system,
 454                                                id_map=id_map)
 455        net_charge_residues=create_charge_map(label="residue_map",
 456                                                espresso_system=espresso_system,
 457                                                id_map=id_map)       
 458        if dimensionless:
 459            mean_charge=np.mean(np.array(list(net_charge_molecules.values())))
 460        else:
 461            mean_charge=np.mean(np.array([value.magnitude for value in net_charge_molecules.values()]))*self.units.Quantity(1,'reduced_charge')
 462        return {"mean": mean_charge, "molecules": net_charge_molecules, "residues": net_charge_residues}
 463
 464    def center_molecule_in_simulation_box(self, molecule_id, espresso_system):
 465        """
 466        Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`.
 467        
 468        Args:
 469            molecule_id(`int`): Id of the molecule to be centered.
 470            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 471        """
 472        if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0:
 473            raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.")      
 474        center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system)
 475        box_center = [espresso_system.box_l[0]/2.0,
 476                      espresso_system.box_l[1]/2.0,
 477                      espresso_system.box_l[2]/2.0]
 478        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
 479        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
 480        for pid in particle_id_list:
 481            es_pos = espresso_system.part.by_id(pid).pos
 482            espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center
 483        return 
 484
 485    def check_aminoacid_key(self, key):
 486        """
 487        Checks if `key` corresponds to a valid aminoacid letter code.
 488
 489        Args:
 490            key(`str`): key to be checked.
 491
 492        Returns:
 493            `bool`: True if `key` is a valid aminoacid letter code, False otherwise.
 494        """
 495        valid_AA_keys=['V', #'VAL'
 496                       'I', #'ILE'
 497                       'L', #'LEU'
 498                       'E', #'GLU'
 499                       'Q', #'GLN'
 500                       'D', #'ASP'
 501                       'N', #'ASN'
 502                       'H', #'HIS'
 503                       'W', #'TRP'
 504                       'F', #'PHE'
 505                       'Y', #'TYR'
 506                       'R', #'ARG' 
 507                       'K', #'LYS'
 508                       'S', #'SER'
 509                       'T', #'THR'
 510                       'M', #'MET'
 511                       'A', #'ALA'
 512                       'G', #'GLY'
 513                       'P', #'PRO'
 514                       'C'] #'CYS'
 515        if key in valid_AA_keys:
 516            return True
 517        else:
 518            return False
 519
 520    def check_dimensionality(self, variable, expected_dimensionality):
 521        """
 522        Checks if the dimensionality of `variable` matches `expected_dimensionality`.
 523
 524        Args:
 525            variable(`pint.Quantity`): Quantity to be checked.
 526            expected_dimensionality(`str`): Expected dimension of the variable.
 527
 528        Returns:
 529            (`bool`): `True` if the variable if of the expected dimensionality, `False` otherwise.
 530
 531        Note:
 532            - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality).
 533            - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"`    
 534        """
 535        correct_dimensionality=variable.check(f"{expected_dimensionality}")      
 536        if not correct_dimensionality:
 537            raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}")
 538        return correct_dimensionality
 539
 540    def check_if_df_cell_has_a_value(self, index,key):
 541        """
 542        Checks if a cell in the `pmb.df` at the specified index and column has a value.
 543
 544        Args:
 545            index(`int`): Index of the row to check.
 546            key(`str`): Column label to check.
 547
 548        Returns:
 549            `bool`: `True` if the cell has a value, `False` otherwise.
 550        """
 551        idx = pd.IndexSlice
 552        import warnings
 553        with warnings.catch_warnings():
 554            warnings.simplefilter("ignore")
 555            return not pd.isna(self.df.loc[index, idx[key]])
 556
 557    def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined):
 558        """
 559        Checks if `name` is defined in `pmb.df`.
 560
 561        Args:
 562            name(`str`): label to check if defined in `pmb.df`.
 563            pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`.
 564
 565        Returns:
 566            `bool`: `True` for success, `False` otherwise.
 567        """
 568        if name in self.df['name'].unique():
 569            current_object_type = self.df[self.df['name']==name].pmb_type.values[0]
 570            if current_object_type != pmb_type_to_be_defined:
 571                raise ValueError (f"The name {name} is already defined in the df with a pmb_type = {current_object_type}, pymMBE does not support objects with the same name but different pmb_types")
 572            return True            
 573        else:
 574            return False
 575
 576    def check_if_metal_ion(self,key):
 577        """
 578        Checks if `key` corresponds to a label of a supported metal ion.
 579
 580        Args:
 581            key(`str`): key to be checked
 582
 583        Returns:
 584            (`bool`): True if `key`  is a supported metal ion, False otherwise.
 585        """
 586        if key in self.get_metal_ions_charge_number_map().keys():
 587            return True
 588        else:
 589            return False
 590
 591    def check_pka_set(self, pka_set):
 592        """
 593        Checks that `pka_set` has the formatting expected by the pyMBE library.
 594       
 595        Args:
 596            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
 597        """
 598        required_keys=['pka_value','acidity']
 599        for required_key in required_keys:
 600            for pka_name, pka_entry in pka_set.items():
 601                if required_key not in pka_entry:
 602                    raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")')
 603        return
 604
 605    def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")):
 606        """
 607        Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a pd.NA value.
 608
 609        Args:
 610            index(`int`): Index of the row to clean.
 611            columns_keys_to_clean(`list` of `str`, optional): List with the column keys to be cleaned. Defaults to [`particle_id`, `particle_id2`, `residue_id`, `molecule_id`].
 612        """   
 613        for column_key in columns_keys_to_clean:
 614            self.add_value_to_df(key=(column_key,''),index=index,new_value=pd.NA)
 615        self.df.fillna(pd.NA, inplace=True)
 616        return
 617
 618    def convert_columns_to_original_format (self, df):
 619        """
 620        Converts the columns of the Dataframe to the original format in pyMBE.
 621        
 622        Args:
 623            df(`DataFrame`): dataframe with pyMBE information as a string  
 624        
 625        """
 626
 627        columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', ('state_one','es_type'),('state_two','es_type'),('state_one','z'),('state_two','z') ]  
 628
 629        columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset']
 630
 631        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence', 'chain_map', 'node_map']
 632
 633        for column_name in columns_dtype_int:
 634            df[column_name] = df[column_name].astype(pd.Int64Dtype())
 635            
 636        for column_name in columns_with_list_or_dict:
 637            if df[column_name].isnull().all():
 638                df[column_name] = df[column_name].astype(object)
 639            else:
 640                df[column_name] = df[column_name].apply(lambda x: json.loads(x) if pd.notnull(x) else x)
 641
 642        for column_name in columns_with_units:
 643            df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x)
 644
 645        df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x)
 646        df["l0"] = df["l0"].astype(object)
 647        df["pka"] = df["pka"].astype(object)
 648        return df
 649    
 650    def convert_str_to_bond_object (self, bond_str):
 651        """
 652        Convert a row read as a `str` to the corresponding ESPResSo bond object. 
 653
 654        Args:
 655            bond_str(`str`): string with the information of a bond object.
 656
 657        Returns:
 658            bond_object(`obj`): ESPResSo bond object.
 659
 660        Note:
 661            Current supported bonds are: HarmonicBond and FeneBond
 662        """
 663        import espressomd.interactions
 664
 665        supported_bonds = ['HarmonicBond', 'FeneBond']
 666        m = re.search(r'^([A-Za-z0-9_]+)\((\{.+\})\)$', bond_str)
 667        if m is None:
 668            raise ValueError(f'Cannot parse bond "{bond_str}"')
 669        bond = m.group(1)
 670        if bond not in supported_bonds:
 671            raise NotImplementedError(f"Bond type '{bond}' currently not implemented in pyMBE, accepted types are {supported_bonds}")
 672        params = json.loads(m.group(2))
 673        bond_id = params.pop("bond_id")
 674        bond_object = getattr(espressomd.interactions, bond)(**params)
 675        bond_object._bond_id = bond_id
 676        return bond_object
 677
 678    def copy_df_entry(self, name, column_name, number_of_copies):
 679        '''
 680        Creates 'number_of_copies' of a given 'name' in `pymbe.df`.
 681
 682        Args:
 683            name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df`
 684            column_name(`str`): Column name to use as a filter. 
 685            number_of_copies(`int`): number of copies of `name` to be created.
 686        
 687        Note:
 688            - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 
 689        '''
 690
 691        valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ]
 692        if column_name not in valid_column_names:
 693            raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}")
 694        df_by_name = self.df.loc[self.df.name == name]
 695        if number_of_copies != 1:           
 696            if df_by_name[column_name].isnull().values.any():       
 697                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True)
 698            else:
 699                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
 700                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
 701                df_by_name_repeated[column_name] = pd.NA
 702            # Concatenate the new particle rows to  `df`
 703            self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
 704        else:
 705            if not df_by_name[column_name].isnull().values.any():     
 706                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
 707                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
 708                df_by_name_repeated[column_name] = pd.NA
 709                self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
 710        return
 711
 712    def create_added_salt (self, espresso_system, cation_name, anion_name, c_salt, verbose=True):    
 713        """
 714        Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`.
 715
 716        Args:
 717            espresso_system(`espressomd.system.System`): instance of an espresso system object.
 718            cation_name(`str`): `name` of a particle with a positive charge.
 719            anion_name(`str`): `name` of a particle with a negative charge.
 720            c_salt(`float`): Salt concentration.
 721            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
 722            
 723        Returns:
 724            c_salt_calculated(`float`): Calculated salt concentration added to `espresso_system`.
 725        """
 726        cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
 727        anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
 728        if cation_name_charge <= 0:
 729            raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge)
 730        if anion_name_charge >= 0:
 731            raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge)
 732        # Calculate the number of ions in the simulation box
 733        volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3')
 734        if c_salt.check('[substance] [length]**-3'):
 735            N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude)
 736            c_salt_calculated=N_ions/(volume*self.N_A)
 737        elif c_salt.check('[length]**-3'):
 738            N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude)
 739            c_salt_calculated=N_ions/volume
 740        else:
 741            raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt)
 742        N_cation = N_ions*abs(anion_name_charge)
 743        N_anion = N_ions*abs(cation_name_charge)
 744        self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation)
 745        self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion)
 746        if verbose:
 747            if c_salt_calculated.check('[substance] [length]**-3'):
 748                print(f"\n Added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions")
 749            elif c_salt_calculated.check('[length]**-3'):
 750                print(f"\n Added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions")
 751        return c_salt_calculated
 752
 753    def create_bond_in_espresso(self, bond_type, bond_parameters):
 754        '''
 755        Creates either a harmonic or a FENE bond in ESPResSo
 756
 757        Args:
 758            bond_type(`str`): label to identify the potential to model the bond.
 759            bond_parameters(`dict`): parameters of the potential of the bond.
 760
 761        Note:
 762            Currently, only HARMONIC and FENE bonds are supported.
 763
 764            For a HARMONIC bond the dictionary must contain:
 765
 766                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
 767                using the `pmb.units` UnitRegistry.
 768                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
 769                the `pmb.units` UnitRegistry.
 770           
 771            For a FENE bond the dictionary must additionally contain:
 772                
 773                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
 774                units of length using the `pmb.units` UnitRegistry. Default 'None'.
 775
 776        Returns:
 777              bond_object (`obj`): an ESPResSo bond object
 778        '''
 779        from espressomd import interactions
 780
 781        valid_bond_types   = ["harmonic", "FENE"]
 782        
 783        if 'k' in bond_parameters:
 784            bond_magnitude     = bond_parameters['k'].to('reduced_energy / reduced_length**2')
 785        else:
 786            raise ValueError("Magnitude of the potential (k) is missing")
 787        
 788        if bond_type == 'harmonic':
 789            if 'r_0' in bond_parameters:
 790                bond_length        = bond_parameters['r_0'].to('reduced_length')
 791            else:
 792                raise ValueError("Equilibrium bond length (r_0) is missing")
 793            bond_object    = interactions.HarmonicBond(k   = bond_magnitude.magnitude,
 794                                                       r_0 = bond_length.magnitude)
 795        elif bond_type == 'FENE':
 796            if 'r_0' in bond_parameters:
 797                bond_length        = bond_parameters['r_0'].to('reduced_length').magnitude
 798            else:
 799                print("WARNING: No value provided for r_0. Defaulting to r_0 = 0")
 800                bond_length=0
 801            if 'd_r_max' in bond_parameters:
 802                max_bond_stret = bond_parameters['d_r_max'].to('reduced_length')
 803            else:
 804                raise ValueError("Maximal stretching length (d_r_max) is missing")
 805            bond_object    = interactions.FeneBond(r_0     = bond_length, 
 806                                                   k       = bond_magnitude.magnitude,
 807                                                   d_r_max = max_bond_stret.magnitude)
 808        else:
 809            raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}")
 810        return bond_object
 811
 812
 813    def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True):
 814        """
 815        Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`.
 816        
 817        Args:
 818            object_name(`str`): `name` of a pymbe object.
 819            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 820            cation_name(`str`): `name` of a particle with a positive charge.
 821            anion_name(`str`): `name` of a particle with a negative charge.
 822            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
 823
 824        Returns: 
 825            counterion_number(`dict`): {"name": number}
 826        """
 827        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0]
 828        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0]
 829        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
 830        counterion_number={}
 831        object_charge={}
 832        for name in ['positive', 'negative']:
 833            object_charge[name]=0
 834        for id in object_ids:
 835            if espresso_system.part.by_id(id).q > 0:
 836                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 837            elif espresso_system.part.by_id(id).q < 0:
 838                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 839        if object_charge['positive'] % abs(anion_charge) == 0:
 840            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
 841        else:
 842            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
 843        if object_charge['negative'] % abs(cation_charge) == 0:
 844            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
 845        else:
 846            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
 847        if counterion_number[cation_name] > 0: 
 848            self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name])
 849        else:
 850            counterion_number[cation_name]=0
 851        if counterion_number[anion_name] > 0:
 852            self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name])
 853        else:
 854            counterion_number[anion_name] = 0
 855        if verbose:
 856            print('The following counter-ions have been created: ')
 857            for name in counterion_number.keys():
 858                print(f'Ion type: {name} created number: {counterion_number[name]}')
 859        return counterion_number
 860
 861    def create_hydrogel(self, name, espresso_system):
 862        """ 
 863        creates the hydrogel `name` in espresso_system
 864        Args:
 865            name(`str`): Label of the hydrogel to be created. `name` must be defined in the `pmb.df`
 866            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 867
 868        Returns:
 869            hydrogel_info(`dict`):  {"name":hydrogel_name, "chains": {chain_id1: {residue_id1: {'central_bead_id': central_bead_id, 'side_chain_ids': [particle_id1,...]},...,"node_start":node_start,"node_end":node_end}, chain_id2: {...},...}, "nodes":{node1:[node1_id],...}}     
 870        """
 871        if not self.check_if_name_is_defined_in_df(name=name, pmb_type_to_be_defined='hydrogel'):
 872            raise ValueError(f"Hydrogel with name '{name}' is not defined in the DataFrame.")
 873        hydrogel_info={"name":name, "chains":{}, "nodes":{}}
 874        # placing nodes
 875        node_positions = {}
 876        node_topology = self.df[self.df["name"]==name]["node_map"].iloc[0]
 877        for node_info in node_topology:
 878            node_index = node_info["lattice_index"]
 879            node_name = node_info["particle_name"]
 880            node_pos, node_id = self.create_hydrogel_node(self.format_node(node_index), node_name, espresso_system)
 881            hydrogel_info["nodes"][self.format_node(node_index)]=node_id
 882            node_positions[node_id[0]]=node_pos
 883        
 884        # Placing chains between nodes
 885        # Looping over all the 16 chains
 886        chain_topology = self.df[self.df["name"]==name]["chain_map"].iloc[0]
 887        for chain_info in chain_topology:
 888            node_s = chain_info["node_start"]
 889            node_e = chain_info["node_end"]
 890            molecule_info = self.create_hydrogel_chain(node_s, node_e, node_positions, espresso_system)
 891            for molecule_id in molecule_info:
 892                hydrogel_info["chains"][molecule_id] = molecule_info[molecule_id]
 893                hydrogel_info["chains"][molecule_id]["node_start"]=node_s
 894                hydrogel_info["chains"][molecule_id]["node_end"]=node_e
 895        return hydrogel_info
 896
 897    def create_hydrogel_chain(self, node_start, node_end, node_positions, espresso_system):
 898        """
 899        Creates a chain between two nodes of a hydrogel.
 900
 901        Args:
 902            node_start(`str`): name of the starting node particle at which the first residue of the chain will be attached. 
 903            node_end(`str`): name of the ending node particle at which the last residue of the chain will be attached.
 904            node_positions(`dict`): dictionary with the positions of the nodes. The keys are the node names and the values are the positions.
 905            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 906
 907        Note:
 908            For example, 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.
 909            The chain will be placed in the direction of the vector between `node_start` and `node_end`. 
 910        """
 911        if self.lattice_builder is None:
 912            raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.")
 913
 914        molecule_name = "chain_"+node_start+"_"+node_end
 915        sequence = self.df[self.df['name']==molecule_name].residue_list.values [0]
 916        assert len(sequence) != 0 and not isinstance(sequence, str)
 917        assert len(sequence) == self.lattice_builder.mpc
 918
 919        key, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end)
 920        assert node_start != node_end or sequence == sequence[::-1], \
 921            (f"chain cannot be defined between '{node_start}' and '{node_end}' since it "
 922            "would form a loop with a non-symmetric sequence (under-defined stereocenter)")
 923
 924        if reverse:
 925            sequence = sequence[::-1]
 926
 927        node_start_pos = np.array(list(int(x) for x in node_start.strip('[]').split()))*0.25*self.lattice_builder.box_l
 928        node_end_pos = np.array(list(int(x) for x in node_end.strip('[]').split()))*0.25*self.lattice_builder.box_l
 929        node1 = espresso_system.part.select(lambda p: (p.pos == node_start_pos).all()).id
 930        node2 = espresso_system.part.select(lambda p: (p.pos == node_end_pos).all()).id
 931
 932        if not node1[0] in node_positions or not node2[0] in node_positions:
 933            raise ValueError("Set node position before placing a chain between them")
 934         
 935        # Finding a backbone vector between node_start and node_end
 936        vec_between_nodes = np.array(node_positions[node2[0]]) - np.array(node_positions[node1[0]])
 937        vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l)
 938        backbone_vector = list(vec_between_nodes/(self.lattice_builder.mpc + 1))
 939        node_start_name = self.df[(self.df["particle_id"]==node1[0]) & (self.df["pmb_type"]=="particle")]["name"].values[0]
 940        first_res_name = self.df[(self.df["pmb_type"]=="residue") & (self.df["name"]==sequence[0])]["central_bead"].values[0]
 941        l0 = self.get_bond_length(node_start_name, first_res_name, hard_check=True)
 942        chain_molecule_info = self.create_molecule(
 943                name=molecule_name,  # Use the name defined earlier
 944                number_of_molecules=1,  # Creating one chain
 945                espresso_system=espresso_system,
 946                list_of_first_residue_positions=[list(np.array(node_positions[node1[0]]) + np.array(backbone_vector))],#Start at the first node
 947                backbone_vector=np.array(backbone_vector)/l0,
 948                use_default_bond=False  # Use defaut bonds between monomers
 949            )
 950        # Collecting ids of beads of the chain/molecule
 951        chain_ids = []
 952        residue_ids = []
 953        for molecule_id in chain_molecule_info:
 954            for residue_id in chain_molecule_info[molecule_id]:
 955                residue_ids.append(residue_id)
 956                bead_id = chain_molecule_info[molecule_id][residue_id]['central_bead_id']
 957                chain_ids.append(bead_id)
 958
 959        self.lattice_builder.chains[key] = sequence
 960        # Search bonds between nodes and chain ends
 961        BeadType_near_to_node_start = self.df[(self.df["residue_id"] == residue_ids[0]) & (self.df["central_bead"].notnull())]["central_bead"].drop_duplicates().iloc[0]
 962        BeadType_near_to_node_end = self.df[(self.df["residue_id"] == residue_ids[-1]) & (self.df["central_bead"].notnull())]["central_bead"].drop_duplicates().iloc[0]
 963        bond_node1_first_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_start],
 964                                                    particle_name2 = BeadType_near_to_node_start,
 965                                                    hard_check=False,
 966                                                    use_default_bond=False)
 967        bond_node2_last_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_end],
 968                                                    particle_name2 = BeadType_near_to_node_end,
 969                                                    hard_check=False,
 970                                                    use_default_bond=False)
 971
 972        espresso_system.part.by_id(node1[0]).add_bond((bond_node1_first_monomer, chain_ids[0]))
 973        espresso_system.part.by_id(node2[0]).add_bond((bond_node2_last_monomer, chain_ids[-1]))
 974        # Add bonds to data frame
 975        bond_index1 = self.add_bond_in_df(particle_id1=node1[0],
 976                                    particle_id2=chain_ids[0],
 977                                    use_default_bond=False)
 978        self.add_value_to_df(key=('molecule_id',''),
 979                        index=int(bond_index1),
 980                        new_value=molecule_id,
 981                        overwrite=True)
 982        self.add_value_to_df(key=('residue_id',''),
 983                            index=int (bond_index1),
 984                            new_value=residue_ids[0],
 985                            overwrite=True)
 986        bond_index2 = self.add_bond_in_df(particle_id1=node2[0],
 987                                    particle_id2=chain_ids[-1],
 988                                    use_default_bond=False)
 989        self.add_value_to_df(key=('molecule_id',''),
 990                        index=int(bond_index2),
 991                        new_value=molecule_id,
 992                        overwrite=True)
 993        self.add_value_to_df(key=('residue_id',''),
 994                        index=int (bond_index2),
 995                        new_value=residue_ids[-1],
 996                        overwrite=True)
 997        return chain_molecule_info
 998    
 999    def create_hydrogel_node(self, node_index, node_name, espresso_system):
1000        """
1001        Set a node residue type.
1002        
1003        Args:
1004            node_index(`str`): Lattice node index in the form of a string, e.g. "[0 0 0]".
1005            node_name(`str`): name of the node particle defined in pyMBE.
1006        Returns:
1007            node_position(`list`): Position of the node in the lattice.
1008            p_id(`int`): Particle ID of the node.
1009        """
1010        if self.lattice_builder is None:
1011            raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.")
1012
1013        node_position = np.array(list(int(x) for x in node_index.strip('[]').split()))*0.25*self.lattice_builder.box_l
1014        p_id = self.create_particle(name = node_name,
1015                         espresso_system=espresso_system,
1016                         number_of_particles=1,
1017                         position = [node_position])
1018        key = self.lattice_builder._get_node_by_label(node_index)
1019        self.lattice_builder.nodes[key] = node_name
1020
1021        return node_position.tolist(), p_id
1022
1023    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False):
1024        """
1025        Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1026
1027        Args:
1028            name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df`
1029            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1030            number_of_molecules(`int`): Number of molecules of type `name` to be created.
1031            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.
1032            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. 
1033            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df`
1034
1035        Returns:
1036            molecules_info(`dict`):  {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 
1037        """
1038        if list_of_first_residue_positions is not None:
1039            for item in list_of_first_residue_positions:
1040                if not isinstance(item, list):
1041                    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.")
1042                elif len(item) != 3:
1043                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
1044
1045            if len(list_of_first_residue_positions) != number_of_molecules:
1046                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
1047        if number_of_molecules <= 0:
1048            return 0
1049        # Generate an arbitrary random unit vector
1050        if backbone_vector is None:
1051            backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
1052                                                    radius=1, 
1053                                                    n_samples=1,
1054                                                    on_surface=True)[0]
1055        else:
1056            backbone_vector = np.array(backbone_vector)
1057        self.check_if_name_is_defined_in_df(name=name,
1058                                            pmb_type_to_be_defined='molecule')
1059            
1060        first_residue = True
1061        molecules_info = {}
1062        residue_list = self.df[self.df['name']==name].residue_list.values [0]
1063        self.copy_df_entry(name=name,
1064                        column_name='molecule_id',
1065                        number_of_copies=number_of_molecules)
1066
1067        molecules_index = np.where(self.df['name']==name)
1068        molecule_index_list =list(molecules_index[0])[-number_of_molecules:]
1069        pos_index = 0 
1070        for molecule_index in molecule_index_list:        
1071            molecule_id = self.assign_molecule_id(molecule_index=molecule_index)
1072            molecules_info[molecule_id] = {}
1073            for residue in residue_list:
1074                if first_residue:
1075                    if list_of_first_residue_positions is None:
1076                        residue_position = None
1077                    else:
1078                        for item in list_of_first_residue_positions:
1079                            residue_position = [np.array(list_of_first_residue_positions[pos_index])]
1080                            
1081                    residues_info = self.create_residue(name=residue,
1082                                                        espresso_system=espresso_system, 
1083                                                        central_bead_position=residue_position,  
1084                                                        use_default_bond= use_default_bond, 
1085                                                        backbone_vector=backbone_vector)
1086                    residue_id = next(iter(residues_info))
1087                    # Add the correct molecule_id to all particles in the residue
1088                    for index in self.df[self.df['residue_id']==residue_id].index:
1089                        self.add_value_to_df(key=('molecule_id',''),
1090                                            index=int (index),
1091                                            new_value=molecule_id,
1092                                            overwrite=True)
1093                    central_bead_id = residues_info[residue_id]['central_bead_id']
1094                    previous_residue = residue
1095                    residue_position = espresso_system.part.by_id(central_bead_id).pos
1096                    previous_residue_id = central_bead_id
1097                    first_residue = False          
1098                else:                    
1099                    previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0]
1100                    new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0]
1101                    bond = self.search_bond(particle_name1=previous_central_bead_name, 
1102                                            particle_name2=new_central_bead_name, 
1103                                            hard_check=True, 
1104                                            use_default_bond=use_default_bond)                
1105                    l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 
1106                                            particle_name2=new_central_bead_name, 
1107                                            hard_check=True, 
1108                                            use_default_bond=use_default_bond)                
1109                    
1110                    residue_position = residue_position+backbone_vector*l0
1111                    residues_info = self.create_residue(name=residue, 
1112                                                        espresso_system=espresso_system, 
1113                                                        central_bead_position=[residue_position],
1114                                                        use_default_bond= use_default_bond, 
1115                                                        backbone_vector=backbone_vector)
1116                    residue_id = next(iter(residues_info))      
1117                    for index in self.df[self.df['residue_id']==residue_id].index:
1118                        self.add_value_to_df(key=('molecule_id',''),
1119                                            index=int (index),
1120                                            new_value=molecule_id,
1121                                            overwrite=True)            
1122                    central_bead_id = residues_info[residue_id]['central_bead_id']
1123                    espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id))
1124                    bond_index = self.add_bond_in_df(particle_id1=central_bead_id,
1125                                        particle_id2=previous_residue_id,
1126                                        use_default_bond=use_default_bond) 
1127                    self.add_value_to_df(key=('molecule_id',''),
1128                                            index=int (bond_index),
1129                                            new_value=molecule_id,
1130                                            overwrite=True)           
1131                    previous_residue_id = central_bead_id
1132                    previous_residue = residue                    
1133                molecules_info[molecule_id][residue_id] = residues_info[residue_id]
1134            first_residue = True
1135            pos_index+=1
1136        
1137        return molecules_info
1138    
1139    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
1140        """
1141        Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`.
1142        
1143        Args:
1144            name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 
1145            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1146            number_of_particles(`int`): Number of particles to be created.
1147            position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
1148            fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
1149        Returns:
1150            created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`.
1151        """       
1152        if number_of_particles <=0:
1153            return []
1154        self.check_if_name_is_defined_in_df(name=name,
1155                                       pmb_type_to_be_defined='particle')
1156        # Copy the data of the particle `number_of_particles` times in the `df`
1157        self.copy_df_entry(name=name,
1158                          column_name='particle_id',
1159                          number_of_copies=number_of_particles)
1160        # Get information from the particle type `name` from the df     
1161        z = self.df.loc[self.df['name']==name].state_one.z.values[0]
1162        z = 0. if z is None else z
1163        es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
1164        # Get a list of the index in `df` corresponding to the new particles to be created
1165        index = np.where(self.df['name']==name)
1166        index_list =list(index[0])[-number_of_particles:]
1167        # Create the new particles into  `espresso_system`
1168        created_pid_list=[]
1169        for index in range (number_of_particles):
1170            df_index=int (index_list[index])
1171            self.clean_df_row(index=df_index)
1172            if position is None:
1173                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1174            else:
1175                particle_position = position[index]
1176            if len(espresso_system.part.all()) == 0:
1177                bead_id = 0
1178            else:
1179                bead_id = max (espresso_system.part.all().id) + 1
1180            created_pid_list.append(bead_id)
1181            kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z)
1182            if fix:
1183                kwargs["fix"] = 3 * [fix]
1184            espresso_system.part.add(**kwargs)
1185            self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id)                  
1186        return created_pid_list
1187
1188    def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None):
1189        """
1190        Creates all `particle`s associated to `pmb object` into  `espresso` a number of times equal to `number_of_objects`.
1191        
1192        Args:
1193            name(`str`): Unique label of the `pmb object` to be created. 
1194            number_of_objects(`int`): Number of `pmb object`s to be created.
1195            espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library.
1196            position(`list`): Coordinates where the particles should be created.
1197            use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`.
1198            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. 
1199
1200        Note:
1201            - If no `position` is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length. 
1202        """
1203        allowed_objects=['particle','residue','molecule']
1204        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1205        if pmb_type not in allowed_objects:
1206            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1207        if pmb_type == 'particle':
1208            self.create_particle(name=name, 
1209                                number_of_particles=number_of_objects, 
1210                                espresso_system=espresso_system, 
1211                                position=position)
1212        elif pmb_type == 'residue':
1213            self.create_residue(name=name,  
1214                                espresso_system=espresso_system, 
1215                                central_bead_position=position,
1216                                use_default_bond=use_default_bond,
1217                                backbone_vector=backbone_vector)
1218        elif pmb_type == 'molecule':
1219            self.create_molecule(name=name, 
1220                                number_of_molecules=number_of_objects, 
1221                                espresso_system=espresso_system, 
1222                                use_default_bond=use_default_bond, 
1223                                list_of_first_residue_positions=position,
1224                                backbone_vector=backbone_vector)
1225        return
1226
1227    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1228        """
1229        Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions`
1230
1231        Args:
1232            name(`str`): Label of the protein to be created. 
1233            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1234            number_of_proteins(`int`): Number of proteins to be created.
1235            positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}}
1236        """
1237
1238        if number_of_proteins <=0:
1239            return
1240        self.check_if_name_is_defined_in_df(name=name,
1241                                        pmb_type_to_be_defined='protein')
1242        self.copy_df_entry(name=name,
1243                            column_name='molecule_id',
1244                            number_of_copies=number_of_proteins)
1245        protein_index = np.where(self.df['name']==name)
1246        protein_index_list =list(protein_index[0])[-number_of_proteins:]
1247        
1248        box_half=espresso_system.box_l[0]/2.0
1249
1250        for molecule_index in protein_index_list:     
1251
1252            molecule_id = self.assign_molecule_id(molecule_index=molecule_index)
1253
1254            protein_center = self.generate_coordinates_outside_sphere(radius = 1, 
1255                                                                        max_dist=box_half, 
1256                                                                        n_samples=1, 
1257                                                                        center=[box_half]*3)[0]
1258   
1259            for residue in topology_dict.keys():
1260
1261                residue_name = re.split(r'\d+', residue)[0]
1262                residue_number = re.split(r'(\d+)', residue)[1]
1263                residue_position = topology_dict[residue]['initial_pos']
1264                position = residue_position + protein_center
1265
1266                particle_id = self.create_particle(name=residue_name,
1267                                                            espresso_system=espresso_system,
1268                                                            number_of_particles=1,
1269                                                            position=[position], 
1270                                                            fix = True)
1271                
1272                index = self.df[self.df['particle_id']==particle_id[0]].index.values[0]
1273                self.add_value_to_df(key=('residue_id',''),
1274                                            index=int (index),
1275                                            new_value=int(residue_number),
1276                                            overwrite=True)
1277
1278                self.add_value_to_df(key=('molecule_id',''),
1279                                        index=int (index),
1280                                        new_value=molecule_id,
1281                                        overwrite=True)
1282
1283        return
1284
1285    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None):
1286        """
1287        Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1288
1289        Args:
1290            name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df`
1291            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1292            central_bead_position(`list` of `float`): Position of the central bead.
1293            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 `pmb.df`
1294            backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`.
1295
1296        Returns:
1297            residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1298        """
1299        self.check_if_name_is_defined_in_df(name=name,
1300                                            pmb_type_to_be_defined='residue')
1301        # Copy the data of a residue in the `df
1302        self.copy_df_entry(name=name,
1303                            column_name='residue_id',
1304                            number_of_copies=1)
1305        residues_index = np.where(self.df['name']==name)
1306        residue_index_list =list(residues_index[0])[-1:]
1307        # search for defined particle and residue names
1308        particle_and_residue_df = self.df.loc[(self.df['pmb_type']== "particle") | (self.df['pmb_type']== "residue")]
1309        particle_and_residue_names = particle_and_residue_df["name"].tolist()
1310        for residue_index in residue_index_list:  
1311            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1312            for side_chain_element in side_chain_list:
1313                if side_chain_element not in particle_and_residue_names:              
1314                    raise ValueError (f"{side_chain_element} is not defined")
1315        # Internal bookkepping of the residue info (important for side-chain residues)
1316        # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1317        residues_info={}
1318        for residue_index in residue_index_list:     
1319            self.clean_df_row(index=int(residue_index))
1320            # Assign a residue_id
1321            if self.df['residue_id'].isnull().all():
1322                residue_id=0
1323            else:
1324                residue_id = self.df['residue_id'].max() + 1
1325            self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id)
1326            # create the principal bead   
1327            central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0]            
1328            central_bead_id = self.create_particle(name=central_bead_name,
1329                                                                espresso_system=espresso_system,
1330                                                                position=central_bead_position,
1331                                                                number_of_particles = 1)[0]
1332            central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1333            #assigns same residue_id to the central_bead particle created.
1334            index = self.df[self.df['particle_id']==central_bead_id].index.values[0]
1335            self.df.at [index,'residue_id'] = residue_id
1336            # Internal bookkeeping of the central bead id
1337            residues_info[residue_id]={}
1338            residues_info[residue_id]['central_bead_id']=central_bead_id
1339            # create the lateral beads  
1340            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1341            side_chain_beads_ids = []
1342            for side_chain_element in side_chain_list:  
1343                pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 
1344                if pmb_type == 'particle':
1345                    bond = self.search_bond(particle_name1=central_bead_name, 
1346                                            particle_name2=side_chain_element, 
1347                                            hard_check=True, 
1348                                            use_default_bond=use_default_bond)
1349                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1350                                              particle_name2=side_chain_element, 
1351                                              hard_check=True, 
1352                                              use_default_bond=use_default_bond)
1353
1354                    if backbone_vector is None:
1355                        bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1356                                                                 radius=l0, 
1357                                                                 n_samples=1,
1358                                                                 on_surface=True)[0]
1359                    else:
1360                        bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector),
1361                                                                                                    magnitude=l0)
1362                     
1363                    side_bead_id = self.create_particle(name=side_chain_element, 
1364                                                                    espresso_system=espresso_system,
1365                                                                    position=[bead_position], 
1366                                                                    number_of_particles=1)[0]
1367                    index = self.df[self.df['particle_id']==side_bead_id].index.values[0]
1368                    self.add_value_to_df(key=('residue_id',''),
1369                                        index=int (index),
1370                                        new_value=residue_id, 
1371                                        overwrite=True)
1372                    side_chain_beads_ids.append(side_bead_id)
1373                    espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id))
1374                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1375                                        particle_id2=side_bead_id,
1376                                        use_default_bond=use_default_bond)
1377                    self.add_value_to_df(key=('residue_id',''),
1378                                        index=int (index),
1379                                        new_value=residue_id, 
1380                                        overwrite=True)
1381
1382                elif pmb_type == 'residue':
1383                    central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0]
1384                    bond = self.search_bond(particle_name1=central_bead_name, 
1385                                            particle_name2=central_bead_side_chain, 
1386                                            hard_check=True, 
1387                                            use_default_bond=use_default_bond)
1388                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1389                                              particle_name2=central_bead_side_chain, 
1390                                              hard_check=True, 
1391                                              use_default_bond=use_default_bond)
1392                    if backbone_vector is None:
1393                        residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1394                                                                    radius=l0, 
1395                                                                    n_samples=1,
1396                                                                    on_surface=True)[0]
1397                    else:
1398                        residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1399                                                                                                        magnitude=l0)
1400                    lateral_residue_info = self.create_residue(name=side_chain_element, 
1401                                                                espresso_system=espresso_system,
1402                                                                central_bead_position=[residue_position],
1403                                                                use_default_bond=use_default_bond)
1404                    lateral_residue_dict=list(lateral_residue_info.values())[0]
1405                    central_bead_side_chain_id=lateral_residue_dict['central_bead_id']
1406                    lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids']
1407                    residue_id_side_chain=list(lateral_residue_info.keys())[0]
1408                    # Change the residue_id of the residue in the side chain to the one of the bigger residue
1409                    index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0]
1410                    self.add_value_to_df(key=('residue_id',''),
1411                                        index=int(index),
1412                                        new_value=residue_id, 
1413                                        overwrite=True)
1414                    # Change the residue_id of the particles in the residue in the side chain
1415                    side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids
1416                    for particle_id in side_chain_beads_ids:
1417                        index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0]
1418                        self.add_value_to_df(key=('residue_id',''),
1419                                            index=int (index),
1420                                            new_value=residue_id, 
1421                                            overwrite=True)
1422                    espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id))
1423                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1424                                        particle_id2=central_bead_side_chain_id,
1425                                        use_default_bond=use_default_bond)
1426                    self.add_value_to_df(key=('residue_id',''),
1427                                        index=int (index),
1428                                        new_value=residue_id, 
1429                                        overwrite=True)
1430                    # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue
1431                    for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index:        
1432                        self.add_value_to_df(key=('residue_id',''),
1433                                            index=int(index),
1434                                            new_value=residue_id, 
1435                                            overwrite=True)
1436            # Internal bookkeeping of the side chain beads ids
1437            residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids
1438        return  residues_info
1439
1440    def create_variable_with_units(self, variable):
1441        """
1442        Returns a pint object with the value and units defined in `variable`.
1443
1444        Args:
1445            variable(`dict` or `str`): {'value': value, 'units': units}
1446        Returns:
1447            variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry.
1448        """        
1449        if isinstance(variable, dict):
1450            value=variable.pop('value')
1451            units=variable.pop('units')
1452        elif isinstance(variable, str):
1453            value = float(re.split(r'\s+', variable)[0])
1454            units = re.split(r'\s+', variable)[1]
1455        variable_with_units=value*self.units(units)
1456        return variable_with_units 
1457
1458    def define_AA_residues(self, sequence, model):
1459        """
1460        Defines in `pmb.df` all the different residues in `sequence`.
1461
1462        Args:
1463            sequence(`lst`):  Sequence of the peptide or protein.
1464            model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1465
1466        Returns:
1467            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.             
1468        """
1469        residue_list = []
1470        for residue_name in sequence:
1471            if model == '1beadAA':
1472                central_bead = residue_name
1473                side_chains = []
1474            elif model == '2beadAA':
1475                if residue_name in ['c','n', 'G']: 
1476                    central_bead = residue_name
1477                    side_chains = []
1478                else:
1479                    central_bead = 'CA'              
1480                    side_chains = [residue_name]
1481            if residue_name not in residue_list:   
1482                self.define_residue(name = 'AA-'+residue_name, 
1483                                    central_bead = central_bead,
1484                                    side_chains = side_chains)              
1485            residue_list.append('AA-'+residue_name)
1486        return residue_list
1487
1488    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1489        
1490        '''
1491        Defines a pmb object of type `bond` in `pymbe.df`.
1492
1493        Args:
1494            bond_type(`str`): label to identify the potential to model the bond.
1495            bond_parameters(`dict`): parameters of the potential of the bond.
1496            particle_pairs(`lst`): list of the `names` of the `particles` to be bonded.
1497
1498        Note:
1499            Currently, only HARMONIC and FENE bonds are supported.
1500
1501            For a HARMONIC bond the dictionary must contain the following parameters:
1502
1503                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
1504                using the `pmb.units` UnitRegistry.
1505                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
1506                the `pmb.units` UnitRegistry.
1507           
1508            For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:
1509                
1510                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
1511                units of length using the `pmb.units` UnitRegistry. Default 'None'.
1512        '''
1513
1514        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1515        for particle_name1, particle_name2 in particle_pairs:
1516
1517            lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1,
1518                                                 particle_name2 = particle_name2,
1519                                                 combining_rule = 'Lorentz-Berthelot')
1520      
1521            l0 = self.calculate_initial_bond_length(bond_object = bond_object,
1522                                                    bond_type   = bond_type,
1523                                                    epsilon     = lj_parameters["epsilon"],
1524                                                    sigma       = lj_parameters["sigma"],
1525                                                    cutoff      = lj_parameters["cutoff"],
1526                                                    offset      = lj_parameters["offset"],)
1527            index = len(self.df)
1528            for label in [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']:
1529                self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond")
1530            self.df.at [index,'name']= f'{particle_name1}-{particle_name2}'
1531            self.df.at [index,'bond_object'] = bond_object
1532            self.df.at [index,'l0'] = l0
1533            self.add_value_to_df(index=index,
1534                                    key=('pmb_type',''),
1535                                    new_value='bond')
1536            self.add_value_to_df(index=index,
1537                                    key=('parameters_of_the_potential',''),
1538                                    new_value=bond_object.get_params(),
1539                                    non_standard_value=True)
1540        self.df.fillna(pd.NA, inplace=True)
1541        return
1542
1543    
1544    def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None):
1545        """
1546        Asigns `bond` in `pmb.df` as the default bond.
1547        The LJ parameters can be optionally provided to calculate the initial bond length
1548
1549        Args:
1550            bond_type(`str`): label to identify the potential to model the bond.
1551            bond_parameters(`dict`): parameters of the potential of the bond.
1552            sigma(`float`, optional): LJ sigma of the interaction between the particles.
1553            epsilon(`float`, optional): LJ epsilon for the interaction between the particles.
1554            cutoff(`float`, optional): cutoff-radius of the LJ interaction.
1555            offset(`float`, optional): offset of the LJ interaction.
1556
1557        Note:
1558            - Currently, only harmonic and FENE bonds are supported. 
1559        """
1560        
1561        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1562
1563        if epsilon is None:
1564            epsilon=1*self.units('reduced_energy')
1565        if sigma is None:
1566            sigma=1*self.units('reduced_length')
1567        if cutoff is None:
1568            cutoff=2**(1.0/6.0)*self.units('reduced_length')
1569        if offset is None:
1570            offset=0*self.units('reduced_length')
1571        l0 = self.calculate_initial_bond_length(bond_object = bond_object, 
1572                                                bond_type   = bond_type, 
1573                                                epsilon     = epsilon,
1574                                                sigma       = sigma,
1575                                                cutoff      = cutoff,
1576                                                offset      = offset)
1577
1578        if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'):
1579            return
1580        if len(self.df.index) != 0:
1581            index = max(self.df.index)+1
1582        else:
1583            index = 0
1584        self.df.at [index,'name']        = 'default'
1585        self.df.at [index,'bond_object'] = bond_object
1586        self.df.at [index,'l0']          = l0
1587        self.add_value_to_df(index       = index,
1588                            key          = ('pmb_type',''),
1589                            new_value    = 'bond')
1590        self.add_value_to_df(index       = index,
1591                            key          = ('parameters_of_the_potential',''),
1592                            new_value    = bond_object.get_params(),
1593                            non_standard_value=True)
1594        self.df.fillna(pd.NA, inplace=True)
1595        return
1596    
1597    def define_hydrogel(self, name, node_map, chain_map):
1598        """
1599        Defines a pyMBE object of type `hydrogel` in `pymbe.df`.
1600
1601        Args:
1602            name(`str`): Unique label that identifies the `hydrogel`.
1603            node_map(`list of ict`): [{"particle_name": , "lattice_index": }, ... ]
1604            chain_map(`list of dict`): [{"node_start": , "node_end": , "residue_list": , ... ]
1605        """
1606        node_indices = {tuple(entry['lattice_index']) for entry in node_map}
1607        diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices}
1608        if node_indices != diamond_indices:
1609            raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ")
1610        
1611        chain_map_connectivity = set()
1612        for entry in chain_map:
1613            start = self.lattice_builder.node_labels[entry['node_start']]
1614            end = self.lattice_builder.node_labels[entry['node_end']]
1615            chain_map_connectivity.add((start,end))
1616
1617        if self.lattice_builder.lattice.connectivity != chain_map_connectivity:
1618            raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs")
1619
1620        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='hydrogel'):
1621            return
1622
1623        index = len(self.df)
1624        self.df.at [index, "name"] = name
1625        self.df.at [index, "pmb_type"] = "hydrogel"
1626        self.add_value_to_df(index = index,
1627                        key = ('node_map',''),
1628                        new_value = node_map,
1629                        non_standard_value=True)
1630        self.add_value_to_df(index = index,
1631                        key = ('chain_map',''),
1632                        new_value = chain_map,
1633                        non_standard_value=True)
1634        for chain_label in chain_map:
1635            node_start = chain_label["node_start"]
1636            node_end = chain_label["node_end"]
1637            residue_list = chain_label['residue_list']
1638            # Molecule name
1639            molecule_name = "chain_"+node_start+"_"+node_end
1640            self.define_molecule(name=molecule_name, residue_list=residue_list)
1641        return
1642
1643    def define_molecule(self, name, residue_list):
1644        """
1645        Defines a pyMBE object of type `molecule` in `pymbe.df`.
1646
1647        Args:
1648            name(`str`): Unique label that identifies the `molecule`.
1649            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.  
1650        """
1651        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'):
1652            return
1653        index = len(self.df)
1654        self.df.at [index,'name'] = name
1655        self.df.at [index,'pmb_type'] = 'molecule'
1656        self.df.at [index,('residue_list','')] = residue_list
1657        self.df.fillna(pd.NA, inplace=True)
1658        return
1659
1660    def define_particle_entry_in_df(self,name):
1661        """
1662        Defines a particle entry in pmb.df.
1663
1664        Args:
1665            name(`str`): Unique label that identifies this particle type.
1666
1667        Returns:
1668            index(`int`): Index of the particle in pmb.df  
1669        """
1670
1671        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'):
1672            index = self.df[self.df['name']==name].index[0]                                   
1673        else:
1674            index = len(self.df)
1675            self.df.at [index, 'name'] = name
1676            self.df.at [index,'pmb_type'] = 'particle'
1677        self.df.fillna(pd.NA, inplace=True)
1678        return index
1679
1680    def define_particle(self, name, z=0, acidity=pd.NA, pka=pd.NA, sigma=pd.NA, epsilon=pd.NA, cutoff=pd.NA, offset=pd.NA,overwrite=False):
1681        """
1682        Defines the properties of a particle object.
1683
1684        Args:
1685            name(`str`): Unique label that identifies this particle type.  
1686            z(`int`, optional): Permanent charge number of this particle type. Defaults to 0.
1687            acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to pd.NA.
1688            pka(`float`, optional): If `particle` is an acid or a base, it defines its  pka-value. Defaults to pd.NA.
1689            sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1690            cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1691            offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1692            epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA.
1693            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
1694
1695        Note:
1696            - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units.
1697            - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units.
1698            - `cutoff` defaults to `2**(1./6.) reduced_length`. 
1699            - `offset` defaults to 0.
1700            - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions.
1701            - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`.
1702        """ 
1703        index=self.define_particle_entry_in_df(name=name)
1704        
1705        # If `cutoff` and `offset` are not defined, default them to the following values
1706        if pd.isna(cutoff):
1707            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1708        if pd.isna(offset):
1709            offset=self.units.Quantity(0, "reduced_length")
1710
1711        # Define LJ parameters
1712        parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"},
1713                                        "cutoff":{"value": cutoff, "dimensionality": "[length]"},
1714                                        "offset":{"value": offset, "dimensionality": "[length]"},
1715                                        "epsilon":{"value": epsilon, "dimensionality": "[energy]"},}
1716
1717        for parameter_key in parameters_with_dimensionality.keys():
1718            if not pd.isna(parameters_with_dimensionality[parameter_key]["value"]):
1719                self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 
1720                                          expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"])
1721                self.add_value_to_df(key=(parameter_key,''),
1722                                    index=index,
1723                                    new_value=parameters_with_dimensionality[parameter_key]["value"],
1724                                    overwrite=overwrite)
1725
1726        # Define particle acid/base properties
1727        self.set_particle_acidity(name=name, 
1728                                acidity=acidity, 
1729                                default_charge_number=z, 
1730                                pka=pka,
1731                                overwrite=overwrite)
1732        self.df.fillna(pd.NA, inplace=True)
1733        return 
1734    
1735    def define_particles(self, parameters, overwrite=False):
1736        '''
1737        Defines a particle object in pyMBE for each particle name in `particle_names`
1738
1739        Args:
1740            parameters(`dict`):  dictionary with the particle parameters. 
1741            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1742
1743        Note:
1744            - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},}
1745        '''
1746        if not parameters:
1747            return 0
1748        for particle_name in parameters.keys():
1749            parameters[particle_name]["overwrite"]=overwrite
1750            self.define_particle(**parameters[particle_name])
1751        return
1752    
1753    def define_peptide(self, name, sequence, model):
1754        """
1755        Defines a pyMBE object of type `peptide` in the `pymbe.df`.
1756
1757        Args:
1758            name (`str`): Unique label that identifies the `peptide`.
1759            sequence (`string`): Sequence of the `peptide`.
1760            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1761        """
1762        if self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'):
1763            return
1764        valid_keys = ['1beadAA','2beadAA']
1765        if model not in valid_keys:
1766            raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
1767        clean_sequence = self.protein_sequence_parser(sequence=sequence)    
1768        residue_list = self.define_AA_residues(sequence=clean_sequence,
1769                                                model=model)
1770        self.define_molecule(name = name, residue_list=residue_list)
1771        index = self.df.loc[self.df['name'] == name].index.item() 
1772        self.df.at [index,'model'] = model
1773        self.df.at [index,('sequence','')] = clean_sequence
1774        self.df.fillna(pd.NA, inplace=True)
1775        
1776    
1777    def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False):
1778        """
1779        Defines a globular protein pyMBE object  in `pymbe.df`.
1780
1781        Args:
1782            name (`str`): Unique label that identifies the protein.
1783            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1784            topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value}
1785            lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca".
1786            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1787
1788        Note:
1789            - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential.
1790        """
1791
1792        # Sanity checks
1793        valid_model_keys = ['1beadAA','2beadAA']
1794        valid_lj_setups = ["wca"]
1795        if model not in valid_model_keys:
1796            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
1797        if lj_setup_mode not in valid_lj_setups:
1798            raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}')
1799        if lj_setup_mode == "wca":
1800            sigma = 1*self.units.Quantity("reduced_length")
1801            epsilon = 1*self.units.Quantity("reduced_energy")
1802        part_dict={}
1803        sequence=[]
1804        metal_ions_charge_number_map=self.get_metal_ions_charge_number_map()
1805        for particle in topology_dict.keys():
1806            particle_name = re.split(r'\d+', particle)[0] 
1807            if particle_name not in part_dict.keys():
1808                if lj_setup_mode == "wca":
1809                    part_dict[particle_name]={"sigma": sigma,
1810                                        "offset": topology_dict[particle]['radius']*2-sigma,
1811                                        "epsilon": epsilon,
1812                                        "name": particle_name}
1813                if self.check_if_metal_ion(key=particle_name):
1814                    z=metal_ions_charge_number_map[particle_name]
1815                else:
1816                    z=0
1817                part_dict[particle_name]["z"]=z
1818            
1819            if self.check_aminoacid_key(key=particle_name):
1820                sequence.append(particle_name) 
1821            
1822        self.define_particles(parameters=part_dict,
1823                            overwrite=overwrite)
1824        residue_list = self.define_AA_residues(sequence=sequence, 
1825                                               model=model)
1826        index = len(self.df)
1827        self.df.at [index,'name'] = name
1828        self.df.at [index,'pmb_type'] = 'protein'
1829        self.df.at [index,'model'] = model
1830        self.df.at [index,('sequence','')] = sequence  
1831        self.df.at [index,('residue_list','')] = residue_list
1832        self.df.fillna(pd.NA, inplace=True)    
1833        return 
1834    
1835    def define_residue(self, name, central_bead, side_chains):
1836        """
1837        Defines a pyMBE object of type `residue` in `pymbe.df`.
1838
1839        Args:
1840            name(`str`): Unique label that identifies the `residue`.
1841            central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`.
1842            side_chains(`list` of `str`): List of `name`s of the pmb_objects to be placed as side_chains of the `residue`. Currently, only pmb_objects of type `particle`s or `residue`s are supported.
1843        """
1844        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'):
1845            return
1846        index = len(self.df)
1847        self.df.at [index, 'name'] = name
1848        self.df.at [index,'pmb_type'] = 'residue'
1849        self.df.at [index,'central_bead'] = central_bead
1850        self.df.at [index,('side_chains','')] = side_chains
1851        self.df.fillna(pd.NA, inplace=True)
1852        return 
1853
1854    def delete_entries_in_df(self, entry_name):
1855        """
1856        Deletes entries with name `entry_name` from the DataFrame if it exists.
1857
1858        Args:
1859            entry_name (`str`): The name of the entry in the dataframe to delete.
1860
1861        """
1862        if entry_name in self.df["name"].values:
1863            self.df = self.df[self.df["name"] != entry_name].reset_index(drop=True)
1864
1865    def destroy_pmb_object_in_system(self, name, espresso_system):
1866        """
1867        Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 
1868
1869        Args:
1870            name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`.
1871            espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library.
1872
1873        Note:
1874            - If `name`  is a object_type=`particle`, only the matching particles that are not part of bigger objects (e.g. `residue`, `molecule`) will be destroyed. To destroy particles in such objects, destroy the bigger object instead.
1875        """
1876        allowed_objects = ['particle','residue','molecule']
1877        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1878        if pmb_type not in allowed_objects:
1879            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1880        if pmb_type == 'particle':
1881            particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 
1882                                                  & (self.df['molecule_id'].isna())]
1883            particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna())
1884                                                 & (self.df['molecule_id'].isna())].particle_id.tolist()
1885            for particle_id in particle_ids_list:
1886                espresso_system.part.by_id(particle_id).remove()
1887            self.df = self.df.drop(index=particles_index)
1888        if pmb_type == 'residue':
1889            residues_id = self.df.loc[self.df['name']== name].residue_id.to_list()
1890            for residue_id in residues_id:
1891                molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0]
1892                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1893                self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index)
1894                for particle_id in particle_ids_list:
1895                    espresso_system.part.by_id(particle_id).remove()
1896                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)    
1897        if pmb_type == 'molecule':
1898            molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list()
1899            for molecule_id in molecules_id:
1900                molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']=="molecule")].name.values[0]
1901                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1902                self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index)
1903                for particle_id in particle_ids_list:
1904                    espresso_system.part.by_id(particle_id).remove()   
1905                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)             
1906        
1907        self.df.reset_index(drop=True,inplace=True)
1908
1909        return
1910
1911    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
1912        """
1913        Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration.
1914        To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an
1915        ion pair and the concentrations of the various species is solved numerically using a self-consistent approach.
1916        More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831).
1917        This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237).
1918
1919        Args:
1920            pH_res('float'): pH-value in the reservoir.
1921            c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
1922            activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
1923
1924        Returns:
1925            cH_res('pint.Quantity'): Concentration of H+ ions.
1926            cOH_res('pint.Quantity'): Concentration of OH- ions.
1927            cNa_res('pint.Quantity'): Concentration of Na+ ions.
1928            cCl_res('pint.Quantity'): Concentration of Cl- ions.
1929        """
1930
1931        self_consistent_run = 0
1932        cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+
1933
1934        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
1935            #Calculate and initial guess for the concentrations of various species based on ideal gas estimate
1936            cOH_res = self.Kw / cH_res 
1937            cNa_res = None
1938            cCl_res = None
1939            if cOH_res>=cH_res:
1940                #adjust the concentration of sodium if there is excess OH- in the reservoir:
1941                cNa_res = c_salt_res + (cOH_res-cH_res)
1942                cCl_res = c_salt_res
1943            else:
1944                # adjust the concentration of chloride if there is excess H+ in the reservoir
1945                cCl_res = c_salt_res + (cH_res-cOH_res)
1946                cNa_res = c_salt_res
1947                
1948            def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res):
1949                nonlocal max_number_sc_runs, self_consistent_run
1950                if self_consistent_run<max_number_sc_runs:
1951                    self_consistent_run+=1
1952                    ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1953                    cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
1954                    if cOH_res>=cH_res:
1955                        #adjust the concentration of sodium if there is excess OH- in the reservoir:
1956                        cNa_res = c_salt_res + (cOH_res-cH_res)
1957                        cCl_res = c_salt_res
1958                    else:
1959                        # adjust the concentration of chloride if there is excess H+ in the reservoir
1960                        cCl_res = c_salt_res + (cH_res-cOH_res)
1961                        cNa_res = c_salt_res
1962                    return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1963                else:
1964                    return cH_res, cOH_res, cNa_res, cCl_res
1965            return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1966
1967        cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1968        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1969        determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1970
1971        while abs(determined_pH-pH_res)>1e-6:
1972            if determined_pH > pH_res:
1973                cH_res *= 1.005
1974            else:
1975                cH_res /= 1.003
1976            cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1977            ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1978            determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1979            self_consistent_run=0
1980
1981        return cH_res, cOH_res, cNa_res, cCl_res
1982
1983    def enable_motion_of_rigid_object(self, name, espresso_system):
1984        '''
1985        Enables the motion of the rigid object `name` in the `espresso_system`.
1986
1987        Args:
1988            name(`str`): Label of the object.
1989            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1990
1991        Note:
1992            - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"].
1993        '''
1994        print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
1995        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1996        if pmb_type != 'protein':
1997            raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein')
1998        molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list()
1999        for molecule_id in molecule_ids_list:    
2000            particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list()
2001            center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system)
2002            rigid_object_center = espresso_system.part.add(pos=center_of_mass,
2003                                                           rotation=[True,True,True], 
2004                                                           type=self.propose_unused_type())
2005            
2006            rigid_object_center.mass = len(particle_ids_list)
2007            momI = 0
2008            for pid in particle_ids_list:
2009                momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
2010
2011            rigid_object_center.rinertia = np.ones(3) * momI
2012            
2013            for particle_id in particle_ids_list:
2014                pid = espresso_system.part.by_id(particle_id)
2015                pid.vs_auto_relate_to(rigid_object_center.id)
2016        return
2017
2018    def filter_df(self, pmb_type):
2019        """
2020        Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns.
2021        
2022        Args:
2023            pmb_type(`str`): pmb_object_type to filter in `pmb.df`.
2024
2025        Returns:
2026            pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`.
2027        """
2028        pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type]
2029        pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1)
2030        return pmb_type_df
2031
2032    def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False):
2033        """
2034        Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2035        
2036        Args:
2037            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2038            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2039            use_default_bond(`bool`, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to 'False'. 
2040
2041        Returns:
2042            bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` if a matching bond exists
2043
2044        Note:
2045            - If `use_default_bond`=`True`, it returns "default" if no key is found.
2046        """
2047        bond_keys = [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']
2048        bond_defined=False
2049        for bond_key in bond_keys:
2050            if bond_key in self.df["name"].values:
2051                bond_defined=True
2052                correct_key=bond_key
2053                break
2054        if bond_defined:
2055            return correct_key
2056        elif use_default_bond:
2057            return 'default'
2058        else:
2059            return 
2060
2061    def find_value_from_es_type(self, es_type, column_name):
2062        """
2063        Finds a value in `pmb.df` for a `column_name` and `es_type` pair.
2064
2065        Args:
2066            es_type(`int`): value of the espresso type
2067            column_name(`str`): name of the column in `pymbe.df`
2068
2069        Returns:
2070            Value in `pymbe.df` matching  `column_name` and `es_type`
2071        """
2072        idx = pd.IndexSlice
2073        for state in ['state_one', 'state_two']:            
2074            index = self.df.loc[self.df[(state, 'es_type')] == es_type].index
2075            if len(index) > 0:
2076                if column_name == 'label':
2077                    label = self.df.loc[idx[index[0]], idx[(state,column_name)]]
2078                    return label
2079                else: 
2080                    column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]]
2081                    return column_name_value
2082        return None
2083
2084    def format_node(self, node_list):
2085        return "[" + " ".join(map(str, node_list)) + "]"
2086
2087
2088    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2089        """
2090        Generates coordinates outside a sphere centered at `center`.
2091
2092        Args:
2093            center(`lst`): Coordinates of the center of the sphere.
2094            radius(`float`): Radius of the sphere.
2095            max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates.
2096            n_samples(`int`): Number of sample points.
2097
2098        Returns:
2099            coord_list(`lst`): Coordinates of the sample points.
2100        """
2101        if not radius > 0: 
2102            raise ValueError (f'The value of {radius} must be a positive value')
2103        if not radius < max_dist:
2104            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
2105        coord_list = []
2106        counter = 0
2107        while counter<n_samples:
2108            coord = self.generate_random_points_in_a_sphere(center=center, 
2109                                            radius=max_dist,
2110                                            n_samples=1)[0]
2111            if np.linalg.norm(coord-np.asarray(center))>=radius:
2112                coord_list.append (coord)
2113                counter += 1
2114        return coord_list
2115    
2116    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2117        """
2118        Uniformly samples points from a hypersphere. If on_surface is set to True, the points are
2119        uniformly sampled from the surface of the hypersphere.
2120        
2121        Args:
2122            center(`lst`): Array with the coordinates of the center of the spheres.
2123            radius(`float`): Radius of the sphere.
2124            n_samples(`int`): Number of sample points to generate inside the sphere.
2125            on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere.
2126
2127        Returns:
2128            samples(`list`): Coordinates of the sample points inside the hypersphere.
2129        """
2130        # initial values
2131        center=np.array(center)
2132        d = center.shape[0]
2133        # sample n_samples points in d dimensions from a standard normal distribution
2134        samples = self.rng.normal(size=(n_samples, d))
2135        # make the samples lie on the surface of the unit hypersphere
2136        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
2137        samples /= normalize_radii
2138        if not on_surface:
2139            # make the samples lie inside the hypersphere with the correct density
2140            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
2141            new_radii = np.power(uniform_points, 1/d)
2142            samples *= new_radii
2143        # scale the points to have the correct radius and center
2144        samples = samples * radius + center
2145        return samples 
2146
2147    def generate_trial_perpendicular_vector(self,vector,magnitude):
2148        """
2149        Generates an orthogonal vector to the input `vector`.
2150
2151        Args:
2152            vector(`lst`): arbitrary vector.
2153            magnitude(`float`): magnitude of the orthogonal vector.
2154            
2155        Returns:
2156            (`lst`): Orthogonal vector with the same magnitude as the input vector.
2157        """ 
2158        np_vec = np.array(vector) 
2159        np_vec /= np.linalg.norm(np_vec) 
2160        if np.all(np_vec == 0):
2161            raise ValueError('Zero vector')
2162        # Generate a random vector 
2163        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
2164                                                                center=[0,0,0],
2165                                                                n_samples=1, 
2166                                                                on_surface=True)[0]
2167        # Project the random vector onto the input vector and subtract the projection
2168        projection = np.dot(random_vector, np_vec) * np_vec
2169        perpendicular_vector = random_vector - projection
2170        # Normalize the perpendicular vector to have the same magnitude as the input vector
2171        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
2172        return perpendicular_vector*magnitude
2173
2174    def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2175        """
2176        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length.
2177        If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead.
2178        If no bond is found, it prints a message and it does not return anything. If `hard_check` is activated, the code stops if no bond is found.
2179
2180        Args:
2181            particle_name1(str): label of the type of the first particle type of the bonded particles.
2182            particle_name2(str): label of the type of the second particle type of the bonded particles.
2183            hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2184            use_default_bond(bool, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to False. 
2185
2186        Returns:
2187            l0(`pint.Quantity`): bond length
2188        
2189        Note:
2190            - If `use_default_bond`=True and no bond is defined between `particle_name1` and `particle_name2`, it returns the default bond defined in `pmb.df`.
2191            - If `hard_check`=`True` stops the code when no bond is found.
2192        """
2193        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2194                                    particle_name2=particle_name2, 
2195                                    use_default_bond=use_default_bond)
2196        if bond_key:
2197            return self.df[self.df['name']==bond_key].l0.values[0]
2198        else:
2199            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2200            if hard_check:
2201                sys.exit(1)
2202            else:
2203                return
2204
2205    def get_charge_number_map(self):
2206        '''
2207        Gets the charge number of each `espresso_type` in `pymbe.df`.
2208        
2209        Returns:
2210            charge_number_map(`dict`): {espresso_type: z}.
2211        '''
2212        if self.df.state_one['es_type'].isnull().values.any():         
2213            df_state_one = self.df.state_one.dropna()     
2214            df_state_two = self.df.state_two.dropna()  
2215        else:    
2216            df_state_one = self.df.state_one
2217            if self.df.state_two['es_type'].isnull().values.any():
2218                df_state_two = self.df.state_two.dropna()   
2219            else:
2220                df_state_two = self.df.state_two
2221        state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values)
2222        state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values)
2223        charge_number_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2224        return charge_number_map
2225
2226    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2227        """
2228        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2229        `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule.
2230
2231        Args:
2232            particle_name1 (str): label of the type of the first particle type
2233            particle_name2 (str): label of the type of the second particle type
2234            combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2235
2236        Returns:
2237            {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2238
2239        Note:
2240            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
2241            - 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.
2242        """
2243        supported_combining_rules=["Lorentz-Berthelot"]
2244        lj_parameters_keys=["sigma","epsilon","offset","cutoff"]
2245        if combining_rule not in supported_combining_rules:
2246            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2247        lj_parameters={}
2248        for key in lj_parameters_keys:
2249            lj_parameters[key]=[]
2250        # Search the LJ parameters of the type pair
2251        for name in [particle_name1,particle_name2]:
2252            for key in lj_parameters_keys:
2253                lj_parameters[key].append(self.df[self.df.name == name][key].values[0])
2254        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2255        if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]):
2256            return {}
2257        # Apply combining rule
2258        if combining_rule == 'Lorentz-Berthelot':
2259            lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2
2260            lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2
2261            lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2
2262            lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1])
2263        return lj_parameters
2264
2265    def get_metal_ions_charge_number_map(self):
2266        """
2267        Gets a map with the charge numbers of all the metal ions supported.
2268
2269        Returns:
2270            metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number}
2271
2272        """
2273        metal_charge_number_map = {"Ca": 2}
2274        return metal_charge_number_map
2275
2276    def get_particle_id_map(self, object_name):
2277        '''
2278        Gets all the ids associated with the object with name `object_name` in `pmb.df`
2279
2280        Args:
2281            object_name(`str`): name of the object
2282        
2283        Returns:
2284            id_map(`dict`): dict of the structure {"all": [all_ids_with_object_name], "residue_map": {res_id: [particle_ids_in_res_id]}, "molecule_map": {mol_id: [particle_ids_in_mol_id]}, }
2285        '''
2286        object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0]
2287        valid_types = ["particle", "molecule", "residue", "protein"]
2288        if object_type not in valid_types:
2289            raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}")
2290        id_list = []
2291        mol_map = {}
2292        res_map = {}
2293        def do_res_map(res_ids):
2294            for res_id in res_ids:
2295                res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2296                res_map[res_id]=res_list
2297            return res_map
2298        if object_type in ['molecule', 'protein']:
2299            mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist()
2300            for mol_id in mol_ids:
2301                res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist())
2302                res_map=do_res_map(res_ids=res_ids)    
2303                mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2304                id_list+=mol_list
2305                mol_map[mol_id]=mol_list
2306        elif object_type == 'residue':     
2307            res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist()
2308            res_map=do_res_map(res_ids=res_ids)
2309            id_list=[]
2310            for res_id_list in res_map.values():
2311                id_list+=res_id_list
2312        elif object_type == 'particle':
2313            id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist()
2314        return {"all": id_list, "molecule_map": mol_map, "residue_map": res_map}
2315
2316    def get_pka_set(self):
2317        '''
2318        Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df`
2319        
2320        Returns:
2321            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
2322        '''
2323        titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna()
2324        pka_set = {}
2325        for index in titratables_AA_df.name.keys():
2326            name = titratables_AA_df.name[index]
2327            pka_value = titratables_AA_df.pka[index]
2328            acidity = titratables_AA_df.acidity[index]   
2329            pka_set[name] = {'pka_value':pka_value,'acidity':acidity}
2330        return pka_set 
2331    
2332    def get_radius_map(self, dimensionless=True):
2333        '''
2334        Gets the effective radius of each `espresso type` in `pmb.df`. 
2335
2336        Args:
2337            dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False.
2338        
2339        Returns:
2340            radius_map(`dict`): {espresso_type: radius}.
2341
2342        Note:
2343            The radius corresponds to (sigma+offset)/2
2344        '''
2345        df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates()
2346        df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates()
2347        state_one = pd.Series((df_state_one.sigma.values+df_state_one.offset.values)/2.0,index=df_state_one.state_one.es_type.values)
2348        state_two = pd.Series((df_state_two.sigma.values+df_state_two.offset.values)/2.0,index=df_state_two.state_two.es_type.values)
2349        radius_map  = pd.concat([state_one,state_two],axis=0).to_dict()  
2350        if dimensionless:
2351            for key in radius_map:
2352                radius_map[key] = radius_map[key].magnitude
2353        return radius_map
2354
2355    def get_reduced_units(self):
2356        """
2357        Returns the  current set of reduced units defined in pyMBE.
2358
2359        Returns:
2360            reduced_units_text(`str`): text with information about the current set of reduced units.
2361
2362        """
2363        unit_length=self.units.Quantity(1,'reduced_length')
2364        unit_energy=self.units.Quantity(1,'reduced_energy')
2365        unit_charge=self.units.Quantity(1,'reduced_charge')
2366        reduced_units_text = "\n".join(["Current set of reduced units:",
2367                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2368                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2369                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2370                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2371                                        ])   
2372        return reduced_units_text
2373
2374    def get_resource(self, path):
2375        '''
2376        Locate a file resource of the pyMBE package.
2377
2378        Args:
2379            path(`str`): Relative path to the resource
2380
2381        Returns:
2382            path(`str`): Absolute path to the resource
2383
2384        '''
2385        import os
2386        return os.path.join(os.path.dirname(__file__), path)
2387
2388    def get_type_map(self):
2389        """
2390        Gets all different espresso types assigned to particles  in `pmb.df`.
2391        
2392        Returns:
2393            type_map(`dict`): {"name": espresso_type}.
2394        """
2395        df_state_one = self.df.state_one.dropna(how='all')     
2396        df_state_two = self.df.state_two.dropna(how='all')  
2397        state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label)
2398        state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label)
2399        type_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2400        return type_map
2401
2402    def initialize_lattice_builder(self, diamond_lattice):
2403        """
2404        Initialize the lattice builder with the DiamondLattice object.
2405
2406        Args:
2407            diamond_lattice(`DiamondLattice`): DiamondLattice object from the `lib/lattice` module to be used in the LatticeBuilder.
2408        """
2409        if not isinstance(diamond_lattice, DiamondLattice):
2410            raise TypeError("Currently only DiamondLattice objects are supported.")
2411        self.lattice_builder = LatticeBuilder(lattice=diamond_lattice)
2412        logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}")
2413        return self.lattice_builder
2414
2415
2416    def load_interaction_parameters(self, filename, overwrite=False):
2417        """
2418        Loads the interaction parameters stored in `filename` into `pmb.df`
2419        
2420        Args:
2421            filename(`str`): name of the file to be read
2422            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2423        """
2424        without_units = ['z','es_type']
2425        with_units = ['sigma','epsilon','offset','cutoff']
2426
2427        with open(filename, 'r') as f:
2428            interaction_data = json.load(f)
2429            interaction_parameter_set = interaction_data["data"]
2430
2431        for key in interaction_parameter_set:
2432            param_dict=interaction_parameter_set[key]
2433            object_type=param_dict.pop('object_type')
2434            if object_type == 'particle':
2435                not_required_attributes={}    
2436                for not_required_key in without_units+with_units:
2437                    if not_required_key in param_dict.keys():
2438                        if not_required_key in with_units:
2439                            not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key))
2440                        elif not_required_key in without_units:
2441                            not_required_attributes[not_required_key]=param_dict.pop(not_required_key)
2442                    else:
2443                        not_required_attributes[not_required_key]=pd.NA
2444                self.define_particle(name=param_dict.pop('name'),
2445                                z=not_required_attributes.pop('z'),
2446                                sigma=not_required_attributes.pop('sigma'),
2447                                offset=not_required_attributes.pop('offset'),
2448                                cutoff=not_required_attributes.pop('cutoff'),
2449                                epsilon=not_required_attributes.pop('epsilon'),
2450                                overwrite=overwrite)
2451            elif object_type == 'residue':
2452                self.define_residue(**param_dict)
2453            elif object_type == 'molecule':
2454                self.define_molecule(**param_dict)
2455            elif object_type == 'peptide':
2456                self.define_peptide(**param_dict)
2457            elif object_type == 'bond':
2458                particle_pairs = param_dict.pop('particle_pairs')
2459                bond_parameters = param_dict.pop('bond_parameters')
2460                bond_type = param_dict.pop('bond_type')
2461                if bond_type == 'harmonic':
2462                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2463                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2464                    bond = {'r_0'    : r_0,
2465                            'k'      : k,
2466                            }
2467
2468                elif bond_type == 'FENE':
2469                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2470                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2471                    d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max'))
2472                    bond =  {'r_0'    : r_0,
2473                             'k'      : k,
2474                             'd_r_max': d_r_max,
2475                             }
2476                else:
2477                    raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds")
2478                self.define_bond(bond_type=bond_type, 
2479                                bond_parameters=bond, 
2480                                particle_pairs=particle_pairs)
2481            else:
2482                raise ValueError(object_type+' is not a known pmb object type')
2483            
2484        return
2485    
2486    def load_pka_set(self, filename, overwrite=True):
2487        """
2488        Loads the pka_set stored in `filename` into `pmb.df`.
2489        
2490        Args:
2491            filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}.
2492            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 
2493        """
2494        with open(filename, 'r') as f:
2495            pka_data = json.load(f)
2496            pka_set = pka_data["data"]
2497
2498        self.check_pka_set(pka_set=pka_set)
2499
2500        for key in pka_set:
2501            acidity = pka_set[key]['acidity']
2502            pka_value = pka_set[key]['pka_value']
2503            self.set_particle_acidity(name=key, 
2504                                      acidity=acidity, 
2505                                      pka=pka_value, 
2506                                      overwrite=overwrite)
2507        return
2508
2509
2510    def propose_unused_type(self):
2511        """
2512        Searches in `pmb.df` all the different particle types defined and returns a new one.
2513
2514        Returns:
2515            unused_type(`int`): unused particle type
2516        """
2517        type_map = self.get_type_map()
2518        if not type_map:  
2519            unused_type = 0
2520        else:
2521            valid_values = [v for v in type_map.values() if pd.notna(v)]  # Filter out pd.NA values
2522            unused_type = max(valid_values) + 1 if valid_values else 0  # Ensure max() doesn't fail if all values are NA
2523        return unused_type
2524
2525    def protein_sequence_parser(self, sequence):
2526        '''
2527        Parses `sequence` to the one letter code for amino acids.
2528        
2529        Args:
2530            sequence(`str` or `lst`): Sequence of the amino acid. 
2531
2532        Returns:
2533            clean_sequence(`lst`): `sequence` using the one letter code.
2534        
2535        Note:
2536            - Accepted formats for `sequence` are:
2537                - `lst` with one letter or three letter code of each aminoacid in each element
2538                - `str` with the sequence using the one letter code
2539                - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-"
2540        
2541        '''
2542        # Aminoacid key
2543        keys={"ALA": "A",
2544                "ARG": "R",
2545                "ASN": "N",
2546                "ASP": "D",
2547                "CYS": "C",
2548                "GLU": "E",
2549                "GLN": "Q",
2550                "GLY": "G",
2551                "HIS": "H",
2552                "ILE": "I",
2553                "LEU": "L",
2554                "LYS": "K",
2555                "MET": "M",
2556                "PHE": "F",
2557                "PRO": "P",
2558                "SER": "S",
2559                "THR": "T",
2560                "TRP": "W",
2561                "TYR": "Y",
2562                "VAL": "V",
2563                "PSER": "J",
2564                "PTHR": "U",
2565                "PTyr": "Z",
2566                "NH2": "n",
2567                "COOH": "c"}
2568        clean_sequence=[]
2569        if isinstance(sequence, str):
2570            if sequence.find("-") != -1:
2571                splited_sequence=sequence.split("-")
2572                for residue in splited_sequence:
2573                    if len(residue) == 1:
2574                        if residue in keys.values():
2575                            residue_ok=residue
2576                        else:
2577                            if residue.upper() in keys.values():
2578                                residue_ok=residue.upper()
2579                            else:
2580                                raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence")
2581                        clean_sequence.append(residue_ok)
2582                    else:
2583                        if residue in keys.keys():
2584                            clean_sequence.append(keys[residue])
2585                        else:
2586                            if residue.upper() in keys.keys():
2587                                clean_sequence.append(keys[residue.upper()])
2588                            else:
2589                                raise ValueError("Unknown  code for a residue: ", residue, " please review the input sequence")
2590            else:
2591                for residue in sequence:
2592                    if residue in keys.values():
2593                        residue_ok=residue
2594                    else:
2595                        if residue.upper() in keys.values():
2596                            residue_ok=residue.upper()
2597                        else:
2598                            raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence")
2599                    clean_sequence.append(residue_ok)
2600        if isinstance(sequence, list):
2601            for residue in sequence:
2602                if residue in keys.values():
2603                    residue_ok=residue
2604                else:
2605                    if residue.upper() in keys.values():
2606                        residue_ok=residue.upper()
2607                    elif (residue.upper() in keys.keys()):
2608                        residue_ok= keys[residue.upper()]
2609                    else:
2610                        raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence")
2611                clean_sequence.append(residue_ok)
2612        return clean_sequence
2613    
2614    def read_pmb_df (self,filename):
2615        """
2616        Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE.
2617
2618        Args:
2619            filename(`str`): path to file.
2620
2621        Note:
2622            This function only accepts files with CSV format. 
2623        """
2624        
2625        if filename.rsplit(".", 1)[1] != "csv":
2626            raise ValueError("Only files with CSV format are supported")
2627        df = pd.read_csv (filename,header=[0, 1], index_col=0)
2628        columns_names = self.setup_df()
2629        
2630        multi_index = pd.MultiIndex.from_tuples(columns_names)
2631        df.columns = multi_index
2632        
2633        self.df = self.convert_columns_to_original_format(df)
2634        self.df.fillna(pd.NA, inplace=True)
2635        return self.df
2636    
2637    def read_protein_vtf_in_df (self,filename,unit_length=None):
2638        """
2639        Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`.
2640
2641        Args:
2642            filename(`str`): Path to the vtf file with the coarse-grained model.
2643            unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None.
2644
2645        Returns:
2646            topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}
2647
2648        Note:
2649            - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom.
2650        """
2651
2652        print (f'Loading protein coarse grain model file: {filename}')
2653
2654        coord_list = []
2655        particles_dict = {}
2656
2657        if unit_length is None:
2658            unit_length = 1 * self.units.angstrom 
2659
2660        with open (filename,'r') as protein_model:
2661            for line in protein_model :
2662                line_split = line.split()
2663                if line_split:
2664                    line_header = line_split[0]
2665                    if line_header == 'atom':
2666                        atom_id  = line_split[1]
2667                        atom_name = line_split[3]
2668                        atom_resname = line_split[5]
2669                        chain_id = line_split[9]
2670                        radius = float(line_split [11])*unit_length 
2671                        particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius]
2672                    elif line_header.isnumeric(): 
2673                        atom_coord = line_split[1:] 
2674                        atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord]
2675                        coord_list.append (atom_coord)
2676
2677        numbered_label = []
2678        i = 0   
2679        
2680        for atom_id in particles_dict.keys():
2681    
2682            if atom_id == 1:
2683                atom_name = particles_dict[atom_id][0]
2684                numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2685                numbered_label.append(numbered_name)
2686
2687            elif atom_id != 1: 
2688            
2689                if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]:
2690                    i += 1                    
2691                    count = 1
2692                    atom_name = particles_dict[atom_id][0]
2693                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2694                    numbered_label.append(numbered_name)
2695                    
2696                elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]:
2697                    if count == 2 or particles_dict[atom_id][1] == 'GLY':
2698                        i +=1  
2699                        count = 0
2700                    atom_name = particles_dict[atom_id][0]
2701                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2702                    numbered_label.append(numbered_name)
2703                    count +=1
2704
2705        topology_dict = {}
2706
2707        for i in range (0, len(numbered_label)):   
2708            topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] ,
2709                                                    'chain_id':numbered_label[i][1],
2710                                                    'radius':numbered_label[i][2] }
2711
2712        return topology_dict
2713
2714    def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2715        """
2716        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2717        If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead.
2718        If no bond is found, it prints a message and it does not return anything. If `hard_check` is activated, the code stops if no bond is found.
2719
2720        Args:
2721            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2722            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2723            hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2724            use_default_bond(`bool`, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to False. 
2725
2726        Returns:
2727            bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library.
2728        
2729        Note:
2730            - If `use_default_bond`=True and no bond is defined between `particle_name1` and `particle_name2`, it returns the default bond defined in `pmb.df`.
2731            - If `hard_check`=`True` stops the code when no bond is found.
2732        """
2733        
2734        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2735                                    particle_name2=particle_name2, 
2736                                    use_default_bond=use_default_bond)
2737        if use_default_bond:
2738            if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'):
2739                raise ValueError(f"use_default_bond is set to {use_default_bond} but no default bond has been defined. Please define a default bond with pmb.define_default_bond")
2740        if bond_key:
2741            return self.df[self.df['name']==bond_key].bond_object.values[0]
2742        else:
2743            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2744            if hard_check:
2745                sys.exit(1)
2746            else:
2747                return
2748
2749    def search_particles_in_residue(self, residue_name):
2750        '''
2751        Searches for all particles in a given residue of name `residue_name`.
2752
2753        Args:
2754            residue_name (`str`): name of the residue to be searched
2755
2756        Returns:
2757            list_of_particles_in_residue (`lst`): list of the names of all particles in the residue
2758
2759        Note:
2760            - The function returns a name per particle in residue, i.e. if there are multiple particles with the same type `list_of_particles_in_residue` will have repeated items.
2761 
2762        '''
2763        index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 
2764        central_bead = self.df.at [index_residue, ('central_bead', '')]
2765        list_of_side_chains = self.df.at [index_residue, ('side_chains', '')]
2766
2767        list_of_particles_in_residue = []
2768        list_of_particles_in_residue.append(central_bead)
2769
2770        for side_chain in list_of_side_chains: 
2771            object_type = self.df[self.df['name']==side_chain].pmb_type.values[0]
2772
2773            if object_type == "residue":
2774                list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain)
2775                list_of_particles_in_residue += list_of_particles_in_side_chain_residue
2776            elif object_type == "particle":
2777                list_of_particles_in_residue.append(side_chain)
2778
2779        return list_of_particles_in_residue
2780
2781         
2782
2783    def set_particle_acidity(self, name, acidity=pd.NA, default_charge_number=0, pka=pd.NA, overwrite=True):
2784        """
2785        Sets the particle acidity including the charges in each of its possible states. 
2786
2787        Args:
2788            name(`str`): Unique label that identifies the `particle`. 
2789            acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None.
2790            default_charge_number (`int`): Charge number of the particle. Defaults to 0.
2791            pka(`float`, optional):  If `particle` is an acid or a base, it defines its pka-value. Defaults to pandas.NA.
2792            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2793     
2794        Note:
2795            - For particles with  `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 
2796        deprotonated states, respectively. 
2797            - For particles without an acidity `acidity = pandas.NA`, only `state_one` is defined.
2798            - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 
2799        """
2800        acidity_valid_keys = ['inert','acidic', 'basic']
2801        if not pd.isna(acidity):
2802            if acidity not in acidity_valid_keys:
2803                raise ValueError(f"Acidity {acidity} provided for particle name  {name} is not supproted. Valid keys are: {acidity_valid_keys}")
2804            if acidity in ['acidic', 'basic'] and pd.isna(pka):
2805                raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.")   
2806            if acidity == "inert":
2807                acidity = pd.NA
2808                print("Deprecation warning: the keyword 'inert' for acidity has been replaced by setting acidity = pd.NA. For backwards compatibility, acidity has been set to pd.NA. Support for `acidity = 'inert'` may be deprecated in future releases of pyMBE")
2809
2810        self.define_particle_entry_in_df(name=name)
2811        
2812        for index in self.df[self.df['name']==name].index:       
2813            if pka is not pd.NA:
2814                self.add_value_to_df(key=('pka',''),
2815                                    index=index,
2816                                    new_value=pka, 
2817                                    overwrite=overwrite)
2818            
2819            self.add_value_to_df(key=('acidity',''),
2820                                 index=index,
2821                                 new_value=acidity, 
2822                                 overwrite=overwrite) 
2823            if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')):
2824                self.add_value_to_df(key=('state_one','es_type'),
2825                                     index=index,
2826                                     new_value=self.propose_unused_type(), 
2827                                     overwrite=overwrite)  
2828            if pd.isna(self.df.loc [self.df['name']  == name].acidity.iloc[0]):
2829                self.add_value_to_df(key=('state_one','z'),
2830                                     index=index,
2831                                     new_value=default_charge_number, 
2832                                     overwrite=overwrite)
2833                self.add_value_to_df(key=('state_one','label'),
2834                                     index=index,
2835                                     new_value=name, 
2836                                    overwrite=overwrite)
2837            else:
2838                protonated_label = f'{name}H'
2839                self.add_value_to_df(key=('state_one','label'),
2840                                     index=index,
2841                                     new_value=protonated_label, 
2842                                    overwrite=overwrite)
2843                self.add_value_to_df(key=('state_two','label'),
2844                                     index=index,
2845                                     new_value=name, 
2846                                    overwrite=overwrite)
2847                if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')):
2848                    self.add_value_to_df(key=('state_two','es_type'),
2849                                         index=index,
2850                                         new_value=self.propose_unused_type(), 
2851                                         overwrite=overwrite)
2852                if self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'acidic':        
2853                    self.add_value_to_df(key=('state_one','z'),
2854                                         index=index,new_value=0, 
2855                                         overwrite=overwrite)
2856                    self.add_value_to_df(key=('state_two','z'),
2857                                         index=index,
2858                                         new_value=-1, 
2859                                         overwrite=overwrite)
2860                elif self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'basic':
2861                    self.add_value_to_df(key=('state_one','z'),
2862                                         index=index,new_value=+1, 
2863                                         overwrite=overwrite)
2864                    self.add_value_to_df(key=('state_two','z'),
2865                                         index=index,
2866                                         new_value=0, 
2867                                         overwrite=overwrite)   
2868        self.df.fillna(pd.NA, inplace=True)
2869        return
2870    
2871    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True):
2872        """
2873        Sets the set of reduced units used by pyMBE.units and it prints it.
2874
2875        Args:
2876            unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 
2877            unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
2878            temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 
2879            Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
2880            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
2881
2882        Note:
2883            - If no `temperature` is given, a value of 298.15 K is assumed by default.
2884            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
2885            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
2886            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2887        """
2888        if unit_length is None:
2889            unit_length= 0.355*self.units.nm
2890        if temperature is None:
2891            temperature = 298.15 * self.units.K
2892        if unit_charge is None:
2893            unit_charge = scipy.constants.e * self.units.C
2894        if Kw is None:
2895            Kw = 1e-14
2896        # Sanity check
2897        variables=[unit_length,temperature,unit_charge]
2898        dimensionalities=["[length]","[temperature]","[charge]"]
2899        for variable,dimensionality in zip(variables,dimensionalities):
2900            self.check_dimensionality(variable,dimensionality)
2901        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2902        self.kT=temperature*self.kB
2903        self.units._build_cache()
2904        self.units.define(f'reduced_energy = {self.kT} ')
2905        self.units.define(f'reduced_length = {unit_length}')
2906        self.units.define(f'reduced_charge = {unit_charge}')
2907        if verbose:        
2908            print(self.get_reduced_units())
2909        return
2910
2911    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
2912        """
2913        Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 
2914
2915        Args:
2916            counter_ion(`str`): `name` of the counter_ion `particle`.
2917            constant_pH(`float`): pH-value.
2918            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
2919            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2920            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
2921
2922        Returns:
2923            RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2924            sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
2925        """
2926        from espressomd import reaction_methods
2927        if exclusion_range is None:
2928            exclusion_range = max(self.get_radius_map().values())*2.0
2929        if pka_set is None:
2930            pka_set=self.get_pka_set()    
2931        self.check_pka_set(pka_set=pka_set)
2932        if use_exclusion_radius_per_type:
2933            exclusion_radius_per_type = self.get_radius_map()
2934        else:
2935            exclusion_radius_per_type = {}
2936        
2937        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2938                                                    exclusion_range=exclusion_range, 
2939                                                    seed=self.seed, 
2940                                                    constant_pH=constant_pH,
2941                                                    exclusion_radius_per_type = exclusion_radius_per_type
2942                                                    )
2943        sucessfull_reactions_labels=[]
2944        charge_number_map = self.get_charge_number_map()
2945        for name in pka_set.keys():
2946            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
2947                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
2948                continue
2949            gamma=10**-pka_set[name]['pka_value']
2950            state_one_type   = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
2951            state_two_type   = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
2952            counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0]
2953            RE.add_reaction(gamma=gamma,
2954                            reactant_types=[state_one_type],
2955                            product_types=[state_two_type, counter_ion_type],
2956                            default_charges={state_one_type: charge_number_map[state_one_type],
2957                            state_two_type: charge_number_map[state_two_type],
2958                            counter_ion_type: charge_number_map[counter_ion_type]})
2959            sucessfull_reactions_labels.append(name)
2960        return RE, sucessfull_reactions_labels
2961
2962    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2963        """
2964        Sets up grand-canonical coupling to a reservoir of salt.
2965        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2966
2967        Args:
2968            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2969            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
2970            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
2971            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
2972            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
2973            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2974
2975        Returns:
2976            RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
2977        """
2978        from espressomd import reaction_methods
2979        if exclusion_range is None:
2980            exclusion_range = max(self.get_radius_map().values())*2.0
2981        if use_exclusion_radius_per_type:
2982            exclusion_radius_per_type = self.get_radius_map()
2983        else:
2984            exclusion_radius_per_type = {}
2985        
2986        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2987                                                    exclusion_range=exclusion_range, 
2988                                                    seed=self.seed, 
2989                                                    exclusion_radius_per_type = exclusion_radius_per_type
2990                                                    )
2991
2992        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
2993        determined_activity_coefficient = activity_coefficient(c_salt_res)
2994        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
2995
2996        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
2997        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
2998
2999        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
3000        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
3001
3002        if salt_cation_charge <= 0:
3003            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3004        if salt_anion_charge >= 0:
3005            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3006
3007        # Grand-canonical coupling to the reservoir
3008        RE.add_reaction(
3009            gamma = K_salt.magnitude,
3010            reactant_types = [],
3011            reactant_coefficients = [],
3012            product_types = [ salt_cation_es_type, salt_anion_es_type ],
3013            product_coefficients = [ 1, 1 ],
3014            default_charges = {
3015                salt_cation_es_type: salt_cation_charge, 
3016                salt_anion_es_type: salt_anion_charge, 
3017            }
3018        )
3019
3020        return RE
3021
3022    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, pka_set=None, use_exclusion_radius_per_type = False):
3023        """
3024        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3025        reservoir of small ions. 
3026        This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
3027
3028        [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3029
3030        Args:
3031            pH_res ('float): pH-value in the reservoir.
3032            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3033            proton_name ('str'): Name of the proton (H+) particle.
3034            hydroxide_name ('str'): Name of the hydroxide (OH-) particle.
3035            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
3036            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
3037            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3038            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
3039            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3040            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
3041
3042        Returns:
3043            RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3044            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3045            ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3046        """
3047        from espressomd import reaction_methods
3048        if exclusion_range is None:
3049            exclusion_range = max(self.get_radius_map().values())*2.0
3050        if pka_set is None:
3051            pka_set=self.get_pka_set()    
3052        self.check_pka_set(pka_set=pka_set)
3053        if use_exclusion_radius_per_type:
3054            exclusion_radius_per_type = self.get_radius_map()
3055        else:
3056            exclusion_radius_per_type = {}
3057        
3058        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3059                                                    exclusion_range=exclusion_range, 
3060                                                    seed=self.seed, 
3061                                                    exclusion_radius_per_type = exclusion_radius_per_type
3062                                                    )
3063
3064        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3065        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3066        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3067        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3068        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3069        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3070        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3071
3072        proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0]
3073        hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0]     
3074        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
3075        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
3076
3077        proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0]
3078        hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0]     
3079        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
3080        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
3081
3082        if proton_charge <= 0:
3083            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
3084        if salt_cation_charge <= 0:
3085            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3086        if hydroxide_charge >= 0:
3087            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
3088        if salt_anion_charge >= 0:
3089            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3090
3091        # Grand-canonical coupling to the reservoir
3092        # 0 = H+ + OH-
3093        RE.add_reaction(
3094            gamma = K_W.magnitude,
3095            reactant_types = [],
3096            reactant_coefficients = [],
3097            product_types = [ proton_es_type, hydroxide_es_type ],
3098            product_coefficients = [ 1, 1 ],
3099            default_charges = {
3100                proton_es_type: proton_charge, 
3101                hydroxide_es_type: hydroxide_charge, 
3102            }
3103        )
3104
3105        # 0 = Na+ + Cl-
3106        RE.add_reaction(
3107            gamma = K_NACL.magnitude,
3108            reactant_types = [],
3109            reactant_coefficients = [],
3110            product_types = [ salt_cation_es_type, salt_anion_es_type ],
3111            product_coefficients = [ 1, 1 ],
3112            default_charges = {
3113                salt_cation_es_type: salt_cation_charge, 
3114                salt_anion_es_type: salt_anion_charge, 
3115            }
3116        )
3117
3118        # 0 = Na+ + OH-
3119        RE.add_reaction(
3120            gamma = (K_NACL * K_W / K_HCL).magnitude,
3121            reactant_types = [],
3122            reactant_coefficients = [],
3123            product_types = [ salt_cation_es_type, hydroxide_es_type ],
3124            product_coefficients = [ 1, 1 ],
3125            default_charges = {
3126                salt_cation_es_type: salt_cation_charge, 
3127                hydroxide_es_type: hydroxide_charge, 
3128            }
3129        )
3130
3131        # 0 = H+ + Cl-
3132        RE.add_reaction(
3133            gamma = K_HCL.magnitude,
3134            reactant_types = [],
3135            reactant_coefficients = [],
3136            product_types = [ proton_es_type, salt_anion_es_type ],
3137            product_coefficients = [ 1, 1 ],
3138            default_charges = {
3139                proton_es_type: proton_charge, 
3140                salt_anion_es_type: salt_anion_charge, 
3141            }
3142        )
3143
3144        # Annealing moves to ensure sufficient sampling
3145        # Cation annealing H+ = Na+
3146        RE.add_reaction(
3147            gamma = (K_NACL / K_HCL).magnitude,
3148            reactant_types = [proton_es_type],
3149            reactant_coefficients = [ 1 ],
3150            product_types = [ salt_cation_es_type ],
3151            product_coefficients = [ 1 ],
3152            default_charges = {
3153                proton_es_type: proton_charge, 
3154                salt_cation_es_type: salt_cation_charge, 
3155            }
3156        )
3157
3158        # Anion annealing OH- = Cl- 
3159        RE.add_reaction(
3160            gamma = (K_HCL / K_W).magnitude,
3161            reactant_types = [hydroxide_es_type],
3162            reactant_coefficients = [ 1 ],
3163            product_types = [ salt_anion_es_type ],
3164            product_coefficients = [ 1 ],
3165            default_charges = {
3166                hydroxide_es_type: hydroxide_charge, 
3167                salt_anion_es_type: salt_anion_charge, 
3168            }
3169        )
3170
3171        sucessful_reactions_labels=[]
3172        charge_number_map = self.get_charge_number_map()
3173        for name in pka_set.keys():
3174            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3175                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3176                continue
3177
3178            Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3179
3180            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3181            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3182
3183            # Reaction in terms of proton: HA = A + H+
3184            RE.add_reaction(gamma=Ka.magnitude,
3185                            reactant_types=[state_one_type],
3186                            reactant_coefficients=[1],
3187                            product_types=[state_two_type, proton_es_type],
3188                            product_coefficients=[1, 1],
3189                            default_charges={state_one_type: charge_number_map[state_one_type],
3190                            state_two_type: charge_number_map[state_two_type],
3191                            proton_es_type: proton_charge})
3192
3193            # Reaction in terms of salt cation: HA = A + Na+
3194            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
3195                            reactant_types=[state_one_type],
3196                            reactant_coefficients=[1],
3197                            product_types=[state_two_type, salt_cation_es_type],
3198                            product_coefficients=[1, 1],
3199                            default_charges={state_one_type: charge_number_map[state_one_type],
3200                            state_two_type: charge_number_map[state_two_type],
3201                            salt_cation_es_type: salt_cation_charge})
3202
3203            # Reaction in terms of hydroxide: OH- + HA = A
3204            RE.add_reaction(gamma=(Ka / K_W).magnitude,
3205                            reactant_types=[state_one_type, hydroxide_es_type],
3206                            reactant_coefficients=[1, 1],
3207                            product_types=[state_two_type],
3208                            product_coefficients=[1],
3209                            default_charges={state_one_type: charge_number_map[state_one_type],
3210                            state_two_type: charge_number_map[state_two_type],
3211                            hydroxide_es_type: hydroxide_charge})
3212
3213            # Reaction in terms of salt anion: Cl- + HA = A
3214            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
3215                            reactant_types=[state_one_type, salt_anion_es_type],
3216                            reactant_coefficients=[1, 1],
3217                            product_types=[state_two_type],
3218                            product_coefficients=[1],
3219                            default_charges={state_one_type: charge_number_map[state_one_type],
3220                            state_two_type: charge_number_map[state_two_type],
3221                            salt_anion_es_type: salt_anion_charge})
3222
3223            sucessful_reactions_labels.append(name)
3224        return RE, sucessful_reactions_labels, ionic_strength_res
3225
3226    def setup_grxmc_unified(self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
3227        """
3228        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3229        reservoir of small ions. 
3230        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-}. 
3231        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'.
3232
3233        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3234        [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3235
3236        Args:
3237            pH_res ('float'): pH-value in the reservoir.
3238            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3239            cation_name ('str'): Name of the cationic particle.
3240            anion_name ('str'): Name of the anionic particle.
3241            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3242            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
3243            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3244            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`.
3245
3246        Returns:
3247            RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3248            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3249            ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3250        """
3251        from espressomd import reaction_methods
3252        if exclusion_range is None:
3253            exclusion_range = max(self.get_radius_map().values())*2.0
3254        if pka_set is None:
3255            pka_set=self.get_pka_set()    
3256        self.check_pka_set(pka_set=pka_set)
3257        if use_exclusion_radius_per_type:
3258            exclusion_radius_per_type = self.get_radius_map()
3259        else:
3260            exclusion_radius_per_type = {}
3261        
3262        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3263                                                    exclusion_range=exclusion_range, 
3264                                                    seed=self.seed, 
3265                                                    exclusion_radius_per_type = exclusion_radius_per_type
3266                                                    )
3267
3268        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3269        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3270        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3271        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3272        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3273        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3274        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3275        K_XX = a_cation * a_anion
3276
3277        cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0]
3278        anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0]     
3279        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
3280        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
3281        if cation_charge <= 0:
3282            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3283        if anion_charge >= 0:
3284            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3285
3286        # Coupling to the reservoir: 0 = X+ + X-
3287        RE.add_reaction(
3288            gamma = K_XX.magnitude,
3289            reactant_types = [],
3290            reactant_coefficients = [],
3291            product_types = [ cation_es_type, anion_es_type ],
3292            product_coefficients = [ 1, 1 ],
3293            default_charges = {
3294                cation_es_type: cation_charge, 
3295                anion_es_type: anion_charge, 
3296            }
3297        )
3298
3299        sucessful_reactions_labels=[]
3300        charge_number_map = self.get_charge_number_map()
3301        for name in pka_set.keys():
3302            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3303                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3304                continue
3305
3306            Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l
3307            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3308
3309            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3310            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3311
3312            # Reaction in terms of small cation: HA = A + X+
3313            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3314                            reactant_types=[state_one_type],
3315                            reactant_coefficients=[1],
3316                            product_types=[state_two_type, cation_es_type],
3317                            product_coefficients=[1, 1],
3318                            default_charges={state_one_type: charge_number_map[state_one_type],
3319                            state_two_type: charge_number_map[state_two_type],
3320                            cation_es_type: cation_charge})
3321
3322            # Reaction in terms of small anion: X- + HA = A
3323            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3324                            reactant_types=[state_one_type, anion_es_type],
3325                            reactant_coefficients=[1, 1],
3326                            product_types=[state_two_type],
3327                            product_coefficients=[1],
3328                            default_charges={state_one_type: charge_number_map[state_one_type],
3329                            state_two_type: charge_number_map[state_two_type],
3330                            anion_es_type: anion_charge})
3331
3332            sucessful_reactions_labels.append(name)
3333        return RE, sucessful_reactions_labels, ionic_strength_res
3334
3335    def setup_df (self):
3336        """
3337        Sets up the pyMBE's dataframe `pymbe.df`.
3338
3339        Returns:
3340            columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe
3341        """
3342        
3343        columns_dtypes = {
3344            'name': {
3345                '': str},
3346            'pmb_type': {
3347                '': str},
3348            'particle_id': {
3349                '': pd.Int64Dtype()},
3350            'particle_id2':  {
3351                '': pd.Int64Dtype()},
3352            'residue_id':  {
3353                '': pd.Int64Dtype()},
3354            'molecule_id':  {
3355                '': pd.Int64Dtype()},
3356            'acidity':  {
3357                '': str},
3358            'pka':  {
3359                '': object},
3360            'central_bead':  {
3361                '': object},
3362            'side_chains': {
3363                '': object},
3364            'residue_list': {
3365                '': object},
3366            'model': {
3367                '': str},
3368            'sigma': {
3369                '': object},
3370            'cutoff': {
3371                '': object},
3372            'offset': {
3373                '': object},
3374            'epsilon': {
3375                '': object},
3376            'state_one': {
3377                'label': str,
3378                'es_type': pd.Int64Dtype(),
3379                'z': pd.Int64Dtype()},
3380            'state_two': {
3381                'label': str,
3382                'es_type': pd.Int64Dtype(),
3383                'z': pd.Int64Dtype()},
3384            'sequence': {
3385                '': object},
3386            'bond_object': {
3387                '': object},
3388            'parameters_of_the_potential':{
3389                '': object},
3390            'l0': {
3391                '': float},
3392            'node_map':{
3393                '':object},
3394            'chain_map':{
3395                '':object}}
3396        
3397        self.df = pd.DataFrame(columns=pd.MultiIndex.from_tuples([(col_main, col_sub) for col_main, sub_cols in columns_dtypes.items() for col_sub in sub_cols.keys()]))
3398        
3399        for level1, sub_dtypes in columns_dtypes.items():
3400            for level2, dtype in sub_dtypes.items():
3401                self.df[level1, level2] = self.df[level1, level2].astype(dtype)
3402                
3403        columns_names = pd.MultiIndex.from_frame(self.df)
3404        columns_names = columns_names.names
3405                
3406        return columns_names
3407
3408    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True):
3409        """
3410        Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`.
3411
3412        Args:
3413            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
3414            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.
3415            combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3416            warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True.
3417
3418        Note:
3419            - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 
3420            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
3421            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3422
3423        """
3424        from itertools import combinations_with_replacement
3425        import warnings
3426        implemented_combining_rules = ['Lorentz-Berthelot']
3427        compulsory_parameters_in_df = ['sigma','epsilon']
3428        # Sanity check
3429        if combining_rule not in implemented_combining_rules:
3430            raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules)
3431        shift=0
3432        if shift_potential:
3433            shift="auto"
3434        # List which particles have sigma and epsilon values defined in pmb.df and which ones don't
3435        particles_types_with_LJ_parameters = []
3436        non_parametrized_labels= []
3437        for particle_type in self.get_type_map().values():
3438            check_list=[]
3439            for key in compulsory_parameters_in_df:
3440                value_in_df=self.find_value_from_es_type(es_type=particle_type,
3441                                                        column_name=key)
3442                check_list.append(pd.isna(value_in_df))
3443            if any(check_list):
3444                non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 
3445                                                                            column_name='label'))
3446            else:
3447                particles_types_with_LJ_parameters.append(particle_type)
3448        # Set up LJ interactions between all particle types
3449        for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 
3450            particle_name1 = self.find_value_from_es_type(es_type=type_pair[0],
3451                                                        column_name="name")
3452            particle_name2 = self.find_value_from_es_type(es_type=type_pair[1],
3453                                                        column_name="name")
3454            lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1,
3455                                                 particle_name2 = particle_name2,
3456                                                 combining_rule = combining_rule)
3457            
3458            # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
3459            if not lj_parameters:
3460                continue
3461            espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 
3462                                                                                    sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 
3463                                                                                    cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude,
3464                                                                                    offset = lj_parameters["offset"].to("reduced_length").magnitude, 
3465                                                                                    shift = shift)                                                                                          
3466            index = len(self.df)
3467            label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label")
3468            label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label")
3469            self.df.at [index, 'name'] = f'LJ: {label1}-{label2}'
3470            lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params()
3471
3472            self.add_value_to_df(index=index,
3473                                key=('pmb_type',''),
3474                                new_value='LennardJones')
3475            
3476            self.add_value_to_df(index=index,
3477                                key=('parameters_of_the_potential',''),
3478                                new_value=lj_params,
3479                                non_standard_value=True)
3480        if non_parametrized_labels and warnings:
3481            warnings.warn(f'The following particles do not have a defined value of sigma or epsilon in pmb.df: {non_parametrized_labels}. No LJ interaction has been added in ESPResSo for those particles.',UserWarning) 
3482        return
3483
3484    def write_pmb_df (self, filename):
3485        '''
3486        Writes the pyMBE dataframe into a csv file
3487        
3488        Args:
3489            filename(`str`): Path to the csv file 
3490        '''
3491
3492        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence']
3493        df = self.df.copy(deep=True)
3494        for column_name in columns_with_list_or_dict:
3495            df[column_name] = df[column_name].apply(lambda x: json.dumps(x) if isinstance(x, (np.ndarray, tuple, list, dict)) or pd.notnull(x) else x)
3496        df['bond_object'] = df['bond_object'].apply(lambda x: f'{x.__class__.__name__}({json.dumps({**x.get_params(), "bond_id": x._bond_id})})' if pd.notnull(x) else x)
3497        df.fillna(pd.NA, inplace=True)
3498        df.to_csv(filename)
3499        return
class pymbe_library:
  31class pymbe_library():
  32    """
  33    The library for the Molecular Builder for ESPResSo (pyMBE)
  34
  35    Attributes:
  36        N_A(`pint.Quantity`): Avogadro number.
  37        Kb(`pint.Quantity`): Boltzmann constant.
  38        e(`pint.Quantity`): Elementary charge.
  39        df(`Pandas.Dataframe`): Dataframe used to bookkeep all the information stored in pyMBE. Typically refered as `pmb.df`. 
  40        kT(`pint.Quantity`): Thermal energy.
  41        Kw(`pint.Quantity`): Ionic product of water. Used in the setup of the G-RxMC method.
  42    """
  43
  44    class NumpyEncoder(json.JSONEncoder):
  45        """
  46        Custom JSON encoder that converts NumPy arrays to Python lists
  47        and NumPy scalars to Python scalars.
  48        """
  49        def default(self, obj):
  50            if isinstance(obj, np.ndarray):
  51                return obj.tolist()
  52            if isinstance(obj, np.generic):
  53                return obj.item()
  54            return super().default(obj)
  55
  56    def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None):
  57        """
  58        Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 
  59        and sets up  the `pmb.df` for bookkeeping.
  60
  61        Args:
  62            temperature(`pint.Quantity`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None.
  63            unit_length(`pint.Quantity`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None.
  64            unit_charge (`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
  65            Kw (`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
  66        
  67        Note:
  68            - If no `temperature` is given, a value of 298.15 K is assumed by default.
  69            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
  70            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
  71            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
  72        """
  73        # Seed and RNG
  74        self.seed=seed
  75        self.rng = np.random.default_rng(seed)
  76        self.units=pint.UnitRegistry()
  77        self.N_A=scipy.constants.N_A / self.units.mol
  78        self.kB=scipy.constants.k * self.units.J / self.units.K
  79        self.e=scipy.constants.e * self.units.C
  80        self.set_reduced_units(unit_length=unit_length, 
  81                               unit_charge=unit_charge,
  82                               temperature=temperature, 
  83                               Kw=Kw, 
  84                               verbose=False)
  85        self.setup_df()
  86        self.lattice_builder = None
  87        return
  88
  89    def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False):
  90        """
  91        Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles.
  92
  93        Args:
  94            particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles
  95            particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles
  96            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False.
  97
  98        Returns:
  99            index(`int`): Row index where the bond information has been added in pmb.df.
 100        """
 101        particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0]
 102        particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0]
 103        
 104        bond_key = self.find_bond_key(particle_name1=particle_name1,
 105                                    particle_name2=particle_name2, 
 106                                    use_default_bond=use_default_bond)
 107        if not bond_key:
 108            return None
 109        self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1)
 110        indexs = np.where(self.df['name']==bond_key)
 111        index_list = list (indexs[0])
 112        used_bond_df = self.df.loc[self.df['particle_id2'].notnull()]
 113        #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict
 114        used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 )
 115        used_bond_index = used_bond_df.index.to_list()
 116        if not index_list:
 117            return None
 118        for index in index_list:
 119            if index not in used_bond_index:
 120                self.clean_df_row(index=int(index))
 121                self.df.at[index,'particle_id'] = particle_id1
 122                self.df.at[index,'particle_id2'] = particle_id2
 123                break
 124        return index
 125
 126    def add_bonds_to_espresso(self, espresso_system) :
 127        """
 128        Adds all bonds defined in `pmb.df` to `espresso_system`.
 129
 130        Args:
 131            espresso_system(`espressomd.system.System`): system object of espressomd library
 132        """
 133
 134        if 'bond' in self.df["pmb_type"].values:
 135            bond_df = self.df.loc[self.df ['pmb_type'] == 'bond']
 136            bond_list = bond_df.bond_object.values.tolist()
 137            for bond in bond_list:
 138                espresso_system.bonded_inter.add(bond)
 139
 140        else:
 141            print ('WARNING: There are no bonds defined in pymbe.df')
 142        
 143        return
 144
 145    def add_value_to_df(self,index,key,new_value, non_standard_value=False, overwrite=False):
 146        """
 147        Adds a value to a cell in the `pmb.df` DataFrame.
 148
 149        Args:
 150            index(`int`): index of the row to add the value to.
 151            key(`str`): the column label to add the value to.
 152            non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False.
 153            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
 154        """
 155
 156        token = "#protected:"
 157
 158        def protect(obj):
 159            if non_standard_value:
 160                return token + json.dumps(obj, cls=self.NumpyEncoder)
 161            return obj
 162
 163        def deprotect(obj):
 164            if non_standard_value and isinstance(obj, str) and obj.startswith(token):
 165                return json.loads(obj.removeprefix(token))
 166            return obj
 167
 168        # Make sure index is a scalar integer value
 169        index = int(index)
 170        assert isinstance(index, int), '`index` should be a scalar integer value.'
 171        idx = pd.IndexSlice
 172        if self.check_if_df_cell_has_a_value(index=index,key=key):
 173            old_value = self.df.loc[index,idx[key]]
 174            if not pd.Series([protect(old_value)]).equals(pd.Series([protect(new_value)])):
 175                name=self.df.loc[index,('name','')]
 176                pmb_type=self.df.loc[index,('pmb_type','')]
 177                logging.debug(f"You are attempting to redefine the properties of {name} of pmb_type {pmb_type}")    
 178                if overwrite:
 179                    logging.info(f'Overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}')
 180                if not overwrite:
 181                    logging.debug(f"pyMBE has preserved of the entry `{key}`: old_value = {old_value}. If you want to overwrite it with new_value = {new_value}, activate the switch overwrite = True ")
 182                    return
 183
 184        self.df.loc[index,idx[key]] = protect(new_value)
 185        if non_standard_value:
 186            self.df[key] = self.df[key].apply(deprotect)
 187        return
 188    
 189    def assign_molecule_id(self, molecule_index):
 190        """
 191        Assigns the `molecule_id` of the pmb object given by `pmb_type`
 192        
 193        Args:
 194            molecule_index(`int`): index of the current `pmb_object_type` to assign the `molecule_id`
 195        Returns:
 196            molecule_id(`int`): Id of the molecule
 197        """
 198
 199        self.clean_df_row(index=int(molecule_index))
 200        
 201        if self.df['molecule_id'].isnull().values.all():
 202            molecule_id = 0        
 203        else:
 204            molecule_id = self.df['molecule_id'].max() +1
 205
 206        self.add_value_to_df (key=('molecule_id',''),
 207                                index=int(molecule_index),
 208                                new_value=molecule_id)
 209
 210        return molecule_id
 211    
 212    def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system):
 213        """
 214        Calculates the center of the molecule with a given molecule_id.
 215
 216        Args:
 217            molecule_id(`int`): Id of the molecule whose center of mass is to be calculated.
 218            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 219        
 220        Returns:
 221            center_of_mass(`lst`): Coordinates of the center of mass.
 222        """
 223        center_of_mass = np.zeros(3)
 224        axis_list = [0,1,2]
 225        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
 226        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
 227        for pid in particle_id_list:
 228            for axis in axis_list:
 229                center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis]
 230        center_of_mass = center_of_mass / len(particle_id_list)
 231        return center_of_mass
 232
 233    def calculate_HH(self, molecule_name, pH_list=None, pka_set=None):
 234        """
 235        Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 
 236        for molecules with the name `molecule_name`.
 237
 238        Args:
 239            molecule_name(`str`): name of the molecule to calculate the ideal charge for
 240            pH_list(`lst`): pH-values to calculate. 
 241            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
 242
 243        Returns:
 244            Z_HH(`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list`
 245
 246        Note:
 247            - This function supports objects with pmb types: "molecule", "peptide" and "protein".
 248            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
 249            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
 250            - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan`  should be used.
 251        """
 252        if pH_list is None:
 253            pH_list=np.linspace(2,12,50)
 254        if pka_set is None:
 255            pka_set=self.get_pka_set() 
 256        self.check_pka_set(pka_set=pka_set)
 257        charge_number_map = self.get_charge_number_map()
 258        Z_HH=[]
 259        for pH_value in pH_list:    
 260            Z=0
 261            index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 
 262            residue_list = self.df.at [index,('residue_list','')]
 263            sequence = self.df.at [index,('sequence','')]
 264            if np.any(pd.isnull(sequence)):
 265                # Molecule has no sequence
 266                for residue in residue_list:
 267                    list_of_particles_in_residue = self.search_particles_in_residue(residue)
 268                    for particle in list_of_particles_in_residue:
 269                        if particle in pka_set.keys():
 270                            if pka_set[particle]['acidity'] == 'acidic':
 271                                psi=-1
 272                            elif pka_set[particle]['acidity']== 'basic':
 273                                psi=+1
 274                            else:
 275                                psi=0
 276                            Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value'])))                      
 277                Z_HH.append(Z)
 278            else:
 279                # Molecule has a sequence
 280                for name in sequence:
 281                    if name in pka_set.keys():
 282                        if pka_set[name]['acidity'] == 'acidic':
 283                            psi=-1
 284                        elif pka_set[name]['acidity']== 'basic':
 285                            psi=+1
 286                        else:
 287                            psi=0
 288                        Z+=psi/(1+10**(psi*(pH_value-pka_set[name]['pka_value'])))
 289                    else:
 290                        state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
 291                        Z+=charge_number_map[state_one_type]
 292                Z_HH.append(Z)
 293
 294        return Z_HH
 295
 296    def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
 297        """
 298        Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning.
 299
 300        Args:
 301            c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 302            c_salt('float'): Salt concentration in the reservoir.
 303            pH_list('lst'): List of pH-values in the reservoir. 
 304            pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}.
 305
 306        Returns:
 307            {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
 308            pH_system_list ('lst'): List of pH_values in the system.
 309            partition_coefficients_list ('lst'): List of partition coefficients of cations.
 310
 311        Note:
 312            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
 313            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
 314            - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result.
 315        """
 316        if pH_list is None:
 317            pH_list=np.linspace(2,12,50)
 318        if pka_set is None:
 319            pka_set=self.get_pka_set() 
 320        self.check_pka_set(pka_set=pka_set)
 321
 322        partition_coefficients_list = []
 323        pH_system_list = []
 324        Z_HH_Donnan={}
 325        for key in c_macro:
 326            Z_HH_Donnan[key] = []
 327
 328        def calc_charges(c_macro, pH):
 329            """
 330            Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation.
 331
 332            Args:
 333                c_macro('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 334                pH('float'): pH-value that is used in the HH equation.
 335
 336            Returns:
 337                charge('dict'): {"molecule_name": charge}
 338            """
 339            charge = {}
 340            for name in c_macro:
 341                charge[name] = self.calculate_HH(name, [pH], pka_set)[0]
 342            return charge
 343
 344        def calc_partition_coefficient(charge, c_macro):
 345            """
 346            Calculates the partition coefficients of positive ions according to the ideal Donnan theory.
 347
 348            Args:
 349                charge('dict'): {"molecule_name": charge}
 350                c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
 351            """
 352            nonlocal ionic_strength_res
 353            charge_density = 0.0
 354            for key in charge:
 355                charge_density += charge[key] * c_macro[key]
 356            return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude
 357
 358        for pH_value in pH_list:    
 359            # calculate the ionic strength of the reservoir
 360            if pH_value <= 7.0:
 361                ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 
 362            elif pH_value > 7.0:
 363                ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt
 364
 365            #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations
 366            #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation.
 367            #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects.
 368            equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro))
 369            logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq")
 370            partition_coefficient = 10**logxi.root
 371
 372            charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient))
 373            for key in c_macro:
 374                Z_HH_Donnan[key].append(charges_temp[key])
 375
 376            pH_system_list.append(pH_value - np.log10(partition_coefficient))
 377            partition_coefficients_list.append(partition_coefficient)
 378
 379        return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
 380
 381    def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset):
 382        """
 383        Calculates the initial bond length that is used when setting up molecules,
 384        based on the minimum of the sum of bonded and short-range (LJ) interactions.
 385        
 386        Args:
 387            bond_object(`espressomd.interactions.BondedInteractions`): instance of a bond object from espressomd library
 388            bond_type(`str`): label identifying the used bonded potential
 389            epsilon(`pint.Quantity`): LJ epsilon of the interaction between the particles
 390            sigma(`pint.Quantity`): LJ sigma of the interaction between the particles
 391            cutoff(`pint.Quantity`): cutoff-radius of the LJ interaction 
 392            offset(`pint.Quantity`): offset of the LJ interaction
 393        """    
 394        def truncated_lj_potential(x, epsilon, sigma, cutoff,offset):
 395            if x>cutoff:
 396                return 0.0
 397            else:
 398                return 4*epsilon*((sigma/(x-offset))**12-(sigma/(x-offset))**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6)
 399
 400        epsilon_red=epsilon.to('reduced_energy').magnitude
 401        sigma_red=sigma.to('reduced_length').magnitude
 402        cutoff_red=cutoff.to('reduced_length').magnitude
 403        offset_red=offset.to('reduced_length').magnitude
 404
 405        if bond_type == "harmonic":
 406            r_0 = bond_object.params.get('r_0')
 407            k = bond_object.params.get('k')
 408            l0 = scipy.optimize.minimize(lambda x: 0.5*k*(x-r_0)**2 + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red, offset_red), x0=r_0).x
 409        elif bond_type == "FENE":
 410            r_0 = bond_object.params.get('r_0')
 411            k = bond_object.params.get('k')
 412            d_r_max = bond_object.params.get('d_r_max')
 413            l0 = scipy.optimize.minimize(lambda x: -0.5*k*(d_r_max**2)*np.log(1-((x-r_0)/d_r_max)**2) + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red,offset_red), x0=1.0).x
 414        return l0
 415
 416    def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False):
 417        '''
 418        Calculates the net charge per molecule of molecules with `name` = molecule_name. 
 419        Returns the net charge per molecule and a maps with the net charge per residue and molecule.
 420
 421        Args:
 422            espresso_system(`espressomd.system.System`): system information 
 423            molecule_name(`str`): name of the molecule to calculate the net charge
 424            dimensionless(`bool'): sets if the charge is returned with a dimension or not
 425
 426        Returns:
 427            (`dict`): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }}
 428
 429        Note:
 430            - The net charge of the molecule is averaged over all molecules of type `name` 
 431            - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name`
 432        '''        
 433        valid_pmb_types = ["molecule", "protein"]
 434        pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0]
 435        if pmb_type not in valid_pmb_types:
 436            raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}")      
 437
 438        id_map = self.get_particle_id_map(object_name=molecule_name)
 439        def create_charge_map(espresso_system,id_map,label):
 440            charge_number_map={}
 441            for super_id in id_map[label].keys():
 442                if dimensionless:
 443                    net_charge=0
 444                else:
 445                    net_charge=0 * self.units.Quantity(1,'reduced_charge')
 446                for pid in id_map[label][super_id]:
 447                    if dimensionless:
 448                        net_charge+=espresso_system.part.by_id(pid).q
 449                    else:
 450                        net_charge+=espresso_system.part.by_id(pid).q * self.units.Quantity(1,'reduced_charge')
 451                charge_number_map[super_id]=net_charge
 452            return charge_number_map
 453        net_charge_molecules=create_charge_map(label="molecule_map",
 454                                                espresso_system=espresso_system,
 455                                                id_map=id_map)
 456        net_charge_residues=create_charge_map(label="residue_map",
 457                                                espresso_system=espresso_system,
 458                                                id_map=id_map)       
 459        if dimensionless:
 460            mean_charge=np.mean(np.array(list(net_charge_molecules.values())))
 461        else:
 462            mean_charge=np.mean(np.array([value.magnitude for value in net_charge_molecules.values()]))*self.units.Quantity(1,'reduced_charge')
 463        return {"mean": mean_charge, "molecules": net_charge_molecules, "residues": net_charge_residues}
 464
 465    def center_molecule_in_simulation_box(self, molecule_id, espresso_system):
 466        """
 467        Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`.
 468        
 469        Args:
 470            molecule_id(`int`): Id of the molecule to be centered.
 471            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 472        """
 473        if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0:
 474            raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.")      
 475        center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system)
 476        box_center = [espresso_system.box_l[0]/2.0,
 477                      espresso_system.box_l[1]/2.0,
 478                      espresso_system.box_l[2]/2.0]
 479        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
 480        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
 481        for pid in particle_id_list:
 482            es_pos = espresso_system.part.by_id(pid).pos
 483            espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center
 484        return 
 485
 486    def check_aminoacid_key(self, key):
 487        """
 488        Checks if `key` corresponds to a valid aminoacid letter code.
 489
 490        Args:
 491            key(`str`): key to be checked.
 492
 493        Returns:
 494            `bool`: True if `key` is a valid aminoacid letter code, False otherwise.
 495        """
 496        valid_AA_keys=['V', #'VAL'
 497                       'I', #'ILE'
 498                       'L', #'LEU'
 499                       'E', #'GLU'
 500                       'Q', #'GLN'
 501                       'D', #'ASP'
 502                       'N', #'ASN'
 503                       'H', #'HIS'
 504                       'W', #'TRP'
 505                       'F', #'PHE'
 506                       'Y', #'TYR'
 507                       'R', #'ARG' 
 508                       'K', #'LYS'
 509                       'S', #'SER'
 510                       'T', #'THR'
 511                       'M', #'MET'
 512                       'A', #'ALA'
 513                       'G', #'GLY'
 514                       'P', #'PRO'
 515                       'C'] #'CYS'
 516        if key in valid_AA_keys:
 517            return True
 518        else:
 519            return False
 520
 521    def check_dimensionality(self, variable, expected_dimensionality):
 522        """
 523        Checks if the dimensionality of `variable` matches `expected_dimensionality`.
 524
 525        Args:
 526            variable(`pint.Quantity`): Quantity to be checked.
 527            expected_dimensionality(`str`): Expected dimension of the variable.
 528
 529        Returns:
 530            (`bool`): `True` if the variable if of the expected dimensionality, `False` otherwise.
 531
 532        Note:
 533            - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality).
 534            - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"`    
 535        """
 536        correct_dimensionality=variable.check(f"{expected_dimensionality}")      
 537        if not correct_dimensionality:
 538            raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}")
 539        return correct_dimensionality
 540
 541    def check_if_df_cell_has_a_value(self, index,key):
 542        """
 543        Checks if a cell in the `pmb.df` at the specified index and column has a value.
 544
 545        Args:
 546            index(`int`): Index of the row to check.
 547            key(`str`): Column label to check.
 548
 549        Returns:
 550            `bool`: `True` if the cell has a value, `False` otherwise.
 551        """
 552        idx = pd.IndexSlice
 553        import warnings
 554        with warnings.catch_warnings():
 555            warnings.simplefilter("ignore")
 556            return not pd.isna(self.df.loc[index, idx[key]])
 557
 558    def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined):
 559        """
 560        Checks if `name` is defined in `pmb.df`.
 561
 562        Args:
 563            name(`str`): label to check if defined in `pmb.df`.
 564            pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`.
 565
 566        Returns:
 567            `bool`: `True` for success, `False` otherwise.
 568        """
 569        if name in self.df['name'].unique():
 570            current_object_type = self.df[self.df['name']==name].pmb_type.values[0]
 571            if current_object_type != pmb_type_to_be_defined:
 572                raise ValueError (f"The name {name} is already defined in the df with a pmb_type = {current_object_type}, pymMBE does not support objects with the same name but different pmb_types")
 573            return True            
 574        else:
 575            return False
 576
 577    def check_if_metal_ion(self,key):
 578        """
 579        Checks if `key` corresponds to a label of a supported metal ion.
 580
 581        Args:
 582            key(`str`): key to be checked
 583
 584        Returns:
 585            (`bool`): True if `key`  is a supported metal ion, False otherwise.
 586        """
 587        if key in self.get_metal_ions_charge_number_map().keys():
 588            return True
 589        else:
 590            return False
 591
 592    def check_pka_set(self, pka_set):
 593        """
 594        Checks that `pka_set` has the formatting expected by the pyMBE library.
 595       
 596        Args:
 597            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
 598        """
 599        required_keys=['pka_value','acidity']
 600        for required_key in required_keys:
 601            for pka_name, pka_entry in pka_set.items():
 602                if required_key not in pka_entry:
 603                    raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")')
 604        return
 605
 606    def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")):
 607        """
 608        Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a pd.NA value.
 609
 610        Args:
 611            index(`int`): Index of the row to clean.
 612            columns_keys_to_clean(`list` of `str`, optional): List with the column keys to be cleaned. Defaults to [`particle_id`, `particle_id2`, `residue_id`, `molecule_id`].
 613        """   
 614        for column_key in columns_keys_to_clean:
 615            self.add_value_to_df(key=(column_key,''),index=index,new_value=pd.NA)
 616        self.df.fillna(pd.NA, inplace=True)
 617        return
 618
 619    def convert_columns_to_original_format (self, df):
 620        """
 621        Converts the columns of the Dataframe to the original format in pyMBE.
 622        
 623        Args:
 624            df(`DataFrame`): dataframe with pyMBE information as a string  
 625        
 626        """
 627
 628        columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', ('state_one','es_type'),('state_two','es_type'),('state_one','z'),('state_two','z') ]  
 629
 630        columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset']
 631
 632        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence', 'chain_map', 'node_map']
 633
 634        for column_name in columns_dtype_int:
 635            df[column_name] = df[column_name].astype(pd.Int64Dtype())
 636            
 637        for column_name in columns_with_list_or_dict:
 638            if df[column_name].isnull().all():
 639                df[column_name] = df[column_name].astype(object)
 640            else:
 641                df[column_name] = df[column_name].apply(lambda x: json.loads(x) if pd.notnull(x) else x)
 642
 643        for column_name in columns_with_units:
 644            df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x)
 645
 646        df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x)
 647        df["l0"] = df["l0"].astype(object)
 648        df["pka"] = df["pka"].astype(object)
 649        return df
 650    
 651    def convert_str_to_bond_object (self, bond_str):
 652        """
 653        Convert a row read as a `str` to the corresponding ESPResSo bond object. 
 654
 655        Args:
 656            bond_str(`str`): string with the information of a bond object.
 657
 658        Returns:
 659            bond_object(`obj`): ESPResSo bond object.
 660
 661        Note:
 662            Current supported bonds are: HarmonicBond and FeneBond
 663        """
 664        import espressomd.interactions
 665
 666        supported_bonds = ['HarmonicBond', 'FeneBond']
 667        m = re.search(r'^([A-Za-z0-9_]+)\((\{.+\})\)$', bond_str)
 668        if m is None:
 669            raise ValueError(f'Cannot parse bond "{bond_str}"')
 670        bond = m.group(1)
 671        if bond not in supported_bonds:
 672            raise NotImplementedError(f"Bond type '{bond}' currently not implemented in pyMBE, accepted types are {supported_bonds}")
 673        params = json.loads(m.group(2))
 674        bond_id = params.pop("bond_id")
 675        bond_object = getattr(espressomd.interactions, bond)(**params)
 676        bond_object._bond_id = bond_id
 677        return bond_object
 678
 679    def copy_df_entry(self, name, column_name, number_of_copies):
 680        '''
 681        Creates 'number_of_copies' of a given 'name' in `pymbe.df`.
 682
 683        Args:
 684            name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df`
 685            column_name(`str`): Column name to use as a filter. 
 686            number_of_copies(`int`): number of copies of `name` to be created.
 687        
 688        Note:
 689            - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 
 690        '''
 691
 692        valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ]
 693        if column_name not in valid_column_names:
 694            raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}")
 695        df_by_name = self.df.loc[self.df.name == name]
 696        if number_of_copies != 1:           
 697            if df_by_name[column_name].isnull().values.any():       
 698                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True)
 699            else:
 700                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
 701                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
 702                df_by_name_repeated[column_name] = pd.NA
 703            # Concatenate the new particle rows to  `df`
 704            self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
 705        else:
 706            if not df_by_name[column_name].isnull().values.any():     
 707                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
 708                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
 709                df_by_name_repeated[column_name] = pd.NA
 710                self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
 711        return
 712
 713    def create_added_salt (self, espresso_system, cation_name, anion_name, c_salt, verbose=True):    
 714        """
 715        Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`.
 716
 717        Args:
 718            espresso_system(`espressomd.system.System`): instance of an espresso system object.
 719            cation_name(`str`): `name` of a particle with a positive charge.
 720            anion_name(`str`): `name` of a particle with a negative charge.
 721            c_salt(`float`): Salt concentration.
 722            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
 723            
 724        Returns:
 725            c_salt_calculated(`float`): Calculated salt concentration added to `espresso_system`.
 726        """
 727        cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
 728        anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
 729        if cation_name_charge <= 0:
 730            raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge)
 731        if anion_name_charge >= 0:
 732            raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge)
 733        # Calculate the number of ions in the simulation box
 734        volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3')
 735        if c_salt.check('[substance] [length]**-3'):
 736            N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude)
 737            c_salt_calculated=N_ions/(volume*self.N_A)
 738        elif c_salt.check('[length]**-3'):
 739            N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude)
 740            c_salt_calculated=N_ions/volume
 741        else:
 742            raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt)
 743        N_cation = N_ions*abs(anion_name_charge)
 744        N_anion = N_ions*abs(cation_name_charge)
 745        self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation)
 746        self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion)
 747        if verbose:
 748            if c_salt_calculated.check('[substance] [length]**-3'):
 749                print(f"\n Added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions")
 750            elif c_salt_calculated.check('[length]**-3'):
 751                print(f"\n Added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions")
 752        return c_salt_calculated
 753
 754    def create_bond_in_espresso(self, bond_type, bond_parameters):
 755        '''
 756        Creates either a harmonic or a FENE bond in ESPResSo
 757
 758        Args:
 759            bond_type(`str`): label to identify the potential to model the bond.
 760            bond_parameters(`dict`): parameters of the potential of the bond.
 761
 762        Note:
 763            Currently, only HARMONIC and FENE bonds are supported.
 764
 765            For a HARMONIC bond the dictionary must contain:
 766
 767                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
 768                using the `pmb.units` UnitRegistry.
 769                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
 770                the `pmb.units` UnitRegistry.
 771           
 772            For a FENE bond the dictionary must additionally contain:
 773                
 774                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
 775                units of length using the `pmb.units` UnitRegistry. Default 'None'.
 776
 777        Returns:
 778              bond_object (`obj`): an ESPResSo bond object
 779        '''
 780        from espressomd import interactions
 781
 782        valid_bond_types   = ["harmonic", "FENE"]
 783        
 784        if 'k' in bond_parameters:
 785            bond_magnitude     = bond_parameters['k'].to('reduced_energy / reduced_length**2')
 786        else:
 787            raise ValueError("Magnitude of the potential (k) is missing")
 788        
 789        if bond_type == 'harmonic':
 790            if 'r_0' in bond_parameters:
 791                bond_length        = bond_parameters['r_0'].to('reduced_length')
 792            else:
 793                raise ValueError("Equilibrium bond length (r_0) is missing")
 794            bond_object    = interactions.HarmonicBond(k   = bond_magnitude.magnitude,
 795                                                       r_0 = bond_length.magnitude)
 796        elif bond_type == 'FENE':
 797            if 'r_0' in bond_parameters:
 798                bond_length        = bond_parameters['r_0'].to('reduced_length').magnitude
 799            else:
 800                print("WARNING: No value provided for r_0. Defaulting to r_0 = 0")
 801                bond_length=0
 802            if 'd_r_max' in bond_parameters:
 803                max_bond_stret = bond_parameters['d_r_max'].to('reduced_length')
 804            else:
 805                raise ValueError("Maximal stretching length (d_r_max) is missing")
 806            bond_object    = interactions.FeneBond(r_0     = bond_length, 
 807                                                   k       = bond_magnitude.magnitude,
 808                                                   d_r_max = max_bond_stret.magnitude)
 809        else:
 810            raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}")
 811        return bond_object
 812
 813
 814    def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True):
 815        """
 816        Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`.
 817        
 818        Args:
 819            object_name(`str`): `name` of a pymbe object.
 820            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 821            cation_name(`str`): `name` of a particle with a positive charge.
 822            anion_name(`str`): `name` of a particle with a negative charge.
 823            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
 824
 825        Returns: 
 826            counterion_number(`dict`): {"name": number}
 827        """
 828        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0]
 829        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0]
 830        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
 831        counterion_number={}
 832        object_charge={}
 833        for name in ['positive', 'negative']:
 834            object_charge[name]=0
 835        for id in object_ids:
 836            if espresso_system.part.by_id(id).q > 0:
 837                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 838            elif espresso_system.part.by_id(id).q < 0:
 839                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
 840        if object_charge['positive'] % abs(anion_charge) == 0:
 841            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
 842        else:
 843            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
 844        if object_charge['negative'] % abs(cation_charge) == 0:
 845            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
 846        else:
 847            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
 848        if counterion_number[cation_name] > 0: 
 849            self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name])
 850        else:
 851            counterion_number[cation_name]=0
 852        if counterion_number[anion_name] > 0:
 853            self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name])
 854        else:
 855            counterion_number[anion_name] = 0
 856        if verbose:
 857            print('The following counter-ions have been created: ')
 858            for name in counterion_number.keys():
 859                print(f'Ion type: {name} created number: {counterion_number[name]}')
 860        return counterion_number
 861
 862    def create_hydrogel(self, name, espresso_system):
 863        """ 
 864        creates the hydrogel `name` in espresso_system
 865        Args:
 866            name(`str`): Label of the hydrogel to be created. `name` must be defined in the `pmb.df`
 867            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 868
 869        Returns:
 870            hydrogel_info(`dict`):  {"name":hydrogel_name, "chains": {chain_id1: {residue_id1: {'central_bead_id': central_bead_id, 'side_chain_ids': [particle_id1,...]},...,"node_start":node_start,"node_end":node_end}, chain_id2: {...},...}, "nodes":{node1:[node1_id],...}}     
 871        """
 872        if not self.check_if_name_is_defined_in_df(name=name, pmb_type_to_be_defined='hydrogel'):
 873            raise ValueError(f"Hydrogel with name '{name}' is not defined in the DataFrame.")
 874        hydrogel_info={"name":name, "chains":{}, "nodes":{}}
 875        # placing nodes
 876        node_positions = {}
 877        node_topology = self.df[self.df["name"]==name]["node_map"].iloc[0]
 878        for node_info in node_topology:
 879            node_index = node_info["lattice_index"]
 880            node_name = node_info["particle_name"]
 881            node_pos, node_id = self.create_hydrogel_node(self.format_node(node_index), node_name, espresso_system)
 882            hydrogel_info["nodes"][self.format_node(node_index)]=node_id
 883            node_positions[node_id[0]]=node_pos
 884        
 885        # Placing chains between nodes
 886        # Looping over all the 16 chains
 887        chain_topology = self.df[self.df["name"]==name]["chain_map"].iloc[0]
 888        for chain_info in chain_topology:
 889            node_s = chain_info["node_start"]
 890            node_e = chain_info["node_end"]
 891            molecule_info = self.create_hydrogel_chain(node_s, node_e, node_positions, espresso_system)
 892            for molecule_id in molecule_info:
 893                hydrogel_info["chains"][molecule_id] = molecule_info[molecule_id]
 894                hydrogel_info["chains"][molecule_id]["node_start"]=node_s
 895                hydrogel_info["chains"][molecule_id]["node_end"]=node_e
 896        return hydrogel_info
 897
 898    def create_hydrogel_chain(self, node_start, node_end, node_positions, espresso_system):
 899        """
 900        Creates a chain between two nodes of a hydrogel.
 901
 902        Args:
 903            node_start(`str`): name of the starting node particle at which the first residue of the chain will be attached. 
 904            node_end(`str`): name of the ending node particle at which the last residue of the chain will be attached.
 905            node_positions(`dict`): dictionary with the positions of the nodes. The keys are the node names and the values are the positions.
 906            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 907
 908        Note:
 909            For example, 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.
 910            The chain will be placed in the direction of the vector between `node_start` and `node_end`. 
 911        """
 912        if self.lattice_builder is None:
 913            raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.")
 914
 915        molecule_name = "chain_"+node_start+"_"+node_end
 916        sequence = self.df[self.df['name']==molecule_name].residue_list.values [0]
 917        assert len(sequence) != 0 and not isinstance(sequence, str)
 918        assert len(sequence) == self.lattice_builder.mpc
 919
 920        key, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end)
 921        assert node_start != node_end or sequence == sequence[::-1], \
 922            (f"chain cannot be defined between '{node_start}' and '{node_end}' since it "
 923            "would form a loop with a non-symmetric sequence (under-defined stereocenter)")
 924
 925        if reverse:
 926            sequence = sequence[::-1]
 927
 928        node_start_pos = np.array(list(int(x) for x in node_start.strip('[]').split()))*0.25*self.lattice_builder.box_l
 929        node_end_pos = np.array(list(int(x) for x in node_end.strip('[]').split()))*0.25*self.lattice_builder.box_l
 930        node1 = espresso_system.part.select(lambda p: (p.pos == node_start_pos).all()).id
 931        node2 = espresso_system.part.select(lambda p: (p.pos == node_end_pos).all()).id
 932
 933        if not node1[0] in node_positions or not node2[0] in node_positions:
 934            raise ValueError("Set node position before placing a chain between them")
 935         
 936        # Finding a backbone vector between node_start and node_end
 937        vec_between_nodes = np.array(node_positions[node2[0]]) - np.array(node_positions[node1[0]])
 938        vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l)
 939        backbone_vector = list(vec_between_nodes/(self.lattice_builder.mpc + 1))
 940        node_start_name = self.df[(self.df["particle_id"]==node1[0]) & (self.df["pmb_type"]=="particle")]["name"].values[0]
 941        first_res_name = self.df[(self.df["pmb_type"]=="residue") & (self.df["name"]==sequence[0])]["central_bead"].values[0]
 942        l0 = self.get_bond_length(node_start_name, first_res_name, hard_check=True)
 943        chain_molecule_info = self.create_molecule(
 944                name=molecule_name,  # Use the name defined earlier
 945                number_of_molecules=1,  # Creating one chain
 946                espresso_system=espresso_system,
 947                list_of_first_residue_positions=[list(np.array(node_positions[node1[0]]) + np.array(backbone_vector))],#Start at the first node
 948                backbone_vector=np.array(backbone_vector)/l0,
 949                use_default_bond=False  # Use defaut bonds between monomers
 950            )
 951        # Collecting ids of beads of the chain/molecule
 952        chain_ids = []
 953        residue_ids = []
 954        for molecule_id in chain_molecule_info:
 955            for residue_id in chain_molecule_info[molecule_id]:
 956                residue_ids.append(residue_id)
 957                bead_id = chain_molecule_info[molecule_id][residue_id]['central_bead_id']
 958                chain_ids.append(bead_id)
 959
 960        self.lattice_builder.chains[key] = sequence
 961        # Search bonds between nodes and chain ends
 962        BeadType_near_to_node_start = self.df[(self.df["residue_id"] == residue_ids[0]) & (self.df["central_bead"].notnull())]["central_bead"].drop_duplicates().iloc[0]
 963        BeadType_near_to_node_end = self.df[(self.df["residue_id"] == residue_ids[-1]) & (self.df["central_bead"].notnull())]["central_bead"].drop_duplicates().iloc[0]
 964        bond_node1_first_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_start],
 965                                                    particle_name2 = BeadType_near_to_node_start,
 966                                                    hard_check=False,
 967                                                    use_default_bond=False)
 968        bond_node2_last_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_end],
 969                                                    particle_name2 = BeadType_near_to_node_end,
 970                                                    hard_check=False,
 971                                                    use_default_bond=False)
 972
 973        espresso_system.part.by_id(node1[0]).add_bond((bond_node1_first_monomer, chain_ids[0]))
 974        espresso_system.part.by_id(node2[0]).add_bond((bond_node2_last_monomer, chain_ids[-1]))
 975        # Add bonds to data frame
 976        bond_index1 = self.add_bond_in_df(particle_id1=node1[0],
 977                                    particle_id2=chain_ids[0],
 978                                    use_default_bond=False)
 979        self.add_value_to_df(key=('molecule_id',''),
 980                        index=int(bond_index1),
 981                        new_value=molecule_id,
 982                        overwrite=True)
 983        self.add_value_to_df(key=('residue_id',''),
 984                            index=int (bond_index1),
 985                            new_value=residue_ids[0],
 986                            overwrite=True)
 987        bond_index2 = self.add_bond_in_df(particle_id1=node2[0],
 988                                    particle_id2=chain_ids[-1],
 989                                    use_default_bond=False)
 990        self.add_value_to_df(key=('molecule_id',''),
 991                        index=int(bond_index2),
 992                        new_value=molecule_id,
 993                        overwrite=True)
 994        self.add_value_to_df(key=('residue_id',''),
 995                        index=int (bond_index2),
 996                        new_value=residue_ids[-1],
 997                        overwrite=True)
 998        return chain_molecule_info
 999    
1000    def create_hydrogel_node(self, node_index, node_name, espresso_system):
1001        """
1002        Set a node residue type.
1003        
1004        Args:
1005            node_index(`str`): Lattice node index in the form of a string, e.g. "[0 0 0]".
1006            node_name(`str`): name of the node particle defined in pyMBE.
1007        Returns:
1008            node_position(`list`): Position of the node in the lattice.
1009            p_id(`int`): Particle ID of the node.
1010        """
1011        if self.lattice_builder is None:
1012            raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.")
1013
1014        node_position = np.array(list(int(x) for x in node_index.strip('[]').split()))*0.25*self.lattice_builder.box_l
1015        p_id = self.create_particle(name = node_name,
1016                         espresso_system=espresso_system,
1017                         number_of_particles=1,
1018                         position = [node_position])
1019        key = self.lattice_builder._get_node_by_label(node_index)
1020        self.lattice_builder.nodes[key] = node_name
1021
1022        return node_position.tolist(), p_id
1023
1024    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False):
1025        """
1026        Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1027
1028        Args:
1029            name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df`
1030            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1031            number_of_molecules(`int`): Number of molecules of type `name` to be created.
1032            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.
1033            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. 
1034            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df`
1035
1036        Returns:
1037            molecules_info(`dict`):  {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 
1038        """
1039        if list_of_first_residue_positions is not None:
1040            for item in list_of_first_residue_positions:
1041                if not isinstance(item, list):
1042                    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.")
1043                elif len(item) != 3:
1044                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
1045
1046            if len(list_of_first_residue_positions) != number_of_molecules:
1047                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
1048        if number_of_molecules <= 0:
1049            return 0
1050        # Generate an arbitrary random unit vector
1051        if backbone_vector is None:
1052            backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
1053                                                    radius=1, 
1054                                                    n_samples=1,
1055                                                    on_surface=True)[0]
1056        else:
1057            backbone_vector = np.array(backbone_vector)
1058        self.check_if_name_is_defined_in_df(name=name,
1059                                            pmb_type_to_be_defined='molecule')
1060            
1061        first_residue = True
1062        molecules_info = {}
1063        residue_list = self.df[self.df['name']==name].residue_list.values [0]
1064        self.copy_df_entry(name=name,
1065                        column_name='molecule_id',
1066                        number_of_copies=number_of_molecules)
1067
1068        molecules_index = np.where(self.df['name']==name)
1069        molecule_index_list =list(molecules_index[0])[-number_of_molecules:]
1070        pos_index = 0 
1071        for molecule_index in molecule_index_list:        
1072            molecule_id = self.assign_molecule_id(molecule_index=molecule_index)
1073            molecules_info[molecule_id] = {}
1074            for residue in residue_list:
1075                if first_residue:
1076                    if list_of_first_residue_positions is None:
1077                        residue_position = None
1078                    else:
1079                        for item in list_of_first_residue_positions:
1080                            residue_position = [np.array(list_of_first_residue_positions[pos_index])]
1081                            
1082                    residues_info = self.create_residue(name=residue,
1083                                                        espresso_system=espresso_system, 
1084                                                        central_bead_position=residue_position,  
1085                                                        use_default_bond= use_default_bond, 
1086                                                        backbone_vector=backbone_vector)
1087                    residue_id = next(iter(residues_info))
1088                    # Add the correct molecule_id to all particles in the residue
1089                    for index in self.df[self.df['residue_id']==residue_id].index:
1090                        self.add_value_to_df(key=('molecule_id',''),
1091                                            index=int (index),
1092                                            new_value=molecule_id,
1093                                            overwrite=True)
1094                    central_bead_id = residues_info[residue_id]['central_bead_id']
1095                    previous_residue = residue
1096                    residue_position = espresso_system.part.by_id(central_bead_id).pos
1097                    previous_residue_id = central_bead_id
1098                    first_residue = False          
1099                else:                    
1100                    previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0]
1101                    new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0]
1102                    bond = self.search_bond(particle_name1=previous_central_bead_name, 
1103                                            particle_name2=new_central_bead_name, 
1104                                            hard_check=True, 
1105                                            use_default_bond=use_default_bond)                
1106                    l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 
1107                                            particle_name2=new_central_bead_name, 
1108                                            hard_check=True, 
1109                                            use_default_bond=use_default_bond)                
1110                    
1111                    residue_position = residue_position+backbone_vector*l0
1112                    residues_info = self.create_residue(name=residue, 
1113                                                        espresso_system=espresso_system, 
1114                                                        central_bead_position=[residue_position],
1115                                                        use_default_bond= use_default_bond, 
1116                                                        backbone_vector=backbone_vector)
1117                    residue_id = next(iter(residues_info))      
1118                    for index in self.df[self.df['residue_id']==residue_id].index:
1119                        self.add_value_to_df(key=('molecule_id',''),
1120                                            index=int (index),
1121                                            new_value=molecule_id,
1122                                            overwrite=True)            
1123                    central_bead_id = residues_info[residue_id]['central_bead_id']
1124                    espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id))
1125                    bond_index = self.add_bond_in_df(particle_id1=central_bead_id,
1126                                        particle_id2=previous_residue_id,
1127                                        use_default_bond=use_default_bond) 
1128                    self.add_value_to_df(key=('molecule_id',''),
1129                                            index=int (bond_index),
1130                                            new_value=molecule_id,
1131                                            overwrite=True)           
1132                    previous_residue_id = central_bead_id
1133                    previous_residue = residue                    
1134                molecules_info[molecule_id][residue_id] = residues_info[residue_id]
1135            first_residue = True
1136            pos_index+=1
1137        
1138        return molecules_info
1139    
1140    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
1141        """
1142        Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`.
1143        
1144        Args:
1145            name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 
1146            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1147            number_of_particles(`int`): Number of particles to be created.
1148            position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
1149            fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
1150        Returns:
1151            created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`.
1152        """       
1153        if number_of_particles <=0:
1154            return []
1155        self.check_if_name_is_defined_in_df(name=name,
1156                                       pmb_type_to_be_defined='particle')
1157        # Copy the data of the particle `number_of_particles` times in the `df`
1158        self.copy_df_entry(name=name,
1159                          column_name='particle_id',
1160                          number_of_copies=number_of_particles)
1161        # Get information from the particle type `name` from the df     
1162        z = self.df.loc[self.df['name']==name].state_one.z.values[0]
1163        z = 0. if z is None else z
1164        es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
1165        # Get a list of the index in `df` corresponding to the new particles to be created
1166        index = np.where(self.df['name']==name)
1167        index_list =list(index[0])[-number_of_particles:]
1168        # Create the new particles into  `espresso_system`
1169        created_pid_list=[]
1170        for index in range (number_of_particles):
1171            df_index=int (index_list[index])
1172            self.clean_df_row(index=df_index)
1173            if position is None:
1174                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1175            else:
1176                particle_position = position[index]
1177            if len(espresso_system.part.all()) == 0:
1178                bead_id = 0
1179            else:
1180                bead_id = max (espresso_system.part.all().id) + 1
1181            created_pid_list.append(bead_id)
1182            kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z)
1183            if fix:
1184                kwargs["fix"] = 3 * [fix]
1185            espresso_system.part.add(**kwargs)
1186            self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id)                  
1187        return created_pid_list
1188
1189    def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None):
1190        """
1191        Creates all `particle`s associated to `pmb object` into  `espresso` a number of times equal to `number_of_objects`.
1192        
1193        Args:
1194            name(`str`): Unique label of the `pmb object` to be created. 
1195            number_of_objects(`int`): Number of `pmb object`s to be created.
1196            espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library.
1197            position(`list`): Coordinates where the particles should be created.
1198            use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`.
1199            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. 
1200
1201        Note:
1202            - If no `position` is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length. 
1203        """
1204        allowed_objects=['particle','residue','molecule']
1205        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1206        if pmb_type not in allowed_objects:
1207            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1208        if pmb_type == 'particle':
1209            self.create_particle(name=name, 
1210                                number_of_particles=number_of_objects, 
1211                                espresso_system=espresso_system, 
1212                                position=position)
1213        elif pmb_type == 'residue':
1214            self.create_residue(name=name,  
1215                                espresso_system=espresso_system, 
1216                                central_bead_position=position,
1217                                use_default_bond=use_default_bond,
1218                                backbone_vector=backbone_vector)
1219        elif pmb_type == 'molecule':
1220            self.create_molecule(name=name, 
1221                                number_of_molecules=number_of_objects, 
1222                                espresso_system=espresso_system, 
1223                                use_default_bond=use_default_bond, 
1224                                list_of_first_residue_positions=position,
1225                                backbone_vector=backbone_vector)
1226        return
1227
1228    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1229        """
1230        Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions`
1231
1232        Args:
1233            name(`str`): Label of the protein to be created. 
1234            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1235            number_of_proteins(`int`): Number of proteins to be created.
1236            positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}}
1237        """
1238
1239        if number_of_proteins <=0:
1240            return
1241        self.check_if_name_is_defined_in_df(name=name,
1242                                        pmb_type_to_be_defined='protein')
1243        self.copy_df_entry(name=name,
1244                            column_name='molecule_id',
1245                            number_of_copies=number_of_proteins)
1246        protein_index = np.where(self.df['name']==name)
1247        protein_index_list =list(protein_index[0])[-number_of_proteins:]
1248        
1249        box_half=espresso_system.box_l[0]/2.0
1250
1251        for molecule_index in protein_index_list:     
1252
1253            molecule_id = self.assign_molecule_id(molecule_index=molecule_index)
1254
1255            protein_center = self.generate_coordinates_outside_sphere(radius = 1, 
1256                                                                        max_dist=box_half, 
1257                                                                        n_samples=1, 
1258                                                                        center=[box_half]*3)[0]
1259   
1260            for residue in topology_dict.keys():
1261
1262                residue_name = re.split(r'\d+', residue)[0]
1263                residue_number = re.split(r'(\d+)', residue)[1]
1264                residue_position = topology_dict[residue]['initial_pos']
1265                position = residue_position + protein_center
1266
1267                particle_id = self.create_particle(name=residue_name,
1268                                                            espresso_system=espresso_system,
1269                                                            number_of_particles=1,
1270                                                            position=[position], 
1271                                                            fix = True)
1272                
1273                index = self.df[self.df['particle_id']==particle_id[0]].index.values[0]
1274                self.add_value_to_df(key=('residue_id',''),
1275                                            index=int (index),
1276                                            new_value=int(residue_number),
1277                                            overwrite=True)
1278
1279                self.add_value_to_df(key=('molecule_id',''),
1280                                        index=int (index),
1281                                        new_value=molecule_id,
1282                                        overwrite=True)
1283
1284        return
1285
1286    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None):
1287        """
1288        Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1289
1290        Args:
1291            name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df`
1292            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1293            central_bead_position(`list` of `float`): Position of the central bead.
1294            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 `pmb.df`
1295            backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`.
1296
1297        Returns:
1298            residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1299        """
1300        self.check_if_name_is_defined_in_df(name=name,
1301                                            pmb_type_to_be_defined='residue')
1302        # Copy the data of a residue in the `df
1303        self.copy_df_entry(name=name,
1304                            column_name='residue_id',
1305                            number_of_copies=1)
1306        residues_index = np.where(self.df['name']==name)
1307        residue_index_list =list(residues_index[0])[-1:]
1308        # search for defined particle and residue names
1309        particle_and_residue_df = self.df.loc[(self.df['pmb_type']== "particle") | (self.df['pmb_type']== "residue")]
1310        particle_and_residue_names = particle_and_residue_df["name"].tolist()
1311        for residue_index in residue_index_list:  
1312            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1313            for side_chain_element in side_chain_list:
1314                if side_chain_element not in particle_and_residue_names:              
1315                    raise ValueError (f"{side_chain_element} is not defined")
1316        # Internal bookkepping of the residue info (important for side-chain residues)
1317        # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1318        residues_info={}
1319        for residue_index in residue_index_list:     
1320            self.clean_df_row(index=int(residue_index))
1321            # Assign a residue_id
1322            if self.df['residue_id'].isnull().all():
1323                residue_id=0
1324            else:
1325                residue_id = self.df['residue_id'].max() + 1
1326            self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id)
1327            # create the principal bead   
1328            central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0]            
1329            central_bead_id = self.create_particle(name=central_bead_name,
1330                                                                espresso_system=espresso_system,
1331                                                                position=central_bead_position,
1332                                                                number_of_particles = 1)[0]
1333            central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1334            #assigns same residue_id to the central_bead particle created.
1335            index = self.df[self.df['particle_id']==central_bead_id].index.values[0]
1336            self.df.at [index,'residue_id'] = residue_id
1337            # Internal bookkeeping of the central bead id
1338            residues_info[residue_id]={}
1339            residues_info[residue_id]['central_bead_id']=central_bead_id
1340            # create the lateral beads  
1341            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1342            side_chain_beads_ids = []
1343            for side_chain_element in side_chain_list:  
1344                pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 
1345                if pmb_type == 'particle':
1346                    bond = self.search_bond(particle_name1=central_bead_name, 
1347                                            particle_name2=side_chain_element, 
1348                                            hard_check=True, 
1349                                            use_default_bond=use_default_bond)
1350                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1351                                              particle_name2=side_chain_element, 
1352                                              hard_check=True, 
1353                                              use_default_bond=use_default_bond)
1354
1355                    if backbone_vector is None:
1356                        bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1357                                                                 radius=l0, 
1358                                                                 n_samples=1,
1359                                                                 on_surface=True)[0]
1360                    else:
1361                        bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector),
1362                                                                                                    magnitude=l0)
1363                     
1364                    side_bead_id = self.create_particle(name=side_chain_element, 
1365                                                                    espresso_system=espresso_system,
1366                                                                    position=[bead_position], 
1367                                                                    number_of_particles=1)[0]
1368                    index = self.df[self.df['particle_id']==side_bead_id].index.values[0]
1369                    self.add_value_to_df(key=('residue_id',''),
1370                                        index=int (index),
1371                                        new_value=residue_id, 
1372                                        overwrite=True)
1373                    side_chain_beads_ids.append(side_bead_id)
1374                    espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id))
1375                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1376                                        particle_id2=side_bead_id,
1377                                        use_default_bond=use_default_bond)
1378                    self.add_value_to_df(key=('residue_id',''),
1379                                        index=int (index),
1380                                        new_value=residue_id, 
1381                                        overwrite=True)
1382
1383                elif pmb_type == 'residue':
1384                    central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0]
1385                    bond = self.search_bond(particle_name1=central_bead_name, 
1386                                            particle_name2=central_bead_side_chain, 
1387                                            hard_check=True, 
1388                                            use_default_bond=use_default_bond)
1389                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1390                                              particle_name2=central_bead_side_chain, 
1391                                              hard_check=True, 
1392                                              use_default_bond=use_default_bond)
1393                    if backbone_vector is None:
1394                        residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1395                                                                    radius=l0, 
1396                                                                    n_samples=1,
1397                                                                    on_surface=True)[0]
1398                    else:
1399                        residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1400                                                                                                        magnitude=l0)
1401                    lateral_residue_info = self.create_residue(name=side_chain_element, 
1402                                                                espresso_system=espresso_system,
1403                                                                central_bead_position=[residue_position],
1404                                                                use_default_bond=use_default_bond)
1405                    lateral_residue_dict=list(lateral_residue_info.values())[0]
1406                    central_bead_side_chain_id=lateral_residue_dict['central_bead_id']
1407                    lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids']
1408                    residue_id_side_chain=list(lateral_residue_info.keys())[0]
1409                    # Change the residue_id of the residue in the side chain to the one of the bigger residue
1410                    index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0]
1411                    self.add_value_to_df(key=('residue_id',''),
1412                                        index=int(index),
1413                                        new_value=residue_id, 
1414                                        overwrite=True)
1415                    # Change the residue_id of the particles in the residue in the side chain
1416                    side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids
1417                    for particle_id in side_chain_beads_ids:
1418                        index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0]
1419                        self.add_value_to_df(key=('residue_id',''),
1420                                            index=int (index),
1421                                            new_value=residue_id, 
1422                                            overwrite=True)
1423                    espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id))
1424                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1425                                        particle_id2=central_bead_side_chain_id,
1426                                        use_default_bond=use_default_bond)
1427                    self.add_value_to_df(key=('residue_id',''),
1428                                        index=int (index),
1429                                        new_value=residue_id, 
1430                                        overwrite=True)
1431                    # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue
1432                    for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index:        
1433                        self.add_value_to_df(key=('residue_id',''),
1434                                            index=int(index),
1435                                            new_value=residue_id, 
1436                                            overwrite=True)
1437            # Internal bookkeeping of the side chain beads ids
1438            residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids
1439        return  residues_info
1440
1441    def create_variable_with_units(self, variable):
1442        """
1443        Returns a pint object with the value and units defined in `variable`.
1444
1445        Args:
1446            variable(`dict` or `str`): {'value': value, 'units': units}
1447        Returns:
1448            variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry.
1449        """        
1450        if isinstance(variable, dict):
1451            value=variable.pop('value')
1452            units=variable.pop('units')
1453        elif isinstance(variable, str):
1454            value = float(re.split(r'\s+', variable)[0])
1455            units = re.split(r'\s+', variable)[1]
1456        variable_with_units=value*self.units(units)
1457        return variable_with_units 
1458
1459    def define_AA_residues(self, sequence, model):
1460        """
1461        Defines in `pmb.df` all the different residues in `sequence`.
1462
1463        Args:
1464            sequence(`lst`):  Sequence of the peptide or protein.
1465            model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1466
1467        Returns:
1468            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.             
1469        """
1470        residue_list = []
1471        for residue_name in sequence:
1472            if model == '1beadAA':
1473                central_bead = residue_name
1474                side_chains = []
1475            elif model == '2beadAA':
1476                if residue_name in ['c','n', 'G']: 
1477                    central_bead = residue_name
1478                    side_chains = []
1479                else:
1480                    central_bead = 'CA'              
1481                    side_chains = [residue_name]
1482            if residue_name not in residue_list:   
1483                self.define_residue(name = 'AA-'+residue_name, 
1484                                    central_bead = central_bead,
1485                                    side_chains = side_chains)              
1486            residue_list.append('AA-'+residue_name)
1487        return residue_list
1488
1489    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1490        
1491        '''
1492        Defines a pmb object of type `bond` in `pymbe.df`.
1493
1494        Args:
1495            bond_type(`str`): label to identify the potential to model the bond.
1496            bond_parameters(`dict`): parameters of the potential of the bond.
1497            particle_pairs(`lst`): list of the `names` of the `particles` to be bonded.
1498
1499        Note:
1500            Currently, only HARMONIC and FENE bonds are supported.
1501
1502            For a HARMONIC bond the dictionary must contain the following parameters:
1503
1504                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
1505                using the `pmb.units` UnitRegistry.
1506                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
1507                the `pmb.units` UnitRegistry.
1508           
1509            For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:
1510                
1511                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
1512                units of length using the `pmb.units` UnitRegistry. Default 'None'.
1513        '''
1514
1515        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1516        for particle_name1, particle_name2 in particle_pairs:
1517
1518            lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1,
1519                                                 particle_name2 = particle_name2,
1520                                                 combining_rule = 'Lorentz-Berthelot')
1521      
1522            l0 = self.calculate_initial_bond_length(bond_object = bond_object,
1523                                                    bond_type   = bond_type,
1524                                                    epsilon     = lj_parameters["epsilon"],
1525                                                    sigma       = lj_parameters["sigma"],
1526                                                    cutoff      = lj_parameters["cutoff"],
1527                                                    offset      = lj_parameters["offset"],)
1528            index = len(self.df)
1529            for label in [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']:
1530                self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond")
1531            self.df.at [index,'name']= f'{particle_name1}-{particle_name2}'
1532            self.df.at [index,'bond_object'] = bond_object
1533            self.df.at [index,'l0'] = l0
1534            self.add_value_to_df(index=index,
1535                                    key=('pmb_type',''),
1536                                    new_value='bond')
1537            self.add_value_to_df(index=index,
1538                                    key=('parameters_of_the_potential',''),
1539                                    new_value=bond_object.get_params(),
1540                                    non_standard_value=True)
1541        self.df.fillna(pd.NA, inplace=True)
1542        return
1543
1544    
1545    def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None):
1546        """
1547        Asigns `bond` in `pmb.df` as the default bond.
1548        The LJ parameters can be optionally provided to calculate the initial bond length
1549
1550        Args:
1551            bond_type(`str`): label to identify the potential to model the bond.
1552            bond_parameters(`dict`): parameters of the potential of the bond.
1553            sigma(`float`, optional): LJ sigma of the interaction between the particles.
1554            epsilon(`float`, optional): LJ epsilon for the interaction between the particles.
1555            cutoff(`float`, optional): cutoff-radius of the LJ interaction.
1556            offset(`float`, optional): offset of the LJ interaction.
1557
1558        Note:
1559            - Currently, only harmonic and FENE bonds are supported. 
1560        """
1561        
1562        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1563
1564        if epsilon is None:
1565            epsilon=1*self.units('reduced_energy')
1566        if sigma is None:
1567            sigma=1*self.units('reduced_length')
1568        if cutoff is None:
1569            cutoff=2**(1.0/6.0)*self.units('reduced_length')
1570        if offset is None:
1571            offset=0*self.units('reduced_length')
1572        l0 = self.calculate_initial_bond_length(bond_object = bond_object, 
1573                                                bond_type   = bond_type, 
1574                                                epsilon     = epsilon,
1575                                                sigma       = sigma,
1576                                                cutoff      = cutoff,
1577                                                offset      = offset)
1578
1579        if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'):
1580            return
1581        if len(self.df.index) != 0:
1582            index = max(self.df.index)+1
1583        else:
1584            index = 0
1585        self.df.at [index,'name']        = 'default'
1586        self.df.at [index,'bond_object'] = bond_object
1587        self.df.at [index,'l0']          = l0
1588        self.add_value_to_df(index       = index,
1589                            key          = ('pmb_type',''),
1590                            new_value    = 'bond')
1591        self.add_value_to_df(index       = index,
1592                            key          = ('parameters_of_the_potential',''),
1593                            new_value    = bond_object.get_params(),
1594                            non_standard_value=True)
1595        self.df.fillna(pd.NA, inplace=True)
1596        return
1597    
1598    def define_hydrogel(self, name, node_map, chain_map):
1599        """
1600        Defines a pyMBE object of type `hydrogel` in `pymbe.df`.
1601
1602        Args:
1603            name(`str`): Unique label that identifies the `hydrogel`.
1604            node_map(`list of ict`): [{"particle_name": , "lattice_index": }, ... ]
1605            chain_map(`list of dict`): [{"node_start": , "node_end": , "residue_list": , ... ]
1606        """
1607        node_indices = {tuple(entry['lattice_index']) for entry in node_map}
1608        diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices}
1609        if node_indices != diamond_indices:
1610            raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ")
1611        
1612        chain_map_connectivity = set()
1613        for entry in chain_map:
1614            start = self.lattice_builder.node_labels[entry['node_start']]
1615            end = self.lattice_builder.node_labels[entry['node_end']]
1616            chain_map_connectivity.add((start,end))
1617
1618        if self.lattice_builder.lattice.connectivity != chain_map_connectivity:
1619            raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs")
1620
1621        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='hydrogel'):
1622            return
1623
1624        index = len(self.df)
1625        self.df.at [index, "name"] = name
1626        self.df.at [index, "pmb_type"] = "hydrogel"
1627        self.add_value_to_df(index = index,
1628                        key = ('node_map',''),
1629                        new_value = node_map,
1630                        non_standard_value=True)
1631        self.add_value_to_df(index = index,
1632                        key = ('chain_map',''),
1633                        new_value = chain_map,
1634                        non_standard_value=True)
1635        for chain_label in chain_map:
1636            node_start = chain_label["node_start"]
1637            node_end = chain_label["node_end"]
1638            residue_list = chain_label['residue_list']
1639            # Molecule name
1640            molecule_name = "chain_"+node_start+"_"+node_end
1641            self.define_molecule(name=molecule_name, residue_list=residue_list)
1642        return
1643
1644    def define_molecule(self, name, residue_list):
1645        """
1646        Defines a pyMBE object of type `molecule` in `pymbe.df`.
1647
1648        Args:
1649            name(`str`): Unique label that identifies the `molecule`.
1650            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.  
1651        """
1652        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'):
1653            return
1654        index = len(self.df)
1655        self.df.at [index,'name'] = name
1656        self.df.at [index,'pmb_type'] = 'molecule'
1657        self.df.at [index,('residue_list','')] = residue_list
1658        self.df.fillna(pd.NA, inplace=True)
1659        return
1660
1661    def define_particle_entry_in_df(self,name):
1662        """
1663        Defines a particle entry in pmb.df.
1664
1665        Args:
1666            name(`str`): Unique label that identifies this particle type.
1667
1668        Returns:
1669            index(`int`): Index of the particle in pmb.df  
1670        """
1671
1672        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'):
1673            index = self.df[self.df['name']==name].index[0]                                   
1674        else:
1675            index = len(self.df)
1676            self.df.at [index, 'name'] = name
1677            self.df.at [index,'pmb_type'] = 'particle'
1678        self.df.fillna(pd.NA, inplace=True)
1679        return index
1680
1681    def define_particle(self, name, z=0, acidity=pd.NA, pka=pd.NA, sigma=pd.NA, epsilon=pd.NA, cutoff=pd.NA, offset=pd.NA,overwrite=False):
1682        """
1683        Defines the properties of a particle object.
1684
1685        Args:
1686            name(`str`): Unique label that identifies this particle type.  
1687            z(`int`, optional): Permanent charge number of this particle type. Defaults to 0.
1688            acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to pd.NA.
1689            pka(`float`, optional): If `particle` is an acid or a base, it defines its  pka-value. Defaults to pd.NA.
1690            sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1691            cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1692            offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1693            epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA.
1694            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
1695
1696        Note:
1697            - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units.
1698            - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units.
1699            - `cutoff` defaults to `2**(1./6.) reduced_length`. 
1700            - `offset` defaults to 0.
1701            - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions.
1702            - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`.
1703        """ 
1704        index=self.define_particle_entry_in_df(name=name)
1705        
1706        # If `cutoff` and `offset` are not defined, default them to the following values
1707        if pd.isna(cutoff):
1708            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1709        if pd.isna(offset):
1710            offset=self.units.Quantity(0, "reduced_length")
1711
1712        # Define LJ parameters
1713        parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"},
1714                                        "cutoff":{"value": cutoff, "dimensionality": "[length]"},
1715                                        "offset":{"value": offset, "dimensionality": "[length]"},
1716                                        "epsilon":{"value": epsilon, "dimensionality": "[energy]"},}
1717
1718        for parameter_key in parameters_with_dimensionality.keys():
1719            if not pd.isna(parameters_with_dimensionality[parameter_key]["value"]):
1720                self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 
1721                                          expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"])
1722                self.add_value_to_df(key=(parameter_key,''),
1723                                    index=index,
1724                                    new_value=parameters_with_dimensionality[parameter_key]["value"],
1725                                    overwrite=overwrite)
1726
1727        # Define particle acid/base properties
1728        self.set_particle_acidity(name=name, 
1729                                acidity=acidity, 
1730                                default_charge_number=z, 
1731                                pka=pka,
1732                                overwrite=overwrite)
1733        self.df.fillna(pd.NA, inplace=True)
1734        return 
1735    
1736    def define_particles(self, parameters, overwrite=False):
1737        '''
1738        Defines a particle object in pyMBE for each particle name in `particle_names`
1739
1740        Args:
1741            parameters(`dict`):  dictionary with the particle parameters. 
1742            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1743
1744        Note:
1745            - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},}
1746        '''
1747        if not parameters:
1748            return 0
1749        for particle_name in parameters.keys():
1750            parameters[particle_name]["overwrite"]=overwrite
1751            self.define_particle(**parameters[particle_name])
1752        return
1753    
1754    def define_peptide(self, name, sequence, model):
1755        """
1756        Defines a pyMBE object of type `peptide` in the `pymbe.df`.
1757
1758        Args:
1759            name (`str`): Unique label that identifies the `peptide`.
1760            sequence (`string`): Sequence of the `peptide`.
1761            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1762        """
1763        if self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'):
1764            return
1765        valid_keys = ['1beadAA','2beadAA']
1766        if model not in valid_keys:
1767            raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
1768        clean_sequence = self.protein_sequence_parser(sequence=sequence)    
1769        residue_list = self.define_AA_residues(sequence=clean_sequence,
1770                                                model=model)
1771        self.define_molecule(name = name, residue_list=residue_list)
1772        index = self.df.loc[self.df['name'] == name].index.item() 
1773        self.df.at [index,'model'] = model
1774        self.df.at [index,('sequence','')] = clean_sequence
1775        self.df.fillna(pd.NA, inplace=True)
1776        
1777    
1778    def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False):
1779        """
1780        Defines a globular protein pyMBE object  in `pymbe.df`.
1781
1782        Args:
1783            name (`str`): Unique label that identifies the protein.
1784            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1785            topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value}
1786            lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca".
1787            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1788
1789        Note:
1790            - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential.
1791        """
1792
1793        # Sanity checks
1794        valid_model_keys = ['1beadAA','2beadAA']
1795        valid_lj_setups = ["wca"]
1796        if model not in valid_model_keys:
1797            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
1798        if lj_setup_mode not in valid_lj_setups:
1799            raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}')
1800        if lj_setup_mode == "wca":
1801            sigma = 1*self.units.Quantity("reduced_length")
1802            epsilon = 1*self.units.Quantity("reduced_energy")
1803        part_dict={}
1804        sequence=[]
1805        metal_ions_charge_number_map=self.get_metal_ions_charge_number_map()
1806        for particle in topology_dict.keys():
1807            particle_name = re.split(r'\d+', particle)[0] 
1808            if particle_name not in part_dict.keys():
1809                if lj_setup_mode == "wca":
1810                    part_dict[particle_name]={"sigma": sigma,
1811                                        "offset": topology_dict[particle]['radius']*2-sigma,
1812                                        "epsilon": epsilon,
1813                                        "name": particle_name}
1814                if self.check_if_metal_ion(key=particle_name):
1815                    z=metal_ions_charge_number_map[particle_name]
1816                else:
1817                    z=0
1818                part_dict[particle_name]["z"]=z
1819            
1820            if self.check_aminoacid_key(key=particle_name):
1821                sequence.append(particle_name) 
1822            
1823        self.define_particles(parameters=part_dict,
1824                            overwrite=overwrite)
1825        residue_list = self.define_AA_residues(sequence=sequence, 
1826                                               model=model)
1827        index = len(self.df)
1828        self.df.at [index,'name'] = name
1829        self.df.at [index,'pmb_type'] = 'protein'
1830        self.df.at [index,'model'] = model
1831        self.df.at [index,('sequence','')] = sequence  
1832        self.df.at [index,('residue_list','')] = residue_list
1833        self.df.fillna(pd.NA, inplace=True)    
1834        return 
1835    
1836    def define_residue(self, name, central_bead, side_chains):
1837        """
1838        Defines a pyMBE object of type `residue` in `pymbe.df`.
1839
1840        Args:
1841            name(`str`): Unique label that identifies the `residue`.
1842            central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`.
1843            side_chains(`list` of `str`): List of `name`s of the pmb_objects to be placed as side_chains of the `residue`. Currently, only pmb_objects of type `particle`s or `residue`s are supported.
1844        """
1845        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'):
1846            return
1847        index = len(self.df)
1848        self.df.at [index, 'name'] = name
1849        self.df.at [index,'pmb_type'] = 'residue'
1850        self.df.at [index,'central_bead'] = central_bead
1851        self.df.at [index,('side_chains','')] = side_chains
1852        self.df.fillna(pd.NA, inplace=True)
1853        return 
1854
1855    def delete_entries_in_df(self, entry_name):
1856        """
1857        Deletes entries with name `entry_name` from the DataFrame if it exists.
1858
1859        Args:
1860            entry_name (`str`): The name of the entry in the dataframe to delete.
1861
1862        """
1863        if entry_name in self.df["name"].values:
1864            self.df = self.df[self.df["name"] != entry_name].reset_index(drop=True)
1865
1866    def destroy_pmb_object_in_system(self, name, espresso_system):
1867        """
1868        Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 
1869
1870        Args:
1871            name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`.
1872            espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library.
1873
1874        Note:
1875            - If `name`  is a object_type=`particle`, only the matching particles that are not part of bigger objects (e.g. `residue`, `molecule`) will be destroyed. To destroy particles in such objects, destroy the bigger object instead.
1876        """
1877        allowed_objects = ['particle','residue','molecule']
1878        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1879        if pmb_type not in allowed_objects:
1880            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1881        if pmb_type == 'particle':
1882            particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 
1883                                                  & (self.df['molecule_id'].isna())]
1884            particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna())
1885                                                 & (self.df['molecule_id'].isna())].particle_id.tolist()
1886            for particle_id in particle_ids_list:
1887                espresso_system.part.by_id(particle_id).remove()
1888            self.df = self.df.drop(index=particles_index)
1889        if pmb_type == 'residue':
1890            residues_id = self.df.loc[self.df['name']== name].residue_id.to_list()
1891            for residue_id in residues_id:
1892                molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0]
1893                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1894                self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index)
1895                for particle_id in particle_ids_list:
1896                    espresso_system.part.by_id(particle_id).remove()
1897                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)    
1898        if pmb_type == 'molecule':
1899            molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list()
1900            for molecule_id in molecules_id:
1901                molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']=="molecule")].name.values[0]
1902                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1903                self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index)
1904                for particle_id in particle_ids_list:
1905                    espresso_system.part.by_id(particle_id).remove()   
1906                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)             
1907        
1908        self.df.reset_index(drop=True,inplace=True)
1909
1910        return
1911
1912    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
1913        """
1914        Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration.
1915        To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an
1916        ion pair and the concentrations of the various species is solved numerically using a self-consistent approach.
1917        More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831).
1918        This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237).
1919
1920        Args:
1921            pH_res('float'): pH-value in the reservoir.
1922            c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
1923            activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
1924
1925        Returns:
1926            cH_res('pint.Quantity'): Concentration of H+ ions.
1927            cOH_res('pint.Quantity'): Concentration of OH- ions.
1928            cNa_res('pint.Quantity'): Concentration of Na+ ions.
1929            cCl_res('pint.Quantity'): Concentration of Cl- ions.
1930        """
1931
1932        self_consistent_run = 0
1933        cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+
1934
1935        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
1936            #Calculate and initial guess for the concentrations of various species based on ideal gas estimate
1937            cOH_res = self.Kw / cH_res 
1938            cNa_res = None
1939            cCl_res = None
1940            if cOH_res>=cH_res:
1941                #adjust the concentration of sodium if there is excess OH- in the reservoir:
1942                cNa_res = c_salt_res + (cOH_res-cH_res)
1943                cCl_res = c_salt_res
1944            else:
1945                # adjust the concentration of chloride if there is excess H+ in the reservoir
1946                cCl_res = c_salt_res + (cH_res-cOH_res)
1947                cNa_res = c_salt_res
1948                
1949            def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res):
1950                nonlocal max_number_sc_runs, self_consistent_run
1951                if self_consistent_run<max_number_sc_runs:
1952                    self_consistent_run+=1
1953                    ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1954                    cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
1955                    if cOH_res>=cH_res:
1956                        #adjust the concentration of sodium if there is excess OH- in the reservoir:
1957                        cNa_res = c_salt_res + (cOH_res-cH_res)
1958                        cCl_res = c_salt_res
1959                    else:
1960                        # adjust the concentration of chloride if there is excess H+ in the reservoir
1961                        cCl_res = c_salt_res + (cH_res-cOH_res)
1962                        cNa_res = c_salt_res
1963                    return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1964                else:
1965                    return cH_res, cOH_res, cNa_res, cCl_res
1966            return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1967
1968        cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1969        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1970        determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1971
1972        while abs(determined_pH-pH_res)>1e-6:
1973            if determined_pH > pH_res:
1974                cH_res *= 1.005
1975            else:
1976                cH_res /= 1.003
1977            cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1978            ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1979            determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1980            self_consistent_run=0
1981
1982        return cH_res, cOH_res, cNa_res, cCl_res
1983
1984    def enable_motion_of_rigid_object(self, name, espresso_system):
1985        '''
1986        Enables the motion of the rigid object `name` in the `espresso_system`.
1987
1988        Args:
1989            name(`str`): Label of the object.
1990            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1991
1992        Note:
1993            - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"].
1994        '''
1995        print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
1996        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1997        if pmb_type != 'protein':
1998            raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein')
1999        molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list()
2000        for molecule_id in molecule_ids_list:    
2001            particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list()
2002            center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system)
2003            rigid_object_center = espresso_system.part.add(pos=center_of_mass,
2004                                                           rotation=[True,True,True], 
2005                                                           type=self.propose_unused_type())
2006            
2007            rigid_object_center.mass = len(particle_ids_list)
2008            momI = 0
2009            for pid in particle_ids_list:
2010                momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
2011
2012            rigid_object_center.rinertia = np.ones(3) * momI
2013            
2014            for particle_id in particle_ids_list:
2015                pid = espresso_system.part.by_id(particle_id)
2016                pid.vs_auto_relate_to(rigid_object_center.id)
2017        return
2018
2019    def filter_df(self, pmb_type):
2020        """
2021        Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns.
2022        
2023        Args:
2024            pmb_type(`str`): pmb_object_type to filter in `pmb.df`.
2025
2026        Returns:
2027            pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`.
2028        """
2029        pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type]
2030        pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1)
2031        return pmb_type_df
2032
2033    def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False):
2034        """
2035        Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2036        
2037        Args:
2038            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2039            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2040            use_default_bond(`bool`, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to 'False'. 
2041
2042        Returns:
2043            bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` if a matching bond exists
2044
2045        Note:
2046            - If `use_default_bond`=`True`, it returns "default" if no key is found.
2047        """
2048        bond_keys = [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']
2049        bond_defined=False
2050        for bond_key in bond_keys:
2051            if bond_key in self.df["name"].values:
2052                bond_defined=True
2053                correct_key=bond_key
2054                break
2055        if bond_defined:
2056            return correct_key
2057        elif use_default_bond:
2058            return 'default'
2059        else:
2060            return 
2061
2062    def find_value_from_es_type(self, es_type, column_name):
2063        """
2064        Finds a value in `pmb.df` for a `column_name` and `es_type` pair.
2065
2066        Args:
2067            es_type(`int`): value of the espresso type
2068            column_name(`str`): name of the column in `pymbe.df`
2069
2070        Returns:
2071            Value in `pymbe.df` matching  `column_name` and `es_type`
2072        """
2073        idx = pd.IndexSlice
2074        for state in ['state_one', 'state_two']:            
2075            index = self.df.loc[self.df[(state, 'es_type')] == es_type].index
2076            if len(index) > 0:
2077                if column_name == 'label':
2078                    label = self.df.loc[idx[index[0]], idx[(state,column_name)]]
2079                    return label
2080                else: 
2081                    column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]]
2082                    return column_name_value
2083        return None
2084
2085    def format_node(self, node_list):
2086        return "[" + " ".join(map(str, node_list)) + "]"
2087
2088
2089    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2090        """
2091        Generates coordinates outside a sphere centered at `center`.
2092
2093        Args:
2094            center(`lst`): Coordinates of the center of the sphere.
2095            radius(`float`): Radius of the sphere.
2096            max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates.
2097            n_samples(`int`): Number of sample points.
2098
2099        Returns:
2100            coord_list(`lst`): Coordinates of the sample points.
2101        """
2102        if not radius > 0: 
2103            raise ValueError (f'The value of {radius} must be a positive value')
2104        if not radius < max_dist:
2105            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
2106        coord_list = []
2107        counter = 0
2108        while counter<n_samples:
2109            coord = self.generate_random_points_in_a_sphere(center=center, 
2110                                            radius=max_dist,
2111                                            n_samples=1)[0]
2112            if np.linalg.norm(coord-np.asarray(center))>=radius:
2113                coord_list.append (coord)
2114                counter += 1
2115        return coord_list
2116    
2117    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2118        """
2119        Uniformly samples points from a hypersphere. If on_surface is set to True, the points are
2120        uniformly sampled from the surface of the hypersphere.
2121        
2122        Args:
2123            center(`lst`): Array with the coordinates of the center of the spheres.
2124            radius(`float`): Radius of the sphere.
2125            n_samples(`int`): Number of sample points to generate inside the sphere.
2126            on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere.
2127
2128        Returns:
2129            samples(`list`): Coordinates of the sample points inside the hypersphere.
2130        """
2131        # initial values
2132        center=np.array(center)
2133        d = center.shape[0]
2134        # sample n_samples points in d dimensions from a standard normal distribution
2135        samples = self.rng.normal(size=(n_samples, d))
2136        # make the samples lie on the surface of the unit hypersphere
2137        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
2138        samples /= normalize_radii
2139        if not on_surface:
2140            # make the samples lie inside the hypersphere with the correct density
2141            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
2142            new_radii = np.power(uniform_points, 1/d)
2143            samples *= new_radii
2144        # scale the points to have the correct radius and center
2145        samples = samples * radius + center
2146        return samples 
2147
2148    def generate_trial_perpendicular_vector(self,vector,magnitude):
2149        """
2150        Generates an orthogonal vector to the input `vector`.
2151
2152        Args:
2153            vector(`lst`): arbitrary vector.
2154            magnitude(`float`): magnitude of the orthogonal vector.
2155            
2156        Returns:
2157            (`lst`): Orthogonal vector with the same magnitude as the input vector.
2158        """ 
2159        np_vec = np.array(vector) 
2160        np_vec /= np.linalg.norm(np_vec) 
2161        if np.all(np_vec == 0):
2162            raise ValueError('Zero vector')
2163        # Generate a random vector 
2164        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
2165                                                                center=[0,0,0],
2166                                                                n_samples=1, 
2167                                                                on_surface=True)[0]
2168        # Project the random vector onto the input vector and subtract the projection
2169        projection = np.dot(random_vector, np_vec) * np_vec
2170        perpendicular_vector = random_vector - projection
2171        # Normalize the perpendicular vector to have the same magnitude as the input vector
2172        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
2173        return perpendicular_vector*magnitude
2174
2175    def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2176        """
2177        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length.
2178        If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead.
2179        If no bond is found, it prints a message and it does not return anything. If `hard_check` is activated, the code stops if no bond is found.
2180
2181        Args:
2182            particle_name1(str): label of the type of the first particle type of the bonded particles.
2183            particle_name2(str): label of the type of the second particle type of the bonded particles.
2184            hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2185            use_default_bond(bool, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to False. 
2186
2187        Returns:
2188            l0(`pint.Quantity`): bond length
2189        
2190        Note:
2191            - If `use_default_bond`=True and no bond is defined between `particle_name1` and `particle_name2`, it returns the default bond defined in `pmb.df`.
2192            - If `hard_check`=`True` stops the code when no bond is found.
2193        """
2194        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2195                                    particle_name2=particle_name2, 
2196                                    use_default_bond=use_default_bond)
2197        if bond_key:
2198            return self.df[self.df['name']==bond_key].l0.values[0]
2199        else:
2200            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2201            if hard_check:
2202                sys.exit(1)
2203            else:
2204                return
2205
2206    def get_charge_number_map(self):
2207        '''
2208        Gets the charge number of each `espresso_type` in `pymbe.df`.
2209        
2210        Returns:
2211            charge_number_map(`dict`): {espresso_type: z}.
2212        '''
2213        if self.df.state_one['es_type'].isnull().values.any():         
2214            df_state_one = self.df.state_one.dropna()     
2215            df_state_two = self.df.state_two.dropna()  
2216        else:    
2217            df_state_one = self.df.state_one
2218            if self.df.state_two['es_type'].isnull().values.any():
2219                df_state_two = self.df.state_two.dropna()   
2220            else:
2221                df_state_two = self.df.state_two
2222        state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values)
2223        state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values)
2224        charge_number_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2225        return charge_number_map
2226
2227    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2228        """
2229        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2230        `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule.
2231
2232        Args:
2233            particle_name1 (str): label of the type of the first particle type
2234            particle_name2 (str): label of the type of the second particle type
2235            combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2236
2237        Returns:
2238            {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2239
2240        Note:
2241            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
2242            - 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.
2243        """
2244        supported_combining_rules=["Lorentz-Berthelot"]
2245        lj_parameters_keys=["sigma","epsilon","offset","cutoff"]
2246        if combining_rule not in supported_combining_rules:
2247            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2248        lj_parameters={}
2249        for key in lj_parameters_keys:
2250            lj_parameters[key]=[]
2251        # Search the LJ parameters of the type pair
2252        for name in [particle_name1,particle_name2]:
2253            for key in lj_parameters_keys:
2254                lj_parameters[key].append(self.df[self.df.name == name][key].values[0])
2255        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2256        if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]):
2257            return {}
2258        # Apply combining rule
2259        if combining_rule == 'Lorentz-Berthelot':
2260            lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2
2261            lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2
2262            lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2
2263            lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1])
2264        return lj_parameters
2265
2266    def get_metal_ions_charge_number_map(self):
2267        """
2268        Gets a map with the charge numbers of all the metal ions supported.
2269
2270        Returns:
2271            metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number}
2272
2273        """
2274        metal_charge_number_map = {"Ca": 2}
2275        return metal_charge_number_map
2276
2277    def get_particle_id_map(self, object_name):
2278        '''
2279        Gets all the ids associated with the object with name `object_name` in `pmb.df`
2280
2281        Args:
2282            object_name(`str`): name of the object
2283        
2284        Returns:
2285            id_map(`dict`): dict of the structure {"all": [all_ids_with_object_name], "residue_map": {res_id: [particle_ids_in_res_id]}, "molecule_map": {mol_id: [particle_ids_in_mol_id]}, }
2286        '''
2287        object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0]
2288        valid_types = ["particle", "molecule", "residue", "protein"]
2289        if object_type not in valid_types:
2290            raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}")
2291        id_list = []
2292        mol_map = {}
2293        res_map = {}
2294        def do_res_map(res_ids):
2295            for res_id in res_ids:
2296                res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2297                res_map[res_id]=res_list
2298            return res_map
2299        if object_type in ['molecule', 'protein']:
2300            mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist()
2301            for mol_id in mol_ids:
2302                res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist())
2303                res_map=do_res_map(res_ids=res_ids)    
2304                mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2305                id_list+=mol_list
2306                mol_map[mol_id]=mol_list
2307        elif object_type == 'residue':     
2308            res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist()
2309            res_map=do_res_map(res_ids=res_ids)
2310            id_list=[]
2311            for res_id_list in res_map.values():
2312                id_list+=res_id_list
2313        elif object_type == 'particle':
2314            id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist()
2315        return {"all": id_list, "molecule_map": mol_map, "residue_map": res_map}
2316
2317    def get_pka_set(self):
2318        '''
2319        Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df`
2320        
2321        Returns:
2322            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
2323        '''
2324        titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna()
2325        pka_set = {}
2326        for index in titratables_AA_df.name.keys():
2327            name = titratables_AA_df.name[index]
2328            pka_value = titratables_AA_df.pka[index]
2329            acidity = titratables_AA_df.acidity[index]   
2330            pka_set[name] = {'pka_value':pka_value,'acidity':acidity}
2331        return pka_set 
2332    
2333    def get_radius_map(self, dimensionless=True):
2334        '''
2335        Gets the effective radius of each `espresso type` in `pmb.df`. 
2336
2337        Args:
2338            dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False.
2339        
2340        Returns:
2341            radius_map(`dict`): {espresso_type: radius}.
2342
2343        Note:
2344            The radius corresponds to (sigma+offset)/2
2345        '''
2346        df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates()
2347        df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates()
2348        state_one = pd.Series((df_state_one.sigma.values+df_state_one.offset.values)/2.0,index=df_state_one.state_one.es_type.values)
2349        state_two = pd.Series((df_state_two.sigma.values+df_state_two.offset.values)/2.0,index=df_state_two.state_two.es_type.values)
2350        radius_map  = pd.concat([state_one,state_two],axis=0).to_dict()  
2351        if dimensionless:
2352            for key in radius_map:
2353                radius_map[key] = radius_map[key].magnitude
2354        return radius_map
2355
2356    def get_reduced_units(self):
2357        """
2358        Returns the  current set of reduced units defined in pyMBE.
2359
2360        Returns:
2361            reduced_units_text(`str`): text with information about the current set of reduced units.
2362
2363        """
2364        unit_length=self.units.Quantity(1,'reduced_length')
2365        unit_energy=self.units.Quantity(1,'reduced_energy')
2366        unit_charge=self.units.Quantity(1,'reduced_charge')
2367        reduced_units_text = "\n".join(["Current set of reduced units:",
2368                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2369                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2370                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2371                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2372                                        ])   
2373        return reduced_units_text
2374
2375    def get_resource(self, path):
2376        '''
2377        Locate a file resource of the pyMBE package.
2378
2379        Args:
2380            path(`str`): Relative path to the resource
2381
2382        Returns:
2383            path(`str`): Absolute path to the resource
2384
2385        '''
2386        import os
2387        return os.path.join(os.path.dirname(__file__), path)
2388
2389    def get_type_map(self):
2390        """
2391        Gets all different espresso types assigned to particles  in `pmb.df`.
2392        
2393        Returns:
2394            type_map(`dict`): {"name": espresso_type}.
2395        """
2396        df_state_one = self.df.state_one.dropna(how='all')     
2397        df_state_two = self.df.state_two.dropna(how='all')  
2398        state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label)
2399        state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label)
2400        type_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2401        return type_map
2402
2403    def initialize_lattice_builder(self, diamond_lattice):
2404        """
2405        Initialize the lattice builder with the DiamondLattice object.
2406
2407        Args:
2408            diamond_lattice(`DiamondLattice`): DiamondLattice object from the `lib/lattice` module to be used in the LatticeBuilder.
2409        """
2410        if not isinstance(diamond_lattice, DiamondLattice):
2411            raise TypeError("Currently only DiamondLattice objects are supported.")
2412        self.lattice_builder = LatticeBuilder(lattice=diamond_lattice)
2413        logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}")
2414        return self.lattice_builder
2415
2416
2417    def load_interaction_parameters(self, filename, overwrite=False):
2418        """
2419        Loads the interaction parameters stored in `filename` into `pmb.df`
2420        
2421        Args:
2422            filename(`str`): name of the file to be read
2423            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2424        """
2425        without_units = ['z','es_type']
2426        with_units = ['sigma','epsilon','offset','cutoff']
2427
2428        with open(filename, 'r') as f:
2429            interaction_data = json.load(f)
2430            interaction_parameter_set = interaction_data["data"]
2431
2432        for key in interaction_parameter_set:
2433            param_dict=interaction_parameter_set[key]
2434            object_type=param_dict.pop('object_type')
2435            if object_type == 'particle':
2436                not_required_attributes={}    
2437                for not_required_key in without_units+with_units:
2438                    if not_required_key in param_dict.keys():
2439                        if not_required_key in with_units:
2440                            not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key))
2441                        elif not_required_key in without_units:
2442                            not_required_attributes[not_required_key]=param_dict.pop(not_required_key)
2443                    else:
2444                        not_required_attributes[not_required_key]=pd.NA
2445                self.define_particle(name=param_dict.pop('name'),
2446                                z=not_required_attributes.pop('z'),
2447                                sigma=not_required_attributes.pop('sigma'),
2448                                offset=not_required_attributes.pop('offset'),
2449                                cutoff=not_required_attributes.pop('cutoff'),
2450                                epsilon=not_required_attributes.pop('epsilon'),
2451                                overwrite=overwrite)
2452            elif object_type == 'residue':
2453                self.define_residue(**param_dict)
2454            elif object_type == 'molecule':
2455                self.define_molecule(**param_dict)
2456            elif object_type == 'peptide':
2457                self.define_peptide(**param_dict)
2458            elif object_type == 'bond':
2459                particle_pairs = param_dict.pop('particle_pairs')
2460                bond_parameters = param_dict.pop('bond_parameters')
2461                bond_type = param_dict.pop('bond_type')
2462                if bond_type == 'harmonic':
2463                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2464                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2465                    bond = {'r_0'    : r_0,
2466                            'k'      : k,
2467                            }
2468
2469                elif bond_type == 'FENE':
2470                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2471                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2472                    d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max'))
2473                    bond =  {'r_0'    : r_0,
2474                             'k'      : k,
2475                             'd_r_max': d_r_max,
2476                             }
2477                else:
2478                    raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds")
2479                self.define_bond(bond_type=bond_type, 
2480                                bond_parameters=bond, 
2481                                particle_pairs=particle_pairs)
2482            else:
2483                raise ValueError(object_type+' is not a known pmb object type')
2484            
2485        return
2486    
2487    def load_pka_set(self, filename, overwrite=True):
2488        """
2489        Loads the pka_set stored in `filename` into `pmb.df`.
2490        
2491        Args:
2492            filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}.
2493            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 
2494        """
2495        with open(filename, 'r') as f:
2496            pka_data = json.load(f)
2497            pka_set = pka_data["data"]
2498
2499        self.check_pka_set(pka_set=pka_set)
2500
2501        for key in pka_set:
2502            acidity = pka_set[key]['acidity']
2503            pka_value = pka_set[key]['pka_value']
2504            self.set_particle_acidity(name=key, 
2505                                      acidity=acidity, 
2506                                      pka=pka_value, 
2507                                      overwrite=overwrite)
2508        return
2509
2510
2511    def propose_unused_type(self):
2512        """
2513        Searches in `pmb.df` all the different particle types defined and returns a new one.
2514
2515        Returns:
2516            unused_type(`int`): unused particle type
2517        """
2518        type_map = self.get_type_map()
2519        if not type_map:  
2520            unused_type = 0
2521        else:
2522            valid_values = [v for v in type_map.values() if pd.notna(v)]  # Filter out pd.NA values
2523            unused_type = max(valid_values) + 1 if valid_values else 0  # Ensure max() doesn't fail if all values are NA
2524        return unused_type
2525
2526    def protein_sequence_parser(self, sequence):
2527        '''
2528        Parses `sequence` to the one letter code for amino acids.
2529        
2530        Args:
2531            sequence(`str` or `lst`): Sequence of the amino acid. 
2532
2533        Returns:
2534            clean_sequence(`lst`): `sequence` using the one letter code.
2535        
2536        Note:
2537            - Accepted formats for `sequence` are:
2538                - `lst` with one letter or three letter code of each aminoacid in each element
2539                - `str` with the sequence using the one letter code
2540                - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-"
2541        
2542        '''
2543        # Aminoacid key
2544        keys={"ALA": "A",
2545                "ARG": "R",
2546                "ASN": "N",
2547                "ASP": "D",
2548                "CYS": "C",
2549                "GLU": "E",
2550                "GLN": "Q",
2551                "GLY": "G",
2552                "HIS": "H",
2553                "ILE": "I",
2554                "LEU": "L",
2555                "LYS": "K",
2556                "MET": "M",
2557                "PHE": "F",
2558                "PRO": "P",
2559                "SER": "S",
2560                "THR": "T",
2561                "TRP": "W",
2562                "TYR": "Y",
2563                "VAL": "V",
2564                "PSER": "J",
2565                "PTHR": "U",
2566                "PTyr": "Z",
2567                "NH2": "n",
2568                "COOH": "c"}
2569        clean_sequence=[]
2570        if isinstance(sequence, str):
2571            if sequence.find("-") != -1:
2572                splited_sequence=sequence.split("-")
2573                for residue in splited_sequence:
2574                    if len(residue) == 1:
2575                        if residue in keys.values():
2576                            residue_ok=residue
2577                        else:
2578                            if residue.upper() in keys.values():
2579                                residue_ok=residue.upper()
2580                            else:
2581                                raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence")
2582                        clean_sequence.append(residue_ok)
2583                    else:
2584                        if residue in keys.keys():
2585                            clean_sequence.append(keys[residue])
2586                        else:
2587                            if residue.upper() in keys.keys():
2588                                clean_sequence.append(keys[residue.upper()])
2589                            else:
2590                                raise ValueError("Unknown  code for a residue: ", residue, " please review the input sequence")
2591            else:
2592                for residue in sequence:
2593                    if residue in keys.values():
2594                        residue_ok=residue
2595                    else:
2596                        if residue.upper() in keys.values():
2597                            residue_ok=residue.upper()
2598                        else:
2599                            raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence")
2600                    clean_sequence.append(residue_ok)
2601        if isinstance(sequence, list):
2602            for residue in sequence:
2603                if residue in keys.values():
2604                    residue_ok=residue
2605                else:
2606                    if residue.upper() in keys.values():
2607                        residue_ok=residue.upper()
2608                    elif (residue.upper() in keys.keys()):
2609                        residue_ok= keys[residue.upper()]
2610                    else:
2611                        raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence")
2612                clean_sequence.append(residue_ok)
2613        return clean_sequence
2614    
2615    def read_pmb_df (self,filename):
2616        """
2617        Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE.
2618
2619        Args:
2620            filename(`str`): path to file.
2621
2622        Note:
2623            This function only accepts files with CSV format. 
2624        """
2625        
2626        if filename.rsplit(".", 1)[1] != "csv":
2627            raise ValueError("Only files with CSV format are supported")
2628        df = pd.read_csv (filename,header=[0, 1], index_col=0)
2629        columns_names = self.setup_df()
2630        
2631        multi_index = pd.MultiIndex.from_tuples(columns_names)
2632        df.columns = multi_index
2633        
2634        self.df = self.convert_columns_to_original_format(df)
2635        self.df.fillna(pd.NA, inplace=True)
2636        return self.df
2637    
2638    def read_protein_vtf_in_df (self,filename,unit_length=None):
2639        """
2640        Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`.
2641
2642        Args:
2643            filename(`str`): Path to the vtf file with the coarse-grained model.
2644            unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None.
2645
2646        Returns:
2647            topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}
2648
2649        Note:
2650            - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom.
2651        """
2652
2653        print (f'Loading protein coarse grain model file: {filename}')
2654
2655        coord_list = []
2656        particles_dict = {}
2657
2658        if unit_length is None:
2659            unit_length = 1 * self.units.angstrom 
2660
2661        with open (filename,'r') as protein_model:
2662            for line in protein_model :
2663                line_split = line.split()
2664                if line_split:
2665                    line_header = line_split[0]
2666                    if line_header == 'atom':
2667                        atom_id  = line_split[1]
2668                        atom_name = line_split[3]
2669                        atom_resname = line_split[5]
2670                        chain_id = line_split[9]
2671                        radius = float(line_split [11])*unit_length 
2672                        particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius]
2673                    elif line_header.isnumeric(): 
2674                        atom_coord = line_split[1:] 
2675                        atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord]
2676                        coord_list.append (atom_coord)
2677
2678        numbered_label = []
2679        i = 0   
2680        
2681        for atom_id in particles_dict.keys():
2682    
2683            if atom_id == 1:
2684                atom_name = particles_dict[atom_id][0]
2685                numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2686                numbered_label.append(numbered_name)
2687
2688            elif atom_id != 1: 
2689            
2690                if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]:
2691                    i += 1                    
2692                    count = 1
2693                    atom_name = particles_dict[atom_id][0]
2694                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2695                    numbered_label.append(numbered_name)
2696                    
2697                elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]:
2698                    if count == 2 or particles_dict[atom_id][1] == 'GLY':
2699                        i +=1  
2700                        count = 0
2701                    atom_name = particles_dict[atom_id][0]
2702                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2703                    numbered_label.append(numbered_name)
2704                    count +=1
2705
2706        topology_dict = {}
2707
2708        for i in range (0, len(numbered_label)):   
2709            topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] ,
2710                                                    'chain_id':numbered_label[i][1],
2711                                                    'radius':numbered_label[i][2] }
2712
2713        return topology_dict
2714
2715    def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2716        """
2717        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2718        If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead.
2719        If no bond is found, it prints a message and it does not return anything. If `hard_check` is activated, the code stops if no bond is found.
2720
2721        Args:
2722            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2723            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2724            hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2725            use_default_bond(`bool`, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to False. 
2726
2727        Returns:
2728            bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library.
2729        
2730        Note:
2731            - If `use_default_bond`=True and no bond is defined between `particle_name1` and `particle_name2`, it returns the default bond defined in `pmb.df`.
2732            - If `hard_check`=`True` stops the code when no bond is found.
2733        """
2734        
2735        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2736                                    particle_name2=particle_name2, 
2737                                    use_default_bond=use_default_bond)
2738        if use_default_bond:
2739            if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'):
2740                raise ValueError(f"use_default_bond is set to {use_default_bond} but no default bond has been defined. Please define a default bond with pmb.define_default_bond")
2741        if bond_key:
2742            return self.df[self.df['name']==bond_key].bond_object.values[0]
2743        else:
2744            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2745            if hard_check:
2746                sys.exit(1)
2747            else:
2748                return
2749
2750    def search_particles_in_residue(self, residue_name):
2751        '''
2752        Searches for all particles in a given residue of name `residue_name`.
2753
2754        Args:
2755            residue_name (`str`): name of the residue to be searched
2756
2757        Returns:
2758            list_of_particles_in_residue (`lst`): list of the names of all particles in the residue
2759
2760        Note:
2761            - The function returns a name per particle in residue, i.e. if there are multiple particles with the same type `list_of_particles_in_residue` will have repeated items.
2762 
2763        '''
2764        index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 
2765        central_bead = self.df.at [index_residue, ('central_bead', '')]
2766        list_of_side_chains = self.df.at [index_residue, ('side_chains', '')]
2767
2768        list_of_particles_in_residue = []
2769        list_of_particles_in_residue.append(central_bead)
2770
2771        for side_chain in list_of_side_chains: 
2772            object_type = self.df[self.df['name']==side_chain].pmb_type.values[0]
2773
2774            if object_type == "residue":
2775                list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain)
2776                list_of_particles_in_residue += list_of_particles_in_side_chain_residue
2777            elif object_type == "particle":
2778                list_of_particles_in_residue.append(side_chain)
2779
2780        return list_of_particles_in_residue
2781
2782         
2783
2784    def set_particle_acidity(self, name, acidity=pd.NA, default_charge_number=0, pka=pd.NA, overwrite=True):
2785        """
2786        Sets the particle acidity including the charges in each of its possible states. 
2787
2788        Args:
2789            name(`str`): Unique label that identifies the `particle`. 
2790            acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None.
2791            default_charge_number (`int`): Charge number of the particle. Defaults to 0.
2792            pka(`float`, optional):  If `particle` is an acid or a base, it defines its pka-value. Defaults to pandas.NA.
2793            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2794     
2795        Note:
2796            - For particles with  `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 
2797        deprotonated states, respectively. 
2798            - For particles without an acidity `acidity = pandas.NA`, only `state_one` is defined.
2799            - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 
2800        """
2801        acidity_valid_keys = ['inert','acidic', 'basic']
2802        if not pd.isna(acidity):
2803            if acidity not in acidity_valid_keys:
2804                raise ValueError(f"Acidity {acidity} provided for particle name  {name} is not supproted. Valid keys are: {acidity_valid_keys}")
2805            if acidity in ['acidic', 'basic'] and pd.isna(pka):
2806                raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.")   
2807            if acidity == "inert":
2808                acidity = pd.NA
2809                print("Deprecation warning: the keyword 'inert' for acidity has been replaced by setting acidity = pd.NA. For backwards compatibility, acidity has been set to pd.NA. Support for `acidity = 'inert'` may be deprecated in future releases of pyMBE")
2810
2811        self.define_particle_entry_in_df(name=name)
2812        
2813        for index in self.df[self.df['name']==name].index:       
2814            if pka is not pd.NA:
2815                self.add_value_to_df(key=('pka',''),
2816                                    index=index,
2817                                    new_value=pka, 
2818                                    overwrite=overwrite)
2819            
2820            self.add_value_to_df(key=('acidity',''),
2821                                 index=index,
2822                                 new_value=acidity, 
2823                                 overwrite=overwrite) 
2824            if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')):
2825                self.add_value_to_df(key=('state_one','es_type'),
2826                                     index=index,
2827                                     new_value=self.propose_unused_type(), 
2828                                     overwrite=overwrite)  
2829            if pd.isna(self.df.loc [self.df['name']  == name].acidity.iloc[0]):
2830                self.add_value_to_df(key=('state_one','z'),
2831                                     index=index,
2832                                     new_value=default_charge_number, 
2833                                     overwrite=overwrite)
2834                self.add_value_to_df(key=('state_one','label'),
2835                                     index=index,
2836                                     new_value=name, 
2837                                    overwrite=overwrite)
2838            else:
2839                protonated_label = f'{name}H'
2840                self.add_value_to_df(key=('state_one','label'),
2841                                     index=index,
2842                                     new_value=protonated_label, 
2843                                    overwrite=overwrite)
2844                self.add_value_to_df(key=('state_two','label'),
2845                                     index=index,
2846                                     new_value=name, 
2847                                    overwrite=overwrite)
2848                if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')):
2849                    self.add_value_to_df(key=('state_two','es_type'),
2850                                         index=index,
2851                                         new_value=self.propose_unused_type(), 
2852                                         overwrite=overwrite)
2853                if self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'acidic':        
2854                    self.add_value_to_df(key=('state_one','z'),
2855                                         index=index,new_value=0, 
2856                                         overwrite=overwrite)
2857                    self.add_value_to_df(key=('state_two','z'),
2858                                         index=index,
2859                                         new_value=-1, 
2860                                         overwrite=overwrite)
2861                elif self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'basic':
2862                    self.add_value_to_df(key=('state_one','z'),
2863                                         index=index,new_value=+1, 
2864                                         overwrite=overwrite)
2865                    self.add_value_to_df(key=('state_two','z'),
2866                                         index=index,
2867                                         new_value=0, 
2868                                         overwrite=overwrite)   
2869        self.df.fillna(pd.NA, inplace=True)
2870        return
2871    
2872    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True):
2873        """
2874        Sets the set of reduced units used by pyMBE.units and it prints it.
2875
2876        Args:
2877            unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 
2878            unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
2879            temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 
2880            Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
2881            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
2882
2883        Note:
2884            - If no `temperature` is given, a value of 298.15 K is assumed by default.
2885            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
2886            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
2887            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2888        """
2889        if unit_length is None:
2890            unit_length= 0.355*self.units.nm
2891        if temperature is None:
2892            temperature = 298.15 * self.units.K
2893        if unit_charge is None:
2894            unit_charge = scipy.constants.e * self.units.C
2895        if Kw is None:
2896            Kw = 1e-14
2897        # Sanity check
2898        variables=[unit_length,temperature,unit_charge]
2899        dimensionalities=["[length]","[temperature]","[charge]"]
2900        for variable,dimensionality in zip(variables,dimensionalities):
2901            self.check_dimensionality(variable,dimensionality)
2902        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2903        self.kT=temperature*self.kB
2904        self.units._build_cache()
2905        self.units.define(f'reduced_energy = {self.kT} ')
2906        self.units.define(f'reduced_length = {unit_length}')
2907        self.units.define(f'reduced_charge = {unit_charge}')
2908        if verbose:        
2909            print(self.get_reduced_units())
2910        return
2911
2912    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
2913        """
2914        Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 
2915
2916        Args:
2917            counter_ion(`str`): `name` of the counter_ion `particle`.
2918            constant_pH(`float`): pH-value.
2919            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
2920            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2921            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
2922
2923        Returns:
2924            RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2925            sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
2926        """
2927        from espressomd import reaction_methods
2928        if exclusion_range is None:
2929            exclusion_range = max(self.get_radius_map().values())*2.0
2930        if pka_set is None:
2931            pka_set=self.get_pka_set()    
2932        self.check_pka_set(pka_set=pka_set)
2933        if use_exclusion_radius_per_type:
2934            exclusion_radius_per_type = self.get_radius_map()
2935        else:
2936            exclusion_radius_per_type = {}
2937        
2938        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2939                                                    exclusion_range=exclusion_range, 
2940                                                    seed=self.seed, 
2941                                                    constant_pH=constant_pH,
2942                                                    exclusion_radius_per_type = exclusion_radius_per_type
2943                                                    )
2944        sucessfull_reactions_labels=[]
2945        charge_number_map = self.get_charge_number_map()
2946        for name in pka_set.keys():
2947            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
2948                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
2949                continue
2950            gamma=10**-pka_set[name]['pka_value']
2951            state_one_type   = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
2952            state_two_type   = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
2953            counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0]
2954            RE.add_reaction(gamma=gamma,
2955                            reactant_types=[state_one_type],
2956                            product_types=[state_two_type, counter_ion_type],
2957                            default_charges={state_one_type: charge_number_map[state_one_type],
2958                            state_two_type: charge_number_map[state_two_type],
2959                            counter_ion_type: charge_number_map[counter_ion_type]})
2960            sucessfull_reactions_labels.append(name)
2961        return RE, sucessfull_reactions_labels
2962
2963    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2964        """
2965        Sets up grand-canonical coupling to a reservoir of salt.
2966        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2967
2968        Args:
2969            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2970            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
2971            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
2972            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
2973            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
2974            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2975
2976        Returns:
2977            RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
2978        """
2979        from espressomd import reaction_methods
2980        if exclusion_range is None:
2981            exclusion_range = max(self.get_radius_map().values())*2.0
2982        if use_exclusion_radius_per_type:
2983            exclusion_radius_per_type = self.get_radius_map()
2984        else:
2985            exclusion_radius_per_type = {}
2986        
2987        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2988                                                    exclusion_range=exclusion_range, 
2989                                                    seed=self.seed, 
2990                                                    exclusion_radius_per_type = exclusion_radius_per_type
2991                                                    )
2992
2993        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
2994        determined_activity_coefficient = activity_coefficient(c_salt_res)
2995        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
2996
2997        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
2998        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
2999
3000        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
3001        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
3002
3003        if salt_cation_charge <= 0:
3004            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3005        if salt_anion_charge >= 0:
3006            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3007
3008        # Grand-canonical coupling to the reservoir
3009        RE.add_reaction(
3010            gamma = K_salt.magnitude,
3011            reactant_types = [],
3012            reactant_coefficients = [],
3013            product_types = [ salt_cation_es_type, salt_anion_es_type ],
3014            product_coefficients = [ 1, 1 ],
3015            default_charges = {
3016                salt_cation_es_type: salt_cation_charge, 
3017                salt_anion_es_type: salt_anion_charge, 
3018            }
3019        )
3020
3021        return RE
3022
3023    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, pka_set=None, use_exclusion_radius_per_type = False):
3024        """
3025        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3026        reservoir of small ions. 
3027        This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
3028
3029        [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3030
3031        Args:
3032            pH_res ('float): pH-value in the reservoir.
3033            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3034            proton_name ('str'): Name of the proton (H+) particle.
3035            hydroxide_name ('str'): Name of the hydroxide (OH-) particle.
3036            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
3037            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
3038            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3039            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
3040            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3041            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
3042
3043        Returns:
3044            RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3045            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3046            ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3047        """
3048        from espressomd import reaction_methods
3049        if exclusion_range is None:
3050            exclusion_range = max(self.get_radius_map().values())*2.0
3051        if pka_set is None:
3052            pka_set=self.get_pka_set()    
3053        self.check_pka_set(pka_set=pka_set)
3054        if use_exclusion_radius_per_type:
3055            exclusion_radius_per_type = self.get_radius_map()
3056        else:
3057            exclusion_radius_per_type = {}
3058        
3059        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3060                                                    exclusion_range=exclusion_range, 
3061                                                    seed=self.seed, 
3062                                                    exclusion_radius_per_type = exclusion_radius_per_type
3063                                                    )
3064
3065        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3066        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3067        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3068        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3069        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3070        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3071        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3072
3073        proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0]
3074        hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0]     
3075        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
3076        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
3077
3078        proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0]
3079        hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0]     
3080        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
3081        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
3082
3083        if proton_charge <= 0:
3084            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
3085        if salt_cation_charge <= 0:
3086            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3087        if hydroxide_charge >= 0:
3088            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
3089        if salt_anion_charge >= 0:
3090            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3091
3092        # Grand-canonical coupling to the reservoir
3093        # 0 = H+ + OH-
3094        RE.add_reaction(
3095            gamma = K_W.magnitude,
3096            reactant_types = [],
3097            reactant_coefficients = [],
3098            product_types = [ proton_es_type, hydroxide_es_type ],
3099            product_coefficients = [ 1, 1 ],
3100            default_charges = {
3101                proton_es_type: proton_charge, 
3102                hydroxide_es_type: hydroxide_charge, 
3103            }
3104        )
3105
3106        # 0 = Na+ + Cl-
3107        RE.add_reaction(
3108            gamma = K_NACL.magnitude,
3109            reactant_types = [],
3110            reactant_coefficients = [],
3111            product_types = [ salt_cation_es_type, salt_anion_es_type ],
3112            product_coefficients = [ 1, 1 ],
3113            default_charges = {
3114                salt_cation_es_type: salt_cation_charge, 
3115                salt_anion_es_type: salt_anion_charge, 
3116            }
3117        )
3118
3119        # 0 = Na+ + OH-
3120        RE.add_reaction(
3121            gamma = (K_NACL * K_W / K_HCL).magnitude,
3122            reactant_types = [],
3123            reactant_coefficients = [],
3124            product_types = [ salt_cation_es_type, hydroxide_es_type ],
3125            product_coefficients = [ 1, 1 ],
3126            default_charges = {
3127                salt_cation_es_type: salt_cation_charge, 
3128                hydroxide_es_type: hydroxide_charge, 
3129            }
3130        )
3131
3132        # 0 = H+ + Cl-
3133        RE.add_reaction(
3134            gamma = K_HCL.magnitude,
3135            reactant_types = [],
3136            reactant_coefficients = [],
3137            product_types = [ proton_es_type, salt_anion_es_type ],
3138            product_coefficients = [ 1, 1 ],
3139            default_charges = {
3140                proton_es_type: proton_charge, 
3141                salt_anion_es_type: salt_anion_charge, 
3142            }
3143        )
3144
3145        # Annealing moves to ensure sufficient sampling
3146        # Cation annealing H+ = Na+
3147        RE.add_reaction(
3148            gamma = (K_NACL / K_HCL).magnitude,
3149            reactant_types = [proton_es_type],
3150            reactant_coefficients = [ 1 ],
3151            product_types = [ salt_cation_es_type ],
3152            product_coefficients = [ 1 ],
3153            default_charges = {
3154                proton_es_type: proton_charge, 
3155                salt_cation_es_type: salt_cation_charge, 
3156            }
3157        )
3158
3159        # Anion annealing OH- = Cl- 
3160        RE.add_reaction(
3161            gamma = (K_HCL / K_W).magnitude,
3162            reactant_types = [hydroxide_es_type],
3163            reactant_coefficients = [ 1 ],
3164            product_types = [ salt_anion_es_type ],
3165            product_coefficients = [ 1 ],
3166            default_charges = {
3167                hydroxide_es_type: hydroxide_charge, 
3168                salt_anion_es_type: salt_anion_charge, 
3169            }
3170        )
3171
3172        sucessful_reactions_labels=[]
3173        charge_number_map = self.get_charge_number_map()
3174        for name in pka_set.keys():
3175            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3176                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3177                continue
3178
3179            Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3180
3181            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3182            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3183
3184            # Reaction in terms of proton: HA = A + H+
3185            RE.add_reaction(gamma=Ka.magnitude,
3186                            reactant_types=[state_one_type],
3187                            reactant_coefficients=[1],
3188                            product_types=[state_two_type, proton_es_type],
3189                            product_coefficients=[1, 1],
3190                            default_charges={state_one_type: charge_number_map[state_one_type],
3191                            state_two_type: charge_number_map[state_two_type],
3192                            proton_es_type: proton_charge})
3193
3194            # Reaction in terms of salt cation: HA = A + Na+
3195            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
3196                            reactant_types=[state_one_type],
3197                            reactant_coefficients=[1],
3198                            product_types=[state_two_type, salt_cation_es_type],
3199                            product_coefficients=[1, 1],
3200                            default_charges={state_one_type: charge_number_map[state_one_type],
3201                            state_two_type: charge_number_map[state_two_type],
3202                            salt_cation_es_type: salt_cation_charge})
3203
3204            # Reaction in terms of hydroxide: OH- + HA = A
3205            RE.add_reaction(gamma=(Ka / K_W).magnitude,
3206                            reactant_types=[state_one_type, hydroxide_es_type],
3207                            reactant_coefficients=[1, 1],
3208                            product_types=[state_two_type],
3209                            product_coefficients=[1],
3210                            default_charges={state_one_type: charge_number_map[state_one_type],
3211                            state_two_type: charge_number_map[state_two_type],
3212                            hydroxide_es_type: hydroxide_charge})
3213
3214            # Reaction in terms of salt anion: Cl- + HA = A
3215            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
3216                            reactant_types=[state_one_type, salt_anion_es_type],
3217                            reactant_coefficients=[1, 1],
3218                            product_types=[state_two_type],
3219                            product_coefficients=[1],
3220                            default_charges={state_one_type: charge_number_map[state_one_type],
3221                            state_two_type: charge_number_map[state_two_type],
3222                            salt_anion_es_type: salt_anion_charge})
3223
3224            sucessful_reactions_labels.append(name)
3225        return RE, sucessful_reactions_labels, ionic_strength_res
3226
3227    def setup_grxmc_unified(self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
3228        """
3229        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3230        reservoir of small ions. 
3231        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-}. 
3232        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'.
3233
3234        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3235        [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3236
3237        Args:
3238            pH_res ('float'): pH-value in the reservoir.
3239            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3240            cation_name ('str'): Name of the cationic particle.
3241            anion_name ('str'): Name of the anionic particle.
3242            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3243            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
3244            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3245            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`.
3246
3247        Returns:
3248            RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3249            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3250            ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3251        """
3252        from espressomd import reaction_methods
3253        if exclusion_range is None:
3254            exclusion_range = max(self.get_radius_map().values())*2.0
3255        if pka_set is None:
3256            pka_set=self.get_pka_set()    
3257        self.check_pka_set(pka_set=pka_set)
3258        if use_exclusion_radius_per_type:
3259            exclusion_radius_per_type = self.get_radius_map()
3260        else:
3261            exclusion_radius_per_type = {}
3262        
3263        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3264                                                    exclusion_range=exclusion_range, 
3265                                                    seed=self.seed, 
3266                                                    exclusion_radius_per_type = exclusion_radius_per_type
3267                                                    )
3268
3269        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3270        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3271        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3272        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3273        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3274        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3275        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3276        K_XX = a_cation * a_anion
3277
3278        cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0]
3279        anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0]     
3280        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
3281        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
3282        if cation_charge <= 0:
3283            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3284        if anion_charge >= 0:
3285            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3286
3287        # Coupling to the reservoir: 0 = X+ + X-
3288        RE.add_reaction(
3289            gamma = K_XX.magnitude,
3290            reactant_types = [],
3291            reactant_coefficients = [],
3292            product_types = [ cation_es_type, anion_es_type ],
3293            product_coefficients = [ 1, 1 ],
3294            default_charges = {
3295                cation_es_type: cation_charge, 
3296                anion_es_type: anion_charge, 
3297            }
3298        )
3299
3300        sucessful_reactions_labels=[]
3301        charge_number_map = self.get_charge_number_map()
3302        for name in pka_set.keys():
3303            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3304                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3305                continue
3306
3307            Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l
3308            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3309
3310            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3311            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3312
3313            # Reaction in terms of small cation: HA = A + X+
3314            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3315                            reactant_types=[state_one_type],
3316                            reactant_coefficients=[1],
3317                            product_types=[state_two_type, cation_es_type],
3318                            product_coefficients=[1, 1],
3319                            default_charges={state_one_type: charge_number_map[state_one_type],
3320                            state_two_type: charge_number_map[state_two_type],
3321                            cation_es_type: cation_charge})
3322
3323            # Reaction in terms of small anion: X- + HA = A
3324            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3325                            reactant_types=[state_one_type, anion_es_type],
3326                            reactant_coefficients=[1, 1],
3327                            product_types=[state_two_type],
3328                            product_coefficients=[1],
3329                            default_charges={state_one_type: charge_number_map[state_one_type],
3330                            state_two_type: charge_number_map[state_two_type],
3331                            anion_es_type: anion_charge})
3332
3333            sucessful_reactions_labels.append(name)
3334        return RE, sucessful_reactions_labels, ionic_strength_res
3335
3336    def setup_df (self):
3337        """
3338        Sets up the pyMBE's dataframe `pymbe.df`.
3339
3340        Returns:
3341            columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe
3342        """
3343        
3344        columns_dtypes = {
3345            'name': {
3346                '': str},
3347            'pmb_type': {
3348                '': str},
3349            'particle_id': {
3350                '': pd.Int64Dtype()},
3351            'particle_id2':  {
3352                '': pd.Int64Dtype()},
3353            'residue_id':  {
3354                '': pd.Int64Dtype()},
3355            'molecule_id':  {
3356                '': pd.Int64Dtype()},
3357            'acidity':  {
3358                '': str},
3359            'pka':  {
3360                '': object},
3361            'central_bead':  {
3362                '': object},
3363            'side_chains': {
3364                '': object},
3365            'residue_list': {
3366                '': object},
3367            'model': {
3368                '': str},
3369            'sigma': {
3370                '': object},
3371            'cutoff': {
3372                '': object},
3373            'offset': {
3374                '': object},
3375            'epsilon': {
3376                '': object},
3377            'state_one': {
3378                'label': str,
3379                'es_type': pd.Int64Dtype(),
3380                'z': pd.Int64Dtype()},
3381            'state_two': {
3382                'label': str,
3383                'es_type': pd.Int64Dtype(),
3384                'z': pd.Int64Dtype()},
3385            'sequence': {
3386                '': object},
3387            'bond_object': {
3388                '': object},
3389            'parameters_of_the_potential':{
3390                '': object},
3391            'l0': {
3392                '': float},
3393            'node_map':{
3394                '':object},
3395            'chain_map':{
3396                '':object}}
3397        
3398        self.df = pd.DataFrame(columns=pd.MultiIndex.from_tuples([(col_main, col_sub) for col_main, sub_cols in columns_dtypes.items() for col_sub in sub_cols.keys()]))
3399        
3400        for level1, sub_dtypes in columns_dtypes.items():
3401            for level2, dtype in sub_dtypes.items():
3402                self.df[level1, level2] = self.df[level1, level2].astype(dtype)
3403                
3404        columns_names = pd.MultiIndex.from_frame(self.df)
3405        columns_names = columns_names.names
3406                
3407        return columns_names
3408
3409    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True):
3410        """
3411        Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`.
3412
3413        Args:
3414            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
3415            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.
3416            combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3417            warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True.
3418
3419        Note:
3420            - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 
3421            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
3422            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3423
3424        """
3425        from itertools import combinations_with_replacement
3426        import warnings
3427        implemented_combining_rules = ['Lorentz-Berthelot']
3428        compulsory_parameters_in_df = ['sigma','epsilon']
3429        # Sanity check
3430        if combining_rule not in implemented_combining_rules:
3431            raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules)
3432        shift=0
3433        if shift_potential:
3434            shift="auto"
3435        # List which particles have sigma and epsilon values defined in pmb.df and which ones don't
3436        particles_types_with_LJ_parameters = []
3437        non_parametrized_labels= []
3438        for particle_type in self.get_type_map().values():
3439            check_list=[]
3440            for key in compulsory_parameters_in_df:
3441                value_in_df=self.find_value_from_es_type(es_type=particle_type,
3442                                                        column_name=key)
3443                check_list.append(pd.isna(value_in_df))
3444            if any(check_list):
3445                non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 
3446                                                                            column_name='label'))
3447            else:
3448                particles_types_with_LJ_parameters.append(particle_type)
3449        # Set up LJ interactions between all particle types
3450        for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 
3451            particle_name1 = self.find_value_from_es_type(es_type=type_pair[0],
3452                                                        column_name="name")
3453            particle_name2 = self.find_value_from_es_type(es_type=type_pair[1],
3454                                                        column_name="name")
3455            lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1,
3456                                                 particle_name2 = particle_name2,
3457                                                 combining_rule = combining_rule)
3458            
3459            # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
3460            if not lj_parameters:
3461                continue
3462            espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 
3463                                                                                    sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 
3464                                                                                    cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude,
3465                                                                                    offset = lj_parameters["offset"].to("reduced_length").magnitude, 
3466                                                                                    shift = shift)                                                                                          
3467            index = len(self.df)
3468            label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label")
3469            label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label")
3470            self.df.at [index, 'name'] = f'LJ: {label1}-{label2}'
3471            lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params()
3472
3473            self.add_value_to_df(index=index,
3474                                key=('pmb_type',''),
3475                                new_value='LennardJones')
3476            
3477            self.add_value_to_df(index=index,
3478                                key=('parameters_of_the_potential',''),
3479                                new_value=lj_params,
3480                                non_standard_value=True)
3481        if non_parametrized_labels and warnings:
3482            warnings.warn(f'The following particles do not have a defined value of sigma or epsilon in pmb.df: {non_parametrized_labels}. No LJ interaction has been added in ESPResSo for those particles.',UserWarning) 
3483        return
3484
3485    def write_pmb_df (self, filename):
3486        '''
3487        Writes the pyMBE dataframe into a csv file
3488        
3489        Args:
3490            filename(`str`): Path to the csv file 
3491        '''
3492
3493        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence']
3494        df = self.df.copy(deep=True)
3495        for column_name in columns_with_list_or_dict:
3496            df[column_name] = df[column_name].apply(lambda x: json.dumps(x) if isinstance(x, (np.ndarray, tuple, list, dict)) or pd.notnull(x) else x)
3497        df['bond_object'] = df['bond_object'].apply(lambda x: f'{x.__class__.__name__}({json.dumps({**x.get_params(), "bond_id": x._bond_id})})' if pd.notnull(x) else x)
3498        df.fillna(pd.NA, inplace=True)
3499        df.to_csv(filename)
3500        return

The library for the Molecular Builder for ESPResSo (pyMBE)

Attributes:
  • N_A(pint.Quantity): Avogadro number.
  • Kb(pint.Quantity): Boltzmann constant.
  • e(pint.Quantity): Elementary charge.
  • df(Pandas.Dataframe): Dataframe used to bookkeep all the information stored in pyMBE. Typically refered as pmb.df.
  • kT(pint.Quantity): Thermal energy.
  • Kw(pint.Quantity): Ionic product of water. Used in the setup of the G-RxMC method.
pymbe_library(seed, temperature=None, unit_length=None, unit_charge=None, Kw=None)
56    def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None):
57        """
58        Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 
59        and sets up  the `pmb.df` for bookkeeping.
60
61        Args:
62            temperature(`pint.Quantity`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None.
63            unit_length(`pint.Quantity`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None.
64            unit_charge (`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
65            Kw (`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
66        
67        Note:
68            - If no `temperature` is given, a value of 298.15 K is assumed by default.
69            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
70            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
71            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
72        """
73        # Seed and RNG
74        self.seed=seed
75        self.rng = np.random.default_rng(seed)
76        self.units=pint.UnitRegistry()
77        self.N_A=scipy.constants.N_A / self.units.mol
78        self.kB=scipy.constants.k * self.units.J / self.units.K
79        self.e=scipy.constants.e * self.units.C
80        self.set_reduced_units(unit_length=unit_length, 
81                               unit_charge=unit_charge,
82                               temperature=temperature, 
83                               Kw=Kw, 
84                               verbose=False)
85        self.setup_df()
86        self.lattice_builder = None
87        return

Initializes the pymbe_library by setting up the reduced unit system with temperature and reduced_length and sets up the pmb.df for bookkeeping.

Arguments:
  • temperature(pint.Quantity,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None.
  • unit_length(pint.Quantity, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None.
  • unit_charge (pint.Quantity,optional): Reduced unit of charge 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.
Note:
  • 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.
seed
rng
units
N_A
kB
e
lattice_builder
def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False):
 89    def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False):
 90        """
 91        Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles.
 92
 93        Args:
 94            particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles
 95            particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles
 96            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False.
 97
 98        Returns:
 99            index(`int`): Row index where the bond information has been added in pmb.df.
100        """
101        particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0]
102        particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0]
103        
104        bond_key = self.find_bond_key(particle_name1=particle_name1,
105                                    particle_name2=particle_name2, 
106                                    use_default_bond=use_default_bond)
107        if not bond_key:
108            return None
109        self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1)
110        indexs = np.where(self.df['name']==bond_key)
111        index_list = list (indexs[0])
112        used_bond_df = self.df.loc[self.df['particle_id2'].notnull()]
113        #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict
114        used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 )
115        used_bond_index = used_bond_df.index.to_list()
116        if not index_list:
117            return None
118        for index in index_list:
119            if index not in used_bond_index:
120                self.clean_df_row(index=int(index))
121                self.df.at[index,'particle_id'] = particle_id1
122                self.df.at[index,'particle_id2'] = particle_id2
123                break
124        return index

Adds a bond entry on the pymbe.df storing the particle_ids of the two bonded particles.

Arguments:
  • particle_id1(int): particle_id of the type of the first particle type of the bonded particles
  • particle_id2(int): particle_id of the type of the second particle type of the bonded particles
  • use_default_bond(bool, optional): Controls if a bond of type default is used to bond particle whose bond types are not defined in pmb.df. Defaults to False.
Returns:

index(int): Row index where the bond information has been added in pmb.df.

def add_bonds_to_espresso(self, espresso_system):
126    def add_bonds_to_espresso(self, espresso_system) :
127        """
128        Adds all bonds defined in `pmb.df` to `espresso_system`.
129
130        Args:
131            espresso_system(`espressomd.system.System`): system object of espressomd library
132        """
133
134        if 'bond' in self.df["pmb_type"].values:
135            bond_df = self.df.loc[self.df ['pmb_type'] == 'bond']
136            bond_list = bond_df.bond_object.values.tolist()
137            for bond in bond_list:
138                espresso_system.bonded_inter.add(bond)
139
140        else:
141            print ('WARNING: There are no bonds defined in pymbe.df')
142        
143        return

Adds all bonds defined in pmb.df to espresso_system.

Arguments:
  • espresso_system(espressomd.system.System): system object of espressomd library
def add_value_to_df( self, index, key, new_value, non_standard_value=False, overwrite=False):
145    def add_value_to_df(self,index,key,new_value, non_standard_value=False, overwrite=False):
146        """
147        Adds a value to a cell in the `pmb.df` DataFrame.
148
149        Args:
150            index(`int`): index of the row to add the value to.
151            key(`str`): the column label to add the value to.
152            non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False.
153            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
154        """
155
156        token = "#protected:"
157
158        def protect(obj):
159            if non_standard_value:
160                return token + json.dumps(obj, cls=self.NumpyEncoder)
161            return obj
162
163        def deprotect(obj):
164            if non_standard_value and isinstance(obj, str) and obj.startswith(token):
165                return json.loads(obj.removeprefix(token))
166            return obj
167
168        # Make sure index is a scalar integer value
169        index = int(index)
170        assert isinstance(index, int), '`index` should be a scalar integer value.'
171        idx = pd.IndexSlice
172        if self.check_if_df_cell_has_a_value(index=index,key=key):
173            old_value = self.df.loc[index,idx[key]]
174            if not pd.Series([protect(old_value)]).equals(pd.Series([protect(new_value)])):
175                name=self.df.loc[index,('name','')]
176                pmb_type=self.df.loc[index,('pmb_type','')]
177                logging.debug(f"You are attempting to redefine the properties of {name} of pmb_type {pmb_type}")    
178                if overwrite:
179                    logging.info(f'Overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}')
180                if not overwrite:
181                    logging.debug(f"pyMBE has preserved of the entry `{key}`: old_value = {old_value}. If you want to overwrite it with new_value = {new_value}, activate the switch overwrite = True ")
182                    return
183
184        self.df.loc[index,idx[key]] = protect(new_value)
185        if non_standard_value:
186            self.df[key] = self.df[key].apply(deprotect)
187        return

Adds a value to a cell in the pmb.df DataFrame.

Arguments:
  • index(int): index of the row to add the value to.
  • key(str): the column label to add the value to.
  • non_standard_value(bool, optional): Switch to enable insertion of non-standard values, such as dict objects. Defaults to False.
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
def assign_molecule_id(self, molecule_index):
189    def assign_molecule_id(self, molecule_index):
190        """
191        Assigns the `molecule_id` of the pmb object given by `pmb_type`
192        
193        Args:
194            molecule_index(`int`): index of the current `pmb_object_type` to assign the `molecule_id`
195        Returns:
196            molecule_id(`int`): Id of the molecule
197        """
198
199        self.clean_df_row(index=int(molecule_index))
200        
201        if self.df['molecule_id'].isnull().values.all():
202            molecule_id = 0        
203        else:
204            molecule_id = self.df['molecule_id'].max() +1
205
206        self.add_value_to_df (key=('molecule_id',''),
207                                index=int(molecule_index),
208                                new_value=molecule_id)
209
210        return molecule_id

Assigns the molecule_id of the pmb object given by pmb_type

Arguments:
  • molecule_index(int): index of the current pmb_object_type to assign the molecule_id
Returns:

molecule_id(int): Id of the molecule

def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system):
212    def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system):
213        """
214        Calculates the center of the molecule with a given molecule_id.
215
216        Args:
217            molecule_id(`int`): Id of the molecule whose center of mass is to be calculated.
218            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
219        
220        Returns:
221            center_of_mass(`lst`): Coordinates of the center of mass.
222        """
223        center_of_mass = np.zeros(3)
224        axis_list = [0,1,2]
225        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
226        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
227        for pid in particle_id_list:
228            for axis in axis_list:
229                center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis]
230        center_of_mass = center_of_mass / len(particle_id_list)
231        return center_of_mass

Calculates the center of the molecule with a given molecule_id.

Arguments:
  • molecule_id(int): Id of the molecule whose center of mass is to be calculated.
  • espresso_system(espressomd.system.System): Instance of a system object from the espressomd library.
Returns:

center_of_mass(lst): Coordinates of the center of mass.

def calculate_HH(self, molecule_name, pH_list=None, pka_set=None):
233    def calculate_HH(self, molecule_name, pH_list=None, pka_set=None):
234        """
235        Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 
236        for molecules with the name `molecule_name`.
237
238        Args:
239            molecule_name(`str`): name of the molecule to calculate the ideal charge for
240            pH_list(`lst`): pH-values to calculate. 
241            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
242
243        Returns:
244            Z_HH(`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list`
245
246        Note:
247            - This function supports objects with pmb types: "molecule", "peptide" and "protein".
248            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
249            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
250            - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan`  should be used.
251        """
252        if pH_list is None:
253            pH_list=np.linspace(2,12,50)
254        if pka_set is None:
255            pka_set=self.get_pka_set() 
256        self.check_pka_set(pka_set=pka_set)
257        charge_number_map = self.get_charge_number_map()
258        Z_HH=[]
259        for pH_value in pH_list:    
260            Z=0
261            index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 
262            residue_list = self.df.at [index,('residue_list','')]
263            sequence = self.df.at [index,('sequence','')]
264            if np.any(pd.isnull(sequence)):
265                # Molecule has no sequence
266                for residue in residue_list:
267                    list_of_particles_in_residue = self.search_particles_in_residue(residue)
268                    for particle in list_of_particles_in_residue:
269                        if particle in pka_set.keys():
270                            if pka_set[particle]['acidity'] == 'acidic':
271                                psi=-1
272                            elif pka_set[particle]['acidity']== 'basic':
273                                psi=+1
274                            else:
275                                psi=0
276                            Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value'])))                      
277                Z_HH.append(Z)
278            else:
279                # Molecule has a sequence
280                for name in sequence:
281                    if name in pka_set.keys():
282                        if pka_set[name]['acidity'] == 'acidic':
283                            psi=-1
284                        elif pka_set[name]['acidity']== 'basic':
285                            psi=+1
286                        else:
287                            psi=0
288                        Z+=psi/(1+10**(psi*(pH_value-pka_set[name]['pka_value'])))
289                    else:
290                        state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
291                        Z+=charge_number_map[state_one_type]
292                Z_HH.append(Z)
293
294        return Z_HH

Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve for molecules with the name molecule_name.

Arguments:
  • molecule_name(str): name of the molecule to calculate the ideal charge for
  • pH_list(lst): pH-values to calculate.
  • pka_set(dict): {"name" : {"pka_value": pka, "acidity": acidity}}
Returns:

Z_HH(lst): Henderson-Hasselbalch prediction of the charge of sequence in pH_list

Note:
  • This function supports objects with pmb types: "molecule", "peptide" and "protein".
  • If no pH_list is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
  • If no pka_set is given, the pKa values are taken from pmb.df
  • This function should only be used for single-phase systems. For two-phase systems pmb.calculate_HH_Donnan should be used.
def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
296    def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
297        """
298        Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning.
299
300        Args:
301            c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
302            c_salt('float'): Salt concentration in the reservoir.
303            pH_list('lst'): List of pH-values in the reservoir. 
304            pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}.
305
306        Returns:
307            {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
308            pH_system_list ('lst'): List of pH_values in the system.
309            partition_coefficients_list ('lst'): List of partition coefficients of cations.
310
311        Note:
312            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
313            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
314            - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result.
315        """
316        if pH_list is None:
317            pH_list=np.linspace(2,12,50)
318        if pka_set is None:
319            pka_set=self.get_pka_set() 
320        self.check_pka_set(pka_set=pka_set)
321
322        partition_coefficients_list = []
323        pH_system_list = []
324        Z_HH_Donnan={}
325        for key in c_macro:
326            Z_HH_Donnan[key] = []
327
328        def calc_charges(c_macro, pH):
329            """
330            Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation.
331
332            Args:
333                c_macro('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
334                pH('float'): pH-value that is used in the HH equation.
335
336            Returns:
337                charge('dict'): {"molecule_name": charge}
338            """
339            charge = {}
340            for name in c_macro:
341                charge[name] = self.calculate_HH(name, [pH], pka_set)[0]
342            return charge
343
344        def calc_partition_coefficient(charge, c_macro):
345            """
346            Calculates the partition coefficients of positive ions according to the ideal Donnan theory.
347
348            Args:
349                charge('dict'): {"molecule_name": charge}
350                c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
351            """
352            nonlocal ionic_strength_res
353            charge_density = 0.0
354            for key in charge:
355                charge_density += charge[key] * c_macro[key]
356            return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude
357
358        for pH_value in pH_list:    
359            # calculate the ionic strength of the reservoir
360            if pH_value <= 7.0:
361                ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 
362            elif pH_value > 7.0:
363                ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt
364
365            #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations
366            #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation.
367            #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects.
368            equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro))
369            logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq")
370            partition_coefficient = 10**logxi.root
371
372            charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient))
373            for key in c_macro:
374                Z_HH_Donnan[key].append(charges_temp[key])
375
376            pH_system_list.append(pH_value - np.log10(partition_coefficient))
377            partition_coefficients_list.append(partition_coefficient)
378
379        return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}

Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning.

Arguments:
  • c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system.
  • c_salt('float'): Salt concentration in the reservoir.
  • pH_list('lst'): List of pH-values in the reservoir.
  • pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}.
Returns:

{"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} pH_system_list ('lst'): List of pH_values in the system. partition_coefficients_list ('lst'): List of partition coefficients of cations.

Note:
  • If no pH_list is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
  • If no pka_set is given, the pKa values are taken from pmb.df
  • If c_macro does not contain all charged molecules in the system, this function is likely to provide the wrong result.
def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset):
381    def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset):
382        """
383        Calculates the initial bond length that is used when setting up molecules,
384        based on the minimum of the sum of bonded and short-range (LJ) interactions.
385        
386        Args:
387            bond_object(`espressomd.interactions.BondedInteractions`): instance of a bond object from espressomd library
388            bond_type(`str`): label identifying the used bonded potential
389            epsilon(`pint.Quantity`): LJ epsilon of the interaction between the particles
390            sigma(`pint.Quantity`): LJ sigma of the interaction between the particles
391            cutoff(`pint.Quantity`): cutoff-radius of the LJ interaction 
392            offset(`pint.Quantity`): offset of the LJ interaction
393        """    
394        def truncated_lj_potential(x, epsilon, sigma, cutoff,offset):
395            if x>cutoff:
396                return 0.0
397            else:
398                return 4*epsilon*((sigma/(x-offset))**12-(sigma/(x-offset))**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6)
399
400        epsilon_red=epsilon.to('reduced_energy').magnitude
401        sigma_red=sigma.to('reduced_length').magnitude
402        cutoff_red=cutoff.to('reduced_length').magnitude
403        offset_red=offset.to('reduced_length').magnitude
404
405        if bond_type == "harmonic":
406            r_0 = bond_object.params.get('r_0')
407            k = bond_object.params.get('k')
408            l0 = scipy.optimize.minimize(lambda x: 0.5*k*(x-r_0)**2 + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red, offset_red), x0=r_0).x
409        elif bond_type == "FENE":
410            r_0 = bond_object.params.get('r_0')
411            k = bond_object.params.get('k')
412            d_r_max = bond_object.params.get('d_r_max')
413            l0 = scipy.optimize.minimize(lambda x: -0.5*k*(d_r_max**2)*np.log(1-((x-r_0)/d_r_max)**2) + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red,offset_red), x0=1.0).x
414        return l0

Calculates the initial bond length that is used when setting up molecules, based on the minimum of the sum of bonded and short-range (LJ) interactions.

Arguments:
  • bond_object(espressomd.interactions.BondedInteractions): instance of a bond object from espressomd library
  • bond_type(str): label identifying the used bonded potential
  • epsilon(pint.Quantity): LJ epsilon of the interaction between the particles
  • sigma(pint.Quantity): LJ sigma of the interaction between the particles
  • cutoff(pint.Quantity): cutoff-radius of the LJ interaction
  • offset(pint.Quantity): offset of the LJ interaction
def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False):
416    def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False):
417        '''
418        Calculates the net charge per molecule of molecules with `name` = molecule_name. 
419        Returns the net charge per molecule and a maps with the net charge per residue and molecule.
420
421        Args:
422            espresso_system(`espressomd.system.System`): system information 
423            molecule_name(`str`): name of the molecule to calculate the net charge
424            dimensionless(`bool'): sets if the charge is returned with a dimension or not
425
426        Returns:
427            (`dict`): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }}
428
429        Note:
430            - The net charge of the molecule is averaged over all molecules of type `name` 
431            - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name`
432        '''        
433        valid_pmb_types = ["molecule", "protein"]
434        pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0]
435        if pmb_type not in valid_pmb_types:
436            raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}")      
437
438        id_map = self.get_particle_id_map(object_name=molecule_name)
439        def create_charge_map(espresso_system,id_map,label):
440            charge_number_map={}
441            for super_id in id_map[label].keys():
442                if dimensionless:
443                    net_charge=0
444                else:
445                    net_charge=0 * self.units.Quantity(1,'reduced_charge')
446                for pid in id_map[label][super_id]:
447                    if dimensionless:
448                        net_charge+=espresso_system.part.by_id(pid).q
449                    else:
450                        net_charge+=espresso_system.part.by_id(pid).q * self.units.Quantity(1,'reduced_charge')
451                charge_number_map[super_id]=net_charge
452            return charge_number_map
453        net_charge_molecules=create_charge_map(label="molecule_map",
454                                                espresso_system=espresso_system,
455                                                id_map=id_map)
456        net_charge_residues=create_charge_map(label="residue_map",
457                                                espresso_system=espresso_system,
458                                                id_map=id_map)       
459        if dimensionless:
460            mean_charge=np.mean(np.array(list(net_charge_molecules.values())))
461        else:
462            mean_charge=np.mean(np.array([value.magnitude for value in net_charge_molecules.values()]))*self.units.Quantity(1,'reduced_charge')
463        return {"mean": mean_charge, "molecules": net_charge_molecules, "residues": net_charge_residues}

Calculates the net charge per molecule of molecules with name = molecule_name. Returns the net charge per molecule and a maps with the net charge per residue and molecule.

Arguments:
  • espresso_system(espressomd.system.System): system information
  • molecule_name(str): name of the molecule to calculate the net charge
  • dimensionless(`bool'): sets if the charge is returned with a dimension or not
Returns:

(dict): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }}

Note:
  • The net charge of the molecule is averaged over all molecules of type name
  • The net charge of each particle type is averaged over all particle of the same type in all molecules of type name
def center_molecule_in_simulation_box(self, molecule_id, espresso_system):
465    def center_molecule_in_simulation_box(self, molecule_id, espresso_system):
466        """
467        Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`.
468        
469        Args:
470            molecule_id(`int`): Id of the molecule to be centered.
471            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
472        """
473        if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0:
474            raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.")      
475        center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system)
476        box_center = [espresso_system.box_l[0]/2.0,
477                      espresso_system.box_l[1]/2.0,
478                      espresso_system.box_l[2]/2.0]
479        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
480        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
481        for pid in particle_id_list:
482            es_pos = espresso_system.part.by_id(pid).pos
483            espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center
484        return 

Centers the pmb object matching molecule_id in the center of the simulation box in espresso_md.

Arguments:
  • molecule_id(int): Id of the molecule to be centered.
  • espresso_system(espressomd.system.System): Instance of a system object from the espressomd library.
def check_aminoacid_key(self, key):
486    def check_aminoacid_key(self, key):
487        """
488        Checks if `key` corresponds to a valid aminoacid letter code.
489
490        Args:
491            key(`str`): key to be checked.
492
493        Returns:
494            `bool`: True if `key` is a valid aminoacid letter code, False otherwise.
495        """
496        valid_AA_keys=['V', #'VAL'
497                       'I', #'ILE'
498                       'L', #'LEU'
499                       'E', #'GLU'
500                       'Q', #'GLN'
501                       'D', #'ASP'
502                       'N', #'ASN'
503                       'H', #'HIS'
504                       'W', #'TRP'
505                       'F', #'PHE'
506                       'Y', #'TYR'
507                       'R', #'ARG' 
508                       'K', #'LYS'
509                       'S', #'SER'
510                       'T', #'THR'
511                       'M', #'MET'
512                       'A', #'ALA'
513                       'G', #'GLY'
514                       'P', #'PRO'
515                       'C'] #'CYS'
516        if key in valid_AA_keys:
517            return True
518        else:
519            return False

Checks if key corresponds to a valid aminoacid letter code.

Arguments:
  • key(str): key to be checked.
Returns:

bool: True if key is a valid aminoacid letter code, False otherwise.

def check_dimensionality(self, variable, expected_dimensionality):
521    def check_dimensionality(self, variable, expected_dimensionality):
522        """
523        Checks if the dimensionality of `variable` matches `expected_dimensionality`.
524
525        Args:
526            variable(`pint.Quantity`): Quantity to be checked.
527            expected_dimensionality(`str`): Expected dimension of the variable.
528
529        Returns:
530            (`bool`): `True` if the variable if of the expected dimensionality, `False` otherwise.
531
532        Note:
533            - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality).
534            - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"`    
535        """
536        correct_dimensionality=variable.check(f"{expected_dimensionality}")      
537        if not correct_dimensionality:
538            raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}")
539        return correct_dimensionality

Checks if the dimensionality of variable matches expected_dimensionality.

Arguments:
  • variable(pint.Quantity): Quantity to be checked.
  • expected_dimensionality(str): Expected dimension of the variable.
Returns:

(bool): True if the variable if of the expected dimensionality, False otherwise.

Note:
  • expected_dimensionality takes dimensionality following the Pint standards docs.
  • For example, to check for a variable corresponding to a velocity expected_dimensionality = "[length]/[time]"
def check_if_df_cell_has_a_value(self, index, key):
541    def check_if_df_cell_has_a_value(self, index,key):
542        """
543        Checks if a cell in the `pmb.df` at the specified index and column has a value.
544
545        Args:
546            index(`int`): Index of the row to check.
547            key(`str`): Column label to check.
548
549        Returns:
550            `bool`: `True` if the cell has a value, `False` otherwise.
551        """
552        idx = pd.IndexSlice
553        import warnings
554        with warnings.catch_warnings():
555            warnings.simplefilter("ignore")
556            return not pd.isna(self.df.loc[index, idx[key]])

Checks if a cell in the pmb.df at the specified index and column has a value.

Arguments:
  • index(int): Index of the row to check.
  • key(str): Column label to check.
Returns:

bool: True if the cell has a value, False otherwise.

def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined):
558    def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined):
559        """
560        Checks if `name` is defined in `pmb.df`.
561
562        Args:
563            name(`str`): label to check if defined in `pmb.df`.
564            pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`.
565
566        Returns:
567            `bool`: `True` for success, `False` otherwise.
568        """
569        if name in self.df['name'].unique():
570            current_object_type = self.df[self.df['name']==name].pmb_type.values[0]
571            if current_object_type != pmb_type_to_be_defined:
572                raise ValueError (f"The name {name} is already defined in the df with a pmb_type = {current_object_type}, pymMBE does not support objects with the same name but different pmb_types")
573            return True            
574        else:
575            return False

Checks if name is defined in pmb.df.

Arguments:
  • name(str): label to check if defined in pmb.df.
  • pmb_type_to_be_defined(str): pmb object type corresponding to name.
Returns:

bool: True for success, False otherwise.

def check_if_metal_ion(self, key):
577    def check_if_metal_ion(self,key):
578        """
579        Checks if `key` corresponds to a label of a supported metal ion.
580
581        Args:
582            key(`str`): key to be checked
583
584        Returns:
585            (`bool`): True if `key`  is a supported metal ion, False otherwise.
586        """
587        if key in self.get_metal_ions_charge_number_map().keys():
588            return True
589        else:
590            return False

Checks if key corresponds to a label of a supported metal ion.

Arguments:
  • key(str): key to be checked
Returns:

(bool): True if key is a supported metal ion, False otherwise.

def check_pka_set(self, pka_set):
592    def check_pka_set(self, pka_set):
593        """
594        Checks that `pka_set` has the formatting expected by the pyMBE library.
595       
596        Args:
597            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
598        """
599        required_keys=['pka_value','acidity']
600        for required_key in required_keys:
601            for pka_name, pka_entry in pka_set.items():
602                if required_key not in pka_entry:
603                    raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")')
604        return

Checks that pka_set has the formatting expected by the pyMBE library.

Arguments:
  • pka_set(dict): {"name" : {"pka_value": pka, "acidity": acidity}}
def clean_df_row( self, index, columns_keys_to_clean=('particle_id', 'particle_id2', 'residue_id', 'molecule_id')):
606    def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")):
607        """
608        Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a pd.NA value.
609
610        Args:
611            index(`int`): Index of the row to clean.
612            columns_keys_to_clean(`list` of `str`, optional): List with the column keys to be cleaned. Defaults to [`particle_id`, `particle_id2`, `residue_id`, `molecule_id`].
613        """   
614        for column_key in columns_keys_to_clean:
615            self.add_value_to_df(key=(column_key,''),index=index,new_value=pd.NA)
616        self.df.fillna(pd.NA, inplace=True)
617        return

Cleans the columns of pmb.df in columns_keys_to_clean of the row with index index by assigning them a pd.NA value.

Arguments:
  • index(int): Index of the row to clean.
  • columns_keys_to_clean(list of str, optional): List with the column keys to be cleaned. Defaults to [particle_id, particle_id2, residue_id, molecule_id].
def convert_columns_to_original_format(self, df):
619    def convert_columns_to_original_format (self, df):
620        """
621        Converts the columns of the Dataframe to the original format in pyMBE.
622        
623        Args:
624            df(`DataFrame`): dataframe with pyMBE information as a string  
625        
626        """
627
628        columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', ('state_one','es_type'),('state_two','es_type'),('state_one','z'),('state_two','z') ]  
629
630        columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset']
631
632        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence', 'chain_map', 'node_map']
633
634        for column_name in columns_dtype_int:
635            df[column_name] = df[column_name].astype(pd.Int64Dtype())
636            
637        for column_name in columns_with_list_or_dict:
638            if df[column_name].isnull().all():
639                df[column_name] = df[column_name].astype(object)
640            else:
641                df[column_name] = df[column_name].apply(lambda x: json.loads(x) if pd.notnull(x) else x)
642
643        for column_name in columns_with_units:
644            df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x)
645
646        df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x)
647        df["l0"] = df["l0"].astype(object)
648        df["pka"] = df["pka"].astype(object)
649        return df

Converts the columns of the Dataframe to the original format in pyMBE.

Arguments:
  • df(DataFrame): dataframe with pyMBE information as a string
def convert_str_to_bond_object(self, bond_str):
651    def convert_str_to_bond_object (self, bond_str):
652        """
653        Convert a row read as a `str` to the corresponding ESPResSo bond object. 
654
655        Args:
656            bond_str(`str`): string with the information of a bond object.
657
658        Returns:
659            bond_object(`obj`): ESPResSo bond object.
660
661        Note:
662            Current supported bonds are: HarmonicBond and FeneBond
663        """
664        import espressomd.interactions
665
666        supported_bonds = ['HarmonicBond', 'FeneBond']
667        m = re.search(r'^([A-Za-z0-9_]+)\((\{.+\})\)$', bond_str)
668        if m is None:
669            raise ValueError(f'Cannot parse bond "{bond_str}"')
670        bond = m.group(1)
671        if bond not in supported_bonds:
672            raise NotImplementedError(f"Bond type '{bond}' currently not implemented in pyMBE, accepted types are {supported_bonds}")
673        params = json.loads(m.group(2))
674        bond_id = params.pop("bond_id")
675        bond_object = getattr(espressomd.interactions, bond)(**params)
676        bond_object._bond_id = bond_id
677        return bond_object

Convert a row read as a str to the corresponding ESPResSo bond object.

Arguments:
  • bond_str(str): string with the information of a bond object.
Returns:

bond_object(obj): ESPResSo bond object.

Note:

Current supported bonds are: HarmonicBond and FeneBond

def copy_df_entry(self, name, column_name, number_of_copies):
679    def copy_df_entry(self, name, column_name, number_of_copies):
680        '''
681        Creates 'number_of_copies' of a given 'name' in `pymbe.df`.
682
683        Args:
684            name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df`
685            column_name(`str`): Column name to use as a filter. 
686            number_of_copies(`int`): number of copies of `name` to be created.
687        
688        Note:
689            - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 
690        '''
691
692        valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ]
693        if column_name not in valid_column_names:
694            raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}")
695        df_by_name = self.df.loc[self.df.name == name]
696        if number_of_copies != 1:           
697            if df_by_name[column_name].isnull().values.any():       
698                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True)
699            else:
700                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
701                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
702                df_by_name_repeated[column_name] = pd.NA
703            # Concatenate the new particle rows to  `df`
704            self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
705        else:
706            if not df_by_name[column_name].isnull().values.any():     
707                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
708                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
709                df_by_name_repeated[column_name] = pd.NA
710                self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
711        return

Creates 'number_of_copies' of a given 'name' in pymbe.df.

Arguments:
  • name(str): Label of the particle/residue/molecule type to be created. name must be defined in pmb.df
  • column_name(str): Column name to use as a filter.
  • number_of_copies(int): number of copies of name to be created.
Note:
  • Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id"
def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt, verbose=True):
713    def create_added_salt (self, espresso_system, cation_name, anion_name, c_salt, verbose=True):    
714        """
715        Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`.
716
717        Args:
718            espresso_system(`espressomd.system.System`): instance of an espresso system object.
719            cation_name(`str`): `name` of a particle with a positive charge.
720            anion_name(`str`): `name` of a particle with a negative charge.
721            c_salt(`float`): Salt concentration.
722            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
723            
724        Returns:
725            c_salt_calculated(`float`): Calculated salt concentration added to `espresso_system`.
726        """
727        cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
728        anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
729        if cation_name_charge <= 0:
730            raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge)
731        if anion_name_charge >= 0:
732            raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge)
733        # Calculate the number of ions in the simulation box
734        volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3')
735        if c_salt.check('[substance] [length]**-3'):
736            N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude)
737            c_salt_calculated=N_ions/(volume*self.N_A)
738        elif c_salt.check('[length]**-3'):
739            N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude)
740            c_salt_calculated=N_ions/volume
741        else:
742            raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt)
743        N_cation = N_ions*abs(anion_name_charge)
744        N_anion = N_ions*abs(cation_name_charge)
745        self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation)
746        self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion)
747        if verbose:
748            if c_salt_calculated.check('[substance] [length]**-3'):
749                print(f"\n Added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions")
750            elif c_salt_calculated.check('[length]**-3'):
751                print(f"\n Added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions")
752        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.
  • verbose(bool): switch to activate/deactivate verbose. Defaults to True.
Returns:

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

def create_bond_in_espresso(self, bond_type, bond_parameters):
754    def create_bond_in_espresso(self, bond_type, bond_parameters):
755        '''
756        Creates either a harmonic or a FENE bond in ESPResSo
757
758        Args:
759            bond_type(`str`): label to identify the potential to model the bond.
760            bond_parameters(`dict`): parameters of the potential of the bond.
761
762        Note:
763            Currently, only HARMONIC and FENE bonds are supported.
764
765            For a HARMONIC bond the dictionary must contain:
766
767                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
768                using the `pmb.units` UnitRegistry.
769                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
770                the `pmb.units` UnitRegistry.
771           
772            For a FENE bond the dictionary must additionally contain:
773                
774                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
775                units of length using the `pmb.units` UnitRegistry. Default 'None'.
776
777        Returns:
778              bond_object (`obj`): an ESPResSo bond object
779        '''
780        from espressomd import interactions
781
782        valid_bond_types   = ["harmonic", "FENE"]
783        
784        if 'k' in bond_parameters:
785            bond_magnitude     = bond_parameters['k'].to('reduced_energy / reduced_length**2')
786        else:
787            raise ValueError("Magnitude of the potential (k) is missing")
788        
789        if bond_type == 'harmonic':
790            if 'r_0' in bond_parameters:
791                bond_length        = bond_parameters['r_0'].to('reduced_length')
792            else:
793                raise ValueError("Equilibrium bond length (r_0) is missing")
794            bond_object    = interactions.HarmonicBond(k   = bond_magnitude.magnitude,
795                                                       r_0 = bond_length.magnitude)
796        elif bond_type == 'FENE':
797            if 'r_0' in bond_parameters:
798                bond_length        = bond_parameters['r_0'].to('reduced_length').magnitude
799            else:
800                print("WARNING: No value provided for r_0. Defaulting to r_0 = 0")
801                bond_length=0
802            if 'd_r_max' in bond_parameters:
803                max_bond_stret = bond_parameters['d_r_max'].to('reduced_length')
804            else:
805                raise ValueError("Maximal stretching length (d_r_max) is missing")
806            bond_object    = interactions.FeneBond(r_0     = bond_length, 
807                                                   k       = bond_magnitude.magnitude,
808                                                   d_r_max = max_bond_stret.magnitude)
809        else:
810            raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}")
811        return bond_object

Creates either a harmonic or a FENE bond in ESPResSo

Arguments:
  • bond_type(str): label to identify the potential to model the bond.
  • bond_parameters(dict): parameters of the potential of the bond.
Note:

Currently, only HARMONIC and FENE bonds are supported.

For a HARMONIC bond the dictionary must contain:

- k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
using the `pmb.units` UnitRegistry.
- r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
the `pmb.units` UnitRegistry.

For a FENE bond the dictionary must additionally contain:

- d_r_max (`obj`): Maximal stretching length for FENE. It should have 
units of length using the `pmb.units` UnitRegistry. Default 'None'.
Returns:

bond_object (obj): an ESPResSo bond object

def create_counterions( self, object_name, cation_name, anion_name, espresso_system, verbose=True):
814    def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True):
815        """
816        Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`.
817        
818        Args:
819            object_name(`str`): `name` of a pymbe object.
820            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
821            cation_name(`str`): `name` of a particle with a positive charge.
822            anion_name(`str`): `name` of a particle with a negative charge.
823            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
824
825        Returns: 
826            counterion_number(`dict`): {"name": number}
827        """
828        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0]
829        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0]
830        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
831        counterion_number={}
832        object_charge={}
833        for name in ['positive', 'negative']:
834            object_charge[name]=0
835        for id in object_ids:
836            if espresso_system.part.by_id(id).q > 0:
837                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
838            elif espresso_system.part.by_id(id).q < 0:
839                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
840        if object_charge['positive'] % abs(anion_charge) == 0:
841            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
842        else:
843            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
844        if object_charge['negative'] % abs(cation_charge) == 0:
845            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
846        else:
847            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
848        if counterion_number[cation_name] > 0: 
849            self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name])
850        else:
851            counterion_number[cation_name]=0
852        if counterion_number[anion_name] > 0:
853            self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name])
854        else:
855            counterion_number[anion_name] = 0
856        if verbose:
857            print('The following counter-ions have been created: ')
858            for name in counterion_number.keys():
859                print(f'Ion type: {name} created number: {counterion_number[name]}')
860        return counterion_number

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

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.
  • verbose(bool): switch to activate/deactivate verbose. Defaults to True.

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

def create_hydrogel(self, name, espresso_system):
862    def create_hydrogel(self, name, espresso_system):
863        """ 
864        creates the hydrogel `name` in espresso_system
865        Args:
866            name(`str`): Label of the hydrogel to be created. `name` must be defined in the `pmb.df`
867            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
868
869        Returns:
870            hydrogel_info(`dict`):  {"name":hydrogel_name, "chains": {chain_id1: {residue_id1: {'central_bead_id': central_bead_id, 'side_chain_ids': [particle_id1,...]},...,"node_start":node_start,"node_end":node_end}, chain_id2: {...},...}, "nodes":{node1:[node1_id],...}}     
871        """
872        if not self.check_if_name_is_defined_in_df(name=name, pmb_type_to_be_defined='hydrogel'):
873            raise ValueError(f"Hydrogel with name '{name}' is not defined in the DataFrame.")
874        hydrogel_info={"name":name, "chains":{}, "nodes":{}}
875        # placing nodes
876        node_positions = {}
877        node_topology = self.df[self.df["name"]==name]["node_map"].iloc[0]
878        for node_info in node_topology:
879            node_index = node_info["lattice_index"]
880            node_name = node_info["particle_name"]
881            node_pos, node_id = self.create_hydrogel_node(self.format_node(node_index), node_name, espresso_system)
882            hydrogel_info["nodes"][self.format_node(node_index)]=node_id
883            node_positions[node_id[0]]=node_pos
884        
885        # Placing chains between nodes
886        # Looping over all the 16 chains
887        chain_topology = self.df[self.df["name"]==name]["chain_map"].iloc[0]
888        for chain_info in chain_topology:
889            node_s = chain_info["node_start"]
890            node_e = chain_info["node_end"]
891            molecule_info = self.create_hydrogel_chain(node_s, node_e, node_positions, espresso_system)
892            for molecule_id in molecule_info:
893                hydrogel_info["chains"][molecule_id] = molecule_info[molecule_id]
894                hydrogel_info["chains"][molecule_id]["node_start"]=node_s
895                hydrogel_info["chains"][molecule_id]["node_end"]=node_e
896        return hydrogel_info

creates the hydrogel name in espresso_system

Arguments:
  • name(str): Label of the hydrogel to be created. name must be defined in the pmb.df
  • espresso_system(espressomd.system.System): Instance of a system object from the espressomd library.
Returns:

hydrogel_info(dict): {"name":hydrogel_name, "chains": {chain_id1: {residue_id1: {'central_bead_id': central_bead_id, 'side_chain_ids': [particle_id1,...]},...,"node_start":node_start,"node_end":node_end}, chain_id2: {...},...}, "nodes":{node1:[node1_id],...}}

def create_hydrogel_chain(self, node_start, node_end, node_positions, espresso_system):
898    def create_hydrogel_chain(self, node_start, node_end, node_positions, espresso_system):
899        """
900        Creates a chain between two nodes of a hydrogel.
901
902        Args:
903            node_start(`str`): name of the starting node particle at which the first residue of the chain will be attached. 
904            node_end(`str`): name of the ending node particle at which the last residue of the chain will be attached.
905            node_positions(`dict`): dictionary with the positions of the nodes. The keys are the node names and the values are the positions.
906            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
907
908        Note:
909            For example, 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.
910            The chain will be placed in the direction of the vector between `node_start` and `node_end`. 
911        """
912        if self.lattice_builder is None:
913            raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.")
914
915        molecule_name = "chain_"+node_start+"_"+node_end
916        sequence = self.df[self.df['name']==molecule_name].residue_list.values [0]
917        assert len(sequence) != 0 and not isinstance(sequence, str)
918        assert len(sequence) == self.lattice_builder.mpc
919
920        key, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end)
921        assert node_start != node_end or sequence == sequence[::-1], \
922            (f"chain cannot be defined between '{node_start}' and '{node_end}' since it "
923            "would form a loop with a non-symmetric sequence (under-defined stereocenter)")
924
925        if reverse:
926            sequence = sequence[::-1]
927
928        node_start_pos = np.array(list(int(x) for x in node_start.strip('[]').split()))*0.25*self.lattice_builder.box_l
929        node_end_pos = np.array(list(int(x) for x in node_end.strip('[]').split()))*0.25*self.lattice_builder.box_l
930        node1 = espresso_system.part.select(lambda p: (p.pos == node_start_pos).all()).id
931        node2 = espresso_system.part.select(lambda p: (p.pos == node_end_pos).all()).id
932
933        if not node1[0] in node_positions or not node2[0] in node_positions:
934            raise ValueError("Set node position before placing a chain between them")
935         
936        # Finding a backbone vector between node_start and node_end
937        vec_between_nodes = np.array(node_positions[node2[0]]) - np.array(node_positions[node1[0]])
938        vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l)
939        backbone_vector = list(vec_between_nodes/(self.lattice_builder.mpc + 1))
940        node_start_name = self.df[(self.df["particle_id"]==node1[0]) & (self.df["pmb_type"]=="particle")]["name"].values[0]
941        first_res_name = self.df[(self.df["pmb_type"]=="residue") & (self.df["name"]==sequence[0])]["central_bead"].values[0]
942        l0 = self.get_bond_length(node_start_name, first_res_name, hard_check=True)
943        chain_molecule_info = self.create_molecule(
944                name=molecule_name,  # Use the name defined earlier
945                number_of_molecules=1,  # Creating one chain
946                espresso_system=espresso_system,
947                list_of_first_residue_positions=[list(np.array(node_positions[node1[0]]) + np.array(backbone_vector))],#Start at the first node
948                backbone_vector=np.array(backbone_vector)/l0,
949                use_default_bond=False  # Use defaut bonds between monomers
950            )
951        # Collecting ids of beads of the chain/molecule
952        chain_ids = []
953        residue_ids = []
954        for molecule_id in chain_molecule_info:
955            for residue_id in chain_molecule_info[molecule_id]:
956                residue_ids.append(residue_id)
957                bead_id = chain_molecule_info[molecule_id][residue_id]['central_bead_id']
958                chain_ids.append(bead_id)
959
960        self.lattice_builder.chains[key] = sequence
961        # Search bonds between nodes and chain ends
962        BeadType_near_to_node_start = self.df[(self.df["residue_id"] == residue_ids[0]) & (self.df["central_bead"].notnull())]["central_bead"].drop_duplicates().iloc[0]
963        BeadType_near_to_node_end = self.df[(self.df["residue_id"] == residue_ids[-1]) & (self.df["central_bead"].notnull())]["central_bead"].drop_duplicates().iloc[0]
964        bond_node1_first_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_start],
965                                                    particle_name2 = BeadType_near_to_node_start,
966                                                    hard_check=False,
967                                                    use_default_bond=False)
968        bond_node2_last_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_end],
969                                                    particle_name2 = BeadType_near_to_node_end,
970                                                    hard_check=False,
971                                                    use_default_bond=False)
972
973        espresso_system.part.by_id(node1[0]).add_bond((bond_node1_first_monomer, chain_ids[0]))
974        espresso_system.part.by_id(node2[0]).add_bond((bond_node2_last_monomer, chain_ids[-1]))
975        # Add bonds to data frame
976        bond_index1 = self.add_bond_in_df(particle_id1=node1[0],
977                                    particle_id2=chain_ids[0],
978                                    use_default_bond=False)
979        self.add_value_to_df(key=('molecule_id',''),
980                        index=int(bond_index1),
981                        new_value=molecule_id,
982                        overwrite=True)
983        self.add_value_to_df(key=('residue_id',''),
984                            index=int (bond_index1),
985                            new_value=residue_ids[0],
986                            overwrite=True)
987        bond_index2 = self.add_bond_in_df(particle_id1=node2[0],
988                                    particle_id2=chain_ids[-1],
989                                    use_default_bond=False)
990        self.add_value_to_df(key=('molecule_id',''),
991                        index=int(bond_index2),
992                        new_value=molecule_id,
993                        overwrite=True)
994        self.add_value_to_df(key=('residue_id',''),
995                        index=int (bond_index2),
996                        new_value=residue_ids[-1],
997                        overwrite=True)
998        return chain_molecule_info

Creates a chain between two nodes of a hydrogel.

Arguments:
  • node_start(str): name of the starting node particle at which the first residue of the chain will be attached.
  • node_end(str): name of the ending node particle at which the last residue of the chain will be attached.
  • node_positions(dict): dictionary with the positions of the nodes. The keys are the node names and the values are the positions.
  • espresso_system(espressomd.system.System): Instance of a system object from the espressomd library.
Note:

For example, 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. The chain will be placed in the direction of the vector between node_start and node_end.

def create_hydrogel_node(self, node_index, node_name, espresso_system):
1000    def create_hydrogel_node(self, node_index, node_name, espresso_system):
1001        """
1002        Set a node residue type.
1003        
1004        Args:
1005            node_index(`str`): Lattice node index in the form of a string, e.g. "[0 0 0]".
1006            node_name(`str`): name of the node particle defined in pyMBE.
1007        Returns:
1008            node_position(`list`): Position of the node in the lattice.
1009            p_id(`int`): Particle ID of the node.
1010        """
1011        if self.lattice_builder is None:
1012            raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.")
1013
1014        node_position = np.array(list(int(x) for x in node_index.strip('[]').split()))*0.25*self.lattice_builder.box_l
1015        p_id = self.create_particle(name = node_name,
1016                         espresso_system=espresso_system,
1017                         number_of_particles=1,
1018                         position = [node_position])
1019        key = self.lattice_builder._get_node_by_label(node_index)
1020        self.lattice_builder.nodes[key] = node_name
1021
1022        return node_position.tolist(), p_id

Set a node residue type.

Arguments:
  • node_index(str): Lattice node index in the form of a string, e.g. "[0 0 0]".
  • node_name(str): name of the node particle defined in pyMBE.
Returns:

node_position(list): Position of the node in the lattice. p_id(int): Particle ID of the node.

def create_molecule( self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False):
1024    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False):
1025        """
1026        Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1027
1028        Args:
1029            name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df`
1030            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1031            number_of_molecules(`int`): Number of molecules of type `name` to be created.
1032            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.
1033            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. 
1034            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df`
1035
1036        Returns:
1037            molecules_info(`dict`):  {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 
1038        """
1039        if list_of_first_residue_positions is not None:
1040            for item in list_of_first_residue_positions:
1041                if not isinstance(item, list):
1042                    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.")
1043                elif len(item) != 3:
1044                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
1045
1046            if len(list_of_first_residue_positions) != number_of_molecules:
1047                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
1048        if number_of_molecules <= 0:
1049            return 0
1050        # Generate an arbitrary random unit vector
1051        if backbone_vector is None:
1052            backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
1053                                                    radius=1, 
1054                                                    n_samples=1,
1055                                                    on_surface=True)[0]
1056        else:
1057            backbone_vector = np.array(backbone_vector)
1058        self.check_if_name_is_defined_in_df(name=name,
1059                                            pmb_type_to_be_defined='molecule')
1060            
1061        first_residue = True
1062        molecules_info = {}
1063        residue_list = self.df[self.df['name']==name].residue_list.values [0]
1064        self.copy_df_entry(name=name,
1065                        column_name='molecule_id',
1066                        number_of_copies=number_of_molecules)
1067
1068        molecules_index = np.where(self.df['name']==name)
1069        molecule_index_list =list(molecules_index[0])[-number_of_molecules:]
1070        pos_index = 0 
1071        for molecule_index in molecule_index_list:        
1072            molecule_id = self.assign_molecule_id(molecule_index=molecule_index)
1073            molecules_info[molecule_id] = {}
1074            for residue in residue_list:
1075                if first_residue:
1076                    if list_of_first_residue_positions is None:
1077                        residue_position = None
1078                    else:
1079                        for item in list_of_first_residue_positions:
1080                            residue_position = [np.array(list_of_first_residue_positions[pos_index])]
1081                            
1082                    residues_info = self.create_residue(name=residue,
1083                                                        espresso_system=espresso_system, 
1084                                                        central_bead_position=residue_position,  
1085                                                        use_default_bond= use_default_bond, 
1086                                                        backbone_vector=backbone_vector)
1087                    residue_id = next(iter(residues_info))
1088                    # Add the correct molecule_id to all particles in the residue
1089                    for index in self.df[self.df['residue_id']==residue_id].index:
1090                        self.add_value_to_df(key=('molecule_id',''),
1091                                            index=int (index),
1092                                            new_value=molecule_id,
1093                                            overwrite=True)
1094                    central_bead_id = residues_info[residue_id]['central_bead_id']
1095                    previous_residue = residue
1096                    residue_position = espresso_system.part.by_id(central_bead_id).pos
1097                    previous_residue_id = central_bead_id
1098                    first_residue = False          
1099                else:                    
1100                    previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0]
1101                    new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0]
1102                    bond = self.search_bond(particle_name1=previous_central_bead_name, 
1103                                            particle_name2=new_central_bead_name, 
1104                                            hard_check=True, 
1105                                            use_default_bond=use_default_bond)                
1106                    l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 
1107                                            particle_name2=new_central_bead_name, 
1108                                            hard_check=True, 
1109                                            use_default_bond=use_default_bond)                
1110                    
1111                    residue_position = residue_position+backbone_vector*l0
1112                    residues_info = self.create_residue(name=residue, 
1113                                                        espresso_system=espresso_system, 
1114                                                        central_bead_position=[residue_position],
1115                                                        use_default_bond= use_default_bond, 
1116                                                        backbone_vector=backbone_vector)
1117                    residue_id = next(iter(residues_info))      
1118                    for index in self.df[self.df['residue_id']==residue_id].index:
1119                        self.add_value_to_df(key=('molecule_id',''),
1120                                            index=int (index),
1121                                            new_value=molecule_id,
1122                                            overwrite=True)            
1123                    central_bead_id = residues_info[residue_id]['central_bead_id']
1124                    espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id))
1125                    bond_index = self.add_bond_in_df(particle_id1=central_bead_id,
1126                                        particle_id2=previous_residue_id,
1127                                        use_default_bond=use_default_bond) 
1128                    self.add_value_to_df(key=('molecule_id',''),
1129                                            index=int (bond_index),
1130                                            new_value=molecule_id,
1131                                            overwrite=True)           
1132                    previous_residue_id = central_bead_id
1133                    previous_residue = residue                    
1134                molecules_info[molecule_id][residue_id] = residues_info[residue_id]
1135            first_residue = True
1136            pos_index+=1
1137        
1138        return molecules_info

Creates number_of_molecules molecule of type name into espresso_system and bookkeeps them into pmb.df.

Arguments:
  • name(str): Label of the molecule type to be created. name must be defined in pmb.df
  • espresso_system(espressomd.system.System): Instance of a system object from espressomd library.
  • number_of_molecules(int): Number of molecules 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 particle with undefined bonds in pymbe.df
Returns:

molecules_info(dict): {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}}

def create_particle( self, name, espresso_system, number_of_particles, position=None, fix=False):
1140    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
1141        """
1142        Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`.
1143        
1144        Args:
1145            name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 
1146            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1147            number_of_particles(`int`): Number of particles to be created.
1148            position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
1149            fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
1150        Returns:
1151            created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`.
1152        """       
1153        if number_of_particles <=0:
1154            return []
1155        self.check_if_name_is_defined_in_df(name=name,
1156                                       pmb_type_to_be_defined='particle')
1157        # Copy the data of the particle `number_of_particles` times in the `df`
1158        self.copy_df_entry(name=name,
1159                          column_name='particle_id',
1160                          number_of_copies=number_of_particles)
1161        # Get information from the particle type `name` from the df     
1162        z = self.df.loc[self.df['name']==name].state_one.z.values[0]
1163        z = 0. if z is None else z
1164        es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
1165        # Get a list of the index in `df` corresponding to the new particles to be created
1166        index = np.where(self.df['name']==name)
1167        index_list =list(index[0])[-number_of_particles:]
1168        # Create the new particles into  `espresso_system`
1169        created_pid_list=[]
1170        for index in range (number_of_particles):
1171            df_index=int (index_list[index])
1172            self.clean_df_row(index=df_index)
1173            if position is None:
1174                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1175            else:
1176                particle_position = position[index]
1177            if len(espresso_system.part.all()) == 0:
1178                bead_id = 0
1179            else:
1180                bead_id = max (espresso_system.part.all().id) + 1
1181            created_pid_list.append(bead_id)
1182            kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z)
1183            if fix:
1184                kwargs["fix"] = 3 * [fix]
1185            espresso_system.part.add(**kwargs)
1186            self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id)                  
1187        return created_pid_list

Creates number_of_particles particles of type name into espresso_system and bookkeeps them into pymbe.df.

Arguments:
  • name(str): Label of the particle type to be created. name must be a particle defined in pmb_df.
  • 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:

created_pid_list(list of float): List with the ids of the particles created into espresso_system.

def create_pmb_object( self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None):
1189    def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None):
1190        """
1191        Creates all `particle`s associated to `pmb object` into  `espresso` a number of times equal to `number_of_objects`.
1192        
1193        Args:
1194            name(`str`): Unique label of the `pmb object` to be created. 
1195            number_of_objects(`int`): Number of `pmb object`s to be created.
1196            espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library.
1197            position(`list`): Coordinates where the particles should be created.
1198            use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`.
1199            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. 
1200
1201        Note:
1202            - If no `position` is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length. 
1203        """
1204        allowed_objects=['particle','residue','molecule']
1205        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1206        if pmb_type not in allowed_objects:
1207            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1208        if pmb_type == 'particle':
1209            self.create_particle(name=name, 
1210                                number_of_particles=number_of_objects, 
1211                                espresso_system=espresso_system, 
1212                                position=position)
1213        elif pmb_type == 'residue':
1214            self.create_residue(name=name,  
1215                                espresso_system=espresso_system, 
1216                                central_bead_position=position,
1217                                use_default_bond=use_default_bond,
1218                                backbone_vector=backbone_vector)
1219        elif pmb_type == 'molecule':
1220            self.create_molecule(name=name, 
1221                                number_of_molecules=number_of_objects, 
1222                                espresso_system=espresso_system, 
1223                                use_default_bond=use_default_bond, 
1224                                list_of_first_residue_positions=position,
1225                                backbone_vector=backbone_vector)
1226        return

Creates all particles associated to pmb object into espresso a number of times equal to number_of_objects.

Arguments:
  • name(str): Unique label of the pmb object to be created.
  • number_of_objects(int): Number of pmb objects to be created.
  • espresso_system(espressomd.system.System): Instance of an espresso system object from espressomd library.
  • position(list): Coordinates where the particles should be created.
  • use_default_bond(bool,optional): Controls if a default bond is used to bond particles with undefined bonds in pmb.df. Defaults to False.
  • 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.
Note:
  • If no position is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length.
def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1228    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1229        """
1230        Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions`
1231
1232        Args:
1233            name(`str`): Label of the protein to be created. 
1234            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1235            number_of_proteins(`int`): Number of proteins to be created.
1236            positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}}
1237        """
1238
1239        if number_of_proteins <=0:
1240            return
1241        self.check_if_name_is_defined_in_df(name=name,
1242                                        pmb_type_to_be_defined='protein')
1243        self.copy_df_entry(name=name,
1244                            column_name='molecule_id',
1245                            number_of_copies=number_of_proteins)
1246        protein_index = np.where(self.df['name']==name)
1247        protein_index_list =list(protein_index[0])[-number_of_proteins:]
1248        
1249        box_half=espresso_system.box_l[0]/2.0
1250
1251        for molecule_index in protein_index_list:     
1252
1253            molecule_id = self.assign_molecule_id(molecule_index=molecule_index)
1254
1255            protein_center = self.generate_coordinates_outside_sphere(radius = 1, 
1256                                                                        max_dist=box_half, 
1257                                                                        n_samples=1, 
1258                                                                        center=[box_half]*3)[0]
1259   
1260            for residue in topology_dict.keys():
1261
1262                residue_name = re.split(r'\d+', residue)[0]
1263                residue_number = re.split(r'(\d+)', residue)[1]
1264                residue_position = topology_dict[residue]['initial_pos']
1265                position = residue_position + protein_center
1266
1267                particle_id = self.create_particle(name=residue_name,
1268                                                            espresso_system=espresso_system,
1269                                                            number_of_particles=1,
1270                                                            position=[position], 
1271                                                            fix = True)
1272                
1273                index = self.df[self.df['particle_id']==particle_id[0]].index.values[0]
1274                self.add_value_to_df(key=('residue_id',''),
1275                                            index=int (index),
1276                                            new_value=int(residue_number),
1277                                            overwrite=True)
1278
1279                self.add_value_to_df(key=('molecule_id',''),
1280                                        index=int (index),
1281                                        new_value=molecule_id,
1282                                        overwrite=True)
1283
1284        return

Creates number_of_proteins molecules of type name into espresso_system at the coordinates in positions

Arguments:
  • name(str): Label of the protein to be created.
  • espresso_system(espressomd.system.System): Instance of a system object from the espressomd library.
  • number_of_proteins(int): Number of proteins to be created.
  • positions(dict): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}}
def create_residue( self, name, espresso_system, central_bead_position=None, use_default_bond=False, backbone_vector=None):
1286    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None):
1287        """
1288        Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1289
1290        Args:
1291            name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df`
1292            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1293            central_bead_position(`list` of `float`): Position of the central bead.
1294            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 `pmb.df`
1295            backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`.
1296
1297        Returns:
1298            residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1299        """
1300        self.check_if_name_is_defined_in_df(name=name,
1301                                            pmb_type_to_be_defined='residue')
1302        # Copy the data of a residue in the `df
1303        self.copy_df_entry(name=name,
1304                            column_name='residue_id',
1305                            number_of_copies=1)
1306        residues_index = np.where(self.df['name']==name)
1307        residue_index_list =list(residues_index[0])[-1:]
1308        # search for defined particle and residue names
1309        particle_and_residue_df = self.df.loc[(self.df['pmb_type']== "particle") | (self.df['pmb_type']== "residue")]
1310        particle_and_residue_names = particle_and_residue_df["name"].tolist()
1311        for residue_index in residue_index_list:  
1312            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1313            for side_chain_element in side_chain_list:
1314                if side_chain_element not in particle_and_residue_names:              
1315                    raise ValueError (f"{side_chain_element} is not defined")
1316        # Internal bookkepping of the residue info (important for side-chain residues)
1317        # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1318        residues_info={}
1319        for residue_index in residue_index_list:     
1320            self.clean_df_row(index=int(residue_index))
1321            # Assign a residue_id
1322            if self.df['residue_id'].isnull().all():
1323                residue_id=0
1324            else:
1325                residue_id = self.df['residue_id'].max() + 1
1326            self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id)
1327            # create the principal bead   
1328            central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0]            
1329            central_bead_id = self.create_particle(name=central_bead_name,
1330                                                                espresso_system=espresso_system,
1331                                                                position=central_bead_position,
1332                                                                number_of_particles = 1)[0]
1333            central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1334            #assigns same residue_id to the central_bead particle created.
1335            index = self.df[self.df['particle_id']==central_bead_id].index.values[0]
1336            self.df.at [index,'residue_id'] = residue_id
1337            # Internal bookkeeping of the central bead id
1338            residues_info[residue_id]={}
1339            residues_info[residue_id]['central_bead_id']=central_bead_id
1340            # create the lateral beads  
1341            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1342            side_chain_beads_ids = []
1343            for side_chain_element in side_chain_list:  
1344                pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 
1345                if pmb_type == 'particle':
1346                    bond = self.search_bond(particle_name1=central_bead_name, 
1347                                            particle_name2=side_chain_element, 
1348                                            hard_check=True, 
1349                                            use_default_bond=use_default_bond)
1350                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1351                                              particle_name2=side_chain_element, 
1352                                              hard_check=True, 
1353                                              use_default_bond=use_default_bond)
1354
1355                    if backbone_vector is None:
1356                        bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1357                                                                 radius=l0, 
1358                                                                 n_samples=1,
1359                                                                 on_surface=True)[0]
1360                    else:
1361                        bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector),
1362                                                                                                    magnitude=l0)
1363                     
1364                    side_bead_id = self.create_particle(name=side_chain_element, 
1365                                                                    espresso_system=espresso_system,
1366                                                                    position=[bead_position], 
1367                                                                    number_of_particles=1)[0]
1368                    index = self.df[self.df['particle_id']==side_bead_id].index.values[0]
1369                    self.add_value_to_df(key=('residue_id',''),
1370                                        index=int (index),
1371                                        new_value=residue_id, 
1372                                        overwrite=True)
1373                    side_chain_beads_ids.append(side_bead_id)
1374                    espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id))
1375                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1376                                        particle_id2=side_bead_id,
1377                                        use_default_bond=use_default_bond)
1378                    self.add_value_to_df(key=('residue_id',''),
1379                                        index=int (index),
1380                                        new_value=residue_id, 
1381                                        overwrite=True)
1382
1383                elif pmb_type == 'residue':
1384                    central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0]
1385                    bond = self.search_bond(particle_name1=central_bead_name, 
1386                                            particle_name2=central_bead_side_chain, 
1387                                            hard_check=True, 
1388                                            use_default_bond=use_default_bond)
1389                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1390                                              particle_name2=central_bead_side_chain, 
1391                                              hard_check=True, 
1392                                              use_default_bond=use_default_bond)
1393                    if backbone_vector is None:
1394                        residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1395                                                                    radius=l0, 
1396                                                                    n_samples=1,
1397                                                                    on_surface=True)[0]
1398                    else:
1399                        residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1400                                                                                                        magnitude=l0)
1401                    lateral_residue_info = self.create_residue(name=side_chain_element, 
1402                                                                espresso_system=espresso_system,
1403                                                                central_bead_position=[residue_position],
1404                                                                use_default_bond=use_default_bond)
1405                    lateral_residue_dict=list(lateral_residue_info.values())[0]
1406                    central_bead_side_chain_id=lateral_residue_dict['central_bead_id']
1407                    lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids']
1408                    residue_id_side_chain=list(lateral_residue_info.keys())[0]
1409                    # Change the residue_id of the residue in the side chain to the one of the bigger residue
1410                    index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0]
1411                    self.add_value_to_df(key=('residue_id',''),
1412                                        index=int(index),
1413                                        new_value=residue_id, 
1414                                        overwrite=True)
1415                    # Change the residue_id of the particles in the residue in the side chain
1416                    side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids
1417                    for particle_id in side_chain_beads_ids:
1418                        index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0]
1419                        self.add_value_to_df(key=('residue_id',''),
1420                                            index=int (index),
1421                                            new_value=residue_id, 
1422                                            overwrite=True)
1423                    espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id))
1424                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1425                                        particle_id2=central_bead_side_chain_id,
1426                                        use_default_bond=use_default_bond)
1427                    self.add_value_to_df(key=('residue_id',''),
1428                                        index=int (index),
1429                                        new_value=residue_id, 
1430                                        overwrite=True)
1431                    # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue
1432                    for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index:        
1433                        self.add_value_to_df(key=('residue_id',''),
1434                                            index=int(index),
1435                                            new_value=residue_id, 
1436                                            overwrite=True)
1437            # Internal bookkeeping of the side chain beads ids
1438            residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids
1439        return  residues_info

Creates a residue of type name into espresso_system and bookkeeps them into pmb.df.

Arguments:
  • name(str): Label of the residue type to be created. name must be defined in pmb.df
  • 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 pmb.df
  • backbone_vector(list of float): Backbone vector of the molecule. All side chains are created perpendicularly to backbone_vector.
Returns:

residues_info(dict): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}

def create_variable_with_units(self, variable):
1441    def create_variable_with_units(self, variable):
1442        """
1443        Returns a pint object with the value and units defined in `variable`.
1444
1445        Args:
1446            variable(`dict` or `str`): {'value': value, 'units': units}
1447        Returns:
1448            variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry.
1449        """        
1450        if isinstance(variable, dict):
1451            value=variable.pop('value')
1452            units=variable.pop('units')
1453        elif isinstance(variable, str):
1454            value = float(re.split(r'\s+', variable)[0])
1455            units = re.split(r'\s+', variable)[1]
1456        variable_with_units=value*self.units(units)
1457        return variable_with_units 

Returns a pint object with the value and units defined in variable.

Arguments:
  • variable(dict or str): {'value': value, 'units': units}
Returns:

variable_with_units(obj): variable with units using the pyMBE UnitRegistry.

def define_AA_residues(self, sequence, model):
1459    def define_AA_residues(self, sequence, model):
1460        """
1461        Defines in `pmb.df` all the different residues in `sequence`.
1462
1463        Args:
1464            sequence(`lst`):  Sequence of the peptide or protein.
1465            model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1466
1467        Returns:
1468            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.             
1469        """
1470        residue_list = []
1471        for residue_name in sequence:
1472            if model == '1beadAA':
1473                central_bead = residue_name
1474                side_chains = []
1475            elif model == '2beadAA':
1476                if residue_name in ['c','n', 'G']: 
1477                    central_bead = residue_name
1478                    side_chains = []
1479                else:
1480                    central_bead = 'CA'              
1481                    side_chains = [residue_name]
1482            if residue_name not in residue_list:   
1483                self.define_residue(name = 'AA-'+residue_name, 
1484                                    central_bead = central_bead,
1485                                    side_chains = side_chains)              
1486            residue_list.append('AA-'+residue_name)
1487        return residue_list

Defines in pmb.df all the different residues in sequence.

Arguments:
  • sequence(lst): Sequence of the peptide or protein.
  • model(string): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
Returns:

residue_list(list of str): List of the names of the residues in the sequence of the molecule.

def define_bond(self, bond_type, bond_parameters, particle_pairs):
1489    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1490        
1491        '''
1492        Defines a pmb object of type `bond` in `pymbe.df`.
1493
1494        Args:
1495            bond_type(`str`): label to identify the potential to model the bond.
1496            bond_parameters(`dict`): parameters of the potential of the bond.
1497            particle_pairs(`lst`): list of the `names` of the `particles` to be bonded.
1498
1499        Note:
1500            Currently, only HARMONIC and FENE bonds are supported.
1501
1502            For a HARMONIC bond the dictionary must contain the following parameters:
1503
1504                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
1505                using the `pmb.units` UnitRegistry.
1506                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
1507                the `pmb.units` UnitRegistry.
1508           
1509            For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:
1510                
1511                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
1512                units of length using the `pmb.units` UnitRegistry. Default 'None'.
1513        '''
1514
1515        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1516        for particle_name1, particle_name2 in particle_pairs:
1517
1518            lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1,
1519                                                 particle_name2 = particle_name2,
1520                                                 combining_rule = 'Lorentz-Berthelot')
1521      
1522            l0 = self.calculate_initial_bond_length(bond_object = bond_object,
1523                                                    bond_type   = bond_type,
1524                                                    epsilon     = lj_parameters["epsilon"],
1525                                                    sigma       = lj_parameters["sigma"],
1526                                                    cutoff      = lj_parameters["cutoff"],
1527                                                    offset      = lj_parameters["offset"],)
1528            index = len(self.df)
1529            for label in [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']:
1530                self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond")
1531            self.df.at [index,'name']= f'{particle_name1}-{particle_name2}'
1532            self.df.at [index,'bond_object'] = bond_object
1533            self.df.at [index,'l0'] = l0
1534            self.add_value_to_df(index=index,
1535                                    key=('pmb_type',''),
1536                                    new_value='bond')
1537            self.add_value_to_df(index=index,
1538                                    key=('parameters_of_the_potential',''),
1539                                    new_value=bond_object.get_params(),
1540                                    non_standard_value=True)
1541        self.df.fillna(pd.NA, inplace=True)
1542        return

Defines a pmb object of type bond in pymbe.df.

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.
Note:

Currently, only HARMONIC and FENE bonds are supported.

For a HARMONIC bond the dictionary must contain the following parameters:

- k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
using the `pmb.units` UnitRegistry.
- r_0 (`obj`)    : 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 (`obj`): 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, epsilon=None, sigma=None, cutoff=None, offset=None):
1545    def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None):
1546        """
1547        Asigns `bond` in `pmb.df` as the default bond.
1548        The LJ parameters can be optionally provided to calculate the initial bond length
1549
1550        Args:
1551            bond_type(`str`): label to identify the potential to model the bond.
1552            bond_parameters(`dict`): parameters of the potential of the bond.
1553            sigma(`float`, optional): LJ sigma of the interaction between the particles.
1554            epsilon(`float`, optional): LJ epsilon for the interaction between the particles.
1555            cutoff(`float`, optional): cutoff-radius of the LJ interaction.
1556            offset(`float`, optional): offset of the LJ interaction.
1557
1558        Note:
1559            - Currently, only harmonic and FENE bonds are supported. 
1560        """
1561        
1562        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1563
1564        if epsilon is None:
1565            epsilon=1*self.units('reduced_energy')
1566        if sigma is None:
1567            sigma=1*self.units('reduced_length')
1568        if cutoff is None:
1569            cutoff=2**(1.0/6.0)*self.units('reduced_length')
1570        if offset is None:
1571            offset=0*self.units('reduced_length')
1572        l0 = self.calculate_initial_bond_length(bond_object = bond_object, 
1573                                                bond_type   = bond_type, 
1574                                                epsilon     = epsilon,
1575                                                sigma       = sigma,
1576                                                cutoff      = cutoff,
1577                                                offset      = offset)
1578
1579        if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'):
1580            return
1581        if len(self.df.index) != 0:
1582            index = max(self.df.index)+1
1583        else:
1584            index = 0
1585        self.df.at [index,'name']        = 'default'
1586        self.df.at [index,'bond_object'] = bond_object
1587        self.df.at [index,'l0']          = l0
1588        self.add_value_to_df(index       = index,
1589                            key          = ('pmb_type',''),
1590                            new_value    = 'bond')
1591        self.add_value_to_df(index       = index,
1592                            key          = ('parameters_of_the_potential',''),
1593                            new_value    = bond_object.get_params(),
1594                            non_standard_value=True)
1595        self.df.fillna(pd.NA, inplace=True)
1596        return

Asigns bond in pmb.df as the default bond. The LJ parameters can be optionally provided to calculate the initial bond length

Arguments:
  • bond_type(str): label to identify the potential to model the bond.
  • bond_parameters(dict): parameters of the potential of the bond.
  • sigma(float, optional): LJ sigma of the interaction between the particles.
  • epsilon(float, optional): LJ epsilon for the interaction between the particles.
  • cutoff(float, optional): cutoff-radius of the LJ interaction.
  • offset(float, optional): offset of the LJ interaction.
Note:
  • Currently, only harmonic and FENE bonds are supported.
def define_hydrogel(self, name, node_map, chain_map):
1598    def define_hydrogel(self, name, node_map, chain_map):
1599        """
1600        Defines a pyMBE object of type `hydrogel` in `pymbe.df`.
1601
1602        Args:
1603            name(`str`): Unique label that identifies the `hydrogel`.
1604            node_map(`list of ict`): [{"particle_name": , "lattice_index": }, ... ]
1605            chain_map(`list of dict`): [{"node_start": , "node_end": , "residue_list": , ... ]
1606        """
1607        node_indices = {tuple(entry['lattice_index']) for entry in node_map}
1608        diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices}
1609        if node_indices != diamond_indices:
1610            raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ")
1611        
1612        chain_map_connectivity = set()
1613        for entry in chain_map:
1614            start = self.lattice_builder.node_labels[entry['node_start']]
1615            end = self.lattice_builder.node_labels[entry['node_end']]
1616            chain_map_connectivity.add((start,end))
1617
1618        if self.lattice_builder.lattice.connectivity != chain_map_connectivity:
1619            raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs")
1620
1621        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='hydrogel'):
1622            return
1623
1624        index = len(self.df)
1625        self.df.at [index, "name"] = name
1626        self.df.at [index, "pmb_type"] = "hydrogel"
1627        self.add_value_to_df(index = index,
1628                        key = ('node_map',''),
1629                        new_value = node_map,
1630                        non_standard_value=True)
1631        self.add_value_to_df(index = index,
1632                        key = ('chain_map',''),
1633                        new_value = chain_map,
1634                        non_standard_value=True)
1635        for chain_label in chain_map:
1636            node_start = chain_label["node_start"]
1637            node_end = chain_label["node_end"]
1638            residue_list = chain_label['residue_list']
1639            # Molecule name
1640            molecule_name = "chain_"+node_start+"_"+node_end
1641            self.define_molecule(name=molecule_name, residue_list=residue_list)
1642        return

Defines a pyMBE object of type hydrogel in pymbe.df.

Arguments:
  • name(str): Unique label that identifies the hydrogel.
  • node_map(list of ict): [{"particle_name": , "lattice_index": }, ... ]
  • chain_map(list of dict): [{"node_start": , "node_end": , "residue_list": , ... ]
def define_molecule(self, name, residue_list):
1644    def define_molecule(self, name, residue_list):
1645        """
1646        Defines a pyMBE object of type `molecule` in `pymbe.df`.
1647
1648        Args:
1649            name(`str`): Unique label that identifies the `molecule`.
1650            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.  
1651        """
1652        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'):
1653            return
1654        index = len(self.df)
1655        self.df.at [index,'name'] = name
1656        self.df.at [index,'pmb_type'] = 'molecule'
1657        self.df.at [index,('residue_list','')] = residue_list
1658        self.df.fillna(pd.NA, inplace=True)
1659        return

Defines a pyMBE object of type molecule in pymbe.df.

Arguments:
  • name(str): Unique label that identifies the molecule.
  • residue_list(list of str): List of the names of the residues in the sequence of the molecule.
def define_particle_entry_in_df(self, name):
1661    def define_particle_entry_in_df(self,name):
1662        """
1663        Defines a particle entry in pmb.df.
1664
1665        Args:
1666            name(`str`): Unique label that identifies this particle type.
1667
1668        Returns:
1669            index(`int`): Index of the particle in pmb.df  
1670        """
1671
1672        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'):
1673            index = self.df[self.df['name']==name].index[0]                                   
1674        else:
1675            index = len(self.df)
1676            self.df.at [index, 'name'] = name
1677            self.df.at [index,'pmb_type'] = 'particle'
1678        self.df.fillna(pd.NA, inplace=True)
1679        return index

Defines a particle entry in pmb.df.

Arguments:
  • name(str): Unique label that identifies this particle type.
Returns:

index(int): Index of the particle in pmb.df

def define_particle( self, name, z=0, acidity=<NA>, pka=<NA>, sigma=<NA>, epsilon=<NA>, cutoff=<NA>, offset=<NA>, overwrite=False):
1681    def define_particle(self, name, z=0, acidity=pd.NA, pka=pd.NA, sigma=pd.NA, epsilon=pd.NA, cutoff=pd.NA, offset=pd.NA,overwrite=False):
1682        """
1683        Defines the properties of a particle object.
1684
1685        Args:
1686            name(`str`): Unique label that identifies this particle type.  
1687            z(`int`, optional): Permanent charge number of this particle type. Defaults to 0.
1688            acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to pd.NA.
1689            pka(`float`, optional): If `particle` is an acid or a base, it defines its  pka-value. Defaults to pd.NA.
1690            sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1691            cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1692            offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA.
1693            epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA.
1694            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
1695
1696        Note:
1697            - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units.
1698            - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units.
1699            - `cutoff` defaults to `2**(1./6.) reduced_length`. 
1700            - `offset` defaults to 0.
1701            - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions.
1702            - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`.
1703        """ 
1704        index=self.define_particle_entry_in_df(name=name)
1705        
1706        # If `cutoff` and `offset` are not defined, default them to the following values
1707        if pd.isna(cutoff):
1708            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1709        if pd.isna(offset):
1710            offset=self.units.Quantity(0, "reduced_length")
1711
1712        # Define LJ parameters
1713        parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"},
1714                                        "cutoff":{"value": cutoff, "dimensionality": "[length]"},
1715                                        "offset":{"value": offset, "dimensionality": "[length]"},
1716                                        "epsilon":{"value": epsilon, "dimensionality": "[energy]"},}
1717
1718        for parameter_key in parameters_with_dimensionality.keys():
1719            if not pd.isna(parameters_with_dimensionality[parameter_key]["value"]):
1720                self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 
1721                                          expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"])
1722                self.add_value_to_df(key=(parameter_key,''),
1723                                    index=index,
1724                                    new_value=parameters_with_dimensionality[parameter_key]["value"],
1725                                    overwrite=overwrite)
1726
1727        # Define particle acid/base properties
1728        self.set_particle_acidity(name=name, 
1729                                acidity=acidity, 
1730                                default_charge_number=z, 
1731                                pka=pka,
1732                                overwrite=overwrite)
1733        self.df.fillna(pd.NA, inplace=True)
1734        return 

Defines the properties of a particle object.

Arguments:
  • name(str): Unique label that identifies this particle type.
  • 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.
  • sigma(pint.Quantity, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. 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.
  • epsilon(pint.Quantity, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA.
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
Note:
  • 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.
  • The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions.
  • For more information on sigma, epsilon, cutoff and offset check pmb.setup_lj_interactions().
def define_particles(self, parameters, overwrite=False):
1736    def define_particles(self, parameters, overwrite=False):
1737        '''
1738        Defines a particle object in pyMBE for each particle name in `particle_names`
1739
1740        Args:
1741            parameters(`dict`):  dictionary with the particle parameters. 
1742            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1743
1744        Note:
1745            - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},}
1746        '''
1747        if not parameters:
1748            return 0
1749        for particle_name in parameters.keys():
1750            parameters[particle_name]["overwrite"]=overwrite
1751            self.define_particle(**parameters[particle_name])
1752        return

Defines a particle object in pyMBE for each particle name in particle_names

Arguments:
  • parameters(dict): dictionary with the particle parameters.
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
Note:
  • parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},}
def define_peptide(self, name, sequence, model):
1754    def define_peptide(self, name, sequence, model):
1755        """
1756        Defines a pyMBE object of type `peptide` in the `pymbe.df`.
1757
1758        Args:
1759            name (`str`): Unique label that identifies the `peptide`.
1760            sequence (`string`): Sequence of the `peptide`.
1761            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1762        """
1763        if self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'):
1764            return
1765        valid_keys = ['1beadAA','2beadAA']
1766        if model not in valid_keys:
1767            raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
1768        clean_sequence = self.protein_sequence_parser(sequence=sequence)    
1769        residue_list = self.define_AA_residues(sequence=clean_sequence,
1770                                                model=model)
1771        self.define_molecule(name = name, residue_list=residue_list)
1772        index = self.df.loc[self.df['name'] == name].index.item() 
1773        self.df.at [index,'model'] = model
1774        self.df.at [index,('sequence','')] = clean_sequence
1775        self.df.fillna(pd.NA, inplace=True)

Defines a pyMBE object of type peptide in the pymbe.df.

Arguments:
  • name (str): Unique label that identifies the peptide.
  • sequence (string): Sequence of the peptide.
  • model (string): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
def define_protein( self, name, model, topology_dict, lj_setup_mode='wca', overwrite=False):
1778    def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False):
1779        """
1780        Defines a globular protein pyMBE object  in `pymbe.df`.
1781
1782        Args:
1783            name (`str`): Unique label that identifies the protein.
1784            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1785            topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value}
1786            lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca".
1787            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1788
1789        Note:
1790            - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential.
1791        """
1792
1793        # Sanity checks
1794        valid_model_keys = ['1beadAA','2beadAA']
1795        valid_lj_setups = ["wca"]
1796        if model not in valid_model_keys:
1797            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
1798        if lj_setup_mode not in valid_lj_setups:
1799            raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}')
1800        if lj_setup_mode == "wca":
1801            sigma = 1*self.units.Quantity("reduced_length")
1802            epsilon = 1*self.units.Quantity("reduced_energy")
1803        part_dict={}
1804        sequence=[]
1805        metal_ions_charge_number_map=self.get_metal_ions_charge_number_map()
1806        for particle in topology_dict.keys():
1807            particle_name = re.split(r'\d+', particle)[0] 
1808            if particle_name not in part_dict.keys():
1809                if lj_setup_mode == "wca":
1810                    part_dict[particle_name]={"sigma": sigma,
1811                                        "offset": topology_dict[particle]['radius']*2-sigma,
1812                                        "epsilon": epsilon,
1813                                        "name": particle_name}
1814                if self.check_if_metal_ion(key=particle_name):
1815                    z=metal_ions_charge_number_map[particle_name]
1816                else:
1817                    z=0
1818                part_dict[particle_name]["z"]=z
1819            
1820            if self.check_aminoacid_key(key=particle_name):
1821                sequence.append(particle_name) 
1822            
1823        self.define_particles(parameters=part_dict,
1824                            overwrite=overwrite)
1825        residue_list = self.define_AA_residues(sequence=sequence, 
1826                                               model=model)
1827        index = len(self.df)
1828        self.df.at [index,'name'] = name
1829        self.df.at [index,'pmb_type'] = 'protein'
1830        self.df.at [index,'model'] = model
1831        self.df.at [index,('sequence','')] = sequence  
1832        self.df.at [index,('residue_list','')] = residue_list
1833        self.df.fillna(pd.NA, inplace=True)    
1834        return 

Defines a globular protein pyMBE object in pymbe.df.

Arguments:
  • name (str): Unique label that identifies the protein.
  • model (string): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
  • topology_dict (dict): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value}
  • lj_setup_mode(str): Key for the setup for the LJ potential. Defaults to "wca".
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
Note:
  • 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):
1836    def define_residue(self, name, central_bead, side_chains):
1837        """
1838        Defines a pyMBE object of type `residue` in `pymbe.df`.
1839
1840        Args:
1841            name(`str`): Unique label that identifies the `residue`.
1842            central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`.
1843            side_chains(`list` of `str`): List of `name`s of the pmb_objects to be placed as side_chains of the `residue`. Currently, only pmb_objects of type `particle`s or `residue`s are supported.
1844        """
1845        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'):
1846            return
1847        index = len(self.df)
1848        self.df.at [index, 'name'] = name
1849        self.df.at [index,'pmb_type'] = 'residue'
1850        self.df.at [index,'central_bead'] = central_bead
1851        self.df.at [index,('side_chains','')] = side_chains
1852        self.df.fillna(pd.NA, inplace=True)
1853        return 

Defines a pyMBE object of type residue in pymbe.df.

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 names of the pmb_objects to be placed as side_chains of the residue. Currently, only pmb_objects of type particles or residues are supported.
def delete_entries_in_df(self, entry_name):
1855    def delete_entries_in_df(self, entry_name):
1856        """
1857        Deletes entries with name `entry_name` from the DataFrame if it exists.
1858
1859        Args:
1860            entry_name (`str`): The name of the entry in the dataframe to delete.
1861
1862        """
1863        if entry_name in self.df["name"].values:
1864            self.df = self.df[self.df["name"] != entry_name].reset_index(drop=True)

Deletes entries with name entry_name from the DataFrame if it exists.

Arguments:
  • entry_name (str): The name of the entry in the dataframe to delete.
def destroy_pmb_object_in_system(self, name, espresso_system):
1866    def destroy_pmb_object_in_system(self, name, espresso_system):
1867        """
1868        Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 
1869
1870        Args:
1871            name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`.
1872            espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library.
1873
1874        Note:
1875            - If `name`  is a object_type=`particle`, only the matching particles that are not part of bigger objects (e.g. `residue`, `molecule`) will be destroyed. To destroy particles in such objects, destroy the bigger object instead.
1876        """
1877        allowed_objects = ['particle','residue','molecule']
1878        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1879        if pmb_type not in allowed_objects:
1880            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1881        if pmb_type == 'particle':
1882            particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 
1883                                                  & (self.df['molecule_id'].isna())]
1884            particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna())
1885                                                 & (self.df['molecule_id'].isna())].particle_id.tolist()
1886            for particle_id in particle_ids_list:
1887                espresso_system.part.by_id(particle_id).remove()
1888            self.df = self.df.drop(index=particles_index)
1889        if pmb_type == 'residue':
1890            residues_id = self.df.loc[self.df['name']== name].residue_id.to_list()
1891            for residue_id in residues_id:
1892                molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0]
1893                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1894                self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index)
1895                for particle_id in particle_ids_list:
1896                    espresso_system.part.by_id(particle_id).remove()
1897                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)    
1898        if pmb_type == 'molecule':
1899            molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list()
1900            for molecule_id in molecules_id:
1901                molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']=="molecule")].name.values[0]
1902                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1903                self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index)
1904                for particle_id in particle_ids_list:
1905                    espresso_system.part.by_id(particle_id).remove()   
1906                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)             
1907        
1908        self.df.reset_index(drop=True,inplace=True)
1909
1910        return

Destroys all particles associated with name in espresso_system amd removes the destroyed pmb_objects from pmb.df

Arguments:
  • name(str): Label of the pmb object to be destroyed. The pmb object must be defined in pymbe.df.
  • espresso_system(espressomd.system.System): Instance of a system class from espressomd library.
Note:
  • If name is a object_type=particle, only the matching particles that are not part of bigger objects (e.g. residue, molecule) will be destroyed. To destroy particles in such objects, destroy the bigger object instead.
def determine_reservoir_concentrations( self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
1912    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
1913        """
1914        Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration.
1915        To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an
1916        ion pair and the concentrations of the various species is solved numerically using a self-consistent approach.
1917        More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831).
1918        This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237).
1919
1920        Args:
1921            pH_res('float'): pH-value in the reservoir.
1922            c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
1923            activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
1924
1925        Returns:
1926            cH_res('pint.Quantity'): Concentration of H+ ions.
1927            cOH_res('pint.Quantity'): Concentration of OH- ions.
1928            cNa_res('pint.Quantity'): Concentration of Na+ ions.
1929            cCl_res('pint.Quantity'): Concentration of Cl- ions.
1930        """
1931
1932        self_consistent_run = 0
1933        cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+
1934
1935        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
1936            #Calculate and initial guess for the concentrations of various species based on ideal gas estimate
1937            cOH_res = self.Kw / cH_res 
1938            cNa_res = None
1939            cCl_res = None
1940            if cOH_res>=cH_res:
1941                #adjust the concentration of sodium if there is excess OH- in the reservoir:
1942                cNa_res = c_salt_res + (cOH_res-cH_res)
1943                cCl_res = c_salt_res
1944            else:
1945                # adjust the concentration of chloride if there is excess H+ in the reservoir
1946                cCl_res = c_salt_res + (cH_res-cOH_res)
1947                cNa_res = c_salt_res
1948                
1949            def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res):
1950                nonlocal max_number_sc_runs, self_consistent_run
1951                if self_consistent_run<max_number_sc_runs:
1952                    self_consistent_run+=1
1953                    ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1954                    cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
1955                    if cOH_res>=cH_res:
1956                        #adjust the concentration of sodium if there is excess OH- in the reservoir:
1957                        cNa_res = c_salt_res + (cOH_res-cH_res)
1958                        cCl_res = c_salt_res
1959                    else:
1960                        # adjust the concentration of chloride if there is excess H+ in the reservoir
1961                        cCl_res = c_salt_res + (cH_res-cOH_res)
1962                        cNa_res = c_salt_res
1963                    return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1964                else:
1965                    return cH_res, cOH_res, cNa_res, cCl_res
1966            return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1967
1968        cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1969        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1970        determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1971
1972        while abs(determined_pH-pH_res)>1e-6:
1973            if determined_pH > pH_res:
1974                cH_res *= 1.005
1975            else:
1976                cH_res /= 1.003
1977            cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1978            ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1979            determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1980            self_consistent_run=0
1981
1982        return cH_res, cOH_res, cNa_res, cCl_res

Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration. To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an ion pair and the concentrations of the various species is solved numerically using a self-consistent approach. More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831). This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237).

Arguments:
  • pH_res('float'): 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', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
Returns:

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.

def enable_motion_of_rigid_object(self, name, espresso_system):
1984    def enable_motion_of_rigid_object(self, name, espresso_system):
1985        '''
1986        Enables the motion of the rigid object `name` in the `espresso_system`.
1987
1988        Args:
1989            name(`str`): Label of the object.
1990            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1991
1992        Note:
1993            - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"].
1994        '''
1995        print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
1996        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1997        if pmb_type != 'protein':
1998            raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein')
1999        molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list()
2000        for molecule_id in molecule_ids_list:    
2001            particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list()
2002            center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system)
2003            rigid_object_center = espresso_system.part.add(pos=center_of_mass,
2004                                                           rotation=[True,True,True], 
2005                                                           type=self.propose_unused_type())
2006            
2007            rigid_object_center.mass = len(particle_ids_list)
2008            momI = 0
2009            for pid in particle_ids_list:
2010                momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
2011
2012            rigid_object_center.rinertia = np.ones(3) * momI
2013            
2014            for particle_id in particle_ids_list:
2015                pid = espresso_system.part.by_id(particle_id)
2016                pid.vs_auto_relate_to(rigid_object_center.id)
2017        return

Enables the motion of the rigid object name in the espresso_system.

Arguments:
  • name(str): Label of the object.
  • espresso_system(espressomd.system.System): Instance of a system object from the espressomd library.
Note:
  • It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"].
def filter_df(self, pmb_type):
2019    def filter_df(self, pmb_type):
2020        """
2021        Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns.
2022        
2023        Args:
2024            pmb_type(`str`): pmb_object_type to filter in `pmb.df`.
2025
2026        Returns:
2027            pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`.
2028        """
2029        pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type]
2030        pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1)
2031        return pmb_type_df

Filters pmb.df and returns a sub-set of it containing only rows with pmb_object_type=pmb_type and non-NaN columns.

Arguments:
  • pmb_type(str): pmb_object_type to filter in pmb.df.
Returns:

pmb_type_df(Pandas.Dataframe): filtered pmb.df.

def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False):
2033    def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False):
2034        """
2035        Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2036        
2037        Args:
2038            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2039            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2040            use_default_bond(`bool`, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to 'False'. 
2041
2042        Returns:
2043            bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` if a matching bond exists
2044
2045        Note:
2046            - If `use_default_bond`=`True`, it returns "default" if no key is found.
2047        """
2048        bond_keys = [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']
2049        bond_defined=False
2050        for bond_key in bond_keys:
2051            if bond_key in self.df["name"].values:
2052                bond_defined=True
2053                correct_key=bond_key
2054                break
2055        if bond_defined:
2056            return correct_key
2057        elif use_default_bond:
2058            return 'default'
2059        else:
2060            return 

Searches for the name of the bond between particle_name1 and particle_name2 in pymbe.df and returns it.

Arguments:
  • particle_name1(str): label of the type of the first particle type of the bonded particles.
  • particle_name2(str): label of the type of the second particle type of the bonded particles.
  • use_default_bond(bool, optional): If it is activated, the "default" bond is returned if no bond is found between particle_name1 and particle_name2. Defaults to 'False'.
Returns:

bond_key (str): name of the bond between particle_name1 and particle_name2 if a matching bond exists

Note:
  • If use_default_bond=True, it returns "default" if no key is found.
def find_value_from_es_type(self, es_type, column_name):
2062    def find_value_from_es_type(self, es_type, column_name):
2063        """
2064        Finds a value in `pmb.df` for a `column_name` and `es_type` pair.
2065
2066        Args:
2067            es_type(`int`): value of the espresso type
2068            column_name(`str`): name of the column in `pymbe.df`
2069
2070        Returns:
2071            Value in `pymbe.df` matching  `column_name` and `es_type`
2072        """
2073        idx = pd.IndexSlice
2074        for state in ['state_one', 'state_two']:            
2075            index = self.df.loc[self.df[(state, 'es_type')] == es_type].index
2076            if len(index) > 0:
2077                if column_name == 'label':
2078                    label = self.df.loc[idx[index[0]], idx[(state,column_name)]]
2079                    return label
2080                else: 
2081                    column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]]
2082                    return column_name_value
2083        return None

Finds a value in pmb.df for a column_name and es_type pair.

Arguments:
  • es_type(int): value of the espresso type
  • column_name(str): name of the column in pymbe.df
Returns:

Value in pymbe.df matching column_name and es_type

def format_node(self, node_list):
2085    def format_node(self, node_list):
2086        return "[" + " ".join(map(str, node_list)) + "]"
def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2089    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
2090        """
2091        Generates coordinates outside a sphere centered at `center`.
2092
2093        Args:
2094            center(`lst`): Coordinates of the center of the sphere.
2095            radius(`float`): Radius of the sphere.
2096            max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates.
2097            n_samples(`int`): Number of sample points.
2098
2099        Returns:
2100            coord_list(`lst`): Coordinates of the sample points.
2101        """
2102        if not radius > 0: 
2103            raise ValueError (f'The value of {radius} must be a positive value')
2104        if not radius < max_dist:
2105            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
2106        coord_list = []
2107        counter = 0
2108        while counter<n_samples:
2109            coord = self.generate_random_points_in_a_sphere(center=center, 
2110                                            radius=max_dist,
2111                                            n_samples=1)[0]
2112            if np.linalg.norm(coord-np.asarray(center))>=radius:
2113                coord_list.append (coord)
2114                counter += 1
2115        return coord_list

Generates coordinates outside a sphere centered at center.

Arguments:
  • center(lst): Coordinates of the center of the sphere.
  • radius(float): Radius of the sphere.
  • max_dist(float): Maximum distance from the center of the spahre to generate coordinates.
  • n_samples(int): Number of sample points.
Returns:

coord_list(lst): Coordinates of the sample points.

def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2117    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
2118        """
2119        Uniformly samples points from a hypersphere. If on_surface is set to True, the points are
2120        uniformly sampled from the surface of the hypersphere.
2121        
2122        Args:
2123            center(`lst`): Array with the coordinates of the center of the spheres.
2124            radius(`float`): Radius of the sphere.
2125            n_samples(`int`): Number of sample points to generate inside the sphere.
2126            on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere.
2127
2128        Returns:
2129            samples(`list`): Coordinates of the sample points inside the hypersphere.
2130        """
2131        # initial values
2132        center=np.array(center)
2133        d = center.shape[0]
2134        # sample n_samples points in d dimensions from a standard normal distribution
2135        samples = self.rng.normal(size=(n_samples, d))
2136        # make the samples lie on the surface of the unit hypersphere
2137        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
2138        samples /= normalize_radii
2139        if not on_surface:
2140            # make the samples lie inside the hypersphere with the correct density
2141            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
2142            new_radii = np.power(uniform_points, 1/d)
2143            samples *= new_radii
2144        # scale the points to have the correct radius and center
2145        samples = samples * radius + center
2146        return samples 

Uniformly samples points from a hypersphere. If on_surface is set to True, the points are uniformly sampled from the surface of the hypersphere.

Arguments:
  • center(lst): Array with the coordinates of the center of the spheres.
  • radius(float): Radius of the sphere.
  • n_samples(int): Number of sample points to generate inside the sphere.
  • on_surface (bool, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere.
Returns:

samples(list): Coordinates of the sample points inside the hypersphere.

def generate_trial_perpendicular_vector(self, vector, magnitude):
2148    def generate_trial_perpendicular_vector(self,vector,magnitude):
2149        """
2150        Generates an orthogonal vector to the input `vector`.
2151
2152        Args:
2153            vector(`lst`): arbitrary vector.
2154            magnitude(`float`): magnitude of the orthogonal vector.
2155            
2156        Returns:
2157            (`lst`): Orthogonal vector with the same magnitude as the input vector.
2158        """ 
2159        np_vec = np.array(vector) 
2160        np_vec /= np.linalg.norm(np_vec) 
2161        if np.all(np_vec == 0):
2162            raise ValueError('Zero vector')
2163        # Generate a random vector 
2164        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
2165                                                                center=[0,0,0],
2166                                                                n_samples=1, 
2167                                                                on_surface=True)[0]
2168        # Project the random vector onto the input vector and subtract the projection
2169        projection = np.dot(random_vector, np_vec) * np_vec
2170        perpendicular_vector = random_vector - projection
2171        # Normalize the perpendicular vector to have the same magnitude as the input vector
2172        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
2173        return perpendicular_vector*magnitude

Generates an orthogonal vector to the input vector.

Arguments:
  • vector(lst): arbitrary vector.
  • magnitude(float): magnitude of the orthogonal vector.
Returns:

(lst): Orthogonal vector with the same magnitude as the input vector.

def get_bond_length( self, particle_name1, particle_name2, hard_check=False, use_default_bond=False):
2175    def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2176        """
2177        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length.
2178        If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead.
2179        If no bond is found, it prints a message and it does not return anything. If `hard_check` is activated, the code stops if no bond is found.
2180
2181        Args:
2182            particle_name1(str): label of the type of the first particle type of the bonded particles.
2183            particle_name2(str): label of the type of the second particle type of the bonded particles.
2184            hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2185            use_default_bond(bool, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to False. 
2186
2187        Returns:
2188            l0(`pint.Quantity`): bond length
2189        
2190        Note:
2191            - If `use_default_bond`=True and no bond is defined between `particle_name1` and `particle_name2`, it returns the default bond defined in `pmb.df`.
2192            - If `hard_check`=`True` stops the code when no bond is found.
2193        """
2194        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2195                                    particle_name2=particle_name2, 
2196                                    use_default_bond=use_default_bond)
2197        if bond_key:
2198            return self.df[self.df['name']==bond_key].l0.values[0]
2199        else:
2200            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2201            if hard_check:
2202                sys.exit(1)
2203            else:
2204                return

Searches for bonds between the particle types given by particle_name1 and particle_name2 in pymbe.df and returns the initial bond length. If use_default_bond is activated and a "default" bond is defined, returns the length of that default bond instead. If no bond is found, it prints a message and it does not return anything. If hard_check is activated, the code stops if no bond is found.

Arguments:
  • particle_name1(str): label of the type of the first particle type of the bonded particles.
  • particle_name2(str): label of the type of the second particle type of the bonded particles.
  • hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False.
  • use_default_bond(bool, optional): If it is activated, the "default" bond is returned if no bond is found between particle_name1 and particle_name2. Defaults to False.
Returns:

l0(pint.Quantity): bond length

Note:
  • If use_default_bond=True and no bond is defined between particle_name1 and particle_name2, it returns the default bond defined in pmb.df.
  • If hard_check=True stops the code when no bond is found.
def get_charge_number_map(self):
2206    def get_charge_number_map(self):
2207        '''
2208        Gets the charge number of each `espresso_type` in `pymbe.df`.
2209        
2210        Returns:
2211            charge_number_map(`dict`): {espresso_type: z}.
2212        '''
2213        if self.df.state_one['es_type'].isnull().values.any():         
2214            df_state_one = self.df.state_one.dropna()     
2215            df_state_two = self.df.state_two.dropna()  
2216        else:    
2217            df_state_one = self.df.state_one
2218            if self.df.state_two['es_type'].isnull().values.any():
2219                df_state_two = self.df.state_two.dropna()   
2220            else:
2221                df_state_two = self.df.state_two
2222        state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values)
2223        state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values)
2224        charge_number_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2225        return charge_number_map

Gets the charge number of each espresso_type in pymbe.df.

Returns:

charge_number_map(dict): {espresso_type: z}.

def get_lj_parameters( self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2227    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2228        """
2229        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2230        `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule.
2231
2232        Args:
2233            particle_name1 (str): label of the type of the first particle type
2234            particle_name2 (str): label of the type of the second particle type
2235            combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2236
2237        Returns:
2238            {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2239
2240        Note:
2241            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
2242            - 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.
2243        """
2244        supported_combining_rules=["Lorentz-Berthelot"]
2245        lj_parameters_keys=["sigma","epsilon","offset","cutoff"]
2246        if combining_rule not in supported_combining_rules:
2247            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2248        lj_parameters={}
2249        for key in lj_parameters_keys:
2250            lj_parameters[key]=[]
2251        # Search the LJ parameters of the type pair
2252        for name in [particle_name1,particle_name2]:
2253            for key in lj_parameters_keys:
2254                lj_parameters[key].append(self.df[self.df.name == name][key].values[0])
2255        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2256        if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]):
2257            return {}
2258        # Apply combining rule
2259        if combining_rule == 'Lorentz-Berthelot':
2260            lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2
2261            lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2
2262            lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2
2263            lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1])
2264        return lj_parameters

Returns the Lennard-Jones parameters for the interaction between the particle types given by particle_name1 and particle_name2 in pymbe.df, 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:

{"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}

Note:
  • 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_metal_ions_charge_number_map(self):
2266    def get_metal_ions_charge_number_map(self):
2267        """
2268        Gets a map with the charge numbers of all the metal ions supported.
2269
2270        Returns:
2271            metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number}
2272
2273        """
2274        metal_charge_number_map = {"Ca": 2}
2275        return metal_charge_number_map

Gets a map with the charge numbers of all the metal ions supported.

Returns:

metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number}

def get_particle_id_map(self, object_name):
2277    def get_particle_id_map(self, object_name):
2278        '''
2279        Gets all the ids associated with the object with name `object_name` in `pmb.df`
2280
2281        Args:
2282            object_name(`str`): name of the object
2283        
2284        Returns:
2285            id_map(`dict`): dict of the structure {"all": [all_ids_with_object_name], "residue_map": {res_id: [particle_ids_in_res_id]}, "molecule_map": {mol_id: [particle_ids_in_mol_id]}, }
2286        '''
2287        object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0]
2288        valid_types = ["particle", "molecule", "residue", "protein"]
2289        if object_type not in valid_types:
2290            raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}")
2291        id_list = []
2292        mol_map = {}
2293        res_map = {}
2294        def do_res_map(res_ids):
2295            for res_id in res_ids:
2296                res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2297                res_map[res_id]=res_list
2298            return res_map
2299        if object_type in ['molecule', 'protein']:
2300            mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist()
2301            for mol_id in mol_ids:
2302                res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist())
2303                res_map=do_res_map(res_ids=res_ids)    
2304                mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2305                id_list+=mol_list
2306                mol_map[mol_id]=mol_list
2307        elif object_type == 'residue':     
2308            res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist()
2309            res_map=do_res_map(res_ids=res_ids)
2310            id_list=[]
2311            for res_id_list in res_map.values():
2312                id_list+=res_id_list
2313        elif object_type == 'particle':
2314            id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist()
2315        return {"all": id_list, "molecule_map": mol_map, "residue_map": res_map}

Gets all the ids associated with the object with name object_name in pmb.df

Arguments:
  • object_name(str): name of the object
Returns:

id_map(dict): dict of the structure {"all": [all_ids_with_object_name], "residue_map": {res_id: [particle_ids_in_res_id]}, "molecule_map": {mol_id: [particle_ids_in_mol_id]}, }

def get_pka_set(self):
2317    def get_pka_set(self):
2318        '''
2319        Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df`
2320        
2321        Returns:
2322            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
2323        '''
2324        titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna()
2325        pka_set = {}
2326        for index in titratables_AA_df.name.keys():
2327            name = titratables_AA_df.name[index]
2328            pka_value = titratables_AA_df.pka[index]
2329            acidity = titratables_AA_df.acidity[index]   
2330            pka_set[name] = {'pka_value':pka_value,'acidity':acidity}
2331        return pka_set 

Gets the pka-values and acidities of the particles with acid/base properties in pmb.df

Returns:

pka_set(dict): {"name" : {"pka_value": pka, "acidity": acidity}}

def get_radius_map(self, dimensionless=True):
2333    def get_radius_map(self, dimensionless=True):
2334        '''
2335        Gets the effective radius of each `espresso type` in `pmb.df`. 
2336
2337        Args:
2338            dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False.
2339        
2340        Returns:
2341            radius_map(`dict`): {espresso_type: radius}.
2342
2343        Note:
2344            The radius corresponds to (sigma+offset)/2
2345        '''
2346        df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates()
2347        df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates()
2348        state_one = pd.Series((df_state_one.sigma.values+df_state_one.offset.values)/2.0,index=df_state_one.state_one.es_type.values)
2349        state_two = pd.Series((df_state_two.sigma.values+df_state_two.offset.values)/2.0,index=df_state_two.state_two.es_type.values)
2350        radius_map  = pd.concat([state_one,state_two],axis=0).to_dict()  
2351        if dimensionless:
2352            for key in radius_map:
2353                radius_map[key] = radius_map[key].magnitude
2354        return radius_map

Gets the effective radius of each espresso type in pmb.df.

Arguments:
  • dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False.
Returns:

radius_map(dict): {espresso_type: radius}.

Note:

The radius corresponds to (sigma+offset)/2

def get_reduced_units(self):
2356    def get_reduced_units(self):
2357        """
2358        Returns the  current set of reduced units defined in pyMBE.
2359
2360        Returns:
2361            reduced_units_text(`str`): text with information about the current set of reduced units.
2362
2363        """
2364        unit_length=self.units.Quantity(1,'reduced_length')
2365        unit_energy=self.units.Quantity(1,'reduced_energy')
2366        unit_charge=self.units.Quantity(1,'reduced_charge')
2367        reduced_units_text = "\n".join(["Current set of reduced units:",
2368                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2369                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2370                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2371                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2372                                        ])   
2373        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_resource(self, path):
2375    def get_resource(self, path):
2376        '''
2377        Locate a file resource of the pyMBE package.
2378
2379        Args:
2380            path(`str`): Relative path to the resource
2381
2382        Returns:
2383            path(`str`): Absolute path to the resource
2384
2385        '''
2386        import os
2387        return os.path.join(os.path.dirname(__file__), path)

Locate a file resource of the pyMBE package.

Arguments:
  • path(str): Relative path to the resource
Returns:

path(str): Absolute path to the resource

def get_type_map(self):
2389    def get_type_map(self):
2390        """
2391        Gets all different espresso types assigned to particles  in `pmb.df`.
2392        
2393        Returns:
2394            type_map(`dict`): {"name": espresso_type}.
2395        """
2396        df_state_one = self.df.state_one.dropna(how='all')     
2397        df_state_two = self.df.state_two.dropna(how='all')  
2398        state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label)
2399        state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label)
2400        type_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2401        return type_map

Gets all different espresso types assigned to particles in pmb.df.

Returns:

type_map(dict): {"name": espresso_type}.

def initialize_lattice_builder(self, diamond_lattice):
2403    def initialize_lattice_builder(self, diamond_lattice):
2404        """
2405        Initialize the lattice builder with the DiamondLattice object.
2406
2407        Args:
2408            diamond_lattice(`DiamondLattice`): DiamondLattice object from the `lib/lattice` module to be used in the LatticeBuilder.
2409        """
2410        if not isinstance(diamond_lattice, DiamondLattice):
2411            raise TypeError("Currently only DiamondLattice objects are supported.")
2412        self.lattice_builder = LatticeBuilder(lattice=diamond_lattice)
2413        logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}")
2414        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_interaction_parameters(self, filename, overwrite=False):
2417    def load_interaction_parameters(self, filename, overwrite=False):
2418        """
2419        Loads the interaction parameters stored in `filename` into `pmb.df`
2420        
2421        Args:
2422            filename(`str`): name of the file to be read
2423            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2424        """
2425        without_units = ['z','es_type']
2426        with_units = ['sigma','epsilon','offset','cutoff']
2427
2428        with open(filename, 'r') as f:
2429            interaction_data = json.load(f)
2430            interaction_parameter_set = interaction_data["data"]
2431
2432        for key in interaction_parameter_set:
2433            param_dict=interaction_parameter_set[key]
2434            object_type=param_dict.pop('object_type')
2435            if object_type == 'particle':
2436                not_required_attributes={}    
2437                for not_required_key in without_units+with_units:
2438                    if not_required_key in param_dict.keys():
2439                        if not_required_key in with_units:
2440                            not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key))
2441                        elif not_required_key in without_units:
2442                            not_required_attributes[not_required_key]=param_dict.pop(not_required_key)
2443                    else:
2444                        not_required_attributes[not_required_key]=pd.NA
2445                self.define_particle(name=param_dict.pop('name'),
2446                                z=not_required_attributes.pop('z'),
2447                                sigma=not_required_attributes.pop('sigma'),
2448                                offset=not_required_attributes.pop('offset'),
2449                                cutoff=not_required_attributes.pop('cutoff'),
2450                                epsilon=not_required_attributes.pop('epsilon'),
2451                                overwrite=overwrite)
2452            elif object_type == 'residue':
2453                self.define_residue(**param_dict)
2454            elif object_type == 'molecule':
2455                self.define_molecule(**param_dict)
2456            elif object_type == 'peptide':
2457                self.define_peptide(**param_dict)
2458            elif object_type == 'bond':
2459                particle_pairs = param_dict.pop('particle_pairs')
2460                bond_parameters = param_dict.pop('bond_parameters')
2461                bond_type = param_dict.pop('bond_type')
2462                if bond_type == 'harmonic':
2463                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2464                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2465                    bond = {'r_0'    : r_0,
2466                            'k'      : k,
2467                            }
2468
2469                elif bond_type == 'FENE':
2470                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2471                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2472                    d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max'))
2473                    bond =  {'r_0'    : r_0,
2474                             'k'      : k,
2475                             'd_r_max': d_r_max,
2476                             }
2477                else:
2478                    raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds")
2479                self.define_bond(bond_type=bond_type, 
2480                                bond_parameters=bond, 
2481                                particle_pairs=particle_pairs)
2482            else:
2483                raise ValueError(object_type+' is not a known pmb object type')
2484            
2485        return

Loads the interaction parameters stored in filename into pmb.df

Arguments:
  • filename(str): name of the file to be read
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
def load_pka_set(self, filename, overwrite=True):
2487    def load_pka_set(self, filename, overwrite=True):
2488        """
2489        Loads the pka_set stored in `filename` into `pmb.df`.
2490        
2491        Args:
2492            filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}.
2493            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 
2494        """
2495        with open(filename, 'r') as f:
2496            pka_data = json.load(f)
2497            pka_set = pka_data["data"]
2498
2499        self.check_pka_set(pka_set=pka_set)
2500
2501        for key in pka_set:
2502            acidity = pka_set[key]['acidity']
2503            pka_value = pka_set[key]['pka_value']
2504            self.set_particle_acidity(name=key, 
2505                                      acidity=acidity, 
2506                                      pka=pka_value, 
2507                                      overwrite=overwrite)
2508        return

Loads the pka_set stored in filename into pmb.df.

Arguments:
  • filename(str): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}.
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True.
def propose_unused_type(self):
2511    def propose_unused_type(self):
2512        """
2513        Searches in `pmb.df` all the different particle types defined and returns a new one.
2514
2515        Returns:
2516            unused_type(`int`): unused particle type
2517        """
2518        type_map = self.get_type_map()
2519        if not type_map:  
2520            unused_type = 0
2521        else:
2522            valid_values = [v for v in type_map.values() if pd.notna(v)]  # Filter out pd.NA values
2523            unused_type = max(valid_values) + 1 if valid_values else 0  # Ensure max() doesn't fail if all values are NA
2524        return unused_type

Searches in pmb.df all the different particle types defined and returns a new one.

Returns:

unused_type(int): unused particle type

def protein_sequence_parser(self, sequence):
2526    def protein_sequence_parser(self, sequence):
2527        '''
2528        Parses `sequence` to the one letter code for amino acids.
2529        
2530        Args:
2531            sequence(`str` or `lst`): Sequence of the amino acid. 
2532
2533        Returns:
2534            clean_sequence(`lst`): `sequence` using the one letter code.
2535        
2536        Note:
2537            - Accepted formats for `sequence` are:
2538                - `lst` with one letter or three letter code of each aminoacid in each element
2539                - `str` with the sequence using the one letter code
2540                - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-"
2541        
2542        '''
2543        # Aminoacid key
2544        keys={"ALA": "A",
2545                "ARG": "R",
2546                "ASN": "N",
2547                "ASP": "D",
2548                "CYS": "C",
2549                "GLU": "E",
2550                "GLN": "Q",
2551                "GLY": "G",
2552                "HIS": "H",
2553                "ILE": "I",
2554                "LEU": "L",
2555                "LYS": "K",
2556                "MET": "M",
2557                "PHE": "F",
2558                "PRO": "P",
2559                "SER": "S",
2560                "THR": "T",
2561                "TRP": "W",
2562                "TYR": "Y",
2563                "VAL": "V",
2564                "PSER": "J",
2565                "PTHR": "U",
2566                "PTyr": "Z",
2567                "NH2": "n",
2568                "COOH": "c"}
2569        clean_sequence=[]
2570        if isinstance(sequence, str):
2571            if sequence.find("-") != -1:
2572                splited_sequence=sequence.split("-")
2573                for residue in splited_sequence:
2574                    if len(residue) == 1:
2575                        if residue in keys.values():
2576                            residue_ok=residue
2577                        else:
2578                            if residue.upper() in keys.values():
2579                                residue_ok=residue.upper()
2580                            else:
2581                                raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence")
2582                        clean_sequence.append(residue_ok)
2583                    else:
2584                        if residue in keys.keys():
2585                            clean_sequence.append(keys[residue])
2586                        else:
2587                            if residue.upper() in keys.keys():
2588                                clean_sequence.append(keys[residue.upper()])
2589                            else:
2590                                raise ValueError("Unknown  code for a residue: ", residue, " please review the input sequence")
2591            else:
2592                for residue in sequence:
2593                    if residue in keys.values():
2594                        residue_ok=residue
2595                    else:
2596                        if residue.upper() in keys.values():
2597                            residue_ok=residue.upper()
2598                        else:
2599                            raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence")
2600                    clean_sequence.append(residue_ok)
2601        if isinstance(sequence, list):
2602            for residue in sequence:
2603                if residue in keys.values():
2604                    residue_ok=residue
2605                else:
2606                    if residue.upper() in keys.values():
2607                        residue_ok=residue.upper()
2608                    elif (residue.upper() in keys.keys()):
2609                        residue_ok= keys[residue.upper()]
2610                    else:
2611                        raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence")
2612                clean_sequence.append(residue_ok)
2613        return clean_sequence

Parses sequence to the one letter code for amino acids.

Arguments:
  • sequence(str or lst): Sequence of the amino acid.
Returns:

clean_sequence(lst): sequence using the one letter code.

Note:
  • Accepted formats for sequence are:
    • lst with one letter or three letter code of each aminoacid in each element
    • str with the sequence using the one letter code
    • str with the squence using the three letter code, each aminoacid must be separated by a hyphen "-"
def read_pmb_df(self, filename):
2615    def read_pmb_df (self,filename):
2616        """
2617        Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE.
2618
2619        Args:
2620            filename(`str`): path to file.
2621
2622        Note:
2623            This function only accepts files with CSV format. 
2624        """
2625        
2626        if filename.rsplit(".", 1)[1] != "csv":
2627            raise ValueError("Only files with CSV format are supported")
2628        df = pd.read_csv (filename,header=[0, 1], index_col=0)
2629        columns_names = self.setup_df()
2630        
2631        multi_index = pd.MultiIndex.from_tuples(columns_names)
2632        df.columns = multi_index
2633        
2634        self.df = self.convert_columns_to_original_format(df)
2635        self.df.fillna(pd.NA, inplace=True)
2636        return self.df

Reads a pyMBE's Dataframe stored in filename and stores the information into pyMBE.

Arguments:
  • filename(str): path to file.
Note:

This function only accepts files with CSV format.

def read_protein_vtf_in_df(self, filename, unit_length=None):
2638    def read_protein_vtf_in_df (self,filename,unit_length=None):
2639        """
2640        Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`.
2641
2642        Args:
2643            filename(`str`): Path to the vtf file with the coarse-grained model.
2644            unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None.
2645
2646        Returns:
2647            topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}
2648
2649        Note:
2650            - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom.
2651        """
2652
2653        print (f'Loading protein coarse grain model file: {filename}')
2654
2655        coord_list = []
2656        particles_dict = {}
2657
2658        if unit_length is None:
2659            unit_length = 1 * self.units.angstrom 
2660
2661        with open (filename,'r') as protein_model:
2662            for line in protein_model :
2663                line_split = line.split()
2664                if line_split:
2665                    line_header = line_split[0]
2666                    if line_header == 'atom':
2667                        atom_id  = line_split[1]
2668                        atom_name = line_split[3]
2669                        atom_resname = line_split[5]
2670                        chain_id = line_split[9]
2671                        radius = float(line_split [11])*unit_length 
2672                        particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius]
2673                    elif line_header.isnumeric(): 
2674                        atom_coord = line_split[1:] 
2675                        atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord]
2676                        coord_list.append (atom_coord)
2677
2678        numbered_label = []
2679        i = 0   
2680        
2681        for atom_id in particles_dict.keys():
2682    
2683            if atom_id == 1:
2684                atom_name = particles_dict[atom_id][0]
2685                numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2686                numbered_label.append(numbered_name)
2687
2688            elif atom_id != 1: 
2689            
2690                if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]:
2691                    i += 1                    
2692                    count = 1
2693                    atom_name = particles_dict[atom_id][0]
2694                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2695                    numbered_label.append(numbered_name)
2696                    
2697                elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]:
2698                    if count == 2 or particles_dict[atom_id][1] == 'GLY':
2699                        i +=1  
2700                        count = 0
2701                    atom_name = particles_dict[atom_id][0]
2702                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2703                    numbered_label.append(numbered_name)
2704                    count +=1
2705
2706        topology_dict = {}
2707
2708        for i in range (0, len(numbered_label)):   
2709            topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] ,
2710                                                    'chain_id':numbered_label[i][1],
2711                                                    'radius':numbered_label[i][2] }
2712
2713        return topology_dict

Loads a coarse-grained protein model in a vtf file filename into the pmb.df and it labels it with name.

Arguments:
  • filename(str): Path to the vtf file with the coarse-grained model.
  • unit_length(obj): unit of length of the the coordinates in filename using the pyMBE UnitRegistry. Defaults to None.
Returns:

topology_dict(dict): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}

Note:
  • If no unit_length is provided, it is assumed that the coordinates are in Angstrom.
def search_bond( self, particle_name1, particle_name2, hard_check=False, use_default_bond=False):
2715    def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2716        """
2717        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2718        If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead.
2719        If no bond is found, it prints a message and it does not return anything. If `hard_check` is activated, the code stops if no bond is found.
2720
2721        Args:
2722            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2723            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2724            hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2725            use_default_bond(`bool`, optional): If it is activated, the "default" bond is returned if no bond is found between `particle_name1` and `particle_name2`. Defaults to False. 
2726
2727        Returns:
2728            bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library.
2729        
2730        Note:
2731            - If `use_default_bond`=True and no bond is defined between `particle_name1` and `particle_name2`, it returns the default bond defined in `pmb.df`.
2732            - If `hard_check`=`True` stops the code when no bond is found.
2733        """
2734        
2735        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2736                                    particle_name2=particle_name2, 
2737                                    use_default_bond=use_default_bond)
2738        if use_default_bond:
2739            if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'):
2740                raise ValueError(f"use_default_bond is set to {use_default_bond} but no default bond has been defined. Please define a default bond with pmb.define_default_bond")
2741        if bond_key:
2742            return self.df[self.df['name']==bond_key].bond_object.values[0]
2743        else:
2744            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2745            if hard_check:
2746                sys.exit(1)
2747            else:
2748                return

Searches for bonds between the particle types given by particle_name1 and particle_name2 in pymbe.df and returns it. If use_default_bond is activated and a "default" bond is defined, returns that default bond instead. If no bond is found, it prints a message and it does not return anything. If hard_check is activated, the code stops if no bond is found.

Arguments:
  • particle_name1(str): label of the type of the first particle type of the bonded particles.
  • particle_name2(str): label of the type of the second particle type of the bonded particles.
  • hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False.
  • use_default_bond(bool, optional): If it is activated, the "default" bond is returned if no bond is found between particle_name1 and particle_name2. Defaults to False.
Returns:

bond(espressomd.interactions.BondedInteractions): bond object from the espressomd library.

Note:
  • If use_default_bond=True and no bond is defined between particle_name1 and particle_name2, it returns the default bond defined in pmb.df.
  • If hard_check=True stops the code when no bond is found.
def search_particles_in_residue(self, residue_name):
2750    def search_particles_in_residue(self, residue_name):
2751        '''
2752        Searches for all particles in a given residue of name `residue_name`.
2753
2754        Args:
2755            residue_name (`str`): name of the residue to be searched
2756
2757        Returns:
2758            list_of_particles_in_residue (`lst`): list of the names of all particles in the residue
2759
2760        Note:
2761            - The function returns a name per particle in residue, i.e. if there are multiple particles with the same type `list_of_particles_in_residue` will have repeated items.
2762 
2763        '''
2764        index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 
2765        central_bead = self.df.at [index_residue, ('central_bead', '')]
2766        list_of_side_chains = self.df.at [index_residue, ('side_chains', '')]
2767
2768        list_of_particles_in_residue = []
2769        list_of_particles_in_residue.append(central_bead)
2770
2771        for side_chain in list_of_side_chains: 
2772            object_type = self.df[self.df['name']==side_chain].pmb_type.values[0]
2773
2774            if object_type == "residue":
2775                list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain)
2776                list_of_particles_in_residue += list_of_particles_in_side_chain_residue
2777            elif object_type == "particle":
2778                list_of_particles_in_residue.append(side_chain)
2779
2780        return list_of_particles_in_residue

Searches for all particles in a given residue of name residue_name.

Arguments:
  • residue_name (str): name of the residue to be searched
Returns:

list_of_particles_in_residue (lst): list of the names of all particles in the residue

Note:
  • The function returns a name per particle in residue, i.e. if there are multiple particles with the same type list_of_particles_in_residue will have repeated items.
def set_particle_acidity( self, name, acidity=<NA>, default_charge_number=0, pka=<NA>, overwrite=True):
2784    def set_particle_acidity(self, name, acidity=pd.NA, default_charge_number=0, pka=pd.NA, overwrite=True):
2785        """
2786        Sets the particle acidity including the charges in each of its possible states. 
2787
2788        Args:
2789            name(`str`): Unique label that identifies the `particle`. 
2790            acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None.
2791            default_charge_number (`int`): Charge number of the particle. Defaults to 0.
2792            pka(`float`, optional):  If `particle` is an acid or a base, it defines its pka-value. Defaults to pandas.NA.
2793            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2794     
2795        Note:
2796            - For particles with  `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 
2797        deprotonated states, respectively. 
2798            - For particles without an acidity `acidity = pandas.NA`, only `state_one` is defined.
2799            - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 
2800        """
2801        acidity_valid_keys = ['inert','acidic', 'basic']
2802        if not pd.isna(acidity):
2803            if acidity not in acidity_valid_keys:
2804                raise ValueError(f"Acidity {acidity} provided for particle name  {name} is not supproted. Valid keys are: {acidity_valid_keys}")
2805            if acidity in ['acidic', 'basic'] and pd.isna(pka):
2806                raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.")   
2807            if acidity == "inert":
2808                acidity = pd.NA
2809                print("Deprecation warning: the keyword 'inert' for acidity has been replaced by setting acidity = pd.NA. For backwards compatibility, acidity has been set to pd.NA. Support for `acidity = 'inert'` may be deprecated in future releases of pyMBE")
2810
2811        self.define_particle_entry_in_df(name=name)
2812        
2813        for index in self.df[self.df['name']==name].index:       
2814            if pka is not pd.NA:
2815                self.add_value_to_df(key=('pka',''),
2816                                    index=index,
2817                                    new_value=pka, 
2818                                    overwrite=overwrite)
2819            
2820            self.add_value_to_df(key=('acidity',''),
2821                                 index=index,
2822                                 new_value=acidity, 
2823                                 overwrite=overwrite) 
2824            if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')):
2825                self.add_value_to_df(key=('state_one','es_type'),
2826                                     index=index,
2827                                     new_value=self.propose_unused_type(), 
2828                                     overwrite=overwrite)  
2829            if pd.isna(self.df.loc [self.df['name']  == name].acidity.iloc[0]):
2830                self.add_value_to_df(key=('state_one','z'),
2831                                     index=index,
2832                                     new_value=default_charge_number, 
2833                                     overwrite=overwrite)
2834                self.add_value_to_df(key=('state_one','label'),
2835                                     index=index,
2836                                     new_value=name, 
2837                                    overwrite=overwrite)
2838            else:
2839                protonated_label = f'{name}H'
2840                self.add_value_to_df(key=('state_one','label'),
2841                                     index=index,
2842                                     new_value=protonated_label, 
2843                                    overwrite=overwrite)
2844                self.add_value_to_df(key=('state_two','label'),
2845                                     index=index,
2846                                     new_value=name, 
2847                                    overwrite=overwrite)
2848                if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')):
2849                    self.add_value_to_df(key=('state_two','es_type'),
2850                                         index=index,
2851                                         new_value=self.propose_unused_type(), 
2852                                         overwrite=overwrite)
2853                if self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'acidic':        
2854                    self.add_value_to_df(key=('state_one','z'),
2855                                         index=index,new_value=0, 
2856                                         overwrite=overwrite)
2857                    self.add_value_to_df(key=('state_two','z'),
2858                                         index=index,
2859                                         new_value=-1, 
2860                                         overwrite=overwrite)
2861                elif self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'basic':
2862                    self.add_value_to_df(key=('state_one','z'),
2863                                         index=index,new_value=+1, 
2864                                         overwrite=overwrite)
2865                    self.add_value_to_df(key=('state_two','z'),
2866                                         index=index,
2867                                         new_value=0, 
2868                                         overwrite=overwrite)   
2869        self.df.fillna(pd.NA, inplace=True)
2870        return

Sets the particle acidity including the charges in each of its possible states.

Arguments:
  • name(str): Unique label that identifies the particle.
  • acidity(str): Identifies whether the particle is acidic or basic, used to setup constant pH simulations. Defaults to None.
  • default_charge_number (int): Charge number of the particle. Defaults to 0.
  • pka(float, optional): If particle is an acid or a base, it defines its pka-value. Defaults to pandas.NA.
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
Note:
  • For particles with acidity = acidic or acidity = basic, state_one and state_two correspond to the protonated and

deprotonated states, respectively. - For particles without an acidity acidity = pandas.NA, only state_one is defined. - Each state has the following properties as sub-indexes: label,charge and es_type.

def set_reduced_units( self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True):
2872    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True):
2873        """
2874        Sets the set of reduced units used by pyMBE.units and it prints it.
2875
2876        Args:
2877            unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 
2878            unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
2879            temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 
2880            Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
2881            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
2882
2883        Note:
2884            - If no `temperature` is given, a value of 298.15 K is assumed by default.
2885            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
2886            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
2887            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2888        """
2889        if unit_length is None:
2890            unit_length= 0.355*self.units.nm
2891        if temperature is None:
2892            temperature = 298.15 * self.units.K
2893        if unit_charge is None:
2894            unit_charge = scipy.constants.e * self.units.C
2895        if Kw is None:
2896            Kw = 1e-14
2897        # Sanity check
2898        variables=[unit_length,temperature,unit_charge]
2899        dimensionalities=["[length]","[temperature]","[charge]"]
2900        for variable,dimensionality in zip(variables,dimensionalities):
2901            self.check_dimensionality(variable,dimensionality)
2902        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2903        self.kT=temperature*self.kB
2904        self.units._build_cache()
2905        self.units.define(f'reduced_energy = {self.kT} ')
2906        self.units.define(f'reduced_length = {unit_length}')
2907        self.units.define(f'reduced_charge = {unit_charge}')
2908        if verbose:        
2909            print(self.get_reduced_units())
2910        return

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.
  • verbose(bool, optional): Switch to activate/deactivate verbose. Defaults to True.
Note:
  • 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, pka_set=None, use_exclusion_radius_per_type=False):
2912    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
2913        """
2914        Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 
2915
2916        Args:
2917            counter_ion(`str`): `name` of the counter_ion `particle`.
2918            constant_pH(`float`): pH-value.
2919            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
2920            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2921            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
2922
2923        Returns:
2924            RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2925            sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
2926        """
2927        from espressomd import reaction_methods
2928        if exclusion_range is None:
2929            exclusion_range = max(self.get_radius_map().values())*2.0
2930        if pka_set is None:
2931            pka_set=self.get_pka_set()    
2932        self.check_pka_set(pka_set=pka_set)
2933        if use_exclusion_radius_per_type:
2934            exclusion_radius_per_type = self.get_radius_map()
2935        else:
2936            exclusion_radius_per_type = {}
2937        
2938        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2939                                                    exclusion_range=exclusion_range, 
2940                                                    seed=self.seed, 
2941                                                    constant_pH=constant_pH,
2942                                                    exclusion_radius_per_type = exclusion_radius_per_type
2943                                                    )
2944        sucessfull_reactions_labels=[]
2945        charge_number_map = self.get_charge_number_map()
2946        for name in pka_set.keys():
2947            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
2948                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
2949                continue
2950            gamma=10**-pka_set[name]['pka_value']
2951            state_one_type   = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
2952            state_two_type   = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
2953            counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0]
2954            RE.add_reaction(gamma=gamma,
2955                            reactant_types=[state_one_type],
2956                            product_types=[state_two_type, counter_ion_type],
2957                            default_charges={state_one_type: charge_number_map[state_one_type],
2958                            state_two_type: charge_number_map[state_two_type],
2959                            counter_ion_type: charge_number_map[counter_ion_type]})
2960            sucessfull_reactions_labels.append(name)
2961        return RE, sucessfull_reactions_labels

Sets up the Acid/Base reactions for acidic/basic particles defined in pmb.df 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.
  • pka_set(dict,optional): Desired pka_set, pka_set(dict): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
Returns:

RE(reaction_methods.ConstantpHEnsemble): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library. sucessfull_reactions_labels(lst): Labels of the reactions set up by pyMBE.

def setup_gcmc( self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type=False):
2963    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2964        """
2965        Sets up grand-canonical coupling to a reservoir of salt.
2966        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2967
2968        Args:
2969            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2970            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
2971            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
2972            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
2973            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
2974            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2975
2976        Returns:
2977            RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
2978        """
2979        from espressomd import reaction_methods
2980        if exclusion_range is None:
2981            exclusion_range = max(self.get_radius_map().values())*2.0
2982        if use_exclusion_radius_per_type:
2983            exclusion_radius_per_type = self.get_radius_map()
2984        else:
2985            exclusion_radius_per_type = {}
2986        
2987        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2988                                                    exclusion_range=exclusion_range, 
2989                                                    seed=self.seed, 
2990                                                    exclusion_radius_per_type = exclusion_radius_per_type
2991                                                    )
2992
2993        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
2994        determined_activity_coefficient = activity_coefficient(c_salt_res)
2995        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
2996
2997        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
2998        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
2999
3000        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
3001        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
3002
3003        if salt_cation_charge <= 0:
3004            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3005        if salt_anion_charge >= 0:
3006            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3007
3008        # Grand-canonical coupling to the reservoir
3009        RE.add_reaction(
3010            gamma = K_salt.magnitude,
3011            reactant_types = [],
3012            reactant_coefficients = [],
3013            product_types = [ salt_cation_es_type, salt_anion_es_type ],
3014            product_coefficients = [ 1, 1 ],
3015            default_charges = {
3016                salt_cation_es_type: salt_cation_charge, 
3017                salt_anion_es_type: salt_anion_charge, 
3018            }
3019        )
3020
3021        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:

RE (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, pka_set=None, use_exclusion_radius_per_type=False):
3023    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, pka_set=None, use_exclusion_radius_per_type = False):
3024        """
3025        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3026        reservoir of small ions. 
3027        This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
3028
3029        [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3030
3031        Args:
3032            pH_res ('float): pH-value in the reservoir.
3033            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3034            proton_name ('str'): Name of the proton (H+) particle.
3035            hydroxide_name ('str'): Name of the hydroxide (OH-) particle.
3036            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
3037            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
3038            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3039            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
3040            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3041            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
3042
3043        Returns:
3044            RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3045            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3046            ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3047        """
3048        from espressomd import reaction_methods
3049        if exclusion_range is None:
3050            exclusion_range = max(self.get_radius_map().values())*2.0
3051        if pka_set is None:
3052            pka_set=self.get_pka_set()    
3053        self.check_pka_set(pka_set=pka_set)
3054        if use_exclusion_radius_per_type:
3055            exclusion_radius_per_type = self.get_radius_map()
3056        else:
3057            exclusion_radius_per_type = {}
3058        
3059        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3060                                                    exclusion_range=exclusion_range, 
3061                                                    seed=self.seed, 
3062                                                    exclusion_radius_per_type = exclusion_radius_per_type
3063                                                    )
3064
3065        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3066        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3067        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3068        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3069        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3070        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3071        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
3072
3073        proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0]
3074        hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0]     
3075        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
3076        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
3077
3078        proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0]
3079        hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0]     
3080        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
3081        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
3082
3083        if proton_charge <= 0:
3084            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
3085        if salt_cation_charge <= 0:
3086            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
3087        if hydroxide_charge >= 0:
3088            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
3089        if salt_anion_charge >= 0:
3090            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
3091
3092        # Grand-canonical coupling to the reservoir
3093        # 0 = H+ + OH-
3094        RE.add_reaction(
3095            gamma = K_W.magnitude,
3096            reactant_types = [],
3097            reactant_coefficients = [],
3098            product_types = [ proton_es_type, hydroxide_es_type ],
3099            product_coefficients = [ 1, 1 ],
3100            default_charges = {
3101                proton_es_type: proton_charge, 
3102                hydroxide_es_type: hydroxide_charge, 
3103            }
3104        )
3105
3106        # 0 = Na+ + Cl-
3107        RE.add_reaction(
3108            gamma = K_NACL.magnitude,
3109            reactant_types = [],
3110            reactant_coefficients = [],
3111            product_types = [ salt_cation_es_type, salt_anion_es_type ],
3112            product_coefficients = [ 1, 1 ],
3113            default_charges = {
3114                salt_cation_es_type: salt_cation_charge, 
3115                salt_anion_es_type: salt_anion_charge, 
3116            }
3117        )
3118
3119        # 0 = Na+ + OH-
3120        RE.add_reaction(
3121            gamma = (K_NACL * K_W / K_HCL).magnitude,
3122            reactant_types = [],
3123            reactant_coefficients = [],
3124            product_types = [ salt_cation_es_type, hydroxide_es_type ],
3125            product_coefficients = [ 1, 1 ],
3126            default_charges = {
3127                salt_cation_es_type: salt_cation_charge, 
3128                hydroxide_es_type: hydroxide_charge, 
3129            }
3130        )
3131
3132        # 0 = H+ + Cl-
3133        RE.add_reaction(
3134            gamma = K_HCL.magnitude,
3135            reactant_types = [],
3136            reactant_coefficients = [],
3137            product_types = [ proton_es_type, salt_anion_es_type ],
3138            product_coefficients = [ 1, 1 ],
3139            default_charges = {
3140                proton_es_type: proton_charge, 
3141                salt_anion_es_type: salt_anion_charge, 
3142            }
3143        )
3144
3145        # Annealing moves to ensure sufficient sampling
3146        # Cation annealing H+ = Na+
3147        RE.add_reaction(
3148            gamma = (K_NACL / K_HCL).magnitude,
3149            reactant_types = [proton_es_type],
3150            reactant_coefficients = [ 1 ],
3151            product_types = [ salt_cation_es_type ],
3152            product_coefficients = [ 1 ],
3153            default_charges = {
3154                proton_es_type: proton_charge, 
3155                salt_cation_es_type: salt_cation_charge, 
3156            }
3157        )
3158
3159        # Anion annealing OH- = Cl- 
3160        RE.add_reaction(
3161            gamma = (K_HCL / K_W).magnitude,
3162            reactant_types = [hydroxide_es_type],
3163            reactant_coefficients = [ 1 ],
3164            product_types = [ salt_anion_es_type ],
3165            product_coefficients = [ 1 ],
3166            default_charges = {
3167                hydroxide_es_type: hydroxide_charge, 
3168                salt_anion_es_type: salt_anion_charge, 
3169            }
3170        )
3171
3172        sucessful_reactions_labels=[]
3173        charge_number_map = self.get_charge_number_map()
3174        for name in pka_set.keys():
3175            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3176                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3177                continue
3178
3179            Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3180
3181            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3182            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3183
3184            # Reaction in terms of proton: HA = A + H+
3185            RE.add_reaction(gamma=Ka.magnitude,
3186                            reactant_types=[state_one_type],
3187                            reactant_coefficients=[1],
3188                            product_types=[state_two_type, proton_es_type],
3189                            product_coefficients=[1, 1],
3190                            default_charges={state_one_type: charge_number_map[state_one_type],
3191                            state_two_type: charge_number_map[state_two_type],
3192                            proton_es_type: proton_charge})
3193
3194            # Reaction in terms of salt cation: HA = A + Na+
3195            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
3196                            reactant_types=[state_one_type],
3197                            reactant_coefficients=[1],
3198                            product_types=[state_two_type, salt_cation_es_type],
3199                            product_coefficients=[1, 1],
3200                            default_charges={state_one_type: charge_number_map[state_one_type],
3201                            state_two_type: charge_number_map[state_two_type],
3202                            salt_cation_es_type: salt_cation_charge})
3203
3204            # Reaction in terms of hydroxide: OH- + HA = A
3205            RE.add_reaction(gamma=(Ka / K_W).magnitude,
3206                            reactant_types=[state_one_type, hydroxide_es_type],
3207                            reactant_coefficients=[1, 1],
3208                            product_types=[state_two_type],
3209                            product_coefficients=[1],
3210                            default_charges={state_one_type: charge_number_map[state_one_type],
3211                            state_two_type: charge_number_map[state_two_type],
3212                            hydroxide_es_type: hydroxide_charge})
3213
3214            # Reaction in terms of salt anion: Cl- + HA = A
3215            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
3216                            reactant_types=[state_one_type, salt_anion_es_type],
3217                            reactant_coefficients=[1, 1],
3218                            product_types=[state_two_type],
3219                            product_coefficients=[1],
3220                            default_charges={state_one_type: charge_number_map[state_one_type],
3221                            state_two_type: charge_number_map[state_two_type],
3222                            salt_anion_es_type: salt_anion_charge})
3223
3224            sucessful_reactions_labels.append(name)
3225        return RE, sucessful_reactions_labels, ionic_strength_res

Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a reservoir of small ions. 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., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.

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.
  • pka_set(dict,optional): Desired pka_set, pka_set(dict): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
  • use_exclusion_radius_per_type(bool,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to False.
Returns:

RE (obj): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. sucessful_reactions_labels(lst): Labels of the reactions set up by pyMBE. ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients).

def setup_grxmc_unified( self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type=False):
3227    def setup_grxmc_unified(self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
3228        """
3229        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3230        reservoir of small ions. 
3231        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-}. 
3232        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'.
3233
3234        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3235        [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.
3236
3237        Args:
3238            pH_res ('float'): pH-value in the reservoir.
3239            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3240            cation_name ('str'): Name of the cationic particle.
3241            anion_name ('str'): Name of the anionic particle.
3242            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3243            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
3244            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3245            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`.
3246
3247        Returns:
3248            RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3249            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3250            ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3251        """
3252        from espressomd import reaction_methods
3253        if exclusion_range is None:
3254            exclusion_range = max(self.get_radius_map().values())*2.0
3255        if pka_set is None:
3256            pka_set=self.get_pka_set()    
3257        self.check_pka_set(pka_set=pka_set)
3258        if use_exclusion_radius_per_type:
3259            exclusion_radius_per_type = self.get_radius_map()
3260        else:
3261            exclusion_radius_per_type = {}
3262        
3263        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3264                                                    exclusion_range=exclusion_range, 
3265                                                    seed=self.seed, 
3266                                                    exclusion_radius_per_type = exclusion_radius_per_type
3267                                                    )
3268
3269        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3270        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3271        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3272        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3273        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3274        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3275        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3276        K_XX = a_cation * a_anion
3277
3278        cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0]
3279        anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0]     
3280        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
3281        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
3282        if cation_charge <= 0:
3283            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3284        if anion_charge >= 0:
3285            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3286
3287        # Coupling to the reservoir: 0 = X+ + X-
3288        RE.add_reaction(
3289            gamma = K_XX.magnitude,
3290            reactant_types = [],
3291            reactant_coefficients = [],
3292            product_types = [ cation_es_type, anion_es_type ],
3293            product_coefficients = [ 1, 1 ],
3294            default_charges = {
3295                cation_es_type: cation_charge, 
3296                anion_es_type: anion_charge, 
3297            }
3298        )
3299
3300        sucessful_reactions_labels=[]
3301        charge_number_map = self.get_charge_number_map()
3302        for name in pka_set.keys():
3303            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3304                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3305                continue
3306
3307            Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l
3308            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3309
3310            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3311            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3312
3313            # Reaction in terms of small cation: HA = A + X+
3314            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3315                            reactant_types=[state_one_type],
3316                            reactant_coefficients=[1],
3317                            product_types=[state_two_type, cation_es_type],
3318                            product_coefficients=[1, 1],
3319                            default_charges={state_one_type: charge_number_map[state_one_type],
3320                            state_two_type: charge_number_map[state_two_type],
3321                            cation_es_type: cation_charge})
3322
3323            # Reaction in terms of small anion: X- + HA = A
3324            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3325                            reactant_types=[state_one_type, anion_es_type],
3326                            reactant_coefficients=[1, 1],
3327                            product_types=[state_two_type],
3328                            product_coefficients=[1],
3329                            default_charges={state_one_type: charge_number_map[state_one_type],
3330                            state_two_type: charge_number_map[state_two_type],
3331                            anion_es_type: anion_charge})
3332
3333            sucessful_reactions_labels.append(name)
3334        return RE, sucessful_reactions_labels, ionic_strength_res

Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a reservoir of small ions. 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., KosÌŒovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020.

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.
  • pka_set(dict,optional): Desired pka_set, pka_set(dict): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
  • use_exclusion_radius_per_type(bool,optional): Controls if one exclusion_radius per each espresso_type. Defaults to False.
Returns:

RE (reaction_ensemble.ReactionEnsemble): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. sucessful_reactions_labels(lst): Labels of the reactions set up by pyMBE. ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).

def setup_df(self):
3336    def setup_df (self):
3337        """
3338        Sets up the pyMBE's dataframe `pymbe.df`.
3339
3340        Returns:
3341            columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe
3342        """
3343        
3344        columns_dtypes = {
3345            'name': {
3346                '': str},
3347            'pmb_type': {
3348                '': str},
3349            'particle_id': {
3350                '': pd.Int64Dtype()},
3351            'particle_id2':  {
3352                '': pd.Int64Dtype()},
3353            'residue_id':  {
3354                '': pd.Int64Dtype()},
3355            'molecule_id':  {
3356                '': pd.Int64Dtype()},
3357            'acidity':  {
3358                '': str},
3359            'pka':  {
3360                '': object},
3361            'central_bead':  {
3362                '': object},
3363            'side_chains': {
3364                '': object},
3365            'residue_list': {
3366                '': object},
3367            'model': {
3368                '': str},
3369            'sigma': {
3370                '': object},
3371            'cutoff': {
3372                '': object},
3373            'offset': {
3374                '': object},
3375            'epsilon': {
3376                '': object},
3377            'state_one': {
3378                'label': str,
3379                'es_type': pd.Int64Dtype(),
3380                'z': pd.Int64Dtype()},
3381            'state_two': {
3382                'label': str,
3383                'es_type': pd.Int64Dtype(),
3384                'z': pd.Int64Dtype()},
3385            'sequence': {
3386                '': object},
3387            'bond_object': {
3388                '': object},
3389            'parameters_of_the_potential':{
3390                '': object},
3391            'l0': {
3392                '': float},
3393            'node_map':{
3394                '':object},
3395            'chain_map':{
3396                '':object}}
3397        
3398        self.df = pd.DataFrame(columns=pd.MultiIndex.from_tuples([(col_main, col_sub) for col_main, sub_cols in columns_dtypes.items() for col_sub in sub_cols.keys()]))
3399        
3400        for level1, sub_dtypes in columns_dtypes.items():
3401            for level2, dtype in sub_dtypes.items():
3402                self.df[level1, level2] = self.df[level1, level2].astype(dtype)
3403                
3404        columns_names = pd.MultiIndex.from_frame(self.df)
3405        columns_names = columns_names.names
3406                
3407        return columns_names

Sets up the pyMBE's dataframe pymbe.df.

Returns:

columns_names(obj): pandas multiindex object with the column names of the pyMBE's dataframe

def setup_lj_interactions( self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True):
3409    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True):
3410        """
3411        Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`.
3412
3413        Args:
3414            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
3415            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.
3416            combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3417            warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True.
3418
3419        Note:
3420            - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 
3421            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
3422            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3423
3424        """
3425        from itertools import combinations_with_replacement
3426        import warnings
3427        implemented_combining_rules = ['Lorentz-Berthelot']
3428        compulsory_parameters_in_df = ['sigma','epsilon']
3429        # Sanity check
3430        if combining_rule not in implemented_combining_rules:
3431            raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules)
3432        shift=0
3433        if shift_potential:
3434            shift="auto"
3435        # List which particles have sigma and epsilon values defined in pmb.df and which ones don't
3436        particles_types_with_LJ_parameters = []
3437        non_parametrized_labels= []
3438        for particle_type in self.get_type_map().values():
3439            check_list=[]
3440            for key in compulsory_parameters_in_df:
3441                value_in_df=self.find_value_from_es_type(es_type=particle_type,
3442                                                        column_name=key)
3443                check_list.append(pd.isna(value_in_df))
3444            if any(check_list):
3445                non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 
3446                                                                            column_name='label'))
3447            else:
3448                particles_types_with_LJ_parameters.append(particle_type)
3449        # Set up LJ interactions between all particle types
3450        for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 
3451            particle_name1 = self.find_value_from_es_type(es_type=type_pair[0],
3452                                                        column_name="name")
3453            particle_name2 = self.find_value_from_es_type(es_type=type_pair[1],
3454                                                        column_name="name")
3455            lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1,
3456                                                 particle_name2 = particle_name2,
3457                                                 combining_rule = combining_rule)
3458            
3459            # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
3460            if not lj_parameters:
3461                continue
3462            espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 
3463                                                                                    sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 
3464                                                                                    cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude,
3465                                                                                    offset = lj_parameters["offset"].to("reduced_length").magnitude, 
3466                                                                                    shift = shift)                                                                                          
3467            index = len(self.df)
3468            label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label")
3469            label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label")
3470            self.df.at [index, 'name'] = f'LJ: {label1}-{label2}'
3471            lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params()
3472
3473            self.add_value_to_df(index=index,
3474                                key=('pmb_type',''),
3475                                new_value='LennardJones')
3476            
3477            self.add_value_to_df(index=index,
3478                                key=('parameters_of_the_potential',''),
3479                                new_value=lj_params,
3480                                non_standard_value=True)
3481        if non_parametrized_labels and warnings:
3482            warnings.warn(f'The following particles do not have a defined value of sigma or epsilon in pmb.df: {non_parametrized_labels}. No LJ interaction has been added in ESPResSo for those particles.',UserWarning) 
3483        return

Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for sigma, offset, and epsilon stored in pymbe.df.

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.
Note:
  • LJ interactions will only be set up between particles with defined values of sigma and epsilon in the pmb.df.
  • Currently, the only combining_rule supported is Lorentz-Berthelot.
  • Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
def write_pmb_df(self, filename):
3485    def write_pmb_df (self, filename):
3486        '''
3487        Writes the pyMBE dataframe into a csv file
3488        
3489        Args:
3490            filename(`str`): Path to the csv file 
3491        '''
3492
3493        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence']
3494        df = self.df.copy(deep=True)
3495        for column_name in columns_with_list_or_dict:
3496            df[column_name] = df[column_name].apply(lambda x: json.dumps(x) if isinstance(x, (np.ndarray, tuple, list, dict)) or pd.notnull(x) else x)
3497        df['bond_object'] = df['bond_object'].apply(lambda x: f'{x.__class__.__name__}({json.dumps({**x.get_params(), "bond_id": x._bond_id})})' if pd.notnull(x) else x)
3498        df.fillna(pd.NA, inplace=True)
3499        df.to_csv(filename)
3500        return

Writes the pyMBE dataframe into a csv file

Arguments:
  • filename(str): Path to the csv file
class pymbe_library.NumpyEncoder(json.encoder.JSONEncoder):
44    class NumpyEncoder(json.JSONEncoder):
45        """
46        Custom JSON encoder that converts NumPy arrays to Python lists
47        and NumPy scalars to Python scalars.
48        """
49        def default(self, obj):
50            if isinstance(obj, np.ndarray):
51                return obj.tolist()
52            if isinstance(obj, np.generic):
53                return obj.item()
54            return super().default(obj)

Custom JSON encoder that converts NumPy arrays to Python lists and NumPy scalars to Python scalars.

def default(self, obj):
49        def default(self, obj):
50            if isinstance(obj, np.ndarray):
51                return obj.tolist()
52            if isinstance(obj, np.generic):
53                return obj.item()
54            return super().default(obj)

Implement this method in a subclass such that it returns a serializable object for o, or calls the base implementation (to raise a TypeError).

For example, to support arbitrary iterators, you could implement default like this::

def default(self, o):
    try:
        iterable = iter(o)
    except TypeError:
        pass
    else:
        return list(iterable)
    # Let the base class default method raise the TypeError
    return JSONEncoder.default(self, o)
Inherited Members
json.encoder.JSONEncoder
JSONEncoder
item_separator
key_separator
skipkeys
ensure_ascii
check_circular
allow_nan
sort_keys
indent
encode
iterencode