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

The library for the Molecular Builder for ESPResSo (pyMBE)

Attributes:
  • N_A(pint.Quantity): Avogadro number.
  • Kb(pint.Quantity): Boltzmann constant.
  • e(pint.Quantity): Elementary charge.
  • df(Pandas.Dataframe): Dataframe used to bookkeep all the information stored in pyMBE. Typically refered as pmb.df.
  • kT(pint.Quantity): Thermal energy.
  • Kw(pint.Quantity): Ionic product of water. Used in the setup of the G-RxMC method.
pymbe_library(seed, temperature=None, unit_length=None, unit_charge=None, Kw=None)
56    def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None):
57        """
58        Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 
59        and sets up  the `pmb.df` for bookkeeping.
60
61        Args:
62            temperature(`pint.Quantity`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None.
63            unit_length(`pint.Quantity`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None.
64            unit_charge (`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
65            Kw (`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
66        
67        Note:
68            - If no `temperature` is given, a value of 298.15 K is assumed by default.
69            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
70            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
71            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
72        """
73        # Seed and RNG
74        self.seed=seed
75        self.rng = np.random.default_rng(seed)
76        self.units=pint.UnitRegistry()
77        self.N_A=scipy.constants.N_A / self.units.mol
78        self.kB=scipy.constants.k * self.units.J / self.units.K
79        self.e=scipy.constants.e * self.units.C
80        self.set_reduced_units(unit_length=unit_length, 
81                               unit_charge=unit_charge,
82                               temperature=temperature, 
83                               Kw=Kw, 
84                               verbose=False)
85        self.setup_df()
86        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):
 88    def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False):
 89        """
 90        Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles.
 91
 92        Args:
 93            particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles
 94            particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles
 95            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False.
 96
 97        Returns:
 98            index(`int`): Row index where the bond information has been added in pmb.df.
 99        """
100        particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0]
101        particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0]
102        
103        bond_key = self.find_bond_key(particle_name1=particle_name1,
104                                    particle_name2=particle_name2, 
105                                    use_default_bond=use_default_bond)
106        if not bond_key:
107            return None
108        self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1)
109        indexs = np.where(self.df['name']==bond_key)
110        index_list = list (indexs[0])
111        used_bond_df = self.df.loc[self.df['particle_id2'].notnull()]
112        #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict
113        used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 )
114        used_bond_index = used_bond_df.index.to_list()
115        if not index_list:
116            return None
117        for index in index_list:
118            if index not in used_bond_index:
119                self.clean_df_row(index=int(index))
120                self.df.at[index,'particle_id'] = particle_id1
121                self.df.at[index,'particle_id2'] = particle_id2
122                break
123        return index

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):
125    def add_bonds_to_espresso(self, espresso_system) :
126        """
127        Adds all bonds defined in `pmb.df` to `espresso_system`.
128
129        Args:
130            espresso_system(`espressomd.system.System`): system object of espressomd library
131        """
132
133        if 'bond' in self.df.values:
134            bond_df = self.df.loc[self.df ['pmb_type'] == 'bond']
135            bond_list = bond_df.bond_object.values.tolist()
136            for bond in bond_list:
137                espresso_system.bonded_inter.add(bond)
138        else:
139            print ('WARNING: There are no bonds defined in pymbe.df')
140        
141        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):
143    def add_value_to_df(self,index,key,new_value, verbose=True, non_standard_value=False, overwrite=False):
144        """
145        Adds a value to a cell in the `pmb.df` DataFrame.
146
147        Args:
148            index(`int`): index of the row to add the value to.
149            key(`str`): the column label to add the value to.
150            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
151            non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False.
152            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
153        """
154
155        token = "#protected:"
156
157        def protect(obj):
158            if non_standard_value:
159                return token + json.dumps(obj, cls=self.NumpyEncoder)
160            return obj
161
162        def deprotect(obj):
163            if non_standard_value and isinstance(obj, str) and obj.startswith(token):
164                return json.loads(obj.removeprefix(token))
165            return obj
166
167        # Make sure index is a scalar integer value
168        index = int (index)
169        assert isinstance(index, int), '`index` should be a scalar integer value.'
170        idx = pd.IndexSlice
171        if self.check_if_df_cell_has_a_value(index=index,key=key):
172            old_value= self.df.loc[index,idx[key]]
173            if protect(old_value) != protect(new_value):
174                name=self.df.loc[index,('name','')]
175                pmb_type=self.df.loc[index,('pmb_type','')]
176                if verbose:
177                    print(f"WARNING: you are attempting to redefine the properties of {name} of pmb_type {pmb_type}")    
178                if overwrite and verbose:
179                    print(f'WARNING: overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}')
180                if not overwrite:
181                    if verbose:
182                        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 ")
183                    return
184        self.df.loc[index,idx[key]] = protect(new_value)
185        if non_standard_value:
186            self.df[key] = self.df[key].apply(deprotect)
187        return

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

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

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

Creates either a harmonic or a FENE bond in ESPResSo

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

Currently, only HARMONIC and FENE bonds are supported.

For a HARMONIC bond the dictionary must contain:

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

For a FENE bond the dictionary must additionally contain:

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

bond_object (obj): an ESPResSo bond object

def create_counterions( self, object_name, cation_name, anion_name, espresso_system, verbose=True):
822    def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True):
823        """
824        Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`.
825        
826        Args:
827            object_name(`str`): `name` of a pymbe object.
828            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
829            cation_name(`str`): `name` of a particle with a positive charge.
830            anion_name(`str`): `name` of a particle with a negative charge.
831            verbose(`bool`): switch to activate/deactivate verbose. Defaults to True.
832
833        Returns: 
834            counterion_number(`dict`): {"name": number}
835        """
836        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0]
837        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0]
838        object_ids = self.get_particle_id_map(object_name=object_name)["all"]
839        counterion_number={}
840        object_charge={}
841        for name in ['positive', 'negative']:
842            object_charge[name]=0
843        for id in object_ids:
844            if espresso_system.part.by_id(id).q > 0:
845                object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q ))
846            elif espresso_system.part.by_id(id).q < 0:
847                object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q ))
848        if object_charge['positive'] % abs(anion_charge) == 0:
849            counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge))
850        else:
851            raise ValueError('The number of positive charges in the pmb_object must be divisible by the  charge of the anion')
852        if object_charge['negative'] % abs(cation_charge) == 0:
853            counterion_number[cation_name]=int(object_charge['negative']/cation_charge)
854        else:
855            raise ValueError('The number of negative charges in the pmb_object must be divisible by the  charge of the cation')
856        if counterion_number[cation_name] > 0: 
857            self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name])
858        else:
859            counterion_number[cation_name]=0
860        if counterion_number[anion_name] > 0:
861            self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name])
862        else:
863            counterion_number[anion_name] = 0
864        if verbose:
865            print('The following counter-ions have been created: ')
866            for name in counterion_number.keys():
867                print(f'Ion type: {name} created number: {counterion_number[name]}')
868        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):
870    def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False):
871        """
872        Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
873
874        Args:
875            name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df`
876            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
877            number_of_molecules(`int`): Number of molecules of type `name` to be created.
878            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.
879            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. 
880            use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df`
881
882        Returns:
883            molecules_info(`dict`):  {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 
884        """
885        if list_of_first_residue_positions is not None:
886            for item in list_of_first_residue_positions:
887                if not isinstance(item, list):
888                    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.")
889                elif len(item) != 3:
890                    raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.")
891
892            if len(list_of_first_residue_positions) != number_of_molecules:
893                raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}")
894        if number_of_molecules <= 0:
895            return 0
896        self.check_if_name_is_defined_in_df(name=name,
897                                            pmb_type_to_be_defined='molecule')
898            
899        first_residue = True
900        molecules_info = {}
901        residue_list = self.df[self.df['name']==name].residue_list.values [0]
902        self.copy_df_entry(name=name,
903                        column_name='molecule_id',
904                        number_of_copies=number_of_molecules)
905
906        molecules_index = np.where(self.df['name']==name)
907        molecule_index_list =list(molecules_index[0])[-number_of_molecules:]
908        used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist()
909        pos_index = 0 
910        for molecule_index in molecule_index_list:        
911            molecule_id = self.assign_molecule_id(name=name,
912                                                pmb_type='molecule',
913                                                used_molecules_id=used_molecules_id,
914                                                molecule_index=molecule_index)
915            molecules_info[molecule_id] = {}
916            for residue in residue_list:
917                if first_residue:
918                    if list_of_first_residue_positions is None:
919                        residue_position = None
920                    else:
921                        for item in list_of_first_residue_positions:
922                            residue_position = [np.array(list_of_first_residue_positions[pos_index])]
923                    # Generate an arbitrary random unit vector
924                    if backbone_vector is None:
925                        backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0],
926                                                                radius=1, 
927                                                                n_samples=1,
928                                                                on_surface=True)[0]
929                    residues_info = self.create_residue(name=residue,
930                                                        espresso_system=espresso_system, 
931                                                        central_bead_position=residue_position,  
932                                                        use_default_bond= use_default_bond, 
933                                                        backbone_vector=backbone_vector)
934                    residue_id = next(iter(residues_info))
935                    # Add the correct molecule_id to all particles in the residue
936                    for index in self.df[self.df['residue_id']==residue_id].index:
937                        self.add_value_to_df(key=('molecule_id',''),
938                                            index=int (index),
939                                            new_value=molecule_id,
940                                            overwrite=True,
941                                            verbose=False)
942                    central_bead_id = residues_info[residue_id]['central_bead_id']
943                    previous_residue = residue
944                    residue_position = espresso_system.part.by_id(central_bead_id).pos
945                    previous_residue_id = central_bead_id
946                    first_residue = False          
947                else:                    
948                    previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0]
949                    new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0]       
950                    bond = self.search_bond(particle_name1=previous_central_bead_name, 
951                                            particle_name2=new_central_bead_name, 
952                                            hard_check=True, 
953                                            use_default_bond=use_default_bond)                
954                    l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 
955                                            particle_name2=new_central_bead_name, 
956                                            hard_check=True, 
957                                            use_default_bond=use_default_bond)                
958                    residue_position = residue_position+backbone_vector*l0
959                    residues_info = self.create_residue(name=residue, 
960                                                        espresso_system=espresso_system, 
961                                                        central_bead_position=[residue_position],
962                                                        use_default_bond= use_default_bond, 
963                                                        backbone_vector=backbone_vector)
964                    residue_id = next(iter(residues_info))      
965                    for index in self.df[self.df['residue_id']==residue_id].index:
966                        self.add_value_to_df(key=('molecule_id',''),
967                                            index=int (index),
968                                            new_value=molecule_id,
969                                            verbose=False,
970                                            overwrite=True)            
971                    central_bead_id = residues_info[residue_id]['central_bead_id']
972                    espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id))
973                    bond_index = self.add_bond_in_df(particle_id1=central_bead_id,
974                                        particle_id2=previous_residue_id,
975                                        use_default_bond=use_default_bond) 
976                    self.add_value_to_df(key=('molecule_id',''),
977                                            index=int (bond_index),
978                                            new_value=molecule_id,
979                                            verbose=False,
980                                            overwrite=True)           
981                    previous_residue_id = central_bead_id
982                    previous_residue = residue                    
983                molecules_info[molecule_id][residue_id] = residues_info[residue_id]
984            first_residue = True
985            pos_index+=1
986        
987        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):
 989    def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False):
 990        """
 991        Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`.
 992        
 993        Args:
 994            name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 
 995            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
 996            number_of_particles(`int`): Number of particles to be created.
 997            position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None.
 998            fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False.
 999        Returns:
1000            created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`.
1001        """       
1002        if number_of_particles <=0:
1003            return []
1004        self.check_if_name_is_defined_in_df(name=name,
1005                                       pmb_type_to_be_defined='particle')
1006        # Copy the data of the particle `number_of_particles` times in the `df`
1007        self.copy_df_entry(name=name,
1008                          column_name='particle_id',
1009                          number_of_copies=number_of_particles)
1010        # Get information from the particle type `name` from the df     
1011        z = self.df.loc[self.df['name']==name].state_one.z.values[0]
1012        z = 0. if z is None else z
1013        es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
1014        # Get a list of the index in `df` corresponding to the new particles to be created
1015        index = np.where(self.df['name']==name)
1016        index_list =list(index[0])[-number_of_particles:]
1017        # Create the new particles into  `espresso_system`
1018        created_pid_list=[]
1019        for index in range (number_of_particles):
1020            df_index=int (index_list[index])
1021            self.clean_df_row(index=df_index)
1022            if position is None:
1023                particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l)
1024            else:
1025                particle_position = position[index]
1026            if len(espresso_system.part.all()) == 0:
1027                bead_id = 0
1028            else:
1029                bead_id = max (espresso_system.part.all().id) + 1
1030            created_pid_list.append(bead_id)
1031            kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z)
1032            if fix:
1033                kwargs["fix"] = 3 * [fix]
1034            espresso_system.part.add(**kwargs)
1035            self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id, verbose=False)                  
1036        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):
1038    def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None):
1039        """
1040        Creates all `particle`s associated to `pmb object` into  `espresso` a number of times equal to `number_of_objects`.
1041        
1042        Args:
1043            name(`str`): Unique label of the `pmb object` to be created. 
1044            number_of_objects(`int`): Number of `pmb object`s to be created.
1045            espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library.
1046            position(`list`): Coordinates where the particles should be created.
1047            use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`.
1048            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. 
1049
1050        Note:
1051            - 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. 
1052        """
1053        allowed_objects=['particle','residue','molecule']
1054        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1055        if pmb_type not in allowed_objects:
1056            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1057        if pmb_type == 'particle':
1058            self.create_particle(name=name, 
1059                                number_of_particles=number_of_objects, 
1060                                espresso_system=espresso_system, 
1061                                position=position)
1062        elif pmb_type == 'residue':
1063            self.create_residue(name=name,  
1064                                espresso_system=espresso_system, 
1065                                central_bead_position=position,
1066                                use_default_bond=use_default_bond,
1067                                backbone_vector=backbone_vector)
1068        elif pmb_type == 'molecule':
1069            self.create_molecule(name=name, 
1070                                number_of_molecules=number_of_objects, 
1071                                espresso_system=espresso_system, 
1072                                use_default_bond=use_default_bond, 
1073                                list_of_first_residue_positions=position,
1074                                backbone_vector=backbone_vector)
1075        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):
1077    def create_protein(self, name, number_of_proteins, espresso_system, topology_dict):
1078        """
1079        Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions`
1080
1081        Args:
1082            name(`str`): Label of the protein to be created. 
1083            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1084            number_of_proteins(`int`): Number of proteins to be created.
1085            positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}}
1086        """
1087
1088        if number_of_proteins <=0:
1089            return
1090        self.check_if_name_is_defined_in_df(name=name,
1091                                        pmb_type_to_be_defined='protein')
1092        self.copy_df_entry(name=name,
1093                            column_name='molecule_id',
1094                            number_of_copies=number_of_proteins)
1095        protein_index = np.where(self.df['name']==name)
1096        protein_index_list =list(protein_index[0])[-number_of_proteins:]
1097        used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist()
1098        
1099        box_half=espresso_system.box_l[0]/2.0
1100
1101        for molecule_index in protein_index_list:     
1102
1103            molecule_id = self.assign_molecule_id (name=name,pmb_type='protein',used_molecules_id=used_molecules_id,molecule_index=molecule_index)
1104
1105            protein_center = self.generate_coordinates_outside_sphere(radius = 1, 
1106                                                                        max_dist=box_half, 
1107                                                                        n_samples=1, 
1108                                                                        center=[box_half]*3)[0]
1109   
1110            for residue in topology_dict.keys():
1111
1112                residue_name = re.split(r'\d+', residue)[0]
1113                residue_number = re.split(r'(\d+)', residue)[1]
1114                residue_position = topology_dict[residue]['initial_pos']
1115                position = residue_position + protein_center
1116
1117                particle_id = self.create_particle(name=residue_name,
1118                                                            espresso_system=espresso_system,
1119                                                            number_of_particles=1,
1120                                                            position=[position], 
1121                                                            fix = True)
1122                
1123                index = self.df[self.df['particle_id']==particle_id[0]].index.values[0]
1124                
1125                self.add_value_to_df(key=('residue_id',''),
1126                                            index=int (index),
1127                                            new_value=residue_number,
1128                                            overwrite=True,
1129                                            verbose=False)
1130
1131                self.add_value_to_df(key=('molecule_id',''),
1132                                        index=int (index),
1133                                        new_value=molecule_id,
1134                                        overwrite=True,
1135                                        verbose=False)
1136
1137        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):
1139    def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None):
1140        """
1141        Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`.
1142
1143        Args:
1144            name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df`
1145            espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library.
1146            central_bead_position(`list` of `float`): Position of the central bead.
1147            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`
1148            backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`.
1149
1150        Returns:
1151            residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1152        """
1153        self.check_if_name_is_defined_in_df(name=name,
1154                                            pmb_type_to_be_defined='residue')
1155        # Copy the data of a residue in the `df
1156        self.copy_df_entry(name=name,
1157                            column_name='residue_id',
1158                            number_of_copies=1)
1159        residues_index = np.where(self.df['name']==name)
1160        residue_index_list =list(residues_index[0])[-1:]
1161        for residue_index in residue_index_list:  
1162            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1163            for side_chain_element in side_chain_list:
1164                if side_chain_element not in self.df.values:              
1165                    raise ValueError (side_chain_element +'is not defined')
1166        # Internal bookkepping of the residue info (important for side-chain residues)
1167        # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}}
1168        residues_info={}
1169        for residue_index in residue_index_list:     
1170            self.clean_df_row(index=int(residue_index))
1171            # Assign a residue_id
1172            if self.df['residue_id'].isnull().all():
1173                residue_id=0
1174            else:
1175                residue_id = self.df['residue_id'].max() + 1
1176            self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id, verbose=False)
1177            # create the principal bead   
1178            central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0]            
1179            central_bead_id = self.create_particle(name=central_bead_name,
1180                                                                espresso_system=espresso_system,
1181                                                                position=central_bead_position,
1182                                                                number_of_particles = 1)[0]
1183            central_bead_position=espresso_system.part.by_id(central_bead_id).pos
1184            #assigns same residue_id to the central_bead particle created.
1185            index = self.df[self.df['particle_id']==central_bead_id].index.values[0]
1186            self.df.at [index,'residue_id'] = residue_id
1187            # Internal bookkeeping of the central bead id
1188            residues_info[residue_id]={}
1189            residues_info[residue_id]['central_bead_id']=central_bead_id
1190            # create the lateral beads  
1191            side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0]
1192            side_chain_beads_ids = []
1193            for side_chain_element in side_chain_list:
1194                
1195                pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 
1196                if pmb_type == 'particle':
1197                    bond = self.search_bond(particle_name1=central_bead_name, 
1198                                            particle_name2=side_chain_element, 
1199                                            hard_check=True, 
1200                                            use_default_bond=use_default_bond)
1201                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1202                                              particle_name2=side_chain_element, 
1203                                              hard_check=True, 
1204                                              use_default_bond=use_default_bond)
1205
1206                    if backbone_vector is None:
1207                        bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1208                                                                 radius=l0, 
1209                                                                 n_samples=1,
1210                                                                 on_surface=True)[0]
1211                    else:
1212                        bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1213                                                                                                    magnitude=l0)
1214                    
1215                    side_bead_id = self.create_particle(name=side_chain_element, 
1216                                                                    espresso_system=espresso_system,
1217                                                                    position=[bead_position], 
1218                                                                    number_of_particles=1)[0]
1219                    index = self.df[self.df['particle_id']==side_bead_id].index.values[0]
1220                    self.add_value_to_df(key=('residue_id',''),
1221                                        index=int (index),
1222                                        new_value=residue_id, 
1223                                        verbose=False,
1224                                        overwrite=True)
1225                    side_chain_beads_ids.append(side_bead_id)
1226                    espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id))
1227                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1228                                        particle_id2=side_bead_id,
1229                                        use_default_bond=use_default_bond)
1230                    self.add_value_to_df(key=('residue_id',''),
1231                                        index=int (index),
1232                                        new_value=residue_id, 
1233                                        verbose=False,
1234                                        overwrite=True)
1235
1236                elif pmb_type == 'residue':
1237                    central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0]
1238                    bond = self.search_bond(particle_name1=central_bead_name, 
1239                                            particle_name2=central_bead_side_chain, 
1240                                            hard_check=True, 
1241                                            use_default_bond=use_default_bond)
1242                    l0 = self.get_bond_length(particle_name1=central_bead_name, 
1243                                              particle_name2=central_bead_side_chain, 
1244                                              hard_check=True, 
1245                                              use_default_bond=use_default_bond)
1246                    if backbone_vector is None:
1247                        residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 
1248                                                                    radius=l0, 
1249                                                                    n_samples=1,
1250                                                                    on_surface=True)[0]
1251                    else:
1252                        residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector,
1253                                                                                                        magnitude=l0)
1254                    lateral_residue_info = self.create_residue(name=side_chain_element, 
1255                                                                espresso_system=espresso_system,
1256                                                                central_bead_position=[residue_position],
1257                                                                use_default_bond=use_default_bond)
1258                    lateral_residue_dict=list(lateral_residue_info.values())[0]
1259                    central_bead_side_chain_id=lateral_residue_dict['central_bead_id']
1260                    lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids']
1261                    residue_id_side_chain=list(lateral_residue_info.keys())[0]
1262                    # Change the residue_id of the residue in the side chain to the one of the bigger residue
1263                    index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0]
1264                    self.add_value_to_df(key=('residue_id',''),
1265                                        index=int(index),
1266                                        new_value=residue_id, 
1267                                        verbose=False,
1268                                        overwrite=True)
1269                    # Change the residue_id of the particles in the residue in the side chain
1270                    side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids
1271                    for particle_id in side_chain_beads_ids:
1272                        index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0]
1273                        self.add_value_to_df(key=('residue_id',''),
1274                                            index=int (index),
1275                                            new_value=residue_id, 
1276                                            verbose=False,
1277                                            overwrite=True)
1278                    espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id))
1279                    index = self.add_bond_in_df(particle_id1=central_bead_id,
1280                                        particle_id2=central_bead_side_chain_id,
1281                                        use_default_bond=use_default_bond)
1282                    self.add_value_to_df(key=('residue_id',''),
1283                                        index=int (index),
1284                                        new_value=residue_id, 
1285                                        verbose=False,
1286                                        overwrite=True)
1287                    # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue
1288                    for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index:        
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            # Internal bookkeeping of the side chain beads ids
1295            residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids
1296        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):
1298    def create_variable_with_units(self, variable):
1299        """
1300        Returns a pint object with the value and units defined in `variable`.
1301
1302        Args:
1303            variable(`dict` or `str`): {'value': value, 'units': units}
1304        Returns:
1305            variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry.
1306        """        
1307        
1308        if isinstance(variable, dict):
1309
1310            value=variable.pop('value')
1311            units=variable.pop('units')
1312
1313        elif isinstance(variable, str):
1314
1315            value = float(re.split(r'\s+', variable)[0])
1316            units = re.split(r'\s+', variable)[1]
1317        
1318        variable_with_units=value*self.units(units)
1319
1320        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):
1322    def define_AA_residues(self, sequence, model):
1323        """
1324        Defines in `pmb.df` all the different residues in `sequence`.
1325
1326        Args:
1327            sequence(`lst`):  Sequence of the peptide or protein.
1328            model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1329
1330        Returns:
1331            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.             
1332        """
1333        residue_list = []
1334        for residue_name in sequence:
1335            if model == '1beadAA':
1336                central_bead = residue_name
1337                side_chains = []
1338            elif model == '2beadAA':
1339                if residue_name in ['c','n', 'G']: 
1340                    central_bead = residue_name
1341                    side_chains = []
1342                else:
1343                    central_bead = 'CA'              
1344                    side_chains = [residue_name]
1345            if residue_name not in residue_list:   
1346                self.define_residue(name = 'AA-'+residue_name, 
1347                                    central_bead = central_bead,
1348                                    side_chains = side_chains)              
1349            residue_list.append('AA-'+residue_name)
1350        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):
1352    def define_bond(self, bond_type, bond_parameters, particle_pairs):
1353        
1354        '''
1355        Defines a pmb object of type `bond` in `pymbe.df`.
1356
1357        Args:
1358            bond_type(`str`): label to identify the potential to model the bond.
1359            bond_parameters(`dict`): parameters of the potential of the bond.
1360            particle_pairs(`lst`): list of the `names` of the `particles` to be bonded.
1361
1362        Note:
1363            Currently, only HARMONIC and FENE bonds are supported.
1364
1365            For a HARMONIC bond the dictionary must contain the following parameters:
1366
1367                - k (`obj`)      : Magnitude of the bond. It should have units of energy/length**2 
1368                using the `pmb.units` UnitRegistry.
1369                - r_0 (`obj`)    : Equilibrium bond length. It should have units of length using 
1370                the `pmb.units` UnitRegistry.
1371           
1372            For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:
1373                
1374                - d_r_max (`obj`): Maximal stretching length for FENE. It should have 
1375                units of length using the `pmb.units` UnitRegistry. Default 'None'.
1376        '''
1377
1378        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1379
1380        for particle_name1, particle_name2 in particle_pairs:
1381
1382            lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1,
1383                                                 particle_name2 = particle_name2,
1384                                                 combining_rule = 'Lorentz-Berthelot')
1385
1386            l0 = self.calculate_initial_bond_length(bond_object = bond_object,
1387                                                    bond_type   = bond_type,
1388                                                    epsilon     = lj_parameters["epsilon"],
1389                                                    sigma       = lj_parameters["sigma"],
1390                                                    cutoff      = lj_parameters["cutoff"],
1391                                                    offset      = lj_parameters["offset"],)
1392            index = len(self.df)
1393            for label in [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']:
1394                self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond")
1395            self.df.at [index,'name']= f'{particle_name1}-{particle_name2}'
1396            self.df.at [index,'bond_object'] = bond_object
1397            self.df.at [index,'l0'] = l0
1398            self.add_value_to_df(index=index,
1399                                    key=('pmb_type',''),
1400                                    new_value='bond')
1401            self.add_value_to_df(index=index,
1402                                    key=('parameters_of_the_potential',''),
1403                                    new_value=bond_object.get_params(),
1404                                    non_standard_value=True)
1405
1406        return

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

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

Currently, only HARMONIC and FENE bonds are supported.

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

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

For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:

- d_r_max (`obj`): Maximal stretching length for FENE. It should have 
units of length using the `pmb.units` UnitRegistry. Default 'None'.
def define_default_bond( self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None):
1409    def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None):
1410        """
1411        Asigns `bond` in `pmb.df` as the default bond.
1412        The LJ parameters can be optionally provided to calculate the initial bond length
1413
1414        Args:
1415            bond_type(`str`): label to identify the potential to model the bond.
1416            bond_parameters(`dict`): parameters of the potential of the bond.
1417            sigma(`float`, optional): LJ sigma of the interaction between the particles.
1418            epsilon(`float`, optional): LJ epsilon for the interaction between the particles.
1419            cutoff(`float`, optional): cutoff-radius of the LJ interaction.
1420            offset(`float`, optional): offset of the LJ interaction.
1421
1422        Note:
1423            - Currently, only harmonic and FENE bonds are supported. 
1424        """
1425        
1426        bond_object=self.create_bond_in_espresso(bond_type, bond_parameters)
1427
1428        if epsilon is None:
1429            epsilon=1*self.units('reduced_energy')
1430        if sigma is None:
1431            sigma=1*self.units('reduced_length')
1432        if cutoff is None:
1433            cutoff=2**(1.0/6.0)*self.units('reduced_length')
1434        if offset is None:
1435            offset=0*self.units('reduced_length')
1436        l0 = self.calculate_initial_bond_length(bond_object = bond_object, 
1437                                                bond_type   = bond_type, 
1438                                                epsilon     = epsilon,
1439                                                sigma       = sigma,
1440                                                cutoff      = cutoff,
1441                                                offset      = offset)
1442
1443        if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'):
1444            return
1445        if len(self.df.index) != 0:
1446            index = max(self.df.index)+1
1447        else:
1448            index = 0
1449        self.df.at [index,'name']        = 'default'
1450        self.df.at [index,'bond_object'] = bond_object
1451        self.df.at [index,'l0']          = l0
1452        self.add_value_to_df(index       = index,
1453                            key          = ('pmb_type',''),
1454                            new_value    = 'bond')
1455        self.add_value_to_df(index       = index,
1456                            key          = ('parameters_of_the_potential',''),
1457                            new_value    = bond_object.get_params(),
1458                            non_standard_value=True)
1459        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):
1461    def define_molecule(self, name, residue_list):
1462        """
1463        Defines a pyMBE object of type `molecule` in `pymbe.df`.
1464
1465        Args:
1466            name(`str`): Unique label that identifies the `molecule`.
1467            residue_list(`list` of `str`): List of the `name`s of the `residue`s  in the sequence of the `molecule`.  
1468        """
1469        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'):
1470            return
1471        index = len(self.df)
1472        self.df.at [index,'name'] = name
1473        self.df.at [index,'pmb_type'] = 'molecule'
1474        self.df.at [index,('residue_list','')] = residue_list
1475        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):
1477    def define_particle_entry_in_df(self,name):
1478        """
1479        Defines a particle entry in pmb.df.
1480
1481        Args:
1482            name(`str`): Unique label that identifies this particle type.
1483
1484        Returns:
1485            index(`int`): Index of the particle in pmb.df  
1486        """
1487
1488        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'):
1489            index = self.df[self.df['name']==name].index[0]                                   
1490        else:
1491            index = len(self.df)
1492            self.df.at [index, 'name'] = name
1493            self.df.at [index,'pmb_type'] = 'particle'
1494        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):
1496    def define_particle(self, name, z=0, acidity=None, pka=None, sigma=None, epsilon=None, cutoff=None, offset=None,verbose=True,overwrite=False):
1497        """
1498        Defines the properties of a particle object.
1499
1500        Args:
1501            name(`str`): Unique label that identifies this particle type.  
1502            z(`int`, optional): Permanent charge number of this particle type. Defaults to 0.
1503            acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None.
1504            pka(`float`, optional): If `particle` is an acid or a base, it defines its  pka-value. Defaults to None.
1505            sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
1506            cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
1507            offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None.
1508            epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to None.
1509            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
1510            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
1511
1512        Note:
1513            - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units.
1514            - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units.
1515            - `cutoff` defaults to `2**(1./6.) reduced_length`. 
1516            - `offset` defaults to 0.
1517            - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions.
1518            - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`.
1519        """ 
1520        index=self.define_particle_entry_in_df(name=name)
1521        
1522        # If `cutoff` and `offset` are not defined, default them to the following values
1523        if cutoff is None:
1524            cutoff=self.units.Quantity(2**(1./6.), "reduced_length")
1525        if offset is None:
1526            offset=self.units.Quantity(0, "reduced_length")
1527
1528        # Define LJ parameters
1529        parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"},
1530                                        "cutoff":{"value": cutoff, "dimensionality": "[length]"},
1531                                        "offset":{"value": offset, "dimensionality": "[length]"},
1532                                        "epsilon":{"value": epsilon, "dimensionality": "[energy]"},}
1533
1534        for parameter_key in parameters_with_dimensionality.keys():
1535            if parameters_with_dimensionality[parameter_key]["value"] is not None:
1536                self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 
1537                                          expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"])
1538                self.add_value_to_df(key=(parameter_key,''),
1539                                    index=index,
1540                                    new_value=parameters_with_dimensionality[parameter_key]["value"],
1541                                    verbose=verbose,
1542                                    overwrite=overwrite)
1543
1544        # Define particle acid/base properties
1545        self.set_particle_acidity(name=name, 
1546                                acidity=acidity, 
1547                                default_charge_number=z, 
1548                                pka=pka,
1549                                verbose=verbose,
1550                                overwrite=overwrite)
1551        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):
1553    def define_particles(self, parameters, overwrite=False, verbose=True):
1554        '''
1555        Defines a particle object in pyMBE for each particle name in `particle_names`
1556
1557        Args:
1558            parameters(`dict`):  dictionary with the particle parameters. 
1559            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1560            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
1561
1562        Note:
1563            - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},}
1564        '''
1565        if not parameters:
1566            return 0
1567        for particle_name in parameters.keys():
1568            parameters[particle_name]["overwrite"]=overwrite
1569            parameters[particle_name]["verbose"]=verbose
1570            self.define_particle(**parameters[particle_name])
1571        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):
1573    def define_peptide(self, name, sequence, model):
1574        """
1575        Defines a pyMBE object of type `peptide` in the `pymbe.df`.
1576
1577        Args:
1578            name (`str`): Unique label that identifies the `peptide`.
1579            sequence (`string`): Sequence of the `peptide`.
1580            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1581        """
1582        if not self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'):
1583            valid_keys = ['1beadAA','2beadAA']
1584            if model not in valid_keys:
1585                raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA')
1586            clean_sequence = self.protein_sequence_parser(sequence=sequence)    
1587            residue_list = self.define_AA_residues(sequence=clean_sequence,
1588                                                    model=model)
1589            self.define_molecule(name = name, residue_list=residue_list)
1590            index = self.df.loc[self.df['name'] == name].index.item() 
1591            self.df.at [index,'model'] = model
1592            self.df.at [index,('sequence','')] = clean_sequence
1593        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):
1595    def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False, verbose=True):
1596        """
1597        Defines a globular protein pyMBE object  in `pymbe.df`.
1598
1599        Args:
1600            name (`str`): Unique label that identifies the protein.
1601            model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
1602            topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value}
1603            lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca".
1604            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
1605            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
1606
1607        Note:
1608            - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential.
1609        """
1610
1611        # Sanity checks
1612        valid_model_keys = ['1beadAA','2beadAA']
1613        valid_lj_setups = ["wca"]
1614        if model not in valid_model_keys:
1615            raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}')
1616        if lj_setup_mode not in valid_lj_setups:
1617            raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}')
1618        if lj_setup_mode == "wca":
1619            sigma = 1*self.units.Quantity("reduced_length")
1620            epsilon = 1*self.units.Quantity("reduced_energy")
1621        part_dict={}
1622        sequence=[]
1623        metal_ions_charge_number_map=self.get_metal_ions_charge_number_map()
1624        for particle in topology_dict.keys():
1625            particle_name = re.split(r'\d+', particle)[0] 
1626            if particle_name not in part_dict.keys():
1627                if lj_setup_mode == "wca":
1628                    part_dict[particle_name]={"sigma": sigma,
1629                                        "offset": topology_dict[particle]['radius']*2-sigma,
1630                                        "epsilon": epsilon,
1631                                        "name": particle_name}
1632                if self.check_if_metal_ion(key=particle_name):
1633                    z=metal_ions_charge_number_map[particle_name]
1634                else:
1635                    z=0
1636                part_dict[particle_name]["z"]=z
1637            
1638            if self.check_aminoacid_key(key=particle_name):
1639                sequence.append(particle_name) 
1640            
1641        self.define_particles(parameters=part_dict,
1642                            overwrite=overwrite,  
1643                            verbose=verbose)
1644        residue_list = self.define_AA_residues(sequence=sequence, 
1645                                               model=model)
1646        index = len(self.df)
1647        self.df.at [index,'name'] = name
1648        self.df.at [index,'pmb_type'] = 'protein'
1649        self.df.at [index,'model'] = model
1650        self.df.at [index,('sequence','')] = sequence  
1651        self.df.at [index,('residue_list','')] = residue_list    
1652        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):
1654    def define_residue(self, name, central_bead, side_chains):
1655        """
1656        Defines a pyMBE object of type `residue` in `pymbe.df`.
1657
1658        Args:
1659            name(`str`): Unique label that identifies the `residue`.
1660            central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`.
1661            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.
1662        """
1663        if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'):
1664            return
1665        index = len(self.df)
1666        self.df.at [index, 'name'] = name
1667        self.df.at [index,'pmb_type'] = 'residue'
1668        self.df.at [index,'central_bead'] = central_bead
1669        self.df.at [index,('side_chains','')] = side_chains
1670        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):
1672    def destroy_pmb_object_in_system(self, name, espresso_system):
1673        """
1674        Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 
1675
1676        Args:
1677            name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`.
1678            espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library.
1679
1680        Note:
1681            - 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.
1682        """
1683        allowed_objects = ['particle','residue','molecule']
1684        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1685        if pmb_type not in allowed_objects:
1686            raise ValueError('Object type not supported, supported types are ', allowed_objects)
1687        if pmb_type == 'particle':
1688            particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 
1689                                                  & (self.df['molecule_id'].isna())]
1690            particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna())
1691                                                 & (self.df['molecule_id'].isna())].particle_id.tolist()
1692            for particle_id in particle_ids_list:
1693                espresso_system.part.by_id(particle_id).remove()
1694            self.df = self.df.drop(index=particles_index)
1695        if pmb_type == 'residue':
1696            residues_id = self.df.loc[self.df['name']== name].residue_id.to_list()
1697            for residue_id in residues_id:
1698                molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0]
1699                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1700                self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index)
1701                for particle_id in particle_ids_list:
1702                    espresso_system.part.by_id(particle_id).remove()
1703                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)    
1704        if pmb_type == 'molecule':
1705            molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list()
1706            for molecule_id in molecules_id:
1707                molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']=="molecule")].name.values[0]
1708                particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"]
1709                self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index)
1710                for particle_id in particle_ids_list:
1711                    espresso_system.part.by_id(particle_id).remove()   
1712                    self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index)             
1713        
1714        self.df.reset_index(drop=True,inplace=True)
1715
1716        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):
1718    def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200):
1719        """
1720        Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration.
1721        To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an
1722        ion pair and the concentrations of the various species is solved numerically using a self-consistent approach.
1723        More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831).
1724        This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237).
1725
1726        Args:
1727            pH_res('float'): pH-value in the reservoir.
1728            c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
1729            activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
1730
1731        Returns:
1732            cH_res('pint.Quantity'): Concentration of H+ ions.
1733            cOH_res('pint.Quantity'): Concentration of OH- ions.
1734            cNa_res('pint.Quantity'): Concentration of Na+ ions.
1735            cCl_res('pint.Quantity'): Concentration of Cl- ions.
1736        """
1737
1738        self_consistent_run = 0
1739        cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+
1740
1741        def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res):
1742            #Calculate and initial guess for the concentrations of various species based on ideal gas estimate
1743            cOH_res = self.Kw / cH_res 
1744            cNa_res = None
1745            cCl_res = None
1746            if cOH_res>=cH_res:
1747                #adjust the concentration of sodium if there is excess OH- in the reservoir:
1748                cNa_res = c_salt_res + (cOH_res-cH_res)
1749                cCl_res = c_salt_res
1750            else:
1751                # adjust the concentration of chloride if there is excess H+ in the reservoir
1752                cCl_res = c_salt_res + (cH_res-cOH_res)
1753                cNa_res = c_salt_res
1754                
1755            def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res):
1756                nonlocal max_number_sc_runs, self_consistent_run
1757                if self_consistent_run<max_number_sc_runs:
1758                    self_consistent_run+=1
1759                    ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1760                    cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res))
1761                    if cOH_res>=cH_res:
1762                        #adjust the concentration of sodium if there is excess OH- in the reservoir:
1763                        cNa_res = c_salt_res + (cOH_res-cH_res)
1764                        cCl_res = c_salt_res
1765                    else:
1766                        # adjust the concentration of chloride if there is excess H+ in the reservoir
1767                        cCl_res = c_salt_res + (cH_res-cOH_res)
1768                        cNa_res = c_salt_res
1769                    return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1770                else:
1771                    return cH_res, cOH_res, cNa_res, cCl_res
1772            return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res)
1773
1774        cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1775        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1776        determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1777
1778        while abs(determined_pH-pH_res)>1e-6:
1779            if determined_pH > pH_res:
1780                cH_res *= 1.005
1781            else:
1782                cH_res /= 1.003
1783            cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res)
1784            ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
1785            determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res)))
1786            self_consistent_run=0
1787
1788        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):
1790    def enable_motion_of_rigid_object(self, name, espresso_system):
1791        '''
1792        Enables the motion of the rigid object `name` in the `espresso_system`.
1793
1794        Args:
1795            name(`str`): Label of the object.
1796            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
1797
1798        Note:
1799            - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"].
1800        '''
1801        print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]')
1802        pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0]
1803        if pmb_type != 'protein':
1804            raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein')
1805        molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list()
1806        for molecule_id in molecule_ids_list:    
1807            particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list()
1808            center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system)
1809            rigid_object_center = espresso_system.part.add(pos=center_of_mass,
1810                                                           rotation=[True,True,True], 
1811                                                           type=self.propose_unused_type())
1812            
1813            rigid_object_center.mass = len(particle_ids_list)
1814            momI = 0
1815            for pid in particle_ids_list:
1816                momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2)
1817
1818            rigid_object_center.rinertia = np.ones(3) * momI
1819            
1820            for particle_id in particle_ids_list:
1821                pid = espresso_system.part.by_id(particle_id)
1822                pid.vs_auto_relate_to(rigid_object_center.id)
1823        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):
1825    def filter_df(self, pmb_type):
1826        """
1827        Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns.
1828        
1829        Args:
1830            pmb_type(`str`): pmb_object_type to filter in `pmb.df`.
1831
1832        Returns:
1833            pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`.
1834        """
1835        pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type]
1836        pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1)
1837        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):
1839    def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False):
1840        """
1841        Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
1842        
1843        Args:
1844            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
1845            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
1846            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'. 
1847
1848        Returns:
1849            bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` if a matching bond exists
1850
1851        Note:
1852            - If `use_default_bond`=`True`, it returns "default" if no key is found.
1853        """
1854        bond_keys = [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']
1855        bond_defined=False
1856        for bond_key in bond_keys:
1857            if bond_key in self.df["name"].values:
1858                bond_defined=True
1859                correct_key=bond_key
1860                break
1861        if bond_defined:
1862            return correct_key
1863        elif use_default_bond:
1864            return 'default'
1865        else:
1866            return 

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

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

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

Note:
  • If use_default_bond=True, it returns "default" if no key is found.
def find_value_from_es_type(self, es_type, column_name):
1868    def find_value_from_es_type(self, es_type, column_name):
1869        """
1870        Finds a value in `pmb.df` for a `column_name` and `es_type` pair.
1871
1872        Args:
1873            es_type(`int`): value of the espresso type
1874            column_name(`str`): name of the column in `pymbe.df`
1875
1876        Returns:
1877            Value in `pymbe.df` matching  `column_name` and `es_type`
1878        """
1879        idx = pd.IndexSlice
1880        for state in ['state_one', 'state_two']:
1881            index = np.where(self.df[(state, 'es_type')] == es_type)[0]
1882            if len(index) > 0:
1883                if column_name == 'label':
1884                    label = self.df.loc[idx[index[0]], idx[(state,column_name)]]
1885                    return label
1886                else: 
1887                    column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]]
1888                    return column_name_value
1889        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):
1891    def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples):
1892        """
1893        Generates coordinates outside a sphere centered at `center`.
1894
1895        Args:
1896            center(`lst`): Coordinates of the center of the sphere.
1897            radius(`float`): Radius of the sphere.
1898            max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates.
1899            n_samples(`int`): Number of sample points.
1900
1901        Returns:
1902            coord_list(`lst`): Coordinates of the sample points.
1903        """
1904        if not radius > 0: 
1905            raise ValueError (f'The value of {radius} must be a positive value')
1906        if not radius < max_dist:
1907            raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))')
1908        coord_list = []
1909        counter = 0
1910        while counter<n_samples:
1911            coord = self.generate_random_points_in_a_sphere(center=center, 
1912                                            radius=max_dist,
1913                                            n_samples=1)[0]
1914            if np.linalg.norm(coord-np.asarray(center))>=radius:
1915                coord_list.append (coord)
1916                counter += 1
1917        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):
1919    def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False):
1920        """
1921        Uniformly samples points from a hypersphere. If on_surface is set to True, the points are
1922        uniformly sampled from the surface of the hypersphere.
1923        
1924        Args:
1925            center(`lst`): Array with the coordinates of the center of the spheres.
1926            radius(`float`): Radius of the sphere.
1927            n_samples(`int`): Number of sample points to generate inside the sphere.
1928            on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere.
1929
1930        Returns:
1931            samples(`list`): Coordinates of the sample points inside the hypersphere.
1932        """
1933        # initial values
1934        center=np.array(center)
1935        d = center.shape[0]
1936        # sample n_samples points in d dimensions from a standard normal distribution
1937        samples = self.rng.normal(size=(n_samples, d))
1938        # make the samples lie on the surface of the unit hypersphere
1939        normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis]
1940        samples /= normalize_radii
1941        if not on_surface:
1942            # make the samples lie inside the hypersphere with the correct density
1943            uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis]
1944            new_radii = np.power(uniform_points, 1/d)
1945            samples *= new_radii
1946        # scale the points to have the correct radius and center
1947        samples = samples * radius + center
1948        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):
1950    def generate_trial_perpendicular_vector(self,vector,magnitude):
1951        """
1952        Generates an orthogonal vector to the input `vector`.
1953
1954        Args:
1955            vector(`lst`): arbitrary vector.
1956            magnitude(`float`): magnitude of the orthogonal vector.
1957            
1958        Returns:
1959            (`lst`): Orthogonal vector with the same magnitude as the input vector.
1960        """ 
1961        np_vec = np.array(vector) 
1962        np_vec /= np.linalg.norm(np_vec) 
1963        if np.all(np_vec == 0):
1964            raise ValueError('Zero vector')
1965        # Generate a random vector 
1966        random_vector = self.generate_random_points_in_a_sphere(radius=1, 
1967                                                                center=[0,0,0],
1968                                                                n_samples=1, 
1969                                                                on_surface=True)[0]
1970        # Project the random vector onto the input vector and subtract the projection
1971        projection = np.dot(random_vector, np_vec) * np_vec
1972        perpendicular_vector = random_vector - projection
1973        # Normalize the perpendicular vector to have the same magnitude as the input vector
1974        perpendicular_vector /= np.linalg.norm(perpendicular_vector) 
1975        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):
1977    def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
1978        """
1979        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length.
1980        If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead.
1981        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.
1982
1983        Args:
1984            particle_name1(str): label of the type of the first particle type of the bonded particles.
1985            particle_name2(str): label of the type of the second particle type of the bonded particles.
1986            hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
1987            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. 
1988
1989        Returns:
1990            l0(`pint.Quantity`): bond length
1991        
1992        Note:
1993            - 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`.
1994            - If `hard_check`=`True` stops the code when no bond is found.
1995        """
1996        bond_key = self.find_bond_key(particle_name1=particle_name1, 
1997                                    particle_name2=particle_name2, 
1998                                    use_default_bond=use_default_bond)
1999        if bond_key:
2000            return self.df[self.df['name']==bond_key].l0.values[0]
2001        else:
2002            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2003            if hard_check:
2004                sys.exit(1)
2005            else:
2006                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):
2008    def get_charge_number_map(self):
2009        '''
2010        Gets the charge number of each `espresso_type` in `pymbe.df`.
2011        
2012        Returns:
2013            charge_number_map(`dict`): {espresso_type: z}.
2014        '''
2015        if self.df.state_one['es_type'].isnull().values.any():         
2016            df_state_one = self.df.state_one.dropna()     
2017            df_state_two = self.df.state_two.dropna()  
2018        else:    
2019            df_state_one = self.df.state_one
2020            if self.df.state_two['es_type'].isnull().values.any():
2021                df_state_two = self.df.state_two.dropna()   
2022            else:
2023                df_state_two = self.df.state_two
2024        state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values)
2025        state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values)
2026        charge_number_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2027        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'):
2029    def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'):
2030        """
2031        Returns the Lennard-Jones parameters for the interaction between the particle types given by 
2032        `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule.
2033
2034        Args:
2035            particle_name1 (str): label of the type of the first particle type
2036            particle_name2 (str): label of the type of the second particle type
2037            combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
2038
2039        Returns:
2040            {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
2041
2042        Note:
2043            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
2044            - 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.
2045        """
2046        supported_combining_rules=["Lorentz-Berthelot"]
2047        lj_parameters_keys=["sigma","epsilon","offset","cutoff"]
2048        if combining_rule not in supported_combining_rules:
2049            raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}")
2050        lj_parameters={}
2051        for key in lj_parameters_keys:
2052            lj_parameters[key]=[]
2053        # Search the LJ parameters of the type pair
2054        for name in [particle_name1,particle_name2]:
2055            for key in lj_parameters_keys:
2056                lj_parameters[key].append(self.df[self.df.name == name][key].values[0])
2057        # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
2058        if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]):
2059            return {}
2060        # Apply combining rule
2061        if combining_rule == 'Lorentz-Berthelot':
2062            lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2
2063            lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2
2064            lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2
2065            lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1])
2066        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):
2068    def get_metal_ions_charge_number_map(self):
2069        """
2070        Gets a map with the charge numbers of all the metal ions supported.
2071
2072        Returns:
2073            metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number}
2074
2075        """
2076        metal_charge_number_map = {"Ca": 2}
2077        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):
2079    def get_particle_id_map(self, object_name):
2080        '''
2081        Gets all the ids associated with the object with name `object_name` in `pmb.df`
2082
2083        Args:
2084            object_name(`str`): name of the object
2085        
2086        Returns:
2087            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]}, }
2088        '''
2089        object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0]
2090        valid_types = ["particle", "molecule", "residue", "protein"]
2091        if object_type not in valid_types:
2092            raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}")
2093        id_list = []
2094        mol_map = {}
2095        res_map = {}
2096        def do_res_map(res_ids):
2097            for res_id in res_ids:
2098                res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2099                res_map[res_id]=res_list
2100            return res_map
2101        if object_type in ['molecule', 'protein']:
2102            mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist()
2103            for mol_id in mol_ids:
2104                res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist())
2105                res_map=do_res_map(res_ids=res_ids)    
2106                mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist()
2107                id_list+=mol_list
2108                mol_map[mol_id]=mol_list
2109        elif object_type == 'residue':     
2110            res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist()
2111            res_map=do_res_map(res_ids=res_ids)
2112            id_list=[]
2113            for res_id_list in res_map.values():
2114                id_list+=res_id_list
2115        elif object_type == 'particle':
2116            id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist()
2117        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):
2119    def get_pka_set(self):
2120        '''
2121        Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df`
2122        
2123        Returns:
2124            pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}
2125        '''
2126        titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna()
2127        pka_set = {}
2128        for index in titratables_AA_df.name.keys():
2129            name = titratables_AA_df.name[index]
2130            pka_value = titratables_AA_df.pka[index]
2131            acidity = titratables_AA_df.acidity[index]   
2132            pka_set[name] = {'pka_value':pka_value,'acidity':acidity}
2133        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):
2135    def get_radius_map(self, dimensionless=True):
2136        '''
2137        Gets the effective radius of each `espresso type` in `pmb.df`. 
2138
2139        Args:
2140            dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False.
2141        
2142        Returns:
2143            radius_map(`dict`): {espresso_type: radius}.
2144
2145        Note:
2146            The radius corresponds to (sigma+offset)/2
2147        '''
2148        df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates()
2149        df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates()
2150        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)
2151        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)
2152        radius_map  = pd.concat([state_one,state_two],axis=0).to_dict()  
2153        if dimensionless:
2154            for key in radius_map:
2155                radius_map[key] = radius_map[key].magnitude
2156        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):
2158    def get_reduced_units(self):
2159        """
2160        Returns the  current set of reduced units defined in pyMBE.
2161
2162        Returns:
2163            reduced_units_text(`str`): text with information about the current set of reduced units.
2164
2165        """
2166        unit_length=self.units.Quantity(1,'reduced_length')
2167        unit_energy=self.units.Quantity(1,'reduced_energy')
2168        unit_charge=self.units.Quantity(1,'reduced_charge')
2169        reduced_units_text = "\n".join(["Current set of reduced units:",
2170                                       f"{unit_length.to('nm'):.5g} = {unit_length}",
2171                                       f"{unit_energy.to('J'):.5g} = {unit_energy}",
2172                                       f"{unit_charge.to('C'):.5g} = {unit_charge}",
2173                                       f"Temperature: {(self.kT/self.kB).to('K'):.5g}"
2174                                        ])   
2175        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):
2177    def get_resource(self, path):
2178        '''
2179        Locate a file resource of the pyMBE package.
2180
2181        Args:
2182            path(`str`): Relative path to the resource
2183
2184        Returns:
2185            path(`str`): Absolute path to the resource
2186
2187        '''
2188        import os
2189        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):
2191    def get_type_map(self):
2192        """
2193        Gets all different espresso types assigned to particles  in `pmb.df`.
2194        
2195        Returns:
2196            type_map(`dict`): {"name": espresso_type}.
2197        """
2198        if self.df.state_one['es_type'].isnull().values.any():         
2199            df_state_one = self.df.state_one.dropna(how='all')     
2200            df_state_two = self.df.state_two.dropna(how='all')  
2201        else:    
2202            df_state_one = self.df.state_one
2203            if self.df.state_two['es_type'].isnull().values.any():
2204                df_state_two = self.df.state_two.dropna(how='all')   
2205            else:
2206                df_state_two = self.df.state_two
2207        state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label)
2208        state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label)
2209        type_map  = pd.concat([state_one,state_two],axis=0).to_dict()
2210        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):
2212    def load_interaction_parameters(self, filename, verbose=False, overwrite=False):
2213        """
2214        Loads the interaction parameters stored in `filename` into `pmb.df`
2215        
2216        Args:
2217            filename(`str`): name of the file to be read
2218            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to False.
2219            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2220        """
2221        without_units = ['z','es_type','acidity']
2222        with_units = ['sigma','epsilon','offset','cutoff']
2223
2224        with open(filename, 'r') as f:
2225            interaction_data = json.load(f)
2226            interaction_parameter_set = interaction_data["data"]
2227
2228        for key in interaction_parameter_set:
2229            param_dict=interaction_parameter_set[key]
2230            object_type=param_dict['object_type']
2231            if object_type == 'particle':
2232                not_required_attributes={}    
2233                for not_required_key in without_units+with_units:
2234                    if not_required_key in param_dict.keys():
2235                        if not_required_key in with_units:
2236                            not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key))
2237                        elif not_required_key in without_units:
2238                            not_required_attributes[not_required_key]=param_dict.pop(not_required_key)
2239                    else:
2240                        not_required_attributes[not_required_key]=None
2241                self.define_particle(name=param_dict.pop('name'),
2242                                z=not_required_attributes.pop('z'),
2243                                sigma=not_required_attributes.pop('sigma'),
2244                                offset=not_required_attributes.pop('offset'),
2245                                cutoff=not_required_attributes.pop('cutoff'),
2246                                acidity=not_required_attributes.pop('acidity'),
2247                                epsilon=not_required_attributes.pop('epsilon'),
2248                                verbose=verbose,
2249                                overwrite=overwrite)
2250            elif object_type == 'residue':
2251                self.define_residue (name = param_dict.pop('name'),
2252                                    central_bead = param_dict.pop('central_bead'),
2253                                    side_chains = param_dict.pop('side_chains'))
2254            elif object_type == 'molecule':
2255                self.define_molecule(name=param_dict.pop('name'),
2256                                    residue_list=param_dict.pop('residue_list'))
2257            elif object_type == 'peptide':
2258                self.define_peptide(name=param_dict.pop('name'),
2259                                    sequence=param_dict.pop('sequence'),
2260                                    model=param_dict.pop('model'))
2261            elif object_type == 'bond':
2262                particle_pairs = param_dict.pop('particle_pairs')
2263                bond_parameters = param_dict.pop('bond_parameters')
2264                bond_type = param_dict.pop('bond_type')
2265                if bond_type == 'harmonic':
2266                    k = self.create_variable_with_units(variable=bond_parameters.pop('k'))
2267                    r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0'))
2268                    bond = {'r_0'    : r_0,
2269                            'k'      : k,
2270                            }
2271
2272                elif bond_type == 'FENE':
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                    d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max'))
2276                    bond =  {'r_0'    : r_0,
2277                             'k'      : k,
2278                             'd_r_max': d_r_max,
2279                             }
2280                else:
2281                    raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds")
2282                self.define_bond(bond_type=bond_type, 
2283                                bond_parameters=bond, 
2284                                particle_pairs=particle_pairs)
2285            else:
2286                raise ValueError(object_type+' is not a known pmb object type')
2287            
2288        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):
2290    def load_pka_set(self, filename, verbose=False, overwrite=True):
2291        """
2292        Loads the pka_set stored in `filename` into `pmb.df`.
2293        
2294        Args:
2295            filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}.
2296            verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to False.
2297            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 
2298        """
2299        with open(filename, 'r') as f:
2300            pka_data = json.load(f)
2301            pka_set = pka_data["data"]
2302
2303        self.check_pka_set(pka_set=pka_set)
2304
2305        for key in pka_set:
2306            acidity = pka_set[key]['acidity']
2307            pka_value = pka_set[key]['pka_value']
2308            self.set_particle_acidity(name=key, 
2309                                      acidity=acidity, 
2310                                      pka=pka_value, 
2311                                      verbose=verbose, 
2312                                      overwrite=overwrite)
2313        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 propose_unused_type(self):
2316    def propose_unused_type(self):
2317        """
2318        Searches in `pmb.df` all the different particle types defined and returns a new one.
2319
2320        Returns:
2321            unused_type(`int`): unused particle type
2322        """
2323        type_map=self.get_type_map()
2324        if type_map == {}:    
2325            unused_type = 0
2326        else:
2327            unused_type=max(type_map.values())+1    
2328        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):
2330    def protein_sequence_parser(self, sequence):
2331        '''
2332        Parses `sequence` to the one letter code for amino acids.
2333        
2334        Args:
2335            sequence(`str` or `lst`): Sequence of the amino acid. 
2336
2337        Returns:
2338            clean_sequence(`lst`): `sequence` using the one letter code.
2339        
2340        Note:
2341            - Accepted formats for `sequence` are:
2342                - `lst` with one letter or three letter code of each aminoacid in each element
2343                - `str` with the sequence using the one letter code
2344                - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-"
2345        
2346        '''
2347        # Aminoacid key
2348        keys={"ALA": "A",
2349                "ARG": "R",
2350                "ASN": "N",
2351                "ASP": "D",
2352                "CYS": "C",
2353                "GLU": "E",
2354                "GLN": "Q",
2355                "GLY": "G",
2356                "HIS": "H",
2357                "ILE": "I",
2358                "LEU": "L",
2359                "LYS": "K",
2360                "MET": "M",
2361                "PHE": "F",
2362                "PRO": "P",
2363                "SER": "S",
2364                "THR": "T",
2365                "TRP": "W",
2366                "TYR": "Y",
2367                "VAL": "V",
2368                "PSER": "J",
2369                "PTHR": "U",
2370                "PTyr": "Z",
2371                "NH2": "n",
2372                "COOH": "c"}
2373        clean_sequence=[]
2374        if isinstance(sequence, str):
2375            if sequence.find("-") != -1:
2376                splited_sequence=sequence.split("-")
2377                for residue in splited_sequence:
2378                    if len(residue) == 1:
2379                        if residue in keys.values():
2380                            residue_ok=residue
2381                        else:
2382                            if residue.upper() in keys.values():
2383                                residue_ok=residue.upper()
2384                            else:
2385                                raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence")
2386                        clean_sequence.append(residue_ok)
2387                    else:
2388                        if residue in keys.keys():
2389                            clean_sequence.append(keys[residue])
2390                        else:
2391                            if residue.upper() in keys.keys():
2392                                clean_sequence.append(keys[residue.upper()])
2393                            else:
2394                                raise ValueError("Unknown  code for a residue: ", residue, " please review the input sequence")
2395            else:
2396                for residue in sequence:
2397                    if residue in keys.values():
2398                        residue_ok=residue
2399                    else:
2400                        if residue.upper() in keys.values():
2401                            residue_ok=residue.upper()
2402                        else:
2403                            raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence")
2404                    clean_sequence.append(residue_ok)
2405        if isinstance(sequence, list):
2406            for residue in sequence:
2407                if residue in keys.values():
2408                    residue_ok=residue
2409                else:
2410                    if residue.upper() in keys.values():
2411                        residue_ok=residue.upper()
2412                    elif (residue.upper() in keys.keys()):
2413                        residue_ok= keys[residue.upper()]
2414                    else:
2415                        raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence")
2416                clean_sequence.append(residue_ok)
2417        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):
2419    def read_pmb_df (self,filename):
2420        """
2421        Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE.
2422
2423        Args:
2424            filename(`str`): path to file.
2425
2426        Note:
2427            This function only accepts files with CSV format. 
2428        """
2429        
2430        if filename.rsplit(".", 1)[1] != "csv":
2431            raise ValueError("Only files with CSV format are supported")
2432        df = pd.read_csv (filename,header=[0, 1], index_col=0)
2433        columns_names = self.setup_df()
2434        
2435        multi_index = pd.MultiIndex.from_tuples(columns_names)
2436        df.columns = multi_index
2437        
2438        self.df = self.convert_columns_to_original_format(df)
2439        
2440        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):
2442    def read_protein_vtf_in_df (self,filename,unit_length=None):
2443        """
2444        Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`.
2445
2446        Args:
2447            filename(`str`): Path to the vtf file with the coarse-grained model.
2448            unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None.
2449
2450        Returns:
2451            topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}
2452
2453        Note:
2454            - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom.
2455        """
2456
2457        print (f'Loading protein coarse grain model file: {filename}')
2458
2459        coord_list = []
2460        particles_dict = {}
2461
2462        if unit_length is None:
2463            unit_length = 1 * self.units.angstrom 
2464
2465        with open (filename,'r') as protein_model:
2466            for line in protein_model :
2467                line_split = line.split()
2468                if line_split:
2469                    line_header = line_split[0]
2470                    if line_header == 'atom':
2471                        atom_id  = line_split[1]
2472                        atom_name = line_split[3]
2473                        atom_resname = line_split[5]
2474                        chain_id = line_split[9]
2475                        radius = float(line_split [11])*unit_length 
2476                        particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius]
2477                    elif line_header.isnumeric(): 
2478                        atom_coord = line_split[1:] 
2479                        atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord]
2480                        coord_list.append (atom_coord)
2481
2482        numbered_label = []
2483        i = 0   
2484        
2485        for atom_id in particles_dict.keys():
2486    
2487            if atom_id == 1:
2488                atom_name = particles_dict[atom_id][0]
2489                numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2490                numbered_label.append(numbered_name)
2491
2492            elif atom_id != 1: 
2493            
2494                if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]:
2495                    i += 1                    
2496                    count = 1
2497                    atom_name = particles_dict[atom_id][0]
2498                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2499                    numbered_label.append(numbered_name)
2500                    
2501                elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]:
2502                    if count == 2 or particles_dict[atom_id][1] == 'GLY':
2503                        i +=1  
2504                        count = 0
2505                    atom_name = particles_dict[atom_id][0]
2506                    numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]]
2507                    numbered_label.append(numbered_name)
2508                    count +=1
2509
2510        topology_dict = {}
2511
2512        for i in range (0, len(numbered_label)):   
2513            topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] ,
2514                                                    'chain_id':numbered_label[i][1],
2515                                                    'radius':numbered_label[i][2] }
2516
2517        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):
2519    def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) :
2520        """
2521        Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it.
2522        If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead.
2523        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.
2524
2525        Args:
2526            particle_name1(`str`): label of the type of the first particle type of the bonded particles.
2527            particle_name2(`str`): label of the type of the second particle type of the bonded particles.
2528            hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 
2529            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. 
2530
2531        Returns:
2532            bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library.
2533        
2534        Note:
2535            - 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`.
2536            - If `hard_check`=`True` stops the code when no bond is found.
2537        """
2538        
2539        bond_key = self.find_bond_key(particle_name1=particle_name1, 
2540                                    particle_name2=particle_name2, 
2541                                    use_default_bond=use_default_bond)
2542        if use_default_bond:
2543            if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'):
2544                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")
2545        if bond_key:
2546            return self.df[self.df['name']==bond_key].bond_object.values[0]
2547        else:
2548            print(f"Bond not defined between particles {particle_name1} and {particle_name2}")
2549            if hard_check:
2550                sys.exit(1)
2551            else:
2552                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):
2554    def search_particles_in_residue(self, residue_name):
2555        '''
2556        Searches for all particles in a given residue of name `residue_name`.
2557
2558        Args:
2559            residue_name (`str`): name of the residue to be searched
2560
2561        Returns:
2562            list_of_particles_in_residue (`lst`): list of the names of all particles in the residue
2563
2564        Note:
2565            - 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.
2566 
2567        '''
2568        index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 
2569        central_bead = self.df.at [index_residue, ('central_bead', '')]
2570        list_of_side_chains = self.df.at [index_residue, ('side_chains', '')]
2571
2572        list_of_particles_in_residue = []
2573        list_of_particles_in_residue.append(central_bead)
2574
2575        for side_chain in list_of_side_chains: 
2576            object_type = self.df[self.df['name']==side_chain].pmb_type.values[0]
2577
2578            if object_type == "residue":
2579                list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain)
2580                list_of_particles_in_residue += list_of_particles_in_side_chain_residue
2581            elif object_type == "particle":
2582                list_of_particles_in_residue.append(side_chain)
2583
2584        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):
2586    def set_particle_acidity(self, name, acidity=None, default_charge_number=0, pka=None, verbose=True, overwrite=True):
2587        """
2588        Sets the particle acidity including the charges in each of its possible states. 
2589
2590        Args:
2591            name(`str`): Unique label that identifies the `particle`. 
2592            acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None.
2593            default_charge_number (`int`): Charge number of the particle. Defaults to 0.
2594            pka(`float`, optional):  If `particle` is an acid or a base, it defines its pka-value. Defaults to None.
2595            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
2596            overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 
2597     
2598        Note:
2599            - For particles with  `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 
2600        deprotonated states, respectively. 
2601            - For particles without an acidity `acidity = None`, only `state_one` is defined.
2602            - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 
2603        """
2604        acidity_valid_keys = ['inert','acidic', 'basic']
2605        if acidity is not None:
2606            if acidity not in acidity_valid_keys:
2607                raise ValueError(f"Acidity {acidity} provided for particle name  {name} is not supproted. Valid keys are: {acidity_valid_keys}")
2608            if acidity in ['acidic', 'basic'] and pka is None:
2609                raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.")   
2610        if acidity == "inert":
2611            acidity = None
2612            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")
2613
2614        self.define_particle_entry_in_df(name=name)
2615        
2616        for index in self.df[self.df['name']==name].index:       
2617            if pka:
2618                self.add_value_to_df(key=('pka',''),
2619                                    index=index,
2620                                    new_value=pka, 
2621                                    verbose=verbose,
2622                                    overwrite=overwrite)
2623           
2624            self.add_value_to_df(key=('acidity',''),
2625                                 index=index,
2626                                 new_value=acidity, 
2627                                 verbose=verbose,
2628                                 overwrite=overwrite) 
2629            if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')):
2630                self.add_value_to_df(key=('state_one','es_type'),
2631                                     index=index,
2632                                     new_value=self.propose_unused_type(), 
2633                                     verbose=verbose,
2634                                     overwrite=overwrite)  
2635            if pd.isna(self.df.loc [self.df['name']  == name].acidity.iloc[0]):
2636                self.add_value_to_df(key=('state_one','z'),
2637                                     index=index,
2638                                     new_value=default_charge_number, 
2639                                     verbose=verbose,
2640                                     overwrite=overwrite)
2641                self.add_value_to_df(key=('state_one','label'),
2642                                     index=index,
2643                                     new_value=name, 
2644                                     verbose=verbose,
2645                                    overwrite=overwrite)
2646            else:
2647                protonated_label = f'{name}H'
2648                self.add_value_to_df(key=('state_one','label'),
2649                                     index=index,
2650                                     new_value=protonated_label, 
2651                                     verbose=verbose,
2652                                    overwrite=overwrite)
2653                self.add_value_to_df(key=('state_two','label'),
2654                                     index=index,
2655                                     new_value=name, 
2656                                     verbose=verbose,
2657                                    overwrite=overwrite)
2658                if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')):
2659                    self.add_value_to_df(key=('state_two','es_type'),
2660                                         index=index,
2661                                         new_value=self.propose_unused_type(), 
2662                                         verbose=verbose,
2663                                         overwrite=overwrite)
2664                if self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'acidic':        
2665                    self.add_value_to_df(key=('state_one','z'),
2666                                         index=index,new_value=0, 
2667                                         verbose=verbose,
2668                                         overwrite=overwrite)
2669                    self.add_value_to_df(key=('state_two','z'),
2670                                         index=index,
2671                                         new_value=-1, 
2672                                         verbose=verbose,
2673                                         overwrite=overwrite)
2674                elif self.df.loc [self.df['name']  == name].acidity.iloc[0] == 'basic':
2675                    self.add_value_to_df(key=('state_one','z'),
2676                                         index=index,new_value=+1, 
2677                                         verbose=verbose,
2678                                         overwrite=overwrite)
2679                    self.add_value_to_df(key=('state_two','z'),
2680                                         index=index,
2681                                         new_value=0, 
2682                                         verbose=verbose,
2683                                         overwrite=overwrite)   
2684        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):
2686    def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True):
2687        """
2688        Sets the set of reduced units used by pyMBE.units and it prints it.
2689
2690        Args:
2691            unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 
2692            unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 
2693            temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 
2694            Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 
2695            verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True.
2696
2697        Note:
2698            - If no `temperature` is given, a value of 298.15 K is assumed by default.
2699            - If no `unit_length` is given, a value of 0.355 nm is assumed by default.
2700            - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 
2701            - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 
2702        """
2703        if unit_length is None:
2704            unit_length= 0.355*self.units.nm
2705        if temperature is None:
2706            temperature = 298.15 * self.units.K
2707        if unit_charge is None:
2708            unit_charge = scipy.constants.e * self.units.C
2709        if Kw is None:
2710            Kw = 1e-14
2711        # Sanity check
2712        variables=[unit_length,temperature,unit_charge]
2713        dimensionalities=["[length]","[temperature]","[charge]"]
2714        for variable,dimensionality in zip(variables,dimensionalities):
2715            self.check_dimensionality(variable,dimensionality)
2716        self.Kw=Kw*self.units.mol**2 / (self.units.l**2)
2717        self.kT=temperature*self.kB
2718        self.units._build_cache()
2719        self.units.define(f'reduced_energy = {self.kT} ')
2720        self.units.define(f'reduced_length = {unit_length}')
2721        self.units.define(f'reduced_charge = {unit_charge}')
2722        if verbose:        
2723            print(self.get_reduced_units())
2724        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):
2726    def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False):
2727        """
2728        Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 
2729
2730        Args:
2731            counter_ion(`str`): `name` of the counter_ion `particle`.
2732            constant_pH(`float`): pH-value.
2733            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
2734            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2735            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
2736
2737        Returns:
2738            RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library.
2739            sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
2740        """
2741        from espressomd import reaction_methods
2742        if exclusion_range is None:
2743            exclusion_range = max(self.get_radius_map().values())*2.0
2744        if pka_set is None:
2745            pka_set=self.get_pka_set()    
2746        self.check_pka_set(pka_set=pka_set)
2747        if use_exclusion_radius_per_type:
2748            exclusion_radius_per_type = self.get_radius_map()
2749        else:
2750            exclusion_radius_per_type = {}
2751        
2752        RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2753                                                    exclusion_range=exclusion_range, 
2754                                                    seed=self.seed, 
2755                                                    constant_pH=constant_pH,
2756                                                    exclusion_radius_per_type = exclusion_radius_per_type
2757                                                    )
2758        sucessfull_reactions_labels=[]
2759        charge_number_map = self.get_charge_number_map()
2760        for name in pka_set.keys():
2761            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
2762                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
2763                continue
2764            gamma=10**-pka_set[name]['pka_value']
2765            state_one_type   = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
2766            state_two_type   = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
2767            counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0]
2768            RE.add_reaction(gamma=gamma,
2769                            reactant_types=[state_one_type],
2770                            product_types=[state_two_type, counter_ion_type],
2771                            default_charges={state_one_type: charge_number_map[state_one_type],
2772                            state_two_type: charge_number_map[state_two_type],
2773                            counter_ion_type: charge_number_map[counter_ion_type]})
2774            sucessfull_reactions_labels.append(name)
2775        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):
2777    def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False):
2778        """
2779        Sets up grand-canonical coupling to a reservoir of salt.
2780        For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead.
2781
2782        Args:
2783            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2784            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
2785            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
2786            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
2787            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
2788            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2789
2790        Returns:
2791            RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library.
2792        """
2793        from espressomd import reaction_methods
2794        if exclusion_range is None:
2795            exclusion_range = max(self.get_radius_map().values())*2.0
2796        if use_exclusion_radius_per_type:
2797            exclusion_radius_per_type = self.get_radius_map()
2798        else:
2799            exclusion_radius_per_type = {}
2800        
2801        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2802                                                    exclusion_range=exclusion_range, 
2803                                                    seed=self.seed, 
2804                                                    exclusion_radius_per_type = exclusion_radius_per_type
2805                                                    )
2806
2807        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
2808        determined_activity_coefficient = activity_coefficient(c_salt_res)
2809        K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient
2810
2811        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
2812        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
2813
2814        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
2815        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
2816
2817        if salt_cation_charge <= 0:
2818            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
2819        if salt_anion_charge >= 0:
2820            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
2821
2822        # Grand-canonical coupling to the reservoir
2823        RE.add_reaction(
2824            gamma = K_salt.magnitude,
2825            reactant_types = [],
2826            reactant_coefficients = [],
2827            product_types = [ salt_cation_es_type, salt_anion_es_type ],
2828            product_coefficients = [ 1, 1 ],
2829            default_charges = {
2830                salt_cation_es_type: salt_cation_charge, 
2831                salt_anion_es_type: salt_anion_charge, 
2832            }
2833        )
2834
2835        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):
2837    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):
2838        """
2839        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
2840        reservoir of small ions. 
2841        This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1].
2842
2843        [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.
2844
2845        Args:
2846            pH_res ('float): pH-value in the reservoir.
2847            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
2848            proton_name ('str'): Name of the proton (H+) particle.
2849            hydroxide_name ('str'): Name of the hydroxide (OH-) particle.
2850            salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle.
2851            salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle.
2852            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
2853            exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted.
2854            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
2855            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`.
2856
2857        Returns:
2858            RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
2859            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
2860            ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients).
2861        """
2862        from espressomd import reaction_methods
2863        if exclusion_range is None:
2864            exclusion_range = max(self.get_radius_map().values())*2.0
2865        if pka_set is None:
2866            pka_set=self.get_pka_set()    
2867        self.check_pka_set(pka_set=pka_set)
2868        if use_exclusion_radius_per_type:
2869            exclusion_radius_per_type = self.get_radius_map()
2870        else:
2871            exclusion_radius_per_type = {}
2872        
2873        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
2874                                                    exclusion_range=exclusion_range, 
2875                                                    seed=self.seed, 
2876                                                    exclusion_radius_per_type = exclusion_radius_per_type
2877                                                    )
2878
2879        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
2880        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
2881        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
2882        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
2883        K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
2884        K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
2885        K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient
2886
2887        proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0]
2888        hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0]     
2889        salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0]
2890        salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0]     
2891
2892        proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0]
2893        hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0]     
2894        salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0]
2895        salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0]     
2896
2897        if proton_charge <= 0:
2898            raise ValueError('ERROR proton charge must be positive, charge ', proton_charge)
2899        if salt_cation_charge <= 0:
2900            raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge)
2901        if hydroxide_charge >= 0:
2902            raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge)
2903        if salt_anion_charge >= 0:
2904            raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge)
2905
2906        # Grand-canonical coupling to the reservoir
2907        # 0 = H+ + OH-
2908        RE.add_reaction(
2909            gamma = K_W.magnitude,
2910            reactant_types = [],
2911            reactant_coefficients = [],
2912            product_types = [ proton_es_type, hydroxide_es_type ],
2913            product_coefficients = [ 1, 1 ],
2914            default_charges = {
2915                proton_es_type: proton_charge, 
2916                hydroxide_es_type: hydroxide_charge, 
2917            }
2918        )
2919
2920        # 0 = Na+ + Cl-
2921        RE.add_reaction(
2922            gamma = K_NACL.magnitude,
2923            reactant_types = [],
2924            reactant_coefficients = [],
2925            product_types = [ salt_cation_es_type, salt_anion_es_type ],
2926            product_coefficients = [ 1, 1 ],
2927            default_charges = {
2928                salt_cation_es_type: salt_cation_charge, 
2929                salt_anion_es_type: salt_anion_charge, 
2930            }
2931        )
2932
2933        # 0 = Na+ + OH-
2934        RE.add_reaction(
2935            gamma = (K_NACL * K_W / K_HCL).magnitude,
2936            reactant_types = [],
2937            reactant_coefficients = [],
2938            product_types = [ salt_cation_es_type, hydroxide_es_type ],
2939            product_coefficients = [ 1, 1 ],
2940            default_charges = {
2941                salt_cation_es_type: salt_cation_charge, 
2942                hydroxide_es_type: hydroxide_charge, 
2943            }
2944        )
2945
2946        # 0 = H+ + Cl-
2947        RE.add_reaction(
2948            gamma = K_HCL.magnitude,
2949            reactant_types = [],
2950            reactant_coefficients = [],
2951            product_types = [ proton_es_type, salt_anion_es_type ],
2952            product_coefficients = [ 1, 1 ],
2953            default_charges = {
2954                proton_es_type: proton_charge, 
2955                salt_anion_es_type: salt_anion_charge, 
2956            }
2957        )
2958
2959        # Annealing moves to ensure sufficient sampling
2960        # Cation annealing H+ = Na+
2961        RE.add_reaction(
2962            gamma = (K_NACL / K_HCL).magnitude,
2963            reactant_types = [proton_es_type],
2964            reactant_coefficients = [ 1 ],
2965            product_types = [ salt_cation_es_type ],
2966            product_coefficients = [ 1 ],
2967            default_charges = {
2968                proton_es_type: proton_charge, 
2969                salt_cation_es_type: salt_cation_charge, 
2970            }
2971        )
2972
2973        # Anion annealing OH- = Cl- 
2974        RE.add_reaction(
2975            gamma = (K_HCL / K_W).magnitude,
2976            reactant_types = [hydroxide_es_type],
2977            reactant_coefficients = [ 1 ],
2978            product_types = [ salt_anion_es_type ],
2979            product_coefficients = [ 1 ],
2980            default_charges = {
2981                hydroxide_es_type: hydroxide_charge, 
2982                salt_anion_es_type: salt_anion_charge, 
2983            }
2984        )
2985
2986        sucessful_reactions_labels=[]
2987        charge_number_map = self.get_charge_number_map()
2988        for name in pka_set.keys():
2989            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
2990                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
2991                continue
2992
2993            Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
2994
2995            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
2996            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
2997
2998            # Reaction in terms of proton: HA = A + H+
2999            RE.add_reaction(gamma=Ka.magnitude,
3000                            reactant_types=[state_one_type],
3001                            reactant_coefficients=[1],
3002                            product_types=[state_two_type, proton_es_type],
3003                            product_coefficients=[1, 1],
3004                            default_charges={state_one_type: charge_number_map[state_one_type],
3005                            state_two_type: charge_number_map[state_two_type],
3006                            proton_es_type: proton_charge})
3007
3008            # Reaction in terms of salt cation: HA = A + Na+
3009            RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude,
3010                            reactant_types=[state_one_type],
3011                            reactant_coefficients=[1],
3012                            product_types=[state_two_type, salt_cation_es_type],
3013                            product_coefficients=[1, 1],
3014                            default_charges={state_one_type: charge_number_map[state_one_type],
3015                            state_two_type: charge_number_map[state_two_type],
3016                            salt_cation_es_type: salt_cation_charge})
3017
3018            # Reaction in terms of hydroxide: OH- + HA = A
3019            RE.add_reaction(gamma=(Ka / K_W).magnitude,
3020                            reactant_types=[state_one_type, hydroxide_es_type],
3021                            reactant_coefficients=[1, 1],
3022                            product_types=[state_two_type],
3023                            product_coefficients=[1],
3024                            default_charges={state_one_type: charge_number_map[state_one_type],
3025                            state_two_type: charge_number_map[state_two_type],
3026                            hydroxide_es_type: hydroxide_charge})
3027
3028            # Reaction in terms of salt anion: Cl- + HA = A
3029            RE.add_reaction(gamma=(Ka / K_HCL).magnitude,
3030                            reactant_types=[state_one_type, salt_anion_es_type],
3031                            reactant_coefficients=[1, 1],
3032                            product_types=[state_two_type],
3033                            product_coefficients=[1],
3034                            default_charges={state_one_type: charge_number_map[state_one_type],
3035                            state_two_type: charge_number_map[state_two_type],
3036                            salt_anion_es_type: salt_anion_charge})
3037
3038            sucessful_reactions_labels.append(name)
3039        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):
3041    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):
3042        """
3043        Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 
3044        reservoir of small ions. 
3045        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-}. 
3046        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'.
3047
3048        [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4).
3049        [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.
3050
3051        Args:
3052            pH_res ('float'): pH-value in the reservoir.
3053            c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir.
3054            cation_name ('str'): Name of the cationic particle.
3055            anion_name ('str'): Name of the anionic particle.
3056            activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
3057            exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted.
3058            pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None.
3059            use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`.
3060
3061        Returns:
3062            RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
3063            sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE.
3064            ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).
3065        """
3066        from espressomd import reaction_methods
3067        if exclusion_range is None:
3068            exclusion_range = max(self.get_radius_map().values())*2.0
3069        if pka_set is None:
3070            pka_set=self.get_pka_set()    
3071        self.check_pka_set(pka_set=pka_set)
3072        if use_exclusion_radius_per_type:
3073            exclusion_radius_per_type = self.get_radius_map()
3074        else:
3075            exclusion_radius_per_type = {}
3076        
3077        RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude,
3078                                                    exclusion_range=exclusion_range, 
3079                                                    seed=self.seed, 
3080                                                    exclusion_radius_per_type = exclusion_radius_per_type
3081                                                    )
3082
3083        # Determine the concentrations of the various species in the reservoir and the equilibrium constants
3084        cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient)
3085        ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res)
3086        determined_activity_coefficient = activity_coefficient(ionic_strength_res)
3087        a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)')
3088        a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3089        a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient)
3090        K_XX = a_cation * a_anion
3091
3092        cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0]
3093        anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0]     
3094        cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0]
3095        anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0]     
3096        if cation_charge <= 0:
3097            raise ValueError('ERROR cation charge must be positive, charge ', cation_charge)
3098        if anion_charge >= 0:
3099            raise ValueError('ERROR anion charge must be negative, charge ', anion_charge)
3100
3101        # Coupling to the reservoir: 0 = X+ + X-
3102        RE.add_reaction(
3103            gamma = K_XX.magnitude,
3104            reactant_types = [],
3105            reactant_coefficients = [],
3106            product_types = [ cation_es_type, anion_es_type ],
3107            product_coefficients = [ 1, 1 ],
3108            default_charges = {
3109                cation_es_type: cation_charge, 
3110                anion_es_type: anion_charge, 
3111            }
3112        )
3113
3114        sucessful_reactions_labels=[]
3115        charge_number_map = self.get_charge_number_map()
3116        for name in pka_set.keys():
3117            if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False :
3118                print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.')
3119                continue
3120
3121            Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l
3122            gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen
3123
3124            state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0]
3125            state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0]
3126
3127            # Reaction in terms of small cation: HA = A + X+
3128            RE.add_reaction(gamma=gamma_K_AX.magnitude,
3129                            reactant_types=[state_one_type],
3130                            reactant_coefficients=[1],
3131                            product_types=[state_two_type, cation_es_type],
3132                            product_coefficients=[1, 1],
3133                            default_charges={state_one_type: charge_number_map[state_one_type],
3134                            state_two_type: charge_number_map[state_two_type],
3135                            cation_es_type: cation_charge})
3136
3137            # Reaction in terms of small anion: X- + HA = A
3138            RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude,
3139                            reactant_types=[state_one_type, anion_es_type],
3140                            reactant_coefficients=[1, 1],
3141                            product_types=[state_two_type],
3142                            product_coefficients=[1],
3143                            default_charges={state_one_type: charge_number_map[state_one_type],
3144                            state_two_type: charge_number_map[state_two_type],
3145                            anion_es_type: anion_charge})
3146
3147            sucessful_reactions_labels.append(name)
3148        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):
3150    def setup_df (self):
3151        """
3152        Sets up the pyMBE's dataframe `pymbe.df`.
3153
3154        Returns:
3155            columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe
3156        """
3157        
3158        columns_dtypes = {
3159            'name': {
3160                '': str},
3161            'pmb_type': {
3162                '': str},
3163            'particle_id': {
3164                '': object},
3165            'particle_id2':  {
3166                '': object},
3167            'residue_id':  {
3168                '': object},
3169            'molecule_id':  {
3170                '': object},
3171            'acidity':  {
3172                '': str},
3173            'pka':  {
3174                '': float},
3175            'central_bead':  {
3176                '': object},
3177            'side_chains': {
3178                '': object},
3179            'residue_list': {
3180                '': object},
3181            'model': {
3182                '': str},
3183            'sigma': {
3184                '': object},
3185            'cutoff': {
3186                '': object},
3187            'offset': {
3188                '': object},
3189            'epsilon': {
3190                '': object},
3191            'state_one': {
3192                'label': str,
3193                'es_type': object,
3194                'z': object },
3195            'state_two': {
3196                'label': str,
3197                'es_type': object,
3198                'z': object },
3199            'sequence': {
3200                '': object},
3201            'bond_object': {
3202                '': object},
3203            'parameters_of_the_potential':{
3204                '': object},
3205            'l0': {
3206                '': float},
3207            }
3208        
3209        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()]))
3210        
3211        for level1, sub_dtypes in columns_dtypes.items():
3212            for level2, dtype in sub_dtypes.items():
3213                self.df[level1, level2] = self.df[level1, level2].astype(dtype)
3214                
3215        columns_names = pd.MultiIndex.from_frame(self.df)
3216        columns_names = columns_names.names
3217                
3218        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):
3220    def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True):
3221        """
3222        Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`.
3223
3224        Args:
3225            espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library.
3226            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.
3227            combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'.
3228            warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True.
3229
3230        Note:
3231            - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 
3232            - Currently, the only `combining_rule` supported is Lorentz-Berthelot.
3233            - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3234
3235        """
3236        from itertools import combinations_with_replacement
3237        implemented_combining_rules = ['Lorentz-Berthelot']
3238        compulsory_parameters_in_df = ['sigma','epsilon']
3239        # Sanity check
3240        if combining_rule not in implemented_combining_rules:
3241            raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules)
3242        if shift_potential:
3243            shift="auto"
3244        else:
3245            shift=0
3246        # List which particles have sigma and epsilon values defined in pmb.df and which ones don't
3247        particles_types_with_LJ_parameters = []
3248        non_parametrized_labels= []
3249        for particle_type in self.get_type_map().values():
3250            check_list=[]
3251            for key in compulsory_parameters_in_df:
3252                value_in_df=self.find_value_from_es_type(es_type=particle_type,
3253                                                        column_name=key)
3254                check_list.append(np.isnan(value_in_df))
3255            if any(check_list):
3256                non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 
3257                                                                            column_name='label'))
3258            else:
3259                particles_types_with_LJ_parameters.append(particle_type)
3260        # Set up LJ interactions between all particle types
3261        for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 
3262            particle_name1 = self.find_value_from_es_type(es_type=type_pair[0],
3263                                                        column_name="name")
3264            particle_name2 = self.find_value_from_es_type(es_type=type_pair[1],
3265                                                        column_name="name")
3266            lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1,
3267                                                 particle_name2 = particle_name2,
3268                                                 combining_rule = combining_rule)
3269            
3270            # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others    
3271            if not lj_parameters:
3272                continue
3273            espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 
3274                                                                                    sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 
3275                                                                                    cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude,
3276                                                                                    offset = lj_parameters["offset"].to("reduced_length").magnitude, 
3277                                                                                    shift = shift)                                                                                          
3278            index = len(self.df)
3279            label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label")
3280            label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label")
3281            self.df.at [index, 'name'] = f'LJ: {label1}-{label2}'
3282            lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params()
3283
3284            self.add_value_to_df(index=index,
3285                                key=('pmb_type',''),
3286                                new_value='LennardJones')
3287            
3288            self.add_value_to_df(index=index,
3289                                key=('parameters_of_the_potential',''),
3290                                new_value=lj_params,
3291                                non_standard_value=True)
3292        if non_parametrized_labels and warnings:
3293            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.')
3294            
3295        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):
3297    def write_pmb_df (self, filename):
3298        '''
3299        Writes the pyMBE dataframe into a csv file
3300        
3301        Args:
3302            filename(`str`): Path to the csv file 
3303        '''
3304
3305        columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence']
3306        df = self.df.copy(deep=True)
3307        for column_name in columns_with_list_or_dict:
3308            df[column_name] = df[column_name].apply(lambda x: json.dumps(x) if isinstance(x, (np.ndarray, tuple, list, dict)) or pd.notnull(x) else x)
3309        df['bond_object'] = df['bond_object'].apply(lambda x: f'{x.__class__.__name__}({json.dumps(x.get_params())})' if pd.notnull(x) else x)
3310        df.to_csv(filename)
3311
3312        return

Writes the pyMBE dataframe into a csv file

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

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

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

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

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

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