pyMBE

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

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

Arguments:
  • temperature(pint.Quantity,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None.
  • unit_length(pint.Quantity, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None.
  • unit_charge (pint.Quantity,optional): Reduced unit of charge defined using the pmb.units UnitRegistry. Defaults to None.
  • Kw (pint.Quantity,optional): Ionic product of water in mol^2/l^2. Defaults to None.
Note:
  • If no temperature is given, a value of 298.15 K is assumed by default.
  • If no unit_length is given, a value of 0.355 nm is assumed by default.
  • If no unit_charge is given, a value of 1 elementary charge is assumed by default.
  • If no Kw is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default.
seed
rng
units
N_A
kB
e
def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False):
 89    def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False):
 90        """
 91        Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles.
 92
 93        Args:
 94            particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles
 95            particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles
 96            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False.
 97
 98        Returns:
 99            index(`int`): Row index where the bond information has been added in pmb.df.
100        """
101        particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0]
102        particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0]
103        
104        bond_key = self.find_bond_key(particle_name1=particle_name1,
105                                    particle_name2=particle_name2, 
106                                    use_default_bond=use_default_bond)
107        if not bond_key:
108            return
109        self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1)
110        indexs = np.where(self.df['name']==bond_key)
111        index_list = list (indexs[0])
112        used_bond_df = self.df.loc[self.df['particle_id2'].notnull()]
113        #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict
114        used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 )
115        used_bond_index = used_bond_df.index.to_list()
116        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.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, verbose=True, non_standard_value=False, overwrite=False):
142    def add_value_to_df(self,index,key,new_value, verbose=True, 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            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
150            non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False.
151            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
152        """
153
154        token = "#protected:"
155
156        def protect(obj):
157            if non_standard_value:
158                return token + json.dumps(obj, cls=self.NumpyEncoder)
159            return obj
160
161        def deprotect(obj):
162            if non_standard_value and isinstance(obj, str) and obj.startswith(token):
163                return json.loads(obj.removeprefix(token))
164            return obj
165
166        # Make sure index is a scalar integer value
167        index = int (index)
168        assert isinstance(index, int), '`index` should be a scalar integer value.'
169        idx = pd.IndexSlice
170        if self.check_if_df_cell_has_a_value(index=index,key=key):
171            old_value= self.df.loc[index,idx[key]]
172            if protect(old_value) != protect(new_value):
173                name=self.df.loc[index,('name','')]
174                pmb_type=self.df.loc[index,('pmb_type','')]
175                if verbose:
176                    print(f"WARNING: you are attempting to redefine the properties of {name} of pmb_type {pmb_type}")    
177                if overwrite and verbose:
178                    print(f'WARNING: overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}')
179                if not overwrite:
180                    if verbose:
181                        print(f"WARNING: pyMBE has preserved of the entry `{key}`: old_value = {old_value}. If you want to overwrite it with new_value = {new_value}, activate the switch overwrite = True ")
182                    return
183        self.df.loc[index,idx[key]] = protect(new_value)
184        if non_standard_value:
185            self.df[key] = self.df[key].apply(deprotect)
186        return

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.
  • verbose(bool, optional): Switch to activate/deactivate verbose. Defaults to True.
  • 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, name, molecule_index, pmb_type, used_molecules_id):
188    def assign_molecule_id(self, name, molecule_index, pmb_type, used_molecules_id):
189        """
190        Assigns the `molecule_id` of the pmb object given by `pmb_type`
191        
192        Args:
193            name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df`
194            pmb_type(`str`): pmb_object_type to assign the `molecule_id` 
195            molecule_index(`int`): index of the current `pmb_object_type` to assign the `molecule_id`
196            used_molecules_id(`lst`): list with the `molecule_id` values already used.
197        
198        Returns:
199            molecule_id(`int`): Id of the molecule
200        """
201
202        self.clean_df_row(index=int(molecule_index))
203        
204        if self.df['molecule_id'].isnull().values.all():
205            molecule_id = 0        
206        else:
207            # check if a residue is part of another molecule
208            check_residue_name = self.df[self.df['residue_list'].astype(str).str.contains(name)]
209            mol_pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
210            if not check_residue_name.empty and mol_pmb_type == pmb_type:
211                for value in check_residue_name.index.to_list():                  
212                    if value not in used_molecules_id:                              
213                        molecule_id = self.df.loc[value].molecule_id.values[0]                    
214                        break
215            else:
216                molecule_id = self.df['molecule_id'].max() +1
217
218        self.add_value_to_df (key=('molecule_id',''),
219                                index=int(molecule_index),
220                                new_value=molecule_id, 
221                                verbose=False)
222
223        return molecule_id

Assigns the molecule_id of the pmb object given by pmb_type

Arguments:
  • name(str): Label of the molecule type to be created. name must be defined in pmb.df
  • pmb_type(str): pmb_object_type to assign the molecule_id
  • molecule_index(int): index of the current pmb_object_type to assign the molecule_id
  • used_molecules_id(lst): list with the molecule_id values already used.
Returns:

molecule_id(int): Id of the molecule

def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system):
225    def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system):
226        """
227        Calculates the center of the molecule with a given molecule_id.
228
229        Args:
230            molecule_id(`int`): Id of the molecule whose center of mass is to be calculated.
231            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
232        
233        Returns:
234            center_of_mass(`lst`): Coordinates of the center of mass.
235        """
236        center_of_mass = np.zeros(3)
237        axis_list = [0,1,2]
238        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
239        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
240        for pid in particle_id_list:
241            for axis in axis_list:
242                center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis]
243        center_of_mass = center_of_mass / len(particle_id_list)
244        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):
246    def calculate_HH(self, molecule_name, pH_list=None, pka_set=None):
247        """
248        Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 
249        for molecules with the name `molecule_name`.
250
251        Args:
252            molecule_name(`str`): name of the molecule to calculate the ideal charge for
253            pH_list(`lst`): pH-values to calculate. 
254            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
255
256        Returns:
257            Z_HH(`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list`
258
259        Note:
260            - This function supports objects with pmb types: "molecule", "peptide" and "protein".
261            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
262            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
263            - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan`  should be used.
264        """
265        if pH_list is None:
266            pH_list=np.linspace(2,12,50)
267        if pka_set is None:
268            pka_set=self.get_pka_set() 
269        self.check_pka_set(pka_set=pka_set)
270        charge_number_map = self.get_charge_number_map()
271        Z_HH=[]
272        for pH_value in pH_list:    
273            Z=0
274            index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 
275            residue_list = self.df.at [index,('residue_list','')]
276            sequence = self.df.at [index,('sequence','')]
277            if np.any(pd.isnull(sequence)):
278                # Molecule has no sequence
279                for residue in residue_list:
280                    list_of_particles_in_residue = self.search_particles_in_residue(residue)
281                    for particle in list_of_particles_in_residue:
282                        if particle in pka_set.keys():
283                            if pka_set[particle]['acidity'] == 'acidic':
284                                psi=-1
285                            elif pka_set[particle]['acidity']== 'basic':
286                                psi=+1
287                            else:
288                                psi=0
289                            Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value'])))                      
290                Z_HH.append(Z)
291            else:
292                # Molecule has a sequence
293                if not isinstance(sequence, list):
294                    # If the df has been read by file, the sequence needs to be parsed.
295                    sequence = self.parse_sequence_from_file(sequence=sequence)
296                for name in sequence:
297                    if name in pka_set.keys():
298                        if pka_set[name]['acidity'] == 'acidic':
299                            psi=-1
300                        elif pka_set[name]['acidity']== 'basic':
301                            psi=+1
302                        else:
303                            psi=0
304                        Z+=psi/(1+10**(psi*(pH_value-pka_set[name]['pka_value'])))
305                    else:
306                        state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
307                        Z+=charge_number_map[state_one_type]
308                Z_HH.append(Z)
309
310        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):
312    def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None):
313        """
314        Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning.
315
316        Args:
317            c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
318            c_salt('float'): Salt concentration in the reservoir.
319            pH_list('lst'): List of pH-values in the reservoir. 
320            pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}.
321
322        Returns:
323            {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
324            pH_system_list ('lst'): List of pH_values in the system.
325            partition_coefficients_list ('lst'): List of partition coefficients of cations.
326
327        Note:
328            - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated
329            - If no `pka_set` is given, the pKa values are taken from `pmb.df`
330            - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result.
331        """
332        if pH_list is None:
333            pH_list=np.linspace(2,12,50)
334        if pka_set is None:
335            pka_set=self.get_pka_set() 
336        self.check_pka_set(pka_set=pka_set)
337
338        partition_coefficients_list = []
339        pH_system_list = []
340        Z_HH_Donnan={}
341        for key in c_macro:
342            Z_HH_Donnan[key] = []
343
344        def calc_charges(c_macro, pH):
345            """
346            Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation.
347
348            Args:
349                c_macro('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
350                pH('float'): pH-value that is used in the HH equation.
351
352            Returns:
353                charge('dict'): {"molecule_name": charge}
354            """
355            charge = {}
356            for name in c_macro:
357                charge[name] = self.calculate_HH(name, [pH], pka_set)[0]
358            return charge
359
360        def calc_partition_coefficient(charge, c_macro):
361            """
362            Calculates the partition coefficients of positive ions according to the ideal Donnan theory.
363
364            Args:
365                charge('dict'): {"molecule_name": charge}
366                c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 
367            """
368            nonlocal ionic_strength_res
369            charge_density = 0.0
370            for key in charge:
371                charge_density += charge[key] * c_macro[key]
372            return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude
373
374        for pH_value in pH_list:    
375            # calculate the ionic strength of the reservoir
376            if pH_value <= 7.0:
377                ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 
378            elif pH_value > 7.0:
379                ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt
380
381            #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations
382            #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation.
383            #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects.
384            equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro))
385            logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq")
386            partition_coefficient = 10**logxi.root
387
388            charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient))
389            for key in c_macro:
390                Z_HH_Donnan[key].append(charges_temp[key])
391
392            pH_system_list.append(pH_value - np.log10(partition_coefficient))
393            partition_coefficients_list.append(partition_coefficient)
394
395        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):
397    def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset):
398        """
399        Calculates the initial bond length that is used when setting up molecules,
400        based on the minimum of the sum of bonded and short-range (LJ) interactions.
401        
402        Args:
403            bond_object(`espressomd.interactions.BondedInteractions`): instance of a bond object from espressomd library
404            bond_type(`str`): label identifying the used bonded potential
405            epsilon(`pint.Quantity`): LJ epsilon of the interaction between the particles
406            sigma(`pint.Quantity`): LJ sigma of the interaction between the particles
407            cutoff(`pint.Quantity`): cutoff-radius of the LJ interaction 
408            offset(`pint.Quantity`): offset of the LJ interaction
409        """    
410        def truncated_lj_potential(x, epsilon, sigma, cutoff,offset):
411            if x>cutoff:
412                return 0.0
413            else:
414                return 4*epsilon*((sigma/(x-offset))**12-(sigma/(x-offset))**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6)
415
416        epsilon_red=epsilon.to('reduced_energy').magnitude
417        sigma_red=sigma.to('reduced_length').magnitude
418        cutoff_red=cutoff.to('reduced_length').magnitude
419        offset_red=offset.to('reduced_length').magnitude
420
421        if bond_type == "harmonic":
422            r_0 = bond_object.params.get('r_0')
423            k = bond_object.params.get('k')
424            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
425        elif bond_type == "FENE":
426            r_0 = bond_object.params.get('r_0')
427            k = bond_object.params.get('k')
428            d_r_max = bond_object.params.get('d_r_max')
429            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
430        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):
432    def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False):
433        '''
434        Calculates the net charge per molecule of molecules with `name` = molecule_name. 
435        Returns the net charge per molecule and a maps with the net charge per residue and molecule.
436
437        Args:
438            espresso_system(`espressomd.system.System`): system information 
439            molecule_name(`str`): name of the molecule to calculate the net charge
440            dimensionless(`bool'): sets if the charge is returned with a dimension or not
441
442        Returns:
443            (`dict`): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }}
444
445        Note:
446            - The net charge of the molecule is averaged over all molecules of type `name` 
447            - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name`
448        '''        
449        valid_pmb_types = ["molecule", "protein"]
450        pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0]
451        if pmb_type not in valid_pmb_types:
452            raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}")      
453
454        id_map = self.get_particle_id_map(object_name=molecule_name)
455        def create_charge_map(espresso_system,id_map,label):
456            charge_number_map={}
457            for super_id in id_map[label].keys():
458                if dimensionless:
459                    net_charge=0
460                else:
461                    net_charge=0 * self.units.Quantity(1,'reduced_charge')
462                for pid in id_map[label][super_id]:
463                    if dimensionless:
464                        net_charge+=espresso_system.part.by_id(pid).q
465                    else:
466                        net_charge+=espresso_system.part.by_id(pid).q * self.units.Quantity(1,'reduced_charge')
467                charge_number_map[super_id]=net_charge
468            return charge_number_map
469        net_charge_molecules=create_charge_map(label="molecule_map",
470                                                espresso_system=espresso_system,
471                                                id_map=id_map)
472        net_charge_residues=create_charge_map(label="residue_map",
473                                                espresso_system=espresso_system,
474                                                id_map=id_map)       
475        if dimensionless:
476            mean_charge=np.mean(np.array(list(net_charge_molecules.values())))
477        else:
478            mean_charge=np.mean(np.array([value.magnitude for value in net_charge_molecules.values()]))*self.units.Quantity(1,'reduced_charge')
479        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):
481    def center_molecule_in_simulation_box(self, molecule_id, espresso_system):
482        """
483        Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`.
484        
485        Args:
486            molecule_id(`int`): Id of the molecule to be centered.
487            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
488        """
489        if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0:
490            raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.")      
491        center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system)
492        box_center = [espresso_system.box_l[0]/2.0,
493                      espresso_system.box_l[1]/2.0,
494                      espresso_system.box_l[2]/2.0]
495        molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0]
496        particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"]
497        for pid in particle_id_list:
498            es_pos = espresso_system.part.by_id(pid).pos
499            espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center
500        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):
502    def check_aminoacid_key(self, key):
503        """
504        Checks if `key` corresponds to a valid aminoacid letter code.
505
506        Args:
507            key(`str`): key to be checked.
508
509        Returns:
510            `bool`: True if `key` is a valid aminoacid letter code, False otherwise.
511        """
512        valid_AA_keys=['V', #'VAL'
513                       'I', #'ILE'
514                       'L', #'LEU'
515                       'E', #'GLU'
516                       'Q', #'GLN'
517                       'D', #'ASP'
518                       'N', #'ASN'
519                       'H', #'HIS'
520                       'W', #'TRP'
521                       'F', #'PHE'
522                       'Y', #'TYR'
523                       'R', #'ARG' 
524                       'K', #'LYS'
525                       'S', #'SER'
526                       'T', #'THR'
527                       'M', #'MET'
528                       'A', #'ALA'
529                       'G', #'GLY'
530                       'P', #'PRO'
531                       'C'] #'CYS'
532        if key in valid_AA_keys:
533            return True
534        else:
535            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):
537    def check_dimensionality(self, variable, expected_dimensionality):
538        """
539        Checks if the dimensionality of `variable` matches `expected_dimensionality`.
540
541        Args:
542            variable(`pint.Quantity`): Quantity to be checked.
543            expected_dimensionality(`str`): Expected dimension of the variable.
544
545        Returns:
546            (`bool`): `True` if the variable if of the expected dimensionality, `False` otherwise.
547
548        Note:
549            - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality).
550            - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"`    
551        """
552        correct_dimensionality=variable.check(f"{expected_dimensionality}")      
553        if not correct_dimensionality:
554            raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}")
555        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):
557    def check_if_df_cell_has_a_value(self, index,key):
558        """
559        Checks if a cell in the `pmb.df` at the specified index and column has a value.
560
561        Args:
562            index(`int`): Index of the row to check.
563            key(`str`): Column label to check.
564
565        Returns:
566            `bool`: `True` if the cell has a value, `False` otherwise.
567        """
568        idx = pd.IndexSlice
569        import warnings
570        with warnings.catch_warnings():
571            warnings.simplefilter("ignore")
572            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):
574    def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined):
575        """
576        Checks if `name` is defined in `pmb.df`.
577
578        Args:
579            name(`str`): label to check if defined in `pmb.df`.
580            pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`.
581
582        Returns:
583            `bool`: `True` for success, `False` otherwise.
584        """
585        if name in self.df['name'].unique():
586            current_object_type = self.df[self.df['name']==name].pmb_type.values[0]
587            if current_object_type != pmb_type_to_be_defined:
588                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")
589            return True            
590        else:
591            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):
593    def check_if_metal_ion(self,key):
594        """
595        Checks if `key` corresponds to a label of a supported metal ion.
596
597        Args:
598            key(`str`): key to be checked
599
600        Returns:
601            (`bool`): True if `key`  is a supported metal ion, False otherwise.
602        """
603        if key in self.get_metal_ions_charge_number_map().keys():
604            return True
605        else:
606            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):
608    def check_pka_set(self, pka_set):
609        """
610        Checks that `pka_set` has the formatting expected by the pyMBE library.
611       
612        Args:
613            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
614        """
615        required_keys=['pka_value','acidity']
616        for required_key in required_keys:
617            for pka_entry in pka_set.values():
618                if required_key not in pka_entry.keys():
619                    raise ValueError(f'missing a required key "{required_key}" in the following entry of pka_set "{pka_entry}"')
620        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')):
622    def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")):
623        """
624        Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a np.nan value.
625
626        Args:
627            index(`int`): Index of the row to clean.
628            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`].
629        """   
630        for column_key in columns_keys_to_clean:
631            self.add_value_to_df(key=(column_key,''),index=index,new_value=np.nan, verbose=False)
632        return

Cleans the columns of pmb.df in columns_keys_to_clean of the row with index index by assigning them a np.nan 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):
634    def convert_columns_to_original_format (self, df):
635        """
636        Converts the columns of the Dataframe to the original format in pyMBE.
637        
638        Args:
639            df(`DataFrame`): dataframe with pyMBE information as a string  
640        
641        """
642
643        columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', 'model',('state_one','es_type'),('state_two','es_type'),('state_one','z'),('state_two','z') ]  
644
645        columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset']
646
647        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence']
648
649        for column_name in columns_dtype_int:
650            df[column_name] = df[column_name].astype(object)
651            
652        for column_name in columns_with_list_or_dict:
653            if df[column_name].isnull().all():
654                df[column_name] = df[column_name].astype(object)
655            else:
656                df[column_name] = df[column_name].apply(lambda x: ast.literal_eval(str(x)) if pd.notnull(x) else x)
657
658        for column_name in columns_with_units:
659            df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x)
660
661        df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x)
662
663        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):
665    def convert_str_to_bond_object (self, bond_str):
666        
667        """
668        Convert a row read as a `str` to the corresponding bond object. There are two supported bonds: HarmonicBond and FeneBond
669
670        Args:
671            bond_str(`str`): string with the information of a bond object
672
673        Returns:
674            bond_object(`obj`): EsPRESSo bond object 
675        """
676        
677        from espressomd.interactions import HarmonicBond
678        from espressomd.interactions import FeneBond
679
680        supported_bonds = ['HarmonicBond', 'FeneBond']
681
682        for bond in supported_bonds:
683            variable = re.subn(f'{bond}', '', bond_str)
684            if variable[1] == 1: 
685                params = ast.literal_eval(variable[0])
686                if bond == 'HarmonicBond':
687                    bond_object = HarmonicBond(r_cut =params['r_cut'], k = params['k'], r_0=params['r_0'])
688                elif bond == 'FeneBond':
689                    bond_object = FeneBond(k = params['k'], d_r_max =params['d_r_max'], r_0=params['r_0'])
690
691        return bond_object 

Convert a row read as a str to the corresponding bond object. There are two supported bonds: HarmonicBond and FeneBond

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

bond_object(obj): EsPRESSo bond object

def copy_df_entry(self, name, column_name, number_of_copies):
693    def copy_df_entry(self, name, column_name, number_of_copies):
694        '''
695        Creates 'number_of_copies' of a given 'name' in `pymbe.df`.
696
697        Args:
698            name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df`
699            column_name(`str`): Column name to use as a filter. 
700            number_of_copies(`int`): number of copies of `name` to be created.
701        
702        Note:
703            - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 
704        '''
705
706        valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ]
707        if column_name not in valid_column_names:
708            raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}")
709        df_by_name = self.df.loc[self.df.name == name]
710        if number_of_copies != 1:           
711            if df_by_name[column_name].isnull().values.any():       
712                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True)
713            else:
714                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
715                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
716                df_by_name_repeated[column_name] =np.nan
717            # Concatenate the new particle rows to  `df`
718            self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
719        else:
720            if not df_by_name[column_name].isnull().values.any():     
721                df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 
722                df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True)
723                df_by_name_repeated[column_name] =np.nan
724                self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True)
725        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):
727    def create_added_salt (self, espresso_system, cation_name, anion_name, c_salt, verbose=True):    
728        """
729        Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`.
730
731        Args:
732            espresso_system(`espressomd.system.System`): instance of an espresso system object.
733            cation_name(`str`): `name` of a particle with a positive charge.
734            anion_name(`str`): `name` of a particle with a negative charge.
735            c_salt(`float`): Salt concentration.
736            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
737            
738        Returns:
739            c_salt_calculated(`float`): Calculated salt concentration added to `espresso_system`.
740        """
741        cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
742        anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
743        if cation_name_charge <= 0:
744            raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge)
745        if anion_name_charge >= 0:
746            raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge)
747        # Calculate the number of ions in the simulation box
748        volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3')
749        if c_salt.check('[substance] [length]**-3'):
750            N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude)
751            c_salt_calculated=N_ions/(volume*self.N_A)
752        elif c_salt.check('[length]**-3'):
753            N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude)
754            c_salt_calculated=N_ions/volume
755        else:
756            raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt)
757        N_cation = N_ions*abs(anion_name_charge)
758        N_anion = N_ions*abs(cation_name_charge)
759        self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation)
760        self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion)
761        if verbose:
762            if c_salt_calculated.check('[substance] [length]**-3'):
763                print(f"\n Added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions")
764            elif c_salt_calculated.check('[length]**-3'):
765                print(f"\n Added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions")
766        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):
768    def create_bond_in_espresso(self, bond_type, bond_parameters):
769        '''
770        Creates either a harmonic or a FENE bond in ESPREesSo
771
772        Args:
773            bond_type(`str`): label to identify the potential to model the bond.
774            bond_parameters(`dict`): parameters of the potential of the bond.
775
776        Note:
777            Currently, only HARMONIC and FENE bonds are supported.
778
779            For a HARMONIC bond the dictionary must contain:
780
781                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
782                using the `pmb.units` UnitRegistry.
783                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
784                the `pmb.units` UnitRegistry.
785           
786            For a FENE bond the dictionay must additionally contain:
787                
788                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
789                units of length using the `pmb.units` UnitRegistry. Default 'None'.
790
791        Returns:
792              bond_object (`obj`): a harmonic or a FENE bond object in ESPREesSo
793        '''
794        from espressomd import interactions
795
796        valid_bond_types   = ["harmonic", "FENE"]
797        
798        if 'k' in bond_parameters:
799            bond_magnitude     = bond_parameters['k'].to('reduced_energy / reduced_length**2')
800        else:
801            raise ValueError("Magnitud of the potential (k) is missing")
802        
803        if bond_type == 'harmonic':
804            if 'r_0' in bond_parameters:
805                bond_length        = bond_parameters['r_0'].to('reduced_length')
806            else:
807                raise ValueError("Equilibrium bond length (r_0) is missing")
808            bond_object    = interactions.HarmonicBond(k   = bond_magnitude.magnitude,
809                                                       r_0 = bond_length.magnitude)
810        elif bond_type == 'FENE':
811            if 'r_0' in bond_parameters:
812                bond_length        = bond_parameters['r_0'].to('reduced_length').magnitude
813            else:
814                print("WARNING: No value provided for r_0. Defaulting to r_0 = 0")
815                bond_length=0
816            if 'd_r_max' in bond_parameters:
817                max_bond_stret = bond_parameters['d_r_max'].to('reduced_length')
818            else:
819                raise ValueError("Maximal stretching length (d_r_max) is missing")
820            bond_object    = interactions.FeneBond(r_0     = bond_length, 
821                                                   k       = bond_magnitude.magnitude,
822                                                   d_r_max = max_bond_stret.magnitude)
823        else:
824            raise ValueError(f"Bond type {bond_type} currently not implemented in pyMBE, valid types are {valid_bond_types}")
825
826        return bond_object

Creates either a harmonic or a FENE bond in ESPREesSo

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 dictionay 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): a harmonic or a FENE bond object in ESPREesSo

def create_counterions( self, object_name, cation_name, anion_name, espresso_system, verbose=True):
829    def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True):
830        """
831        Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`.
832        
833        Args:
834            object_name(`str`): `name` of a pymbe object.
835            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
836            cation_name(`str`): `name` of a particle with a positive charge.
837            anion_name(`str`): `name` of a particle with a negative charge.
838            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
839
840        Returns: 
841            counterion_number(`dict`): {"name": number}
842        """
843        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0]
844        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0]
845        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
846        counterion_number={}
847        object_charge={}
848        for name in ['positive', 'negative']:
849            object_charge[name]=0
850        for id in object_ids:
851            if espresso_system.part.by_id(id).q > 0:
852                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
853            elif espresso_system.part.by_id(id).q < 0:
854                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
855        if object_charge['positive'] % abs(anion_charge) == 0:
856            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
857        else:
858            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
859        if object_charge['negative'] % abs(cation_charge) == 0:
860            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
861        else:
862            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
863        if counterion_number[cation_name] > 0: 
864            self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name])
865        else:
866            counterion_number[cation_name]=0
867        if counterion_number[anion_name] > 0:
868            self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name])
869        else:
870            counterion_number[anion_name] = 0
871        if verbose:
872            print('The following counter-ions have been created: ')
873            for name in counterion_number.keys():
874                print(f'Ion type: {name} created number: {counterion_number[name]}')
875        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):
877    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False):
878        """
879        Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
880
881        Args:
882            name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df`
883            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
884            number_of_molecules(`int`): Number of molecules of type `name` to be created.
885            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.
886            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. 
887            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df`
888
889        Returns:
890            molecules_info(`dict`):  {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 
891        """
892        if list_of_first_residue_positions is not None:
893            for item in list_of_first_residue_positions:
894                if not isinstance(item, list):
895                    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.")
896                elif len(item) != 3:
897                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
898
899            if len(list_of_first_residue_positions) != number_of_molecules:
900                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
901        if number_of_molecules <= 0:
902            return 0
903        self.check_if_name_is_defined_in_df(name=name,
904                                            pmb_type_to_be_defined='molecule')
905            
906        first_residue = True
907        molecules_info = {}
908        residue_list = self.df[self.df['name']==name].residue_list.values [0]
909        self.copy_df_entry(name=name,
910                        column_name='molecule_id',
911                        number_of_copies=number_of_molecules)
912
913        molecules_index = np.where(self.df['name']==name)
914        molecule_index_list =list(molecules_index[0])[-number_of_molecules:]
915        used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist()
916        pos_index = 0 
917        for molecule_index in molecule_index_list:        
918            molecule_id = self.assign_molecule_id(name=name,
919                                                pmb_type='molecule',
920                                                used_molecules_id=used_molecules_id,
921                                                molecule_index=molecule_index)
922            molecules_info[molecule_id] = {}
923            for residue in residue_list:
924                if first_residue:
925                    if list_of_first_residue_positions is None:
926                        residue_position = None
927                    else:
928                        for item in list_of_first_residue_positions:
929                            residue_position = [np.array(list_of_first_residue_positions[pos_index])]
930                    # Generate an arbitrary random unit vector
931                    if backbone_vector is None:
932                        backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
933                                                                radius=1, 
934                                                                n_samples=1,
935                                                                on_surface=True)[0]
936                    residues_info = self.create_residue(name=residue,
937                                                        espresso_system=espresso_system, 
938                                                        central_bead_position=residue_position,  
939                                                        use_default_bond= use_default_bond, 
940                                                        backbone_vector=backbone_vector)
941                    residue_id = next(iter(residues_info))
942                    # Add the correct molecule_id to all particles in the residue
943                    for index in self.df[self.df['residue_id']==residue_id].index:
944                        self.add_value_to_df(key=('molecule_id',''),
945                                            index=int (index),
946                                            new_value=molecule_id,
947                                            overwrite=True,
948                                            verbose=False)
949                    central_bead_id = residues_info[residue_id]['central_bead_id']
950                    previous_residue = residue
951                    residue_position = espresso_system.part.by_id(central_bead_id).pos
952                    previous_residue_id = central_bead_id
953                    first_residue = False          
954                else:                    
955                    previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0]
956                    new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0]       
957                    bond = self.search_bond(particle_name1=previous_central_bead_name, 
958                                            particle_name2=new_central_bead_name, 
959                                            hard_check=True, 
960                                            use_default_bond=use_default_bond)                
961                    l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 
962                                            particle_name2=new_central_bead_name, 
963                                            hard_check=True, 
964                                            use_default_bond=use_default_bond)                
965                    residue_position = residue_position+backbone_vector*l0
966                    residues_info = self.create_residue(name=residue, 
967                                                        espresso_system=espresso_system, 
968                                                        central_bead_position=[residue_position],
969                                                        use_default_bond= use_default_bond, 
970                                                        backbone_vector=backbone_vector)
971                    residue_id = next(iter(residues_info))      
972                    for index in self.df[self.df['residue_id']==residue_id].index:
973                        self.add_value_to_df(key=('molecule_id',''),
974                                            index=int (index),
975                                            new_value=molecule_id,
976                                            verbose=False,
977                                            overwrite=True)            
978                    central_bead_id = residues_info[residue_id]['central_bead_id']
979                    espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id))
980                    bond_index = self.add_bond_in_df(particle_id1=central_bead_id,
981                                        particle_id2=previous_residue_id,
982                                        use_default_bond=use_default_bond) 
983                    self.add_value_to_df(key=('molecule_id',''),
984                                            index=int (bond_index),
985                                            new_value=molecule_id,
986                                            verbose=False,
987                                            overwrite=True)           
988                    previous_residue_id = central_bead_id
989                    previous_residue = residue                    
990                molecules_info[molecule_id][residue_id] = residues_info[residue_id]
991            first_residue = True
992            pos_index+=1
993        
994        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):
 996    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
 997        """
 998        Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`.
 999        
1000        Args:
1001            name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 
1002            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1003            number_of_particles(`int`): Number of particles to be created.
1004            position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
1005            fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
1006        Returns:
1007            created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`.
1008        """       
1009        if number_of_particles <=0:
1010            return []
1011        self.check_if_name_is_defined_in_df(name=name,
1012                                       pmb_type_to_be_defined='particle')
1013        # Copy the data of the particle `number_of_particles` times in the `df`
1014        self.copy_df_entry(name=name,
1015                          column_name='particle_id',
1016                          number_of_copies=number_of_particles)
1017        # Get information from the particle type `name` from the df     
1018        z = self.df.loc[self.df['name']==name].state_one.z.values[0]
1019        z = 0. if z is None else z
1020        es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
1021        # Get a list of the index in `df` corresponding to the new particles to be created
1022        index = np.where(self.df['name']==name)
1023        index_list =list(index[0])[-number_of_particles:]
1024        # Create the new particles into  `espresso_system`
1025        created_pid_list=[]
1026        for index in range (number_of_particles):
1027            df_index=int (index_list[index])
1028            self.clean_df_row(index=df_index)
1029            if position is None:
1030                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1031            else:
1032                particle_position = position[index]
1033            if len(espresso_system.part.all()) == 0:
1034                bead_id = 0
1035            else:
1036                bead_id = max (espresso_system.part.all().id) + 1
1037            created_pid_list.append(bead_id)
1038            kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z)
1039            if fix:
1040                kwargs["fix"] = 3 * [fix]
1041            espresso_system.part.add(**kwargs)
1042            self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id, verbose=False)                  
1043        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):
1045    def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None):
1046        """
1047        Creates all `particle`s associated to `pmb object` into  `espresso` a number of times equal to `number_of_objects`.
1048        
1049        Args:
1050            name(`str`): Unique label of the `pmb object` to be created. 
1051            number_of_objects(`int`): Number of `pmb object`s to be created.
1052            espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library.
1053            position(`list`): Coordinates where the particles should be created.
1054            use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`.
1055            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. 
1056
1057        Note:
1058            - 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. 
1059        """
1060        allowed_objects=['particle','residue','molecule']
1061        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1062        if pmb_type not in allowed_objects:
1063            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1064        if pmb_type == 'particle':
1065            self.create_particle(name=name, 
1066                                number_of_particles=number_of_objects, 
1067                                espresso_system=espresso_system, 
1068                                position=position)
1069        elif pmb_type == 'residue':
1070            self.create_residue(name=name,  
1071                                espresso_system=espresso_system, 
1072                                central_bead_position=position,
1073                                use_default_bond=use_default_bond,
1074                                backbone_vector=backbone_vector)
1075        elif pmb_type == 'molecule':
1076            self.create_molecule(name=name, 
1077                                number_of_molecules=number_of_objects, 
1078                                espresso_system=espresso_system, 
1079                                use_default_bond=use_default_bond, 
1080                                list_of_first_residue_positions=position,
1081                                backbone_vector=backbone_vector)
1082        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):
1084    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1085        """
1086        Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions`
1087
1088        Args:
1089            name(`str`): Label of the protein to be created. 
1090            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1091            number_of_proteins(`int`): Number of proteins to be created.
1092            positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}}
1093        """
1094
1095        if number_of_proteins <=0:
1096            return
1097        self.check_if_name_is_defined_in_df(name=name,
1098                                        pmb_type_to_be_defined='protein')
1099        self.copy_df_entry(name=name,
1100                            column_name='molecule_id',
1101                            number_of_copies=number_of_proteins)
1102        protein_index = np.where(self.df['name']==name)
1103        protein_index_list =list(protein_index[0])[-number_of_proteins:]
1104        used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist()
1105        
1106        box_half=espresso_system.box_l[0]/2.0
1107
1108        for molecule_index in protein_index_list:     
1109
1110            molecule_id = self.assign_molecule_id (name=name,pmb_type='protein',used_molecules_id=used_molecules_id,molecule_index=molecule_index)
1111
1112            protein_center = self.generate_coordinates_outside_sphere(radius = 1, 
1113                                                                        max_dist=box_half, 
1114                                                                        n_samples=1, 
1115                                                                        center=[box_half]*3)[0]
1116   
1117            for residue in topology_dict.keys():
1118
1119                residue_name = re.split(r'\d+', residue)[0]
1120                residue_number = re.split(r'(\d+)', residue)[1]
1121                residue_position = topology_dict[residue]['initial_pos']
1122                position = residue_position + protein_center
1123
1124                particle_id = self.create_particle(name=residue_name,
1125                                                            espresso_system=espresso_system,
1126                                                            number_of_particles=1,
1127                                                            position=[position], 
1128                                                            fix = True)
1129                
1130                index = self.df[self.df['particle_id']==particle_id[0]].index.values[0]
1131                
1132                self.add_value_to_df(key=('residue_id',''),
1133                                            index=int (index),
1134                                            new_value=residue_number,
1135                                            overwrite=True,
1136                                            verbose=False)
1137
1138                self.add_value_to_df(key=('molecule_id',''),
1139                                        index=int (index),
1140                                        new_value=molecule_id,
1141                                        overwrite=True,
1142                                        verbose=False)
1143
1144        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):
1146    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None):
1147        """
1148        Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1149
1150        Args:
1151            name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df`
1152            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1153            central_bead_position(`list` of `float`): Position of the central bead.
1154            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`
1155            backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`.
1156
1157        Returns:
1158            residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1159        """
1160        self.check_if_name_is_defined_in_df(name=name,
1161                                            pmb_type_to_be_defined='residue')
1162        # Copy the data of a residue in the `df
1163        self.copy_df_entry(name=name,
1164                            column_name='residue_id',
1165                            number_of_copies=1)
1166        residues_index = np.where(self.df['name']==name)
1167        residue_index_list =list(residues_index[0])[-1:]
1168        for residue_index in residue_index_list:  
1169            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1170            for side_chain_element in side_chain_list:
1171                if side_chain_element not in self.df.values:              
1172                    raise ValueError (side_chain_element +'is not defined')
1173        # Internal bookkepping of the residue info (important for side-chain residues)
1174        # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1175        residues_info={}
1176        for residue_index in residue_index_list:     
1177            self.clean_df_row(index=int(residue_index))
1178            # Assign a residue_id
1179            if self.df['residue_id'].isnull().all():
1180                residue_id=0
1181            else:
1182                residue_id = self.df['residue_id'].max() + 1
1183            self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id, verbose=False)
1184            # create the principal bead   
1185            central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0]            
1186            central_bead_id = self.create_particle(name=central_bead_name,
1187                                                                espresso_system=espresso_system,
1188                                                                position=central_bead_position,
1189                                                                number_of_particles = 1)[0]
1190            central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1191            #assigns same residue_id to the central_bead particle created.
1192            index = self.df[self.df['particle_id']==central_bead_id].index.values[0]
1193            self.df.at [index,'residue_id'] = residue_id
1194            # Internal bookkeeping of the central bead id
1195            residues_info[residue_id]={}
1196            residues_info[residue_id]['central_bead_id']=central_bead_id
1197            # create the lateral beads  
1198            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1199            side_chain_beads_ids = []
1200            for side_chain_element in side_chain_list:
1201                
1202                pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 
1203                if pmb_type == 'particle':
1204                    bond = self.search_bond(particle_name1=central_bead_name, 
1205                                            particle_name2=side_chain_element, 
1206                                            hard_check=True, 
1207                                            use_default_bond=use_default_bond)
1208                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1209                                              particle_name2=side_chain_element, 
1210                                              hard_check=True, 
1211                                              use_default_bond=use_default_bond)
1212
1213                    if backbone_vector is None:
1214                        bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1215                                                                 radius=l0, 
1216                                                                 n_samples=1,
1217                                                                 on_surface=True)[0]
1218                    else:
1219                        bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1220                                                                                                    magnitude=l0)
1221                    
1222                    side_bead_id = self.create_particle(name=side_chain_element, 
1223                                                                    espresso_system=espresso_system,
1224                                                                    position=[bead_position], 
1225                                                                    number_of_particles=1)[0]
1226                    index = self.df[self.df['particle_id']==side_bead_id].index.values[0]
1227                    self.add_value_to_df(key=('residue_id',''),
1228                                        index=int (index),
1229                                        new_value=residue_id, 
1230                                        verbose=False,
1231                                        overwrite=True)
1232                    side_chain_beads_ids.append(side_bead_id)
1233                    espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id))
1234                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1235                                        particle_id2=side_bead_id,
1236                                        use_default_bond=use_default_bond)
1237                    self.add_value_to_df(key=('residue_id',''),
1238                                        index=int (index),
1239                                        new_value=residue_id, 
1240                                        verbose=False,
1241                                        overwrite=True)
1242
1243                elif pmb_type == 'residue':
1244                    central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0]
1245                    bond = self.search_bond(particle_name1=central_bead_name, 
1246                                            particle_name2=central_bead_side_chain, 
1247                                            hard_check=True, 
1248                                            use_default_bond=use_default_bond)
1249                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1250                                              particle_name2=central_bead_side_chain, 
1251                                              hard_check=True, 
1252                                              use_default_bond=use_default_bond)
1253                    if backbone_vector is None:
1254                        residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1255                                                                    radius=l0, 
1256                                                                    n_samples=1,
1257                                                                    on_surface=True)[0]
1258                    else:
1259                        residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1260                                                                                                        magnitude=l0)
1261                    lateral_residue_info = self.create_residue(name=side_chain_element, 
1262                                                                espresso_system=espresso_system,
1263                                                                central_bead_position=[residue_position],
1264                                                                use_default_bond=use_default_bond)
1265                    lateral_residue_dict=list(lateral_residue_info.values())[0]
1266                    central_bead_side_chain_id=lateral_residue_dict['central_bead_id']
1267                    lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids']
1268                    residue_id_side_chain=list(lateral_residue_info.keys())[0]
1269                    # Change the residue_id of the residue in the side chain to the one of the bigger residue
1270                    index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0]
1271                    self.add_value_to_df(key=('residue_id',''),
1272                                        index=int(index),
1273                                        new_value=residue_id, 
1274                                        verbose=False,
1275                                        overwrite=True)
1276                    # Change the residue_id of the particles in the residue in the side chain
1277                    side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids
1278                    for particle_id in side_chain_beads_ids:
1279                        index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0]
1280                        self.add_value_to_df(key=('residue_id',''),
1281                                            index=int (index),
1282                                            new_value=residue_id, 
1283                                            verbose=False,
1284                                            overwrite=True)
1285                    espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id))
1286                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1287                                        particle_id2=central_bead_side_chain_id,
1288                                        use_default_bond=use_default_bond)
1289                    self.add_value_to_df(key=('residue_id',''),
1290                                        index=int (index),
1291                                        new_value=residue_id, 
1292                                        verbose=False,
1293                                        overwrite=True)
1294                    # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue
1295                    for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index:        
1296                        self.add_value_to_df(key=('residue_id',''),
1297                                            index=int(index),
1298                                            new_value=residue_id, 
1299                                            verbose=False,
1300                                            overwrite=True)
1301            # Internal bookkeeping of the side chain beads ids
1302            residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids
1303        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):
1305    def create_variable_with_units(self, variable):
1306        """
1307        Returns a pint object with the value and units defined in `variable`.
1308
1309        Args:
1310            variable(`dict` or `str`): {'value': value, 'units': units}
1311        Returns:
1312            variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry.
1313        """        
1314        
1315        if isinstance(variable, dict):
1316
1317            value=variable.pop('value')
1318            units=variable.pop('units')
1319
1320        elif isinstance(variable, str):
1321
1322            value = float(re.split(r'\s+', variable)[0])
1323            units = re.split(r'\s+', variable)[1]
1324        
1325        variable_with_units=value*self.units(units)
1326
1327        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):
1329    def define_AA_residues(self, sequence, model):
1330        """
1331        Defines in `pmb.df` all the different residues in `sequence`.
1332
1333        Args:
1334            sequence(`lst`):  Sequence of the peptide or protein.
1335            model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1336
1337        Returns:
1338            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.             
1339        """
1340        residue_list = []
1341        for residue_name in sequence:
1342            if model == '1beadAA':
1343                central_bead = residue_name
1344                side_chains = []
1345            elif model == '2beadAA':
1346                if residue_name in ['c','n', 'G']: 
1347                    central_bead = residue_name
1348                    side_chains = []
1349                else:
1350                    central_bead = 'CA'              
1351                    side_chains = [residue_name]
1352            if residue_name not in residue_list:   
1353                self.define_residue(name = 'AA-'+residue_name, 
1354                                    central_bead = central_bead,
1355                                    side_chains = side_chains)              
1356            residue_list.append('AA-'+residue_name)
1357        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):
1359    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1360        
1361        '''
1362        Defines a pmb object of type `bond` in `pymbe.df`.
1363
1364        Args:
1365            bond_type(`str`): label to identify the potential to model the bond.
1366            bond_parameters(`dict`): parameters of the potential of the bond.
1367            particle_pairs(`lst`): list of the `names` of the `particles` to be bonded.
1368
1369        Note:
1370            Currently, only HARMONIC and FENE bonds are supported.
1371
1372            For a HARMONIC bond the dictionary must contain the following parameters:
1373
1374                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
1375                using the `pmb.units` UnitRegistry.
1376                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
1377                the `pmb.units` UnitRegistry.
1378           
1379            For a FENE bond the dictionay must contain the same parameters as for a HARMONIC bond and:
1380                
1381                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
1382                units of length using the `pmb.units` UnitRegistry. Default 'None'.
1383        '''
1384
1385        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1386
1387        for particle_name1, particle_name2 in particle_pairs:
1388
1389            lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1,
1390                                                 particle_name2 = particle_name2,
1391                                                 combining_rule = 'Lorentz-Berthelot')
1392
1393            l0 = self.calculate_initial_bond_length(bond_object = bond_object,
1394                                                    bond_type   = bond_type,
1395                                                    epsilon     = lj_parameters["epsilon"],
1396                                                    sigma       = lj_parameters["sigma"],
1397                                                    cutoff      = lj_parameters["cutoff"],
1398                                                    offset      = lj_parameters["offset"],)
1399            index = len(self.df)
1400            for label in [particle_name1+'-'+particle_name2,particle_name2+'-'+particle_name1]:
1401                self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond")
1402            self.df.at [index,'name']= particle_name1+'-'+particle_name2
1403            self.df.at [index,'bond_object'] = bond_object
1404            self.df.at [index,'l0'] = l0
1405            self.add_value_to_df(index=index,
1406                                    key=('pmb_type',''),
1407                                    new_value='bond')
1408            self.add_value_to_df(index=index,
1409                                    key=('parameters_of_the_potential',''),
1410                                    new_value=bond_object.get_params(),
1411                                    non_standard_value=True)
1412
1413        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 dictionay 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):
1416    def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None):
1417        """
1418        Asigns `bond` in `pmb.df` as the default bond.
1419        The LJ parameters can be optionally provided to calculate the initial bond length
1420
1421        Args:
1422            bond_type(`str`): label to identify the potential to model the bond.
1423            bond_parameters(`dict`): parameters of the potential of the bond.
1424            sigma(`float`, optional): LJ sigma of the interaction between the particles.
1425            epsilon(`float`, optional): LJ epsilon for the interaction between the particles.
1426            cutoff(`float`, optional): cutoff-radius of the LJ interaction.
1427            offset(`float`, optional): offset of the LJ interaction.
1428
1429        Note:
1430            - Currently, only harmonic and FENE bonds are supported. 
1431        """
1432        
1433        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1434
1435        if epsilon is None:
1436            epsilon=1*self.units('reduced_energy')
1437        if sigma is None:
1438            sigma=1*self.units('reduced_length')
1439        if cutoff is None:
1440            cutoff=2**(1.0/6.0)*self.units('reduced_length')
1441        if offset is None:
1442            offset=0*self.units('reduced_length')
1443        l0 = self.calculate_initial_bond_length(bond_object = bond_object, 
1444                                                bond_type   = bond_type, 
1445                                                epsilon     = epsilon,
1446                                                sigma       = sigma,
1447                                                cutoff      = cutoff,
1448                                                offset      = offset)
1449
1450        if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'):
1451            return
1452        if len(self.df.index) != 0:
1453            index = max(self.df.index)+1
1454        else:
1455            index = 0
1456        self.df.at [index,'name']        = 'default'
1457        self.df.at [index,'bond_object'] = bond_object
1458        self.df.at [index,'l0']          = l0
1459        self.add_value_to_df(index       = index,
1460                            key          = ('pmb_type',''),
1461                            new_value    = 'bond')
1462        self.add_value_to_df(index       = index,
1463                            key          = ('parameters_of_the_potential',''),
1464                            new_value    = bond_object.get_params(),
1465                            non_standard_value=True)
1466        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):
1468    def define_molecule(self, name, residue_list):
1469        """
1470        Defines a pyMBE object of type `molecule` in `pymbe.df`.
1471
1472        Args:
1473            name(`str`): Unique label that identifies the `molecule`.
1474            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.  
1475        """
1476        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'):
1477            return
1478        index = len(self.df)
1479        self.df.at [index,'name'] = name
1480        self.df.at [index,'pmb_type'] = 'molecule'
1481        self.df.at [index,('residue_list','')] = residue_list
1482        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):
1484    def define_particle_entry_in_df(self,name):
1485        """
1486        Defines a particle entry in pmb.df.
1487
1488        Args:
1489            name(`str`): Unique label that identifies this particle type.
1490
1491        Returns:
1492            index(`int`): Index of the particle in pmb.df  
1493        """
1494
1495        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'):
1496            index = self.df[self.df['name']==name].index[0]                                   
1497        else:
1498            index = len(self.df)
1499            self.df.at [index, 'name'] = name
1500            self.df.at [index,'pmb_type'] = 'particle'
1501        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=None, pka=None, sigma=None, epsilon=None, cutoff=None, offset=None, verbose=True, overwrite=False):
1503    def define_particle(self, name, z=0, acidity=None, pka=None, sigma=None, epsilon=None, cutoff=None, offset=None,verbose=True,overwrite=False):
1504        """
1505        Defines the properties of a particle object.
1506
1507        Args:
1508            name(`str`): Unique label that identifies this particle type.  
1509            z(`int`, optional): Permanent charge number of this particle type. Defaults to 0.
1510            acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None.
1511            pka(`float`, optional): If `particle` is an acid or a base, it defines its  pka-value. Defaults to None.
1512            sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
1513            cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
1514            offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
1515            epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to None.
1516            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
1517            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
1518
1519        Note:
1520            - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units.
1521            - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units.
1522            - `cutoff` defaults to `2**(1./6.) reduced_length`. 
1523            - `offset` defaults to 0.
1524            - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions.
1525            - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`.
1526        """ 
1527        index=self.define_particle_entry_in_df(name=name)
1528        
1529        # If `cutoff` and `offset` are not defined, default them to the following values
1530        if cutoff is None:
1531            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1532        if offset is None:
1533            offset=self.units.Quantity(0, "reduced_length")
1534
1535        # Define LJ parameters
1536        parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"},
1537                                        "cutoff":{"value": cutoff, "dimensionality": "[length]"},
1538                                        "offset":{"value": offset, "dimensionality": "[length]"},
1539                                        "epsilon":{"value": epsilon, "dimensionality": "[energy]"},}
1540
1541        for parameter_key in parameters_with_dimensionality.keys():
1542            if parameters_with_dimensionality[parameter_key]["value"] is not None:
1543                self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 
1544                                          expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"])
1545                self.add_value_to_df(key=(parameter_key,''),
1546                                    index=index,
1547                                    new_value=parameters_with_dimensionality[parameter_key]["value"],
1548                                    verbose=verbose,
1549                                    overwrite=overwrite)
1550
1551        # Define particle acid/base properties
1552        self.set_particle_acidity(name=name, 
1553                                acidity=acidity, 
1554                                default_charge_number=z, 
1555                                pka=pka,
1556                                verbose=verbose,
1557                                overwrite=overwrite)
1558        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 None.
  • pka(float, optional): If particle is an acid or a base, it defines its pka-value. Defaults to None.
  • sigma(pint.Quantity, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
  • cutoff(pint.Quantity, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
  • offset(pint.Quantity, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
  • epsilon(pint.Quantity, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to None.
  • verbose(bool, optional): Switch to activate/deactivate verbose. Defaults to True.
  • 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, verbose=True):
1560    def define_particles(self, parameters, overwrite=False, verbose=True):
1561        '''
1562        Defines a particle object in pyMBE for each particle name in `particle_names`
1563
1564        Args:
1565            parameters(`dict`):  dictionary with the particle parameters. 
1566            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1567            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
1568
1569        Note:
1570            - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},}
1571        '''
1572        if not parameters:
1573            return 0
1574        for particle_name in parameters.keys():
1575            parameters[particle_name]["overwrite"]=overwrite
1576            parameters[particle_name]["verbose"]=verbose
1577            self.define_particle(**parameters[particle_name])
1578        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.
  • verbose (bool, optional): Switch to activate/deactivate verbose. Defaults to True.
Note:
  • parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},}
def define_peptide(self, name, sequence, model):
1580    def define_peptide(self, name, sequence, model):
1581        """
1582        Defines a pyMBE object of type `peptide` in the `pymbe.df`.
1583
1584        Args:
1585            name (`str`): Unique label that identifies the `peptide`.
1586            sequence (`string`): Sequence of the `peptide`.
1587            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1588        """
1589        if not self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'):
1590            valid_keys = ['1beadAA','2beadAA']
1591            if model not in valid_keys:
1592                raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
1593            clean_sequence = self.protein_sequence_parser(sequence=sequence)    
1594            residue_list = self.define_AA_residues(sequence=clean_sequence,
1595                                                    model=model)
1596            self.define_molecule(name = name, residue_list=residue_list)
1597            index = self.df.loc[self.df['name'] == name].index.item() 
1598            self.df.at [index,'model'] = model
1599            self.df.at [index,('sequence','')] = clean_sequence
1600        return

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, verbose=True):
1602    def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False, verbose=True):
1603        """
1604        Defines a globular protein pyMBE object  in `pymbe.df`.
1605
1606        Args:
1607            name (`str`): Unique label that identifies the protein.
1608            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1609            topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value}
1610            lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca".
1611            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1612            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
1613
1614        Note:
1615            - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential.
1616        """
1617
1618        # Sanity checks
1619        valid_model_keys = ['1beadAA','2beadAA']
1620        valid_lj_setups = ["wca"]
1621        if model not in valid_model_keys:
1622            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
1623        if lj_setup_mode not in valid_lj_setups:
1624            raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}')
1625        if lj_setup_mode == "wca":
1626            sigma = 1*self.units.Quantity("reduced_length")
1627            epsilon = 1*self.units.Quantity("reduced_energy")
1628        part_dict={}
1629        sequence=[]
1630        metal_ions_charge_number_map=self.get_metal_ions_charge_number_map()
1631        for particle in topology_dict.keys():
1632            particle_name = re.split(r'\d+', particle)[0] 
1633            if particle_name not in part_dict.keys():
1634                if lj_setup_mode == "wca":
1635                    part_dict[particle_name]={"sigma": sigma,
1636                                        "offset": topology_dict[particle]['radius']*2-sigma,
1637                                        "epsilon": epsilon,
1638                                        "name": particle_name}
1639                if self.check_if_metal_ion(key=particle_name):
1640                    z=metal_ions_charge_number_map[particle_name]
1641                else:
1642                    z=0
1643                part_dict[particle_name]["z"]=z
1644            
1645            if self.check_aminoacid_key(key=particle_name):
1646                sequence.append(particle_name) 
1647            
1648        self.define_particles(parameters=part_dict,
1649                            overwrite=overwrite,  
1650                            verbose=verbose)
1651        residue_list = self.define_AA_residues(sequence=sequence, 
1652                                               model=model)
1653        index = len(self.df)
1654        self.df.at [index,'name'] = name
1655        self.df.at [index,'pmb_type'] = 'protein'
1656        self.df.at [index,'model'] = model
1657        self.df.at [index,('sequence','')] = sequence  
1658        self.df.at [index,('residue_list','')] = residue_list    
1659        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.
  • verbose (bool, optional): Switch to activate/deactivate verbose. Defaults to True.
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):
1661    def define_residue(self, name, central_bead, side_chains):
1662        """
1663        Defines a pyMBE object of type `residue` in `pymbe.df`.
1664
1665        Args:
1666            name(`str`): Unique label that identifies the `residue`.
1667            central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`.
1668            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.
1669        """
1670        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'):
1671            return
1672        index = len(self.df)
1673        self.df.at [index, 'name'] = name
1674        self.df.at [index,'pmb_type'] = 'residue'
1675        self.df.at [index,'central_bead'] = central_bead
1676        self.df.at [index,('side_chains','')] = side_chains
1677        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 destroy_pmb_object_in_system(self, name, espresso_system):
1679    def destroy_pmb_object_in_system(self, name, espresso_system):
1680        """
1681        Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 
1682
1683        Args:
1684            name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`.
1685            espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library.
1686
1687        Note:
1688            - 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.
1689        """
1690        allowed_objects = ['particle','residue','molecule']
1691        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1692        if pmb_type not in allowed_objects:
1693            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1694        if pmb_type == 'particle':
1695            particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 
1696                                                  & (self.df['molecule_id'].isna())]
1697            particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna())
1698                                                 & (self.df['molecule_id'].isna())].particle_id.tolist()
1699            for particle_id in particle_ids_list:
1700                espresso_system.part.by_id(particle_id).remove()
1701            self.df = self.df.drop(index=particles_index)
1702        if pmb_type == 'residue':
1703            residues_id = self.df.loc[self.df['name']== name].residue_id.to_list()
1704            for residue_id in residues_id:
1705                molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0]
1706                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1707                self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index)
1708                for particle_id in particle_ids_list:
1709                    espresso_system.part.by_id(particle_id).remove()
1710                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)    
1711        if pmb_type == 'molecule':
1712            molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list()
1713            for molecule_id in molecules_id:
1714                molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']=="molecule")].name.values[0]
1715                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1716                self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index)
1717                for particle_id in particle_ids_list:
1718                    espresso_system.part.by_id(particle_id).remove()   
1719                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)             
1720        
1721        self.df.reset_index(drop=True,inplace=True)
1722
1723        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):
1725    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
1726        """
1727        Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration.
1728        To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an
1729        ion pair and the concentrations of the various species is solved numerically using a self-consistent approach.
1730        More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831).
1731        This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237).
1732
1733        Args:
1734            pH_res('float'): pH-value in the reservoir.
1735            c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
1736            activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
1737
1738        Returns:
1739            cH_res('pint.Quantity'): Concentration of H+ ions.
1740            cOH_res('pint.Quantity'): Concentration of OH- ions.
1741            cNa_res('pint.Quantity'): Concentration of Na+ ions.
1742            cCl_res('pint.Quantity'): Concentration of Cl- ions.
1743        """
1744
1745        self_consistent_run = 0
1746        cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+
1747
1748        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
1749            #Calculate and initial guess for the concentrations of various species based on ideal gas estimate
1750            cOH_res = self.Kw / cH_res 
1751            cNa_res = None
1752            cCl_res = None
1753            if cOH_res>=cH_res:
1754                #adjust the concentration of sodium if there is excess OH- in the reservoir:
1755                cNa_res = c_salt_res + (cOH_res-cH_res)
1756                cCl_res = c_salt_res
1757            else:
1758                # adjust the concentration of chloride if there is excess H+ in the reservoir
1759                cCl_res = c_salt_res + (cH_res-cOH_res)
1760                cNa_res = c_salt_res
1761                
1762            def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res):
1763                nonlocal max_number_sc_runs, self_consistent_run
1764                if self_consistent_run<max_number_sc_runs:
1765                    self_consistent_run+=1
1766                    ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1767                    cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
1768                    if cOH_res>=cH_res:
1769                        #adjust the concentration of sodium if there is excess OH- in the reservoir:
1770                        cNa_res = c_salt_res + (cOH_res-cH_res)
1771                        cCl_res = c_salt_res
1772                    else:
1773                        # adjust the concentration of chloride if there is excess H+ in the reservoir
1774                        cCl_res = c_salt_res + (cH_res-cOH_res)
1775                        cNa_res = c_salt_res
1776                    return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1777                else:
1778                    return cH_res, cOH_res, cNa_res, cCl_res
1779            return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1780
1781        cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1782        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1783        determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1784
1785        while abs(determined_pH-pH_res)>1e-6:
1786            if determined_pH > pH_res:
1787                cH_res *= 1.005
1788            else:
1789                cH_res /= 1.003
1790            cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1791            ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1792            determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1793            self_consistent_run=0
1794
1795        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):
1797    def enable_motion_of_rigid_object(self, name, espresso_system):
1798        '''
1799        Enables the motion of the rigid object `name` in the `espresso_system`.
1800
1801        Args:
1802            name(`str`): Label of the object.
1803            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1804
1805        Note:
1806            - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"].
1807        '''
1808        print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
1809        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1810        if pmb_type != 'protein':
1811            raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein')
1812        molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list()
1813        for molecule_id in molecule_ids_list:    
1814            particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list()
1815            center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system)
1816            rigid_object_center = espresso_system.part.add(pos=center_of_mass,
1817                                                           rotation=[True,True,True], 
1818                                                           type=self.propose_unused_type())
1819            
1820            rigid_object_center.mass = len(particle_ids_list)
1821            momI = 0
1822            for pid in particle_ids_list:
1823                momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
1824
1825            rigid_object_center.rinertia = np.ones(3) * momI
1826            
1827            for particle_id in particle_ids_list:
1828                pid = espresso_system.part.by_id(particle_id)
1829                pid.vs_auto_relate_to(rigid_object_center.id)
1830        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):
1832    def filter_df(self, pmb_type):
1833        """
1834        Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns.
1835        
1836        Args:
1837            pmb_type(`str`): pmb_object_type to filter in `pmb.df`.
1838
1839        Returns:
1840            pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`.
1841        """
1842        pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type]
1843        pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1)
1844        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):
1846    def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False):
1847        """
1848        Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
1849        
1850        Args:
1851            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
1852            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
1853            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'. 
1854
1855        Returns:
1856            bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` 
1857
1858        Note:
1859            - If `use_default_bond`=`True`, it returns "default" if no key is found.
1860        """
1861        bond_keys = [particle_name1 +'-'+ particle_name2, particle_name2 +'-'+ particle_name1 ]
1862        bond_defined=False
1863        for bond_key in bond_keys:
1864            if bond_key in self.df["name"].values:
1865                bond_defined=True
1866                correct_key=bond_key
1867                break
1868        if bond_defined:
1869            return correct_key
1870        elif use_default_bond:
1871            return 'default'
1872        else:
1873            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

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):
1875    def find_value_from_es_type(self, es_type, column_name):
1876        """
1877        Finds a value in `pmb.df` for a `column_name` and `es_type` pair.
1878
1879        Args:
1880            es_type(`int`): value of the espresso type
1881            column_name(`str`): name of the column in `pymbe.df`
1882
1883        Returns:
1884            Value in `pymbe.df` matching  `column_name` and `es_type`
1885        """
1886        idx = pd.IndexSlice
1887        for state in ['state_one', 'state_two']:
1888            index = np.where(self.df[(state, 'es_type')] == es_type)[0]
1889            if len(index) > 0:
1890                if column_name == 'label':
1891                    label = self.df.loc[idx[index[0]], idx[(state,column_name)]]
1892                    return label
1893                else: 
1894                    column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]]
1895                    return column_name_value
1896        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):
1898    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
1899        """
1900        Generates coordinates outside a sphere centered at `center`.
1901
1902        Args:
1903            center(`lst`): Coordinates of the center of the sphere.
1904            radius(`float`): Radius of the sphere.
1905            max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates.
1906            n_samples(`int`): Number of sample points.
1907
1908        Returns:
1909            coord_list(`lst`): Coordinates of the sample points.
1910        """
1911        if not radius > 0: 
1912            raise ValueError (f'The value of {radius} must be a positive value')
1913        if not radius < max_dist:
1914            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
1915        coord_list = []
1916        counter = 0
1917        while counter<n_samples:
1918            coord = self.generate_random_points_in_a_sphere(center=center, 
1919                                            radius=max_dist,
1920                                            n_samples=1)[0]
1921            if np.linalg.norm(coord-np.asarray(center))>=radius:
1922                coord_list.append (coord)
1923                counter += 1
1924        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):
1926    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
1927        """
1928        Uniformly samples points from a hypersphere. If on_surface is set to True, the points are
1929        uniformly sampled from the surface of the hypersphere.
1930        
1931        Args:
1932            center(`lst`): Array with the coordinates of the center of the spheres.
1933            radius(`float`): Radius of the sphere.
1934            n_samples(`int`): Number of sample points to generate inside the sphere.
1935            on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere.
1936
1937        Returns:
1938            samples(`list`): Coordinates of the sample points inside the hypersphere.
1939        """
1940        # initial values
1941        center=np.array(center)
1942        d = center.shape[0]
1943        # sample n_samples points in d dimensions from a standard normal distribution
1944        samples = self.rng.normal(size=(n_samples, d))
1945        # make the samples lie on the surface of the unit hypersphere
1946        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
1947        samples /= normalize_radii
1948        if not on_surface:
1949            # make the samples lie inside the hypersphere with the correct density
1950            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
1951            new_radii = np.power(uniform_points, 1/d)
1952            samples *= new_radii
1953        # scale the points to have the correct radius and center
1954        samples = samples * radius + center
1955        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):
1957    def generate_trial_perpendicular_vector(self,vector,magnitude):
1958        """
1959        Generates an orthogonal vector to the input `vector`.
1960
1961        Args:
1962            vector(`lst`): arbitrary vector.
1963            magnitude(`float`): magnitude of the orthogonal vector.
1964            
1965        Returns:
1966            (`lst`): Orthogonal vector with the same magnitude as the input vector.
1967        """ 
1968        np_vec = np.array(vector) 
1969        np_vec /= np.linalg.norm(np_vec) 
1970        if np.all(np_vec == 0):
1971            raise ValueError('Zero vector')
1972        # Generate a random vector 
1973        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
1974                                                                center=[0,0,0],
1975                                                                n_samples=1, 
1976                                                                on_surface=True)[0]
1977        # Project the random vector onto the input vector and subtract the projection
1978        projection = np.dot(random_vector, np_vec) * np_vec
1979        perpendicular_vector = random_vector - projection
1980        # Normalize the perpendicular vector to have the same magnitude as the input vector
1981        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
1982        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):
1984    def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
1985        """
1986        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length.
1987        If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead.
1988        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.
1989
1990        Args:
1991            particle_name1(str): label of the type of the first particle type of the bonded particles.
1992            particle_name2(str): label of the type of the second particle type of the bonded particles.
1993            hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
1994            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. 
1995
1996        Returns:
1997            l0(`pint.Quantity`): bond length
1998        
1999        Note:
2000            - 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`.
2001            - If `hard_check`=`True` stops the code when no bond is found.
2002        """
2003        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2004                                    particle_name2=particle_name2, 
2005                                    use_default_bond=use_default_bond)
2006        if bond_key:
2007            return self.df[self.df['name']==bond_key].l0.values[0]
2008        else:
2009            print("Bond not defined between particles ", particle_name1, " and ", particle_name2)    
2010            if hard_check:
2011                sys.exit(1)
2012            else:
2013                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):
2015    def get_charge_number_map(self):
2016        '''
2017        Gets the charge number of each `espresso_type` in `pymbe.df`.
2018        
2019        Returns:
2020            charge_number_map(`dict`): {espresso_type: z}.
2021        '''
2022        if self.df.state_one['es_type'].isnull().values.any():         
2023            df_state_one = self.df.state_one.dropna()     
2024            df_state_two = self.df.state_two.dropna()  
2025        else:    
2026            df_state_one = self.df.state_one
2027            if self.df.state_two['es_type'].isnull().values.any():
2028                df_state_two = self.df.state_two.dropna()   
2029            else:
2030                df_state_two = self.df.state_two
2031        state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values)
2032        state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values)
2033        charge_number_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2034        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'):
2036    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2037        """
2038        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2039        `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule.
2040
2041        Args:
2042            particle_name1 (str): label of the type of the first particle type
2043            particle_name2 (str): label of the type of the second particle type
2044            combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2045
2046        Returns:
2047            {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2048
2049        Note:
2050            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
2051            - 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.
2052        """
2053        supported_combining_rules=["Lorentz-Berthelot"]
2054        lj_parameters_keys=["sigma","epsilon","offset","cutoff"]
2055        if combining_rule not in supported_combining_rules:
2056            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2057        lj_parameters={}
2058        for key in lj_parameters_keys:
2059            lj_parameters[key]=[]
2060        # Search the LJ parameters of the type pair
2061        for name in [particle_name1,particle_name2]:
2062            for key in lj_parameters_keys:
2063                lj_parameters[key].append(self.df[self.df.name == name][key].values[0])
2064        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2065        if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]):
2066            return {}
2067        # Apply combining rule
2068        if combining_rule == 'Lorentz-Berthelot':
2069            lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2
2070            lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2
2071            lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2
2072            lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1])
2073        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):
2075    def get_metal_ions_charge_number_map(self):
2076        """
2077        Gets a map with the charge numbers of all the metal ions supported.
2078
2079        Returns:
2080            metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number}
2081
2082        """
2083        metal_charge_number_map = {"Ca": 2}
2084        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):
2086    def get_particle_id_map(self, object_name):
2087        '''
2088        Gets all the ids associated with the object with name `object_name` in `pmb.df`
2089
2090        Args:
2091            object_name(`str`): name of the object
2092        
2093        Returns:
2094            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]}, }
2095        '''
2096        object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0]
2097        valid_types = ["particle", "molecule", "residue", "protein"]
2098        if object_type not in valid_types:
2099            raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}")
2100        id_list = []
2101        mol_map = {}
2102        res_map = {}
2103        def do_res_map(res_ids):
2104            for res_id in res_ids:
2105                res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2106                res_map[res_id]=res_list
2107            return res_map
2108        if object_type in ['molecule', 'protein']:
2109            mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist()
2110            for mol_id in mol_ids:
2111                res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist())
2112                res_map=do_res_map(res_ids=res_ids)    
2113                mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2114                id_list+=mol_list
2115                mol_map[mol_id]=mol_list
2116        elif object_type == 'residue':     
2117            res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist()
2118            res_map=do_res_map(res_ids=res_ids)
2119            id_list=[]
2120            for res_id_list in res_map.values():
2121                id_list+=res_id_list
2122        elif object_type == 'particle':
2123            id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist()
2124        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):
2126    def get_pka_set(self):
2127        '''
2128        Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df`
2129        
2130        Returns:
2131            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
2132        '''
2133        titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna()
2134        pka_set = {}
2135        for index in titratables_AA_df.name.keys():
2136            name = titratables_AA_df.name[index]
2137            pka_value = titratables_AA_df.pka[index]
2138            acidity = titratables_AA_df.acidity[index]   
2139            pka_set[name] = {'pka_value':pka_value,'acidity':acidity}
2140        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):
2142    def get_radius_map(self, dimensionless=True):
2143        '''
2144        Gets the effective radius of each `espresso type` in `pmb.df`. 
2145
2146        Args:
2147            dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False.
2148        
2149        Returns:
2150            radius_map(`dict`): {espresso_type: radius}.
2151
2152        Note:
2153            The radius corresponds to (sigma+offset)/2
2154        '''
2155        df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates()
2156        df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates()
2157        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)
2158        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)
2159        radius_map  = pd.concat([state_one,state_two],axis=0).to_dict()  
2160        if dimensionless:
2161            for key in radius_map:
2162                radius_map[key] = radius_map[key].magnitude
2163        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):
2165    def get_reduced_units(self):
2166        """
2167        Returns the  current set of reduced units defined in pyMBE.
2168
2169        Returns:
2170            reduced_units_text(`str`): text with information about the current set of reduced units.
2171
2172        """
2173        unit_length=self.units.Quantity(1,'reduced_length')
2174        unit_energy=self.units.Quantity(1,'reduced_energy')
2175        unit_charge=self.units.Quantity(1,'reduced_charge')
2176        reduced_units_text = "\n".join(["Current set of reduced units:",
2177                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2178                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2179                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2180                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2181                                        ])   
2182        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):
2184    def get_resource(self, path):
2185        '''
2186        Locate a file resource of the pyMBE package.
2187
2188        Args:
2189            path(`str`): Relative path to the resource
2190
2191        Returns:
2192            path(`str`): Absolute path to the resource
2193
2194        '''
2195        import os
2196        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):
2198    def get_type_map(self):
2199        """
2200        Gets all different espresso types assigned to particles  in `pmb.df`.
2201        
2202        Returns:
2203            type_map(`dict`): {"name": espresso_type}.
2204        """
2205        if self.df.state_one['es_type'].isnull().values.any():         
2206            df_state_one = self.df.state_one.dropna(how='all')     
2207            df_state_two = self.df.state_two.dropna(how='all')  
2208        else:    
2209            df_state_one = self.df.state_one
2210            if self.df.state_two['es_type'].isnull().values.any():
2211                df_state_two = self.df.state_two.dropna(how='all')   
2212            else:
2213                df_state_two = self.df.state_two
2214        state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label)
2215        state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label)
2216        type_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2217        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, verbose=False, overwrite=False):
2219    def load_interaction_parameters(self, filename, verbose=False, overwrite=False):
2220        """
2221        Loads the interaction parameters stored in `filename` into `pmb.df`
2222        
2223        Args:
2224            filename(`str`): name of the file to be read
2225            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to False.
2226            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2227        """
2228        without_units = ['z','es_type','acidity']
2229        with_units = ['sigma','epsilon','offset','cutoff']
2230
2231        with open(filename, 'r') as f:
2232            interaction_data = json.load(f)
2233            interaction_parameter_set = interaction_data["data"]
2234
2235        for key in interaction_parameter_set:
2236            param_dict=interaction_parameter_set[key]
2237            object_type=param_dict['object_type']
2238            if object_type == 'particle':
2239                not_required_attributes={}    
2240                for not_required_key in without_units+with_units:
2241                    if not_required_key in param_dict.keys():
2242                        if not_required_key in with_units:
2243                            not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key))
2244                        elif not_required_key in without_units:
2245                            not_required_attributes[not_required_key]=param_dict.pop(not_required_key)
2246                    else:
2247                        not_required_attributes[not_required_key]=None
2248                self.define_particle(name=param_dict.pop('name'),
2249                                z=not_required_attributes.pop('z'),
2250                                sigma=not_required_attributes.pop('sigma'),
2251                                offset=not_required_attributes.pop('offset'),
2252                                cutoff=not_required_attributes.pop('cutoff'),
2253                                acidity=not_required_attributes.pop('acidity'),
2254                                epsilon=not_required_attributes.pop('epsilon'),
2255                                verbose=verbose,
2256                                overwrite=overwrite)
2257            elif object_type == 'residue':
2258                self.define_residue (name = param_dict.pop('name'),
2259                                    central_bead = param_dict.pop('central_bead'),
2260                                    side_chains = param_dict.pop('side_chains'))
2261            elif object_type == 'molecule':
2262                self.define_molecule(name=param_dict.pop('name'),
2263                                    residue_list=param_dict.pop('residue_list'))
2264            elif object_type == 'peptide':
2265                self.define_peptide(name=param_dict.pop('name'),
2266                                    sequence=param_dict.pop('sequence'),
2267                                    model=param_dict.pop('model'))
2268            elif object_type == 'bond':
2269                particle_pairs = param_dict.pop('particle_pairs')
2270                bond_parameters = param_dict.pop('bond_parameters')
2271                bond_type = param_dict.pop('bond_type')
2272                if bond_type == 'harmonic':
2273                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2274                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2275                    bond = {'r_0'    : r_0,
2276                            'k'      : k,
2277                            }
2278
2279                elif bond_type == 'FENE':
2280                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2281                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2282                    d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max'))
2283                    bond =  {'r_0'    : r_0,
2284                             'k'      : k,
2285                             'd_r_max': d_r_max,
2286                             }
2287                else:
2288                    raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds")
2289                self.define_bond(bond_type=bond_type, 
2290                                bond_parameters=bond, 
2291                                particle_pairs=particle_pairs)
2292            else:
2293                raise ValueError(object_type+' is not a known pmb object type')
2294            
2295        return

Loads the interaction parameters stored in filename into pmb.df

Arguments:
  • filename(str): name of the file to be read
  • verbose (bool, optional): Switch to activate/deactivate verbose. Defaults to False.
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
def load_pka_set(self, filename, verbose=False, overwrite=True):
2297    def load_pka_set(self, filename, verbose=False, overwrite=True):
2298        """
2299        Loads the pka_set stored in `filename` into `pmb.df`.
2300        
2301        Args:
2302            filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}.
2303            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to False.
2304            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 
2305        """
2306        with open(filename, 'r') as f:
2307            pka_data = json.load(f)
2308            pka_set = pka_data["data"]
2309
2310        self.check_pka_set(pka_set=pka_set)
2311
2312        for key in pka_set:
2313            acidity = pka_set[key]['acidity']
2314            pka_value = pka_set[key]['pka_value']
2315            self.set_particle_acidity(name=key, 
2316                                      acidity=acidity, 
2317                                      pka=pka_value, 
2318                                      verbose=verbose, 
2319                                      overwrite=overwrite)
2320        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}}.
  • verbose (bool, optional): Switch to activate/deactivate verbose. Defaults to False.
  • overwrite(bool, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True.
def parse_sequence_from_file(self, sequence):
2322    def parse_sequence_from_file(self,sequence):
2323        """
2324        Parses the given sequence such that it can be used in pyMBE. This function has to be used if the df was read from a file.
2325
2326        Args:
2327            sequence(`str`): sequence to be parsed
2328
2329        Returns:
2330            sequence(`lst`): parsed sequence
2331        """
2332        sequence = sequence.replace(' ', '')
2333        sequence = sequence.replace("'", '')
2334        sequence = sequence.split(",")[1:-1]
2335        return sequence

Parses the given sequence such that it can be used in pyMBE. This function has to be used if the df was read from a file.

Arguments:
  • sequence(str): sequence to be parsed
Returns:

sequence(lst): parsed sequence

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

Writes the pyMBE dataframe into a csv file

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

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

def default(self, obj):
50        def default(self, obj):
51            if isinstance(obj, np.ndarray):
52                return obj.tolist()
53            if isinstance(obj, np.generic):
54                return obj.item()
55            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