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

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

Returns:

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

def load_interaction_parameters(self, filename, overwrite=False):
2195    def load_interaction_parameters(self, filename, overwrite=False):
2196        """
2197        Loads the interaction parameters stored in `filename` into `pmb.df`
2198        
2199        Args:
2200            filename(`str`): name of the file to be read
2201            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2202        """
2203        without_units = ['z','es_type']
2204        with_units = ['sigma','epsilon','offset','cutoff']
2205
2206        with open(filename, 'r') as f:
2207            interaction_data = json.load(f)
2208            interaction_parameter_set = interaction_data["data"]
2209
2210        for key in interaction_parameter_set:
2211            param_dict=interaction_parameter_set[key]
2212            object_type=param_dict.pop('object_type')
2213            if object_type == 'particle':
2214                not_required_attributes={}    
2215                for not_required_key in without_units+with_units:
2216                    if not_required_key in param_dict.keys():
2217                        if not_required_key in with_units:
2218                            not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key))
2219                        elif not_required_key in without_units:
2220                            not_required_attributes[not_required_key]=param_dict.pop(not_required_key)
2221                    else:
2222                        not_required_attributes[not_required_key]=pd.NA
2223                self.define_particle(name=param_dict.pop('name'),
2224                                z=not_required_attributes.pop('z'),
2225                                sigma=not_required_attributes.pop('sigma'),
2226                                offset=not_required_attributes.pop('offset'),
2227                                cutoff=not_required_attributes.pop('cutoff'),
2228                                epsilon=not_required_attributes.pop('epsilon'),
2229                                overwrite=overwrite)
2230            elif object_type == 'residue':
2231                self.define_residue(**param_dict)
2232            elif object_type == 'molecule':
2233                self.define_molecule(**param_dict)
2234            elif object_type == 'peptide':
2235                self.define_peptide(**param_dict)
2236            elif object_type == 'bond':
2237                particle_pairs = param_dict.pop('particle_pairs')
2238                bond_parameters = param_dict.pop('bond_parameters')
2239                bond_type = param_dict.pop('bond_type')
2240                if bond_type == 'harmonic':
2241                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2242                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2243                    bond = {'r_0'    : r_0,
2244                            'k'      : k,
2245                            }
2246
2247                elif bond_type == 'FENE':
2248                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2249                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2250                    d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max'))
2251                    bond =  {'r_0'    : r_0,
2252                             'k'      : k,
2253                             'd_r_max': d_r_max,
2254                             }
2255                else:
2256                    raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds")
2257                self.define_bond(bond_type=bond_type, 
2258                                bond_parameters=bond, 
2259                                particle_pairs=particle_pairs)
2260            else:
2261                raise ValueError(object_type+' is not a known pmb object type')
2262            
2263        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):
2265    def load_pka_set(self, filename, overwrite=True):
2266        """
2267        Loads the pka_set stored in `filename` into `pmb.df`.
2268        
2269        Args:
2270            filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}.
2271            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 
2272        """
2273        with open(filename, 'r') as f:
2274            pka_data = json.load(f)
2275            pka_set = pka_data["data"]
2276
2277        self.check_pka_set(pka_set=pka_set)
2278
2279        for key in pka_set:
2280            acidity = pka_set[key]['acidity']
2281            pka_value = pka_set[key]['pka_value']
2282            self.set_particle_acidity(name=key, 
2283                                      acidity=acidity, 
2284                                      pka=pka_value, 
2285                                      overwrite=overwrite)
2286        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):
2289    def propose_unused_type(self):
2290        """
2291        Searches in `pmb.df` all the different particle types defined and returns a new one.
2292
2293        Returns:
2294            unused_type(`int`): unused particle type
2295        """
2296        type_map = self.get_type_map()
2297        if not type_map:  
2298            unused_type = 0
2299        else:
2300            valid_values = [v for v in type_map.values() if pd.notna(v)]  # Filter out pd.NA values
2301            unused_type = max(valid_values) + 1 if valid_values else 0  # Ensure max() doesn't fail if all values are NA
2302        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):
2304    def protein_sequence_parser(self, sequence):
2305        '''
2306        Parses `sequence` to the one letter code for amino acids.
2307        
2308        Args:
2309            sequence(`str` or `lst`): Sequence of the amino acid. 
2310
2311        Returns:
2312            clean_sequence(`lst`): `sequence` using the one letter code.
2313        
2314        Note:
2315            - Accepted formats for `sequence` are:
2316                - `lst` with one letter or three letter code of each aminoacid in each element
2317                - `str` with the sequence using the one letter code
2318                - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-"
2319        
2320        '''
2321        # Aminoacid key
2322        keys={"ALA": "A",
2323                "ARG": "R",
2324                "ASN": "N",
2325                "ASP": "D",
2326                "CYS": "C",
2327                "GLU": "E",
2328                "GLN": "Q",
2329                "GLY": "G",
2330                "HIS": "H",
2331                "ILE": "I",
2332                "LEU": "L",
2333                "LYS": "K",
2334                "MET": "M",
2335                "PHE": "F",
2336                "PRO": "P",
2337                "SER": "S",
2338                "THR": "T",
2339                "TRP": "W",
2340                "TYR": "Y",
2341                "VAL": "V",
2342                "PSER": "J",
2343                "PTHR": "U",
2344                "PTyr": "Z",
2345                "NH2": "n",
2346                "COOH": "c"}
2347        clean_sequence=[]
2348        if isinstance(sequence, str):
2349            if sequence.find("-") != -1:
2350                splited_sequence=sequence.split("-")
2351                for residue in splited_sequence:
2352                    if len(residue) == 1:
2353                        if residue in keys.values():
2354                            residue_ok=residue
2355                        else:
2356                            if residue.upper() in keys.values():
2357                                residue_ok=residue.upper()
2358                            else:
2359                                raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence")
2360                        clean_sequence.append(residue_ok)
2361                    else:
2362                        if residue in keys.keys():
2363                            clean_sequence.append(keys[residue])
2364                        else:
2365                            if residue.upper() in keys.keys():
2366                                clean_sequence.append(keys[residue.upper()])
2367                            else:
2368                                raise ValueError("Unknown  code for a residue: ", residue, " please review the input sequence")
2369            else:
2370                for residue in sequence:
2371                    if residue in keys.values():
2372                        residue_ok=residue
2373                    else:
2374                        if residue.upper() in keys.values():
2375                            residue_ok=residue.upper()
2376                        else:
2377                            raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence")
2378                    clean_sequence.append(residue_ok)
2379        if isinstance(sequence, list):
2380            for residue in sequence:
2381                if residue in keys.values():
2382                    residue_ok=residue
2383                else:
2384                    if residue.upper() in keys.values():
2385                        residue_ok=residue.upper()
2386                    elif (residue.upper() in keys.keys()):
2387                        residue_ok= keys[residue.upper()]
2388                    else:
2389                        raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence")
2390                clean_sequence.append(residue_ok)
2391        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):
2393    def read_pmb_df (self,filename):
2394        """
2395        Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE.
2396
2397        Args:
2398            filename(`str`): path to file.
2399
2400        Note:
2401            This function only accepts files with CSV format. 
2402        """
2403        
2404        if filename.rsplit(".", 1)[1] != "csv":
2405            raise ValueError("Only files with CSV format are supported")
2406        df = pd.read_csv (filename,header=[0, 1], index_col=0)
2407        columns_names = self.setup_df()
2408        
2409        multi_index = pd.MultiIndex.from_tuples(columns_names)
2410        df.columns = multi_index
2411        
2412        self.df = self.convert_columns_to_original_format(df)
2413        self.df.fillna(pd.NA, inplace=True)
2414        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):
2416    def read_protein_vtf_in_df (self,filename,unit_length=None):
2417        """
2418        Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`.
2419
2420        Args:
2421            filename(`str`): Path to the vtf file with the coarse-grained model.
2422            unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None.
2423
2424        Returns:
2425            topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}
2426
2427        Note:
2428            - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom.
2429        """
2430
2431        print (f'Loading protein coarse grain model file: {filename}')
2432
2433        coord_list = []
2434        particles_dict = {}
2435
2436        if unit_length is None:
2437            unit_length = 1 * self.units.angstrom 
2438
2439        with open (filename,'r') as protein_model:
2440            for line in protein_model :
2441                line_split = line.split()
2442                if line_split:
2443                    line_header = line_split[0]
2444                    if line_header == 'atom':
2445                        atom_id  = line_split[1]
2446                        atom_name = line_split[3]
2447                        atom_resname = line_split[5]
2448                        chain_id = line_split[9]
2449                        radius = float(line_split [11])*unit_length 
2450                        particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius]
2451                    elif line_header.isnumeric(): 
2452                        atom_coord = line_split[1:] 
2453                        atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord]
2454                        coord_list.append (atom_coord)
2455
2456        numbered_label = []
2457        i = 0   
2458        
2459        for atom_id in particles_dict.keys():
2460    
2461            if atom_id == 1:
2462                atom_name = particles_dict[atom_id][0]
2463                numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2464                numbered_label.append(numbered_name)
2465
2466            elif atom_id != 1: 
2467            
2468                if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]:
2469                    i += 1                    
2470                    count = 1
2471                    atom_name = particles_dict[atom_id][0]
2472                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2473                    numbered_label.append(numbered_name)
2474                    
2475                elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]:
2476                    if count == 2 or particles_dict[atom_id][1] == 'GLY':
2477                        i +=1  
2478                        count = 0
2479                    atom_name = particles_dict[atom_id][0]
2480                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2481                    numbered_label.append(numbered_name)
2482                    count +=1
2483
2484        topology_dict = {}
2485
2486        for i in range (0, len(numbered_label)):   
2487            topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] ,
2488                                                    'chain_id':numbered_label[i][1],
2489                                                    'radius':numbered_label[i][2] }
2490
2491        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):
2493    def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2494        """
2495        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2496        If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead.
2497        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.
2498
2499        Args:
2500            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2501            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2502            hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2503            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. 
2504
2505        Returns:
2506            bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library.
2507        
2508        Note:
2509            - 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`.
2510            - If `hard_check`=`True` stops the code when no bond is found.
2511        """
2512        
2513        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2514                                    particle_name2=particle_name2, 
2515                                    use_default_bond=use_default_bond)
2516        if use_default_bond:
2517            if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'):
2518                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")
2519        if bond_key:
2520            return self.df[self.df['name']==bond_key].bond_object.values[0]
2521        else:
2522            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2523            if hard_check:
2524                sys.exit(1)
2525            else:
2526                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):
2528    def search_particles_in_residue(self, residue_name):
2529        '''
2530        Searches for all particles in a given residue of name `residue_name`.
2531
2532        Args:
2533            residue_name (`str`): name of the residue to be searched
2534
2535        Returns:
2536            list_of_particles_in_residue (`lst`): list of the names of all particles in the residue
2537
2538        Note:
2539            - 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.
2540 
2541        '''
2542        index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 
2543        central_bead = self.df.at [index_residue, ('central_bead', '')]
2544        list_of_side_chains = self.df.at [index_residue, ('side_chains', '')]
2545
2546        list_of_particles_in_residue = []
2547        list_of_particles_in_residue.append(central_bead)
2548
2549        for side_chain in list_of_side_chains: 
2550            object_type = self.df[self.df['name']==side_chain].pmb_type.values[0]
2551
2552            if object_type == "residue":
2553                list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain)
2554                list_of_particles_in_residue += list_of_particles_in_side_chain_residue
2555            elif object_type == "particle":
2556                list_of_particles_in_residue.append(side_chain)
2557
2558        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):
2560    def set_particle_acidity(self, name, acidity=pd.NA, default_charge_number=0, pka=pd.NA, overwrite=True):
2561        """
2562        Sets the particle acidity including the charges in each of its possible states. 
2563
2564        Args:
2565            name(`str`): Unique label that identifies the `particle`. 
2566            acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None.
2567            default_charge_number (`int`): Charge number of the particle. Defaults to 0.
2568            pka(`float`, optional):  If `particle` is an acid or a base, it defines its pka-value. Defaults to pandas.NA.
2569            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2570     
2571        Note:
2572            - For particles with  `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 
2573        deprotonated states, respectively. 
2574            - For particles without an acidity `acidity = pandas.NA`, only `state_one` is defined.
2575            - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 
2576        """
2577        acidity_valid_keys = ['inert','acidic', 'basic']
2578        if not pd.isna(acidity):
2579            if acidity not in acidity_valid_keys:
2580                raise ValueError(f"Acidity {acidity} provided for particle name  {name} is not supproted. Valid keys are: {acidity_valid_keys}")
2581            if acidity in ['acidic', 'basic'] and pd.isna(pka):
2582                raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.")   
2583            if acidity == "inert":
2584                acidity = pd.NA
2585                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")
2586
2587        self.define_particle_entry_in_df(name=name)
2588        
2589        for index in self.df[self.df['name']==name].index:       
2590            if pka is not pd.NA:
2591                self.add_value_to_df(key=('pka',''),
2592                                    index=index,
2593                                    new_value=pka, 
2594                                    overwrite=overwrite)
2595            
2596            self.add_value_to_df(key=('acidity',''),
2597                                 index=index,
2598                                 new_value=acidity, 
2599                                 overwrite=overwrite) 
2600            if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')):
2601                self.add_value_to_df(key=('state_one','es_type'),
2602                                     index=index,
2603                                     new_value=self.propose_unused_type(), 
2604                                     overwrite=overwrite)  
2605            if pd.isna(self.df.loc [self.df['name']  == name].acidity.iloc[0]):
2606                self.add_value_to_df(key=('state_one','z'),
2607                                     index=index,
2608                                     new_value=default_charge_number, 
2609                                     overwrite=overwrite)
2610                self.add_value_to_df(key=('state_one','label'),
2611                                     index=index,
2612                                     new_value=name, 
2613                                    overwrite=overwrite)
2614            else:
2615                protonated_label = f'{name}H'
2616                self.add_value_to_df(key=('state_one','label'),
2617                                     index=index,
2618                                     new_value=protonated_label, 
2619                                    overwrite=overwrite)
2620                self.add_value_to_df(key=('state_two','label'),
2621                                     index=index,
2622                                     new_value=name, 
2623                                    overwrite=overwrite)
2624                if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')):
2625                    self.add_value_to_df(key=('state_two','es_type'),
2626                                         index=index,
2627                                         new_value=self.propose_unused_type(), 
2628                                         overwrite=overwrite)
2629                if self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'acidic':        
2630                    self.add_value_to_df(key=('state_one','z'),
2631                                         index=index,new_value=0, 
2632                                         overwrite=overwrite)
2633                    self.add_value_to_df(key=('state_two','z'),
2634                                         index=index,
2635                                         new_value=-1, 
2636                                         overwrite=overwrite)
2637                elif self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'basic':
2638                    self.add_value_to_df(key=('state_one','z'),
2639                                         index=index,new_value=+1, 
2640                                         overwrite=overwrite)
2641                    self.add_value_to_df(key=('state_two','z'),
2642                                         index=index,
2643                                         new_value=0, 
2644                                         overwrite=overwrite)   
2645        self.df.fillna(pd.NA, inplace=True)
2646        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):
2648    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True):
2649        """
2650        Sets the set of reduced units used by pyMBE.units and it prints it.
2651
2652        Args:
2653            unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 
2654            unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
2655            temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 
2656            Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
2657            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
2658
2659        Note:
2660            - If no `temperature` is given, a value of 298.15 K is assumed by default.
2661            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
2662            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
2663            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2664        """
2665        if unit_length is None:
2666            unit_length= 0.355*self.units.nm
2667        if temperature is None:
2668            temperature = 298.15 * self.units.K
2669        if unit_charge is None:
2670            unit_charge = scipy.constants.e * self.units.C
2671        if Kw is None:
2672            Kw = 1e-14
2673        # Sanity check
2674        variables=[unit_length,temperature,unit_charge]
2675        dimensionalities=["[length]","[temperature]","[charge]"]
2676        for variable,dimensionality in zip(variables,dimensionalities):
2677            self.check_dimensionality(variable,dimensionality)
2678        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2679        self.kT=temperature*self.kB
2680        self.units._build_cache()
2681        self.units.define(f'reduced_energy = {self.kT} ')
2682        self.units.define(f'reduced_length = {unit_length}')
2683        self.units.define(f'reduced_charge = {unit_charge}')
2684        if verbose:        
2685            print(self.get_reduced_units())
2686        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):
2688    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
2689        """
2690        Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 
2691
2692        Args:
2693            counter_ion(`str`): `name` of the counter_ion `particle`.
2694            constant_pH(`float`): pH-value.
2695            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
2696            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2697            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
2698
2699        Returns:
2700            RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2701            sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
2702        """
2703        from espressomd import reaction_methods
2704        if exclusion_range is None:
2705            exclusion_range = max(self.get_radius_map().values())*2.0
2706        if pka_set is None:
2707            pka_set=self.get_pka_set()    
2708        self.check_pka_set(pka_set=pka_set)
2709        if use_exclusion_radius_per_type:
2710            exclusion_radius_per_type = self.get_radius_map()
2711        else:
2712            exclusion_radius_per_type = {}
2713        
2714        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2715                                                    exclusion_range=exclusion_range, 
2716                                                    seed=self.seed, 
2717                                                    constant_pH=constant_pH,
2718                                                    exclusion_radius_per_type = exclusion_radius_per_type
2719                                                    )
2720        sucessfull_reactions_labels=[]
2721        charge_number_map = self.get_charge_number_map()
2722        for name in pka_set.keys():
2723            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
2724                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
2725                continue
2726            gamma=10**-pka_set[name]['pka_value']
2727            state_one_type   = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
2728            state_two_type   = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
2729            counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0]
2730            RE.add_reaction(gamma=gamma,
2731                            reactant_types=[state_one_type],
2732                            product_types=[state_two_type, counter_ion_type],
2733                            default_charges={state_one_type: charge_number_map[state_one_type],
2734                            state_two_type: charge_number_map[state_two_type],
2735                            counter_ion_type: charge_number_map[counter_ion_type]})
2736            sucessfull_reactions_labels.append(name)
2737        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):
2739    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2740        """
2741        Sets up grand-canonical coupling to a reservoir of salt.
2742        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2743
2744        Args:
2745            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2746            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
2747            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
2748            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
2749            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
2750            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2751
2752        Returns:
2753            RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
2754        """
2755        from espressomd import reaction_methods
2756        if exclusion_range is None:
2757            exclusion_range = max(self.get_radius_map().values())*2.0
2758        if use_exclusion_radius_per_type:
2759            exclusion_radius_per_type = self.get_radius_map()
2760        else:
2761            exclusion_radius_per_type = {}
2762        
2763        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2764                                                    exclusion_range=exclusion_range, 
2765                                                    seed=self.seed, 
2766                                                    exclusion_radius_per_type = exclusion_radius_per_type
2767                                                    )
2768
2769        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
2770        determined_activity_coefficient = activity_coefficient(c_salt_res)
2771        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
2772
2773        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
2774        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
2775
2776        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
2777        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
2778
2779        if salt_cation_charge <= 0:
2780            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
2781        if salt_anion_charge >= 0:
2782            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
2783
2784        # Grand-canonical coupling to the reservoir
2785        RE.add_reaction(
2786            gamma = K_salt.magnitude,
2787            reactant_types = [],
2788            reactant_coefficients = [],
2789            product_types = [ salt_cation_es_type, salt_anion_es_type ],
2790            product_coefficients = [ 1, 1 ],
2791            default_charges = {
2792                salt_cation_es_type: salt_cation_charge, 
2793                salt_anion_es_type: salt_anion_charge, 
2794            }
2795        )
2796
2797        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):
2799    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):
2800        """
2801        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
2802        reservoir of small ions. 
2803        This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
2804
2805        [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.
2806
2807        Args:
2808            pH_res ('float): pH-value in the reservoir.
2809            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2810            proton_name ('str'): Name of the proton (H+) particle.
2811            hydroxide_name ('str'): Name of the hydroxide (OH-) particle.
2812            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
2813            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
2814            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
2815            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
2816            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
2817            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2818
2819        Returns:
2820            RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
2821            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
2822            ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients).
2823        """
2824        from espressomd import reaction_methods
2825        if exclusion_range is None:
2826            exclusion_range = max(self.get_radius_map().values())*2.0
2827        if pka_set is None:
2828            pka_set=self.get_pka_set()    
2829        self.check_pka_set(pka_set=pka_set)
2830        if use_exclusion_radius_per_type:
2831            exclusion_radius_per_type = self.get_radius_map()
2832        else:
2833            exclusion_radius_per_type = {}
2834        
2835        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2836                                                    exclusion_range=exclusion_range, 
2837                                                    seed=self.seed, 
2838                                                    exclusion_radius_per_type = exclusion_radius_per_type
2839                                                    )
2840
2841        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
2842        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
2843        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
2844        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
2845        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
2846        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
2847        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
2848
2849        proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0]
2850        hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0]     
2851        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
2852        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
2853
2854        proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0]
2855        hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0]     
2856        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
2857        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
2858
2859        if proton_charge <= 0:
2860            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
2861        if salt_cation_charge <= 0:
2862            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
2863        if hydroxide_charge >= 0:
2864            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
2865        if salt_anion_charge >= 0:
2866            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
2867
2868        # Grand-canonical coupling to the reservoir
2869        # 0 = H+ + OH-
2870        RE.add_reaction(
2871            gamma = K_W.magnitude,
2872            reactant_types = [],
2873            reactant_coefficients = [],
2874            product_types = [ proton_es_type, hydroxide_es_type ],
2875            product_coefficients = [ 1, 1 ],
2876            default_charges = {
2877                proton_es_type: proton_charge, 
2878                hydroxide_es_type: hydroxide_charge, 
2879            }
2880        )
2881
2882        # 0 = Na+ + Cl-
2883        RE.add_reaction(
2884            gamma = K_NACL.magnitude,
2885            reactant_types = [],
2886            reactant_coefficients = [],
2887            product_types = [ salt_cation_es_type, salt_anion_es_type ],
2888            product_coefficients = [ 1, 1 ],
2889            default_charges = {
2890                salt_cation_es_type: salt_cation_charge, 
2891                salt_anion_es_type: salt_anion_charge, 
2892            }
2893        )
2894
2895        # 0 = Na+ + OH-
2896        RE.add_reaction(
2897            gamma = (K_NACL * K_W / K_HCL).magnitude,
2898            reactant_types = [],
2899            reactant_coefficients = [],
2900            product_types = [ salt_cation_es_type, hydroxide_es_type ],
2901            product_coefficients = [ 1, 1 ],
2902            default_charges = {
2903                salt_cation_es_type: salt_cation_charge, 
2904                hydroxide_es_type: hydroxide_charge, 
2905            }
2906        )
2907
2908        # 0 = H+ + Cl-
2909        RE.add_reaction(
2910            gamma = K_HCL.magnitude,
2911            reactant_types = [],
2912            reactant_coefficients = [],
2913            product_types = [ proton_es_type, salt_anion_es_type ],
2914            product_coefficients = [ 1, 1 ],
2915            default_charges = {
2916                proton_es_type: proton_charge, 
2917                salt_anion_es_type: salt_anion_charge, 
2918            }
2919        )
2920
2921        # Annealing moves to ensure sufficient sampling
2922        # Cation annealing H+ = Na+
2923        RE.add_reaction(
2924            gamma = (K_NACL / K_HCL).magnitude,
2925            reactant_types = [proton_es_type],
2926            reactant_coefficients = [ 1 ],
2927            product_types = [ salt_cation_es_type ],
2928            product_coefficients = [ 1 ],
2929            default_charges = {
2930                proton_es_type: proton_charge, 
2931                salt_cation_es_type: salt_cation_charge, 
2932            }
2933        )
2934
2935        # Anion annealing OH- = Cl- 
2936        RE.add_reaction(
2937            gamma = (K_HCL / K_W).magnitude,
2938            reactant_types = [hydroxide_es_type],
2939            reactant_coefficients = [ 1 ],
2940            product_types = [ salt_anion_es_type ],
2941            product_coefficients = [ 1 ],
2942            default_charges = {
2943                hydroxide_es_type: hydroxide_charge, 
2944                salt_anion_es_type: salt_anion_charge, 
2945            }
2946        )
2947
2948        sucessful_reactions_labels=[]
2949        charge_number_map = self.get_charge_number_map()
2950        for name in pka_set.keys():
2951            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
2952                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
2953                continue
2954
2955            Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
2956
2957            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
2958            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
2959
2960            # Reaction in terms of proton: HA = A + H+
2961            RE.add_reaction(gamma=Ka.magnitude,
2962                            reactant_types=[state_one_type],
2963                            reactant_coefficients=[1],
2964                            product_types=[state_two_type, proton_es_type],
2965                            product_coefficients=[1, 1],
2966                            default_charges={state_one_type: charge_number_map[state_one_type],
2967                            state_two_type: charge_number_map[state_two_type],
2968                            proton_es_type: proton_charge})
2969
2970            # Reaction in terms of salt cation: HA = A + Na+
2971            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
2972                            reactant_types=[state_one_type],
2973                            reactant_coefficients=[1],
2974                            product_types=[state_two_type, salt_cation_es_type],
2975                            product_coefficients=[1, 1],
2976                            default_charges={state_one_type: charge_number_map[state_one_type],
2977                            state_two_type: charge_number_map[state_two_type],
2978                            salt_cation_es_type: salt_cation_charge})
2979
2980            # Reaction in terms of hydroxide: OH- + HA = A
2981            RE.add_reaction(gamma=(Ka / K_W).magnitude,
2982                            reactant_types=[state_one_type, hydroxide_es_type],
2983                            reactant_coefficients=[1, 1],
2984                            product_types=[state_two_type],
2985                            product_coefficients=[1],
2986                            default_charges={state_one_type: charge_number_map[state_one_type],
2987                            state_two_type: charge_number_map[state_two_type],
2988                            hydroxide_es_type: hydroxide_charge})
2989
2990            # Reaction in terms of salt anion: Cl- + HA = A
2991            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
2992                            reactant_types=[state_one_type, salt_anion_es_type],
2993                            reactant_coefficients=[1, 1],
2994                            product_types=[state_two_type],
2995                            product_coefficients=[1],
2996                            default_charges={state_one_type: charge_number_map[state_one_type],
2997                            state_two_type: charge_number_map[state_two_type],
2998                            salt_anion_es_type: salt_anion_charge})
2999
3000            sucessful_reactions_labels.append(name)
3001        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):
3003    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):
3004        """
3005        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3006        reservoir of small ions. 
3007        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-}. 
3008        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'.
3009
3010        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3011        [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.
3012
3013        Args:
3014            pH_res ('float'): pH-value in the reservoir.
3015            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3016            cation_name ('str'): Name of the cationic particle.
3017            anion_name ('str'): Name of the anionic particle.
3018            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3019            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
3020            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3021            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`.
3022
3023        Returns:
3024            RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3025            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3026            ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3027        """
3028        from espressomd import reaction_methods
3029        if exclusion_range is None:
3030            exclusion_range = max(self.get_radius_map().values())*2.0
3031        if pka_set is None:
3032            pka_set=self.get_pka_set()    
3033        self.check_pka_set(pka_set=pka_set)
3034        if use_exclusion_radius_per_type:
3035            exclusion_radius_per_type = self.get_radius_map()
3036        else:
3037            exclusion_radius_per_type = {}
3038        
3039        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3040                                                    exclusion_range=exclusion_range, 
3041                                                    seed=self.seed, 
3042                                                    exclusion_radius_per_type = exclusion_radius_per_type
3043                                                    )
3044
3045        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3046        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3047        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3048        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3049        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3050        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3051        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3052        K_XX = a_cation * a_anion
3053
3054        cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0]
3055        anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0]     
3056        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
3057        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
3058        if cation_charge <= 0:
3059            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3060        if anion_charge >= 0:
3061            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3062
3063        # Coupling to the reservoir: 0 = X+ + X-
3064        RE.add_reaction(
3065            gamma = K_XX.magnitude,
3066            reactant_types = [],
3067            reactant_coefficients = [],
3068            product_types = [ cation_es_type, anion_es_type ],
3069            product_coefficients = [ 1, 1 ],
3070            default_charges = {
3071                cation_es_type: cation_charge, 
3072                anion_es_type: anion_charge, 
3073            }
3074        )
3075
3076        sucessful_reactions_labels=[]
3077        charge_number_map = self.get_charge_number_map()
3078        for name in pka_set.keys():
3079            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3080                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3081                continue
3082
3083            Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l
3084            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3085
3086            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3087            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3088
3089            # Reaction in terms of small cation: HA = A + X+
3090            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3091                            reactant_types=[state_one_type],
3092                            reactant_coefficients=[1],
3093                            product_types=[state_two_type, cation_es_type],
3094                            product_coefficients=[1, 1],
3095                            default_charges={state_one_type: charge_number_map[state_one_type],
3096                            state_two_type: charge_number_map[state_two_type],
3097                            cation_es_type: cation_charge})
3098
3099            # Reaction in terms of small anion: X- + HA = A
3100            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3101                            reactant_types=[state_one_type, anion_es_type],
3102                            reactant_coefficients=[1, 1],
3103                            product_types=[state_two_type],
3104                            product_coefficients=[1],
3105                            default_charges={state_one_type: charge_number_map[state_one_type],
3106                            state_two_type: charge_number_map[state_two_type],
3107                            anion_es_type: anion_charge})
3108
3109            sucessful_reactions_labels.append(name)
3110        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):
3112    def setup_df (self):
3113        """
3114        Sets up the pyMBE's dataframe `pymbe.df`.
3115
3116        Returns:
3117            columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe
3118        """
3119        
3120        columns_dtypes = {
3121            'name': {
3122                '': str},
3123            'pmb_type': {
3124                '': str},
3125            'particle_id': {
3126                '': pd.Int64Dtype()},
3127            'particle_id2':  {
3128                '': pd.Int64Dtype()},
3129            'residue_id':  {
3130                '': pd.Int64Dtype()},
3131            'molecule_id':  {
3132                '': pd.Int64Dtype()},
3133            'acidity':  {
3134                '': str},
3135            'pka':  {
3136                '': object},
3137            'central_bead':  {
3138                '': object},
3139            'side_chains': {
3140                '': object},
3141            'residue_list': {
3142                '': object},
3143            'model': {
3144                '': str},
3145            'sigma': {
3146                '': object},
3147            'cutoff': {
3148                '': object},
3149            'offset': {
3150                '': object},
3151            'epsilon': {
3152                '': object},
3153            'state_one': {
3154                'label': str,
3155                'es_type': pd.Int64Dtype(),
3156                'z': pd.Int64Dtype()},
3157            'state_two': {
3158                'label': str,
3159                'es_type': pd.Int64Dtype(),
3160                'z': pd.Int64Dtype()},
3161            'sequence': {
3162                '': object},
3163            'bond_object': {
3164                '': object},
3165            'parameters_of_the_potential':{
3166                '': object},
3167            'l0': {
3168                '': object},
3169            }
3170        
3171        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()]))
3172        
3173        for level1, sub_dtypes in columns_dtypes.items():
3174            for level2, dtype in sub_dtypes.items():
3175                self.df[level1, level2] = self.df[level1, level2].astype(dtype)
3176                
3177        columns_names = pd.MultiIndex.from_frame(self.df)
3178        columns_names = columns_names.names
3179                
3180        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):
3182    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True):
3183        """
3184        Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`.
3185
3186        Args:
3187            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
3188            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.
3189            combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3190            warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True.
3191
3192        Note:
3193            - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 
3194            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
3195            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3196
3197        """
3198        from itertools import combinations_with_replacement
3199        import warnings
3200        implemented_combining_rules = ['Lorentz-Berthelot']
3201        compulsory_parameters_in_df = ['sigma','epsilon']
3202        # Sanity check
3203        if combining_rule not in implemented_combining_rules:
3204            raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules)
3205        shift=0
3206        if shift_potential:
3207            shift="auto"
3208        # List which particles have sigma and epsilon values defined in pmb.df and which ones don't
3209        particles_types_with_LJ_parameters = []
3210        non_parametrized_labels= []
3211        for particle_type in self.get_type_map().values():
3212            check_list=[]
3213            for key in compulsory_parameters_in_df:
3214                value_in_df=self.find_value_from_es_type(es_type=particle_type,
3215                                                        column_name=key)
3216                check_list.append(pd.isna(value_in_df))
3217            if any(check_list):
3218                non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 
3219                                                                            column_name='label'))
3220            else:
3221                particles_types_with_LJ_parameters.append(particle_type)
3222        # Set up LJ interactions between all particle types
3223        for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 
3224            particle_name1 = self.find_value_from_es_type(es_type=type_pair[0],
3225                                                        column_name="name")
3226            particle_name2 = self.find_value_from_es_type(es_type=type_pair[1],
3227                                                        column_name="name")
3228            lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1,
3229                                                 particle_name2 = particle_name2,
3230                                                 combining_rule = combining_rule)
3231            
3232            # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
3233            if not lj_parameters:
3234                continue
3235            espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 
3236                                                                                    sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 
3237                                                                                    cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude,
3238                                                                                    offset = lj_parameters["offset"].to("reduced_length").magnitude, 
3239                                                                                    shift = shift)                                                                                          
3240            index = len(self.df)
3241            label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label")
3242            label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label")
3243            self.df.at [index, 'name'] = f'LJ: {label1}-{label2}'
3244            lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params()
3245
3246            self.add_value_to_df(index=index,
3247                                key=('pmb_type',''),
3248                                new_value='LennardJones')
3249            
3250            self.add_value_to_df(index=index,
3251                                key=('parameters_of_the_potential',''),
3252                                new_value=lj_params,
3253                                non_standard_value=True)
3254        if non_parametrized_labels and warnings:
3255            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) 
3256 
3257        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):
3259    def write_pmb_df (self, filename):
3260        '''
3261        Writes the pyMBE dataframe into a csv file
3262        
3263        Args:
3264            filename(`str`): Path to the csv file 
3265        '''
3266
3267        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence']
3268        df = self.df.copy(deep=True)
3269        for column_name in columns_with_list_or_dict:
3270            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)
3271        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)
3272        df.fillna(pd.NA, inplace=True)
3273        df.to_csv(filename)
3274        return

Writes the pyMBE dataframe into a csv file

Arguments:
  • filename(str): Path to the csv file
class pymbe_library.NumpyEncoder(json.encoder.JSONEncoder):
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)

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

def default(self, obj):
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)

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