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