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 28 29class pymbe_library(): 30 """ 31 The library for the Molecular Builder for ESPResSo (pyMBE) 32 33 Attributes: 34 N_A(`pint.Quantity`): Avogadro number. 35 Kb(`pint.Quantity`): Boltzmann constant. 36 e(`pint.Quantity`): Elementary charge. 37 df(`Pandas.Dataframe`): Dataframe used to bookkeep all the information stored in pyMBE. Typically refered as `pmb.df`. 38 kT(`pint.Quantity`): Thermal energy. 39 Kw(`pint.Quantity`): Ionic product of water. Used in the setup of the G-RxMC method. 40 """ 41 42 class NumpyEncoder(json.JSONEncoder): 43 """ 44 Custom JSON encoder that converts NumPy arrays to Python lists 45 and NumPy scalars to Python scalars. 46 """ 47 def default(self, obj): 48 if isinstance(obj, np.ndarray): 49 return obj.tolist() 50 if isinstance(obj, np.generic): 51 return obj.item() 52 return super().default(obj) 53 54 def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None): 55 """ 56 Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 57 and sets up the `pmb.df` for bookkeeping. 58 59 Args: 60 temperature(`pint.Quantity`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. 61 unit_length(`pint.Quantity`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. 62 unit_charge (`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 63 Kw (`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 64 65 Note: 66 - If no `temperature` is given, a value of 298.15 K is assumed by default. 67 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 68 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 69 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 70 """ 71 # Seed and RNG 72 self.seed=seed 73 self.rng = np.random.default_rng(seed) 74 self.units=pint.UnitRegistry() 75 self.N_A=scipy.constants.N_A / self.units.mol 76 self.kB=scipy.constants.k * self.units.J / self.units.K 77 self.e=scipy.constants.e * self.units.C 78 self.set_reduced_units(unit_length=unit_length, 79 unit_charge=unit_charge, 80 temperature=temperature, 81 Kw=Kw) 82 self.setup_df() 83 self.lattice_builder = None 84 self.root = importlib.resources.files(__package__) 85 return 86 87 def _check_if_name_is_defined_in_df(self, name): 88 """ 89 Checks if `name` is defined in `pmb.df`. 90 91 Args: 92 name(`str`): label to check if defined in `pmb.df`. 93 94 Returns: 95 `bool`: `True` for success, `False` otherwise. 96 """ 97 return name in self.df['name'].unique() 98 99 def _check_if_multiple_pmb_types_for_name(self, name, pmb_type_to_be_defined): 100 """ 101 Checks if `name` is defined in `pmb.df` with multiple pmb_types. 102 103 Args: 104 name(`str`): label to check if defined in `pmb.df`. 105 pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`. 106 107 Returns: 108 `bool`: `True` for success, `False` otherwise. 109 """ 110 if name in self.df['name'].unique(): 111 current_object_type = self.df[self.df['name']==name].pmb_type.values[0] 112 if current_object_type != pmb_type_to_be_defined: 113 raise ValueError (f"The name {name} is already defined in the df with a pmb_type = {current_object_type}, pymMBE does not support objects with the same name but different pmb_types") 114 115 116 def _check_supported_molecule(self, molecule_name,valid_pmb_types): 117 """ 118 Checks if the molecule name `molecule_name` is supported by a method of pyMBE. 119 120 Args: 121 molecule_name(`str`): pmb object type to be checked. 122 valid_pmb_types(`list` of `str`): List of valid pmb types supported by the method. 123 124 Returns: 125 pmb_type(`str`): pmb_type of the molecule. 126 """ 127 pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0] 128 if pmb_type not in valid_pmb_types: 129 raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}") 130 return pmb_type 131 132 def _check_if_name_has_right_type(self, name, expected_pmb_type, hard_check=True): 133 """ 134 Checks if `name` is of the expected pmb type. 135 136 Args: 137 name(`str`): label to check if defined in `pmb.df`. 138 expected_pmb_type(`str`): pmb object type corresponding to `name`. 139 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`. 140 141 Returns: 142 `bool`: `True` for success, `False` otherwise. 143 """ 144 pmb_type=self.df.loc[self.df['name']==name].pmb_type.values[0] 145 if pmb_type == expected_pmb_type: 146 return True 147 else: 148 if hard_check: 149 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}") 150 return False 151 152 def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False): 153 """ 154 Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles. 155 156 Args: 157 particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles 158 particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles 159 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False. 160 161 Returns: 162 index(`int`): Row index where the bond information has been added in pmb.df. 163 """ 164 particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0] 165 particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0] 166 167 bond_key = self.find_bond_key(particle_name1=particle_name1, 168 particle_name2=particle_name2, 169 use_default_bond=use_default_bond) 170 if not bond_key: 171 return None 172 self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1) 173 indexs = np.where(self.df['name']==bond_key) 174 index_list = list (indexs[0]) 175 used_bond_df = self.df.loc[self.df['particle_id2'].notnull()] 176 #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict 177 used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 ) 178 used_bond_index = used_bond_df.index.to_list() 179 if not index_list: 180 return None 181 for index in index_list: 182 if index not in used_bond_index: 183 self.clean_df_row(index=int(index)) 184 self.df.at[index,'particle_id'] = particle_id1 185 self.df.at[index,'particle_id2'] = particle_id2 186 break 187 return index 188 189 def add_bonds_to_espresso(self, espresso_system) : 190 """ 191 Adds all bonds defined in `pmb.df` to `espresso_system`. 192 193 Args: 194 espresso_system(`espressomd.system.System`): system object of espressomd library 195 """ 196 197 if 'bond' in self.df["pmb_type"].values: 198 bond_df = self.df.loc[self.df ['pmb_type'] == 'bond'] 199 bond_list = bond_df.bond_object.values.tolist() 200 for bond in bond_list: 201 espresso_system.bonded_inter.add(bond) 202 else: 203 logging.warning('there are no bonds defined in pymbe.df') 204 return 205 206 def add_value_to_df(self,index,key,new_value, non_standard_value=False, overwrite=False): 207 """ 208 Adds a value to a cell in the `pmb.df` DataFrame. 209 210 Args: 211 index(`int`): index of the row to add the value to. 212 key(`str`): the column label to add the value to. 213 non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False. 214 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 215 """ 216 217 token = "#protected:" 218 219 def protect(obj): 220 if non_standard_value: 221 return token + json.dumps(obj, cls=self.NumpyEncoder) 222 return obj 223 224 def deprotect(obj): 225 if non_standard_value and isinstance(obj, str) and obj.startswith(token): 226 return json.loads(obj.removeprefix(token)) 227 return obj 228 229 # Make sure index is a scalar integer value 230 index = int(index) 231 assert isinstance(index, int), '`index` should be a scalar integer value.' 232 idx = pd.IndexSlice 233 if self.check_if_df_cell_has_a_value(index=index,key=key): 234 old_value = self.df.loc[index,idx[key]] 235 if not pd.Series([protect(old_value)]).equals(pd.Series([protect(new_value)])): 236 name=self.df.loc[index,('name','')] 237 pmb_type=self.df.loc[index,('pmb_type','')] 238 logging.debug(f"You are attempting to redefine the properties of {name} of pmb_type {pmb_type}") 239 if overwrite: 240 logging.info(f'Overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}') 241 if not overwrite: 242 logging.debug(f"pyMBE has preserved of the entry `{key}`: old_value = {old_value}. If you want to overwrite it with new_value = {new_value}, activate the switch overwrite = True ") 243 return 244 245 self.df.loc[index,idx[key]] = protect(new_value) 246 if non_standard_value: 247 self.df[key] = self.df[key].apply(deprotect) 248 return 249 250 def assign_molecule_id(self, molecule_index): 251 """ 252 Assigns the `molecule_id` of the pmb object given by `pmb_type` 253 254 Args: 255 molecule_index(`int`): index of the current `pmb_object_type` to assign the `molecule_id` 256 Returns: 257 molecule_id(`int`): Id of the molecule 258 """ 259 260 self.clean_df_row(index=int(molecule_index)) 261 262 if self.df['molecule_id'].isnull().values.all(): 263 molecule_id = 0 264 else: 265 molecule_id = self.df['molecule_id'].max() +1 266 267 self.add_value_to_df (key=('molecule_id',''), 268 index=int(molecule_index), 269 new_value=molecule_id) 270 271 return molecule_id 272 273 def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system): 274 """ 275 Calculates the center of the molecule with a given molecule_id. 276 277 Args: 278 molecule_id(`int`): Id of the molecule whose center of mass is to be calculated. 279 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 280 281 Returns: 282 center_of_mass(`lst`): Coordinates of the center of mass. 283 """ 284 center_of_mass = np.zeros(3) 285 axis_list = [0,1,2] 286 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0] 287 particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"] 288 for pid in particle_id_list: 289 for axis in axis_list: 290 center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis] 291 center_of_mass = center_of_mass / len(particle_id_list) 292 return center_of_mass 293 294 def calculate_HH(self, molecule_name, pH_list=None, pka_set=None): 295 """ 296 Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 297 for molecules with the name `molecule_name`. 298 299 Args: 300 molecule_name(`str`): name of the molecule to calculate the ideal charge for 301 pH_list(`lst`): pH-values to calculate. 302 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 303 304 Returns: 305 Z_HH(`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list` 306 307 Note: 308 - This function supports objects with pmb types: "molecule", "peptide" and "protein". 309 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 310 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 311 - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan` should be used. 312 """ 313 self._check_if_name_is_defined_in_df(name=molecule_name) 314 self._check_supported_molecule(molecule_name=molecule_name, 315 valid_pmb_types=["molecule","peptide","protein"]) 316 if pH_list is None: 317 pH_list=np.linspace(2,12,50) 318 if pka_set is None: 319 pka_set=self.get_pka_set() 320 index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 321 residue_list = self.df.at [index,('residue_list','')].copy() 322 particles_in_molecule = [] 323 for residue in residue_list: 324 list_of_particles_in_residue = self.search_particles_in_residue(residue) 325 if len(list_of_particles_in_residue) == 0: 326 logging.warning(f"The residue {residue} has no particles defined in the pyMBE DataFrame, it will be ignored.") 327 continue 328 particles_in_molecule += list_of_particles_in_residue 329 if len(particles_in_molecule) == 0: 330 return [None]*len(pH_list) 331 self.check_pka_set(pka_set=pka_set) 332 charge_number_map = self.get_charge_number_map() 333 Z_HH=[] 334 for pH_value in pH_list: 335 Z=0 336 for particle in particles_in_molecule: 337 if particle in pka_set.keys(): 338 if pka_set[particle]['acidity'] == 'acidic': 339 psi=-1 340 elif pka_set[particle]['acidity']== 'basic': 341 psi=+1 342 Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value']))) 343 else: 344 state_one_type = self.df.loc[self.df['name']==particle].state_one.es_type.values[0] 345 Z+=charge_number_map[state_one_type] 346 Z_HH.append(Z) 347 return Z_HH 348 349 def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None): 350 """ 351 Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning. 352 353 Args: 354 c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 355 c_salt('float'): Salt concentration in the reservoir. 356 pH_list('lst'): List of pH-values in the reservoir. 357 pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}. 358 359 Returns: 360 {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 361 pH_system_list ('lst'): List of pH_values in the system. 362 partition_coefficients_list ('lst'): List of partition coefficients of cations. 363 364 Note: 365 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 366 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 367 - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result. 368 """ 369 if pH_list is None: 370 pH_list=np.linspace(2,12,50) 371 if pka_set is None: 372 pka_set=self.get_pka_set() 373 self.check_pka_set(pka_set=pka_set) 374 375 partition_coefficients_list = [] 376 pH_system_list = [] 377 Z_HH_Donnan={} 378 for key in c_macro: 379 Z_HH_Donnan[key] = [] 380 381 def calc_charges(c_macro, pH): 382 """ 383 Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation. 384 385 Args: 386 c_macro('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 387 pH('float'): pH-value that is used in the HH equation. 388 389 Returns: 390 charge('dict'): {"molecule_name": charge} 391 """ 392 charge = {} 393 for name in c_macro: 394 charge[name] = self.calculate_HH(name, [pH], pka_set)[0] 395 return charge 396 397 def calc_partition_coefficient(charge, c_macro): 398 """ 399 Calculates the partition coefficients of positive ions according to the ideal Donnan theory. 400 401 Args: 402 charge('dict'): {"molecule_name": charge} 403 c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 404 """ 405 nonlocal ionic_strength_res 406 charge_density = 0.0 407 for key in charge: 408 charge_density += charge[key] * c_macro[key] 409 return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude 410 411 for pH_value in pH_list: 412 # calculate the ionic strength of the reservoir 413 if pH_value <= 7.0: 414 ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 415 elif pH_value > 7.0: 416 ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt 417 418 #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations 419 #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation. 420 #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects. 421 equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro)) 422 logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq") 423 partition_coefficient = 10**logxi.root 424 425 charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient)) 426 for key in c_macro: 427 Z_HH_Donnan[key].append(charges_temp[key]) 428 429 pH_system_list.append(pH_value - np.log10(partition_coefficient)) 430 partition_coefficients_list.append(partition_coefficient) 431 432 return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 433 434 def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset): 435 """ 436 Calculates the initial bond length that is used when setting up molecules, 437 based on the minimum of the sum of bonded and short-range (LJ) interactions. 438 439 Args: 440 bond_object(`espressomd.interactions.BondedInteractions`): instance of a bond object from espressomd library 441 bond_type(`str`): label identifying the used bonded potential 442 epsilon(`pint.Quantity`): LJ epsilon of the interaction between the particles 443 sigma(`pint.Quantity`): LJ sigma of the interaction between the particles 444 cutoff(`pint.Quantity`): cutoff-radius of the LJ interaction 445 offset(`pint.Quantity`): offset of the LJ interaction 446 """ 447 def truncated_lj_potential(x, epsilon, sigma, cutoff,offset): 448 if x>cutoff: 449 return 0.0 450 else: 451 return 4*epsilon*((sigma/(x-offset))**12-(sigma/(x-offset))**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6) 452 453 epsilon_red=epsilon.to('reduced_energy').magnitude 454 sigma_red=sigma.to('reduced_length').magnitude 455 cutoff_red=cutoff.to('reduced_length').magnitude 456 offset_red=offset.to('reduced_length').magnitude 457 458 if bond_type == "harmonic": 459 r_0 = bond_object.params.get('r_0') 460 k = bond_object.params.get('k') 461 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 462 elif bond_type == "FENE": 463 r_0 = bond_object.params.get('r_0') 464 k = bond_object.params.get('k') 465 d_r_max = bond_object.params.get('d_r_max') 466 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 467 return l0 468 469 def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False): 470 ''' 471 Calculates the net charge per molecule of molecules with `name` = molecule_name. 472 Returns the net charge per molecule and a maps with the net charge per residue and molecule. 473 474 Args: 475 espresso_system(`espressomd.system.System`): system information 476 molecule_name(`str`): name of the molecule to calculate the net charge 477 dimensionless(`bool'): sets if the charge is returned with a dimension or not 478 479 Returns: 480 (`dict`): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }} 481 482 Note: 483 - The net charge of the molecule is averaged over all molecules of type `name` 484 - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name` 485 ''' 486 self._check_supported_molecule(molecule_name=molecule_name, 487 valid_pmb_types=["molecule","protein","peptide"]) 488 489 id_map = self.get_particle_id_map(object_name=molecule_name) 490 def create_charge_map(espresso_system,id_map,label): 491 charge_number_map={} 492 for super_id in id_map[label].keys(): 493 if dimensionless: 494 net_charge=0 495 else: 496 net_charge=0 * self.units.Quantity(1,'reduced_charge') 497 for pid in id_map[label][super_id]: 498 if dimensionless: 499 net_charge+=espresso_system.part.by_id(pid).q 500 else: 501 net_charge+=espresso_system.part.by_id(pid).q * self.units.Quantity(1,'reduced_charge') 502 charge_number_map[super_id]=net_charge 503 return charge_number_map 504 net_charge_molecules=create_charge_map(label="molecule_map", 505 espresso_system=espresso_system, 506 id_map=id_map) 507 net_charge_residues=create_charge_map(label="residue_map", 508 espresso_system=espresso_system, 509 id_map=id_map) 510 if dimensionless: 511 mean_charge=np.mean(np.array(list(net_charge_molecules.values()))) 512 else: 513 mean_charge=np.mean(np.array([value.magnitude for value in net_charge_molecules.values()]))*self.units.Quantity(1,'reduced_charge') 514 return {"mean": mean_charge, "molecules": net_charge_molecules, "residues": net_charge_residues} 515 516 def center_molecule_in_simulation_box(self, molecule_id, espresso_system): 517 """ 518 Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`. 519 520 Args: 521 molecule_id(`int`): Id of the molecule to be centered. 522 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 523 """ 524 if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0: 525 raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.") 526 center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system) 527 box_center = [espresso_system.box_l[0]/2.0, 528 espresso_system.box_l[1]/2.0, 529 espresso_system.box_l[2]/2.0] 530 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0] 531 particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"] 532 for pid in particle_id_list: 533 es_pos = espresso_system.part.by_id(pid).pos 534 espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center 535 return 536 537 def check_aminoacid_key(self, key): 538 """ 539 Checks if `key` corresponds to a valid aminoacid letter code. 540 541 Args: 542 key(`str`): key to be checked. 543 544 Returns: 545 `bool`: True if `key` is a valid aminoacid letter code, False otherwise. 546 """ 547 valid_AA_keys=['V', #'VAL' 548 'I', #'ILE' 549 'L', #'LEU' 550 'E', #'GLU' 551 'Q', #'GLN' 552 'D', #'ASP' 553 'N', #'ASN' 554 'H', #'HIS' 555 'W', #'TRP' 556 'F', #'PHE' 557 'Y', #'TYR' 558 'R', #'ARG' 559 'K', #'LYS' 560 'S', #'SER' 561 'T', #'THR' 562 'M', #'MET' 563 'A', #'ALA' 564 'G', #'GLY' 565 'P', #'PRO' 566 'C'] #'CYS' 567 if key in valid_AA_keys: 568 return True 569 else: 570 return False 571 572 def check_dimensionality(self, variable, expected_dimensionality): 573 """ 574 Checks if the dimensionality of `variable` matches `expected_dimensionality`. 575 576 Args: 577 variable(`pint.Quantity`): Quantity to be checked. 578 expected_dimensionality(`str`): Expected dimension of the variable. 579 580 Returns: 581 (`bool`): `True` if the variable if of the expected dimensionality, `False` otherwise. 582 583 Note: 584 - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality). 585 - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"` 586 """ 587 correct_dimensionality=variable.check(f"{expected_dimensionality}") 588 if not correct_dimensionality: 589 raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}") 590 return correct_dimensionality 591 592 def check_if_df_cell_has_a_value(self, index,key): 593 """ 594 Checks if a cell in the `pmb.df` at the specified index and column has a value. 595 596 Args: 597 index(`int`): Index of the row to check. 598 key(`str`): Column label to check. 599 600 Returns: 601 `bool`: `True` if the cell has a value, `False` otherwise. 602 """ 603 idx = pd.IndexSlice 604 import warnings 605 with warnings.catch_warnings(): 606 warnings.simplefilter("ignore") 607 return not pd.isna(self.df.loc[index, idx[key]]) 608 609 610 611 def check_if_metal_ion(self,key): 612 """ 613 Checks if `key` corresponds to a label of a supported metal ion. 614 615 Args: 616 key(`str`): key to be checked 617 618 Returns: 619 (`bool`): True if `key` is a supported metal ion, False otherwise. 620 """ 621 if key in self.get_metal_ions_charge_number_map().keys(): 622 return True 623 else: 624 return False 625 626 def check_pka_set(self, pka_set): 627 """ 628 Checks that `pka_set` has the formatting expected by the pyMBE library. 629 630 Args: 631 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 632 """ 633 required_keys=['pka_value','acidity'] 634 for required_key in required_keys: 635 for pka_name, pka_entry in pka_set.items(): 636 if required_key not in pka_entry: 637 raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")') 638 return 639 640 def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")): 641 """ 642 Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a pd.NA value. 643 644 Args: 645 index(`int`): Index of the row to clean. 646 columns_keys_to_clean(`list` of `str`, optional): List with the column keys to be cleaned. Defaults to [`particle_id`, `particle_id2`, `residue_id`, `molecule_id`]. 647 """ 648 for column_key in columns_keys_to_clean: 649 self.add_value_to_df(key=(column_key,''),index=index,new_value=pd.NA) 650 self.df.fillna(pd.NA, inplace=True) 651 return 652 653 def convert_columns_to_original_format (self, df): 654 """ 655 Converts the columns of the Dataframe to the original format in pyMBE. 656 657 Args: 658 df(`DataFrame`): dataframe with pyMBE information as a string 659 660 """ 661 662 columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', ('state_one','es_type'),('state_two','es_type'),('state_one','z'),('state_two','z') ] 663 664 columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset'] 665 666 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence', 'chain_map', 'node_map'] 667 668 for column_name in columns_dtype_int: 669 df[column_name] = df[column_name].astype(pd.Int64Dtype()) 670 671 for column_name in columns_with_list_or_dict: 672 if df[column_name].isnull().all(): 673 df[column_name] = df[column_name].astype(object) 674 else: 675 df[column_name] = df[column_name].apply(lambda x: json.loads(x) if pd.notnull(x) else x) 676 677 for column_name in columns_with_units: 678 df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x) 679 680 df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x) 681 df["l0"] = df["l0"].astype(object) 682 df["pka"] = df["pka"].astype(object) 683 return df 684 685 def convert_str_to_bond_object (self, bond_str): 686 """ 687 Convert a row read as a `str` to the corresponding ESPResSo bond object. 688 689 Args: 690 bond_str(`str`): string with the information of a bond object. 691 692 Returns: 693 bond_object(`obj`): ESPResSo bond object. 694 695 Note: 696 Current supported bonds are: HarmonicBond and FeneBond 697 """ 698 import espressomd.interactions 699 700 supported_bonds = ['HarmonicBond', 'FeneBond'] 701 m = re.search(r'^([A-Za-z0-9_]+)\((\{.+\})\)$', bond_str) 702 if m is None: 703 raise ValueError(f'Cannot parse bond "{bond_str}"') 704 bond = m.group(1) 705 if bond not in supported_bonds: 706 raise NotImplementedError(f"Bond type '{bond}' currently not implemented in pyMBE, accepted types are {supported_bonds}") 707 params = json.loads(m.group(2)) 708 bond_id = params.pop("bond_id") 709 bond_object = getattr(espressomd.interactions, bond)(**params) 710 bond_object._bond_id = bond_id 711 return bond_object 712 713 def copy_df_entry(self, name, column_name, number_of_copies): 714 ''' 715 Creates 'number_of_copies' of a given 'name' in `pymbe.df`. 716 717 Args: 718 name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df` 719 column_name(`str`): Column name to use as a filter. 720 number_of_copies(`int`): number of copies of `name` to be created. 721 722 Note: 723 - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 724 ''' 725 726 valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ] 727 if column_name not in valid_column_names: 728 raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}") 729 df_by_name = self.df.loc[self.df.name == name] 730 if number_of_copies != 1: 731 if df_by_name[column_name].isnull().values.any(): 732 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True) 733 else: 734 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 735 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 736 df_by_name_repeated[column_name] = pd.NA 737 # Concatenate the new particle rows to `df` 738 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 739 else: 740 if not df_by_name[column_name].isnull().values.any(): 741 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 742 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 743 df_by_name_repeated[column_name] = pd.NA 744 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 745 return 746 747 def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt): 748 """ 749 Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`. 750 751 Args: 752 espresso_system(`espressomd.system.System`): instance of an espresso system object. 753 cation_name(`str`): `name` of a particle with a positive charge. 754 anion_name(`str`): `name` of a particle with a negative charge. 755 c_salt(`float`): Salt concentration. 756 757 Returns: 758 c_salt_calculated(`float`): Calculated salt concentration added to `espresso_system`. 759 """ 760 for name in [cation_name, anion_name]: 761 if not self._check_if_name_is_defined_in_df(name=name): 762 logging.warning(f"Object with name '{name}' is not defined in the DataFrame, no ions will be created.") 763 return 764 self._check_if_name_has_right_type(name=cation_name, 765 expected_pmb_type="particle") 766 self._check_if_name_has_right_type(name=anion_name, 767 expected_pmb_type="particle") 768 cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0] 769 anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0] 770 if cation_name_charge <= 0: 771 raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge) 772 if anion_name_charge >= 0: 773 raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge) 774 # Calculate the number of ions in the simulation box 775 volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3') 776 if c_salt.check('[substance] [length]**-3'): 777 N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude) 778 c_salt_calculated=N_ions/(volume*self.N_A) 779 elif c_salt.check('[length]**-3'): 780 N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude) 781 c_salt_calculated=N_ions/volume 782 else: 783 raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt) 784 N_cation = N_ions*abs(anion_name_charge) 785 N_anion = N_ions*abs(cation_name_charge) 786 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation) 787 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion) 788 if c_salt_calculated.check('[substance] [length]**-3'): 789 logging.info(f"added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions") 790 elif c_salt_calculated.check('[length]**-3'): 791 logging.info(f"added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions") 792 return c_salt_calculated 793 794 def create_bond_in_espresso(self, bond_type, bond_parameters): 795 ''' 796 Creates either a harmonic or a FENE bond in ESPResSo 797 798 Args: 799 bond_type(`str`): label to identify the potential to model the bond. 800 bond_parameters(`dict`): parameters of the potential of the bond. 801 802 Note: 803 Currently, only HARMONIC and FENE bonds are supported. 804 805 For a HARMONIC bond the dictionary must contain: 806 807 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 808 using the `pmb.units` UnitRegistry. 809 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 810 the `pmb.units` UnitRegistry. 811 812 For a FENE bond the dictionary must additionally contain: 813 814 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 815 units of length using the `pmb.units` UnitRegistry. Default 'None'. 816 817 Returns: 818 bond_object (`obj`): an ESPResSo bond object 819 ''' 820 from espressomd import interactions 821 822 valid_bond_types = ["harmonic", "FENE"] 823 824 if 'k' in bond_parameters: 825 bond_magnitude = bond_parameters['k'].to('reduced_energy / reduced_length**2') 826 else: 827 raise ValueError("Magnitude of the potential (k) is missing") 828 829 if bond_type == 'harmonic': 830 if 'r_0' in bond_parameters: 831 bond_length = bond_parameters['r_0'].to('reduced_length') 832 else: 833 raise ValueError("Equilibrium bond length (r_0) is missing") 834 bond_object = interactions.HarmonicBond(k = bond_magnitude.magnitude, 835 r_0 = bond_length.magnitude) 836 elif bond_type == 'FENE': 837 if 'r_0' in bond_parameters: 838 bond_length = bond_parameters['r_0'].to('reduced_length').magnitude 839 else: 840 logging.warning("no value provided for r_0. Defaulting to r_0 = 0") 841 bond_length=0 842 if 'd_r_max' in bond_parameters: 843 max_bond_stret = bond_parameters['d_r_max'].to('reduced_length') 844 else: 845 raise ValueError("Maximal stretching length (d_r_max) is missing") 846 bond_object = interactions.FeneBond(r_0 = bond_length, 847 k = bond_magnitude.magnitude, 848 d_r_max = max_bond_stret.magnitude) 849 else: 850 raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}") 851 return bond_object 852 853 854 def create_counterions(self, object_name, cation_name, anion_name, espresso_system): 855 """ 856 Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`. 857 858 Args: 859 object_name(`str`): `name` of a pymbe object. 860 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 861 cation_name(`str`): `name` of a particle with a positive charge. 862 anion_name(`str`): `name` of a particle with a negative charge. 863 864 Returns: 865 counterion_number(`dict`): {"name": number} 866 867 Note: 868 This function currently does not support the creation of counterions for hydrogels. 869 """ 870 for name in [object_name, cation_name, anion_name]: 871 if not self._check_if_name_is_defined_in_df(name=name): 872 logging.warning(f"Object with name '{name}' is not defined in the DataFrame, no counterions will be created.") 873 return 874 for name in [cation_name, anion_name]: 875 self._check_if_name_has_right_type(name=name, expected_pmb_type="particle") 876 self._check_supported_molecule(molecule_name=object_name, 877 valid_pmb_types=["molecule","peptide","protein"]) 878 879 880 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0] 881 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0] 882 object_ids = self.get_particle_id_map(object_name=object_name)["all"] 883 counterion_number={} 884 object_charge={} 885 for name in ['positive', 'negative']: 886 object_charge[name]=0 887 for id in object_ids: 888 if espresso_system.part.by_id(id).q > 0: 889 object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q )) 890 elif espresso_system.part.by_id(id).q < 0: 891 object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q )) 892 if object_charge['positive'] % abs(anion_charge) == 0: 893 counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge)) 894 else: 895 raise ValueError('The number of positive charges in the pmb_object must be divisible by the charge of the anion') 896 if object_charge['negative'] % abs(cation_charge) == 0: 897 counterion_number[cation_name]=int(object_charge['negative']/cation_charge) 898 else: 899 raise ValueError('The number of negative charges in the pmb_object must be divisible by the charge of the cation') 900 if counterion_number[cation_name] > 0: 901 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name]) 902 else: 903 counterion_number[cation_name]=0 904 if counterion_number[anion_name] > 0: 905 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name]) 906 else: 907 counterion_number[anion_name] = 0 908 logging.info('the following counter-ions have been created: ') 909 for name in counterion_number.keys(): 910 logging.info(f'Ion type: {name} created number: {counterion_number[name]}') 911 return counterion_number 912 913 def create_hydrogel(self, name, espresso_system): 914 """ 915 creates the hydrogel `name` in espresso_system 916 Args: 917 name(`str`): Label of the hydrogel to be created. `name` must be defined in the `pmb.df` 918 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 919 920 Returns: 921 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],...}} 922 """ 923 if not self._check_if_name_is_defined_in_df(name=name): 924 logging.warning(f"Hydrogel with name '{name}' is not defined in the DataFrame, no hydrogel will be created.") 925 return 926 self._check_if_name_has_right_type(name=name, 927 expected_pmb_type="hydrogel") 928 hydrogel_info={"name":name, "chains":{}, "nodes":{}} 929 # placing nodes 930 node_positions = {} 931 node_topology = self.df[self.df["name"]==name]["node_map"].iloc[0] 932 for node_info in node_topology: 933 node_index = node_info["lattice_index"] 934 node_name = node_info["particle_name"] 935 node_pos, node_id = self.create_hydrogel_node(self.format_node(node_index), node_name, espresso_system) 936 hydrogel_info["nodes"][self.format_node(node_index)]=node_id 937 node_positions[node_id[0]]=node_pos 938 939 # Placing chains between nodes 940 # Looping over all the 16 chains 941 chain_topology = self.df[self.df["name"]==name]["chain_map"].iloc[0] 942 for chain_info in chain_topology: 943 node_s = chain_info["node_start"] 944 node_e = chain_info["node_end"] 945 molecule_info = self.create_hydrogel_chain(node_s, node_e, node_positions, espresso_system) 946 for molecule_id in molecule_info: 947 hydrogel_info["chains"][molecule_id] = molecule_info[molecule_id] 948 hydrogel_info["chains"][molecule_id]["node_start"]=node_s 949 hydrogel_info["chains"][molecule_id]["node_end"]=node_e 950 return hydrogel_info 951 952 def create_hydrogel_chain(self, node_start, node_end, node_positions, espresso_system): 953 """ 954 Creates a chain between two nodes of a hydrogel. 955 956 Args: 957 node_start(`str`): name of the starting node particle at which the first residue of the chain will be attached. 958 node_end(`str`): name of the ending node particle at which the last residue of the chain will be attached. 959 node_positions(`dict`): dictionary with the positions of the nodes. The keys are the node names and the values are the positions. 960 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 961 962 Note: 963 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. 964 The chain will be placed in the direction of the vector between `node_start` and `node_end`. 965 """ 966 if self.lattice_builder is None: 967 raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.") 968 969 molecule_name = "chain_"+node_start+"_"+node_end 970 sequence = self.df[self.df['name']==molecule_name].residue_list.values [0] 971 assert len(sequence) != 0 and not isinstance(sequence, str) 972 assert len(sequence) == self.lattice_builder.mpc 973 974 key, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end) 975 assert node_start != node_end or sequence == sequence[::-1], \ 976 (f"chain cannot be defined between '{node_start}' and '{node_end}' since it " 977 "would form a loop with a non-symmetric sequence (under-defined stereocenter)") 978 979 if reverse: 980 sequence = sequence[::-1] 981 982 node_start_pos = np.array(list(int(x) for x in node_start.strip('[]').split()))*0.25*self.lattice_builder.box_l 983 node_end_pos = np.array(list(int(x) for x in node_end.strip('[]').split()))*0.25*self.lattice_builder.box_l 984 node1 = espresso_system.part.select(lambda p: (p.pos == node_start_pos).all()).id 985 node2 = espresso_system.part.select(lambda p: (p.pos == node_end_pos).all()).id 986 987 if not node1[0] in node_positions or not node2[0] in node_positions: 988 raise ValueError("Set node position before placing a chain between them") 989 990 # Finding a backbone vector between node_start and node_end 991 vec_between_nodes = np.array(node_positions[node2[0]]) - np.array(node_positions[node1[0]]) 992 vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l) 993 backbone_vector = list(vec_between_nodes/(self.lattice_builder.mpc + 1)) 994 node_start_name = self.df[(self.df["particle_id"]==node1[0]) & (self.df["pmb_type"]=="particle")]["name"].values[0] 995 first_res_name = self.df[(self.df["pmb_type"]=="residue") & (self.df["name"]==sequence[0])]["central_bead"].values[0] 996 l0 = self.get_bond_length(node_start_name, first_res_name, hard_check=True) 997 chain_molecule_info = self.create_molecule( 998 name=molecule_name, # Use the name defined earlier 999 number_of_molecules=1, # Creating one chain 1000 espresso_system=espresso_system, 1001 list_of_first_residue_positions=[list(np.array(node_positions[node1[0]]) + np.array(backbone_vector))],#Start at the first node 1002 backbone_vector=np.array(backbone_vector)/l0, 1003 use_default_bond=False # Use defaut bonds between monomers 1004 ) 1005 # Collecting ids of beads of the chain/molecule 1006 chain_ids = [] 1007 residue_ids = [] 1008 for molecule_id in chain_molecule_info: 1009 for residue_id in chain_molecule_info[molecule_id]: 1010 residue_ids.append(residue_id) 1011 bead_id = chain_molecule_info[molecule_id][residue_id]['central_bead_id'] 1012 chain_ids.append(bead_id) 1013 1014 self.lattice_builder.chains[key] = sequence 1015 # Search bonds between nodes and chain ends 1016 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] 1017 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] 1018 bond_node1_first_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_start], 1019 particle_name2 = BeadType_near_to_node_start, 1020 hard_check=False, 1021 use_default_bond=False) 1022 bond_node2_last_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_end], 1023 particle_name2 = BeadType_near_to_node_end, 1024 hard_check=False, 1025 use_default_bond=False) 1026 1027 espresso_system.part.by_id(node1[0]).add_bond((bond_node1_first_monomer, chain_ids[0])) 1028 espresso_system.part.by_id(node2[0]).add_bond((bond_node2_last_monomer, chain_ids[-1])) 1029 # Add bonds to data frame 1030 bond_index1 = self.add_bond_in_df(particle_id1=node1[0], 1031 particle_id2=chain_ids[0], 1032 use_default_bond=False) 1033 self.add_value_to_df(key=('molecule_id',''), 1034 index=int(bond_index1), 1035 new_value=molecule_id, 1036 overwrite=True) 1037 self.add_value_to_df(key=('residue_id',''), 1038 index=int (bond_index1), 1039 new_value=residue_ids[0], 1040 overwrite=True) 1041 bond_index2 = self.add_bond_in_df(particle_id1=node2[0], 1042 particle_id2=chain_ids[-1], 1043 use_default_bond=False) 1044 self.add_value_to_df(key=('molecule_id',''), 1045 index=int(bond_index2), 1046 new_value=molecule_id, 1047 overwrite=True) 1048 self.add_value_to_df(key=('residue_id',''), 1049 index=int (bond_index2), 1050 new_value=residue_ids[-1], 1051 overwrite=True) 1052 return chain_molecule_info 1053 1054 def create_hydrogel_node(self, node_index, node_name, espresso_system): 1055 """ 1056 Set a node residue type. 1057 1058 Args: 1059 node_index(`str`): Lattice node index in the form of a string, e.g. "[0 0 0]". 1060 node_name(`str`): name of the node particle defined in pyMBE. 1061 Returns: 1062 node_position(`list`): Position of the node in the lattice. 1063 p_id(`int`): Particle ID of the node. 1064 """ 1065 if self.lattice_builder is None: 1066 raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.") 1067 1068 node_position = np.array(list(int(x) for x in node_index.strip('[]').split()))*0.25*self.lattice_builder.box_l 1069 p_id = self.create_particle(name = node_name, 1070 espresso_system=espresso_system, 1071 number_of_particles=1, 1072 position = [node_position]) 1073 key = self.lattice_builder._get_node_by_label(node_index) 1074 self.lattice_builder.nodes[key] = node_name 1075 1076 return node_position.tolist(), p_id 1077 1078 def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False): 1079 """ 1080 Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1081 1082 Args: 1083 name(`str`): Label of the molecule 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 number_of_molecules(`int`): Number of molecules of type `name` to be created. 1086 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. 1087 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. 1088 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df` 1089 1090 Returns: 1091 molecules_info(`dict`): {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 1092 1093 Note: 1094 Despite its name, this function can be used to create both molecules and peptides. 1095 """ 1096 if not self._check_if_name_is_defined_in_df(name=name): 1097 logging.warning(f"Molecule with name '{name}' is not defined in the pyMBE DataFrame, no molecule will be created.") 1098 return {} 1099 if number_of_molecules <= 0: 1100 return {} 1101 if list_of_first_residue_positions is not None: 1102 for item in list_of_first_residue_positions: 1103 if not isinstance(item, list): 1104 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.") 1105 elif len(item) != 3: 1106 raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.") 1107 1108 if len(list_of_first_residue_positions) != number_of_molecules: 1109 raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}") 1110 1111 # This function works for both molecules and peptides 1112 if not self._check_if_name_has_right_type(name=name, expected_pmb_type="molecule", hard_check=False): 1113 self._check_if_name_has_right_type(name=name, expected_pmb_type="peptide") 1114 1115 # Generate an arbitrary random unit vector 1116 if backbone_vector is None: 1117 backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0], 1118 radius=1, 1119 n_samples=1, 1120 on_surface=True)[0] 1121 else: 1122 backbone_vector = np.array(backbone_vector) 1123 1124 1125 1126 first_residue = True 1127 molecules_info = {} 1128 residue_list = self.df[self.df['name']==name].residue_list.values [0] 1129 self.copy_df_entry(name=name, 1130 column_name='molecule_id', 1131 number_of_copies=number_of_molecules) 1132 1133 molecules_index = np.where(self.df['name']==name) 1134 molecule_index_list =list(molecules_index[0])[-number_of_molecules:] 1135 pos_index = 0 1136 for molecule_index in molecule_index_list: 1137 molecule_id = self.assign_molecule_id(molecule_index=molecule_index) 1138 molecules_info[molecule_id] = {} 1139 for residue in residue_list: 1140 if first_residue: 1141 if list_of_first_residue_positions is None: 1142 residue_position = None 1143 else: 1144 for item in list_of_first_residue_positions: 1145 residue_position = [np.array(list_of_first_residue_positions[pos_index])] 1146 1147 residues_info = self.create_residue(name=residue, 1148 espresso_system=espresso_system, 1149 central_bead_position=residue_position, 1150 use_default_bond= use_default_bond, 1151 backbone_vector=backbone_vector) 1152 residue_id = next(iter(residues_info)) 1153 # Add the correct molecule_id to all particles in the residue 1154 for index in self.df[self.df['residue_id']==residue_id].index: 1155 self.add_value_to_df(key=('molecule_id',''), 1156 index=int (index), 1157 new_value=molecule_id, 1158 overwrite=True) 1159 central_bead_id = residues_info[residue_id]['central_bead_id'] 1160 previous_residue = residue 1161 residue_position = espresso_system.part.by_id(central_bead_id).pos 1162 previous_residue_id = central_bead_id 1163 first_residue = False 1164 else: 1165 previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0] 1166 new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0] 1167 bond = self.search_bond(particle_name1=previous_central_bead_name, 1168 particle_name2=new_central_bead_name, 1169 hard_check=True, 1170 use_default_bond=use_default_bond) 1171 l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 1172 particle_name2=new_central_bead_name, 1173 hard_check=True, 1174 use_default_bond=use_default_bond) 1175 1176 residue_position = residue_position+backbone_vector*l0 1177 residues_info = self.create_residue(name=residue, 1178 espresso_system=espresso_system, 1179 central_bead_position=[residue_position], 1180 use_default_bond= use_default_bond, 1181 backbone_vector=backbone_vector) 1182 residue_id = next(iter(residues_info)) 1183 for index in self.df[self.df['residue_id']==residue_id].index: 1184 self.add_value_to_df(key=('molecule_id',''), 1185 index=int (index), 1186 new_value=molecule_id, 1187 overwrite=True) 1188 central_bead_id = residues_info[residue_id]['central_bead_id'] 1189 espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id)) 1190 bond_index = self.add_bond_in_df(particle_id1=central_bead_id, 1191 particle_id2=previous_residue_id, 1192 use_default_bond=use_default_bond) 1193 self.add_value_to_df(key=('molecule_id',''), 1194 index=int (bond_index), 1195 new_value=molecule_id, 1196 overwrite=True) 1197 previous_residue_id = central_bead_id 1198 previous_residue = residue 1199 molecules_info[molecule_id][residue_id] = residues_info[residue_id] 1200 first_residue = True 1201 pos_index+=1 1202 1203 return molecules_info 1204 1205 def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False): 1206 """ 1207 Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`. 1208 1209 Args: 1210 name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 1211 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 1212 number_of_particles(`int`): Number of particles to be created. 1213 position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None. 1214 fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False. 1215 Returns: 1216 created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`. 1217 """ 1218 if number_of_particles <=0: 1219 return [] 1220 if not self._check_if_name_is_defined_in_df(name=name): 1221 logging.warning(f"Particle with name '{name}' is not defined in the pyMBE DataFrame, no particle will be created.") 1222 return [] 1223 self._check_if_name_has_right_type(name=name, 1224 expected_pmb_type="particle") 1225 # Copy the data of the particle `number_of_particles` times in the `df` 1226 self.copy_df_entry(name=name, 1227 column_name='particle_id', 1228 number_of_copies=number_of_particles) 1229 # Get information from the particle type `name` from the df 1230 z = self.df.loc[self.df['name']==name].state_one.z.values[0] 1231 z = 0. if z is None else z 1232 es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 1233 # Get a list of the index in `df` corresponding to the new particles to be created 1234 index = np.where(self.df['name']==name) 1235 index_list =list(index[0])[-number_of_particles:] 1236 # Create the new particles into `espresso_system` 1237 created_pid_list=[] 1238 for index in range (number_of_particles): 1239 df_index=int (index_list[index]) 1240 self.clean_df_row(index=df_index) 1241 if position is None: 1242 particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l) 1243 else: 1244 particle_position = position[index] 1245 if len(espresso_system.part.all()) == 0: 1246 bead_id = 0 1247 else: 1248 bead_id = max (espresso_system.part.all().id) + 1 1249 created_pid_list.append(bead_id) 1250 kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z) 1251 if fix: 1252 kwargs["fix"] = 3 * [fix] 1253 espresso_system.part.add(**kwargs) 1254 self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id) 1255 return created_pid_list 1256 1257 def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None): 1258 """ 1259 Creates all `particle`s associated to `pmb object` into `espresso` a number of times equal to `number_of_objects`. 1260 1261 Args: 1262 name(`str`): Unique label of the `pmb object` to be created. 1263 number_of_objects(`int`): Number of `pmb object`s to be created. 1264 espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library. 1265 position(`list`): Coordinates where the particles should be created. 1266 use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`. 1267 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. 1268 1269 Note: 1270 - If no `position` is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length. 1271 """ 1272 pmb_type = self._check_supported_molecule(molecule_name=name, 1273 valid_pmb_types=["molecule","particle","residue","peptide"]) 1274 1275 if pmb_type == 'particle': 1276 self.create_particle(name=name, 1277 number_of_particles=number_of_objects, 1278 espresso_system=espresso_system, 1279 position=position) 1280 elif pmb_type == 'residue': 1281 self.create_residue(name=name, 1282 espresso_system=espresso_system, 1283 central_bead_position=position, 1284 use_default_bond=use_default_bond, 1285 backbone_vector=backbone_vector) 1286 elif pmb_type in ['molecule','peptide']: 1287 self.create_molecule(name=name, 1288 number_of_molecules=number_of_objects, 1289 espresso_system=espresso_system, 1290 use_default_bond=use_default_bond, 1291 list_of_first_residue_positions=position, 1292 backbone_vector=backbone_vector) 1293 return 1294 1295 def create_protein(self, name, number_of_proteins, espresso_system, topology_dict): 1296 """ 1297 Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions` 1298 1299 Args: 1300 name(`str`): Label of the protein to be created. 1301 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 1302 number_of_proteins(`int`): Number of proteins to be created. 1303 positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}} 1304 """ 1305 1306 if number_of_proteins <=0: 1307 return 1308 if not self._check_if_name_is_defined_in_df(name=name): 1309 logging.warning(f"Protein with name '{name}' is not defined in the pyMBE DataFrame, no protein will be created.") 1310 return 1311 self._check_if_name_has_right_type(name=name, 1312 expected_pmb_type="protein") 1313 self.copy_df_entry(name=name, 1314 column_name='molecule_id', 1315 number_of_copies=number_of_proteins) 1316 protein_index = np.where(self.df['name']==name) 1317 protein_index_list =list(protein_index[0])[-number_of_proteins:] 1318 1319 box_half=espresso_system.box_l[0]/2.0 1320 1321 for molecule_index in protein_index_list: 1322 1323 molecule_id = self.assign_molecule_id(molecule_index=molecule_index) 1324 1325 protein_center = self.generate_coordinates_outside_sphere(radius = 1, 1326 max_dist=box_half, 1327 n_samples=1, 1328 center=[box_half]*3)[0] 1329 1330 for residue in topology_dict.keys(): 1331 1332 residue_name = re.split(r'\d+', residue)[0] 1333 residue_number = re.split(r'(\d+)', residue)[1] 1334 residue_position = topology_dict[residue]['initial_pos'] 1335 position = residue_position + protein_center 1336 1337 particle_id = self.create_particle(name=residue_name, 1338 espresso_system=espresso_system, 1339 number_of_particles=1, 1340 position=[position], 1341 fix = True) 1342 1343 index = self.df[self.df['particle_id']==particle_id[0]].index.values[0] 1344 self.add_value_to_df(key=('residue_id',''), 1345 index=int (index), 1346 new_value=int(residue_number), 1347 overwrite=True) 1348 1349 self.add_value_to_df(key=('molecule_id',''), 1350 index=int (index), 1351 new_value=molecule_id, 1352 overwrite=True) 1353 1354 return 1355 1356 def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None): 1357 """ 1358 Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1359 1360 Args: 1361 name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df` 1362 espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library. 1363 central_bead_position(`list` of `float`): Position of the central bead. 1364 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` 1365 backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`. 1366 1367 Returns: 1368 residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1369 """ 1370 if not self._check_if_name_is_defined_in_df(name=name): 1371 logging.warning(f"Residue with name '{name}' is not defined in the pyMBE DataFrame, no residue will be created.") 1372 return 1373 self._check_if_name_has_right_type(name=name, 1374 expected_pmb_type="residue") 1375 1376 # Copy the data of a residue in the `df 1377 self.copy_df_entry(name=name, 1378 column_name='residue_id', 1379 number_of_copies=1) 1380 residues_index = np.where(self.df['name']==name) 1381 residue_index_list =list(residues_index[0])[-1:] 1382 # search for defined particle and residue names 1383 particle_and_residue_df = self.df.loc[(self.df['pmb_type']== "particle") | (self.df['pmb_type']== "residue")] 1384 particle_and_residue_names = particle_and_residue_df["name"].tolist() 1385 for residue_index in residue_index_list: 1386 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1387 for side_chain_element in side_chain_list: 1388 if side_chain_element not in particle_and_residue_names: 1389 raise ValueError (f"{side_chain_element} is not defined") 1390 # Internal bookkepping of the residue info (important for side-chain residues) 1391 # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1392 residues_info={} 1393 for residue_index in residue_index_list: 1394 self.clean_df_row(index=int(residue_index)) 1395 # Assign a residue_id 1396 if self.df['residue_id'].isnull().all(): 1397 residue_id=0 1398 else: 1399 residue_id = self.df['residue_id'].max() + 1 1400 self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id) 1401 # create the principal bead 1402 central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0] 1403 central_bead_id = self.create_particle(name=central_bead_name, 1404 espresso_system=espresso_system, 1405 position=central_bead_position, 1406 number_of_particles = 1)[0] 1407 central_bead_position=espresso_system.part.by_id(central_bead_id).pos 1408 #assigns same residue_id to the central_bead particle created. 1409 index = self.df[self.df['particle_id']==central_bead_id].index.values[0] 1410 self.df.at [index,'residue_id'] = residue_id 1411 # Internal bookkeeping of the central bead id 1412 residues_info[residue_id]={} 1413 residues_info[residue_id]['central_bead_id']=central_bead_id 1414 # create the lateral beads 1415 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1416 side_chain_beads_ids = [] 1417 for side_chain_element in side_chain_list: 1418 pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 1419 if pmb_type == 'particle': 1420 bond = self.search_bond(particle_name1=central_bead_name, 1421 particle_name2=side_chain_element, 1422 hard_check=True, 1423 use_default_bond=use_default_bond) 1424 l0 = self.get_bond_length(particle_name1=central_bead_name, 1425 particle_name2=side_chain_element, 1426 hard_check=True, 1427 use_default_bond=use_default_bond) 1428 1429 if backbone_vector is None: 1430 bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1431 radius=l0, 1432 n_samples=1, 1433 on_surface=True)[0] 1434 else: 1435 bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector), 1436 magnitude=l0) 1437 1438 side_bead_id = self.create_particle(name=side_chain_element, 1439 espresso_system=espresso_system, 1440 position=[bead_position], 1441 number_of_particles=1)[0] 1442 index = self.df[self.df['particle_id']==side_bead_id].index.values[0] 1443 self.add_value_to_df(key=('residue_id',''), 1444 index=int (index), 1445 new_value=residue_id, 1446 overwrite=True) 1447 side_chain_beads_ids.append(side_bead_id) 1448 espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id)) 1449 index = self.add_bond_in_df(particle_id1=central_bead_id, 1450 particle_id2=side_bead_id, 1451 use_default_bond=use_default_bond) 1452 self.add_value_to_df(key=('residue_id',''), 1453 index=int (index), 1454 new_value=residue_id, 1455 overwrite=True) 1456 1457 elif pmb_type == 'residue': 1458 central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0] 1459 bond = self.search_bond(particle_name1=central_bead_name, 1460 particle_name2=central_bead_side_chain, 1461 hard_check=True, 1462 use_default_bond=use_default_bond) 1463 l0 = self.get_bond_length(particle_name1=central_bead_name, 1464 particle_name2=central_bead_side_chain, 1465 hard_check=True, 1466 use_default_bond=use_default_bond) 1467 if backbone_vector is None: 1468 residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1469 radius=l0, 1470 n_samples=1, 1471 on_surface=True)[0] 1472 else: 1473 residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1474 magnitude=l0) 1475 lateral_residue_info = self.create_residue(name=side_chain_element, 1476 espresso_system=espresso_system, 1477 central_bead_position=[residue_position], 1478 use_default_bond=use_default_bond) 1479 lateral_residue_dict=list(lateral_residue_info.values())[0] 1480 central_bead_side_chain_id=lateral_residue_dict['central_bead_id'] 1481 lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids'] 1482 residue_id_side_chain=list(lateral_residue_info.keys())[0] 1483 # Change the residue_id of the residue in the side chain to the one of the bigger residue 1484 index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0] 1485 self.add_value_to_df(key=('residue_id',''), 1486 index=int(index), 1487 new_value=residue_id, 1488 overwrite=True) 1489 # Change the residue_id of the particles in the residue in the side chain 1490 side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids 1491 for particle_id in side_chain_beads_ids: 1492 index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0] 1493 self.add_value_to_df(key=('residue_id',''), 1494 index=int (index), 1495 new_value=residue_id, 1496 overwrite=True) 1497 espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id)) 1498 index = self.add_bond_in_df(particle_id1=central_bead_id, 1499 particle_id2=central_bead_side_chain_id, 1500 use_default_bond=use_default_bond) 1501 self.add_value_to_df(key=('residue_id',''), 1502 index=int (index), 1503 new_value=residue_id, 1504 overwrite=True) 1505 # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue 1506 for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index: 1507 self.add_value_to_df(key=('residue_id',''), 1508 index=int(index), 1509 new_value=residue_id, 1510 overwrite=True) 1511 # Internal bookkeeping of the side chain beads ids 1512 residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids 1513 return residues_info 1514 1515 def create_variable_with_units(self, variable): 1516 """ 1517 Returns a pint object with the value and units defined in `variable`. 1518 1519 Args: 1520 variable(`dict` or `str`): {'value': value, 'units': units} 1521 Returns: 1522 variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry. 1523 """ 1524 if isinstance(variable, dict): 1525 value=variable.pop('value') 1526 units=variable.pop('units') 1527 elif isinstance(variable, str): 1528 value = float(re.split(r'\s+', variable)[0]) 1529 units = re.split(r'\s+', variable)[1] 1530 variable_with_units=value*self.units(units) 1531 return variable_with_units 1532 1533 def define_AA_residues(self, sequence, model): 1534 """ 1535 Defines in `pmb.df` all the different residues in `sequence`. 1536 1537 Args: 1538 sequence(`lst`): Sequence of the peptide or protein. 1539 model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1540 1541 Returns: 1542 residue_list(`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1543 """ 1544 residue_list = [] 1545 for residue_name in sequence: 1546 if model == '1beadAA': 1547 central_bead = residue_name 1548 side_chains = [] 1549 elif model == '2beadAA': 1550 if residue_name in ['c','n', 'G']: 1551 central_bead = residue_name 1552 side_chains = [] 1553 else: 1554 central_bead = 'CA' 1555 side_chains = [residue_name] 1556 if residue_name not in residue_list: 1557 self.define_residue(name = 'AA-'+residue_name, 1558 central_bead = central_bead, 1559 side_chains = side_chains) 1560 residue_list.append('AA-'+residue_name) 1561 return residue_list 1562 1563 def define_bond(self, bond_type, bond_parameters, particle_pairs): 1564 1565 ''' 1566 Defines a pmb object of type `bond` in `pymbe.df`. 1567 1568 Args: 1569 bond_type(`str`): label to identify the potential to model the bond. 1570 bond_parameters(`dict`): parameters of the potential of the bond. 1571 particle_pairs(`lst`): list of the `names` of the `particles` to be bonded. 1572 1573 Note: 1574 Currently, only HARMONIC and FENE bonds are supported. 1575 1576 For a HARMONIC bond the dictionary must contain the following parameters: 1577 1578 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 1579 using the `pmb.units` UnitRegistry. 1580 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 1581 the `pmb.units` UnitRegistry. 1582 1583 For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and: 1584 1585 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 1586 units of length using the `pmb.units` UnitRegistry. Default 'None'. 1587 ''' 1588 1589 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1590 for particle_name1, particle_name2 in particle_pairs: 1591 1592 lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1, 1593 particle_name2 = particle_name2, 1594 combining_rule = 'Lorentz-Berthelot') 1595 1596 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1597 bond_type = bond_type, 1598 epsilon = lj_parameters["epsilon"], 1599 sigma = lj_parameters["sigma"], 1600 cutoff = lj_parameters["cutoff"], 1601 offset = lj_parameters["offset"],) 1602 index = len(self.df) 1603 for label in [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']: 1604 self._check_if_multiple_pmb_types_for_name(name=label, pmb_type_to_be_defined="bond") 1605 name=f'{particle_name1}-{particle_name2}' 1606 self._check_if_multiple_pmb_types_for_name(name=name, pmb_type_to_be_defined="bond") 1607 self.df.at [index,'name']= name 1608 self.df.at [index,'bond_object'] = bond_object 1609 self.df.at [index,'l0'] = l0 1610 self.add_value_to_df(index=index, 1611 key=('pmb_type',''), 1612 new_value='bond') 1613 self.add_value_to_df(index=index, 1614 key=('parameters_of_the_potential',''), 1615 new_value=bond_object.get_params(), 1616 non_standard_value=True) 1617 self.df.fillna(pd.NA, inplace=True) 1618 return 1619 1620 1621 def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None): 1622 """ 1623 Asigns `bond` in `pmb.df` as the default bond. 1624 The LJ parameters can be optionally provided to calculate the initial bond length 1625 1626 Args: 1627 bond_type(`str`): label to identify the potential to model the bond. 1628 bond_parameters(`dict`): parameters of the potential of the bond. 1629 sigma(`float`, optional): LJ sigma of the interaction between the particles. 1630 epsilon(`float`, optional): LJ epsilon for the interaction between the particles. 1631 cutoff(`float`, optional): cutoff-radius of the LJ interaction. 1632 offset(`float`, optional): offset of the LJ interaction. 1633 1634 Note: 1635 - Currently, only harmonic and FENE bonds are supported. 1636 """ 1637 1638 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1639 1640 if epsilon is None: 1641 epsilon=1*self.units('reduced_energy') 1642 if sigma is None: 1643 sigma=1*self.units('reduced_length') 1644 if cutoff is None: 1645 cutoff=2**(1.0/6.0)*self.units('reduced_length') 1646 if offset is None: 1647 offset=0*self.units('reduced_length') 1648 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1649 bond_type = bond_type, 1650 epsilon = epsilon, 1651 sigma = sigma, 1652 cutoff = cutoff, 1653 offset = offset) 1654 1655 self._check_if_multiple_pmb_types_for_name(name='default', 1656 pmb_type_to_be_defined='bond') 1657 1658 if len(self.df.index) != 0: 1659 index = max(self.df.index)+1 1660 else: 1661 index = 0 1662 self.df.at [index,'name'] = 'default' 1663 self.df.at [index,'bond_object'] = bond_object 1664 self.df.at [index,'l0'] = l0 1665 self.add_value_to_df(index = index, 1666 key = ('pmb_type',''), 1667 new_value = 'bond') 1668 self.add_value_to_df(index = index, 1669 key = ('parameters_of_the_potential',''), 1670 new_value = bond_object.get_params(), 1671 non_standard_value=True) 1672 self.df.fillna(pd.NA, inplace=True) 1673 return 1674 1675 def define_hydrogel(self, name, node_map, chain_map): 1676 """ 1677 Defines a pyMBE object of type `hydrogel` in `pymbe.df`. 1678 1679 Args: 1680 name(`str`): Unique label that identifies the `hydrogel`. 1681 node_map(`list of ict`): [{"particle_name": , "lattice_index": }, ... ] 1682 chain_map(`list of dict`): [{"node_start": , "node_end": , "residue_list": , ... ] 1683 """ 1684 node_indices = {tuple(entry['lattice_index']) for entry in node_map} 1685 diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices} 1686 if node_indices != diamond_indices: 1687 raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ") 1688 1689 chain_map_connectivity = set() 1690 for entry in chain_map: 1691 start = self.lattice_builder.node_labels[entry['node_start']] 1692 end = self.lattice_builder.node_labels[entry['node_end']] 1693 chain_map_connectivity.add((start,end)) 1694 1695 if self.lattice_builder.lattice.connectivity != chain_map_connectivity: 1696 raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs") 1697 1698 self._check_if_multiple_pmb_types_for_name(name=name, 1699 pmb_type_to_be_defined='hydrogel') 1700 1701 1702 index = len(self.df) 1703 self.df.at [index, "name"] = name 1704 self.df.at [index, "pmb_type"] = "hydrogel" 1705 self.add_value_to_df(index = index, 1706 key = ('node_map',''), 1707 new_value = node_map, 1708 non_standard_value=True) 1709 self.add_value_to_df(index = index, 1710 key = ('chain_map',''), 1711 new_value = chain_map, 1712 non_standard_value=True) 1713 for chain_label in chain_map: 1714 node_start = chain_label["node_start"] 1715 node_end = chain_label["node_end"] 1716 residue_list = chain_label['residue_list'] 1717 # Molecule name 1718 molecule_name = "chain_"+node_start+"_"+node_end 1719 self.define_molecule(name=molecule_name, residue_list=residue_list) 1720 return 1721 1722 def define_molecule(self, name, residue_list): 1723 """ 1724 Defines a pyMBE object of type `molecule` in `pymbe.df`. 1725 1726 Args: 1727 name(`str`): Unique label that identifies the `molecule`. 1728 residue_list(`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1729 """ 1730 self._check_if_multiple_pmb_types_for_name(name=name, 1731 pmb_type_to_be_defined='molecule') 1732 1733 index = len(self.df) 1734 self.df.at [index,'name'] = name 1735 self.df.at [index,'pmb_type'] = 'molecule' 1736 self.df.at [index,('residue_list','')] = residue_list 1737 self.df.fillna(pd.NA, inplace=True) 1738 return 1739 1740 def define_particle_entry_in_df(self,name): 1741 """ 1742 Defines a particle entry in pmb.df. 1743 1744 Args: 1745 name(`str`): Unique label that identifies this particle type. 1746 1747 Returns: 1748 index(`int`): Index of the particle in pmb.df 1749 """ 1750 1751 if self._check_if_name_is_defined_in_df(name=name): 1752 index = self.df[self.df['name']==name].index[0] 1753 else: 1754 index = len(self.df) 1755 self.df.at [index, 'name'] = name 1756 self.df.at [index,'pmb_type'] = 'particle' 1757 self.df.fillna(pd.NA, inplace=True) 1758 return index 1759 1760 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): 1761 """ 1762 Defines the properties of a particle object. 1763 1764 Args: 1765 name(`str`): Unique label that identifies this particle type. 1766 z(`int`, optional): Permanent charge number of this particle type. Defaults to 0. 1767 acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to pd.NA. 1768 pka(`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to pd.NA. 1769 sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1770 cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1771 offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1772 epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA. 1773 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1774 1775 Note: 1776 - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units. 1777 - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units. 1778 - `cutoff` defaults to `2**(1./6.) reduced_length`. 1779 - `offset` defaults to 0. 1780 - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions. 1781 - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`. 1782 """ 1783 index=self.define_particle_entry_in_df(name=name) 1784 self._check_if_multiple_pmb_types_for_name(name=name, 1785 pmb_type_to_be_defined='particle') 1786 1787 # If `cutoff` and `offset` are not defined, default them to the following values 1788 if pd.isna(cutoff): 1789 cutoff=self.units.Quantity(2**(1./6.), "reduced_length") 1790 if pd.isna(offset): 1791 offset=self.units.Quantity(0, "reduced_length") 1792 1793 # Define LJ parameters 1794 parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"}, 1795 "cutoff":{"value": cutoff, "dimensionality": "[length]"}, 1796 "offset":{"value": offset, "dimensionality": "[length]"}, 1797 "epsilon":{"value": epsilon, "dimensionality": "[energy]"},} 1798 1799 for parameter_key in parameters_with_dimensionality.keys(): 1800 if not pd.isna(parameters_with_dimensionality[parameter_key]["value"]): 1801 self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 1802 expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"]) 1803 self.add_value_to_df(key=(parameter_key,''), 1804 index=index, 1805 new_value=parameters_with_dimensionality[parameter_key]["value"], 1806 overwrite=overwrite) 1807 1808 # Define particle acid/base properties 1809 self.set_particle_acidity(name=name, 1810 acidity=acidity, 1811 default_charge_number=z, 1812 pka=pka, 1813 overwrite=overwrite) 1814 self.df.fillna(pd.NA, inplace=True) 1815 return 1816 1817 def define_particles(self, parameters, overwrite=False): 1818 ''' 1819 Defines a particle object in pyMBE for each particle name in `particle_names` 1820 1821 Args: 1822 parameters(`dict`): dictionary with the particle parameters. 1823 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1824 1825 Note: 1826 - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},} 1827 ''' 1828 if not parameters: 1829 return 0 1830 for particle_name in parameters.keys(): 1831 parameters[particle_name]["overwrite"]=overwrite 1832 self.define_particle(**parameters[particle_name]) 1833 return 1834 1835 def define_peptide(self, name, sequence, model): 1836 """ 1837 Defines a pyMBE object of type `peptide` in the `pymbe.df`. 1838 1839 Args: 1840 name (`str`): Unique label that identifies the `peptide`. 1841 sequence (`string`): Sequence of the `peptide`. 1842 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1843 """ 1844 self._check_if_multiple_pmb_types_for_name(name = name, 1845 pmb_type_to_be_defined='peptide') 1846 1847 valid_keys = ['1beadAA','2beadAA'] 1848 if model not in valid_keys: 1849 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1850 clean_sequence = self.protein_sequence_parser(sequence=sequence) 1851 residue_list = self.define_AA_residues(sequence=clean_sequence, 1852 model=model) 1853 self.define_molecule(name = name, residue_list=residue_list) 1854 index = self.df.loc[self.df['name'] == name].index.item() 1855 self.df.at [index,'model'] = model 1856 self.df.at [index,('sequence','')] = clean_sequence 1857 self.df.at [index,'pmb_type'] = "peptide" 1858 self.df.fillna(pd.NA, inplace=True) 1859 1860 1861 def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False): 1862 """ 1863 Defines a globular protein pyMBE object in `pymbe.df`. 1864 1865 Args: 1866 name (`str`): Unique label that identifies the protein. 1867 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1868 topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value} 1869 lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca". 1870 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1871 1872 Note: 1873 - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential. 1874 """ 1875 1876 # Sanity checks 1877 self._check_if_multiple_pmb_types_for_name(name = name, 1878 pmb_type_to_be_defined='protein') 1879 valid_model_keys = ['1beadAA','2beadAA'] 1880 valid_lj_setups = ["wca"] 1881 if model not in valid_model_keys: 1882 raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}') 1883 if lj_setup_mode not in valid_lj_setups: 1884 raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}') 1885 if lj_setup_mode == "wca": 1886 sigma = 1*self.units.Quantity("reduced_length") 1887 epsilon = 1*self.units.Quantity("reduced_energy") 1888 part_dict={} 1889 sequence=[] 1890 metal_ions_charge_number_map=self.get_metal_ions_charge_number_map() 1891 for particle in topology_dict.keys(): 1892 particle_name = re.split(r'\d+', particle)[0] 1893 if particle_name not in part_dict.keys(): 1894 if lj_setup_mode == "wca": 1895 part_dict[particle_name]={"sigma": sigma, 1896 "offset": topology_dict[particle]['radius']*2-sigma, 1897 "epsilon": epsilon, 1898 "name": particle_name} 1899 if self.check_if_metal_ion(key=particle_name): 1900 z=metal_ions_charge_number_map[particle_name] 1901 else: 1902 z=0 1903 part_dict[particle_name]["z"]=z 1904 1905 if self.check_aminoacid_key(key=particle_name): 1906 sequence.append(particle_name) 1907 1908 self.define_particles(parameters=part_dict, 1909 overwrite=overwrite) 1910 residue_list = self.define_AA_residues(sequence=sequence, 1911 model=model) 1912 index = len(self.df) 1913 self.df.at [index,'name'] = name 1914 self.df.at [index,'pmb_type'] = 'protein' 1915 self.df.at [index,'model'] = model 1916 self.df.at [index,('sequence','')] = sequence 1917 self.df.at [index,('residue_list','')] = residue_list 1918 self.df.fillna(pd.NA, inplace=True) 1919 return 1920 1921 def define_residue(self, name, central_bead, side_chains): 1922 """ 1923 Defines a pyMBE object of type `residue` in `pymbe.df`. 1924 1925 Args: 1926 name(`str`): Unique label that identifies the `residue`. 1927 central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`. 1928 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. 1929 """ 1930 self._check_if_multiple_pmb_types_for_name(name=name, 1931 pmb_type_to_be_defined='residue') 1932 1933 index = len(self.df) 1934 self.df.at [index, 'name'] = name 1935 self.df.at [index,'pmb_type'] = 'residue' 1936 self.df.at [index,'central_bead'] = central_bead 1937 self.df.at [index,('side_chains','')] = side_chains 1938 self.df.fillna(pd.NA, inplace=True) 1939 return 1940 1941 def delete_entries_in_df(self, entry_name): 1942 """ 1943 Deletes entries with name `entry_name` from the DataFrame if it exists. 1944 1945 Args: 1946 entry_name (`str`): The name of the entry in the dataframe to delete. 1947 1948 """ 1949 if entry_name in self.df["name"].values: 1950 self.df = self.df[self.df["name"] != entry_name].reset_index(drop=True) 1951 1952 def destroy_pmb_object_in_system(self, name, espresso_system): 1953 """ 1954 Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 1955 1956 Args: 1957 name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`. 1958 espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library. 1959 1960 Note: 1961 - If `name` is a object_type=`particle`, only the matching particles that are not part of bigger objects (e.g. `residue`, `molecule`) will be destroyed. To destroy particles in such objects, destroy the bigger object instead. 1962 """ 1963 pmb_type = self._check_supported_molecule(molecule_name=name, 1964 valid_pmb_types= ['particle','residue','molecule',"peptide"]) 1965 1966 if pmb_type == 'particle': 1967 particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1968 & (self.df['molecule_id'].isna())] 1969 particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1970 & (self.df['molecule_id'].isna())].particle_id.tolist() 1971 for particle_id in particle_ids_list: 1972 espresso_system.part.by_id(particle_id).remove() 1973 self.df = self.df.drop(index=particles_index) 1974 if pmb_type == 'residue': 1975 residues_id = self.df.loc[self.df['name']== name].residue_id.to_list() 1976 for residue_id in residues_id: 1977 molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0] 1978 particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"] 1979 self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index) 1980 for particle_id in particle_ids_list: 1981 espresso_system.part.by_id(particle_id).remove() 1982 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1983 if pmb_type in ['molecule',"peptide"]: 1984 molecules_id = self.df.loc[self.df['name']== name].molecule_id.dropna().to_list() 1985 for molecule_id in molecules_id: 1986 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']==pmb_type)].name.values[0] 1987 particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"] 1988 self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index) 1989 for particle_id in particle_ids_list: 1990 espresso_system.part.by_id(particle_id).remove() 1991 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1992 1993 self.df.reset_index(drop=True,inplace=True) 1994 1995 return 1996 1997 def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200): 1998 """ 1999 Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration. 2000 To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an 2001 ion pair and the concentrations of the various species is solved numerically using a self-consistent approach. 2002 More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831). 2003 This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237). 2004 2005 Args: 2006 pH_res('float'): pH-value in the reservoir. 2007 c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2008 activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2009 2010 Returns: 2011 cH_res('pint.Quantity'): Concentration of H+ ions. 2012 cOH_res('pint.Quantity'): Concentration of OH- ions. 2013 cNa_res('pint.Quantity'): Concentration of Na+ ions. 2014 cCl_res('pint.Quantity'): Concentration of Cl- ions. 2015 """ 2016 2017 self_consistent_run = 0 2018 cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+ 2019 2020 def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res): 2021 #Calculate and initial guess for the concentrations of various species based on ideal gas estimate 2022 cOH_res = self.Kw / cH_res 2023 cNa_res = None 2024 cCl_res = None 2025 if cOH_res>=cH_res: 2026 #adjust the concentration of sodium if there is excess OH- in the reservoir: 2027 cNa_res = c_salt_res + (cOH_res-cH_res) 2028 cCl_res = c_salt_res 2029 else: 2030 # adjust the concentration of chloride if there is excess H+ in the reservoir 2031 cCl_res = c_salt_res + (cH_res-cOH_res) 2032 cNa_res = c_salt_res 2033 2034 def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res): 2035 nonlocal max_number_sc_runs, self_consistent_run 2036 if self_consistent_run<max_number_sc_runs: 2037 self_consistent_run+=1 2038 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2039 cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res)) 2040 if cOH_res>=cH_res: 2041 #adjust the concentration of sodium if there is excess OH- in the reservoir: 2042 cNa_res = c_salt_res + (cOH_res-cH_res) 2043 cCl_res = c_salt_res 2044 else: 2045 # adjust the concentration of chloride if there is excess H+ in the reservoir 2046 cCl_res = c_salt_res + (cH_res-cOH_res) 2047 cNa_res = c_salt_res 2048 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 2049 else: 2050 return cH_res, cOH_res, cNa_res, cCl_res 2051 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 2052 2053 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 2054 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2055 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 2056 2057 while abs(determined_pH-pH_res)>1e-6: 2058 if determined_pH > pH_res: 2059 cH_res *= 1.005 2060 else: 2061 cH_res /= 1.003 2062 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 2063 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2064 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 2065 self_consistent_run=0 2066 2067 return cH_res, cOH_res, cNa_res, cCl_res 2068 2069 def enable_motion_of_rigid_object(self, name, espresso_system): 2070 ''' 2071 Enables the motion of the rigid object `name` in the `espresso_system`. 2072 2073 Args: 2074 name(`str`): Label of the object. 2075 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 2076 2077 Note: 2078 - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]. 2079 ''' 2080 logging.info('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]') 2081 self._check_supported_molecule(molecule_name=name, 2082 valid_pmb_types= ['protein']) 2083 molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list() 2084 for molecule_id in molecule_ids_list: 2085 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 2086 center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system) 2087 rigid_object_center = espresso_system.part.add(pos=center_of_mass, 2088 rotation=[True,True,True], 2089 type=self.propose_unused_type()) 2090 2091 rigid_object_center.mass = len(particle_ids_list) 2092 momI = 0 2093 for pid in particle_ids_list: 2094 momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2) 2095 2096 rigid_object_center.rinertia = np.ones(3) * momI 2097 2098 for particle_id in particle_ids_list: 2099 pid = espresso_system.part.by_id(particle_id) 2100 pid.vs_auto_relate_to(rigid_object_center.id) 2101 return 2102 2103 def filter_df(self, pmb_type): 2104 """ 2105 Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns. 2106 2107 Args: 2108 pmb_type(`str`): pmb_object_type to filter in `pmb.df`. 2109 2110 Returns: 2111 pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`. 2112 """ 2113 pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type] 2114 pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1) 2115 return pmb_type_df 2116 2117 def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False): 2118 """ 2119 Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2120 2121 Args: 2122 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 2123 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 2124 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'. 2125 2126 Returns: 2127 bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` if a matching bond exists 2128 2129 Note: 2130 - If `use_default_bond`=`True`, it returns "default" if no key is found. 2131 """ 2132 bond_keys = [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}'] 2133 bond_defined=False 2134 for bond_key in bond_keys: 2135 if bond_key in self.df["name"].values: 2136 bond_defined=True 2137 correct_key=bond_key 2138 break 2139 if bond_defined: 2140 return correct_key 2141 elif use_default_bond: 2142 return 'default' 2143 else: 2144 return 2145 2146 def find_value_from_es_type(self, es_type, column_name): 2147 """ 2148 Finds a value in `pmb.df` for a `column_name` and `es_type` pair. 2149 2150 Args: 2151 es_type(`int`): value of the espresso type 2152 column_name(`str`): name of the column in `pymbe.df` 2153 2154 Returns: 2155 Value in `pymbe.df` matching `column_name` and `es_type` 2156 """ 2157 idx = pd.IndexSlice 2158 for state in ['state_one', 'state_two']: 2159 index = self.df.loc[self.df[(state, 'es_type')] == es_type].index 2160 if len(index) > 0: 2161 if column_name == 'label': 2162 label = self.df.loc[idx[index[0]], idx[(state,column_name)]] 2163 return label 2164 else: 2165 column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]] 2166 return column_name_value 2167 return None 2168 2169 def format_node(self, node_list): 2170 return "[" + " ".join(map(str, node_list)) + "]" 2171 2172 2173 def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples): 2174 """ 2175 Generates coordinates outside a sphere centered at `center`. 2176 2177 Args: 2178 center(`lst`): Coordinates of the center of the sphere. 2179 radius(`float`): Radius of the sphere. 2180 max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates. 2181 n_samples(`int`): Number of sample points. 2182 2183 Returns: 2184 coord_list(`lst`): Coordinates of the sample points. 2185 """ 2186 if not radius > 0: 2187 raise ValueError (f'The value of {radius} must be a positive value') 2188 if not radius < max_dist: 2189 raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))') 2190 coord_list = [] 2191 counter = 0 2192 while counter<n_samples: 2193 coord = self.generate_random_points_in_a_sphere(center=center, 2194 radius=max_dist, 2195 n_samples=1)[0] 2196 if np.linalg.norm(coord-np.asarray(center))>=radius: 2197 coord_list.append (coord) 2198 counter += 1 2199 return coord_list 2200 2201 def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False): 2202 """ 2203 Uniformly samples points from a hypersphere. If on_surface is set to True, the points are 2204 uniformly sampled from the surface of the hypersphere. 2205 2206 Args: 2207 center(`lst`): Array with the coordinates of the center of the spheres. 2208 radius(`float`): Radius of the sphere. 2209 n_samples(`int`): Number of sample points to generate inside the sphere. 2210 on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere. 2211 2212 Returns: 2213 samples(`list`): Coordinates of the sample points inside the hypersphere. 2214 """ 2215 # initial values 2216 center=np.array(center) 2217 d = center.shape[0] 2218 # sample n_samples points in d dimensions from a standard normal distribution 2219 samples = self.rng.normal(size=(n_samples, d)) 2220 # make the samples lie on the surface of the unit hypersphere 2221 normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis] 2222 samples /= normalize_radii 2223 if not on_surface: 2224 # make the samples lie inside the hypersphere with the correct density 2225 uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis] 2226 new_radii = np.power(uniform_points, 1/d) 2227 samples *= new_radii 2228 # scale the points to have the correct radius and center 2229 samples = samples * radius + center 2230 return samples 2231 2232 def generate_trial_perpendicular_vector(self,vector,magnitude): 2233 """ 2234 Generates an orthogonal vector to the input `vector`. 2235 2236 Args: 2237 vector(`lst`): arbitrary vector. 2238 magnitude(`float`): magnitude of the orthogonal vector. 2239 2240 Returns: 2241 (`lst`): Orthogonal vector with the same magnitude as the input vector. 2242 """ 2243 np_vec = np.array(vector) 2244 if np.all(np_vec == 0): 2245 raise ValueError('Zero vector') 2246 np_vec /= np.linalg.norm(np_vec) 2247 # Generate a random vector 2248 random_vector = self.generate_random_points_in_a_sphere(radius=1, 2249 center=[0,0,0], 2250 n_samples=1, 2251 on_surface=True)[0] 2252 # Project the random vector onto the input vector and subtract the projection 2253 projection = np.dot(random_vector, np_vec) * np_vec 2254 perpendicular_vector = random_vector - projection 2255 # Normalize the perpendicular vector to have the same magnitude as the input vector 2256 perpendicular_vector /= np.linalg.norm(perpendicular_vector) 2257 return perpendicular_vector*magnitude 2258 2259 def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2260 """ 2261 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length. 2262 If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead. 2263 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. 2264 2265 Args: 2266 particle_name1(str): label of the type of the first particle type of the bonded particles. 2267 particle_name2(str): label of the type of the second particle type of the bonded particles. 2268 hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2269 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. 2270 2271 Returns: 2272 l0(`pint.Quantity`): bond length 2273 2274 Note: 2275 - 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`. 2276 - If `hard_check`=`True` stops the code when no bond is found. 2277 """ 2278 bond_key = self.find_bond_key(particle_name1=particle_name1, 2279 particle_name2=particle_name2, 2280 use_default_bond=use_default_bond) 2281 if bond_key: 2282 return self.df[self.df['name']==bond_key].l0.values[0] 2283 else: 2284 msg = f"Bond not defined between particles {particle_name1} and {particle_name2}" 2285 if hard_check: 2286 raise ValueError(msg) 2287 else: 2288 logging.warning(msg) 2289 return 2290 2291 def get_charge_number_map(self): 2292 ''' 2293 Gets the charge number of each `espresso_type` in `pymbe.df`. 2294 2295 Returns: 2296 charge_number_map(`dict`): {espresso_type: z}. 2297 ''' 2298 if self.df.state_one['es_type'].isnull().values.any(): 2299 df_state_one = self.df.state_one.dropna() 2300 df_state_two = self.df.state_two.dropna() 2301 else: 2302 df_state_one = self.df.state_one 2303 if self.df.state_two['es_type'].isnull().values.any(): 2304 df_state_two = self.df.state_two.dropna() 2305 else: 2306 df_state_two = self.df.state_two 2307 state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values) 2308 state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values) 2309 charge_number_map = pd.concat([state_one,state_two],axis=0).to_dict() 2310 return charge_number_map 2311 2312 def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'): 2313 """ 2314 Returns the Lennard-Jones parameters for the interaction between the particle types given by 2315 `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule. 2316 2317 Args: 2318 particle_name1 (str): label of the type of the first particle type 2319 particle_name2 (str): label of the type of the second particle type 2320 combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'. 2321 2322 Returns: 2323 {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value} 2324 2325 Note: 2326 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 2327 - 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. 2328 """ 2329 supported_combining_rules=["Lorentz-Berthelot"] 2330 lj_parameters_keys=["sigma","epsilon","offset","cutoff"] 2331 if combining_rule not in supported_combining_rules: 2332 raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}") 2333 lj_parameters={} 2334 for key in lj_parameters_keys: 2335 lj_parameters[key]=[] 2336 # Search the LJ parameters of the type pair 2337 for name in [particle_name1,particle_name2]: 2338 for key in lj_parameters_keys: 2339 lj_parameters[key].append(self.df[self.df.name == name][key].values[0]) 2340 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 2341 if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]): 2342 return {} 2343 # Apply combining rule 2344 if combining_rule == 'Lorentz-Berthelot': 2345 lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2 2346 lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2 2347 lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2 2348 lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1]) 2349 return lj_parameters 2350 2351 def get_metal_ions_charge_number_map(self): 2352 """ 2353 Gets a map with the charge numbers of all the metal ions supported. 2354 2355 Returns: 2356 metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number} 2357 2358 """ 2359 metal_charge_number_map = {"Ca": 2} 2360 return metal_charge_number_map 2361 2362 def get_particle_id_map(self, object_name): 2363 ''' 2364 Gets all the ids associated with the object with name `object_name` in `pmb.df` 2365 2366 Args: 2367 object_name(`str`): name of the object 2368 2369 Returns: 2370 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]}, } 2371 ''' 2372 object_type=self._check_supported_molecule(molecule_name=object_name, 2373 valid_pmb_types= ['particle','residue','molecule',"peptide","protein"]) 2374 id_list = [] 2375 mol_map = {} 2376 res_map = {} 2377 def do_res_map(res_ids): 2378 for res_id in res_ids: 2379 res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist() 2380 res_map[res_id]=res_list 2381 return res_map 2382 if object_type in ['molecule', 'protein', 'peptide']: 2383 mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist() 2384 for mol_id in mol_ids: 2385 res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist()) 2386 res_map=do_res_map(res_ids=res_ids) 2387 mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist() 2388 id_list+=mol_list 2389 mol_map[mol_id]=mol_list 2390 elif object_type == 'residue': 2391 res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist() 2392 res_map=do_res_map(res_ids=res_ids) 2393 id_list=[] 2394 for res_id_list in res_map.values(): 2395 id_list+=res_id_list 2396 elif object_type == 'particle': 2397 id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist() 2398 return {"all": id_list, "molecule_map": mol_map, "residue_map": res_map} 2399 2400 def get_pka_set(self): 2401 ''' 2402 Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df` 2403 2404 Returns: 2405 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 2406 ''' 2407 titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna() 2408 pka_set = {} 2409 for index in titratables_AA_df.name.keys(): 2410 name = titratables_AA_df.name[index] 2411 pka_value = titratables_AA_df.pka[index] 2412 acidity = titratables_AA_df.acidity[index] 2413 pka_set[name] = {'pka_value':pka_value,'acidity':acidity} 2414 return pka_set 2415 2416 def get_radius_map(self, dimensionless=True): 2417 ''' 2418 Gets the effective radius of each `espresso type` in `pmb.df`. 2419 2420 Args: 2421 dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False. 2422 2423 Returns: 2424 radius_map(`dict`): {espresso_type: radius}. 2425 2426 Note: 2427 The radius corresponds to (sigma+offset)/2 2428 ''' 2429 df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates() 2430 df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates() 2431 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) 2432 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) 2433 radius_map = pd.concat([state_one,state_two],axis=0).to_dict() 2434 if dimensionless: 2435 for key in radius_map: 2436 radius_map[key] = radius_map[key].magnitude 2437 return radius_map 2438 2439 def get_reduced_units(self): 2440 """ 2441 Returns the current set of reduced units defined in pyMBE. 2442 2443 Returns: 2444 reduced_units_text(`str`): text with information about the current set of reduced units. 2445 2446 """ 2447 unit_length=self.units.Quantity(1,'reduced_length') 2448 unit_energy=self.units.Quantity(1,'reduced_energy') 2449 unit_charge=self.units.Quantity(1,'reduced_charge') 2450 reduced_units_text = "\n".join(["Current set of reduced units:", 2451 f"{unit_length.to('nm'):.5g} = {unit_length}", 2452 f"{unit_energy.to('J'):.5g} = {unit_energy}", 2453 f"{unit_charge.to('C'):.5g} = {unit_charge}", 2454 f"Temperature: {(self.kT/self.kB).to('K'):.5g}" 2455 ]) 2456 return reduced_units_text 2457 2458 def get_type_map(self): 2459 """ 2460 Gets all different espresso types assigned to particles in `pmb.df`. 2461 2462 Returns: 2463 type_map(`dict`): {"name": espresso_type}. 2464 """ 2465 df_state_one = self.df.state_one.dropna(how='all') 2466 df_state_two = self.df.state_two.dropna(how='all') 2467 state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label) 2468 state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label) 2469 type_map = pd.concat([state_one,state_two],axis=0).to_dict() 2470 return type_map 2471 2472 def initialize_lattice_builder(self, diamond_lattice): 2473 """ 2474 Initialize the lattice builder with the DiamondLattice object. 2475 2476 Args: 2477 diamond_lattice(`DiamondLattice`): DiamondLattice object from the `lib/lattice` module to be used in the LatticeBuilder. 2478 """ 2479 from .lib.lattice import LatticeBuilder, DiamondLattice 2480 if not isinstance(diamond_lattice, DiamondLattice): 2481 raise TypeError("Currently only DiamondLattice objects are supported.") 2482 self.lattice_builder = LatticeBuilder(lattice=diamond_lattice) 2483 logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}") 2484 return self.lattice_builder 2485 2486 2487 def load_interaction_parameters(self, filename, overwrite=False): 2488 """ 2489 Loads the interaction parameters stored in `filename` into `pmb.df` 2490 2491 Args: 2492 filename(`str`): name of the file to be read 2493 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 2494 """ 2495 without_units = ['z','es_type'] 2496 with_units = ['sigma','epsilon','offset','cutoff'] 2497 2498 with open(filename, 'r') as f: 2499 interaction_data = json.load(f) 2500 interaction_parameter_set = interaction_data["data"] 2501 2502 for key in interaction_parameter_set: 2503 param_dict=interaction_parameter_set[key] 2504 object_type=param_dict.pop('object_type') 2505 if object_type == 'particle': 2506 not_required_attributes={} 2507 for not_required_key in without_units+with_units: 2508 if not_required_key in param_dict.keys(): 2509 if not_required_key in with_units: 2510 not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key)) 2511 elif not_required_key in without_units: 2512 not_required_attributes[not_required_key]=param_dict.pop(not_required_key) 2513 else: 2514 not_required_attributes[not_required_key]=pd.NA 2515 self.define_particle(name=param_dict.pop('name'), 2516 z=not_required_attributes.pop('z'), 2517 sigma=not_required_attributes.pop('sigma'), 2518 offset=not_required_attributes.pop('offset'), 2519 cutoff=not_required_attributes.pop('cutoff'), 2520 epsilon=not_required_attributes.pop('epsilon'), 2521 overwrite=overwrite) 2522 elif object_type == 'residue': 2523 self.define_residue(**param_dict) 2524 elif object_type == 'molecule': 2525 self.define_molecule(**param_dict) 2526 elif object_type == 'peptide': 2527 self.define_peptide(**param_dict) 2528 elif object_type == 'bond': 2529 particle_pairs = param_dict.pop('particle_pairs') 2530 bond_parameters = param_dict.pop('bond_parameters') 2531 bond_type = param_dict.pop('bond_type') 2532 if bond_type == 'harmonic': 2533 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2534 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2535 bond = {'r_0' : r_0, 2536 'k' : k, 2537 } 2538 2539 elif bond_type == 'FENE': 2540 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2541 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2542 d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max')) 2543 bond = {'r_0' : r_0, 2544 'k' : k, 2545 'd_r_max': d_r_max, 2546 } 2547 else: 2548 raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds") 2549 self.define_bond(bond_type=bond_type, 2550 bond_parameters=bond, 2551 particle_pairs=particle_pairs) 2552 else: 2553 raise ValueError(object_type+' is not a known pmb object type') 2554 2555 return 2556 2557 def load_pka_set(self, filename, overwrite=True): 2558 """ 2559 Loads the pka_set stored in `filename` into `pmb.df`. 2560 2561 Args: 2562 filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}. 2563 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 2564 """ 2565 with open(filename, 'r') as f: 2566 pka_data = json.load(f) 2567 pka_set = pka_data["data"] 2568 2569 self.check_pka_set(pka_set=pka_set) 2570 2571 for key in pka_set: 2572 acidity = pka_set[key]['acidity'] 2573 pka_value = pka_set[key]['pka_value'] 2574 self.set_particle_acidity(name=key, 2575 acidity=acidity, 2576 pka=pka_value, 2577 overwrite=overwrite) 2578 return 2579 2580 2581 def propose_unused_type(self): 2582 """ 2583 Searches in `pmb.df` all the different particle types defined and returns a new one. 2584 2585 Returns: 2586 unused_type(`int`): unused particle type 2587 """ 2588 type_map = self.get_type_map() 2589 if not type_map: 2590 unused_type = 0 2591 else: 2592 valid_values = [v for v in type_map.values() if pd.notna(v)] # Filter out pd.NA values 2593 unused_type = max(valid_values) + 1 if valid_values else 0 # Ensure max() doesn't fail if all values are NA 2594 return unused_type 2595 2596 def protein_sequence_parser(self, sequence): 2597 ''' 2598 Parses `sequence` to the one letter code for amino acids. 2599 2600 Args: 2601 sequence(`str` or `lst`): Sequence of the amino acid. 2602 2603 Returns: 2604 clean_sequence(`lst`): `sequence` using the one letter code. 2605 2606 Note: 2607 - Accepted formats for `sequence` are: 2608 - `lst` with one letter or three letter code of each aminoacid in each element 2609 - `str` with the sequence using the one letter code 2610 - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-" 2611 2612 ''' 2613 # Aminoacid key 2614 keys={"ALA": "A", 2615 "ARG": "R", 2616 "ASN": "N", 2617 "ASP": "D", 2618 "CYS": "C", 2619 "GLU": "E", 2620 "GLN": "Q", 2621 "GLY": "G", 2622 "HIS": "H", 2623 "ILE": "I", 2624 "LEU": "L", 2625 "LYS": "K", 2626 "MET": "M", 2627 "PHE": "F", 2628 "PRO": "P", 2629 "SER": "S", 2630 "THR": "T", 2631 "TRP": "W", 2632 "TYR": "Y", 2633 "VAL": "V", 2634 "PSER": "J", 2635 "PTHR": "U", 2636 "PTyr": "Z", 2637 "NH2": "n", 2638 "COOH": "c"} 2639 clean_sequence=[] 2640 if isinstance(sequence, str): 2641 if sequence.find("-") != -1: 2642 splited_sequence=sequence.split("-") 2643 for residue in splited_sequence: 2644 if len(residue) == 1: 2645 if residue in keys.values(): 2646 residue_ok=residue 2647 else: 2648 if residue.upper() in keys.values(): 2649 residue_ok=residue.upper() 2650 else: 2651 raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence") 2652 clean_sequence.append(residue_ok) 2653 else: 2654 if residue in keys.keys(): 2655 clean_sequence.append(keys[residue]) 2656 else: 2657 if residue.upper() in keys.keys(): 2658 clean_sequence.append(keys[residue.upper()]) 2659 else: 2660 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2661 else: 2662 for residue in sequence: 2663 if residue in keys.values(): 2664 residue_ok=residue 2665 else: 2666 if residue.upper() in keys.values(): 2667 residue_ok=residue.upper() 2668 else: 2669 raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence") 2670 clean_sequence.append(residue_ok) 2671 if isinstance(sequence, list): 2672 for residue in sequence: 2673 if residue in keys.values(): 2674 residue_ok=residue 2675 else: 2676 if residue.upper() in keys.values(): 2677 residue_ok=residue.upper() 2678 elif (residue.upper() in keys.keys()): 2679 residue_ok= keys[residue.upper()] 2680 else: 2681 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2682 clean_sequence.append(residue_ok) 2683 return clean_sequence 2684 2685 def read_pmb_df (self,filename): 2686 """ 2687 Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE. 2688 2689 Args: 2690 filename(`str`): path to file. 2691 2692 Note: 2693 This function only accepts files with CSV format. 2694 """ 2695 2696 if filename.rsplit(".", 1)[1] != "csv": 2697 raise ValueError("Only files with CSV format are supported") 2698 df = pd.read_csv (filename,header=[0, 1], index_col=0) 2699 columns_names = self.setup_df() 2700 2701 multi_index = pd.MultiIndex.from_tuples(columns_names) 2702 df.columns = multi_index 2703 2704 self.df = self.convert_columns_to_original_format(df) 2705 self.df.fillna(pd.NA, inplace=True) 2706 return self.df 2707 2708 def read_protein_vtf_in_df (self,filename,unit_length=None): 2709 """ 2710 Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`. 2711 2712 Args: 2713 filename(`str`): Path to the vtf file with the coarse-grained model. 2714 unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None. 2715 2716 Returns: 2717 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 2718 2719 Note: 2720 - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom. 2721 """ 2722 2723 logging.info(f'Loading protein coarse grain model file: {filename}') 2724 2725 coord_list = [] 2726 particles_dict = {} 2727 2728 if unit_length is None: 2729 unit_length = 1 * self.units.angstrom 2730 2731 with open (filename,'r') as protein_model: 2732 for line in protein_model : 2733 line_split = line.split() 2734 if line_split: 2735 line_header = line_split[0] 2736 if line_header == 'atom': 2737 atom_id = line_split[1] 2738 atom_name = line_split[3] 2739 atom_resname = line_split[5] 2740 chain_id = line_split[9] 2741 radius = float(line_split [11])*unit_length 2742 particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius] 2743 elif line_header.isnumeric(): 2744 atom_coord = line_split[1:] 2745 atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord] 2746 coord_list.append (atom_coord) 2747 2748 numbered_label = [] 2749 i = 0 2750 2751 for atom_id in particles_dict.keys(): 2752 2753 if atom_id == 1: 2754 atom_name = particles_dict[atom_id][0] 2755 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2756 numbered_label.append(numbered_name) 2757 2758 elif atom_id != 1: 2759 2760 if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]: 2761 i += 1 2762 count = 1 2763 atom_name = particles_dict[atom_id][0] 2764 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2765 numbered_label.append(numbered_name) 2766 2767 elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]: 2768 if count == 2 or particles_dict[atom_id][1] == 'GLY': 2769 i +=1 2770 count = 0 2771 atom_name = particles_dict[atom_id][0] 2772 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2773 numbered_label.append(numbered_name) 2774 count +=1 2775 2776 topology_dict = {} 2777 2778 for i in range (0, len(numbered_label)): 2779 topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] , 2780 'chain_id':numbered_label[i][1], 2781 'radius':numbered_label[i][2] } 2782 2783 return topology_dict 2784 2785 def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2786 """ 2787 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2788 If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead. 2789 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. 2790 2791 Args: 2792 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 2793 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 2794 hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2795 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. 2796 2797 Returns: 2798 bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library. 2799 2800 Note: 2801 - 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`. 2802 - If `hard_check`=`True` stops the code when no bond is found. 2803 """ 2804 2805 bond_key = self.find_bond_key(particle_name1=particle_name1, 2806 particle_name2=particle_name2, 2807 use_default_bond=use_default_bond) 2808 if use_default_bond: 2809 if not self._check_if_name_is_defined_in_df(name="default"): 2810 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") 2811 if bond_key: 2812 return self.df[self.df['name']==bond_key].bond_object.values[0] 2813 else: 2814 msg= f"Bond not defined between particles {particle_name1} and {particle_name2}" 2815 if hard_check: 2816 raise ValueError(msg) 2817 else: 2818 logging.warning(msg) 2819 return None 2820 def search_particles_in_residue(self, residue_name): 2821 ''' 2822 Searches for all particles in a given residue of name `residue_name`. 2823 2824 Args: 2825 residue_name (`str`): name of the residue to be searched 2826 2827 Returns: 2828 list_of_particles_in_residue (`lst`): list of the names of all particles in the residue 2829 2830 Note: 2831 - 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. 2832 - The function will return an empty list if the residue is not defined in `pmb.df`. 2833 - The function will return an empty list if the particles are not defined in the pyMBE DataFrame. 2834 ''' 2835 if not self._check_if_name_is_defined_in_df(name=residue_name): 2836 logging.warning(f"Residue {residue_name} not defined in pmb.df") 2837 return [] 2838 self._check_if_name_has_right_type(name=residue_name, expected_pmb_type="residue") 2839 index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 2840 central_bead = self.df.at [index_residue, ('central_bead', '')] 2841 list_of_side_chains = self.df.at[index_residue, ('side_chains', '')] 2842 list_of_particles_in_residue = [] 2843 if central_bead is not pd.NA: 2844 if self._check_if_name_is_defined_in_df(name=central_bead): 2845 if self._check_if_name_has_right_type(name=central_bead, expected_pmb_type="particle", hard_check=False): 2846 list_of_particles_in_residue.append(central_bead) 2847 if list_of_side_chains is not pd.NA: 2848 for side_chain in list_of_side_chains: 2849 if self._check_if_name_is_defined_in_df(name=side_chain): 2850 object_type = self.df[self.df['name']==side_chain].pmb_type.values[0] 2851 else: 2852 continue 2853 if object_type == "residue": 2854 list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain) 2855 list_of_particles_in_residue += list_of_particles_in_side_chain_residue 2856 elif object_type == "particle": 2857 if side_chain is not pd.NA: 2858 list_of_particles_in_residue.append(side_chain) 2859 return list_of_particles_in_residue 2860 2861 def set_particle_acidity(self, name, acidity=pd.NA, default_charge_number=0, pka=pd.NA, overwrite=True): 2862 """ 2863 Sets the particle acidity including the charges in each of its possible states. 2864 2865 Args: 2866 name(`str`): Unique label that identifies the `particle`. 2867 acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None. 2868 default_charge_number (`int`): Charge number of the particle. Defaults to 0. 2869 pka(`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to pandas.NA. 2870 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 2871 2872 Note: 2873 - For particles with `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 2874 deprotonated states, respectively. 2875 - For particles without an acidity `acidity = pandas.NA`, only `state_one` is defined. 2876 - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 2877 """ 2878 acidity_valid_keys = ['inert','acidic', 'basic'] 2879 if not pd.isna(acidity): 2880 if acidity not in acidity_valid_keys: 2881 raise ValueError(f"Acidity {acidity} provided for particle name {name} is not supproted. Valid keys are: {acidity_valid_keys}") 2882 if acidity in ['acidic', 'basic'] and pd.isna(pka): 2883 raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.") 2884 if acidity == "inert": 2885 acidity = pd.NA 2886 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") 2887 2888 self.define_particle_entry_in_df(name=name) 2889 2890 for index in self.df[self.df['name']==name].index: 2891 if pka is not pd.NA: 2892 self.add_value_to_df(key=('pka',''), 2893 index=index, 2894 new_value=pka, 2895 overwrite=overwrite) 2896 2897 self.add_value_to_df(key=('acidity',''), 2898 index=index, 2899 new_value=acidity, 2900 overwrite=overwrite) 2901 if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')): 2902 self.add_value_to_df(key=('state_one','es_type'), 2903 index=index, 2904 new_value=self.propose_unused_type(), 2905 overwrite=overwrite) 2906 if pd.isna(self.df.loc [self.df['name'] == name].acidity.iloc[0]): 2907 self.add_value_to_df(key=('state_one','z'), 2908 index=index, 2909 new_value=default_charge_number, 2910 overwrite=overwrite) 2911 self.add_value_to_df(key=('state_one','label'), 2912 index=index, 2913 new_value=name, 2914 overwrite=overwrite) 2915 else: 2916 protonated_label = f'{name}H' 2917 self.add_value_to_df(key=('state_one','label'), 2918 index=index, 2919 new_value=protonated_label, 2920 overwrite=overwrite) 2921 self.add_value_to_df(key=('state_two','label'), 2922 index=index, 2923 new_value=name, 2924 overwrite=overwrite) 2925 if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')): 2926 self.add_value_to_df(key=('state_two','es_type'), 2927 index=index, 2928 new_value=self.propose_unused_type(), 2929 overwrite=overwrite) 2930 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'acidic': 2931 self.add_value_to_df(key=('state_one','z'), 2932 index=index,new_value=0, 2933 overwrite=overwrite) 2934 self.add_value_to_df(key=('state_two','z'), 2935 index=index, 2936 new_value=-1, 2937 overwrite=overwrite) 2938 elif self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'basic': 2939 self.add_value_to_df(key=('state_one','z'), 2940 index=index,new_value=+1, 2941 overwrite=overwrite) 2942 self.add_value_to_df(key=('state_two','z'), 2943 index=index, 2944 new_value=0, 2945 overwrite=overwrite) 2946 self.df.fillna(pd.NA, inplace=True) 2947 return 2948 2949 def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None): 2950 """ 2951 Sets the set of reduced units used by pyMBE.units and it prints it. 2952 2953 Args: 2954 unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 2955 unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 2956 temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 2957 Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 2958 2959 Note: 2960 - If no `temperature` is given, a value of 298.15 K is assumed by default. 2961 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 2962 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 2963 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 2964 """ 2965 if unit_length is None: 2966 unit_length= 0.355*self.units.nm 2967 if temperature is None: 2968 temperature = 298.15 * self.units.K 2969 if unit_charge is None: 2970 unit_charge = scipy.constants.e * self.units.C 2971 if Kw is None: 2972 Kw = 1e-14 2973 # Sanity check 2974 variables=[unit_length,temperature,unit_charge] 2975 dimensionalities=["[length]","[temperature]","[charge]"] 2976 for variable,dimensionality in zip(variables,dimensionalities): 2977 self.check_dimensionality(variable,dimensionality) 2978 self.Kw=Kw*self.units.mol**2 / (self.units.l**2) 2979 self.kT=temperature*self.kB 2980 self.units._build_cache() 2981 self.units.define(f'reduced_energy = {self.kT} ') 2982 self.units.define(f'reduced_length = {unit_length}') 2983 self.units.define(f'reduced_charge = {unit_charge}') 2984 logging.info(self.get_reduced_units()) 2985 return 2986 2987 def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2988 """ 2989 Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 2990 2991 Args: 2992 counter_ion(`str`): `name` of the counter_ion `particle`. 2993 constant_pH(`float`): pH-value. 2994 exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted. 2995 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2996 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2997 2998 Returns: 2999 RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library. 3000 sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3001 """ 3002 from espressomd import reaction_methods 3003 if exclusion_range is None: 3004 exclusion_range = max(self.get_radius_map().values())*2.0 3005 if pka_set is None: 3006 pka_set=self.get_pka_set() 3007 self.check_pka_set(pka_set=pka_set) 3008 if use_exclusion_radius_per_type: 3009 exclusion_radius_per_type = self.get_radius_map() 3010 else: 3011 exclusion_radius_per_type = {} 3012 3013 RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3014 exclusion_range=exclusion_range, 3015 seed=self.seed, 3016 constant_pH=constant_pH, 3017 exclusion_radius_per_type = exclusion_radius_per_type 3018 ) 3019 sucessfull_reactions_labels=[] 3020 charge_number_map = self.get_charge_number_map() 3021 for name in pka_set.keys(): 3022 if not self._check_if_name_is_defined_in_df(name): 3023 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.') 3024 continue 3025 gamma=10**-pka_set[name]['pka_value'] 3026 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3027 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3028 counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0] 3029 RE.add_reaction(gamma=gamma, 3030 reactant_types=[state_one_type], 3031 product_types=[state_two_type, counter_ion_type], 3032 default_charges={state_one_type: charge_number_map[state_one_type], 3033 state_two_type: charge_number_map[state_two_type], 3034 counter_ion_type: charge_number_map[counter_ion_type]}) 3035 sucessfull_reactions_labels.append(name) 3036 return RE, sucessfull_reactions_labels 3037 3038 def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False): 3039 """ 3040 Sets up grand-canonical coupling to a reservoir of salt. 3041 For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead. 3042 3043 Args: 3044 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3045 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 3046 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 3047 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3048 exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted. 3049 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 3050 3051 Returns: 3052 RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library. 3053 """ 3054 from espressomd import reaction_methods 3055 if exclusion_range is None: 3056 exclusion_range = max(self.get_radius_map().values())*2.0 3057 if use_exclusion_radius_per_type: 3058 exclusion_radius_per_type = self.get_radius_map() 3059 else: 3060 exclusion_radius_per_type = {} 3061 3062 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3063 exclusion_range=exclusion_range, 3064 seed=self.seed, 3065 exclusion_radius_per_type = exclusion_radius_per_type 3066 ) 3067 3068 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3069 determined_activity_coefficient = activity_coefficient(c_salt_res) 3070 K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient 3071 3072 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 3073 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 3074 3075 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0] 3076 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0] 3077 3078 if salt_cation_charge <= 0: 3079 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 3080 if salt_anion_charge >= 0: 3081 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 3082 3083 # Grand-canonical coupling to the reservoir 3084 RE.add_reaction( 3085 gamma = K_salt.magnitude, 3086 reactant_types = [], 3087 reactant_coefficients = [], 3088 product_types = [ salt_cation_es_type, salt_anion_es_type ], 3089 product_coefficients = [ 1, 1 ], 3090 default_charges = { 3091 salt_cation_es_type: salt_cation_charge, 3092 salt_anion_es_type: salt_anion_charge, 3093 } 3094 ) 3095 3096 return RE 3097 3098 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): 3099 """ 3100 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 3101 reservoir of small ions. 3102 This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1]. 3103 3104 [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020. 3105 3106 Args: 3107 pH_res ('float): pH-value in the reservoir. 3108 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3109 proton_name ('str'): Name of the proton (H+) particle. 3110 hydroxide_name ('str'): Name of the hydroxide (OH-) particle. 3111 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 3112 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 3113 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3114 exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted. 3115 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 3116 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 3117 3118 Returns: 3119 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 3120 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3121 ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients). 3122 """ 3123 from espressomd import reaction_methods 3124 if exclusion_range is None: 3125 exclusion_range = max(self.get_radius_map().values())*2.0 3126 if pka_set is None: 3127 pka_set=self.get_pka_set() 3128 self.check_pka_set(pka_set=pka_set) 3129 if use_exclusion_radius_per_type: 3130 exclusion_radius_per_type = self.get_radius_map() 3131 else: 3132 exclusion_radius_per_type = {} 3133 3134 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3135 exclusion_range=exclusion_range, 3136 seed=self.seed, 3137 exclusion_radius_per_type = exclusion_radius_per_type 3138 ) 3139 3140 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3141 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 3142 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 3143 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 3144 K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3145 K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3146 K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3147 3148 proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0] 3149 hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0] 3150 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 3151 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 3152 3153 proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0] 3154 hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0] 3155 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0] 3156 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0] 3157 3158 if proton_charge <= 0: 3159 raise ValueError('ERROR proton charge must be positive, charge ', proton_charge) 3160 if salt_cation_charge <= 0: 3161 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 3162 if hydroxide_charge >= 0: 3163 raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge) 3164 if salt_anion_charge >= 0: 3165 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 3166 3167 # Grand-canonical coupling to the reservoir 3168 # 0 = H+ + OH- 3169 RE.add_reaction( 3170 gamma = K_W.magnitude, 3171 reactant_types = [], 3172 reactant_coefficients = [], 3173 product_types = [ proton_es_type, hydroxide_es_type ], 3174 product_coefficients = [ 1, 1 ], 3175 default_charges = { 3176 proton_es_type: proton_charge, 3177 hydroxide_es_type: hydroxide_charge, 3178 } 3179 ) 3180 3181 # 0 = Na+ + Cl- 3182 RE.add_reaction( 3183 gamma = K_NACL.magnitude, 3184 reactant_types = [], 3185 reactant_coefficients = [], 3186 product_types = [ salt_cation_es_type, salt_anion_es_type ], 3187 product_coefficients = [ 1, 1 ], 3188 default_charges = { 3189 salt_cation_es_type: salt_cation_charge, 3190 salt_anion_es_type: salt_anion_charge, 3191 } 3192 ) 3193 3194 # 0 = Na+ + OH- 3195 RE.add_reaction( 3196 gamma = (K_NACL * K_W / K_HCL).magnitude, 3197 reactant_types = [], 3198 reactant_coefficients = [], 3199 product_types = [ salt_cation_es_type, hydroxide_es_type ], 3200 product_coefficients = [ 1, 1 ], 3201 default_charges = { 3202 salt_cation_es_type: salt_cation_charge, 3203 hydroxide_es_type: hydroxide_charge, 3204 } 3205 ) 3206 3207 # 0 = H+ + Cl- 3208 RE.add_reaction( 3209 gamma = K_HCL.magnitude, 3210 reactant_types = [], 3211 reactant_coefficients = [], 3212 product_types = [ proton_es_type, salt_anion_es_type ], 3213 product_coefficients = [ 1, 1 ], 3214 default_charges = { 3215 proton_es_type: proton_charge, 3216 salt_anion_es_type: salt_anion_charge, 3217 } 3218 ) 3219 3220 # Annealing moves to ensure sufficient sampling 3221 # Cation annealing H+ = Na+ 3222 RE.add_reaction( 3223 gamma = (K_NACL / K_HCL).magnitude, 3224 reactant_types = [proton_es_type], 3225 reactant_coefficients = [ 1 ], 3226 product_types = [ salt_cation_es_type ], 3227 product_coefficients = [ 1 ], 3228 default_charges = { 3229 proton_es_type: proton_charge, 3230 salt_cation_es_type: salt_cation_charge, 3231 } 3232 ) 3233 3234 # Anion annealing OH- = Cl- 3235 RE.add_reaction( 3236 gamma = (K_HCL / K_W).magnitude, 3237 reactant_types = [hydroxide_es_type], 3238 reactant_coefficients = [ 1 ], 3239 product_types = [ salt_anion_es_type ], 3240 product_coefficients = [ 1 ], 3241 default_charges = { 3242 hydroxide_es_type: hydroxide_charge, 3243 salt_anion_es_type: salt_anion_charge, 3244 } 3245 ) 3246 3247 sucessful_reactions_labels=[] 3248 charge_number_map = self.get_charge_number_map() 3249 for name in pka_set.keys(): 3250 if not self._check_if_name_is_defined_in_df(name): 3251 logging.warning(f'The acid-base reaction of {name} has not been set up because its particle type is not defined in the dataframe.') 3252 continue 3253 3254 Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 3255 3256 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3257 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3258 3259 # Reaction in terms of proton: HA = A + H+ 3260 RE.add_reaction(gamma=Ka.magnitude, 3261 reactant_types=[state_one_type], 3262 reactant_coefficients=[1], 3263 product_types=[state_two_type, proton_es_type], 3264 product_coefficients=[1, 1], 3265 default_charges={state_one_type: charge_number_map[state_one_type], 3266 state_two_type: charge_number_map[state_two_type], 3267 proton_es_type: proton_charge}) 3268 3269 # Reaction in terms of salt cation: HA = A + Na+ 3270 RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude, 3271 reactant_types=[state_one_type], 3272 reactant_coefficients=[1], 3273 product_types=[state_two_type, salt_cation_es_type], 3274 product_coefficients=[1, 1], 3275 default_charges={state_one_type: charge_number_map[state_one_type], 3276 state_two_type: charge_number_map[state_two_type], 3277 salt_cation_es_type: salt_cation_charge}) 3278 3279 # Reaction in terms of hydroxide: OH- + HA = A 3280 RE.add_reaction(gamma=(Ka / K_W).magnitude, 3281 reactant_types=[state_one_type, hydroxide_es_type], 3282 reactant_coefficients=[1, 1], 3283 product_types=[state_two_type], 3284 product_coefficients=[1], 3285 default_charges={state_one_type: charge_number_map[state_one_type], 3286 state_two_type: charge_number_map[state_two_type], 3287 hydroxide_es_type: hydroxide_charge}) 3288 3289 # Reaction in terms of salt anion: Cl- + HA = A 3290 RE.add_reaction(gamma=(Ka / K_HCL).magnitude, 3291 reactant_types=[state_one_type, salt_anion_es_type], 3292 reactant_coefficients=[1, 1], 3293 product_types=[state_two_type], 3294 product_coefficients=[1], 3295 default_charges={state_one_type: charge_number_map[state_one_type], 3296 state_two_type: charge_number_map[state_two_type], 3297 salt_anion_es_type: salt_anion_charge}) 3298 3299 sucessful_reactions_labels.append(name) 3300 return RE, sucessful_reactions_labels, ionic_strength_res 3301 3302 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): 3303 """ 3304 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 3305 reservoir of small ions. 3306 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-}. 3307 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'. 3308 3309 [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4). 3310 [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020. 3311 3312 Args: 3313 pH_res ('float'): pH-value in the reservoir. 3314 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3315 cation_name ('str'): Name of the cationic particle. 3316 anion_name ('str'): Name of the anionic particle. 3317 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3318 exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted. 3319 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 3320 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`. 3321 3322 Returns: 3323 RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 3324 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3325 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 3326 """ 3327 from espressomd import reaction_methods 3328 if exclusion_range is None: 3329 exclusion_range = max(self.get_radius_map().values())*2.0 3330 if pka_set is None: 3331 pka_set=self.get_pka_set() 3332 self.check_pka_set(pka_set=pka_set) 3333 if use_exclusion_radius_per_type: 3334 exclusion_radius_per_type = self.get_radius_map() 3335 else: 3336 exclusion_radius_per_type = {} 3337 3338 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3339 exclusion_range=exclusion_range, 3340 seed=self.seed, 3341 exclusion_radius_per_type = exclusion_radius_per_type 3342 ) 3343 3344 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3345 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 3346 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 3347 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 3348 a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 3349 a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 3350 a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 3351 K_XX = a_cation * a_anion 3352 3353 cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0] 3354 anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0] 3355 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0] 3356 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0] 3357 if cation_charge <= 0: 3358 raise ValueError('ERROR cation charge must be positive, charge ', cation_charge) 3359 if anion_charge >= 0: 3360 raise ValueError('ERROR anion charge must be negative, charge ', anion_charge) 3361 3362 # Coupling to the reservoir: 0 = X+ + X- 3363 RE.add_reaction( 3364 gamma = K_XX.magnitude, 3365 reactant_types = [], 3366 reactant_coefficients = [], 3367 product_types = [ cation_es_type, anion_es_type ], 3368 product_coefficients = [ 1, 1 ], 3369 default_charges = { 3370 cation_es_type: cation_charge, 3371 anion_es_type: anion_charge, 3372 } 3373 ) 3374 3375 sucessful_reactions_labels=[] 3376 charge_number_map = self.get_charge_number_map() 3377 for name in pka_set.keys(): 3378 if not self._check_if_name_is_defined_in_df(name): 3379 logging.warning(f'The acid-base reaction of {name} has not been set up because its particle type is not defined in the dataframe.') 3380 continue 3381 3382 Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l 3383 gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen 3384 3385 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3386 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3387 3388 # Reaction in terms of small cation: HA = A + X+ 3389 RE.add_reaction(gamma=gamma_K_AX.magnitude, 3390 reactant_types=[state_one_type], 3391 reactant_coefficients=[1], 3392 product_types=[state_two_type, cation_es_type], 3393 product_coefficients=[1, 1], 3394 default_charges={state_one_type: charge_number_map[state_one_type], 3395 state_two_type: charge_number_map[state_two_type], 3396 cation_es_type: cation_charge}) 3397 3398 # Reaction in terms of small anion: X- + HA = A 3399 RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude, 3400 reactant_types=[state_one_type, anion_es_type], 3401 reactant_coefficients=[1, 1], 3402 product_types=[state_two_type], 3403 product_coefficients=[1], 3404 default_charges={state_one_type: charge_number_map[state_one_type], 3405 state_two_type: charge_number_map[state_two_type], 3406 anion_es_type: anion_charge}) 3407 3408 sucessful_reactions_labels.append(name) 3409 return RE, sucessful_reactions_labels, ionic_strength_res 3410 3411 def setup_df (self): 3412 """ 3413 Sets up the pyMBE's dataframe `pymbe.df`. 3414 3415 Returns: 3416 columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe 3417 """ 3418 3419 columns_dtypes = { 3420 'name': { 3421 '': str}, 3422 'pmb_type': { 3423 '': str}, 3424 'particle_id': { 3425 '': pd.Int64Dtype()}, 3426 'particle_id2': { 3427 '': pd.Int64Dtype()}, 3428 'residue_id': { 3429 '': pd.Int64Dtype()}, 3430 'molecule_id': { 3431 '': pd.Int64Dtype()}, 3432 'acidity': { 3433 '': str}, 3434 'pka': { 3435 '': object}, 3436 'central_bead': { 3437 '': object}, 3438 'side_chains': { 3439 '': object}, 3440 'residue_list': { 3441 '': object}, 3442 'model': { 3443 '': str}, 3444 'sigma': { 3445 '': object}, 3446 'cutoff': { 3447 '': object}, 3448 'offset': { 3449 '': object}, 3450 'epsilon': { 3451 '': object}, 3452 'state_one': { 3453 'label': str, 3454 'es_type': pd.Int64Dtype(), 3455 'z': pd.Int64Dtype()}, 3456 'state_two': { 3457 'label': str, 3458 'es_type': pd.Int64Dtype(), 3459 'z': pd.Int64Dtype()}, 3460 'sequence': { 3461 '': object}, 3462 'bond_object': { 3463 '': object}, 3464 'parameters_of_the_potential':{ 3465 '': object}, 3466 'l0': { 3467 '': float}, 3468 'node_map':{ 3469 '':object}, 3470 'chain_map':{ 3471 '':object}} 3472 3473 self.df = pd.DataFrame(columns=pd.MultiIndex.from_tuples([(col_main, col_sub) for col_main, sub_cols in columns_dtypes.items() for col_sub in sub_cols.keys()])) 3474 3475 for level1, sub_dtypes in columns_dtypes.items(): 3476 for level2, dtype in sub_dtypes.items(): 3477 self.df[level1, level2] = self.df[level1, level2].astype(dtype) 3478 3479 columns_names = pd.MultiIndex.from_frame(self.df) 3480 columns_names = columns_names.names 3481 3482 return columns_names 3483 3484 def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot'): 3485 """ 3486 Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`. 3487 3488 Args: 3489 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 3490 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. 3491 combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'. 3492 warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True. 3493 3494 Note: 3495 - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 3496 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 3497 - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html 3498 3499 """ 3500 from itertools import combinations_with_replacement 3501 compulsory_parameters_in_df = ['sigma','epsilon'] 3502 shift=0 3503 if shift_potential: 3504 shift="auto" 3505 # List which particles have sigma and epsilon values defined in pmb.df and which ones don't 3506 particles_types_with_LJ_parameters = [] 3507 non_parametrized_labels= [] 3508 for particle_type in self.get_type_map().values(): 3509 check_list=[] 3510 for key in compulsory_parameters_in_df: 3511 value_in_df=self.find_value_from_es_type(es_type=particle_type, 3512 column_name=key) 3513 check_list.append(pd.isna(value_in_df)) 3514 if any(check_list): 3515 non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 3516 column_name='label')) 3517 else: 3518 particles_types_with_LJ_parameters.append(particle_type) 3519 # Set up LJ interactions between all particle types 3520 for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 3521 particle_name1 = self.find_value_from_es_type(es_type=type_pair[0], 3522 column_name="name") 3523 particle_name2 = self.find_value_from_es_type(es_type=type_pair[1], 3524 column_name="name") 3525 lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1, 3526 particle_name2 = particle_name2, 3527 combining_rule = combining_rule) 3528 3529 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 3530 if not lj_parameters: 3531 continue 3532 espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 3533 sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 3534 cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude, 3535 offset = lj_parameters["offset"].to("reduced_length").magnitude, 3536 shift = shift) 3537 index = len(self.df) 3538 label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label") 3539 label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label") 3540 self.df.at [index, 'name'] = f'LJ: {label1}-{label2}' 3541 lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params() 3542 3543 self.add_value_to_df(index=index, 3544 key=('pmb_type',''), 3545 new_value='LennardJones') 3546 3547 self.add_value_to_df(index=index, 3548 key=('parameters_of_the_potential',''), 3549 new_value=lj_params, 3550 non_standard_value=True) 3551 if non_parametrized_labels: 3552 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.') 3553 return 3554 3555 def write_pmb_df (self, filename): 3556 ''' 3557 Writes the pyMBE dataframe into a csv file 3558 3559 Args: 3560 filename(`str`): Path to the csv file 3561 ''' 3562 3563 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence'] 3564 df = self.df.copy(deep=True) 3565 for column_name in columns_with_list_or_dict: 3566 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) 3567 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) 3568 df.fillna(pd.NA, inplace=True) 3569 df.to_csv(filename) 3570 return
30class pymbe_library(): 31 """ 32 The library for the Molecular Builder for ESPResSo (pyMBE) 33 34 Attributes: 35 N_A(`pint.Quantity`): Avogadro number. 36 Kb(`pint.Quantity`): Boltzmann constant. 37 e(`pint.Quantity`): Elementary charge. 38 df(`Pandas.Dataframe`): Dataframe used to bookkeep all the information stored in pyMBE. Typically refered as `pmb.df`. 39 kT(`pint.Quantity`): Thermal energy. 40 Kw(`pint.Quantity`): Ionic product of water. Used in the setup of the G-RxMC method. 41 """ 42 43 class NumpyEncoder(json.JSONEncoder): 44 """ 45 Custom JSON encoder that converts NumPy arrays to Python lists 46 and NumPy scalars to Python scalars. 47 """ 48 def default(self, obj): 49 if isinstance(obj, np.ndarray): 50 return obj.tolist() 51 if isinstance(obj, np.generic): 52 return obj.item() 53 return super().default(obj) 54 55 def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None): 56 """ 57 Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 58 and sets up the `pmb.df` for bookkeeping. 59 60 Args: 61 temperature(`pint.Quantity`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. 62 unit_length(`pint.Quantity`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. 63 unit_charge (`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 64 Kw (`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 65 66 Note: 67 - If no `temperature` is given, a value of 298.15 K is assumed by default. 68 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 69 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 70 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 71 """ 72 # Seed and RNG 73 self.seed=seed 74 self.rng = np.random.default_rng(seed) 75 self.units=pint.UnitRegistry() 76 self.N_A=scipy.constants.N_A / self.units.mol 77 self.kB=scipy.constants.k * self.units.J / self.units.K 78 self.e=scipy.constants.e * self.units.C 79 self.set_reduced_units(unit_length=unit_length, 80 unit_charge=unit_charge, 81 temperature=temperature, 82 Kw=Kw) 83 self.setup_df() 84 self.lattice_builder = None 85 self.root = importlib.resources.files(__package__) 86 return 87 88 def _check_if_name_is_defined_in_df(self, name): 89 """ 90 Checks if `name` is defined in `pmb.df`. 91 92 Args: 93 name(`str`): label to check if defined in `pmb.df`. 94 95 Returns: 96 `bool`: `True` for success, `False` otherwise. 97 """ 98 return name in self.df['name'].unique() 99 100 def _check_if_multiple_pmb_types_for_name(self, name, pmb_type_to_be_defined): 101 """ 102 Checks if `name` is defined in `pmb.df` with multiple pmb_types. 103 104 Args: 105 name(`str`): label to check if defined in `pmb.df`. 106 pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`. 107 108 Returns: 109 `bool`: `True` for success, `False` otherwise. 110 """ 111 if name in self.df['name'].unique(): 112 current_object_type = self.df[self.df['name']==name].pmb_type.values[0] 113 if current_object_type != pmb_type_to_be_defined: 114 raise ValueError (f"The name {name} is already defined in the df with a pmb_type = {current_object_type}, pymMBE does not support objects with the same name but different pmb_types") 115 116 117 def _check_supported_molecule(self, molecule_name,valid_pmb_types): 118 """ 119 Checks if the molecule name `molecule_name` is supported by a method of pyMBE. 120 121 Args: 122 molecule_name(`str`): pmb object type to be checked. 123 valid_pmb_types(`list` of `str`): List of valid pmb types supported by the method. 124 125 Returns: 126 pmb_type(`str`): pmb_type of the molecule. 127 """ 128 pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0] 129 if pmb_type not in valid_pmb_types: 130 raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}") 131 return pmb_type 132 133 def _check_if_name_has_right_type(self, name, expected_pmb_type, hard_check=True): 134 """ 135 Checks if `name` is of the expected pmb type. 136 137 Args: 138 name(`str`): label to check if defined in `pmb.df`. 139 expected_pmb_type(`str`): pmb object type corresponding to `name`. 140 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`. 141 142 Returns: 143 `bool`: `True` for success, `False` otherwise. 144 """ 145 pmb_type=self.df.loc[self.df['name']==name].pmb_type.values[0] 146 if pmb_type == expected_pmb_type: 147 return True 148 else: 149 if hard_check: 150 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}") 151 return False 152 153 def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False): 154 """ 155 Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles. 156 157 Args: 158 particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles 159 particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles 160 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False. 161 162 Returns: 163 index(`int`): Row index where the bond information has been added in pmb.df. 164 """ 165 particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0] 166 particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0] 167 168 bond_key = self.find_bond_key(particle_name1=particle_name1, 169 particle_name2=particle_name2, 170 use_default_bond=use_default_bond) 171 if not bond_key: 172 return None 173 self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1) 174 indexs = np.where(self.df['name']==bond_key) 175 index_list = list (indexs[0]) 176 used_bond_df = self.df.loc[self.df['particle_id2'].notnull()] 177 #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict 178 used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 ) 179 used_bond_index = used_bond_df.index.to_list() 180 if not index_list: 181 return None 182 for index in index_list: 183 if index not in used_bond_index: 184 self.clean_df_row(index=int(index)) 185 self.df.at[index,'particle_id'] = particle_id1 186 self.df.at[index,'particle_id2'] = particle_id2 187 break 188 return index 189 190 def add_bonds_to_espresso(self, espresso_system) : 191 """ 192 Adds all bonds defined in `pmb.df` to `espresso_system`. 193 194 Args: 195 espresso_system(`espressomd.system.System`): system object of espressomd library 196 """ 197 198 if 'bond' in self.df["pmb_type"].values: 199 bond_df = self.df.loc[self.df ['pmb_type'] == 'bond'] 200 bond_list = bond_df.bond_object.values.tolist() 201 for bond in bond_list: 202 espresso_system.bonded_inter.add(bond) 203 else: 204 logging.warning('there are no bonds defined in pymbe.df') 205 return 206 207 def add_value_to_df(self,index,key,new_value, non_standard_value=False, overwrite=False): 208 """ 209 Adds a value to a cell in the `pmb.df` DataFrame. 210 211 Args: 212 index(`int`): index of the row to add the value to. 213 key(`str`): the column label to add the value to. 214 non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False. 215 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 216 """ 217 218 token = "#protected:" 219 220 def protect(obj): 221 if non_standard_value: 222 return token + json.dumps(obj, cls=self.NumpyEncoder) 223 return obj 224 225 def deprotect(obj): 226 if non_standard_value and isinstance(obj, str) and obj.startswith(token): 227 return json.loads(obj.removeprefix(token)) 228 return obj 229 230 # Make sure index is a scalar integer value 231 index = int(index) 232 assert isinstance(index, int), '`index` should be a scalar integer value.' 233 idx = pd.IndexSlice 234 if self.check_if_df_cell_has_a_value(index=index,key=key): 235 old_value = self.df.loc[index,idx[key]] 236 if not pd.Series([protect(old_value)]).equals(pd.Series([protect(new_value)])): 237 name=self.df.loc[index,('name','')] 238 pmb_type=self.df.loc[index,('pmb_type','')] 239 logging.debug(f"You are attempting to redefine the properties of {name} of pmb_type {pmb_type}") 240 if overwrite: 241 logging.info(f'Overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}') 242 if not overwrite: 243 logging.debug(f"pyMBE has preserved of the entry `{key}`: old_value = {old_value}. If you want to overwrite it with new_value = {new_value}, activate the switch overwrite = True ") 244 return 245 246 self.df.loc[index,idx[key]] = protect(new_value) 247 if non_standard_value: 248 self.df[key] = self.df[key].apply(deprotect) 249 return 250 251 def assign_molecule_id(self, molecule_index): 252 """ 253 Assigns the `molecule_id` of the pmb object given by `pmb_type` 254 255 Args: 256 molecule_index(`int`): index of the current `pmb_object_type` to assign the `molecule_id` 257 Returns: 258 molecule_id(`int`): Id of the molecule 259 """ 260 261 self.clean_df_row(index=int(molecule_index)) 262 263 if self.df['molecule_id'].isnull().values.all(): 264 molecule_id = 0 265 else: 266 molecule_id = self.df['molecule_id'].max() +1 267 268 self.add_value_to_df (key=('molecule_id',''), 269 index=int(molecule_index), 270 new_value=molecule_id) 271 272 return molecule_id 273 274 def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system): 275 """ 276 Calculates the center of the molecule with a given molecule_id. 277 278 Args: 279 molecule_id(`int`): Id of the molecule whose center of mass is to be calculated. 280 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 281 282 Returns: 283 center_of_mass(`lst`): Coordinates of the center of mass. 284 """ 285 center_of_mass = np.zeros(3) 286 axis_list = [0,1,2] 287 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0] 288 particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"] 289 for pid in particle_id_list: 290 for axis in axis_list: 291 center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis] 292 center_of_mass = center_of_mass / len(particle_id_list) 293 return center_of_mass 294 295 def calculate_HH(self, molecule_name, pH_list=None, pka_set=None): 296 """ 297 Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 298 for molecules with the name `molecule_name`. 299 300 Args: 301 molecule_name(`str`): name of the molecule to calculate the ideal charge for 302 pH_list(`lst`): pH-values to calculate. 303 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 304 305 Returns: 306 Z_HH(`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list` 307 308 Note: 309 - This function supports objects with pmb types: "molecule", "peptide" and "protein". 310 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 311 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 312 - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan` should be used. 313 """ 314 self._check_if_name_is_defined_in_df(name=molecule_name) 315 self._check_supported_molecule(molecule_name=molecule_name, 316 valid_pmb_types=["molecule","peptide","protein"]) 317 if pH_list is None: 318 pH_list=np.linspace(2,12,50) 319 if pka_set is None: 320 pka_set=self.get_pka_set() 321 index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 322 residue_list = self.df.at [index,('residue_list','')].copy() 323 particles_in_molecule = [] 324 for residue in residue_list: 325 list_of_particles_in_residue = self.search_particles_in_residue(residue) 326 if len(list_of_particles_in_residue) == 0: 327 logging.warning(f"The residue {residue} has no particles defined in the pyMBE DataFrame, it will be ignored.") 328 continue 329 particles_in_molecule += list_of_particles_in_residue 330 if len(particles_in_molecule) == 0: 331 return [None]*len(pH_list) 332 self.check_pka_set(pka_set=pka_set) 333 charge_number_map = self.get_charge_number_map() 334 Z_HH=[] 335 for pH_value in pH_list: 336 Z=0 337 for particle in particles_in_molecule: 338 if particle in pka_set.keys(): 339 if pka_set[particle]['acidity'] == 'acidic': 340 psi=-1 341 elif pka_set[particle]['acidity']== 'basic': 342 psi=+1 343 Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value']))) 344 else: 345 state_one_type = self.df.loc[self.df['name']==particle].state_one.es_type.values[0] 346 Z+=charge_number_map[state_one_type] 347 Z_HH.append(Z) 348 return Z_HH 349 350 def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None): 351 """ 352 Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning. 353 354 Args: 355 c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 356 c_salt('float'): Salt concentration in the reservoir. 357 pH_list('lst'): List of pH-values in the reservoir. 358 pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}. 359 360 Returns: 361 {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 362 pH_system_list ('lst'): List of pH_values in the system. 363 partition_coefficients_list ('lst'): List of partition coefficients of cations. 364 365 Note: 366 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 367 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 368 - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result. 369 """ 370 if pH_list is None: 371 pH_list=np.linspace(2,12,50) 372 if pka_set is None: 373 pka_set=self.get_pka_set() 374 self.check_pka_set(pka_set=pka_set) 375 376 partition_coefficients_list = [] 377 pH_system_list = [] 378 Z_HH_Donnan={} 379 for key in c_macro: 380 Z_HH_Donnan[key] = [] 381 382 def calc_charges(c_macro, pH): 383 """ 384 Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation. 385 386 Args: 387 c_macro('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 388 pH('float'): pH-value that is used in the HH equation. 389 390 Returns: 391 charge('dict'): {"molecule_name": charge} 392 """ 393 charge = {} 394 for name in c_macro: 395 charge[name] = self.calculate_HH(name, [pH], pka_set)[0] 396 return charge 397 398 def calc_partition_coefficient(charge, c_macro): 399 """ 400 Calculates the partition coefficients of positive ions according to the ideal Donnan theory. 401 402 Args: 403 charge('dict'): {"molecule_name": charge} 404 c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 405 """ 406 nonlocal ionic_strength_res 407 charge_density = 0.0 408 for key in charge: 409 charge_density += charge[key] * c_macro[key] 410 return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude 411 412 for pH_value in pH_list: 413 # calculate the ionic strength of the reservoir 414 if pH_value <= 7.0: 415 ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 416 elif pH_value > 7.0: 417 ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt 418 419 #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations 420 #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation. 421 #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects. 422 equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro)) 423 logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq") 424 partition_coefficient = 10**logxi.root 425 426 charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient)) 427 for key in c_macro: 428 Z_HH_Donnan[key].append(charges_temp[key]) 429 430 pH_system_list.append(pH_value - np.log10(partition_coefficient)) 431 partition_coefficients_list.append(partition_coefficient) 432 433 return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 434 435 def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset): 436 """ 437 Calculates the initial bond length that is used when setting up molecules, 438 based on the minimum of the sum of bonded and short-range (LJ) interactions. 439 440 Args: 441 bond_object(`espressomd.interactions.BondedInteractions`): instance of a bond object from espressomd library 442 bond_type(`str`): label identifying the used bonded potential 443 epsilon(`pint.Quantity`): LJ epsilon of the interaction between the particles 444 sigma(`pint.Quantity`): LJ sigma of the interaction between the particles 445 cutoff(`pint.Quantity`): cutoff-radius of the LJ interaction 446 offset(`pint.Quantity`): offset of the LJ interaction 447 """ 448 def truncated_lj_potential(x, epsilon, sigma, cutoff,offset): 449 if x>cutoff: 450 return 0.0 451 else: 452 return 4*epsilon*((sigma/(x-offset))**12-(sigma/(x-offset))**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6) 453 454 epsilon_red=epsilon.to('reduced_energy').magnitude 455 sigma_red=sigma.to('reduced_length').magnitude 456 cutoff_red=cutoff.to('reduced_length').magnitude 457 offset_red=offset.to('reduced_length').magnitude 458 459 if bond_type == "harmonic": 460 r_0 = bond_object.params.get('r_0') 461 k = bond_object.params.get('k') 462 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 463 elif bond_type == "FENE": 464 r_0 = bond_object.params.get('r_0') 465 k = bond_object.params.get('k') 466 d_r_max = bond_object.params.get('d_r_max') 467 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 468 return l0 469 470 def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False): 471 ''' 472 Calculates the net charge per molecule of molecules with `name` = molecule_name. 473 Returns the net charge per molecule and a maps with the net charge per residue and molecule. 474 475 Args: 476 espresso_system(`espressomd.system.System`): system information 477 molecule_name(`str`): name of the molecule to calculate the net charge 478 dimensionless(`bool'): sets if the charge is returned with a dimension or not 479 480 Returns: 481 (`dict`): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }} 482 483 Note: 484 - The net charge of the molecule is averaged over all molecules of type `name` 485 - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name` 486 ''' 487 self._check_supported_molecule(molecule_name=molecule_name, 488 valid_pmb_types=["molecule","protein","peptide"]) 489 490 id_map = self.get_particle_id_map(object_name=molecule_name) 491 def create_charge_map(espresso_system,id_map,label): 492 charge_number_map={} 493 for super_id in id_map[label].keys(): 494 if dimensionless: 495 net_charge=0 496 else: 497 net_charge=0 * self.units.Quantity(1,'reduced_charge') 498 for pid in id_map[label][super_id]: 499 if dimensionless: 500 net_charge+=espresso_system.part.by_id(pid).q 501 else: 502 net_charge+=espresso_system.part.by_id(pid).q * self.units.Quantity(1,'reduced_charge') 503 charge_number_map[super_id]=net_charge 504 return charge_number_map 505 net_charge_molecules=create_charge_map(label="molecule_map", 506 espresso_system=espresso_system, 507 id_map=id_map) 508 net_charge_residues=create_charge_map(label="residue_map", 509 espresso_system=espresso_system, 510 id_map=id_map) 511 if dimensionless: 512 mean_charge=np.mean(np.array(list(net_charge_molecules.values()))) 513 else: 514 mean_charge=np.mean(np.array([value.magnitude for value in net_charge_molecules.values()]))*self.units.Quantity(1,'reduced_charge') 515 return {"mean": mean_charge, "molecules": net_charge_molecules, "residues": net_charge_residues} 516 517 def center_molecule_in_simulation_box(self, molecule_id, espresso_system): 518 """ 519 Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`. 520 521 Args: 522 molecule_id(`int`): Id of the molecule to be centered. 523 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 524 """ 525 if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0: 526 raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.") 527 center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system) 528 box_center = [espresso_system.box_l[0]/2.0, 529 espresso_system.box_l[1]/2.0, 530 espresso_system.box_l[2]/2.0] 531 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0] 532 particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"] 533 for pid in particle_id_list: 534 es_pos = espresso_system.part.by_id(pid).pos 535 espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center 536 return 537 538 def check_aminoacid_key(self, key): 539 """ 540 Checks if `key` corresponds to a valid aminoacid letter code. 541 542 Args: 543 key(`str`): key to be checked. 544 545 Returns: 546 `bool`: True if `key` is a valid aminoacid letter code, False otherwise. 547 """ 548 valid_AA_keys=['V', #'VAL' 549 'I', #'ILE' 550 'L', #'LEU' 551 'E', #'GLU' 552 'Q', #'GLN' 553 'D', #'ASP' 554 'N', #'ASN' 555 'H', #'HIS' 556 'W', #'TRP' 557 'F', #'PHE' 558 'Y', #'TYR' 559 'R', #'ARG' 560 'K', #'LYS' 561 'S', #'SER' 562 'T', #'THR' 563 'M', #'MET' 564 'A', #'ALA' 565 'G', #'GLY' 566 'P', #'PRO' 567 'C'] #'CYS' 568 if key in valid_AA_keys: 569 return True 570 else: 571 return False 572 573 def check_dimensionality(self, variable, expected_dimensionality): 574 """ 575 Checks if the dimensionality of `variable` matches `expected_dimensionality`. 576 577 Args: 578 variable(`pint.Quantity`): Quantity to be checked. 579 expected_dimensionality(`str`): Expected dimension of the variable. 580 581 Returns: 582 (`bool`): `True` if the variable if of the expected dimensionality, `False` otherwise. 583 584 Note: 585 - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality). 586 - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"` 587 """ 588 correct_dimensionality=variable.check(f"{expected_dimensionality}") 589 if not correct_dimensionality: 590 raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}") 591 return correct_dimensionality 592 593 def check_if_df_cell_has_a_value(self, index,key): 594 """ 595 Checks if a cell in the `pmb.df` at the specified index and column has a value. 596 597 Args: 598 index(`int`): Index of the row to check. 599 key(`str`): Column label to check. 600 601 Returns: 602 `bool`: `True` if the cell has a value, `False` otherwise. 603 """ 604 idx = pd.IndexSlice 605 import warnings 606 with warnings.catch_warnings(): 607 warnings.simplefilter("ignore") 608 return not pd.isna(self.df.loc[index, idx[key]]) 609 610 611 612 def check_if_metal_ion(self,key): 613 """ 614 Checks if `key` corresponds to a label of a supported metal ion. 615 616 Args: 617 key(`str`): key to be checked 618 619 Returns: 620 (`bool`): True if `key` is a supported metal ion, False otherwise. 621 """ 622 if key in self.get_metal_ions_charge_number_map().keys(): 623 return True 624 else: 625 return False 626 627 def check_pka_set(self, pka_set): 628 """ 629 Checks that `pka_set` has the formatting expected by the pyMBE library. 630 631 Args: 632 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 633 """ 634 required_keys=['pka_value','acidity'] 635 for required_key in required_keys: 636 for pka_name, pka_entry in pka_set.items(): 637 if required_key not in pka_entry: 638 raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")') 639 return 640 641 def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")): 642 """ 643 Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a pd.NA value. 644 645 Args: 646 index(`int`): Index of the row to clean. 647 columns_keys_to_clean(`list` of `str`, optional): List with the column keys to be cleaned. Defaults to [`particle_id`, `particle_id2`, `residue_id`, `molecule_id`]. 648 """ 649 for column_key in columns_keys_to_clean: 650 self.add_value_to_df(key=(column_key,''),index=index,new_value=pd.NA) 651 self.df.fillna(pd.NA, inplace=True) 652 return 653 654 def convert_columns_to_original_format (self, df): 655 """ 656 Converts the columns of the Dataframe to the original format in pyMBE. 657 658 Args: 659 df(`DataFrame`): dataframe with pyMBE information as a string 660 661 """ 662 663 columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', ('state_one','es_type'),('state_two','es_type'),('state_one','z'),('state_two','z') ] 664 665 columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset'] 666 667 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence', 'chain_map', 'node_map'] 668 669 for column_name in columns_dtype_int: 670 df[column_name] = df[column_name].astype(pd.Int64Dtype()) 671 672 for column_name in columns_with_list_or_dict: 673 if df[column_name].isnull().all(): 674 df[column_name] = df[column_name].astype(object) 675 else: 676 df[column_name] = df[column_name].apply(lambda x: json.loads(x) if pd.notnull(x) else x) 677 678 for column_name in columns_with_units: 679 df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x) 680 681 df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x) 682 df["l0"] = df["l0"].astype(object) 683 df["pka"] = df["pka"].astype(object) 684 return df 685 686 def convert_str_to_bond_object (self, bond_str): 687 """ 688 Convert a row read as a `str` to the corresponding ESPResSo bond object. 689 690 Args: 691 bond_str(`str`): string with the information of a bond object. 692 693 Returns: 694 bond_object(`obj`): ESPResSo bond object. 695 696 Note: 697 Current supported bonds are: HarmonicBond and FeneBond 698 """ 699 import espressomd.interactions 700 701 supported_bonds = ['HarmonicBond', 'FeneBond'] 702 m = re.search(r'^([A-Za-z0-9_]+)\((\{.+\})\)$', bond_str) 703 if m is None: 704 raise ValueError(f'Cannot parse bond "{bond_str}"') 705 bond = m.group(1) 706 if bond not in supported_bonds: 707 raise NotImplementedError(f"Bond type '{bond}' currently not implemented in pyMBE, accepted types are {supported_bonds}") 708 params = json.loads(m.group(2)) 709 bond_id = params.pop("bond_id") 710 bond_object = getattr(espressomd.interactions, bond)(**params) 711 bond_object._bond_id = bond_id 712 return bond_object 713 714 def copy_df_entry(self, name, column_name, number_of_copies): 715 ''' 716 Creates 'number_of_copies' of a given 'name' in `pymbe.df`. 717 718 Args: 719 name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df` 720 column_name(`str`): Column name to use as a filter. 721 number_of_copies(`int`): number of copies of `name` to be created. 722 723 Note: 724 - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 725 ''' 726 727 valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ] 728 if column_name not in valid_column_names: 729 raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}") 730 df_by_name = self.df.loc[self.df.name == name] 731 if number_of_copies != 1: 732 if df_by_name[column_name].isnull().values.any(): 733 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True) 734 else: 735 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 736 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 737 df_by_name_repeated[column_name] = pd.NA 738 # Concatenate the new particle rows to `df` 739 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 740 else: 741 if not df_by_name[column_name].isnull().values.any(): 742 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 743 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 744 df_by_name_repeated[column_name] = pd.NA 745 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 746 return 747 748 def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt): 749 """ 750 Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`. 751 752 Args: 753 espresso_system(`espressomd.system.System`): instance of an espresso system object. 754 cation_name(`str`): `name` of a particle with a positive charge. 755 anion_name(`str`): `name` of a particle with a negative charge. 756 c_salt(`float`): Salt concentration. 757 758 Returns: 759 c_salt_calculated(`float`): Calculated salt concentration added to `espresso_system`. 760 """ 761 for name in [cation_name, anion_name]: 762 if not self._check_if_name_is_defined_in_df(name=name): 763 logging.warning(f"Object with name '{name}' is not defined in the DataFrame, no ions will be created.") 764 return 765 self._check_if_name_has_right_type(name=cation_name, 766 expected_pmb_type="particle") 767 self._check_if_name_has_right_type(name=anion_name, 768 expected_pmb_type="particle") 769 cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0] 770 anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0] 771 if cation_name_charge <= 0: 772 raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge) 773 if anion_name_charge >= 0: 774 raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge) 775 # Calculate the number of ions in the simulation box 776 volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3') 777 if c_salt.check('[substance] [length]**-3'): 778 N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude) 779 c_salt_calculated=N_ions/(volume*self.N_A) 780 elif c_salt.check('[length]**-3'): 781 N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude) 782 c_salt_calculated=N_ions/volume 783 else: 784 raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt) 785 N_cation = N_ions*abs(anion_name_charge) 786 N_anion = N_ions*abs(cation_name_charge) 787 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation) 788 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion) 789 if c_salt_calculated.check('[substance] [length]**-3'): 790 logging.info(f"added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions") 791 elif c_salt_calculated.check('[length]**-3'): 792 logging.info(f"added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions") 793 return c_salt_calculated 794 795 def create_bond_in_espresso(self, bond_type, bond_parameters): 796 ''' 797 Creates either a harmonic or a FENE bond in ESPResSo 798 799 Args: 800 bond_type(`str`): label to identify the potential to model the bond. 801 bond_parameters(`dict`): parameters of the potential of the bond. 802 803 Note: 804 Currently, only HARMONIC and FENE bonds are supported. 805 806 For a HARMONIC bond the dictionary must contain: 807 808 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 809 using the `pmb.units` UnitRegistry. 810 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 811 the `pmb.units` UnitRegistry. 812 813 For a FENE bond the dictionary must additionally contain: 814 815 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 816 units of length using the `pmb.units` UnitRegistry. Default 'None'. 817 818 Returns: 819 bond_object (`obj`): an ESPResSo bond object 820 ''' 821 from espressomd import interactions 822 823 valid_bond_types = ["harmonic", "FENE"] 824 825 if 'k' in bond_parameters: 826 bond_magnitude = bond_parameters['k'].to('reduced_energy / reduced_length**2') 827 else: 828 raise ValueError("Magnitude of the potential (k) is missing") 829 830 if bond_type == 'harmonic': 831 if 'r_0' in bond_parameters: 832 bond_length = bond_parameters['r_0'].to('reduced_length') 833 else: 834 raise ValueError("Equilibrium bond length (r_0) is missing") 835 bond_object = interactions.HarmonicBond(k = bond_magnitude.magnitude, 836 r_0 = bond_length.magnitude) 837 elif bond_type == 'FENE': 838 if 'r_0' in bond_parameters: 839 bond_length = bond_parameters['r_0'].to('reduced_length').magnitude 840 else: 841 logging.warning("no value provided for r_0. Defaulting to r_0 = 0") 842 bond_length=0 843 if 'd_r_max' in bond_parameters: 844 max_bond_stret = bond_parameters['d_r_max'].to('reduced_length') 845 else: 846 raise ValueError("Maximal stretching length (d_r_max) is missing") 847 bond_object = interactions.FeneBond(r_0 = bond_length, 848 k = bond_magnitude.magnitude, 849 d_r_max = max_bond_stret.magnitude) 850 else: 851 raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}") 852 return bond_object 853 854 855 def create_counterions(self, object_name, cation_name, anion_name, espresso_system): 856 """ 857 Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`. 858 859 Args: 860 object_name(`str`): `name` of a pymbe object. 861 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 862 cation_name(`str`): `name` of a particle with a positive charge. 863 anion_name(`str`): `name` of a particle with a negative charge. 864 865 Returns: 866 counterion_number(`dict`): {"name": number} 867 868 Note: 869 This function currently does not support the creation of counterions for hydrogels. 870 """ 871 for name in [object_name, cation_name, anion_name]: 872 if not self._check_if_name_is_defined_in_df(name=name): 873 logging.warning(f"Object with name '{name}' is not defined in the DataFrame, no counterions will be created.") 874 return 875 for name in [cation_name, anion_name]: 876 self._check_if_name_has_right_type(name=name, expected_pmb_type="particle") 877 self._check_supported_molecule(molecule_name=object_name, 878 valid_pmb_types=["molecule","peptide","protein"]) 879 880 881 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0] 882 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0] 883 object_ids = self.get_particle_id_map(object_name=object_name)["all"] 884 counterion_number={} 885 object_charge={} 886 for name in ['positive', 'negative']: 887 object_charge[name]=0 888 for id in object_ids: 889 if espresso_system.part.by_id(id).q > 0: 890 object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q )) 891 elif espresso_system.part.by_id(id).q < 0: 892 object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q )) 893 if object_charge['positive'] % abs(anion_charge) == 0: 894 counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge)) 895 else: 896 raise ValueError('The number of positive charges in the pmb_object must be divisible by the charge of the anion') 897 if object_charge['negative'] % abs(cation_charge) == 0: 898 counterion_number[cation_name]=int(object_charge['negative']/cation_charge) 899 else: 900 raise ValueError('The number of negative charges in the pmb_object must be divisible by the charge of the cation') 901 if counterion_number[cation_name] > 0: 902 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name]) 903 else: 904 counterion_number[cation_name]=0 905 if counterion_number[anion_name] > 0: 906 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name]) 907 else: 908 counterion_number[anion_name] = 0 909 logging.info('the following counter-ions have been created: ') 910 for name in counterion_number.keys(): 911 logging.info(f'Ion type: {name} created number: {counterion_number[name]}') 912 return counterion_number 913 914 def create_hydrogel(self, name, espresso_system): 915 """ 916 creates the hydrogel `name` in espresso_system 917 Args: 918 name(`str`): Label of the hydrogel to be created. `name` must be defined in the `pmb.df` 919 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 920 921 Returns: 922 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],...}} 923 """ 924 if not self._check_if_name_is_defined_in_df(name=name): 925 logging.warning(f"Hydrogel with name '{name}' is not defined in the DataFrame, no hydrogel will be created.") 926 return 927 self._check_if_name_has_right_type(name=name, 928 expected_pmb_type="hydrogel") 929 hydrogel_info={"name":name, "chains":{}, "nodes":{}} 930 # placing nodes 931 node_positions = {} 932 node_topology = self.df[self.df["name"]==name]["node_map"].iloc[0] 933 for node_info in node_topology: 934 node_index = node_info["lattice_index"] 935 node_name = node_info["particle_name"] 936 node_pos, node_id = self.create_hydrogel_node(self.format_node(node_index), node_name, espresso_system) 937 hydrogel_info["nodes"][self.format_node(node_index)]=node_id 938 node_positions[node_id[0]]=node_pos 939 940 # Placing chains between nodes 941 # Looping over all the 16 chains 942 chain_topology = self.df[self.df["name"]==name]["chain_map"].iloc[0] 943 for chain_info in chain_topology: 944 node_s = chain_info["node_start"] 945 node_e = chain_info["node_end"] 946 molecule_info = self.create_hydrogel_chain(node_s, node_e, node_positions, espresso_system) 947 for molecule_id in molecule_info: 948 hydrogel_info["chains"][molecule_id] = molecule_info[molecule_id] 949 hydrogel_info["chains"][molecule_id]["node_start"]=node_s 950 hydrogel_info["chains"][molecule_id]["node_end"]=node_e 951 return hydrogel_info 952 953 def create_hydrogel_chain(self, node_start, node_end, node_positions, espresso_system): 954 """ 955 Creates a chain between two nodes of a hydrogel. 956 957 Args: 958 node_start(`str`): name of the starting node particle at which the first residue of the chain will be attached. 959 node_end(`str`): name of the ending node particle at which the last residue of the chain will be attached. 960 node_positions(`dict`): dictionary with the positions of the nodes. The keys are the node names and the values are the positions. 961 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 962 963 Note: 964 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. 965 The chain will be placed in the direction of the vector between `node_start` and `node_end`. 966 """ 967 if self.lattice_builder is None: 968 raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.") 969 970 molecule_name = "chain_"+node_start+"_"+node_end 971 sequence = self.df[self.df['name']==molecule_name].residue_list.values [0] 972 assert len(sequence) != 0 and not isinstance(sequence, str) 973 assert len(sequence) == self.lattice_builder.mpc 974 975 key, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end) 976 assert node_start != node_end or sequence == sequence[::-1], \ 977 (f"chain cannot be defined between '{node_start}' and '{node_end}' since it " 978 "would form a loop with a non-symmetric sequence (under-defined stereocenter)") 979 980 if reverse: 981 sequence = sequence[::-1] 982 983 node_start_pos = np.array(list(int(x) for x in node_start.strip('[]').split()))*0.25*self.lattice_builder.box_l 984 node_end_pos = np.array(list(int(x) for x in node_end.strip('[]').split()))*0.25*self.lattice_builder.box_l 985 node1 = espresso_system.part.select(lambda p: (p.pos == node_start_pos).all()).id 986 node2 = espresso_system.part.select(lambda p: (p.pos == node_end_pos).all()).id 987 988 if not node1[0] in node_positions or not node2[0] in node_positions: 989 raise ValueError("Set node position before placing a chain between them") 990 991 # Finding a backbone vector between node_start and node_end 992 vec_between_nodes = np.array(node_positions[node2[0]]) - np.array(node_positions[node1[0]]) 993 vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l) 994 backbone_vector = list(vec_between_nodes/(self.lattice_builder.mpc + 1)) 995 node_start_name = self.df[(self.df["particle_id"]==node1[0]) & (self.df["pmb_type"]=="particle")]["name"].values[0] 996 first_res_name = self.df[(self.df["pmb_type"]=="residue") & (self.df["name"]==sequence[0])]["central_bead"].values[0] 997 l0 = self.get_bond_length(node_start_name, first_res_name, hard_check=True) 998 chain_molecule_info = self.create_molecule( 999 name=molecule_name, # Use the name defined earlier 1000 number_of_molecules=1, # Creating one chain 1001 espresso_system=espresso_system, 1002 list_of_first_residue_positions=[list(np.array(node_positions[node1[0]]) + np.array(backbone_vector))],#Start at the first node 1003 backbone_vector=np.array(backbone_vector)/l0, 1004 use_default_bond=False # Use defaut bonds between monomers 1005 ) 1006 # Collecting ids of beads of the chain/molecule 1007 chain_ids = [] 1008 residue_ids = [] 1009 for molecule_id in chain_molecule_info: 1010 for residue_id in chain_molecule_info[molecule_id]: 1011 residue_ids.append(residue_id) 1012 bead_id = chain_molecule_info[molecule_id][residue_id]['central_bead_id'] 1013 chain_ids.append(bead_id) 1014 1015 self.lattice_builder.chains[key] = sequence 1016 # Search bonds between nodes and chain ends 1017 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] 1018 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] 1019 bond_node1_first_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_start], 1020 particle_name2 = BeadType_near_to_node_start, 1021 hard_check=False, 1022 use_default_bond=False) 1023 bond_node2_last_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_end], 1024 particle_name2 = BeadType_near_to_node_end, 1025 hard_check=False, 1026 use_default_bond=False) 1027 1028 espresso_system.part.by_id(node1[0]).add_bond((bond_node1_first_monomer, chain_ids[0])) 1029 espresso_system.part.by_id(node2[0]).add_bond((bond_node2_last_monomer, chain_ids[-1])) 1030 # Add bonds to data frame 1031 bond_index1 = self.add_bond_in_df(particle_id1=node1[0], 1032 particle_id2=chain_ids[0], 1033 use_default_bond=False) 1034 self.add_value_to_df(key=('molecule_id',''), 1035 index=int(bond_index1), 1036 new_value=molecule_id, 1037 overwrite=True) 1038 self.add_value_to_df(key=('residue_id',''), 1039 index=int (bond_index1), 1040 new_value=residue_ids[0], 1041 overwrite=True) 1042 bond_index2 = self.add_bond_in_df(particle_id1=node2[0], 1043 particle_id2=chain_ids[-1], 1044 use_default_bond=False) 1045 self.add_value_to_df(key=('molecule_id',''), 1046 index=int(bond_index2), 1047 new_value=molecule_id, 1048 overwrite=True) 1049 self.add_value_to_df(key=('residue_id',''), 1050 index=int (bond_index2), 1051 new_value=residue_ids[-1], 1052 overwrite=True) 1053 return chain_molecule_info 1054 1055 def create_hydrogel_node(self, node_index, node_name, espresso_system): 1056 """ 1057 Set a node residue type. 1058 1059 Args: 1060 node_index(`str`): Lattice node index in the form of a string, e.g. "[0 0 0]". 1061 node_name(`str`): name of the node particle defined in pyMBE. 1062 Returns: 1063 node_position(`list`): Position of the node in the lattice. 1064 p_id(`int`): Particle ID of the node. 1065 """ 1066 if self.lattice_builder is None: 1067 raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.") 1068 1069 node_position = np.array(list(int(x) for x in node_index.strip('[]').split()))*0.25*self.lattice_builder.box_l 1070 p_id = self.create_particle(name = node_name, 1071 espresso_system=espresso_system, 1072 number_of_particles=1, 1073 position = [node_position]) 1074 key = self.lattice_builder._get_node_by_label(node_index) 1075 self.lattice_builder.nodes[key] = node_name 1076 1077 return node_position.tolist(), p_id 1078 1079 def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False): 1080 """ 1081 Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1082 1083 Args: 1084 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 1085 espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library. 1086 number_of_molecules(`int`): Number of molecules of type `name` to be created. 1087 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. 1088 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. 1089 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df` 1090 1091 Returns: 1092 molecules_info(`dict`): {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 1093 1094 Note: 1095 Despite its name, this function can be used to create both molecules and peptides. 1096 """ 1097 if not self._check_if_name_is_defined_in_df(name=name): 1098 logging.warning(f"Molecule with name '{name}' is not defined in the pyMBE DataFrame, no molecule will be created.") 1099 return {} 1100 if number_of_molecules <= 0: 1101 return {} 1102 if list_of_first_residue_positions is not None: 1103 for item in list_of_first_residue_positions: 1104 if not isinstance(item, list): 1105 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.") 1106 elif len(item) != 3: 1107 raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.") 1108 1109 if len(list_of_first_residue_positions) != number_of_molecules: 1110 raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}") 1111 1112 # This function works for both molecules and peptides 1113 if not self._check_if_name_has_right_type(name=name, expected_pmb_type="molecule", hard_check=False): 1114 self._check_if_name_has_right_type(name=name, expected_pmb_type="peptide") 1115 1116 # Generate an arbitrary random unit vector 1117 if backbone_vector is None: 1118 backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0], 1119 radius=1, 1120 n_samples=1, 1121 on_surface=True)[0] 1122 else: 1123 backbone_vector = np.array(backbone_vector) 1124 1125 1126 1127 first_residue = True 1128 molecules_info = {} 1129 residue_list = self.df[self.df['name']==name].residue_list.values [0] 1130 self.copy_df_entry(name=name, 1131 column_name='molecule_id', 1132 number_of_copies=number_of_molecules) 1133 1134 molecules_index = np.where(self.df['name']==name) 1135 molecule_index_list =list(molecules_index[0])[-number_of_molecules:] 1136 pos_index = 0 1137 for molecule_index in molecule_index_list: 1138 molecule_id = self.assign_molecule_id(molecule_index=molecule_index) 1139 molecules_info[molecule_id] = {} 1140 for residue in residue_list: 1141 if first_residue: 1142 if list_of_first_residue_positions is None: 1143 residue_position = None 1144 else: 1145 for item in list_of_first_residue_positions: 1146 residue_position = [np.array(list_of_first_residue_positions[pos_index])] 1147 1148 residues_info = self.create_residue(name=residue, 1149 espresso_system=espresso_system, 1150 central_bead_position=residue_position, 1151 use_default_bond= use_default_bond, 1152 backbone_vector=backbone_vector) 1153 residue_id = next(iter(residues_info)) 1154 # Add the correct molecule_id to all particles in the residue 1155 for index in self.df[self.df['residue_id']==residue_id].index: 1156 self.add_value_to_df(key=('molecule_id',''), 1157 index=int (index), 1158 new_value=molecule_id, 1159 overwrite=True) 1160 central_bead_id = residues_info[residue_id]['central_bead_id'] 1161 previous_residue = residue 1162 residue_position = espresso_system.part.by_id(central_bead_id).pos 1163 previous_residue_id = central_bead_id 1164 first_residue = False 1165 else: 1166 previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0] 1167 new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0] 1168 bond = self.search_bond(particle_name1=previous_central_bead_name, 1169 particle_name2=new_central_bead_name, 1170 hard_check=True, 1171 use_default_bond=use_default_bond) 1172 l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 1173 particle_name2=new_central_bead_name, 1174 hard_check=True, 1175 use_default_bond=use_default_bond) 1176 1177 residue_position = residue_position+backbone_vector*l0 1178 residues_info = self.create_residue(name=residue, 1179 espresso_system=espresso_system, 1180 central_bead_position=[residue_position], 1181 use_default_bond= use_default_bond, 1182 backbone_vector=backbone_vector) 1183 residue_id = next(iter(residues_info)) 1184 for index in self.df[self.df['residue_id']==residue_id].index: 1185 self.add_value_to_df(key=('molecule_id',''), 1186 index=int (index), 1187 new_value=molecule_id, 1188 overwrite=True) 1189 central_bead_id = residues_info[residue_id]['central_bead_id'] 1190 espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id)) 1191 bond_index = self.add_bond_in_df(particle_id1=central_bead_id, 1192 particle_id2=previous_residue_id, 1193 use_default_bond=use_default_bond) 1194 self.add_value_to_df(key=('molecule_id',''), 1195 index=int (bond_index), 1196 new_value=molecule_id, 1197 overwrite=True) 1198 previous_residue_id = central_bead_id 1199 previous_residue = residue 1200 molecules_info[molecule_id][residue_id] = residues_info[residue_id] 1201 first_residue = True 1202 pos_index+=1 1203 1204 return molecules_info 1205 1206 def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False): 1207 """ 1208 Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`. 1209 1210 Args: 1211 name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 1212 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 1213 number_of_particles(`int`): Number of particles to be created. 1214 position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None. 1215 fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False. 1216 Returns: 1217 created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`. 1218 """ 1219 if number_of_particles <=0: 1220 return [] 1221 if not self._check_if_name_is_defined_in_df(name=name): 1222 logging.warning(f"Particle with name '{name}' is not defined in the pyMBE DataFrame, no particle will be created.") 1223 return [] 1224 self._check_if_name_has_right_type(name=name, 1225 expected_pmb_type="particle") 1226 # Copy the data of the particle `number_of_particles` times in the `df` 1227 self.copy_df_entry(name=name, 1228 column_name='particle_id', 1229 number_of_copies=number_of_particles) 1230 # Get information from the particle type `name` from the df 1231 z = self.df.loc[self.df['name']==name].state_one.z.values[0] 1232 z = 0. if z is None else z 1233 es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 1234 # Get a list of the index in `df` corresponding to the new particles to be created 1235 index = np.where(self.df['name']==name) 1236 index_list =list(index[0])[-number_of_particles:] 1237 # Create the new particles into `espresso_system` 1238 created_pid_list=[] 1239 for index in range (number_of_particles): 1240 df_index=int (index_list[index]) 1241 self.clean_df_row(index=df_index) 1242 if position is None: 1243 particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l) 1244 else: 1245 particle_position = position[index] 1246 if len(espresso_system.part.all()) == 0: 1247 bead_id = 0 1248 else: 1249 bead_id = max (espresso_system.part.all().id) + 1 1250 created_pid_list.append(bead_id) 1251 kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z) 1252 if fix: 1253 kwargs["fix"] = 3 * [fix] 1254 espresso_system.part.add(**kwargs) 1255 self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id) 1256 return created_pid_list 1257 1258 def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None): 1259 """ 1260 Creates all `particle`s associated to `pmb object` into `espresso` a number of times equal to `number_of_objects`. 1261 1262 Args: 1263 name(`str`): Unique label of the `pmb object` to be created. 1264 number_of_objects(`int`): Number of `pmb object`s to be created. 1265 espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library. 1266 position(`list`): Coordinates where the particles should be created. 1267 use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`. 1268 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. 1269 1270 Note: 1271 - If no `position` is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length. 1272 """ 1273 pmb_type = self._check_supported_molecule(molecule_name=name, 1274 valid_pmb_types=["molecule","particle","residue","peptide"]) 1275 1276 if pmb_type == 'particle': 1277 self.create_particle(name=name, 1278 number_of_particles=number_of_objects, 1279 espresso_system=espresso_system, 1280 position=position) 1281 elif pmb_type == 'residue': 1282 self.create_residue(name=name, 1283 espresso_system=espresso_system, 1284 central_bead_position=position, 1285 use_default_bond=use_default_bond, 1286 backbone_vector=backbone_vector) 1287 elif pmb_type in ['molecule','peptide']: 1288 self.create_molecule(name=name, 1289 number_of_molecules=number_of_objects, 1290 espresso_system=espresso_system, 1291 use_default_bond=use_default_bond, 1292 list_of_first_residue_positions=position, 1293 backbone_vector=backbone_vector) 1294 return 1295 1296 def create_protein(self, name, number_of_proteins, espresso_system, topology_dict): 1297 """ 1298 Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions` 1299 1300 Args: 1301 name(`str`): Label of the protein to be created. 1302 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 1303 number_of_proteins(`int`): Number of proteins to be created. 1304 positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}} 1305 """ 1306 1307 if number_of_proteins <=0: 1308 return 1309 if not self._check_if_name_is_defined_in_df(name=name): 1310 logging.warning(f"Protein with name '{name}' is not defined in the pyMBE DataFrame, no protein will be created.") 1311 return 1312 self._check_if_name_has_right_type(name=name, 1313 expected_pmb_type="protein") 1314 self.copy_df_entry(name=name, 1315 column_name='molecule_id', 1316 number_of_copies=number_of_proteins) 1317 protein_index = np.where(self.df['name']==name) 1318 protein_index_list =list(protein_index[0])[-number_of_proteins:] 1319 1320 box_half=espresso_system.box_l[0]/2.0 1321 1322 for molecule_index in protein_index_list: 1323 1324 molecule_id = self.assign_molecule_id(molecule_index=molecule_index) 1325 1326 protein_center = self.generate_coordinates_outside_sphere(radius = 1, 1327 max_dist=box_half, 1328 n_samples=1, 1329 center=[box_half]*3)[0] 1330 1331 for residue in topology_dict.keys(): 1332 1333 residue_name = re.split(r'\d+', residue)[0] 1334 residue_number = re.split(r'(\d+)', residue)[1] 1335 residue_position = topology_dict[residue]['initial_pos'] 1336 position = residue_position + protein_center 1337 1338 particle_id = self.create_particle(name=residue_name, 1339 espresso_system=espresso_system, 1340 number_of_particles=1, 1341 position=[position], 1342 fix = True) 1343 1344 index = self.df[self.df['particle_id']==particle_id[0]].index.values[0] 1345 self.add_value_to_df(key=('residue_id',''), 1346 index=int (index), 1347 new_value=int(residue_number), 1348 overwrite=True) 1349 1350 self.add_value_to_df(key=('molecule_id',''), 1351 index=int (index), 1352 new_value=molecule_id, 1353 overwrite=True) 1354 1355 return 1356 1357 def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None): 1358 """ 1359 Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1360 1361 Args: 1362 name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df` 1363 espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library. 1364 central_bead_position(`list` of `float`): Position of the central bead. 1365 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` 1366 backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`. 1367 1368 Returns: 1369 residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1370 """ 1371 if not self._check_if_name_is_defined_in_df(name=name): 1372 logging.warning(f"Residue with name '{name}' is not defined in the pyMBE DataFrame, no residue will be created.") 1373 return 1374 self._check_if_name_has_right_type(name=name, 1375 expected_pmb_type="residue") 1376 1377 # Copy the data of a residue in the `df 1378 self.copy_df_entry(name=name, 1379 column_name='residue_id', 1380 number_of_copies=1) 1381 residues_index = np.where(self.df['name']==name) 1382 residue_index_list =list(residues_index[0])[-1:] 1383 # search for defined particle and residue names 1384 particle_and_residue_df = self.df.loc[(self.df['pmb_type']== "particle") | (self.df['pmb_type']== "residue")] 1385 particle_and_residue_names = particle_and_residue_df["name"].tolist() 1386 for residue_index in residue_index_list: 1387 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1388 for side_chain_element in side_chain_list: 1389 if side_chain_element not in particle_and_residue_names: 1390 raise ValueError (f"{side_chain_element} is not defined") 1391 # Internal bookkepping of the residue info (important for side-chain residues) 1392 # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1393 residues_info={} 1394 for residue_index in residue_index_list: 1395 self.clean_df_row(index=int(residue_index)) 1396 # Assign a residue_id 1397 if self.df['residue_id'].isnull().all(): 1398 residue_id=0 1399 else: 1400 residue_id = self.df['residue_id'].max() + 1 1401 self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id) 1402 # create the principal bead 1403 central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0] 1404 central_bead_id = self.create_particle(name=central_bead_name, 1405 espresso_system=espresso_system, 1406 position=central_bead_position, 1407 number_of_particles = 1)[0] 1408 central_bead_position=espresso_system.part.by_id(central_bead_id).pos 1409 #assigns same residue_id to the central_bead particle created. 1410 index = self.df[self.df['particle_id']==central_bead_id].index.values[0] 1411 self.df.at [index,'residue_id'] = residue_id 1412 # Internal bookkeeping of the central bead id 1413 residues_info[residue_id]={} 1414 residues_info[residue_id]['central_bead_id']=central_bead_id 1415 # create the lateral beads 1416 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1417 side_chain_beads_ids = [] 1418 for side_chain_element in side_chain_list: 1419 pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 1420 if pmb_type == 'particle': 1421 bond = self.search_bond(particle_name1=central_bead_name, 1422 particle_name2=side_chain_element, 1423 hard_check=True, 1424 use_default_bond=use_default_bond) 1425 l0 = self.get_bond_length(particle_name1=central_bead_name, 1426 particle_name2=side_chain_element, 1427 hard_check=True, 1428 use_default_bond=use_default_bond) 1429 1430 if backbone_vector is None: 1431 bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1432 radius=l0, 1433 n_samples=1, 1434 on_surface=True)[0] 1435 else: 1436 bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector), 1437 magnitude=l0) 1438 1439 side_bead_id = self.create_particle(name=side_chain_element, 1440 espresso_system=espresso_system, 1441 position=[bead_position], 1442 number_of_particles=1)[0] 1443 index = self.df[self.df['particle_id']==side_bead_id].index.values[0] 1444 self.add_value_to_df(key=('residue_id',''), 1445 index=int (index), 1446 new_value=residue_id, 1447 overwrite=True) 1448 side_chain_beads_ids.append(side_bead_id) 1449 espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id)) 1450 index = self.add_bond_in_df(particle_id1=central_bead_id, 1451 particle_id2=side_bead_id, 1452 use_default_bond=use_default_bond) 1453 self.add_value_to_df(key=('residue_id',''), 1454 index=int (index), 1455 new_value=residue_id, 1456 overwrite=True) 1457 1458 elif pmb_type == 'residue': 1459 central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0] 1460 bond = self.search_bond(particle_name1=central_bead_name, 1461 particle_name2=central_bead_side_chain, 1462 hard_check=True, 1463 use_default_bond=use_default_bond) 1464 l0 = self.get_bond_length(particle_name1=central_bead_name, 1465 particle_name2=central_bead_side_chain, 1466 hard_check=True, 1467 use_default_bond=use_default_bond) 1468 if backbone_vector is None: 1469 residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1470 radius=l0, 1471 n_samples=1, 1472 on_surface=True)[0] 1473 else: 1474 residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1475 magnitude=l0) 1476 lateral_residue_info = self.create_residue(name=side_chain_element, 1477 espresso_system=espresso_system, 1478 central_bead_position=[residue_position], 1479 use_default_bond=use_default_bond) 1480 lateral_residue_dict=list(lateral_residue_info.values())[0] 1481 central_bead_side_chain_id=lateral_residue_dict['central_bead_id'] 1482 lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids'] 1483 residue_id_side_chain=list(lateral_residue_info.keys())[0] 1484 # Change the residue_id of the residue in the side chain to the one of the bigger residue 1485 index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0] 1486 self.add_value_to_df(key=('residue_id',''), 1487 index=int(index), 1488 new_value=residue_id, 1489 overwrite=True) 1490 # Change the residue_id of the particles in the residue in the side chain 1491 side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids 1492 for particle_id in side_chain_beads_ids: 1493 index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0] 1494 self.add_value_to_df(key=('residue_id',''), 1495 index=int (index), 1496 new_value=residue_id, 1497 overwrite=True) 1498 espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id)) 1499 index = self.add_bond_in_df(particle_id1=central_bead_id, 1500 particle_id2=central_bead_side_chain_id, 1501 use_default_bond=use_default_bond) 1502 self.add_value_to_df(key=('residue_id',''), 1503 index=int (index), 1504 new_value=residue_id, 1505 overwrite=True) 1506 # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue 1507 for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index: 1508 self.add_value_to_df(key=('residue_id',''), 1509 index=int(index), 1510 new_value=residue_id, 1511 overwrite=True) 1512 # Internal bookkeeping of the side chain beads ids 1513 residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids 1514 return residues_info 1515 1516 def create_variable_with_units(self, variable): 1517 """ 1518 Returns a pint object with the value and units defined in `variable`. 1519 1520 Args: 1521 variable(`dict` or `str`): {'value': value, 'units': units} 1522 Returns: 1523 variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry. 1524 """ 1525 if isinstance(variable, dict): 1526 value=variable.pop('value') 1527 units=variable.pop('units') 1528 elif isinstance(variable, str): 1529 value = float(re.split(r'\s+', variable)[0]) 1530 units = re.split(r'\s+', variable)[1] 1531 variable_with_units=value*self.units(units) 1532 return variable_with_units 1533 1534 def define_AA_residues(self, sequence, model): 1535 """ 1536 Defines in `pmb.df` all the different residues in `sequence`. 1537 1538 Args: 1539 sequence(`lst`): Sequence of the peptide or protein. 1540 model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1541 1542 Returns: 1543 residue_list(`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1544 """ 1545 residue_list = [] 1546 for residue_name in sequence: 1547 if model == '1beadAA': 1548 central_bead = residue_name 1549 side_chains = [] 1550 elif model == '2beadAA': 1551 if residue_name in ['c','n', 'G']: 1552 central_bead = residue_name 1553 side_chains = [] 1554 else: 1555 central_bead = 'CA' 1556 side_chains = [residue_name] 1557 if residue_name not in residue_list: 1558 self.define_residue(name = 'AA-'+residue_name, 1559 central_bead = central_bead, 1560 side_chains = side_chains) 1561 residue_list.append('AA-'+residue_name) 1562 return residue_list 1563 1564 def define_bond(self, bond_type, bond_parameters, particle_pairs): 1565 1566 ''' 1567 Defines a pmb object of type `bond` in `pymbe.df`. 1568 1569 Args: 1570 bond_type(`str`): label to identify the potential to model the bond. 1571 bond_parameters(`dict`): parameters of the potential of the bond. 1572 particle_pairs(`lst`): list of the `names` of the `particles` to be bonded. 1573 1574 Note: 1575 Currently, only HARMONIC and FENE bonds are supported. 1576 1577 For a HARMONIC bond the dictionary must contain the following parameters: 1578 1579 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 1580 using the `pmb.units` UnitRegistry. 1581 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 1582 the `pmb.units` UnitRegistry. 1583 1584 For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and: 1585 1586 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 1587 units of length using the `pmb.units` UnitRegistry. Default 'None'. 1588 ''' 1589 1590 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1591 for particle_name1, particle_name2 in particle_pairs: 1592 1593 lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1, 1594 particle_name2 = particle_name2, 1595 combining_rule = 'Lorentz-Berthelot') 1596 1597 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1598 bond_type = bond_type, 1599 epsilon = lj_parameters["epsilon"], 1600 sigma = lj_parameters["sigma"], 1601 cutoff = lj_parameters["cutoff"], 1602 offset = lj_parameters["offset"],) 1603 index = len(self.df) 1604 for label in [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']: 1605 self._check_if_multiple_pmb_types_for_name(name=label, pmb_type_to_be_defined="bond") 1606 name=f'{particle_name1}-{particle_name2}' 1607 self._check_if_multiple_pmb_types_for_name(name=name, pmb_type_to_be_defined="bond") 1608 self.df.at [index,'name']= name 1609 self.df.at [index,'bond_object'] = bond_object 1610 self.df.at [index,'l0'] = l0 1611 self.add_value_to_df(index=index, 1612 key=('pmb_type',''), 1613 new_value='bond') 1614 self.add_value_to_df(index=index, 1615 key=('parameters_of_the_potential',''), 1616 new_value=bond_object.get_params(), 1617 non_standard_value=True) 1618 self.df.fillna(pd.NA, inplace=True) 1619 return 1620 1621 1622 def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None): 1623 """ 1624 Asigns `bond` in `pmb.df` as the default bond. 1625 The LJ parameters can be optionally provided to calculate the initial bond length 1626 1627 Args: 1628 bond_type(`str`): label to identify the potential to model the bond. 1629 bond_parameters(`dict`): parameters of the potential of the bond. 1630 sigma(`float`, optional): LJ sigma of the interaction between the particles. 1631 epsilon(`float`, optional): LJ epsilon for the interaction between the particles. 1632 cutoff(`float`, optional): cutoff-radius of the LJ interaction. 1633 offset(`float`, optional): offset of the LJ interaction. 1634 1635 Note: 1636 - Currently, only harmonic and FENE bonds are supported. 1637 """ 1638 1639 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1640 1641 if epsilon is None: 1642 epsilon=1*self.units('reduced_energy') 1643 if sigma is None: 1644 sigma=1*self.units('reduced_length') 1645 if cutoff is None: 1646 cutoff=2**(1.0/6.0)*self.units('reduced_length') 1647 if offset is None: 1648 offset=0*self.units('reduced_length') 1649 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1650 bond_type = bond_type, 1651 epsilon = epsilon, 1652 sigma = sigma, 1653 cutoff = cutoff, 1654 offset = offset) 1655 1656 self._check_if_multiple_pmb_types_for_name(name='default', 1657 pmb_type_to_be_defined='bond') 1658 1659 if len(self.df.index) != 0: 1660 index = max(self.df.index)+1 1661 else: 1662 index = 0 1663 self.df.at [index,'name'] = 'default' 1664 self.df.at [index,'bond_object'] = bond_object 1665 self.df.at [index,'l0'] = l0 1666 self.add_value_to_df(index = index, 1667 key = ('pmb_type',''), 1668 new_value = 'bond') 1669 self.add_value_to_df(index = index, 1670 key = ('parameters_of_the_potential',''), 1671 new_value = bond_object.get_params(), 1672 non_standard_value=True) 1673 self.df.fillna(pd.NA, inplace=True) 1674 return 1675 1676 def define_hydrogel(self, name, node_map, chain_map): 1677 """ 1678 Defines a pyMBE object of type `hydrogel` in `pymbe.df`. 1679 1680 Args: 1681 name(`str`): Unique label that identifies the `hydrogel`. 1682 node_map(`list of ict`): [{"particle_name": , "lattice_index": }, ... ] 1683 chain_map(`list of dict`): [{"node_start": , "node_end": , "residue_list": , ... ] 1684 """ 1685 node_indices = {tuple(entry['lattice_index']) for entry in node_map} 1686 diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices} 1687 if node_indices != diamond_indices: 1688 raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ") 1689 1690 chain_map_connectivity = set() 1691 for entry in chain_map: 1692 start = self.lattice_builder.node_labels[entry['node_start']] 1693 end = self.lattice_builder.node_labels[entry['node_end']] 1694 chain_map_connectivity.add((start,end)) 1695 1696 if self.lattice_builder.lattice.connectivity != chain_map_connectivity: 1697 raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs") 1698 1699 self._check_if_multiple_pmb_types_for_name(name=name, 1700 pmb_type_to_be_defined='hydrogel') 1701 1702 1703 index = len(self.df) 1704 self.df.at [index, "name"] = name 1705 self.df.at [index, "pmb_type"] = "hydrogel" 1706 self.add_value_to_df(index = index, 1707 key = ('node_map',''), 1708 new_value = node_map, 1709 non_standard_value=True) 1710 self.add_value_to_df(index = index, 1711 key = ('chain_map',''), 1712 new_value = chain_map, 1713 non_standard_value=True) 1714 for chain_label in chain_map: 1715 node_start = chain_label["node_start"] 1716 node_end = chain_label["node_end"] 1717 residue_list = chain_label['residue_list'] 1718 # Molecule name 1719 molecule_name = "chain_"+node_start+"_"+node_end 1720 self.define_molecule(name=molecule_name, residue_list=residue_list) 1721 return 1722 1723 def define_molecule(self, name, residue_list): 1724 """ 1725 Defines a pyMBE object of type `molecule` in `pymbe.df`. 1726 1727 Args: 1728 name(`str`): Unique label that identifies the `molecule`. 1729 residue_list(`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1730 """ 1731 self._check_if_multiple_pmb_types_for_name(name=name, 1732 pmb_type_to_be_defined='molecule') 1733 1734 index = len(self.df) 1735 self.df.at [index,'name'] = name 1736 self.df.at [index,'pmb_type'] = 'molecule' 1737 self.df.at [index,('residue_list','')] = residue_list 1738 self.df.fillna(pd.NA, inplace=True) 1739 return 1740 1741 def define_particle_entry_in_df(self,name): 1742 """ 1743 Defines a particle entry in pmb.df. 1744 1745 Args: 1746 name(`str`): Unique label that identifies this particle type. 1747 1748 Returns: 1749 index(`int`): Index of the particle in pmb.df 1750 """ 1751 1752 if self._check_if_name_is_defined_in_df(name=name): 1753 index = self.df[self.df['name']==name].index[0] 1754 else: 1755 index = len(self.df) 1756 self.df.at [index, 'name'] = name 1757 self.df.at [index,'pmb_type'] = 'particle' 1758 self.df.fillna(pd.NA, inplace=True) 1759 return index 1760 1761 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): 1762 """ 1763 Defines the properties of a particle object. 1764 1765 Args: 1766 name(`str`): Unique label that identifies this particle type. 1767 z(`int`, optional): Permanent charge number of this particle type. Defaults to 0. 1768 acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to pd.NA. 1769 pka(`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to pd.NA. 1770 sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1771 cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1772 offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1773 epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA. 1774 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1775 1776 Note: 1777 - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units. 1778 - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units. 1779 - `cutoff` defaults to `2**(1./6.) reduced_length`. 1780 - `offset` defaults to 0. 1781 - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions. 1782 - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`. 1783 """ 1784 index=self.define_particle_entry_in_df(name=name) 1785 self._check_if_multiple_pmb_types_for_name(name=name, 1786 pmb_type_to_be_defined='particle') 1787 1788 # If `cutoff` and `offset` are not defined, default them to the following values 1789 if pd.isna(cutoff): 1790 cutoff=self.units.Quantity(2**(1./6.), "reduced_length") 1791 if pd.isna(offset): 1792 offset=self.units.Quantity(0, "reduced_length") 1793 1794 # Define LJ parameters 1795 parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"}, 1796 "cutoff":{"value": cutoff, "dimensionality": "[length]"}, 1797 "offset":{"value": offset, "dimensionality": "[length]"}, 1798 "epsilon":{"value": epsilon, "dimensionality": "[energy]"},} 1799 1800 for parameter_key in parameters_with_dimensionality.keys(): 1801 if not pd.isna(parameters_with_dimensionality[parameter_key]["value"]): 1802 self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 1803 expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"]) 1804 self.add_value_to_df(key=(parameter_key,''), 1805 index=index, 1806 new_value=parameters_with_dimensionality[parameter_key]["value"], 1807 overwrite=overwrite) 1808 1809 # Define particle acid/base properties 1810 self.set_particle_acidity(name=name, 1811 acidity=acidity, 1812 default_charge_number=z, 1813 pka=pka, 1814 overwrite=overwrite) 1815 self.df.fillna(pd.NA, inplace=True) 1816 return 1817 1818 def define_particles(self, parameters, overwrite=False): 1819 ''' 1820 Defines a particle object in pyMBE for each particle name in `particle_names` 1821 1822 Args: 1823 parameters(`dict`): dictionary with the particle parameters. 1824 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1825 1826 Note: 1827 - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},} 1828 ''' 1829 if not parameters: 1830 return 0 1831 for particle_name in parameters.keys(): 1832 parameters[particle_name]["overwrite"]=overwrite 1833 self.define_particle(**parameters[particle_name]) 1834 return 1835 1836 def define_peptide(self, name, sequence, model): 1837 """ 1838 Defines a pyMBE object of type `peptide` in the `pymbe.df`. 1839 1840 Args: 1841 name (`str`): Unique label that identifies the `peptide`. 1842 sequence (`string`): Sequence of the `peptide`. 1843 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1844 """ 1845 self._check_if_multiple_pmb_types_for_name(name = name, 1846 pmb_type_to_be_defined='peptide') 1847 1848 valid_keys = ['1beadAA','2beadAA'] 1849 if model not in valid_keys: 1850 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1851 clean_sequence = self.protein_sequence_parser(sequence=sequence) 1852 residue_list = self.define_AA_residues(sequence=clean_sequence, 1853 model=model) 1854 self.define_molecule(name = name, residue_list=residue_list) 1855 index = self.df.loc[self.df['name'] == name].index.item() 1856 self.df.at [index,'model'] = model 1857 self.df.at [index,('sequence','')] = clean_sequence 1858 self.df.at [index,'pmb_type'] = "peptide" 1859 self.df.fillna(pd.NA, inplace=True) 1860 1861 1862 def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False): 1863 """ 1864 Defines a globular protein pyMBE object in `pymbe.df`. 1865 1866 Args: 1867 name (`str`): Unique label that identifies the protein. 1868 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1869 topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value} 1870 lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca". 1871 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1872 1873 Note: 1874 - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential. 1875 """ 1876 1877 # Sanity checks 1878 self._check_if_multiple_pmb_types_for_name(name = name, 1879 pmb_type_to_be_defined='protein') 1880 valid_model_keys = ['1beadAA','2beadAA'] 1881 valid_lj_setups = ["wca"] 1882 if model not in valid_model_keys: 1883 raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}') 1884 if lj_setup_mode not in valid_lj_setups: 1885 raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}') 1886 if lj_setup_mode == "wca": 1887 sigma = 1*self.units.Quantity("reduced_length") 1888 epsilon = 1*self.units.Quantity("reduced_energy") 1889 part_dict={} 1890 sequence=[] 1891 metal_ions_charge_number_map=self.get_metal_ions_charge_number_map() 1892 for particle in topology_dict.keys(): 1893 particle_name = re.split(r'\d+', particle)[0] 1894 if particle_name not in part_dict.keys(): 1895 if lj_setup_mode == "wca": 1896 part_dict[particle_name]={"sigma": sigma, 1897 "offset": topology_dict[particle]['radius']*2-sigma, 1898 "epsilon": epsilon, 1899 "name": particle_name} 1900 if self.check_if_metal_ion(key=particle_name): 1901 z=metal_ions_charge_number_map[particle_name] 1902 else: 1903 z=0 1904 part_dict[particle_name]["z"]=z 1905 1906 if self.check_aminoacid_key(key=particle_name): 1907 sequence.append(particle_name) 1908 1909 self.define_particles(parameters=part_dict, 1910 overwrite=overwrite) 1911 residue_list = self.define_AA_residues(sequence=sequence, 1912 model=model) 1913 index = len(self.df) 1914 self.df.at [index,'name'] = name 1915 self.df.at [index,'pmb_type'] = 'protein' 1916 self.df.at [index,'model'] = model 1917 self.df.at [index,('sequence','')] = sequence 1918 self.df.at [index,('residue_list','')] = residue_list 1919 self.df.fillna(pd.NA, inplace=True) 1920 return 1921 1922 def define_residue(self, name, central_bead, side_chains): 1923 """ 1924 Defines a pyMBE object of type `residue` in `pymbe.df`. 1925 1926 Args: 1927 name(`str`): Unique label that identifies the `residue`. 1928 central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`. 1929 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. 1930 """ 1931 self._check_if_multiple_pmb_types_for_name(name=name, 1932 pmb_type_to_be_defined='residue') 1933 1934 index = len(self.df) 1935 self.df.at [index, 'name'] = name 1936 self.df.at [index,'pmb_type'] = 'residue' 1937 self.df.at [index,'central_bead'] = central_bead 1938 self.df.at [index,('side_chains','')] = side_chains 1939 self.df.fillna(pd.NA, inplace=True) 1940 return 1941 1942 def delete_entries_in_df(self, entry_name): 1943 """ 1944 Deletes entries with name `entry_name` from the DataFrame if it exists. 1945 1946 Args: 1947 entry_name (`str`): The name of the entry in the dataframe to delete. 1948 1949 """ 1950 if entry_name in self.df["name"].values: 1951 self.df = self.df[self.df["name"] != entry_name].reset_index(drop=True) 1952 1953 def destroy_pmb_object_in_system(self, name, espresso_system): 1954 """ 1955 Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 1956 1957 Args: 1958 name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`. 1959 espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library. 1960 1961 Note: 1962 - If `name` is a object_type=`particle`, only the matching particles that are not part of bigger objects (e.g. `residue`, `molecule`) will be destroyed. To destroy particles in such objects, destroy the bigger object instead. 1963 """ 1964 pmb_type = self._check_supported_molecule(molecule_name=name, 1965 valid_pmb_types= ['particle','residue','molecule',"peptide"]) 1966 1967 if pmb_type == 'particle': 1968 particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1969 & (self.df['molecule_id'].isna())] 1970 particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1971 & (self.df['molecule_id'].isna())].particle_id.tolist() 1972 for particle_id in particle_ids_list: 1973 espresso_system.part.by_id(particle_id).remove() 1974 self.df = self.df.drop(index=particles_index) 1975 if pmb_type == 'residue': 1976 residues_id = self.df.loc[self.df['name']== name].residue_id.to_list() 1977 for residue_id in residues_id: 1978 molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0] 1979 particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"] 1980 self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index) 1981 for particle_id in particle_ids_list: 1982 espresso_system.part.by_id(particle_id).remove() 1983 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1984 if pmb_type in ['molecule',"peptide"]: 1985 molecules_id = self.df.loc[self.df['name']== name].molecule_id.dropna().to_list() 1986 for molecule_id in molecules_id: 1987 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']==pmb_type)].name.values[0] 1988 particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"] 1989 self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index) 1990 for particle_id in particle_ids_list: 1991 espresso_system.part.by_id(particle_id).remove() 1992 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1993 1994 self.df.reset_index(drop=True,inplace=True) 1995 1996 return 1997 1998 def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200): 1999 """ 2000 Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration. 2001 To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an 2002 ion pair and the concentrations of the various species is solved numerically using a self-consistent approach. 2003 More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831). 2004 This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237). 2005 2006 Args: 2007 pH_res('float'): pH-value in the reservoir. 2008 c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2009 activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2010 2011 Returns: 2012 cH_res('pint.Quantity'): Concentration of H+ ions. 2013 cOH_res('pint.Quantity'): Concentration of OH- ions. 2014 cNa_res('pint.Quantity'): Concentration of Na+ ions. 2015 cCl_res('pint.Quantity'): Concentration of Cl- ions. 2016 """ 2017 2018 self_consistent_run = 0 2019 cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+ 2020 2021 def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res): 2022 #Calculate and initial guess for the concentrations of various species based on ideal gas estimate 2023 cOH_res = self.Kw / cH_res 2024 cNa_res = None 2025 cCl_res = None 2026 if cOH_res>=cH_res: 2027 #adjust the concentration of sodium if there is excess OH- in the reservoir: 2028 cNa_res = c_salt_res + (cOH_res-cH_res) 2029 cCl_res = c_salt_res 2030 else: 2031 # adjust the concentration of chloride if there is excess H+ in the reservoir 2032 cCl_res = c_salt_res + (cH_res-cOH_res) 2033 cNa_res = c_salt_res 2034 2035 def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res): 2036 nonlocal max_number_sc_runs, self_consistent_run 2037 if self_consistent_run<max_number_sc_runs: 2038 self_consistent_run+=1 2039 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2040 cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res)) 2041 if cOH_res>=cH_res: 2042 #adjust the concentration of sodium if there is excess OH- in the reservoir: 2043 cNa_res = c_salt_res + (cOH_res-cH_res) 2044 cCl_res = c_salt_res 2045 else: 2046 # adjust the concentration of chloride if there is excess H+ in the reservoir 2047 cCl_res = c_salt_res + (cH_res-cOH_res) 2048 cNa_res = c_salt_res 2049 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 2050 else: 2051 return cH_res, cOH_res, cNa_res, cCl_res 2052 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 2053 2054 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 2055 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2056 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 2057 2058 while abs(determined_pH-pH_res)>1e-6: 2059 if determined_pH > pH_res: 2060 cH_res *= 1.005 2061 else: 2062 cH_res /= 1.003 2063 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 2064 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2065 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 2066 self_consistent_run=0 2067 2068 return cH_res, cOH_res, cNa_res, cCl_res 2069 2070 def enable_motion_of_rigid_object(self, name, espresso_system): 2071 ''' 2072 Enables the motion of the rigid object `name` in the `espresso_system`. 2073 2074 Args: 2075 name(`str`): Label of the object. 2076 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 2077 2078 Note: 2079 - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]. 2080 ''' 2081 logging.info('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]') 2082 self._check_supported_molecule(molecule_name=name, 2083 valid_pmb_types= ['protein']) 2084 molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list() 2085 for molecule_id in molecule_ids_list: 2086 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 2087 center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system) 2088 rigid_object_center = espresso_system.part.add(pos=center_of_mass, 2089 rotation=[True,True,True], 2090 type=self.propose_unused_type()) 2091 2092 rigid_object_center.mass = len(particle_ids_list) 2093 momI = 0 2094 for pid in particle_ids_list: 2095 momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2) 2096 2097 rigid_object_center.rinertia = np.ones(3) * momI 2098 2099 for particle_id in particle_ids_list: 2100 pid = espresso_system.part.by_id(particle_id) 2101 pid.vs_auto_relate_to(rigid_object_center.id) 2102 return 2103 2104 def filter_df(self, pmb_type): 2105 """ 2106 Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns. 2107 2108 Args: 2109 pmb_type(`str`): pmb_object_type to filter in `pmb.df`. 2110 2111 Returns: 2112 pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`. 2113 """ 2114 pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type] 2115 pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1) 2116 return pmb_type_df 2117 2118 def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False): 2119 """ 2120 Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2121 2122 Args: 2123 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 2124 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 2125 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'. 2126 2127 Returns: 2128 bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` if a matching bond exists 2129 2130 Note: 2131 - If `use_default_bond`=`True`, it returns "default" if no key is found. 2132 """ 2133 bond_keys = [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}'] 2134 bond_defined=False 2135 for bond_key in bond_keys: 2136 if bond_key in self.df["name"].values: 2137 bond_defined=True 2138 correct_key=bond_key 2139 break 2140 if bond_defined: 2141 return correct_key 2142 elif use_default_bond: 2143 return 'default' 2144 else: 2145 return 2146 2147 def find_value_from_es_type(self, es_type, column_name): 2148 """ 2149 Finds a value in `pmb.df` for a `column_name` and `es_type` pair. 2150 2151 Args: 2152 es_type(`int`): value of the espresso type 2153 column_name(`str`): name of the column in `pymbe.df` 2154 2155 Returns: 2156 Value in `pymbe.df` matching `column_name` and `es_type` 2157 """ 2158 idx = pd.IndexSlice 2159 for state in ['state_one', 'state_two']: 2160 index = self.df.loc[self.df[(state, 'es_type')] == es_type].index 2161 if len(index) > 0: 2162 if column_name == 'label': 2163 label = self.df.loc[idx[index[0]], idx[(state,column_name)]] 2164 return label 2165 else: 2166 column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]] 2167 return column_name_value 2168 return None 2169 2170 def format_node(self, node_list): 2171 return "[" + " ".join(map(str, node_list)) + "]" 2172 2173 2174 def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples): 2175 """ 2176 Generates coordinates outside a sphere centered at `center`. 2177 2178 Args: 2179 center(`lst`): Coordinates of the center of the sphere. 2180 radius(`float`): Radius of the sphere. 2181 max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates. 2182 n_samples(`int`): Number of sample points. 2183 2184 Returns: 2185 coord_list(`lst`): Coordinates of the sample points. 2186 """ 2187 if not radius > 0: 2188 raise ValueError (f'The value of {radius} must be a positive value') 2189 if not radius < max_dist: 2190 raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))') 2191 coord_list = [] 2192 counter = 0 2193 while counter<n_samples: 2194 coord = self.generate_random_points_in_a_sphere(center=center, 2195 radius=max_dist, 2196 n_samples=1)[0] 2197 if np.linalg.norm(coord-np.asarray(center))>=radius: 2198 coord_list.append (coord) 2199 counter += 1 2200 return coord_list 2201 2202 def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False): 2203 """ 2204 Uniformly samples points from a hypersphere. If on_surface is set to True, the points are 2205 uniformly sampled from the surface of the hypersphere. 2206 2207 Args: 2208 center(`lst`): Array with the coordinates of the center of the spheres. 2209 radius(`float`): Radius of the sphere. 2210 n_samples(`int`): Number of sample points to generate inside the sphere. 2211 on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere. 2212 2213 Returns: 2214 samples(`list`): Coordinates of the sample points inside the hypersphere. 2215 """ 2216 # initial values 2217 center=np.array(center) 2218 d = center.shape[0] 2219 # sample n_samples points in d dimensions from a standard normal distribution 2220 samples = self.rng.normal(size=(n_samples, d)) 2221 # make the samples lie on the surface of the unit hypersphere 2222 normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis] 2223 samples /= normalize_radii 2224 if not on_surface: 2225 # make the samples lie inside the hypersphere with the correct density 2226 uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis] 2227 new_radii = np.power(uniform_points, 1/d) 2228 samples *= new_radii 2229 # scale the points to have the correct radius and center 2230 samples = samples * radius + center 2231 return samples 2232 2233 def generate_trial_perpendicular_vector(self,vector,magnitude): 2234 """ 2235 Generates an orthogonal vector to the input `vector`. 2236 2237 Args: 2238 vector(`lst`): arbitrary vector. 2239 magnitude(`float`): magnitude of the orthogonal vector. 2240 2241 Returns: 2242 (`lst`): Orthogonal vector with the same magnitude as the input vector. 2243 """ 2244 np_vec = np.array(vector) 2245 if np.all(np_vec == 0): 2246 raise ValueError('Zero vector') 2247 np_vec /= np.linalg.norm(np_vec) 2248 # Generate a random vector 2249 random_vector = self.generate_random_points_in_a_sphere(radius=1, 2250 center=[0,0,0], 2251 n_samples=1, 2252 on_surface=True)[0] 2253 # Project the random vector onto the input vector and subtract the projection 2254 projection = np.dot(random_vector, np_vec) * np_vec 2255 perpendicular_vector = random_vector - projection 2256 # Normalize the perpendicular vector to have the same magnitude as the input vector 2257 perpendicular_vector /= np.linalg.norm(perpendicular_vector) 2258 return perpendicular_vector*magnitude 2259 2260 def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2261 """ 2262 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length. 2263 If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead. 2264 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. 2265 2266 Args: 2267 particle_name1(str): label of the type of the first particle type of the bonded particles. 2268 particle_name2(str): label of the type of the second particle type of the bonded particles. 2269 hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2270 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. 2271 2272 Returns: 2273 l0(`pint.Quantity`): bond length 2274 2275 Note: 2276 - 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`. 2277 - If `hard_check`=`True` stops the code when no bond is found. 2278 """ 2279 bond_key = self.find_bond_key(particle_name1=particle_name1, 2280 particle_name2=particle_name2, 2281 use_default_bond=use_default_bond) 2282 if bond_key: 2283 return self.df[self.df['name']==bond_key].l0.values[0] 2284 else: 2285 msg = f"Bond not defined between particles {particle_name1} and {particle_name2}" 2286 if hard_check: 2287 raise ValueError(msg) 2288 else: 2289 logging.warning(msg) 2290 return 2291 2292 def get_charge_number_map(self): 2293 ''' 2294 Gets the charge number of each `espresso_type` in `pymbe.df`. 2295 2296 Returns: 2297 charge_number_map(`dict`): {espresso_type: z}. 2298 ''' 2299 if self.df.state_one['es_type'].isnull().values.any(): 2300 df_state_one = self.df.state_one.dropna() 2301 df_state_two = self.df.state_two.dropna() 2302 else: 2303 df_state_one = self.df.state_one 2304 if self.df.state_two['es_type'].isnull().values.any(): 2305 df_state_two = self.df.state_two.dropna() 2306 else: 2307 df_state_two = self.df.state_two 2308 state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values) 2309 state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values) 2310 charge_number_map = pd.concat([state_one,state_two],axis=0).to_dict() 2311 return charge_number_map 2312 2313 def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'): 2314 """ 2315 Returns the Lennard-Jones parameters for the interaction between the particle types given by 2316 `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule. 2317 2318 Args: 2319 particle_name1 (str): label of the type of the first particle type 2320 particle_name2 (str): label of the type of the second particle type 2321 combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'. 2322 2323 Returns: 2324 {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value} 2325 2326 Note: 2327 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 2328 - 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. 2329 """ 2330 supported_combining_rules=["Lorentz-Berthelot"] 2331 lj_parameters_keys=["sigma","epsilon","offset","cutoff"] 2332 if combining_rule not in supported_combining_rules: 2333 raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}") 2334 lj_parameters={} 2335 for key in lj_parameters_keys: 2336 lj_parameters[key]=[] 2337 # Search the LJ parameters of the type pair 2338 for name in [particle_name1,particle_name2]: 2339 for key in lj_parameters_keys: 2340 lj_parameters[key].append(self.df[self.df.name == name][key].values[0]) 2341 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 2342 if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]): 2343 return {} 2344 # Apply combining rule 2345 if combining_rule == 'Lorentz-Berthelot': 2346 lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2 2347 lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2 2348 lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2 2349 lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1]) 2350 return lj_parameters 2351 2352 def get_metal_ions_charge_number_map(self): 2353 """ 2354 Gets a map with the charge numbers of all the metal ions supported. 2355 2356 Returns: 2357 metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number} 2358 2359 """ 2360 metal_charge_number_map = {"Ca": 2} 2361 return metal_charge_number_map 2362 2363 def get_particle_id_map(self, object_name): 2364 ''' 2365 Gets all the ids associated with the object with name `object_name` in `pmb.df` 2366 2367 Args: 2368 object_name(`str`): name of the object 2369 2370 Returns: 2371 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]}, } 2372 ''' 2373 object_type=self._check_supported_molecule(molecule_name=object_name, 2374 valid_pmb_types= ['particle','residue','molecule',"peptide","protein"]) 2375 id_list = [] 2376 mol_map = {} 2377 res_map = {} 2378 def do_res_map(res_ids): 2379 for res_id in res_ids: 2380 res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist() 2381 res_map[res_id]=res_list 2382 return res_map 2383 if object_type in ['molecule', 'protein', 'peptide']: 2384 mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist() 2385 for mol_id in mol_ids: 2386 res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist()) 2387 res_map=do_res_map(res_ids=res_ids) 2388 mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist() 2389 id_list+=mol_list 2390 mol_map[mol_id]=mol_list 2391 elif object_type == 'residue': 2392 res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist() 2393 res_map=do_res_map(res_ids=res_ids) 2394 id_list=[] 2395 for res_id_list in res_map.values(): 2396 id_list+=res_id_list 2397 elif object_type == 'particle': 2398 id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist() 2399 return {"all": id_list, "molecule_map": mol_map, "residue_map": res_map} 2400 2401 def get_pka_set(self): 2402 ''' 2403 Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df` 2404 2405 Returns: 2406 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 2407 ''' 2408 titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna() 2409 pka_set = {} 2410 for index in titratables_AA_df.name.keys(): 2411 name = titratables_AA_df.name[index] 2412 pka_value = titratables_AA_df.pka[index] 2413 acidity = titratables_AA_df.acidity[index] 2414 pka_set[name] = {'pka_value':pka_value,'acidity':acidity} 2415 return pka_set 2416 2417 def get_radius_map(self, dimensionless=True): 2418 ''' 2419 Gets the effective radius of each `espresso type` in `pmb.df`. 2420 2421 Args: 2422 dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False. 2423 2424 Returns: 2425 radius_map(`dict`): {espresso_type: radius}. 2426 2427 Note: 2428 The radius corresponds to (sigma+offset)/2 2429 ''' 2430 df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates() 2431 df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates() 2432 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) 2433 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) 2434 radius_map = pd.concat([state_one,state_two],axis=0).to_dict() 2435 if dimensionless: 2436 for key in radius_map: 2437 radius_map[key] = radius_map[key].magnitude 2438 return radius_map 2439 2440 def get_reduced_units(self): 2441 """ 2442 Returns the current set of reduced units defined in pyMBE. 2443 2444 Returns: 2445 reduced_units_text(`str`): text with information about the current set of reduced units. 2446 2447 """ 2448 unit_length=self.units.Quantity(1,'reduced_length') 2449 unit_energy=self.units.Quantity(1,'reduced_energy') 2450 unit_charge=self.units.Quantity(1,'reduced_charge') 2451 reduced_units_text = "\n".join(["Current set of reduced units:", 2452 f"{unit_length.to('nm'):.5g} = {unit_length}", 2453 f"{unit_energy.to('J'):.5g} = {unit_energy}", 2454 f"{unit_charge.to('C'):.5g} = {unit_charge}", 2455 f"Temperature: {(self.kT/self.kB).to('K'):.5g}" 2456 ]) 2457 return reduced_units_text 2458 2459 def get_type_map(self): 2460 """ 2461 Gets all different espresso types assigned to particles in `pmb.df`. 2462 2463 Returns: 2464 type_map(`dict`): {"name": espresso_type}. 2465 """ 2466 df_state_one = self.df.state_one.dropna(how='all') 2467 df_state_two = self.df.state_two.dropna(how='all') 2468 state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label) 2469 state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label) 2470 type_map = pd.concat([state_one,state_two],axis=0).to_dict() 2471 return type_map 2472 2473 def initialize_lattice_builder(self, diamond_lattice): 2474 """ 2475 Initialize the lattice builder with the DiamondLattice object. 2476 2477 Args: 2478 diamond_lattice(`DiamondLattice`): DiamondLattice object from the `lib/lattice` module to be used in the LatticeBuilder. 2479 """ 2480 from .lib.lattice import LatticeBuilder, DiamondLattice 2481 if not isinstance(diamond_lattice, DiamondLattice): 2482 raise TypeError("Currently only DiamondLattice objects are supported.") 2483 self.lattice_builder = LatticeBuilder(lattice=diamond_lattice) 2484 logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}") 2485 return self.lattice_builder 2486 2487 2488 def load_interaction_parameters(self, filename, overwrite=False): 2489 """ 2490 Loads the interaction parameters stored in `filename` into `pmb.df` 2491 2492 Args: 2493 filename(`str`): name of the file to be read 2494 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 2495 """ 2496 without_units = ['z','es_type'] 2497 with_units = ['sigma','epsilon','offset','cutoff'] 2498 2499 with open(filename, 'r') as f: 2500 interaction_data = json.load(f) 2501 interaction_parameter_set = interaction_data["data"] 2502 2503 for key in interaction_parameter_set: 2504 param_dict=interaction_parameter_set[key] 2505 object_type=param_dict.pop('object_type') 2506 if object_type == 'particle': 2507 not_required_attributes={} 2508 for not_required_key in without_units+with_units: 2509 if not_required_key in param_dict.keys(): 2510 if not_required_key in with_units: 2511 not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key)) 2512 elif not_required_key in without_units: 2513 not_required_attributes[not_required_key]=param_dict.pop(not_required_key) 2514 else: 2515 not_required_attributes[not_required_key]=pd.NA 2516 self.define_particle(name=param_dict.pop('name'), 2517 z=not_required_attributes.pop('z'), 2518 sigma=not_required_attributes.pop('sigma'), 2519 offset=not_required_attributes.pop('offset'), 2520 cutoff=not_required_attributes.pop('cutoff'), 2521 epsilon=not_required_attributes.pop('epsilon'), 2522 overwrite=overwrite) 2523 elif object_type == 'residue': 2524 self.define_residue(**param_dict) 2525 elif object_type == 'molecule': 2526 self.define_molecule(**param_dict) 2527 elif object_type == 'peptide': 2528 self.define_peptide(**param_dict) 2529 elif object_type == 'bond': 2530 particle_pairs = param_dict.pop('particle_pairs') 2531 bond_parameters = param_dict.pop('bond_parameters') 2532 bond_type = param_dict.pop('bond_type') 2533 if bond_type == 'harmonic': 2534 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2535 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2536 bond = {'r_0' : r_0, 2537 'k' : k, 2538 } 2539 2540 elif bond_type == 'FENE': 2541 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2542 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2543 d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max')) 2544 bond = {'r_0' : r_0, 2545 'k' : k, 2546 'd_r_max': d_r_max, 2547 } 2548 else: 2549 raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds") 2550 self.define_bond(bond_type=bond_type, 2551 bond_parameters=bond, 2552 particle_pairs=particle_pairs) 2553 else: 2554 raise ValueError(object_type+' is not a known pmb object type') 2555 2556 return 2557 2558 def load_pka_set(self, filename, overwrite=True): 2559 """ 2560 Loads the pka_set stored in `filename` into `pmb.df`. 2561 2562 Args: 2563 filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}. 2564 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 2565 """ 2566 with open(filename, 'r') as f: 2567 pka_data = json.load(f) 2568 pka_set = pka_data["data"] 2569 2570 self.check_pka_set(pka_set=pka_set) 2571 2572 for key in pka_set: 2573 acidity = pka_set[key]['acidity'] 2574 pka_value = pka_set[key]['pka_value'] 2575 self.set_particle_acidity(name=key, 2576 acidity=acidity, 2577 pka=pka_value, 2578 overwrite=overwrite) 2579 return 2580 2581 2582 def propose_unused_type(self): 2583 """ 2584 Searches in `pmb.df` all the different particle types defined and returns a new one. 2585 2586 Returns: 2587 unused_type(`int`): unused particle type 2588 """ 2589 type_map = self.get_type_map() 2590 if not type_map: 2591 unused_type = 0 2592 else: 2593 valid_values = [v for v in type_map.values() if pd.notna(v)] # Filter out pd.NA values 2594 unused_type = max(valid_values) + 1 if valid_values else 0 # Ensure max() doesn't fail if all values are NA 2595 return unused_type 2596 2597 def protein_sequence_parser(self, sequence): 2598 ''' 2599 Parses `sequence` to the one letter code for amino acids. 2600 2601 Args: 2602 sequence(`str` or `lst`): Sequence of the amino acid. 2603 2604 Returns: 2605 clean_sequence(`lst`): `sequence` using the one letter code. 2606 2607 Note: 2608 - Accepted formats for `sequence` are: 2609 - `lst` with one letter or three letter code of each aminoacid in each element 2610 - `str` with the sequence using the one letter code 2611 - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-" 2612 2613 ''' 2614 # Aminoacid key 2615 keys={"ALA": "A", 2616 "ARG": "R", 2617 "ASN": "N", 2618 "ASP": "D", 2619 "CYS": "C", 2620 "GLU": "E", 2621 "GLN": "Q", 2622 "GLY": "G", 2623 "HIS": "H", 2624 "ILE": "I", 2625 "LEU": "L", 2626 "LYS": "K", 2627 "MET": "M", 2628 "PHE": "F", 2629 "PRO": "P", 2630 "SER": "S", 2631 "THR": "T", 2632 "TRP": "W", 2633 "TYR": "Y", 2634 "VAL": "V", 2635 "PSER": "J", 2636 "PTHR": "U", 2637 "PTyr": "Z", 2638 "NH2": "n", 2639 "COOH": "c"} 2640 clean_sequence=[] 2641 if isinstance(sequence, str): 2642 if sequence.find("-") != -1: 2643 splited_sequence=sequence.split("-") 2644 for residue in splited_sequence: 2645 if len(residue) == 1: 2646 if residue in keys.values(): 2647 residue_ok=residue 2648 else: 2649 if residue.upper() in keys.values(): 2650 residue_ok=residue.upper() 2651 else: 2652 raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence") 2653 clean_sequence.append(residue_ok) 2654 else: 2655 if residue in keys.keys(): 2656 clean_sequence.append(keys[residue]) 2657 else: 2658 if residue.upper() in keys.keys(): 2659 clean_sequence.append(keys[residue.upper()]) 2660 else: 2661 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2662 else: 2663 for residue in sequence: 2664 if residue in keys.values(): 2665 residue_ok=residue 2666 else: 2667 if residue.upper() in keys.values(): 2668 residue_ok=residue.upper() 2669 else: 2670 raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence") 2671 clean_sequence.append(residue_ok) 2672 if isinstance(sequence, list): 2673 for residue in sequence: 2674 if residue in keys.values(): 2675 residue_ok=residue 2676 else: 2677 if residue.upper() in keys.values(): 2678 residue_ok=residue.upper() 2679 elif (residue.upper() in keys.keys()): 2680 residue_ok= keys[residue.upper()] 2681 else: 2682 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2683 clean_sequence.append(residue_ok) 2684 return clean_sequence 2685 2686 def read_pmb_df (self,filename): 2687 """ 2688 Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE. 2689 2690 Args: 2691 filename(`str`): path to file. 2692 2693 Note: 2694 This function only accepts files with CSV format. 2695 """ 2696 2697 if filename.rsplit(".", 1)[1] != "csv": 2698 raise ValueError("Only files with CSV format are supported") 2699 df = pd.read_csv (filename,header=[0, 1], index_col=0) 2700 columns_names = self.setup_df() 2701 2702 multi_index = pd.MultiIndex.from_tuples(columns_names) 2703 df.columns = multi_index 2704 2705 self.df = self.convert_columns_to_original_format(df) 2706 self.df.fillna(pd.NA, inplace=True) 2707 return self.df 2708 2709 def read_protein_vtf_in_df (self,filename,unit_length=None): 2710 """ 2711 Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`. 2712 2713 Args: 2714 filename(`str`): Path to the vtf file with the coarse-grained model. 2715 unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None. 2716 2717 Returns: 2718 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 2719 2720 Note: 2721 - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom. 2722 """ 2723 2724 logging.info(f'Loading protein coarse grain model file: {filename}') 2725 2726 coord_list = [] 2727 particles_dict = {} 2728 2729 if unit_length is None: 2730 unit_length = 1 * self.units.angstrom 2731 2732 with open (filename,'r') as protein_model: 2733 for line in protein_model : 2734 line_split = line.split() 2735 if line_split: 2736 line_header = line_split[0] 2737 if line_header == 'atom': 2738 atom_id = line_split[1] 2739 atom_name = line_split[3] 2740 atom_resname = line_split[5] 2741 chain_id = line_split[9] 2742 radius = float(line_split [11])*unit_length 2743 particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius] 2744 elif line_header.isnumeric(): 2745 atom_coord = line_split[1:] 2746 atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord] 2747 coord_list.append (atom_coord) 2748 2749 numbered_label = [] 2750 i = 0 2751 2752 for atom_id in particles_dict.keys(): 2753 2754 if atom_id == 1: 2755 atom_name = particles_dict[atom_id][0] 2756 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2757 numbered_label.append(numbered_name) 2758 2759 elif atom_id != 1: 2760 2761 if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]: 2762 i += 1 2763 count = 1 2764 atom_name = particles_dict[atom_id][0] 2765 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2766 numbered_label.append(numbered_name) 2767 2768 elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]: 2769 if count == 2 or particles_dict[atom_id][1] == 'GLY': 2770 i +=1 2771 count = 0 2772 atom_name = particles_dict[atom_id][0] 2773 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2774 numbered_label.append(numbered_name) 2775 count +=1 2776 2777 topology_dict = {} 2778 2779 for i in range (0, len(numbered_label)): 2780 topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] , 2781 'chain_id':numbered_label[i][1], 2782 'radius':numbered_label[i][2] } 2783 2784 return topology_dict 2785 2786 def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2787 """ 2788 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2789 If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead. 2790 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. 2791 2792 Args: 2793 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 2794 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 2795 hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2796 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. 2797 2798 Returns: 2799 bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library. 2800 2801 Note: 2802 - 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`. 2803 - If `hard_check`=`True` stops the code when no bond is found. 2804 """ 2805 2806 bond_key = self.find_bond_key(particle_name1=particle_name1, 2807 particle_name2=particle_name2, 2808 use_default_bond=use_default_bond) 2809 if use_default_bond: 2810 if not self._check_if_name_is_defined_in_df(name="default"): 2811 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") 2812 if bond_key: 2813 return self.df[self.df['name']==bond_key].bond_object.values[0] 2814 else: 2815 msg= f"Bond not defined between particles {particle_name1} and {particle_name2}" 2816 if hard_check: 2817 raise ValueError(msg) 2818 else: 2819 logging.warning(msg) 2820 return None 2821 def search_particles_in_residue(self, residue_name): 2822 ''' 2823 Searches for all particles in a given residue of name `residue_name`. 2824 2825 Args: 2826 residue_name (`str`): name of the residue to be searched 2827 2828 Returns: 2829 list_of_particles_in_residue (`lst`): list of the names of all particles in the residue 2830 2831 Note: 2832 - 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. 2833 - The function will return an empty list if the residue is not defined in `pmb.df`. 2834 - The function will return an empty list if the particles are not defined in the pyMBE DataFrame. 2835 ''' 2836 if not self._check_if_name_is_defined_in_df(name=residue_name): 2837 logging.warning(f"Residue {residue_name} not defined in pmb.df") 2838 return [] 2839 self._check_if_name_has_right_type(name=residue_name, expected_pmb_type="residue") 2840 index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 2841 central_bead = self.df.at [index_residue, ('central_bead', '')] 2842 list_of_side_chains = self.df.at[index_residue, ('side_chains', '')] 2843 list_of_particles_in_residue = [] 2844 if central_bead is not pd.NA: 2845 if self._check_if_name_is_defined_in_df(name=central_bead): 2846 if self._check_if_name_has_right_type(name=central_bead, expected_pmb_type="particle", hard_check=False): 2847 list_of_particles_in_residue.append(central_bead) 2848 if list_of_side_chains is not pd.NA: 2849 for side_chain in list_of_side_chains: 2850 if self._check_if_name_is_defined_in_df(name=side_chain): 2851 object_type = self.df[self.df['name']==side_chain].pmb_type.values[0] 2852 else: 2853 continue 2854 if object_type == "residue": 2855 list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain) 2856 list_of_particles_in_residue += list_of_particles_in_side_chain_residue 2857 elif object_type == "particle": 2858 if side_chain is not pd.NA: 2859 list_of_particles_in_residue.append(side_chain) 2860 return list_of_particles_in_residue 2861 2862 def set_particle_acidity(self, name, acidity=pd.NA, default_charge_number=0, pka=pd.NA, overwrite=True): 2863 """ 2864 Sets the particle acidity including the charges in each of its possible states. 2865 2866 Args: 2867 name(`str`): Unique label that identifies the `particle`. 2868 acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None. 2869 default_charge_number (`int`): Charge number of the particle. Defaults to 0. 2870 pka(`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to pandas.NA. 2871 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 2872 2873 Note: 2874 - For particles with `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 2875 deprotonated states, respectively. 2876 - For particles without an acidity `acidity = pandas.NA`, only `state_one` is defined. 2877 - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 2878 """ 2879 acidity_valid_keys = ['inert','acidic', 'basic'] 2880 if not pd.isna(acidity): 2881 if acidity not in acidity_valid_keys: 2882 raise ValueError(f"Acidity {acidity} provided for particle name {name} is not supproted. Valid keys are: {acidity_valid_keys}") 2883 if acidity in ['acidic', 'basic'] and pd.isna(pka): 2884 raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.") 2885 if acidity == "inert": 2886 acidity = pd.NA 2887 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") 2888 2889 self.define_particle_entry_in_df(name=name) 2890 2891 for index in self.df[self.df['name']==name].index: 2892 if pka is not pd.NA: 2893 self.add_value_to_df(key=('pka',''), 2894 index=index, 2895 new_value=pka, 2896 overwrite=overwrite) 2897 2898 self.add_value_to_df(key=('acidity',''), 2899 index=index, 2900 new_value=acidity, 2901 overwrite=overwrite) 2902 if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')): 2903 self.add_value_to_df(key=('state_one','es_type'), 2904 index=index, 2905 new_value=self.propose_unused_type(), 2906 overwrite=overwrite) 2907 if pd.isna(self.df.loc [self.df['name'] == name].acidity.iloc[0]): 2908 self.add_value_to_df(key=('state_one','z'), 2909 index=index, 2910 new_value=default_charge_number, 2911 overwrite=overwrite) 2912 self.add_value_to_df(key=('state_one','label'), 2913 index=index, 2914 new_value=name, 2915 overwrite=overwrite) 2916 else: 2917 protonated_label = f'{name}H' 2918 self.add_value_to_df(key=('state_one','label'), 2919 index=index, 2920 new_value=protonated_label, 2921 overwrite=overwrite) 2922 self.add_value_to_df(key=('state_two','label'), 2923 index=index, 2924 new_value=name, 2925 overwrite=overwrite) 2926 if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')): 2927 self.add_value_to_df(key=('state_two','es_type'), 2928 index=index, 2929 new_value=self.propose_unused_type(), 2930 overwrite=overwrite) 2931 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'acidic': 2932 self.add_value_to_df(key=('state_one','z'), 2933 index=index,new_value=0, 2934 overwrite=overwrite) 2935 self.add_value_to_df(key=('state_two','z'), 2936 index=index, 2937 new_value=-1, 2938 overwrite=overwrite) 2939 elif self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'basic': 2940 self.add_value_to_df(key=('state_one','z'), 2941 index=index,new_value=+1, 2942 overwrite=overwrite) 2943 self.add_value_to_df(key=('state_two','z'), 2944 index=index, 2945 new_value=0, 2946 overwrite=overwrite) 2947 self.df.fillna(pd.NA, inplace=True) 2948 return 2949 2950 def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None): 2951 """ 2952 Sets the set of reduced units used by pyMBE.units and it prints it. 2953 2954 Args: 2955 unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 2956 unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 2957 temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 2958 Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 2959 2960 Note: 2961 - If no `temperature` is given, a value of 298.15 K is assumed by default. 2962 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 2963 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 2964 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 2965 """ 2966 if unit_length is None: 2967 unit_length= 0.355*self.units.nm 2968 if temperature is None: 2969 temperature = 298.15 * self.units.K 2970 if unit_charge is None: 2971 unit_charge = scipy.constants.e * self.units.C 2972 if Kw is None: 2973 Kw = 1e-14 2974 # Sanity check 2975 variables=[unit_length,temperature,unit_charge] 2976 dimensionalities=["[length]","[temperature]","[charge]"] 2977 for variable,dimensionality in zip(variables,dimensionalities): 2978 self.check_dimensionality(variable,dimensionality) 2979 self.Kw=Kw*self.units.mol**2 / (self.units.l**2) 2980 self.kT=temperature*self.kB 2981 self.units._build_cache() 2982 self.units.define(f'reduced_energy = {self.kT} ') 2983 self.units.define(f'reduced_length = {unit_length}') 2984 self.units.define(f'reduced_charge = {unit_charge}') 2985 logging.info(self.get_reduced_units()) 2986 return 2987 2988 def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2989 """ 2990 Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 2991 2992 Args: 2993 counter_ion(`str`): `name` of the counter_ion `particle`. 2994 constant_pH(`float`): pH-value. 2995 exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted. 2996 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2997 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2998 2999 Returns: 3000 RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library. 3001 sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3002 """ 3003 from espressomd import reaction_methods 3004 if exclusion_range is None: 3005 exclusion_range = max(self.get_radius_map().values())*2.0 3006 if pka_set is None: 3007 pka_set=self.get_pka_set() 3008 self.check_pka_set(pka_set=pka_set) 3009 if use_exclusion_radius_per_type: 3010 exclusion_radius_per_type = self.get_radius_map() 3011 else: 3012 exclusion_radius_per_type = {} 3013 3014 RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3015 exclusion_range=exclusion_range, 3016 seed=self.seed, 3017 constant_pH=constant_pH, 3018 exclusion_radius_per_type = exclusion_radius_per_type 3019 ) 3020 sucessfull_reactions_labels=[] 3021 charge_number_map = self.get_charge_number_map() 3022 for name in pka_set.keys(): 3023 if not self._check_if_name_is_defined_in_df(name): 3024 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.') 3025 continue 3026 gamma=10**-pka_set[name]['pka_value'] 3027 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3028 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3029 counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0] 3030 RE.add_reaction(gamma=gamma, 3031 reactant_types=[state_one_type], 3032 product_types=[state_two_type, counter_ion_type], 3033 default_charges={state_one_type: charge_number_map[state_one_type], 3034 state_two_type: charge_number_map[state_two_type], 3035 counter_ion_type: charge_number_map[counter_ion_type]}) 3036 sucessfull_reactions_labels.append(name) 3037 return RE, sucessfull_reactions_labels 3038 3039 def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False): 3040 """ 3041 Sets up grand-canonical coupling to a reservoir of salt. 3042 For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead. 3043 3044 Args: 3045 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3046 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 3047 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 3048 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3049 exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted. 3050 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 3051 3052 Returns: 3053 RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library. 3054 """ 3055 from espressomd import reaction_methods 3056 if exclusion_range is None: 3057 exclusion_range = max(self.get_radius_map().values())*2.0 3058 if use_exclusion_radius_per_type: 3059 exclusion_radius_per_type = self.get_radius_map() 3060 else: 3061 exclusion_radius_per_type = {} 3062 3063 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3064 exclusion_range=exclusion_range, 3065 seed=self.seed, 3066 exclusion_radius_per_type = exclusion_radius_per_type 3067 ) 3068 3069 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3070 determined_activity_coefficient = activity_coefficient(c_salt_res) 3071 K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient 3072 3073 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 3074 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 3075 3076 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0] 3077 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0] 3078 3079 if salt_cation_charge <= 0: 3080 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 3081 if salt_anion_charge >= 0: 3082 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 3083 3084 # Grand-canonical coupling to the reservoir 3085 RE.add_reaction( 3086 gamma = K_salt.magnitude, 3087 reactant_types = [], 3088 reactant_coefficients = [], 3089 product_types = [ salt_cation_es_type, salt_anion_es_type ], 3090 product_coefficients = [ 1, 1 ], 3091 default_charges = { 3092 salt_cation_es_type: salt_cation_charge, 3093 salt_anion_es_type: salt_anion_charge, 3094 } 3095 ) 3096 3097 return RE 3098 3099 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): 3100 """ 3101 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 3102 reservoir of small ions. 3103 This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1]. 3104 3105 [1] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020. 3106 3107 Args: 3108 pH_res ('float): pH-value in the reservoir. 3109 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3110 proton_name ('str'): Name of the proton (H+) particle. 3111 hydroxide_name ('str'): Name of the hydroxide (OH-) particle. 3112 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 3113 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 3114 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3115 exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted. 3116 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 3117 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 3118 3119 Returns: 3120 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 3121 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3122 ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients). 3123 """ 3124 from espressomd import reaction_methods 3125 if exclusion_range is None: 3126 exclusion_range = max(self.get_radius_map().values())*2.0 3127 if pka_set is None: 3128 pka_set=self.get_pka_set() 3129 self.check_pka_set(pka_set=pka_set) 3130 if use_exclusion_radius_per_type: 3131 exclusion_radius_per_type = self.get_radius_map() 3132 else: 3133 exclusion_radius_per_type = {} 3134 3135 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3136 exclusion_range=exclusion_range, 3137 seed=self.seed, 3138 exclusion_radius_per_type = exclusion_radius_per_type 3139 ) 3140 3141 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3142 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 3143 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 3144 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 3145 K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3146 K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3147 K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3148 3149 proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0] 3150 hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0] 3151 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 3152 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 3153 3154 proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0] 3155 hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0] 3156 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0] 3157 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0] 3158 3159 if proton_charge <= 0: 3160 raise ValueError('ERROR proton charge must be positive, charge ', proton_charge) 3161 if salt_cation_charge <= 0: 3162 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 3163 if hydroxide_charge >= 0: 3164 raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge) 3165 if salt_anion_charge >= 0: 3166 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 3167 3168 # Grand-canonical coupling to the reservoir 3169 # 0 = H+ + OH- 3170 RE.add_reaction( 3171 gamma = K_W.magnitude, 3172 reactant_types = [], 3173 reactant_coefficients = [], 3174 product_types = [ proton_es_type, hydroxide_es_type ], 3175 product_coefficients = [ 1, 1 ], 3176 default_charges = { 3177 proton_es_type: proton_charge, 3178 hydroxide_es_type: hydroxide_charge, 3179 } 3180 ) 3181 3182 # 0 = Na+ + Cl- 3183 RE.add_reaction( 3184 gamma = K_NACL.magnitude, 3185 reactant_types = [], 3186 reactant_coefficients = [], 3187 product_types = [ salt_cation_es_type, salt_anion_es_type ], 3188 product_coefficients = [ 1, 1 ], 3189 default_charges = { 3190 salt_cation_es_type: salt_cation_charge, 3191 salt_anion_es_type: salt_anion_charge, 3192 } 3193 ) 3194 3195 # 0 = Na+ + OH- 3196 RE.add_reaction( 3197 gamma = (K_NACL * K_W / K_HCL).magnitude, 3198 reactant_types = [], 3199 reactant_coefficients = [], 3200 product_types = [ salt_cation_es_type, hydroxide_es_type ], 3201 product_coefficients = [ 1, 1 ], 3202 default_charges = { 3203 salt_cation_es_type: salt_cation_charge, 3204 hydroxide_es_type: hydroxide_charge, 3205 } 3206 ) 3207 3208 # 0 = H+ + Cl- 3209 RE.add_reaction( 3210 gamma = K_HCL.magnitude, 3211 reactant_types = [], 3212 reactant_coefficients = [], 3213 product_types = [ proton_es_type, salt_anion_es_type ], 3214 product_coefficients = [ 1, 1 ], 3215 default_charges = { 3216 proton_es_type: proton_charge, 3217 salt_anion_es_type: salt_anion_charge, 3218 } 3219 ) 3220 3221 # Annealing moves to ensure sufficient sampling 3222 # Cation annealing H+ = Na+ 3223 RE.add_reaction( 3224 gamma = (K_NACL / K_HCL).magnitude, 3225 reactant_types = [proton_es_type], 3226 reactant_coefficients = [ 1 ], 3227 product_types = [ salt_cation_es_type ], 3228 product_coefficients = [ 1 ], 3229 default_charges = { 3230 proton_es_type: proton_charge, 3231 salt_cation_es_type: salt_cation_charge, 3232 } 3233 ) 3234 3235 # Anion annealing OH- = Cl- 3236 RE.add_reaction( 3237 gamma = (K_HCL / K_W).magnitude, 3238 reactant_types = [hydroxide_es_type], 3239 reactant_coefficients = [ 1 ], 3240 product_types = [ salt_anion_es_type ], 3241 product_coefficients = [ 1 ], 3242 default_charges = { 3243 hydroxide_es_type: hydroxide_charge, 3244 salt_anion_es_type: salt_anion_charge, 3245 } 3246 ) 3247 3248 sucessful_reactions_labels=[] 3249 charge_number_map = self.get_charge_number_map() 3250 for name in pka_set.keys(): 3251 if not self._check_if_name_is_defined_in_df(name): 3252 logging.warning(f'The acid-base reaction of {name} has not been set up because its particle type is not defined in the dataframe.') 3253 continue 3254 3255 Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 3256 3257 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3258 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3259 3260 # Reaction in terms of proton: HA = A + H+ 3261 RE.add_reaction(gamma=Ka.magnitude, 3262 reactant_types=[state_one_type], 3263 reactant_coefficients=[1], 3264 product_types=[state_two_type, proton_es_type], 3265 product_coefficients=[1, 1], 3266 default_charges={state_one_type: charge_number_map[state_one_type], 3267 state_two_type: charge_number_map[state_two_type], 3268 proton_es_type: proton_charge}) 3269 3270 # Reaction in terms of salt cation: HA = A + Na+ 3271 RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude, 3272 reactant_types=[state_one_type], 3273 reactant_coefficients=[1], 3274 product_types=[state_two_type, salt_cation_es_type], 3275 product_coefficients=[1, 1], 3276 default_charges={state_one_type: charge_number_map[state_one_type], 3277 state_two_type: charge_number_map[state_two_type], 3278 salt_cation_es_type: salt_cation_charge}) 3279 3280 # Reaction in terms of hydroxide: OH- + HA = A 3281 RE.add_reaction(gamma=(Ka / K_W).magnitude, 3282 reactant_types=[state_one_type, hydroxide_es_type], 3283 reactant_coefficients=[1, 1], 3284 product_types=[state_two_type], 3285 product_coefficients=[1], 3286 default_charges={state_one_type: charge_number_map[state_one_type], 3287 state_two_type: charge_number_map[state_two_type], 3288 hydroxide_es_type: hydroxide_charge}) 3289 3290 # Reaction in terms of salt anion: Cl- + HA = A 3291 RE.add_reaction(gamma=(Ka / K_HCL).magnitude, 3292 reactant_types=[state_one_type, salt_anion_es_type], 3293 reactant_coefficients=[1, 1], 3294 product_types=[state_two_type], 3295 product_coefficients=[1], 3296 default_charges={state_one_type: charge_number_map[state_one_type], 3297 state_two_type: charge_number_map[state_two_type], 3298 salt_anion_es_type: salt_anion_charge}) 3299 3300 sucessful_reactions_labels.append(name) 3301 return RE, sucessful_reactions_labels, ionic_strength_res 3302 3303 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): 3304 """ 3305 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 3306 reservoir of small ions. 3307 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-}. 3308 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'. 3309 3310 [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4). 3311 [2] Landsgesell, J., Hebbeker, P., Rud, O., Lunkad, R., Košovan, P., & Holm, C. (2020). Grand-reaction method for simulations of ionization equilibria coupled to ion partitioning. Macromolecules, 53(8), 3007-3020. 3312 3313 Args: 3314 pH_res ('float'): pH-value in the reservoir. 3315 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3316 cation_name ('str'): Name of the cationic particle. 3317 anion_name ('str'): Name of the anionic particle. 3318 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3319 exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted. 3320 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 3321 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`. 3322 3323 Returns: 3324 RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 3325 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3326 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 3327 """ 3328 from espressomd import reaction_methods 3329 if exclusion_range is None: 3330 exclusion_range = max(self.get_radius_map().values())*2.0 3331 if pka_set is None: 3332 pka_set=self.get_pka_set() 3333 self.check_pka_set(pka_set=pka_set) 3334 if use_exclusion_radius_per_type: 3335 exclusion_radius_per_type = self.get_radius_map() 3336 else: 3337 exclusion_radius_per_type = {} 3338 3339 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3340 exclusion_range=exclusion_range, 3341 seed=self.seed, 3342 exclusion_radius_per_type = exclusion_radius_per_type 3343 ) 3344 3345 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3346 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 3347 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 3348 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 3349 a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 3350 a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 3351 a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 3352 K_XX = a_cation * a_anion 3353 3354 cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0] 3355 anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0] 3356 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0] 3357 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0] 3358 if cation_charge <= 0: 3359 raise ValueError('ERROR cation charge must be positive, charge ', cation_charge) 3360 if anion_charge >= 0: 3361 raise ValueError('ERROR anion charge must be negative, charge ', anion_charge) 3362 3363 # Coupling to the reservoir: 0 = X+ + X- 3364 RE.add_reaction( 3365 gamma = K_XX.magnitude, 3366 reactant_types = [], 3367 reactant_coefficients = [], 3368 product_types = [ cation_es_type, anion_es_type ], 3369 product_coefficients = [ 1, 1 ], 3370 default_charges = { 3371 cation_es_type: cation_charge, 3372 anion_es_type: anion_charge, 3373 } 3374 ) 3375 3376 sucessful_reactions_labels=[] 3377 charge_number_map = self.get_charge_number_map() 3378 for name in pka_set.keys(): 3379 if not self._check_if_name_is_defined_in_df(name): 3380 logging.warning(f'The acid-base reaction of {name} has not been set up because its particle type is not defined in the dataframe.') 3381 continue 3382 3383 Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l 3384 gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen 3385 3386 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3387 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3388 3389 # Reaction in terms of small cation: HA = A + X+ 3390 RE.add_reaction(gamma=gamma_K_AX.magnitude, 3391 reactant_types=[state_one_type], 3392 reactant_coefficients=[1], 3393 product_types=[state_two_type, cation_es_type], 3394 product_coefficients=[1, 1], 3395 default_charges={state_one_type: charge_number_map[state_one_type], 3396 state_two_type: charge_number_map[state_two_type], 3397 cation_es_type: cation_charge}) 3398 3399 # Reaction in terms of small anion: X- + HA = A 3400 RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude, 3401 reactant_types=[state_one_type, anion_es_type], 3402 reactant_coefficients=[1, 1], 3403 product_types=[state_two_type], 3404 product_coefficients=[1], 3405 default_charges={state_one_type: charge_number_map[state_one_type], 3406 state_two_type: charge_number_map[state_two_type], 3407 anion_es_type: anion_charge}) 3408 3409 sucessful_reactions_labels.append(name) 3410 return RE, sucessful_reactions_labels, ionic_strength_res 3411 3412 def setup_df (self): 3413 """ 3414 Sets up the pyMBE's dataframe `pymbe.df`. 3415 3416 Returns: 3417 columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe 3418 """ 3419 3420 columns_dtypes = { 3421 'name': { 3422 '': str}, 3423 'pmb_type': { 3424 '': str}, 3425 'particle_id': { 3426 '': pd.Int64Dtype()}, 3427 'particle_id2': { 3428 '': pd.Int64Dtype()}, 3429 'residue_id': { 3430 '': pd.Int64Dtype()}, 3431 'molecule_id': { 3432 '': pd.Int64Dtype()}, 3433 'acidity': { 3434 '': str}, 3435 'pka': { 3436 '': object}, 3437 'central_bead': { 3438 '': object}, 3439 'side_chains': { 3440 '': object}, 3441 'residue_list': { 3442 '': object}, 3443 'model': { 3444 '': str}, 3445 'sigma': { 3446 '': object}, 3447 'cutoff': { 3448 '': object}, 3449 'offset': { 3450 '': object}, 3451 'epsilon': { 3452 '': object}, 3453 'state_one': { 3454 'label': str, 3455 'es_type': pd.Int64Dtype(), 3456 'z': pd.Int64Dtype()}, 3457 'state_two': { 3458 'label': str, 3459 'es_type': pd.Int64Dtype(), 3460 'z': pd.Int64Dtype()}, 3461 'sequence': { 3462 '': object}, 3463 'bond_object': { 3464 '': object}, 3465 'parameters_of_the_potential':{ 3466 '': object}, 3467 'l0': { 3468 '': float}, 3469 'node_map':{ 3470 '':object}, 3471 'chain_map':{ 3472 '':object}} 3473 3474 self.df = pd.DataFrame(columns=pd.MultiIndex.from_tuples([(col_main, col_sub) for col_main, sub_cols in columns_dtypes.items() for col_sub in sub_cols.keys()])) 3475 3476 for level1, sub_dtypes in columns_dtypes.items(): 3477 for level2, dtype in sub_dtypes.items(): 3478 self.df[level1, level2] = self.df[level1, level2].astype(dtype) 3479 3480 columns_names = pd.MultiIndex.from_frame(self.df) 3481 columns_names = columns_names.names 3482 3483 return columns_names 3484 3485 def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot'): 3486 """ 3487 Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`. 3488 3489 Args: 3490 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 3491 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. 3492 combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'. 3493 warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True. 3494 3495 Note: 3496 - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 3497 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 3498 - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html 3499 3500 """ 3501 from itertools import combinations_with_replacement 3502 compulsory_parameters_in_df = ['sigma','epsilon'] 3503 shift=0 3504 if shift_potential: 3505 shift="auto" 3506 # List which particles have sigma and epsilon values defined in pmb.df and which ones don't 3507 particles_types_with_LJ_parameters = [] 3508 non_parametrized_labels= [] 3509 for particle_type in self.get_type_map().values(): 3510 check_list=[] 3511 for key in compulsory_parameters_in_df: 3512 value_in_df=self.find_value_from_es_type(es_type=particle_type, 3513 column_name=key) 3514 check_list.append(pd.isna(value_in_df)) 3515 if any(check_list): 3516 non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 3517 column_name='label')) 3518 else: 3519 particles_types_with_LJ_parameters.append(particle_type) 3520 # Set up LJ interactions between all particle types 3521 for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 3522 particle_name1 = self.find_value_from_es_type(es_type=type_pair[0], 3523 column_name="name") 3524 particle_name2 = self.find_value_from_es_type(es_type=type_pair[1], 3525 column_name="name") 3526 lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1, 3527 particle_name2 = particle_name2, 3528 combining_rule = combining_rule) 3529 3530 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 3531 if not lj_parameters: 3532 continue 3533 espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 3534 sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 3535 cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude, 3536 offset = lj_parameters["offset"].to("reduced_length").magnitude, 3537 shift = shift) 3538 index = len(self.df) 3539 label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label") 3540 label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label") 3541 self.df.at [index, 'name'] = f'LJ: {label1}-{label2}' 3542 lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params() 3543 3544 self.add_value_to_df(index=index, 3545 key=('pmb_type',''), 3546 new_value='LennardJones') 3547 3548 self.add_value_to_df(index=index, 3549 key=('parameters_of_the_potential',''), 3550 new_value=lj_params, 3551 non_standard_value=True) 3552 if non_parametrized_labels: 3553 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.') 3554 return 3555 3556 def write_pmb_df (self, filename): 3557 ''' 3558 Writes the pyMBE dataframe into a csv file 3559 3560 Args: 3561 filename(`str`): Path to the csv file 3562 ''' 3563 3564 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence'] 3565 df = self.df.copy(deep=True) 3566 for column_name in columns_with_list_or_dict: 3567 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) 3568 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) 3569 df.fillna(pd.NA, inplace=True) 3570 df.to_csv(filename) 3571 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.
55 def __init__(self, seed, temperature=None, unit_length=None, unit_charge=None, Kw=None): 56 """ 57 Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 58 and sets up the `pmb.df` for bookkeeping. 59 60 Args: 61 temperature(`pint.Quantity`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. 62 unit_length(`pint.Quantity`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. 63 unit_charge (`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 64 Kw (`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 65 66 Note: 67 - If no `temperature` is given, a value of 298.15 K is assumed by default. 68 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 69 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 70 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 71 """ 72 # Seed and RNG 73 self.seed=seed 74 self.rng = np.random.default_rng(seed) 75 self.units=pint.UnitRegistry() 76 self.N_A=scipy.constants.N_A / self.units.mol 77 self.kB=scipy.constants.k * self.units.J / self.units.K 78 self.e=scipy.constants.e * self.units.C 79 self.set_reduced_units(unit_length=unit_length, 80 unit_charge=unit_charge, 81 temperature=temperature, 82 Kw=Kw) 83 self.setup_df() 84 self.lattice_builder = None 85 self.root = importlib.resources.files(__package__) 86 return
Initializes the pymbe_library by setting up the reduced unit system with temperature
and reduced_length
and sets up the pmb.df
for bookkeeping.
Arguments:
- temperature(
pint.Quantity
,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. - unit_length(
pint.Quantity
, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. - unit_charge (
pint.Quantity
,optional): Reduced unit of charge defined using thepmb.units
UnitRegistry. Defaults to None. - Kw (
pint.Quantity
,optional): Ionic product of water in mol^2/l^2. Defaults to None.
Note:
- If no
temperature
is given, a value of 298.15 K is assumed by default.- If no
unit_length
is given, a value of 0.355 nm is assumed by default.- If no
unit_charge
is given, a value of 1 elementary charge is assumed by default.- If no
Kw
is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default.
153 def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False): 154 """ 155 Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles. 156 157 Args: 158 particle_id1(`int`): particle_id of the type of the first particle type of the bonded particles 159 particle_id2(`int`): particle_id of the type of the second particle type of the bonded particles 160 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle whose bond types are not defined in `pmb.df`. Defaults to False. 161 162 Returns: 163 index(`int`): Row index where the bond information has been added in pmb.df. 164 """ 165 particle_name1 = self.df.loc[(self.df['particle_id']==particle_id1) & (self.df['pmb_type']=="particle")].name.values[0] 166 particle_name2 = self.df.loc[(self.df['particle_id']==particle_id2) & (self.df['pmb_type']=="particle")].name.values[0] 167 168 bond_key = self.find_bond_key(particle_name1=particle_name1, 169 particle_name2=particle_name2, 170 use_default_bond=use_default_bond) 171 if not bond_key: 172 return None 173 self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1) 174 indexs = np.where(self.df['name']==bond_key) 175 index_list = list (indexs[0]) 176 used_bond_df = self.df.loc[self.df['particle_id2'].notnull()] 177 #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict 178 used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 ) 179 used_bond_index = used_bond_df.index.to_list() 180 if not index_list: 181 return None 182 for index in index_list: 183 if index not in used_bond_index: 184 self.clean_df_row(index=int(index)) 185 self.df.at[index,'particle_id'] = particle_id1 186 self.df.at[index,'particle_id2'] = particle_id2 187 break 188 return index
Adds a bond entry on the pymbe.df
storing the particle_ids of the two bonded particles.
Arguments:
- particle_id1(
int
): particle_id of the type of the first particle type of the bonded particles - particle_id2(
int
): particle_id of the type of the second particle type of the bonded particles - use_default_bond(
bool
, optional): Controls if a bond of typedefault
is used to bond particle whose bond types are not defined inpmb.df
. Defaults to False.
Returns:
index(
int
): Row index where the bond information has been added in pmb.df.
190 def add_bonds_to_espresso(self, espresso_system) : 191 """ 192 Adds all bonds defined in `pmb.df` to `espresso_system`. 193 194 Args: 195 espresso_system(`espressomd.system.System`): system object of espressomd library 196 """ 197 198 if 'bond' in self.df["pmb_type"].values: 199 bond_df = self.df.loc[self.df ['pmb_type'] == 'bond'] 200 bond_list = bond_df.bond_object.values.tolist() 201 for bond in bond_list: 202 espresso_system.bonded_inter.add(bond) 203 else: 204 logging.warning('there are no bonds defined in pymbe.df') 205 return
Adds all bonds defined in pmb.df
to espresso_system
.
Arguments:
- espresso_system(
espressomd.system.System
): system object of espressomd library
207 def add_value_to_df(self,index,key,new_value, non_standard_value=False, overwrite=False): 208 """ 209 Adds a value to a cell in the `pmb.df` DataFrame. 210 211 Args: 212 index(`int`): index of the row to add the value to. 213 key(`str`): the column label to add the value to. 214 non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False. 215 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 216 """ 217 218 token = "#protected:" 219 220 def protect(obj): 221 if non_standard_value: 222 return token + json.dumps(obj, cls=self.NumpyEncoder) 223 return obj 224 225 def deprotect(obj): 226 if non_standard_value and isinstance(obj, str) and obj.startswith(token): 227 return json.loads(obj.removeprefix(token)) 228 return obj 229 230 # Make sure index is a scalar integer value 231 index = int(index) 232 assert isinstance(index, int), '`index` should be a scalar integer value.' 233 idx = pd.IndexSlice 234 if self.check_if_df_cell_has_a_value(index=index,key=key): 235 old_value = self.df.loc[index,idx[key]] 236 if not pd.Series([protect(old_value)]).equals(pd.Series([protect(new_value)])): 237 name=self.df.loc[index,('name','')] 238 pmb_type=self.df.loc[index,('pmb_type','')] 239 logging.debug(f"You are attempting to redefine the properties of {name} of pmb_type {pmb_type}") 240 if overwrite: 241 logging.info(f'Overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}') 242 if not overwrite: 243 logging.debug(f"pyMBE has preserved of the entry `{key}`: old_value = {old_value}. If you want to overwrite it with new_value = {new_value}, activate the switch overwrite = True ") 244 return 245 246 self.df.loc[index,idx[key]] = protect(new_value) 247 if non_standard_value: 248 self.df[key] = self.df[key].apply(deprotect) 249 return
Adds a value to a cell in the pmb.df
DataFrame.
Arguments:
- index(
int
): index of the row to add the value to. - key(
str
): the column label to add the value to. - non_standard_value(
bool
, optional): Switch to enable insertion of non-standard values, such asdict
objects. Defaults to False. - overwrite(
bool
, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
251 def assign_molecule_id(self, molecule_index): 252 """ 253 Assigns the `molecule_id` of the pmb object given by `pmb_type` 254 255 Args: 256 molecule_index(`int`): index of the current `pmb_object_type` to assign the `molecule_id` 257 Returns: 258 molecule_id(`int`): Id of the molecule 259 """ 260 261 self.clean_df_row(index=int(molecule_index)) 262 263 if self.df['molecule_id'].isnull().values.all(): 264 molecule_id = 0 265 else: 266 molecule_id = self.df['molecule_id'].max() +1 267 268 self.add_value_to_df (key=('molecule_id',''), 269 index=int(molecule_index), 270 new_value=molecule_id) 271 272 return molecule_id
Assigns the molecule_id
of the pmb object given by pmb_type
Arguments:
- molecule_index(
int
): index of the currentpmb_object_type
to assign themolecule_id
Returns:
molecule_id(
int
): Id of the molecule
274 def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system): 275 """ 276 Calculates the center of the molecule with a given molecule_id. 277 278 Args: 279 molecule_id(`int`): Id of the molecule whose center of mass is to be calculated. 280 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 281 282 Returns: 283 center_of_mass(`lst`): Coordinates of the center of mass. 284 """ 285 center_of_mass = np.zeros(3) 286 axis_list = [0,1,2] 287 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0] 288 particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"] 289 for pid in particle_id_list: 290 for axis in axis_list: 291 center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis] 292 center_of_mass = center_of_mass / len(particle_id_list) 293 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.
295 def calculate_HH(self, molecule_name, pH_list=None, pka_set=None): 296 """ 297 Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 298 for molecules with the name `molecule_name`. 299 300 Args: 301 molecule_name(`str`): name of the molecule to calculate the ideal charge for 302 pH_list(`lst`): pH-values to calculate. 303 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 304 305 Returns: 306 Z_HH(`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list` 307 308 Note: 309 - This function supports objects with pmb types: "molecule", "peptide" and "protein". 310 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 311 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 312 - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan` should be used. 313 """ 314 self._check_if_name_is_defined_in_df(name=molecule_name) 315 self._check_supported_molecule(molecule_name=molecule_name, 316 valid_pmb_types=["molecule","peptide","protein"]) 317 if pH_list is None: 318 pH_list=np.linspace(2,12,50) 319 if pka_set is None: 320 pka_set=self.get_pka_set() 321 index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 322 residue_list = self.df.at [index,('residue_list','')].copy() 323 particles_in_molecule = [] 324 for residue in residue_list: 325 list_of_particles_in_residue = self.search_particles_in_residue(residue) 326 if len(list_of_particles_in_residue) == 0: 327 logging.warning(f"The residue {residue} has no particles defined in the pyMBE DataFrame, it will be ignored.") 328 continue 329 particles_in_molecule += list_of_particles_in_residue 330 if len(particles_in_molecule) == 0: 331 return [None]*len(pH_list) 332 self.check_pka_set(pka_set=pka_set) 333 charge_number_map = self.get_charge_number_map() 334 Z_HH=[] 335 for pH_value in pH_list: 336 Z=0 337 for particle in particles_in_molecule: 338 if particle in pka_set.keys(): 339 if pka_set[particle]['acidity'] == 'acidic': 340 psi=-1 341 elif pka_set[particle]['acidity']== 'basic': 342 psi=+1 343 Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value']))) 344 else: 345 state_one_type = self.df.loc[self.df['name']==particle].state_one.es_type.values[0] 346 Z+=charge_number_map[state_one_type] 347 Z_HH.append(Z) 348 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 ofsequence
inpH_list
Note:
- This function supports objects with pmb types: "molecule", "peptide" and "protein".
- If no
pH_list
is given, 50 equispaced pH-values ranging from 2 to 12 are calculated- If no
pka_set
is given, the pKa values are taken frompmb.df
- This function should only be used for single-phase systems. For two-phase systems
pmb.calculate_HH_Donnan
should be used.
350 def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None): 351 """ 352 Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning. 353 354 Args: 355 c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 356 c_salt('float'): Salt concentration in the reservoir. 357 pH_list('lst'): List of pH-values in the reservoir. 358 pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}. 359 360 Returns: 361 {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 362 pH_system_list ('lst'): List of pH_values in the system. 363 partition_coefficients_list ('lst'): List of partition coefficients of cations. 364 365 Note: 366 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 367 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 368 - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result. 369 """ 370 if pH_list is None: 371 pH_list=np.linspace(2,12,50) 372 if pka_set is None: 373 pka_set=self.get_pka_set() 374 self.check_pka_set(pka_set=pka_set) 375 376 partition_coefficients_list = [] 377 pH_system_list = [] 378 Z_HH_Donnan={} 379 for key in c_macro: 380 Z_HH_Donnan[key] = [] 381 382 def calc_charges(c_macro, pH): 383 """ 384 Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation. 385 386 Args: 387 c_macro('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 388 pH('float'): pH-value that is used in the HH equation. 389 390 Returns: 391 charge('dict'): {"molecule_name": charge} 392 """ 393 charge = {} 394 for name in c_macro: 395 charge[name] = self.calculate_HH(name, [pH], pka_set)[0] 396 return charge 397 398 def calc_partition_coefficient(charge, c_macro): 399 """ 400 Calculates the partition coefficients of positive ions according to the ideal Donnan theory. 401 402 Args: 403 charge('dict'): {"molecule_name": charge} 404 c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 405 """ 406 nonlocal ionic_strength_res 407 charge_density = 0.0 408 for key in charge: 409 charge_density += charge[key] * c_macro[key] 410 return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude 411 412 for pH_value in pH_list: 413 # calculate the ionic strength of the reservoir 414 if pH_value <= 7.0: 415 ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 416 elif pH_value > 7.0: 417 ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt 418 419 #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations 420 #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation. 421 #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects. 422 equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro)) 423 logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq") 424 partition_coefficient = 10**logxi.root 425 426 charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient)) 427 for key in c_macro: 428 Z_HH_Donnan[key].append(charges_temp[key]) 429 430 pH_system_list.append(pH_value - np.log10(partition_coefficient)) 431 partition_coefficients_list.append(partition_coefficient) 432 433 return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list}
Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning.
Arguments:
- c_macro('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system.
- c_salt('float'): Salt concentration in the reservoir.
- pH_list('lst'): List of pH-values in the reservoir.
- pka_set('dict'): {"name": {"pka_value": pka, "acidity": acidity}}.
Returns:
{"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} pH_system_list ('lst'): List of pH_values in the system. partition_coefficients_list ('lst'): List of partition coefficients of cations.
Note:
- If no
pH_list
is given, 50 equispaced pH-values ranging from 2 to 12 are calculated- If no
pka_set
is given, the pKa values are taken frompmb.df
- If
c_macro
does not contain all charged molecules in the system, this function is likely to provide the wrong result.
435 def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff, offset): 436 """ 437 Calculates the initial bond length that is used when setting up molecules, 438 based on the minimum of the sum of bonded and short-range (LJ) interactions. 439 440 Args: 441 bond_object(`espressomd.interactions.BondedInteractions`): instance of a bond object from espressomd library 442 bond_type(`str`): label identifying the used bonded potential 443 epsilon(`pint.Quantity`): LJ epsilon of the interaction between the particles 444 sigma(`pint.Quantity`): LJ sigma of the interaction between the particles 445 cutoff(`pint.Quantity`): cutoff-radius of the LJ interaction 446 offset(`pint.Quantity`): offset of the LJ interaction 447 """ 448 def truncated_lj_potential(x, epsilon, sigma, cutoff,offset): 449 if x>cutoff: 450 return 0.0 451 else: 452 return 4*epsilon*((sigma/(x-offset))**12-(sigma/(x-offset))**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6) 453 454 epsilon_red=epsilon.to('reduced_energy').magnitude 455 sigma_red=sigma.to('reduced_length').magnitude 456 cutoff_red=cutoff.to('reduced_length').magnitude 457 offset_red=offset.to('reduced_length').magnitude 458 459 if bond_type == "harmonic": 460 r_0 = bond_object.params.get('r_0') 461 k = bond_object.params.get('k') 462 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 463 elif bond_type == "FENE": 464 r_0 = bond_object.params.get('r_0') 465 k = bond_object.params.get('k') 466 d_r_max = bond_object.params.get('d_r_max') 467 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 468 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
470 def calculate_net_charge(self, espresso_system, molecule_name, dimensionless=False): 471 ''' 472 Calculates the net charge per molecule of molecules with `name` = molecule_name. 473 Returns the net charge per molecule and a maps with the net charge per residue and molecule. 474 475 Args: 476 espresso_system(`espressomd.system.System`): system information 477 molecule_name(`str`): name of the molecule to calculate the net charge 478 dimensionless(`bool'): sets if the charge is returned with a dimension or not 479 480 Returns: 481 (`dict`): {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }} 482 483 Note: 484 - The net charge of the molecule is averaged over all molecules of type `name` 485 - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name` 486 ''' 487 self._check_supported_molecule(molecule_name=molecule_name, 488 valid_pmb_types=["molecule","protein","peptide"]) 489 490 id_map = self.get_particle_id_map(object_name=molecule_name) 491 def create_charge_map(espresso_system,id_map,label): 492 charge_number_map={} 493 for super_id in id_map[label].keys(): 494 if dimensionless: 495 net_charge=0 496 else: 497 net_charge=0 * self.units.Quantity(1,'reduced_charge') 498 for pid in id_map[label][super_id]: 499 if dimensionless: 500 net_charge+=espresso_system.part.by_id(pid).q 501 else: 502 net_charge+=espresso_system.part.by_id(pid).q * self.units.Quantity(1,'reduced_charge') 503 charge_number_map[super_id]=net_charge 504 return charge_number_map 505 net_charge_molecules=create_charge_map(label="molecule_map", 506 espresso_system=espresso_system, 507 id_map=id_map) 508 net_charge_residues=create_charge_map(label="residue_map", 509 espresso_system=espresso_system, 510 id_map=id_map) 511 if dimensionless: 512 mean_charge=np.mean(np.array(list(net_charge_molecules.values()))) 513 else: 514 mean_charge=np.mean(np.array([value.magnitude for value in net_charge_molecules.values()]))*self.units.Quantity(1,'reduced_charge') 515 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
517 def center_molecule_in_simulation_box(self, molecule_id, espresso_system): 518 """ 519 Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`. 520 521 Args: 522 molecule_id(`int`): Id of the molecule to be centered. 523 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 524 """ 525 if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0: 526 raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.") 527 center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system) 528 box_center = [espresso_system.box_l[0]/2.0, 529 espresso_system.box_l[1]/2.0, 530 espresso_system.box_l[2]/2.0] 531 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type'].isin(["molecule","protein"]))].name.values[0] 532 particle_id_list = self.get_particle_id_map(object_name=molecule_name)["all"] 533 for pid in particle_id_list: 534 es_pos = espresso_system.part.by_id(pid).pos 535 espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center 536 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.
538 def check_aminoacid_key(self, key): 539 """ 540 Checks if `key` corresponds to a valid aminoacid letter code. 541 542 Args: 543 key(`str`): key to be checked. 544 545 Returns: 546 `bool`: True if `key` is a valid aminoacid letter code, False otherwise. 547 """ 548 valid_AA_keys=['V', #'VAL' 549 'I', #'ILE' 550 'L', #'LEU' 551 'E', #'GLU' 552 'Q', #'GLN' 553 'D', #'ASP' 554 'N', #'ASN' 555 'H', #'HIS' 556 'W', #'TRP' 557 'F', #'PHE' 558 'Y', #'TYR' 559 'R', #'ARG' 560 'K', #'LYS' 561 'S', #'SER' 562 'T', #'THR' 563 'M', #'MET' 564 'A', #'ALA' 565 'G', #'GLY' 566 'P', #'PRO' 567 'C'] #'CYS' 568 if key in valid_AA_keys: 569 return True 570 else: 571 return False
Checks if key
corresponds to a valid aminoacid letter code.
Arguments:
- key(
str
): key to be checked.
Returns:
bool
: True ifkey
is a valid aminoacid letter code, False otherwise.
573 def check_dimensionality(self, variable, expected_dimensionality): 574 """ 575 Checks if the dimensionality of `variable` matches `expected_dimensionality`. 576 577 Args: 578 variable(`pint.Quantity`): Quantity to be checked. 579 expected_dimensionality(`str`): Expected dimension of the variable. 580 581 Returns: 582 (`bool`): `True` if the variable if of the expected dimensionality, `False` otherwise. 583 584 Note: 585 - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality). 586 - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"` 587 """ 588 correct_dimensionality=variable.check(f"{expected_dimensionality}") 589 if not correct_dimensionality: 590 raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}") 591 return correct_dimensionality
Checks if the dimensionality of variable
matches expected_dimensionality
.
Arguments:
- variable(
pint.Quantity
): Quantity to be checked. - expected_dimensionality(
str
): Expected dimension of the variable.
Returns:
(
bool
):True
if the variable if of the expected dimensionality,False
otherwise.
Note:
expected_dimensionality
takes dimensionality following the Pint standards docs.- For example, to check for a variable corresponding to a velocity
expected_dimensionality = "[length]/[time]"
593 def check_if_df_cell_has_a_value(self, index,key): 594 """ 595 Checks if a cell in the `pmb.df` at the specified index and column has a value. 596 597 Args: 598 index(`int`): Index of the row to check. 599 key(`str`): Column label to check. 600 601 Returns: 602 `bool`: `True` if the cell has a value, `False` otherwise. 603 """ 604 idx = pd.IndexSlice 605 import warnings 606 with warnings.catch_warnings(): 607 warnings.simplefilter("ignore") 608 return not pd.isna(self.df.loc[index, idx[key]])
Checks if a cell in the pmb.df
at the specified index and column has a value.
Arguments:
- index(
int
): Index of the row to check. - key(
str
): Column label to check.
Returns:
bool
:True
if the cell has a value,False
otherwise.
612 def check_if_metal_ion(self,key): 613 """ 614 Checks if `key` corresponds to a label of a supported metal ion. 615 616 Args: 617 key(`str`): key to be checked 618 619 Returns: 620 (`bool`): True if `key` is a supported metal ion, False otherwise. 621 """ 622 if key in self.get_metal_ions_charge_number_map().keys(): 623 return True 624 else: 625 return False
Checks if key
corresponds to a label of a supported metal ion.
Arguments:
- key(
str
): key to be checked
Returns:
(
bool
): True ifkey
is a supported metal ion, False otherwise.
627 def check_pka_set(self, pka_set): 628 """ 629 Checks that `pka_set` has the formatting expected by the pyMBE library. 630 631 Args: 632 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 633 """ 634 required_keys=['pka_value','acidity'] 635 for required_key in required_keys: 636 for pka_name, pka_entry in pka_set.items(): 637 if required_key not in pka_entry: 638 raise ValueError(f'missing a required key "{required_key}" in entry "{pka_name}" of pka_set ("{pka_entry}")') 639 return
Checks that pka_set
has the formatting expected by the pyMBE library.
Arguments:
- pka_set(
dict
): {"name" : {"pka_value": pka, "acidity": acidity}}
641 def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")): 642 """ 643 Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a pd.NA value. 644 645 Args: 646 index(`int`): Index of the row to clean. 647 columns_keys_to_clean(`list` of `str`, optional): List with the column keys to be cleaned. Defaults to [`particle_id`, `particle_id2`, `residue_id`, `molecule_id`]. 648 """ 649 for column_key in columns_keys_to_clean: 650 self.add_value_to_df(key=(column_key,''),index=index,new_value=pd.NA) 651 self.df.fillna(pd.NA, inplace=True) 652 return
Cleans the columns of pmb.df
in columns_keys_to_clean
of the row with index index
by assigning them a pd.NA value.
Arguments:
- index(
int
): Index of the row to clean. - columns_keys_to_clean(
list
ofstr
, optional): List with the column keys to be cleaned. Defaults to [particle_id
,particle_id2
,residue_id
,molecule_id
].
654 def convert_columns_to_original_format (self, df): 655 """ 656 Converts the columns of the Dataframe to the original format in pyMBE. 657 658 Args: 659 df(`DataFrame`): dataframe with pyMBE information as a string 660 661 """ 662 663 columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', ('state_one','es_type'),('state_two','es_type'),('state_one','z'),('state_two','z') ] 664 665 columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset'] 666 667 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence', 'chain_map', 'node_map'] 668 669 for column_name in columns_dtype_int: 670 df[column_name] = df[column_name].astype(pd.Int64Dtype()) 671 672 for column_name in columns_with_list_or_dict: 673 if df[column_name].isnull().all(): 674 df[column_name] = df[column_name].astype(object) 675 else: 676 df[column_name] = df[column_name].apply(lambda x: json.loads(x) if pd.notnull(x) else x) 677 678 for column_name in columns_with_units: 679 df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x) 680 681 df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x) 682 df["l0"] = df["l0"].astype(object) 683 df["pka"] = df["pka"].astype(object) 684 return df
Converts the columns of the Dataframe to the original format in pyMBE.
Arguments:
- df(
DataFrame
): dataframe with pyMBE information as a string
686 def convert_str_to_bond_object (self, bond_str): 687 """ 688 Convert a row read as a `str` to the corresponding ESPResSo bond object. 689 690 Args: 691 bond_str(`str`): string with the information of a bond object. 692 693 Returns: 694 bond_object(`obj`): ESPResSo bond object. 695 696 Note: 697 Current supported bonds are: HarmonicBond and FeneBond 698 """ 699 import espressomd.interactions 700 701 supported_bonds = ['HarmonicBond', 'FeneBond'] 702 m = re.search(r'^([A-Za-z0-9_]+)\((\{.+\})\)$', bond_str) 703 if m is None: 704 raise ValueError(f'Cannot parse bond "{bond_str}"') 705 bond = m.group(1) 706 if bond not in supported_bonds: 707 raise NotImplementedError(f"Bond type '{bond}' currently not implemented in pyMBE, accepted types are {supported_bonds}") 708 params = json.loads(m.group(2)) 709 bond_id = params.pop("bond_id") 710 bond_object = getattr(espressomd.interactions, bond)(**params) 711 bond_object._bond_id = bond_id 712 return bond_object
Convert a row read as a str
to the corresponding ESPResSo bond object.
Arguments:
- bond_str(
str
): string with the information of a bond object.
Returns:
bond_object(
obj
): ESPResSo bond object.
Note:
Current supported bonds are: HarmonicBond and FeneBond
714 def copy_df_entry(self, name, column_name, number_of_copies): 715 ''' 716 Creates 'number_of_copies' of a given 'name' in `pymbe.df`. 717 718 Args: 719 name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df` 720 column_name(`str`): Column name to use as a filter. 721 number_of_copies(`int`): number of copies of `name` to be created. 722 723 Note: 724 - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 725 ''' 726 727 valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ] 728 if column_name not in valid_column_names: 729 raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}") 730 df_by_name = self.df.loc[self.df.name == name] 731 if number_of_copies != 1: 732 if df_by_name[column_name].isnull().values.any(): 733 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True) 734 else: 735 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 736 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 737 df_by_name_repeated[column_name] = pd.NA 738 # Concatenate the new particle rows to `df` 739 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 740 else: 741 if not df_by_name[column_name].isnull().values.any(): 742 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 743 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 744 df_by_name_repeated[column_name] = pd.NA 745 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 746 return
Creates 'number_of_copies' of a given 'name' in pymbe.df
.
Arguments:
- name(
str
): Label of the particle/residue/molecule type to be created.name
must be defined inpmb.df
- column_name(
str
): Column name to use as a filter. - number_of_copies(
int
): number of copies ofname
to be created.
Note:
- Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id"
748 def create_added_salt(self, espresso_system, cation_name, anion_name, c_salt): 749 """ 750 Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`. 751 752 Args: 753 espresso_system(`espressomd.system.System`): instance of an espresso system object. 754 cation_name(`str`): `name` of a particle with a positive charge. 755 anion_name(`str`): `name` of a particle with a negative charge. 756 c_salt(`float`): Salt concentration. 757 758 Returns: 759 c_salt_calculated(`float`): Calculated salt concentration added to `espresso_system`. 760 """ 761 for name in [cation_name, anion_name]: 762 if not self._check_if_name_is_defined_in_df(name=name): 763 logging.warning(f"Object with name '{name}' is not defined in the DataFrame, no ions will be created.") 764 return 765 self._check_if_name_has_right_type(name=cation_name, 766 expected_pmb_type="particle") 767 self._check_if_name_has_right_type(name=anion_name, 768 expected_pmb_type="particle") 769 cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0] 770 anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0] 771 if cation_name_charge <= 0: 772 raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge) 773 if anion_name_charge >= 0: 774 raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge) 775 # Calculate the number of ions in the simulation box 776 volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3') 777 if c_salt.check('[substance] [length]**-3'): 778 N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude) 779 c_salt_calculated=N_ions/(volume*self.N_A) 780 elif c_salt.check('[length]**-3'): 781 N_ions= int((volume*c_salt.to('reduced_length**-3')).magnitude) 782 c_salt_calculated=N_ions/volume 783 else: 784 raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt) 785 N_cation = N_ions*abs(anion_name_charge) 786 N_anion = N_ions*abs(cation_name_charge) 787 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation) 788 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion) 789 if c_salt_calculated.check('[substance] [length]**-3'): 790 logging.info(f"added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions") 791 elif c_salt_calculated.check('[length]**-3'): 792 logging.info(f"added salt concentration of {c_salt_calculated.to('reduced_length**-3')} given by {N_cation} cations and {N_anion} anions") 793 return c_salt_calculated
Creates a c_salt
concentration of cation_name
and anion_name
ions into the espresso_system
.
Arguments:
- espresso_system(
espressomd.system.System
): instance of an espresso system object. - cation_name(
str
):name
of a particle with a positive charge. - anion_name(
str
):name
of a particle with a negative charge. - c_salt(
float
): Salt concentration.
Returns:
c_salt_calculated(
float
): Calculated salt concentration added toespresso_system
.
795 def create_bond_in_espresso(self, bond_type, bond_parameters): 796 ''' 797 Creates either a harmonic or a FENE bond in ESPResSo 798 799 Args: 800 bond_type(`str`): label to identify the potential to model the bond. 801 bond_parameters(`dict`): parameters of the potential of the bond. 802 803 Note: 804 Currently, only HARMONIC and FENE bonds are supported. 805 806 For a HARMONIC bond the dictionary must contain: 807 808 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 809 using the `pmb.units` UnitRegistry. 810 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 811 the `pmb.units` UnitRegistry. 812 813 For a FENE bond the dictionary must additionally contain: 814 815 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 816 units of length using the `pmb.units` UnitRegistry. Default 'None'. 817 818 Returns: 819 bond_object (`obj`): an ESPResSo bond object 820 ''' 821 from espressomd import interactions 822 823 valid_bond_types = ["harmonic", "FENE"] 824 825 if 'k' in bond_parameters: 826 bond_magnitude = bond_parameters['k'].to('reduced_energy / reduced_length**2') 827 else: 828 raise ValueError("Magnitude of the potential (k) is missing") 829 830 if bond_type == 'harmonic': 831 if 'r_0' in bond_parameters: 832 bond_length = bond_parameters['r_0'].to('reduced_length') 833 else: 834 raise ValueError("Equilibrium bond length (r_0) is missing") 835 bond_object = interactions.HarmonicBond(k = bond_magnitude.magnitude, 836 r_0 = bond_length.magnitude) 837 elif bond_type == 'FENE': 838 if 'r_0' in bond_parameters: 839 bond_length = bond_parameters['r_0'].to('reduced_length').magnitude 840 else: 841 logging.warning("no value provided for r_0. Defaulting to r_0 = 0") 842 bond_length=0 843 if 'd_r_max' in bond_parameters: 844 max_bond_stret = bond_parameters['d_r_max'].to('reduced_length') 845 else: 846 raise ValueError("Maximal stretching length (d_r_max) is missing") 847 bond_object = interactions.FeneBond(r_0 = bond_length, 848 k = bond_magnitude.magnitude, 849 d_r_max = max_bond_stret.magnitude) 850 else: 851 raise NotImplementedError(f"Bond type '{bond_type}' currently not implemented in pyMBE, accepted types are {valid_bond_types}") 852 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
855 def create_counterions(self, object_name, cation_name, anion_name, espresso_system): 856 """ 857 Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`. 858 859 Args: 860 object_name(`str`): `name` of a pymbe object. 861 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 862 cation_name(`str`): `name` of a particle with a positive charge. 863 anion_name(`str`): `name` of a particle with a negative charge. 864 865 Returns: 866 counterion_number(`dict`): {"name": number} 867 868 Note: 869 This function currently does not support the creation of counterions for hydrogels. 870 """ 871 for name in [object_name, cation_name, anion_name]: 872 if not self._check_if_name_is_defined_in_df(name=name): 873 logging.warning(f"Object with name '{name}' is not defined in the DataFrame, no counterions will be created.") 874 return 875 for name in [cation_name, anion_name]: 876 self._check_if_name_has_right_type(name=name, expected_pmb_type="particle") 877 self._check_supported_molecule(molecule_name=object_name, 878 valid_pmb_types=["molecule","peptide","protein"]) 879 880 881 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.iloc[0] 882 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.iloc[0] 883 object_ids = self.get_particle_id_map(object_name=object_name)["all"] 884 counterion_number={} 885 object_charge={} 886 for name in ['positive', 'negative']: 887 object_charge[name]=0 888 for id in object_ids: 889 if espresso_system.part.by_id(id).q > 0: 890 object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q )) 891 elif espresso_system.part.by_id(id).q < 0: 892 object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q )) 893 if object_charge['positive'] % abs(anion_charge) == 0: 894 counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge)) 895 else: 896 raise ValueError('The number of positive charges in the pmb_object must be divisible by the charge of the anion') 897 if object_charge['negative'] % abs(cation_charge) == 0: 898 counterion_number[cation_name]=int(object_charge['negative']/cation_charge) 899 else: 900 raise ValueError('The number of negative charges in the pmb_object must be divisible by the charge of the cation') 901 if counterion_number[cation_name] > 0: 902 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name]) 903 else: 904 counterion_number[cation_name]=0 905 if counterion_number[anion_name] > 0: 906 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name]) 907 else: 908 counterion_number[anion_name] = 0 909 logging.info('the following counter-ions have been created: ') 910 for name in counterion_number.keys(): 911 logging.info(f'Ion type: {name} created number: {counterion_number[name]}') 912 return counterion_number
Creates particles of cation_name
and anion_name
in espresso_system
to counter the net charge of pmb_object
.
Arguments:
- object_name(
str
):name
of a pymbe object. - espresso_system(
espressomd.system.System
): Instance of a system object from the espressomd library. - cation_name(
str
):name
of a particle with a positive charge. - anion_name(
str
):name
of a particle with a negative charge.
Returns:
counterion_number(dict
): {"name": number}
Note:
This function currently does not support the creation of counterions for hydrogels.
914 def create_hydrogel(self, name, espresso_system): 915 """ 916 creates the hydrogel `name` in espresso_system 917 Args: 918 name(`str`): Label of the hydrogel to be created. `name` must be defined in the `pmb.df` 919 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 920 921 Returns: 922 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],...}} 923 """ 924 if not self._check_if_name_is_defined_in_df(name=name): 925 logging.warning(f"Hydrogel with name '{name}' is not defined in the DataFrame, no hydrogel will be created.") 926 return 927 self._check_if_name_has_right_type(name=name, 928 expected_pmb_type="hydrogel") 929 hydrogel_info={"name":name, "chains":{}, "nodes":{}} 930 # placing nodes 931 node_positions = {} 932 node_topology = self.df[self.df["name"]==name]["node_map"].iloc[0] 933 for node_info in node_topology: 934 node_index = node_info["lattice_index"] 935 node_name = node_info["particle_name"] 936 node_pos, node_id = self.create_hydrogel_node(self.format_node(node_index), node_name, espresso_system) 937 hydrogel_info["nodes"][self.format_node(node_index)]=node_id 938 node_positions[node_id[0]]=node_pos 939 940 # Placing chains between nodes 941 # Looping over all the 16 chains 942 chain_topology = self.df[self.df["name"]==name]["chain_map"].iloc[0] 943 for chain_info in chain_topology: 944 node_s = chain_info["node_start"] 945 node_e = chain_info["node_end"] 946 molecule_info = self.create_hydrogel_chain(node_s, node_e, node_positions, espresso_system) 947 for molecule_id in molecule_info: 948 hydrogel_info["chains"][molecule_id] = molecule_info[molecule_id] 949 hydrogel_info["chains"][molecule_id]["node_start"]=node_s 950 hydrogel_info["chains"][molecule_id]["node_end"]=node_e 951 return hydrogel_info
creates the hydrogel name
in espresso_system
Arguments:
- name(
str
): Label of the hydrogel to be created.name
must be defined in 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],...}}
953 def create_hydrogel_chain(self, node_start, node_end, node_positions, espresso_system): 954 """ 955 Creates a chain between two nodes of a hydrogel. 956 957 Args: 958 node_start(`str`): name of the starting node particle at which the first residue of the chain will be attached. 959 node_end(`str`): name of the ending node particle at which the last residue of the chain will be attached. 960 node_positions(`dict`): dictionary with the positions of the nodes. The keys are the node names and the values are the positions. 961 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 962 963 Note: 964 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. 965 The chain will be placed in the direction of the vector between `node_start` and `node_end`. 966 """ 967 if self.lattice_builder is None: 968 raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.") 969 970 molecule_name = "chain_"+node_start+"_"+node_end 971 sequence = self.df[self.df['name']==molecule_name].residue_list.values [0] 972 assert len(sequence) != 0 and not isinstance(sequence, str) 973 assert len(sequence) == self.lattice_builder.mpc 974 975 key, reverse = self.lattice_builder._get_node_vector_pair(node_start, node_end) 976 assert node_start != node_end or sequence == sequence[::-1], \ 977 (f"chain cannot be defined between '{node_start}' and '{node_end}' since it " 978 "would form a loop with a non-symmetric sequence (under-defined stereocenter)") 979 980 if reverse: 981 sequence = sequence[::-1] 982 983 node_start_pos = np.array(list(int(x) for x in node_start.strip('[]').split()))*0.25*self.lattice_builder.box_l 984 node_end_pos = np.array(list(int(x) for x in node_end.strip('[]').split()))*0.25*self.lattice_builder.box_l 985 node1 = espresso_system.part.select(lambda p: (p.pos == node_start_pos).all()).id 986 node2 = espresso_system.part.select(lambda p: (p.pos == node_end_pos).all()).id 987 988 if not node1[0] in node_positions or not node2[0] in node_positions: 989 raise ValueError("Set node position before placing a chain between them") 990 991 # Finding a backbone vector between node_start and node_end 992 vec_between_nodes = np.array(node_positions[node2[0]]) - np.array(node_positions[node1[0]]) 993 vec_between_nodes = vec_between_nodes - self.lattice_builder.box_l * np.round(vec_between_nodes/self.lattice_builder.box_l) 994 backbone_vector = list(vec_between_nodes/(self.lattice_builder.mpc + 1)) 995 node_start_name = self.df[(self.df["particle_id"]==node1[0]) & (self.df["pmb_type"]=="particle")]["name"].values[0] 996 first_res_name = self.df[(self.df["pmb_type"]=="residue") & (self.df["name"]==sequence[0])]["central_bead"].values[0] 997 l0 = self.get_bond_length(node_start_name, first_res_name, hard_check=True) 998 chain_molecule_info = self.create_molecule( 999 name=molecule_name, # Use the name defined earlier 1000 number_of_molecules=1, # Creating one chain 1001 espresso_system=espresso_system, 1002 list_of_first_residue_positions=[list(np.array(node_positions[node1[0]]) + np.array(backbone_vector))],#Start at the first node 1003 backbone_vector=np.array(backbone_vector)/l0, 1004 use_default_bond=False # Use defaut bonds between monomers 1005 ) 1006 # Collecting ids of beads of the chain/molecule 1007 chain_ids = [] 1008 residue_ids = [] 1009 for molecule_id in chain_molecule_info: 1010 for residue_id in chain_molecule_info[molecule_id]: 1011 residue_ids.append(residue_id) 1012 bead_id = chain_molecule_info[molecule_id][residue_id]['central_bead_id'] 1013 chain_ids.append(bead_id) 1014 1015 self.lattice_builder.chains[key] = sequence 1016 # Search bonds between nodes and chain ends 1017 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] 1018 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] 1019 bond_node1_first_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_start], 1020 particle_name2 = BeadType_near_to_node_start, 1021 hard_check=False, 1022 use_default_bond=False) 1023 bond_node2_last_monomer = self.search_bond(particle_name1 = self.lattice_builder.nodes[node_end], 1024 particle_name2 = BeadType_near_to_node_end, 1025 hard_check=False, 1026 use_default_bond=False) 1027 1028 espresso_system.part.by_id(node1[0]).add_bond((bond_node1_first_monomer, chain_ids[0])) 1029 espresso_system.part.by_id(node2[0]).add_bond((bond_node2_last_monomer, chain_ids[-1])) 1030 # Add bonds to data frame 1031 bond_index1 = self.add_bond_in_df(particle_id1=node1[0], 1032 particle_id2=chain_ids[0], 1033 use_default_bond=False) 1034 self.add_value_to_df(key=('molecule_id',''), 1035 index=int(bond_index1), 1036 new_value=molecule_id, 1037 overwrite=True) 1038 self.add_value_to_df(key=('residue_id',''), 1039 index=int (bond_index1), 1040 new_value=residue_ids[0], 1041 overwrite=True) 1042 bond_index2 = self.add_bond_in_df(particle_id1=node2[0], 1043 particle_id2=chain_ids[-1], 1044 use_default_bond=False) 1045 self.add_value_to_df(key=('molecule_id',''), 1046 index=int(bond_index2), 1047 new_value=molecule_id, 1048 overwrite=True) 1049 self.add_value_to_df(key=('residue_id',''), 1050 index=int (bond_index2), 1051 new_value=residue_ids[-1], 1052 overwrite=True) 1053 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_start
andnode_end
.
1055 def create_hydrogel_node(self, node_index, node_name, espresso_system): 1056 """ 1057 Set a node residue type. 1058 1059 Args: 1060 node_index(`str`): Lattice node index in the form of a string, e.g. "[0 0 0]". 1061 node_name(`str`): name of the node particle defined in pyMBE. 1062 Returns: 1063 node_position(`list`): Position of the node in the lattice. 1064 p_id(`int`): Particle ID of the node. 1065 """ 1066 if self.lattice_builder is None: 1067 raise ValueError("LatticeBuilder is not initialized. Use `initialize_lattice_builder` first.") 1068 1069 node_position = np.array(list(int(x) for x in node_index.strip('[]').split()))*0.25*self.lattice_builder.box_l 1070 p_id = self.create_particle(name = node_name, 1071 espresso_system=espresso_system, 1072 number_of_particles=1, 1073 position = [node_position]) 1074 key = self.lattice_builder._get_node_by_label(node_index) 1075 self.lattice_builder.nodes[key] = node_name 1076 1077 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.
1079 def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, backbone_vector=None, use_default_bond=False): 1080 """ 1081 Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1082 1083 Args: 1084 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 1085 espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library. 1086 number_of_molecules(`int`): Number of molecules of type `name` to be created. 1087 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. 1088 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. 1089 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df` 1090 1091 Returns: 1092 molecules_info(`dict`): {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 1093 1094 Note: 1095 Despite its name, this function can be used to create both molecules and peptides. 1096 """ 1097 if not self._check_if_name_is_defined_in_df(name=name): 1098 logging.warning(f"Molecule with name '{name}' is not defined in the pyMBE DataFrame, no molecule will be created.") 1099 return {} 1100 if number_of_molecules <= 0: 1101 return {} 1102 if list_of_first_residue_positions is not None: 1103 for item in list_of_first_residue_positions: 1104 if not isinstance(item, list): 1105 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.") 1106 elif len(item) != 3: 1107 raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.") 1108 1109 if len(list_of_first_residue_positions) != number_of_molecules: 1110 raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}") 1111 1112 # This function works for both molecules and peptides 1113 if not self._check_if_name_has_right_type(name=name, expected_pmb_type="molecule", hard_check=False): 1114 self._check_if_name_has_right_type(name=name, expected_pmb_type="peptide") 1115 1116 # Generate an arbitrary random unit vector 1117 if backbone_vector is None: 1118 backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0], 1119 radius=1, 1120 n_samples=1, 1121 on_surface=True)[0] 1122 else: 1123 backbone_vector = np.array(backbone_vector) 1124 1125 1126 1127 first_residue = True 1128 molecules_info = {} 1129 residue_list = self.df[self.df['name']==name].residue_list.values [0] 1130 self.copy_df_entry(name=name, 1131 column_name='molecule_id', 1132 number_of_copies=number_of_molecules) 1133 1134 molecules_index = np.where(self.df['name']==name) 1135 molecule_index_list =list(molecules_index[0])[-number_of_molecules:] 1136 pos_index = 0 1137 for molecule_index in molecule_index_list: 1138 molecule_id = self.assign_molecule_id(molecule_index=molecule_index) 1139 molecules_info[molecule_id] = {} 1140 for residue in residue_list: 1141 if first_residue: 1142 if list_of_first_residue_positions is None: 1143 residue_position = None 1144 else: 1145 for item in list_of_first_residue_positions: 1146 residue_position = [np.array(list_of_first_residue_positions[pos_index])] 1147 1148 residues_info = self.create_residue(name=residue, 1149 espresso_system=espresso_system, 1150 central_bead_position=residue_position, 1151 use_default_bond= use_default_bond, 1152 backbone_vector=backbone_vector) 1153 residue_id = next(iter(residues_info)) 1154 # Add the correct molecule_id to all particles in the residue 1155 for index in self.df[self.df['residue_id']==residue_id].index: 1156 self.add_value_to_df(key=('molecule_id',''), 1157 index=int (index), 1158 new_value=molecule_id, 1159 overwrite=True) 1160 central_bead_id = residues_info[residue_id]['central_bead_id'] 1161 previous_residue = residue 1162 residue_position = espresso_system.part.by_id(central_bead_id).pos 1163 previous_residue_id = central_bead_id 1164 first_residue = False 1165 else: 1166 previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0] 1167 new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0] 1168 bond = self.search_bond(particle_name1=previous_central_bead_name, 1169 particle_name2=new_central_bead_name, 1170 hard_check=True, 1171 use_default_bond=use_default_bond) 1172 l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 1173 particle_name2=new_central_bead_name, 1174 hard_check=True, 1175 use_default_bond=use_default_bond) 1176 1177 residue_position = residue_position+backbone_vector*l0 1178 residues_info = self.create_residue(name=residue, 1179 espresso_system=espresso_system, 1180 central_bead_position=[residue_position], 1181 use_default_bond= use_default_bond, 1182 backbone_vector=backbone_vector) 1183 residue_id = next(iter(residues_info)) 1184 for index in self.df[self.df['residue_id']==residue_id].index: 1185 self.add_value_to_df(key=('molecule_id',''), 1186 index=int (index), 1187 new_value=molecule_id, 1188 overwrite=True) 1189 central_bead_id = residues_info[residue_id]['central_bead_id'] 1190 espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id)) 1191 bond_index = self.add_bond_in_df(particle_id1=central_bead_id, 1192 particle_id2=previous_residue_id, 1193 use_default_bond=use_default_bond) 1194 self.add_value_to_df(key=('molecule_id',''), 1195 index=int (bond_index), 1196 new_value=molecule_id, 1197 overwrite=True) 1198 previous_residue_id = central_bead_id 1199 previous_residue = residue 1200 molecules_info[molecule_id][residue_id] = residues_info[residue_id] 1201 first_residue = True 1202 pos_index+=1 1203 1204 return molecules_info
Creates number_of_molecules
molecule of type name
into espresso_system
and bookkeeps them into pmb.df
.
Arguments:
- name(
str
): Label of the molecule type to be created.name
must be defined inpmb.df
- espresso_system(
espressomd.system.System
): Instance of a system object from espressomd library. - number_of_molecules(
int
): Number of molecules of typename
to be created. - list_of_first_residue_positions(
list
, optional): List of coordinates where the central bead of the first_residue_position will be created, random by default. - backbone_vector(
list
offloat
): Backbone vector of the molecule, random by default. Central beads of the residues in theresidue_list
are placed along this vector. - use_default_bond(
bool
, optional): Controls if a bond of typedefault
is 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.
1206 def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False): 1207 """ 1208 Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`. 1209 1210 Args: 1211 name(`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 1212 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 1213 number_of_particles(`int`): Number of particles to be created. 1214 position(list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None. 1215 fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False. 1216 Returns: 1217 created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`. 1218 """ 1219 if number_of_particles <=0: 1220 return [] 1221 if not self._check_if_name_is_defined_in_df(name=name): 1222 logging.warning(f"Particle with name '{name}' is not defined in the pyMBE DataFrame, no particle will be created.") 1223 return [] 1224 self._check_if_name_has_right_type(name=name, 1225 expected_pmb_type="particle") 1226 # Copy the data of the particle `number_of_particles` times in the `df` 1227 self.copy_df_entry(name=name, 1228 column_name='particle_id', 1229 number_of_copies=number_of_particles) 1230 # Get information from the particle type `name` from the df 1231 z = self.df.loc[self.df['name']==name].state_one.z.values[0] 1232 z = 0. if z is None else z 1233 es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 1234 # Get a list of the index in `df` corresponding to the new particles to be created 1235 index = np.where(self.df['name']==name) 1236 index_list =list(index[0])[-number_of_particles:] 1237 # Create the new particles into `espresso_system` 1238 created_pid_list=[] 1239 for index in range (number_of_particles): 1240 df_index=int (index_list[index]) 1241 self.clean_df_row(index=df_index) 1242 if position is None: 1243 particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l) 1244 else: 1245 particle_position = position[index] 1246 if len(espresso_system.part.all()) == 0: 1247 bead_id = 0 1248 else: 1249 bead_id = max (espresso_system.part.all().id) + 1 1250 created_pid_list.append(bead_id) 1251 kwargs = dict(id=bead_id, pos=particle_position, type=es_type, q=z) 1252 if fix: 1253 kwargs["fix"] = 3 * [fix] 1254 espresso_system.part.add(**kwargs) 1255 self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id) 1256 return created_pid_list
Creates number_of_particles
particles of type name
into espresso_system
and bookkeeps them into pymbe.df
.
Arguments:
- name(
str
): Label of the particle type to be created.name
must be aparticle
defined 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(
list
offloat
): List with the ids of the particles created intoespresso_system
.
1258 def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False, backbone_vector=None): 1259 """ 1260 Creates all `particle`s associated to `pmb object` into `espresso` a number of times equal to `number_of_objects`. 1261 1262 Args: 1263 name(`str`): Unique label of the `pmb object` to be created. 1264 number_of_objects(`int`): Number of `pmb object`s to be created. 1265 espresso_system(`espressomd.system.System`): Instance of an espresso system object from espressomd library. 1266 position(`list`): Coordinates where the particles should be created. 1267 use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`. 1268 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. 1269 1270 Note: 1271 - If no `position` is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length. 1272 """ 1273 pmb_type = self._check_supported_molecule(molecule_name=name, 1274 valid_pmb_types=["molecule","particle","residue","peptide"]) 1275 1276 if pmb_type == 'particle': 1277 self.create_particle(name=name, 1278 number_of_particles=number_of_objects, 1279 espresso_system=espresso_system, 1280 position=position) 1281 elif pmb_type == 'residue': 1282 self.create_residue(name=name, 1283 espresso_system=espresso_system, 1284 central_bead_position=position, 1285 use_default_bond=use_default_bond, 1286 backbone_vector=backbone_vector) 1287 elif pmb_type in ['molecule','peptide']: 1288 self.create_molecule(name=name, 1289 number_of_molecules=number_of_objects, 1290 espresso_system=espresso_system, 1291 use_default_bond=use_default_bond, 1292 list_of_first_residue_positions=position, 1293 backbone_vector=backbone_vector) 1294 return
Creates all particle
s associated to pmb object
into espresso
a number of times equal to number_of_objects
.
Arguments:
- name(
str
): Unique label of thepmb object
to be created. - number_of_objects(
int
): Number ofpmb object
s to be created. - espresso_system(
espressomd.system.System
): Instance of an espresso system object from espressomd library. - position(
list
): Coordinates where the particles should be created. - use_default_bond(
bool
,optional): Controls if adefault
bond is used to bond particles with undefined bonds inpmb.df
. Defaults toFalse
. - backbone_vector(
list
offloat
): Backbone vector of the molecule, random by default. Central beads of the residues in theresidue_list
are placed along this vector.
Note:
- If no
position
is given, particles will be created in random positions. For bonded particles, they will be created at a distance equal to the bond length.
1296 def create_protein(self, name, number_of_proteins, espresso_system, topology_dict): 1297 """ 1298 Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions` 1299 1300 Args: 1301 name(`str`): Label of the protein to be created. 1302 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 1303 number_of_proteins(`int`): Number of proteins to be created. 1304 positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}} 1305 """ 1306 1307 if number_of_proteins <=0: 1308 return 1309 if not self._check_if_name_is_defined_in_df(name=name): 1310 logging.warning(f"Protein with name '{name}' is not defined in the pyMBE DataFrame, no protein will be created.") 1311 return 1312 self._check_if_name_has_right_type(name=name, 1313 expected_pmb_type="protein") 1314 self.copy_df_entry(name=name, 1315 column_name='molecule_id', 1316 number_of_copies=number_of_proteins) 1317 protein_index = np.where(self.df['name']==name) 1318 protein_index_list =list(protein_index[0])[-number_of_proteins:] 1319 1320 box_half=espresso_system.box_l[0]/2.0 1321 1322 for molecule_index in protein_index_list: 1323 1324 molecule_id = self.assign_molecule_id(molecule_index=molecule_index) 1325 1326 protein_center = self.generate_coordinates_outside_sphere(radius = 1, 1327 max_dist=box_half, 1328 n_samples=1, 1329 center=[box_half]*3)[0] 1330 1331 for residue in topology_dict.keys(): 1332 1333 residue_name = re.split(r'\d+', residue)[0] 1334 residue_number = re.split(r'(\d+)', residue)[1] 1335 residue_position = topology_dict[residue]['initial_pos'] 1336 position = residue_position + protein_center 1337 1338 particle_id = self.create_particle(name=residue_name, 1339 espresso_system=espresso_system, 1340 number_of_particles=1, 1341 position=[position], 1342 fix = True) 1343 1344 index = self.df[self.df['particle_id']==particle_id[0]].index.values[0] 1345 self.add_value_to_df(key=('residue_id',''), 1346 index=int (index), 1347 new_value=int(residue_number), 1348 overwrite=True) 1349 1350 self.add_value_to_df(key=('molecule_id',''), 1351 index=int (index), 1352 new_value=molecule_id, 1353 overwrite=True) 1354 1355 return
Creates number_of_proteins
molecules of type name
into espresso_system
at the coordinates in positions
Arguments:
- name(
str
): Label of the protein to be created. - espresso_system(
espressomd.system.System
): Instance of a system object from the espressomd library. - number_of_proteins(
int
): Number of proteins to be created. - positions(
dict
): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}}
1357 def create_residue(self, name, espresso_system, central_bead_position=None,use_default_bond=False, backbone_vector=None): 1358 """ 1359 Creates a residue of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1360 1361 Args: 1362 name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df` 1363 espresso_system(`espressomd.system.System`): Instance of a system object from espressomd library. 1364 central_bead_position(`list` of `float`): Position of the central bead. 1365 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` 1366 backbone_vector(`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`. 1367 1368 Returns: 1369 residues_info(`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1370 """ 1371 if not self._check_if_name_is_defined_in_df(name=name): 1372 logging.warning(f"Residue with name '{name}' is not defined in the pyMBE DataFrame, no residue will be created.") 1373 return 1374 self._check_if_name_has_right_type(name=name, 1375 expected_pmb_type="residue") 1376 1377 # Copy the data of a residue in the `df 1378 self.copy_df_entry(name=name, 1379 column_name='residue_id', 1380 number_of_copies=1) 1381 residues_index = np.where(self.df['name']==name) 1382 residue_index_list =list(residues_index[0])[-1:] 1383 # search for defined particle and residue names 1384 particle_and_residue_df = self.df.loc[(self.df['pmb_type']== "particle") | (self.df['pmb_type']== "residue")] 1385 particle_and_residue_names = particle_and_residue_df["name"].tolist() 1386 for residue_index in residue_index_list: 1387 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1388 for side_chain_element in side_chain_list: 1389 if side_chain_element not in particle_and_residue_names: 1390 raise ValueError (f"{side_chain_element} is not defined") 1391 # Internal bookkepping of the residue info (important for side-chain residues) 1392 # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1393 residues_info={} 1394 for residue_index in residue_index_list: 1395 self.clean_df_row(index=int(residue_index)) 1396 # Assign a residue_id 1397 if self.df['residue_id'].isnull().all(): 1398 residue_id=0 1399 else: 1400 residue_id = self.df['residue_id'].max() + 1 1401 self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id) 1402 # create the principal bead 1403 central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0] 1404 central_bead_id = self.create_particle(name=central_bead_name, 1405 espresso_system=espresso_system, 1406 position=central_bead_position, 1407 number_of_particles = 1)[0] 1408 central_bead_position=espresso_system.part.by_id(central_bead_id).pos 1409 #assigns same residue_id to the central_bead particle created. 1410 index = self.df[self.df['particle_id']==central_bead_id].index.values[0] 1411 self.df.at [index,'residue_id'] = residue_id 1412 # Internal bookkeeping of the central bead id 1413 residues_info[residue_id]={} 1414 residues_info[residue_id]['central_bead_id']=central_bead_id 1415 # create the lateral beads 1416 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1417 side_chain_beads_ids = [] 1418 for side_chain_element in side_chain_list: 1419 pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 1420 if pmb_type == 'particle': 1421 bond = self.search_bond(particle_name1=central_bead_name, 1422 particle_name2=side_chain_element, 1423 hard_check=True, 1424 use_default_bond=use_default_bond) 1425 l0 = self.get_bond_length(particle_name1=central_bead_name, 1426 particle_name2=side_chain_element, 1427 hard_check=True, 1428 use_default_bond=use_default_bond) 1429 1430 if backbone_vector is None: 1431 bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1432 radius=l0, 1433 n_samples=1, 1434 on_surface=True)[0] 1435 else: 1436 bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=np.array(backbone_vector), 1437 magnitude=l0) 1438 1439 side_bead_id = self.create_particle(name=side_chain_element, 1440 espresso_system=espresso_system, 1441 position=[bead_position], 1442 number_of_particles=1)[0] 1443 index = self.df[self.df['particle_id']==side_bead_id].index.values[0] 1444 self.add_value_to_df(key=('residue_id',''), 1445 index=int (index), 1446 new_value=residue_id, 1447 overwrite=True) 1448 side_chain_beads_ids.append(side_bead_id) 1449 espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id)) 1450 index = self.add_bond_in_df(particle_id1=central_bead_id, 1451 particle_id2=side_bead_id, 1452 use_default_bond=use_default_bond) 1453 self.add_value_to_df(key=('residue_id',''), 1454 index=int (index), 1455 new_value=residue_id, 1456 overwrite=True) 1457 1458 elif pmb_type == 'residue': 1459 central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0] 1460 bond = self.search_bond(particle_name1=central_bead_name, 1461 particle_name2=central_bead_side_chain, 1462 hard_check=True, 1463 use_default_bond=use_default_bond) 1464 l0 = self.get_bond_length(particle_name1=central_bead_name, 1465 particle_name2=central_bead_side_chain, 1466 hard_check=True, 1467 use_default_bond=use_default_bond) 1468 if backbone_vector is None: 1469 residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1470 radius=l0, 1471 n_samples=1, 1472 on_surface=True)[0] 1473 else: 1474 residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1475 magnitude=l0) 1476 lateral_residue_info = self.create_residue(name=side_chain_element, 1477 espresso_system=espresso_system, 1478 central_bead_position=[residue_position], 1479 use_default_bond=use_default_bond) 1480 lateral_residue_dict=list(lateral_residue_info.values())[0] 1481 central_bead_side_chain_id=lateral_residue_dict['central_bead_id'] 1482 lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids'] 1483 residue_id_side_chain=list(lateral_residue_info.keys())[0] 1484 # Change the residue_id of the residue in the side chain to the one of the bigger residue 1485 index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0] 1486 self.add_value_to_df(key=('residue_id',''), 1487 index=int(index), 1488 new_value=residue_id, 1489 overwrite=True) 1490 # Change the residue_id of the particles in the residue in the side chain 1491 side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids 1492 for particle_id in side_chain_beads_ids: 1493 index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0] 1494 self.add_value_to_df(key=('residue_id',''), 1495 index=int (index), 1496 new_value=residue_id, 1497 overwrite=True) 1498 espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id)) 1499 index = self.add_bond_in_df(particle_id1=central_bead_id, 1500 particle_id2=central_bead_side_chain_id, 1501 use_default_bond=use_default_bond) 1502 self.add_value_to_df(key=('residue_id',''), 1503 index=int (index), 1504 new_value=residue_id, 1505 overwrite=True) 1506 # Change the residue_id of the bonds in the residues in the side chain to the one of the bigger residue 1507 for index in self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='bond') ].index: 1508 self.add_value_to_df(key=('residue_id',''), 1509 index=int(index), 1510 new_value=residue_id, 1511 overwrite=True) 1512 # Internal bookkeeping of the side chain beads ids 1513 residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids 1514 return residues_info
Creates a residue of type name
into espresso_system
and bookkeeps them into pmb.df
.
Arguments:
- name(
str
): Label of the residue type to be created.name
must be defined inpmb.df
- espresso_system(
espressomd.system.System
): Instance of a system object from espressomd library. - central_bead_position(
list
offloat
): Position of the central bead. - use_default_bond(
bool
): Switch to control if a bond of typedefault
is used to bond a particle whose bonds types are not defined inpmb.df
- backbone_vector(
list
offloat
): 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, ...]}}
1516 def create_variable_with_units(self, variable): 1517 """ 1518 Returns a pint object with the value and units defined in `variable`. 1519 1520 Args: 1521 variable(`dict` or `str`): {'value': value, 'units': units} 1522 Returns: 1523 variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry. 1524 """ 1525 if isinstance(variable, dict): 1526 value=variable.pop('value') 1527 units=variable.pop('units') 1528 elif isinstance(variable, str): 1529 value = float(re.split(r'\s+', variable)[0]) 1530 units = re.split(r'\s+', variable)[1] 1531 variable_with_units=value*self.units(units) 1532 return variable_with_units
Returns a pint object with the value and units defined in variable
.
Arguments:
- variable(
dict
orstr
): {'value': value, 'units': units}
Returns:
variable_with_units(
obj
): variable with units using the pyMBE UnitRegistry.
1534 def define_AA_residues(self, sequence, model): 1535 """ 1536 Defines in `pmb.df` all the different residues in `sequence`. 1537 1538 Args: 1539 sequence(`lst`): Sequence of the peptide or protein. 1540 model(`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1541 1542 Returns: 1543 residue_list(`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1544 """ 1545 residue_list = [] 1546 for residue_name in sequence: 1547 if model == '1beadAA': 1548 central_bead = residue_name 1549 side_chains = [] 1550 elif model == '2beadAA': 1551 if residue_name in ['c','n', 'G']: 1552 central_bead = residue_name 1553 side_chains = [] 1554 else: 1555 central_bead = 'CA' 1556 side_chains = [residue_name] 1557 if residue_name not in residue_list: 1558 self.define_residue(name = 'AA-'+residue_name, 1559 central_bead = central_bead, 1560 side_chains = side_chains) 1561 residue_list.append('AA-'+residue_name) 1562 return residue_list
Defines in pmb.df
all the different residues in sequence
.
Arguments:
- sequence(
lst
): Sequence of the peptide or protein. - model(
string
): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported.
Returns:
residue_list(
list
ofstr
): List of thename
s of theresidue
s in the sequence of themolecule
.
1564 def define_bond(self, bond_type, bond_parameters, particle_pairs): 1565 1566 ''' 1567 Defines a pmb object of type `bond` in `pymbe.df`. 1568 1569 Args: 1570 bond_type(`str`): label to identify the potential to model the bond. 1571 bond_parameters(`dict`): parameters of the potential of the bond. 1572 particle_pairs(`lst`): list of the `names` of the `particles` to be bonded. 1573 1574 Note: 1575 Currently, only HARMONIC and FENE bonds are supported. 1576 1577 For a HARMONIC bond the dictionary must contain the following parameters: 1578 1579 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 1580 using the `pmb.units` UnitRegistry. 1581 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 1582 the `pmb.units` UnitRegistry. 1583 1584 For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and: 1585 1586 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 1587 units of length using the `pmb.units` UnitRegistry. Default 'None'. 1588 ''' 1589 1590 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1591 for particle_name1, particle_name2 in particle_pairs: 1592 1593 lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1, 1594 particle_name2 = particle_name2, 1595 combining_rule = 'Lorentz-Berthelot') 1596 1597 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1598 bond_type = bond_type, 1599 epsilon = lj_parameters["epsilon"], 1600 sigma = lj_parameters["sigma"], 1601 cutoff = lj_parameters["cutoff"], 1602 offset = lj_parameters["offset"],) 1603 index = len(self.df) 1604 for label in [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}']: 1605 self._check_if_multiple_pmb_types_for_name(name=label, pmb_type_to_be_defined="bond") 1606 name=f'{particle_name1}-{particle_name2}' 1607 self._check_if_multiple_pmb_types_for_name(name=name, pmb_type_to_be_defined="bond") 1608 self.df.at [index,'name']= name 1609 self.df.at [index,'bond_object'] = bond_object 1610 self.df.at [index,'l0'] = l0 1611 self.add_value_to_df(index=index, 1612 key=('pmb_type',''), 1613 new_value='bond') 1614 self.add_value_to_df(index=index, 1615 key=('parameters_of_the_potential',''), 1616 new_value=bond_object.get_params(), 1617 non_standard_value=True) 1618 self.df.fillna(pd.NA, inplace=True) 1619 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 thenames
of theparticles
to be bonded.
Note:
Currently, only HARMONIC and FENE bonds are supported.
For a HARMONIC bond the dictionary must contain the following parameters:
- k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 using the `pmb.units` UnitRegistry. - r_0 (`obj`) : Equilibrium bond length. It should have units of length using the `pmb.units` UnitRegistry.
For a FENE bond the dictionary must contain the same parameters as for a HARMONIC bond and:
- d_r_max (`obj`): Maximal stretching length for FENE. It should have units of length using the `pmb.units` UnitRegistry. Default 'None'.
1622 def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None, offset=None): 1623 """ 1624 Asigns `bond` in `pmb.df` as the default bond. 1625 The LJ parameters can be optionally provided to calculate the initial bond length 1626 1627 Args: 1628 bond_type(`str`): label to identify the potential to model the bond. 1629 bond_parameters(`dict`): parameters of the potential of the bond. 1630 sigma(`float`, optional): LJ sigma of the interaction between the particles. 1631 epsilon(`float`, optional): LJ epsilon for the interaction between the particles. 1632 cutoff(`float`, optional): cutoff-radius of the LJ interaction. 1633 offset(`float`, optional): offset of the LJ interaction. 1634 1635 Note: 1636 - Currently, only harmonic and FENE bonds are supported. 1637 """ 1638 1639 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1640 1641 if epsilon is None: 1642 epsilon=1*self.units('reduced_energy') 1643 if sigma is None: 1644 sigma=1*self.units('reduced_length') 1645 if cutoff is None: 1646 cutoff=2**(1.0/6.0)*self.units('reduced_length') 1647 if offset is None: 1648 offset=0*self.units('reduced_length') 1649 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1650 bond_type = bond_type, 1651 epsilon = epsilon, 1652 sigma = sigma, 1653 cutoff = cutoff, 1654 offset = offset) 1655 1656 self._check_if_multiple_pmb_types_for_name(name='default', 1657 pmb_type_to_be_defined='bond') 1658 1659 if len(self.df.index) != 0: 1660 index = max(self.df.index)+1 1661 else: 1662 index = 0 1663 self.df.at [index,'name'] = 'default' 1664 self.df.at [index,'bond_object'] = bond_object 1665 self.df.at [index,'l0'] = l0 1666 self.add_value_to_df(index = index, 1667 key = ('pmb_type',''), 1668 new_value = 'bond') 1669 self.add_value_to_df(index = index, 1670 key = ('parameters_of_the_potential',''), 1671 new_value = bond_object.get_params(), 1672 non_standard_value=True) 1673 self.df.fillna(pd.NA, inplace=True) 1674 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.
1676 def define_hydrogel(self, name, node_map, chain_map): 1677 """ 1678 Defines a pyMBE object of type `hydrogel` in `pymbe.df`. 1679 1680 Args: 1681 name(`str`): Unique label that identifies the `hydrogel`. 1682 node_map(`list of ict`): [{"particle_name": , "lattice_index": }, ... ] 1683 chain_map(`list of dict`): [{"node_start": , "node_end": , "residue_list": , ... ] 1684 """ 1685 node_indices = {tuple(entry['lattice_index']) for entry in node_map} 1686 diamond_indices = {tuple(row) for row in self.lattice_builder.lattice.indices} 1687 if node_indices != diamond_indices: 1688 raise ValueError(f"Incomplete hydrogel: A diamond lattice must contain exactly 8 lattice indices, {diamond_indices} ") 1689 1690 chain_map_connectivity = set() 1691 for entry in chain_map: 1692 start = self.lattice_builder.node_labels[entry['node_start']] 1693 end = self.lattice_builder.node_labels[entry['node_end']] 1694 chain_map_connectivity.add((start,end)) 1695 1696 if self.lattice_builder.lattice.connectivity != chain_map_connectivity: 1697 raise ValueError("Incomplete hydrogel: A diamond lattice must contain correct 16 lattice index pairs") 1698 1699 self._check_if_multiple_pmb_types_for_name(name=name, 1700 pmb_type_to_be_defined='hydrogel') 1701 1702 1703 index = len(self.df) 1704 self.df.at [index, "name"] = name 1705 self.df.at [index, "pmb_type"] = "hydrogel" 1706 self.add_value_to_df(index = index, 1707 key = ('node_map',''), 1708 new_value = node_map, 1709 non_standard_value=True) 1710 self.add_value_to_df(index = index, 1711 key = ('chain_map',''), 1712 new_value = chain_map, 1713 non_standard_value=True) 1714 for chain_label in chain_map: 1715 node_start = chain_label["node_start"] 1716 node_end = chain_label["node_end"] 1717 residue_list = chain_label['residue_list'] 1718 # Molecule name 1719 molecule_name = "chain_"+node_start+"_"+node_end 1720 self.define_molecule(name=molecule_name, residue_list=residue_list) 1721 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": , ... ]
1723 def define_molecule(self, name, residue_list): 1724 """ 1725 Defines a pyMBE object of type `molecule` in `pymbe.df`. 1726 1727 Args: 1728 name(`str`): Unique label that identifies the `molecule`. 1729 residue_list(`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1730 """ 1731 self._check_if_multiple_pmb_types_for_name(name=name, 1732 pmb_type_to_be_defined='molecule') 1733 1734 index = len(self.df) 1735 self.df.at [index,'name'] = name 1736 self.df.at [index,'pmb_type'] = 'molecule' 1737 self.df.at [index,('residue_list','')] = residue_list 1738 self.df.fillna(pd.NA, inplace=True) 1739 return
Defines a pyMBE object of type molecule
in pymbe.df
.
Arguments:
- name(
str
): Unique label that identifies themolecule
. - residue_list(
list
ofstr
): List of thename
s of theresidue
s in the sequence of themolecule
.
1741 def define_particle_entry_in_df(self,name): 1742 """ 1743 Defines a particle entry in pmb.df. 1744 1745 Args: 1746 name(`str`): Unique label that identifies this particle type. 1747 1748 Returns: 1749 index(`int`): Index of the particle in pmb.df 1750 """ 1751 1752 if self._check_if_name_is_defined_in_df(name=name): 1753 index = self.df[self.df['name']==name].index[0] 1754 else: 1755 index = len(self.df) 1756 self.df.at [index, 'name'] = name 1757 self.df.at [index,'pmb_type'] = 'particle' 1758 self.df.fillna(pd.NA, inplace=True) 1759 return index
Defines a particle entry in pmb.df.
Arguments:
- name(
str
): Unique label that identifies this particle type.
Returns:
index(
int
): Index of the particle in pmb.df
1761 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): 1762 """ 1763 Defines the properties of a particle object. 1764 1765 Args: 1766 name(`str`): Unique label that identifies this particle type. 1767 z(`int`, optional): Permanent charge number of this particle type. Defaults to 0. 1768 acidity(`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to pd.NA. 1769 pka(`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to pd.NA. 1770 sigma(`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1771 cutoff(`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1772 offset(`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. 1773 epsilon(`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA. 1774 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1775 1776 Note: 1777 - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units. 1778 - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units. 1779 - `cutoff` defaults to `2**(1./6.) reduced_length`. 1780 - `offset` defaults to 0. 1781 - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions. 1782 - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`. 1783 """ 1784 index=self.define_particle_entry_in_df(name=name) 1785 self._check_if_multiple_pmb_types_for_name(name=name, 1786 pmb_type_to_be_defined='particle') 1787 1788 # If `cutoff` and `offset` are not defined, default them to the following values 1789 if pd.isna(cutoff): 1790 cutoff=self.units.Quantity(2**(1./6.), "reduced_length") 1791 if pd.isna(offset): 1792 offset=self.units.Quantity(0, "reduced_length") 1793 1794 # Define LJ parameters 1795 parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"}, 1796 "cutoff":{"value": cutoff, "dimensionality": "[length]"}, 1797 "offset":{"value": offset, "dimensionality": "[length]"}, 1798 "epsilon":{"value": epsilon, "dimensionality": "[energy]"},} 1799 1800 for parameter_key in parameters_with_dimensionality.keys(): 1801 if not pd.isna(parameters_with_dimensionality[parameter_key]["value"]): 1802 self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 1803 expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"]) 1804 self.add_value_to_df(key=(parameter_key,''), 1805 index=index, 1806 new_value=parameters_with_dimensionality[parameter_key]["value"], 1807 overwrite=overwrite) 1808 1809 # Define particle acid/base properties 1810 self.set_particle_acidity(name=name, 1811 acidity=acidity, 1812 default_charge_number=z, 1813 pka=pka, 1814 overwrite=overwrite) 1815 self.df.fillna(pd.NA, inplace=True) 1816 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 isacidic
orbasic
, used to setup constant pH simulations. Defaults to pd.NA. - pka(
float
, optional): Ifparticle
is an acid or a base, it defines its pka-value. Defaults to pd.NA. - sigma(
pint.Quantity
, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. - cutoff(
pint.Quantity
, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. - offset(
pint.Quantity
, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to pd.NA. - epsilon(
pint.Quantity
, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to pd.NA. - overwrite(
bool
, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
Note:
sigma
,cutoff
andoffset
must have a dimensitonality of[length]
and should be defined using pmb.units.epsilon
must have a dimensitonality of[energy]
and should be defined using pmb.units.cutoff
defaults to2**(1./6.) reduced_length
.offset
defaults to 0.- The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions.
- For more information on
sigma
,epsilon
,cutoff
andoffset
checkpmb.setup_lj_interactions()
.
1818 def define_particles(self, parameters, overwrite=False): 1819 ''' 1820 Defines a particle object in pyMBE for each particle name in `particle_names` 1821 1822 Args: 1823 parameters(`dict`): dictionary with the particle parameters. 1824 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1825 1826 Note: 1827 - parameters = {"particle_name1: {"sigma": sigma_value, "epsilon": epsilon_value, ...}, particle_name2: {...},} 1828 ''' 1829 if not parameters: 1830 return 0 1831 for particle_name in parameters.keys(): 1832 parameters[particle_name]["overwrite"]=overwrite 1833 self.define_particle(**parameters[particle_name]) 1834 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: {...},}
1836 def define_peptide(self, name, sequence, model): 1837 """ 1838 Defines a pyMBE object of type `peptide` in the `pymbe.df`. 1839 1840 Args: 1841 name (`str`): Unique label that identifies the `peptide`. 1842 sequence (`string`): Sequence of the `peptide`. 1843 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1844 """ 1845 self._check_if_multiple_pmb_types_for_name(name = name, 1846 pmb_type_to_be_defined='peptide') 1847 1848 valid_keys = ['1beadAA','2beadAA'] 1849 if model not in valid_keys: 1850 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1851 clean_sequence = self.protein_sequence_parser(sequence=sequence) 1852 residue_list = self.define_AA_residues(sequence=clean_sequence, 1853 model=model) 1854 self.define_molecule(name = name, residue_list=residue_list) 1855 index = self.df.loc[self.df['name'] == name].index.item() 1856 self.df.at [index,'model'] = model 1857 self.df.at [index,('sequence','')] = clean_sequence 1858 self.df.at [index,'pmb_type'] = "peptide" 1859 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.
1862 def define_protein(self, name,model, topology_dict, lj_setup_mode="wca", overwrite=False): 1863 """ 1864 Defines a globular protein pyMBE object in `pymbe.df`. 1865 1866 Args: 1867 name (`str`): Unique label that identifies the protein. 1868 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1869 topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'radius': radius_value} 1870 lj_setup_mode(`str`): Key for the setup for the LJ potential. Defaults to "wca". 1871 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 1872 1873 Note: 1874 - Currently, only `lj_setup_mode="wca"` is supported. This corresponds to setting up the WCA potential. 1875 """ 1876 1877 # Sanity checks 1878 self._check_if_multiple_pmb_types_for_name(name = name, 1879 pmb_type_to_be_defined='protein') 1880 valid_model_keys = ['1beadAA','2beadAA'] 1881 valid_lj_setups = ["wca"] 1882 if model not in valid_model_keys: 1883 raise ValueError('Invalid key for the protein model, supported models are {valid_model_keys}') 1884 if lj_setup_mode not in valid_lj_setups: 1885 raise ValueError('Invalid key for the lj setup, supported setup modes are {valid_lj_setups}') 1886 if lj_setup_mode == "wca": 1887 sigma = 1*self.units.Quantity("reduced_length") 1888 epsilon = 1*self.units.Quantity("reduced_energy") 1889 part_dict={} 1890 sequence=[] 1891 metal_ions_charge_number_map=self.get_metal_ions_charge_number_map() 1892 for particle in topology_dict.keys(): 1893 particle_name = re.split(r'\d+', particle)[0] 1894 if particle_name not in part_dict.keys(): 1895 if lj_setup_mode == "wca": 1896 part_dict[particle_name]={"sigma": sigma, 1897 "offset": topology_dict[particle]['radius']*2-sigma, 1898 "epsilon": epsilon, 1899 "name": particle_name} 1900 if self.check_if_metal_ion(key=particle_name): 1901 z=metal_ions_charge_number_map[particle_name] 1902 else: 1903 z=0 1904 part_dict[particle_name]["z"]=z 1905 1906 if self.check_aminoacid_key(key=particle_name): 1907 sequence.append(particle_name) 1908 1909 self.define_particles(parameters=part_dict, 1910 overwrite=overwrite) 1911 residue_list = self.define_AA_residues(sequence=sequence, 1912 model=model) 1913 index = len(self.df) 1914 self.df.at [index,'name'] = name 1915 self.df.at [index,'pmb_type'] = 'protein' 1916 self.df.at [index,'model'] = model 1917 self.df.at [index,('sequence','')] = sequence 1918 self.df.at [index,('residue_list','')] = residue_list 1919 self.df.fillna(pd.NA, inplace=True) 1920 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.
1922 def define_residue(self, name, central_bead, side_chains): 1923 """ 1924 Defines a pyMBE object of type `residue` in `pymbe.df`. 1925 1926 Args: 1927 name(`str`): Unique label that identifies the `residue`. 1928 central_bead(`str`): `name` of the `particle` to be placed as central_bead of the `residue`. 1929 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. 1930 """ 1931 self._check_if_multiple_pmb_types_for_name(name=name, 1932 pmb_type_to_be_defined='residue') 1933 1934 index = len(self.df) 1935 self.df.at [index, 'name'] = name 1936 self.df.at [index,'pmb_type'] = 'residue' 1937 self.df.at [index,'central_bead'] = central_bead 1938 self.df.at [index,('side_chains','')] = side_chains 1939 self.df.fillna(pd.NA, inplace=True) 1940 return
Defines a pyMBE object of type residue
in pymbe.df
.
Arguments:
- name(
str
): Unique label that identifies theresidue
. - central_bead(
str
):name
of theparticle
to be placed as central_bead of theresidue
. - side_chains(
list
ofstr
): List ofname
s of the pmb_objects to be placed as side_chains of theresidue
. Currently, only pmb_objects of typeparticle
s orresidue
s are supported.
1942 def delete_entries_in_df(self, entry_name): 1943 """ 1944 Deletes entries with name `entry_name` from the DataFrame if it exists. 1945 1946 Args: 1947 entry_name (`str`): The name of the entry in the dataframe to delete. 1948 1949 """ 1950 if entry_name in self.df["name"].values: 1951 self.df = self.df[self.df["name"] != entry_name].reset_index(drop=True)
Deletes entries with name entry_name
from the DataFrame if it exists.
Arguments:
- entry_name (
str
): The name of the entry in the dataframe to delete.
1953 def destroy_pmb_object_in_system(self, name, espresso_system): 1954 """ 1955 Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 1956 1957 Args: 1958 name(`str`): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`. 1959 espresso_system(`espressomd.system.System`): Instance of a system class from espressomd library. 1960 1961 Note: 1962 - If `name` is a object_type=`particle`, only the matching particles that are not part of bigger objects (e.g. `residue`, `molecule`) will be destroyed. To destroy particles in such objects, destroy the bigger object instead. 1963 """ 1964 pmb_type = self._check_supported_molecule(molecule_name=name, 1965 valid_pmb_types= ['particle','residue','molecule',"peptide"]) 1966 1967 if pmb_type == 'particle': 1968 particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1969 & (self.df['molecule_id'].isna())] 1970 particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1971 & (self.df['molecule_id'].isna())].particle_id.tolist() 1972 for particle_id in particle_ids_list: 1973 espresso_system.part.by_id(particle_id).remove() 1974 self.df = self.df.drop(index=particles_index) 1975 if pmb_type == 'residue': 1976 residues_id = self.df.loc[self.df['name']== name].residue_id.to_list() 1977 for residue_id in residues_id: 1978 molecule_name = self.df.loc[(self.df['residue_id']==molecule_id) & (self.df['pmb_type']=="residue")].name.values[0] 1979 particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"] 1980 self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index) 1981 for particle_id in particle_ids_list: 1982 espresso_system.part.by_id(particle_id).remove() 1983 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1984 if pmb_type in ['molecule',"peptide"]: 1985 molecules_id = self.df.loc[self.df['name']== name].molecule_id.dropna().to_list() 1986 for molecule_id in molecules_id: 1987 molecule_name = self.df.loc[(self.df['molecule_id']==molecule_id) & (self.df['pmb_type']==pmb_type)].name.values[0] 1988 particle_ids_list = self.get_particle_id_map(object_name=molecule_name)["all"] 1989 self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index) 1990 for particle_id in particle_ids_list: 1991 espresso_system.part.by_id(particle_id).remove() 1992 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1993 1994 self.df.reset_index(drop=True,inplace=True) 1995 1996 return
Destroys all particles associated with name
in espresso_system
amd removes the destroyed pmb_objects from pmb.df
Arguments:
- name(
str
): Label of the pmb object to be destroyed. The pmb object must be defined inpymbe.df
. - espresso_system(
espressomd.system.System
): Instance of a system class from espressomd library.
Note:
- If
name
is a object_type=particle
, only the matching particles that are not part of bigger objects (e.g.residue
,molecule
) will be destroyed. To destroy particles in such objects, destroy the bigger object instead.
1998 def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200): 1999 """ 2000 Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration. 2001 To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an 2002 ion pair and the concentrations of the various species is solved numerically using a self-consistent approach. 2003 More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831). 2004 This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237). 2005 2006 Args: 2007 pH_res('float'): pH-value in the reservoir. 2008 c_salt_res('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2009 activity_coefficient_monovalent_pair('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2010 2011 Returns: 2012 cH_res('pint.Quantity'): Concentration of H+ ions. 2013 cOH_res('pint.Quantity'): Concentration of OH- ions. 2014 cNa_res('pint.Quantity'): Concentration of Na+ ions. 2015 cCl_res('pint.Quantity'): Concentration of Cl- ions. 2016 """ 2017 2018 self_consistent_run = 0 2019 cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+ 2020 2021 def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res): 2022 #Calculate and initial guess for the concentrations of various species based on ideal gas estimate 2023 cOH_res = self.Kw / cH_res 2024 cNa_res = None 2025 cCl_res = None 2026 if cOH_res>=cH_res: 2027 #adjust the concentration of sodium if there is excess OH- in the reservoir: 2028 cNa_res = c_salt_res + (cOH_res-cH_res) 2029 cCl_res = c_salt_res 2030 else: 2031 # adjust the concentration of chloride if there is excess H+ in the reservoir 2032 cCl_res = c_salt_res + (cH_res-cOH_res) 2033 cNa_res = c_salt_res 2034 2035 def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res): 2036 nonlocal max_number_sc_runs, self_consistent_run 2037 if self_consistent_run<max_number_sc_runs: 2038 self_consistent_run+=1 2039 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2040 cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res)) 2041 if cOH_res>=cH_res: 2042 #adjust the concentration of sodium if there is excess OH- in the reservoir: 2043 cNa_res = c_salt_res + (cOH_res-cH_res) 2044 cCl_res = c_salt_res 2045 else: 2046 # adjust the concentration of chloride if there is excess H+ in the reservoir 2047 cCl_res = c_salt_res + (cH_res-cOH_res) 2048 cNa_res = c_salt_res 2049 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 2050 else: 2051 return cH_res, cOH_res, cNa_res, cCl_res 2052 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 2053 2054 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 2055 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2056 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 2057 2058 while abs(determined_pH-pH_res)>1e-6: 2059 if determined_pH > pH_res: 2060 cH_res *= 1.005 2061 else: 2062 cH_res /= 1.003 2063 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 2064 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2065 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 2066 self_consistent_run=0 2067 2068 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.
2070 def enable_motion_of_rigid_object(self, name, espresso_system): 2071 ''' 2072 Enables the motion of the rigid object `name` in the `espresso_system`. 2073 2074 Args: 2075 name(`str`): Label of the object. 2076 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 2077 2078 Note: 2079 - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]. 2080 ''' 2081 logging.info('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]') 2082 self._check_supported_molecule(molecule_name=name, 2083 valid_pmb_types= ['protein']) 2084 molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list() 2085 for molecule_id in molecule_ids_list: 2086 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 2087 center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system) 2088 rigid_object_center = espresso_system.part.add(pos=center_of_mass, 2089 rotation=[True,True,True], 2090 type=self.propose_unused_type()) 2091 2092 rigid_object_center.mass = len(particle_ids_list) 2093 momI = 0 2094 for pid in particle_ids_list: 2095 momI += np.power(np.linalg.norm(center_of_mass - espresso_system.part.by_id(pid).pos), 2) 2096 2097 rigid_object_center.rinertia = np.ones(3) * momI 2098 2099 for particle_id in particle_ids_list: 2100 pid = espresso_system.part.by_id(particle_id) 2101 pid.vs_auto_relate_to(rigid_object_center.id) 2102 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"].
2104 def filter_df(self, pmb_type): 2105 """ 2106 Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns. 2107 2108 Args: 2109 pmb_type(`str`): pmb_object_type to filter in `pmb.df`. 2110 2111 Returns: 2112 pmb_type_df(`Pandas.Dataframe`): filtered `pmb.df`. 2113 """ 2114 pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type] 2115 pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1) 2116 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
.
2118 def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False): 2119 """ 2120 Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2121 2122 Args: 2123 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 2124 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 2125 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'. 2126 2127 Returns: 2128 bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` if a matching bond exists 2129 2130 Note: 2131 - If `use_default_bond`=`True`, it returns "default" if no key is found. 2132 """ 2133 bond_keys = [f'{particle_name1}-{particle_name2}', f'{particle_name2}-{particle_name1}'] 2134 bond_defined=False 2135 for bond_key in bond_keys: 2136 if bond_key in self.df["name"].values: 2137 bond_defined=True 2138 correct_key=bond_key 2139 break 2140 if bond_defined: 2141 return correct_key 2142 elif use_default_bond: 2143 return 'default' 2144 else: 2145 return
Searches for the name
of the bond between particle_name1
and particle_name2
in pymbe.df
and returns it.
Arguments:
- particle_name1(
str
): label of the type of the first particle type of the bonded particles. - particle_name2(
str
): label of the type of the second particle type of the bonded particles. - use_default_bond(
bool
, optional): If it is activated, the "default" bond is returned if no bond is found betweenparticle_name1
andparticle_name2
. Defaults to 'False'.
Returns:
bond_key (str):
name
of the bond betweenparticle_name1
andparticle_name2
if a matching bond exists
Note:
- If
use_default_bond
=True
, it returns "default" if no key is found.
2147 def find_value_from_es_type(self, es_type, column_name): 2148 """ 2149 Finds a value in `pmb.df` for a `column_name` and `es_type` pair. 2150 2151 Args: 2152 es_type(`int`): value of the espresso type 2153 column_name(`str`): name of the column in `pymbe.df` 2154 2155 Returns: 2156 Value in `pymbe.df` matching `column_name` and `es_type` 2157 """ 2158 idx = pd.IndexSlice 2159 for state in ['state_one', 'state_two']: 2160 index = self.df.loc[self.df[(state, 'es_type')] == es_type].index 2161 if len(index) > 0: 2162 if column_name == 'label': 2163 label = self.df.loc[idx[index[0]], idx[(state,column_name)]] 2164 return label 2165 else: 2166 column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]] 2167 return column_name_value 2168 return None
Finds a value in pmb.df
for a column_name
and es_type
pair.
Arguments:
- es_type(
int
): value of the espresso type - column_name(
str
): name of the column inpymbe.df
Returns:
Value in
pymbe.df
matchingcolumn_name
andes_type
2174 def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples): 2175 """ 2176 Generates coordinates outside a sphere centered at `center`. 2177 2178 Args: 2179 center(`lst`): Coordinates of the center of the sphere. 2180 radius(`float`): Radius of the sphere. 2181 max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates. 2182 n_samples(`int`): Number of sample points. 2183 2184 Returns: 2185 coord_list(`lst`): Coordinates of the sample points. 2186 """ 2187 if not radius > 0: 2188 raise ValueError (f'The value of {radius} must be a positive value') 2189 if not radius < max_dist: 2190 raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))') 2191 coord_list = [] 2192 counter = 0 2193 while counter<n_samples: 2194 coord = self.generate_random_points_in_a_sphere(center=center, 2195 radius=max_dist, 2196 n_samples=1)[0] 2197 if np.linalg.norm(coord-np.asarray(center))>=radius: 2198 coord_list.append (coord) 2199 counter += 1 2200 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.
2202 def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False): 2203 """ 2204 Uniformly samples points from a hypersphere. If on_surface is set to True, the points are 2205 uniformly sampled from the surface of the hypersphere. 2206 2207 Args: 2208 center(`lst`): Array with the coordinates of the center of the spheres. 2209 radius(`float`): Radius of the sphere. 2210 n_samples(`int`): Number of sample points to generate inside the sphere. 2211 on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere. 2212 2213 Returns: 2214 samples(`list`): Coordinates of the sample points inside the hypersphere. 2215 """ 2216 # initial values 2217 center=np.array(center) 2218 d = center.shape[0] 2219 # sample n_samples points in d dimensions from a standard normal distribution 2220 samples = self.rng.normal(size=(n_samples, d)) 2221 # make the samples lie on the surface of the unit hypersphere 2222 normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis] 2223 samples /= normalize_radii 2224 if not on_surface: 2225 # make the samples lie inside the hypersphere with the correct density 2226 uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis] 2227 new_radii = np.power(uniform_points, 1/d) 2228 samples *= new_radii 2229 # scale the points to have the correct radius and center 2230 samples = samples * radius + center 2231 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.
2233 def generate_trial_perpendicular_vector(self,vector,magnitude): 2234 """ 2235 Generates an orthogonal vector to the input `vector`. 2236 2237 Args: 2238 vector(`lst`): arbitrary vector. 2239 magnitude(`float`): magnitude of the orthogonal vector. 2240 2241 Returns: 2242 (`lst`): Orthogonal vector with the same magnitude as the input vector. 2243 """ 2244 np_vec = np.array(vector) 2245 if np.all(np_vec == 0): 2246 raise ValueError('Zero vector') 2247 np_vec /= np.linalg.norm(np_vec) 2248 # Generate a random vector 2249 random_vector = self.generate_random_points_in_a_sphere(radius=1, 2250 center=[0,0,0], 2251 n_samples=1, 2252 on_surface=True)[0] 2253 # Project the random vector onto the input vector and subtract the projection 2254 projection = np.dot(random_vector, np_vec) * np_vec 2255 perpendicular_vector = random_vector - projection 2256 # Normalize the perpendicular vector to have the same magnitude as the input vector 2257 perpendicular_vector /= np.linalg.norm(perpendicular_vector) 2258 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.
2260 def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2261 """ 2262 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length. 2263 If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead. 2264 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. 2265 2266 Args: 2267 particle_name1(str): label of the type of the first particle type of the bonded particles. 2268 particle_name2(str): label of the type of the second particle type of the bonded particles. 2269 hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2270 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. 2271 2272 Returns: 2273 l0(`pint.Quantity`): bond length 2274 2275 Note: 2276 - 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`. 2277 - If `hard_check`=`True` stops the code when no bond is found. 2278 """ 2279 bond_key = self.find_bond_key(particle_name1=particle_name1, 2280 particle_name2=particle_name2, 2281 use_default_bond=use_default_bond) 2282 if bond_key: 2283 return self.df[self.df['name']==bond_key].l0.values[0] 2284 else: 2285 msg = f"Bond not defined between particles {particle_name1} and {particle_name2}" 2286 if hard_check: 2287 raise ValueError(msg) 2288 else: 2289 logging.warning(msg) 2290 return
Searches for bonds between the particle types given by particle_name1
and particle_name2
in pymbe.df
and returns the initial bond length.
If use_default_bond
is activated and a "default" bond is defined, returns the length of that default bond instead.
If no bond is found, it prints a message and it does not return anything. If hard_check
is activated, the code stops if no bond is found.
Arguments:
- particle_name1(str): label of the type of the first particle type of the bonded particles.
- particle_name2(str): label of the type of the second particle type of the bonded particles.
- hard_check(bool, optional): If it is activated, the code stops if no bond is found. Defaults to False.
- use_default_bond(bool, optional): If it is activated, the "default" bond is returned if no bond is found between
particle_name1
andparticle_name2
. Defaults to False.
Returns:
l0(
pint.Quantity
): bond length
Note:
- If
use_default_bond
=True and no bond is defined betweenparticle_name1
andparticle_name2
, it returns the default bond defined inpmb.df
.- If
hard_check
=True
stops the code when no bond is found.
2292 def get_charge_number_map(self): 2293 ''' 2294 Gets the charge number of each `espresso_type` in `pymbe.df`. 2295 2296 Returns: 2297 charge_number_map(`dict`): {espresso_type: z}. 2298 ''' 2299 if self.df.state_one['es_type'].isnull().values.any(): 2300 df_state_one = self.df.state_one.dropna() 2301 df_state_two = self.df.state_two.dropna() 2302 else: 2303 df_state_one = self.df.state_one 2304 if self.df.state_two['es_type'].isnull().values.any(): 2305 df_state_two = self.df.state_two.dropna() 2306 else: 2307 df_state_two = self.df.state_two 2308 state_one = pd.Series (df_state_one.z.values,index=df_state_one.es_type.values) 2309 state_two = pd.Series (df_state_two.z.values,index=df_state_two.es_type.values) 2310 charge_number_map = pd.concat([state_one,state_two],axis=0).to_dict() 2311 return charge_number_map
Gets the charge number of each espresso_type
in pymbe.df
.
Returns:
charge_number_map(
dict
): {espresso_type: z}.
2313 def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'): 2314 """ 2315 Returns the Lennard-Jones parameters for the interaction between the particle types given by 2316 `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule. 2317 2318 Args: 2319 particle_name1 (str): label of the type of the first particle type 2320 particle_name2 (str): label of the type of the second particle type 2321 combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'. 2322 2323 Returns: 2324 {"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value} 2325 2326 Note: 2327 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 2328 - 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. 2329 """ 2330 supported_combining_rules=["Lorentz-Berthelot"] 2331 lj_parameters_keys=["sigma","epsilon","offset","cutoff"] 2332 if combining_rule not in supported_combining_rules: 2333 raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}") 2334 lj_parameters={} 2335 for key in lj_parameters_keys: 2336 lj_parameters[key]=[] 2337 # Search the LJ parameters of the type pair 2338 for name in [particle_name1,particle_name2]: 2339 for key in lj_parameters_keys: 2340 lj_parameters[key].append(self.df[self.df.name == name][key].values[0]) 2341 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 2342 if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]): 2343 return {} 2344 # Apply combining rule 2345 if combining_rule == 'Lorentz-Berthelot': 2346 lj_parameters["sigma"]=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2 2347 lj_parameters["cutoff"]=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2 2348 lj_parameters["offset"]=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2 2349 lj_parameters["epsilon"]=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1]) 2350 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 calculatesigma
andepsilon
for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'.
Returns:
{"epsilon": epsilon_value, "sigma": sigma_value, "offset": offset_value, "cutoff": cutoff_value}
Note:
- Currently, the only
combining_rule
supported is Lorentz-Berthelot.- If the sigma value of
particle_name1
orparticle_name2
is 0, the function will return an empty dictionary. No LJ interactions are set up for particles with sigma = 0.
2352 def get_metal_ions_charge_number_map(self): 2353 """ 2354 Gets a map with the charge numbers of all the metal ions supported. 2355 2356 Returns: 2357 metal_charge_number_map(dict): Has the structure {"metal_name": metal_charge_number} 2358 2359 """ 2360 metal_charge_number_map = {"Ca": 2} 2361 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}
2363 def get_particle_id_map(self, object_name): 2364 ''' 2365 Gets all the ids associated with the object with name `object_name` in `pmb.df` 2366 2367 Args: 2368 object_name(`str`): name of the object 2369 2370 Returns: 2371 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]}, } 2372 ''' 2373 object_type=self._check_supported_molecule(molecule_name=object_name, 2374 valid_pmb_types= ['particle','residue','molecule',"peptide","protein"]) 2375 id_list = [] 2376 mol_map = {} 2377 res_map = {} 2378 def do_res_map(res_ids): 2379 for res_id in res_ids: 2380 res_list=self.df.loc[(self.df['residue_id']== res_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist() 2381 res_map[res_id]=res_list 2382 return res_map 2383 if object_type in ['molecule', 'protein', 'peptide']: 2384 mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist() 2385 for mol_id in mol_ids: 2386 res_ids = set(self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle") ].residue_id.dropna().tolist()) 2387 res_map=do_res_map(res_ids=res_ids) 2388 mol_list=self.df.loc[(self.df['molecule_id']== mol_id) & (self.df['pmb_type']== "particle")].particle_id.dropna().tolist() 2389 id_list+=mol_list 2390 mol_map[mol_id]=mol_list 2391 elif object_type == 'residue': 2392 res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist() 2393 res_map=do_res_map(res_ids=res_ids) 2394 id_list=[] 2395 for res_id_list in res_map.values(): 2396 id_list+=res_id_list 2397 elif object_type == 'particle': 2398 id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist() 2399 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]}, }
2401 def get_pka_set(self): 2402 ''' 2403 Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df` 2404 2405 Returns: 2406 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 2407 ''' 2408 titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna() 2409 pka_set = {} 2410 for index in titratables_AA_df.name.keys(): 2411 name = titratables_AA_df.name[index] 2412 pka_value = titratables_AA_df.pka[index] 2413 acidity = titratables_AA_df.acidity[index] 2414 pka_set[name] = {'pka_value':pka_value,'acidity':acidity} 2415 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}}
2417 def get_radius_map(self, dimensionless=True): 2418 ''' 2419 Gets the effective radius of each `espresso type` in `pmb.df`. 2420 2421 Args: 2422 dimensionless('bool'): controlls if the returned radii have a dimension. Defaults to False. 2423 2424 Returns: 2425 radius_map(`dict`): {espresso_type: radius}. 2426 2427 Note: 2428 The radius corresponds to (sigma+offset)/2 2429 ''' 2430 df_state_one = self.df[[('sigma',''),('offset',''),('state_one','es_type')]].dropna().drop_duplicates() 2431 df_state_two = self.df[[('sigma',''),('offset',''),('state_two','es_type')]].dropna().drop_duplicates() 2432 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) 2433 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) 2434 radius_map = pd.concat([state_one,state_two],axis=0).to_dict() 2435 if dimensionless: 2436 for key in radius_map: 2437 radius_map[key] = radius_map[key].magnitude 2438 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
2440 def get_reduced_units(self): 2441 """ 2442 Returns the current set of reduced units defined in pyMBE. 2443 2444 Returns: 2445 reduced_units_text(`str`): text with information about the current set of reduced units. 2446 2447 """ 2448 unit_length=self.units.Quantity(1,'reduced_length') 2449 unit_energy=self.units.Quantity(1,'reduced_energy') 2450 unit_charge=self.units.Quantity(1,'reduced_charge') 2451 reduced_units_text = "\n".join(["Current set of reduced units:", 2452 f"{unit_length.to('nm'):.5g} = {unit_length}", 2453 f"{unit_energy.to('J'):.5g} = {unit_energy}", 2454 f"{unit_charge.to('C'):.5g} = {unit_charge}", 2455 f"Temperature: {(self.kT/self.kB).to('K'):.5g}" 2456 ]) 2457 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.
2459 def get_type_map(self): 2460 """ 2461 Gets all different espresso types assigned to particles in `pmb.df`. 2462 2463 Returns: 2464 type_map(`dict`): {"name": espresso_type}. 2465 """ 2466 df_state_one = self.df.state_one.dropna(how='all') 2467 df_state_two = self.df.state_two.dropna(how='all') 2468 state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label) 2469 state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label) 2470 type_map = pd.concat([state_one,state_two],axis=0).to_dict() 2471 return type_map
Gets all different espresso types assigned to particles in pmb.df
.
Returns:
type_map(
dict
): {"name": espresso_type}.
2473 def initialize_lattice_builder(self, diamond_lattice): 2474 """ 2475 Initialize the lattice builder with the DiamondLattice object. 2476 2477 Args: 2478 diamond_lattice(`DiamondLattice`): DiamondLattice object from the `lib/lattice` module to be used in the LatticeBuilder. 2479 """ 2480 from .lib.lattice import LatticeBuilder, DiamondLattice 2481 if not isinstance(diamond_lattice, DiamondLattice): 2482 raise TypeError("Currently only DiamondLattice objects are supported.") 2483 self.lattice_builder = LatticeBuilder(lattice=diamond_lattice) 2484 logging.info(f"LatticeBuilder initialized with mpc={diamond_lattice.mpc} and box_l={diamond_lattice.box_l}") 2485 return self.lattice_builder
Initialize the lattice builder with the DiamondLattice object.
Arguments:
- diamond_lattice(
DiamondLattice
): DiamondLattice object from thelib/lattice
module to be used in the LatticeBuilder.
2488 def load_interaction_parameters(self, filename, overwrite=False): 2489 """ 2490 Loads the interaction parameters stored in `filename` into `pmb.df` 2491 2492 Args: 2493 filename(`str`): name of the file to be read 2494 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 2495 """ 2496 without_units = ['z','es_type'] 2497 with_units = ['sigma','epsilon','offset','cutoff'] 2498 2499 with open(filename, 'r') as f: 2500 interaction_data = json.load(f) 2501 interaction_parameter_set = interaction_data["data"] 2502 2503 for key in interaction_parameter_set: 2504 param_dict=interaction_parameter_set[key] 2505 object_type=param_dict.pop('object_type') 2506 if object_type == 'particle': 2507 not_required_attributes={} 2508 for not_required_key in without_units+with_units: 2509 if not_required_key in param_dict.keys(): 2510 if not_required_key in with_units: 2511 not_required_attributes[not_required_key]=self.create_variable_with_units(variable=param_dict.pop(not_required_key)) 2512 elif not_required_key in without_units: 2513 not_required_attributes[not_required_key]=param_dict.pop(not_required_key) 2514 else: 2515 not_required_attributes[not_required_key]=pd.NA 2516 self.define_particle(name=param_dict.pop('name'), 2517 z=not_required_attributes.pop('z'), 2518 sigma=not_required_attributes.pop('sigma'), 2519 offset=not_required_attributes.pop('offset'), 2520 cutoff=not_required_attributes.pop('cutoff'), 2521 epsilon=not_required_attributes.pop('epsilon'), 2522 overwrite=overwrite) 2523 elif object_type == 'residue': 2524 self.define_residue(**param_dict) 2525 elif object_type == 'molecule': 2526 self.define_molecule(**param_dict) 2527 elif object_type == 'peptide': 2528 self.define_peptide(**param_dict) 2529 elif object_type == 'bond': 2530 particle_pairs = param_dict.pop('particle_pairs') 2531 bond_parameters = param_dict.pop('bond_parameters') 2532 bond_type = param_dict.pop('bond_type') 2533 if bond_type == 'harmonic': 2534 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2535 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2536 bond = {'r_0' : r_0, 2537 'k' : k, 2538 } 2539 2540 elif bond_type == 'FENE': 2541 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2542 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2543 d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max')) 2544 bond = {'r_0' : r_0, 2545 'k' : k, 2546 'd_r_max': d_r_max, 2547 } 2548 else: 2549 raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds") 2550 self.define_bond(bond_type=bond_type, 2551 bond_parameters=bond, 2552 particle_pairs=particle_pairs) 2553 else: 2554 raise ValueError(object_type+' is not a known pmb object type') 2555 2556 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.
2558 def load_pka_set(self, filename, overwrite=True): 2559 """ 2560 Loads the pka_set stored in `filename` into `pmb.df`. 2561 2562 Args: 2563 filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}. 2564 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to True. 2565 """ 2566 with open(filename, 'r') as f: 2567 pka_data = json.load(f) 2568 pka_set = pka_data["data"] 2569 2570 self.check_pka_set(pka_set=pka_set) 2571 2572 for key in pka_set: 2573 acidity = pka_set[key]['acidity'] 2574 pka_value = pka_set[key]['pka_value'] 2575 self.set_particle_acidity(name=key, 2576 acidity=acidity, 2577 pka=pka_value, 2578 overwrite=overwrite) 2579 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.
2582 def propose_unused_type(self): 2583 """ 2584 Searches in `pmb.df` all the different particle types defined and returns a new one. 2585 2586 Returns: 2587 unused_type(`int`): unused particle type 2588 """ 2589 type_map = self.get_type_map() 2590 if not type_map: 2591 unused_type = 0 2592 else: 2593 valid_values = [v for v in type_map.values() if pd.notna(v)] # Filter out pd.NA values 2594 unused_type = max(valid_values) + 1 if valid_values else 0 # Ensure max() doesn't fail if all values are NA 2595 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
2597 def protein_sequence_parser(self, sequence): 2598 ''' 2599 Parses `sequence` to the one letter code for amino acids. 2600 2601 Args: 2602 sequence(`str` or `lst`): Sequence of the amino acid. 2603 2604 Returns: 2605 clean_sequence(`lst`): `sequence` using the one letter code. 2606 2607 Note: 2608 - Accepted formats for `sequence` are: 2609 - `lst` with one letter or three letter code of each aminoacid in each element 2610 - `str` with the sequence using the one letter code 2611 - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-" 2612 2613 ''' 2614 # Aminoacid key 2615 keys={"ALA": "A", 2616 "ARG": "R", 2617 "ASN": "N", 2618 "ASP": "D", 2619 "CYS": "C", 2620 "GLU": "E", 2621 "GLN": "Q", 2622 "GLY": "G", 2623 "HIS": "H", 2624 "ILE": "I", 2625 "LEU": "L", 2626 "LYS": "K", 2627 "MET": "M", 2628 "PHE": "F", 2629 "PRO": "P", 2630 "SER": "S", 2631 "THR": "T", 2632 "TRP": "W", 2633 "TYR": "Y", 2634 "VAL": "V", 2635 "PSER": "J", 2636 "PTHR": "U", 2637 "PTyr": "Z", 2638 "NH2": "n", 2639 "COOH": "c"} 2640 clean_sequence=[] 2641 if isinstance(sequence, str): 2642 if sequence.find("-") != -1: 2643 splited_sequence=sequence.split("-") 2644 for residue in splited_sequence: 2645 if len(residue) == 1: 2646 if residue in keys.values(): 2647 residue_ok=residue 2648 else: 2649 if residue.upper() in keys.values(): 2650 residue_ok=residue.upper() 2651 else: 2652 raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence") 2653 clean_sequence.append(residue_ok) 2654 else: 2655 if residue in keys.keys(): 2656 clean_sequence.append(keys[residue]) 2657 else: 2658 if residue.upper() in keys.keys(): 2659 clean_sequence.append(keys[residue.upper()]) 2660 else: 2661 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2662 else: 2663 for residue in sequence: 2664 if residue in keys.values(): 2665 residue_ok=residue 2666 else: 2667 if residue.upper() in keys.values(): 2668 residue_ok=residue.upper() 2669 else: 2670 raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence") 2671 clean_sequence.append(residue_ok) 2672 if isinstance(sequence, list): 2673 for residue in sequence: 2674 if residue in keys.values(): 2675 residue_ok=residue 2676 else: 2677 if residue.upper() in keys.values(): 2678 residue_ok=residue.upper() 2679 elif (residue.upper() in keys.keys()): 2680 residue_ok= keys[residue.upper()] 2681 else: 2682 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2683 clean_sequence.append(residue_ok) 2684 return clean_sequence
Parses sequence
to the one letter code for amino acids.
Arguments:
- sequence(
str
orlst
): Sequence of the amino acid.
Returns:
clean_sequence(
lst
):sequence
using the one letter code.
Note:
- Accepted formats for
sequence
are:
lst
with one letter or three letter code of each aminoacid in each elementstr
with the sequence using the one letter codestr
with the squence using the three letter code, each aminoacid must be separated by a hyphen "-"
2686 def read_pmb_df (self,filename): 2687 """ 2688 Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE. 2689 2690 Args: 2691 filename(`str`): path to file. 2692 2693 Note: 2694 This function only accepts files with CSV format. 2695 """ 2696 2697 if filename.rsplit(".", 1)[1] != "csv": 2698 raise ValueError("Only files with CSV format are supported") 2699 df = pd.read_csv (filename,header=[0, 1], index_col=0) 2700 columns_names = self.setup_df() 2701 2702 multi_index = pd.MultiIndex.from_tuples(columns_names) 2703 df.columns = multi_index 2704 2705 self.df = self.convert_columns_to_original_format(df) 2706 self.df.fillna(pd.NA, inplace=True) 2707 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.
2709 def read_protein_vtf_in_df (self,filename,unit_length=None): 2710 """ 2711 Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`. 2712 2713 Args: 2714 filename(`str`): Path to the vtf file with the coarse-grained model. 2715 unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None. 2716 2717 Returns: 2718 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 2719 2720 Note: 2721 - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom. 2722 """ 2723 2724 logging.info(f'Loading protein coarse grain model file: {filename}') 2725 2726 coord_list = [] 2727 particles_dict = {} 2728 2729 if unit_length is None: 2730 unit_length = 1 * self.units.angstrom 2731 2732 with open (filename,'r') as protein_model: 2733 for line in protein_model : 2734 line_split = line.split() 2735 if line_split: 2736 line_header = line_split[0] 2737 if line_header == 'atom': 2738 atom_id = line_split[1] 2739 atom_name = line_split[3] 2740 atom_resname = line_split[5] 2741 chain_id = line_split[9] 2742 radius = float(line_split [11])*unit_length 2743 particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, radius] 2744 elif line_header.isnumeric(): 2745 atom_coord = line_split[1:] 2746 atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord] 2747 coord_list.append (atom_coord) 2748 2749 numbered_label = [] 2750 i = 0 2751 2752 for atom_id in particles_dict.keys(): 2753 2754 if atom_id == 1: 2755 atom_name = particles_dict[atom_id][0] 2756 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2757 numbered_label.append(numbered_name) 2758 2759 elif atom_id != 1: 2760 2761 if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]: 2762 i += 1 2763 count = 1 2764 atom_name = particles_dict[atom_id][0] 2765 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2766 numbered_label.append(numbered_name) 2767 2768 elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]: 2769 if count == 2 or particles_dict[atom_id][1] == 'GLY': 2770 i +=1 2771 count = 0 2772 atom_name = particles_dict[atom_id][0] 2773 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2774 numbered_label.append(numbered_name) 2775 count +=1 2776 2777 topology_dict = {} 2778 2779 for i in range (0, len(numbered_label)): 2780 topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] , 2781 'chain_id':numbered_label[i][1], 2782 'radius':numbered_label[i][2] } 2783 2784 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 infilename
using the pyMBE UnitRegistry. Defaults to None.
Returns:
topology_dict(
dict
): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}
Note:
- If no
unit_length
is provided, it is assumed that the coordinates are in Angstrom.
2786 def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2787 """ 2788 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2789 If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead. 2790 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. 2791 2792 Args: 2793 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 2794 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 2795 hard_check(`bool`, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2796 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. 2797 2798 Returns: 2799 bond(`espressomd.interactions.BondedInteractions`): bond object from the espressomd library. 2800 2801 Note: 2802 - 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`. 2803 - If `hard_check`=`True` stops the code when no bond is found. 2804 """ 2805 2806 bond_key = self.find_bond_key(particle_name1=particle_name1, 2807 particle_name2=particle_name2, 2808 use_default_bond=use_default_bond) 2809 if use_default_bond: 2810 if not self._check_if_name_is_defined_in_df(name="default"): 2811 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") 2812 if bond_key: 2813 return self.df[self.df['name']==bond_key].bond_object.values[0] 2814 else: 2815 msg= f"Bond not defined between particles {particle_name1} and {particle_name2}" 2816 if hard_check: 2817 raise ValueError(msg) 2818 else: 2819 logging.warning(msg) 2820 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_name1
andparticle_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_name1
andparticle_name2
, it returns the default bond defined inpmb.df
.- If
hard_check
=True
stops the code when no bond is found.
2821 def search_particles_in_residue(self, residue_name): 2822 ''' 2823 Searches for all particles in a given residue of name `residue_name`. 2824 2825 Args: 2826 residue_name (`str`): name of the residue to be searched 2827 2828 Returns: 2829 list_of_particles_in_residue (`lst`): list of the names of all particles in the residue 2830 2831 Note: 2832 - 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. 2833 - The function will return an empty list if the residue is not defined in `pmb.df`. 2834 - The function will return an empty list if the particles are not defined in the pyMBE DataFrame. 2835 ''' 2836 if not self._check_if_name_is_defined_in_df(name=residue_name): 2837 logging.warning(f"Residue {residue_name} not defined in pmb.df") 2838 return [] 2839 self._check_if_name_has_right_type(name=residue_name, expected_pmb_type="residue") 2840 index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 2841 central_bead = self.df.at [index_residue, ('central_bead', '')] 2842 list_of_side_chains = self.df.at[index_residue, ('side_chains', '')] 2843 list_of_particles_in_residue = [] 2844 if central_bead is not pd.NA: 2845 if self._check_if_name_is_defined_in_df(name=central_bead): 2846 if self._check_if_name_has_right_type(name=central_bead, expected_pmb_type="particle", hard_check=False): 2847 list_of_particles_in_residue.append(central_bead) 2848 if list_of_side_chains is not pd.NA: 2849 for side_chain in list_of_side_chains: 2850 if self._check_if_name_is_defined_in_df(name=side_chain): 2851 object_type = self.df[self.df['name']==side_chain].pmb_type.values[0] 2852 else: 2853 continue 2854 if object_type == "residue": 2855 list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain) 2856 list_of_particles_in_residue += list_of_particles_in_side_chain_residue 2857 elif object_type == "particle": 2858 if side_chain is not pd.NA: 2859 list_of_particles_in_residue.append(side_chain) 2860 return list_of_particles_in_residue
Searches for all particles in a given residue of name residue_name
.
Arguments:
- residue_name (
str
): name of the residue to be searched
Returns:
list_of_particles_in_residue (
lst
): list of the names of all particles in the residue
Note:
- The function returns a name per particle in residue, i.e. if there are multiple particles with the same type
list_of_particles_in_residue
will have repeated items.- 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.
2862 def set_particle_acidity(self, name, acidity=pd.NA, default_charge_number=0, pka=pd.NA, overwrite=True): 2863 """ 2864 Sets the particle acidity including the charges in each of its possible states. 2865 2866 Args: 2867 name(`str`): Unique label that identifies the `particle`. 2868 acidity(`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to None. 2869 default_charge_number (`int`): Charge number of the particle. Defaults to 0. 2870 pka(`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to pandas.NA. 2871 overwrite(`bool`, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False. 2872 2873 Note: 2874 - For particles with `acidity = acidic` or `acidity = basic`, `state_one` and `state_two` correspond to the protonated and 2875 deprotonated states, respectively. 2876 - For particles without an acidity `acidity = pandas.NA`, only `state_one` is defined. 2877 - Each state has the following properties as sub-indexes: `label`,`charge` and `es_type`. 2878 """ 2879 acidity_valid_keys = ['inert','acidic', 'basic'] 2880 if not pd.isna(acidity): 2881 if acidity not in acidity_valid_keys: 2882 raise ValueError(f"Acidity {acidity} provided for particle name {name} is not supproted. Valid keys are: {acidity_valid_keys}") 2883 if acidity in ['acidic', 'basic'] and pd.isna(pka): 2884 raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.") 2885 if acidity == "inert": 2886 acidity = pd.NA 2887 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") 2888 2889 self.define_particle_entry_in_df(name=name) 2890 2891 for index in self.df[self.df['name']==name].index: 2892 if pka is not pd.NA: 2893 self.add_value_to_df(key=('pka',''), 2894 index=index, 2895 new_value=pka, 2896 overwrite=overwrite) 2897 2898 self.add_value_to_df(key=('acidity',''), 2899 index=index, 2900 new_value=acidity, 2901 overwrite=overwrite) 2902 if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')): 2903 self.add_value_to_df(key=('state_one','es_type'), 2904 index=index, 2905 new_value=self.propose_unused_type(), 2906 overwrite=overwrite) 2907 if pd.isna(self.df.loc [self.df['name'] == name].acidity.iloc[0]): 2908 self.add_value_to_df(key=('state_one','z'), 2909 index=index, 2910 new_value=default_charge_number, 2911 overwrite=overwrite) 2912 self.add_value_to_df(key=('state_one','label'), 2913 index=index, 2914 new_value=name, 2915 overwrite=overwrite) 2916 else: 2917 protonated_label = f'{name}H' 2918 self.add_value_to_df(key=('state_one','label'), 2919 index=index, 2920 new_value=protonated_label, 2921 overwrite=overwrite) 2922 self.add_value_to_df(key=('state_two','label'), 2923 index=index, 2924 new_value=name, 2925 overwrite=overwrite) 2926 if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')): 2927 self.add_value_to_df(key=('state_two','es_type'), 2928 index=index, 2929 new_value=self.propose_unused_type(), 2930 overwrite=overwrite) 2931 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'acidic': 2932 self.add_value_to_df(key=('state_one','z'), 2933 index=index,new_value=0, 2934 overwrite=overwrite) 2935 self.add_value_to_df(key=('state_two','z'), 2936 index=index, 2937 new_value=-1, 2938 overwrite=overwrite) 2939 elif self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'basic': 2940 self.add_value_to_df(key=('state_one','z'), 2941 index=index,new_value=+1, 2942 overwrite=overwrite) 2943 self.add_value_to_df(key=('state_two','z'), 2944 index=index, 2945 new_value=0, 2946 overwrite=overwrite) 2947 self.df.fillna(pd.NA, inplace=True) 2948 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 isacidic
orbasic
, used to setup constant pH simulations. Defaults to None. - default_charge_number (
int
): Charge number of the particle. Defaults to 0. - pka(
float
, optional): Ifparticle
is an acid or a base, it defines its pka-value. Defaults to pandas.NA. - overwrite(
bool
, optional): Switch to enable overwriting of already existing values in pmb.df. Defaults to False.
Note:
- For particles with
acidity = acidic
oracidity = basic
,state_one
andstate_two
correspond to the protonated and
deprotonated states, respectively.
- For particles without an acidity acidity = pandas.NA
, only state_one
is defined.
- Each state has the following properties as sub-indexes: label
,charge
and es_type
.
2950 def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None): 2951 """ 2952 Sets the set of reduced units used by pyMBE.units and it prints it. 2953 2954 Args: 2955 unit_length(`pint.Quantity`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 2956 unit_charge(`pint.Quantity`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 2957 temperature(`pint.Quantity`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 2958 Kw(`pint.Quantity`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 2959 2960 Note: 2961 - If no `temperature` is given, a value of 298.15 K is assumed by default. 2962 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 2963 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 2964 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 2965 """ 2966 if unit_length is None: 2967 unit_length= 0.355*self.units.nm 2968 if temperature is None: 2969 temperature = 298.15 * self.units.K 2970 if unit_charge is None: 2971 unit_charge = scipy.constants.e * self.units.C 2972 if Kw is None: 2973 Kw = 1e-14 2974 # Sanity check 2975 variables=[unit_length,temperature,unit_charge] 2976 dimensionalities=["[length]","[temperature]","[charge]"] 2977 for variable,dimensionality in zip(variables,dimensionalities): 2978 self.check_dimensionality(variable,dimensionality) 2979 self.Kw=Kw*self.units.mol**2 / (self.units.l**2) 2980 self.kT=temperature*self.kB 2981 self.units._build_cache() 2982 self.units.define(f'reduced_energy = {self.kT} ') 2983 self.units.define(f'reduced_length = {unit_length}') 2984 self.units.define(f'reduced_charge = {unit_charge}') 2985 logging.info(self.get_reduced_units()) 2986 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.units
UnitRegistry. Defaults to None. - unit_charge(
pint.Quantity
,optional): Reduced unit of charge defined using thepmb.units
UnitRegistry. Defaults to None. - temperature(
pint.Quantity
,optional): Temperature of the system, defined using thepmb.units
UnitRegistry. Defaults to None. - Kw(
pint.Quantity
,optional): Ionic product of water in mol^2/l^2. Defaults to None.
Note:
- If no
temperature
is given, a value of 298.15 K is assumed by default.- If no
unit_length
is given, a value of 0.355 nm is assumed by default.- If no
unit_charge
is given, a value of 1 elementary charge is assumed by default.- If no
Kw
is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default.
2988 def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2989 """ 2990 Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 2991 2992 Args: 2993 counter_ion(`str`): `name` of the counter_ion `particle`. 2994 constant_pH(`float`): pH-value. 2995 exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted. 2996 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2997 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2998 2999 Returns: 3000 RE(`reaction_methods.ConstantpHEnsemble`): Instance of a reaction_methods.ConstantpHEnsemble object from the espressomd library. 3001 sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3002 """ 3003 from espressomd import reaction_methods 3004 if exclusion_range is None: 3005 exclusion_range = max(self.get_radius_map().values())*2.0 3006 if pka_set is None: 3007 pka_set=self.get_pka_set() 3008 self.check_pka_set(pka_set=pka_set) 3009 if use_exclusion_radius_per_type: 3010 exclusion_radius_per_type = self.get_radius_map() 3011 else: 3012 exclusion_radius_per_type = {} 3013 3014 RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3015 exclusion_range=exclusion_range, 3016 seed=self.seed, 3017 constant_pH=constant_pH, 3018 exclusion_radius_per_type = exclusion_radius_per_type 3019 ) 3020 sucessfull_reactions_labels=[] 3021 charge_number_map = self.get_charge_number_map() 3022 for name in pka_set.keys(): 3023 if not self._check_if_name_is_defined_in_df(name): 3024 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.') 3025 continue 3026 gamma=10**-pka_set[name]['pka_value'] 3027 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3028 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3029 counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0] 3030 RE.add_reaction(gamma=gamma, 3031 reactant_types=[state_one_type], 3032 product_types=[state_two_type, counter_ion_type], 3033 default_charges={state_one_type: charge_number_map[state_one_type], 3034 state_two_type: charge_number_map[state_two_type], 3035 counter_ion_type: charge_number_map[counter_ion_type]}) 3036 sucessfull_reactions_labels.append(name) 3037 return RE, sucessfull_reactions_labels
Sets up the Acid/Base reactions for acidic/basic particles
defined in pmb.df
to be sampled in the constant pH ensemble.
Arguments:
- counter_ion(
str
):name
of the counter_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.
3039 def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient, exclusion_range=None, use_exclusion_radius_per_type = False): 3040 """ 3041 Sets up grand-canonical coupling to a reservoir of salt. 3042 For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead. 3043 3044 Args: 3045 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3046 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 3047 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 3048 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3049 exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted. 3050 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 3051 3052 Returns: 3053 RE (`reaction_methods.ReactionEnsemble`): Instance of a reaction_methods.ReactionEnsemble object from the espressomd library. 3054 """ 3055 from espressomd import reaction_methods 3056 if exclusion_range is None: 3057 exclusion_range = max(self.get_radius_map().values())*2.0 3058 if use_exclusion_radius_per_type: 3059 exclusion_radius_per_type = self.get_radius_map() 3060 else: 3061 exclusion_radius_per_type = {} 3062 3063 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3064 exclusion_range=exclusion_range, 3065 seed=self.seed, 3066 exclusion_radius_per_type = exclusion_radius_per_type 3067 ) 3068 3069 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3070 determined_activity_coefficient = activity_coefficient(c_salt_res) 3071 K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient 3072 3073 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 3074 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 3075 3076 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0] 3077 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0] 3078 3079 if salt_cation_charge <= 0: 3080 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 3081 if salt_anion_charge >= 0: 3082 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 3083 3084 # Grand-canonical coupling to the reservoir 3085 RE.add_reaction( 3086 gamma = K_salt.magnitude, 3087 reactant_types = [], 3088 reactant_coefficients = [], 3089 product_types = [ salt_cation_es_type, salt_anion_es_type ], 3090 product_coefficients = [ 1, 1 ], 3091 default_charges = { 3092 salt_cation_es_type: salt_cation_charge, 3093 salt_anion_es_type: salt_anion_charge, 3094 } 3095 ) 3096 3097 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.
3099 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): 3100 """ 3101 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 3102 reservoir of small ions. 3103 This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1]. 3104 3105 [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. 3106 3107 Args: 3108 pH_res ('float): pH-value in the reservoir. 3109 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3110 proton_name ('str'): Name of the proton (H+) particle. 3111 hydroxide_name ('str'): Name of the hydroxide (OH-) particle. 3112 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 3113 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 3114 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3115 exclusion_range(`pint.Quantity`, optional): For distances shorter than this value, no particles will be inserted. 3116 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 3117 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 3118 3119 Returns: 3120 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 3121 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3122 ionic_strength_res ('pint.Quantity'): Ionic strength of the reservoir (useful for calculating partition coefficients). 3123 """ 3124 from espressomd import reaction_methods 3125 if exclusion_range is None: 3126 exclusion_range = max(self.get_radius_map().values())*2.0 3127 if pka_set is None: 3128 pka_set=self.get_pka_set() 3129 self.check_pka_set(pka_set=pka_set) 3130 if use_exclusion_radius_per_type: 3131 exclusion_radius_per_type = self.get_radius_map() 3132 else: 3133 exclusion_radius_per_type = {} 3134 3135 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3136 exclusion_range=exclusion_range, 3137 seed=self.seed, 3138 exclusion_radius_per_type = exclusion_radius_per_type 3139 ) 3140 3141 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3142 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 3143 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 3144 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 3145 K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3146 K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3147 K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 3148 3149 proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0] 3150 hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0] 3151 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 3152 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 3153 3154 proton_charge = self.df.loc[self.df['name']==proton_name].state_one.z.values[0] 3155 hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.z.values[0] 3156 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.z.values[0] 3157 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.z.values[0] 3158 3159 if proton_charge <= 0: 3160 raise ValueError('ERROR proton charge must be positive, charge ', proton_charge) 3161 if salt_cation_charge <= 0: 3162 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 3163 if hydroxide_charge >= 0: 3164 raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge) 3165 if salt_anion_charge >= 0: 3166 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 3167 3168 # Grand-canonical coupling to the reservoir 3169 # 0 = H+ + OH- 3170 RE.add_reaction( 3171 gamma = K_W.magnitude, 3172 reactant_types = [], 3173 reactant_coefficients = [], 3174 product_types = [ proton_es_type, hydroxide_es_type ], 3175 product_coefficients = [ 1, 1 ], 3176 default_charges = { 3177 proton_es_type: proton_charge, 3178 hydroxide_es_type: hydroxide_charge, 3179 } 3180 ) 3181 3182 # 0 = Na+ + Cl- 3183 RE.add_reaction( 3184 gamma = K_NACL.magnitude, 3185 reactant_types = [], 3186 reactant_coefficients = [], 3187 product_types = [ salt_cation_es_type, salt_anion_es_type ], 3188 product_coefficients = [ 1, 1 ], 3189 default_charges = { 3190 salt_cation_es_type: salt_cation_charge, 3191 salt_anion_es_type: salt_anion_charge, 3192 } 3193 ) 3194 3195 # 0 = Na+ + OH- 3196 RE.add_reaction( 3197 gamma = (K_NACL * K_W / K_HCL).magnitude, 3198 reactant_types = [], 3199 reactant_coefficients = [], 3200 product_types = [ salt_cation_es_type, hydroxide_es_type ], 3201 product_coefficients = [ 1, 1 ], 3202 default_charges = { 3203 salt_cation_es_type: salt_cation_charge, 3204 hydroxide_es_type: hydroxide_charge, 3205 } 3206 ) 3207 3208 # 0 = H+ + Cl- 3209 RE.add_reaction( 3210 gamma = K_HCL.magnitude, 3211 reactant_types = [], 3212 reactant_coefficients = [], 3213 product_types = [ proton_es_type, salt_anion_es_type ], 3214 product_coefficients = [ 1, 1 ], 3215 default_charges = { 3216 proton_es_type: proton_charge, 3217 salt_anion_es_type: salt_anion_charge, 3218 } 3219 ) 3220 3221 # Annealing moves to ensure sufficient sampling 3222 # Cation annealing H+ = Na+ 3223 RE.add_reaction( 3224 gamma = (K_NACL / K_HCL).magnitude, 3225 reactant_types = [proton_es_type], 3226 reactant_coefficients = [ 1 ], 3227 product_types = [ salt_cation_es_type ], 3228 product_coefficients = [ 1 ], 3229 default_charges = { 3230 proton_es_type: proton_charge, 3231 salt_cation_es_type: salt_cation_charge, 3232 } 3233 ) 3234 3235 # Anion annealing OH- = Cl- 3236 RE.add_reaction( 3237 gamma = (K_HCL / K_W).magnitude, 3238 reactant_types = [hydroxide_es_type], 3239 reactant_coefficients = [ 1 ], 3240 product_types = [ salt_anion_es_type ], 3241 product_coefficients = [ 1 ], 3242 default_charges = { 3243 hydroxide_es_type: hydroxide_charge, 3244 salt_anion_es_type: salt_anion_charge, 3245 } 3246 ) 3247 3248 sucessful_reactions_labels=[] 3249 charge_number_map = self.get_charge_number_map() 3250 for name in pka_set.keys(): 3251 if not self._check_if_name_is_defined_in_df(name): 3252 logging.warning(f'The acid-base reaction of {name} has not been set up because its particle type is not defined in the dataframe.') 3253 continue 3254 3255 Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 3256 3257 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3258 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3259 3260 # Reaction in terms of proton: HA = A + H+ 3261 RE.add_reaction(gamma=Ka.magnitude, 3262 reactant_types=[state_one_type], 3263 reactant_coefficients=[1], 3264 product_types=[state_two_type, proton_es_type], 3265 product_coefficients=[1, 1], 3266 default_charges={state_one_type: charge_number_map[state_one_type], 3267 state_two_type: charge_number_map[state_two_type], 3268 proton_es_type: proton_charge}) 3269 3270 # Reaction in terms of salt cation: HA = A + Na+ 3271 RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude, 3272 reactant_types=[state_one_type], 3273 reactant_coefficients=[1], 3274 product_types=[state_two_type, salt_cation_es_type], 3275 product_coefficients=[1, 1], 3276 default_charges={state_one_type: charge_number_map[state_one_type], 3277 state_two_type: charge_number_map[state_two_type], 3278 salt_cation_es_type: salt_cation_charge}) 3279 3280 # Reaction in terms of hydroxide: OH- + HA = A 3281 RE.add_reaction(gamma=(Ka / K_W).magnitude, 3282 reactant_types=[state_one_type, hydroxide_es_type], 3283 reactant_coefficients=[1, 1], 3284 product_types=[state_two_type], 3285 product_coefficients=[1], 3286 default_charges={state_one_type: charge_number_map[state_one_type], 3287 state_two_type: charge_number_map[state_two_type], 3288 hydroxide_es_type: hydroxide_charge}) 3289 3290 # Reaction in terms of salt anion: Cl- + HA = A 3291 RE.add_reaction(gamma=(Ka / K_HCL).magnitude, 3292 reactant_types=[state_one_type, salt_anion_es_type], 3293 reactant_coefficients=[1, 1], 3294 product_types=[state_two_type], 3295 product_coefficients=[1], 3296 default_charges={state_one_type: charge_number_map[state_one_type], 3297 state_two_type: charge_number_map[state_two_type], 3298 salt_anion_es_type: salt_anion_charge}) 3299 3300 sucessful_reactions_labels.append(name) 3301 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).
3303 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): 3304 """ 3305 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 3306 reservoir of small ions. 3307 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-}. 3308 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'. 3309 3310 [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4). 3311 [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. 3312 3313 Args: 3314 pH_res ('float'): pH-value in the reservoir. 3315 c_salt_res ('pint.Quantity'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 3316 cation_name ('str'): Name of the cationic particle. 3317 anion_name ('str'): Name of the anionic particle. 3318 activity_coefficient ('callable'): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 3319 exclusion_range(`pint.Quantity`, optional): Below this value, no particles will be inserted. 3320 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 3321 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`. 3322 3323 Returns: 3324 RE (`reaction_ensemble.ReactionEnsemble`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 3325 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 3326 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 3327 """ 3328 from espressomd import reaction_methods 3329 if exclusion_range is None: 3330 exclusion_range = max(self.get_radius_map().values())*2.0 3331 if pka_set is None: 3332 pka_set=self.get_pka_set() 3333 self.check_pka_set(pka_set=pka_set) 3334 if use_exclusion_radius_per_type: 3335 exclusion_radius_per_type = self.get_radius_map() 3336 else: 3337 exclusion_radius_per_type = {} 3338 3339 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 3340 exclusion_range=exclusion_range, 3341 seed=self.seed, 3342 exclusion_radius_per_type = exclusion_radius_per_type 3343 ) 3344 3345 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 3346 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 3347 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 3348 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 3349 a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 3350 a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 3351 a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 3352 K_XX = a_cation * a_anion 3353 3354 cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0] 3355 anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0] 3356 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.z.values[0] 3357 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.z.values[0] 3358 if cation_charge <= 0: 3359 raise ValueError('ERROR cation charge must be positive, charge ', cation_charge) 3360 if anion_charge >= 0: 3361 raise ValueError('ERROR anion charge must be negative, charge ', anion_charge) 3362 3363 # Coupling to the reservoir: 0 = X+ + X- 3364 RE.add_reaction( 3365 gamma = K_XX.magnitude, 3366 reactant_types = [], 3367 reactant_coefficients = [], 3368 product_types = [ cation_es_type, anion_es_type ], 3369 product_coefficients = [ 1, 1 ], 3370 default_charges = { 3371 cation_es_type: cation_charge, 3372 anion_es_type: anion_charge, 3373 } 3374 ) 3375 3376 sucessful_reactions_labels=[] 3377 charge_number_map = self.get_charge_number_map() 3378 for name in pka_set.keys(): 3379 if not self._check_if_name_is_defined_in_df(name): 3380 logging.warning(f'The acid-base reaction of {name} has not been set up because its particle type is not defined in the dataframe.') 3381 continue 3382 3383 Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l 3384 gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen 3385 3386 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 3387 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 3388 3389 # Reaction in terms of small cation: HA = A + X+ 3390 RE.add_reaction(gamma=gamma_K_AX.magnitude, 3391 reactant_types=[state_one_type], 3392 reactant_coefficients=[1], 3393 product_types=[state_two_type, cation_es_type], 3394 product_coefficients=[1, 1], 3395 default_charges={state_one_type: charge_number_map[state_one_type], 3396 state_two_type: charge_number_map[state_two_type], 3397 cation_es_type: cation_charge}) 3398 3399 # Reaction in terms of small anion: X- + HA = A 3400 RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude, 3401 reactant_types=[state_one_type, anion_es_type], 3402 reactant_coefficients=[1, 1], 3403 product_types=[state_two_type], 3404 product_coefficients=[1], 3405 default_charges={state_one_type: charge_number_map[state_one_type], 3406 state_two_type: charge_number_map[state_two_type], 3407 anion_es_type: anion_charge}) 3408 3409 sucessful_reactions_labels.append(name) 3410 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).
3412 def setup_df (self): 3413 """ 3414 Sets up the pyMBE's dataframe `pymbe.df`. 3415 3416 Returns: 3417 columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe 3418 """ 3419 3420 columns_dtypes = { 3421 'name': { 3422 '': str}, 3423 'pmb_type': { 3424 '': str}, 3425 'particle_id': { 3426 '': pd.Int64Dtype()}, 3427 'particle_id2': { 3428 '': pd.Int64Dtype()}, 3429 'residue_id': { 3430 '': pd.Int64Dtype()}, 3431 'molecule_id': { 3432 '': pd.Int64Dtype()}, 3433 'acidity': { 3434 '': str}, 3435 'pka': { 3436 '': object}, 3437 'central_bead': { 3438 '': object}, 3439 'side_chains': { 3440 '': object}, 3441 'residue_list': { 3442 '': object}, 3443 'model': { 3444 '': str}, 3445 'sigma': { 3446 '': object}, 3447 'cutoff': { 3448 '': object}, 3449 'offset': { 3450 '': object}, 3451 'epsilon': { 3452 '': object}, 3453 'state_one': { 3454 'label': str, 3455 'es_type': pd.Int64Dtype(), 3456 'z': pd.Int64Dtype()}, 3457 'state_two': { 3458 'label': str, 3459 'es_type': pd.Int64Dtype(), 3460 'z': pd.Int64Dtype()}, 3461 'sequence': { 3462 '': object}, 3463 'bond_object': { 3464 '': object}, 3465 'parameters_of_the_potential':{ 3466 '': object}, 3467 'l0': { 3468 '': float}, 3469 'node_map':{ 3470 '':object}, 3471 'chain_map':{ 3472 '':object}} 3473 3474 self.df = pd.DataFrame(columns=pd.MultiIndex.from_tuples([(col_main, col_sub) for col_main, sub_cols in columns_dtypes.items() for col_sub in sub_cols.keys()])) 3475 3476 for level1, sub_dtypes in columns_dtypes.items(): 3477 for level2, dtype in sub_dtypes.items(): 3478 self.df[level1, level2] = self.df[level1, level2].astype(dtype) 3479 3480 columns_names = pd.MultiIndex.from_frame(self.df) 3481 columns_names = columns_names.names 3482 3483 return columns_names
Sets up the pyMBE's dataframe pymbe.df
.
Returns:
columns_names(
obj
): pandas multiindex object with the column names of the pyMBE's dataframe
3485 def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot'): 3486 """ 3487 Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`. 3488 3489 Args: 3490 espresso_system(`espressomd.system.System`): Instance of a system object from the espressomd library. 3491 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. 3492 combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'. 3493 warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True. 3494 3495 Note: 3496 - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 3497 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 3498 - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html 3499 3500 """ 3501 from itertools import combinations_with_replacement 3502 compulsory_parameters_in_df = ['sigma','epsilon'] 3503 shift=0 3504 if shift_potential: 3505 shift="auto" 3506 # List which particles have sigma and epsilon values defined in pmb.df and which ones don't 3507 particles_types_with_LJ_parameters = [] 3508 non_parametrized_labels= [] 3509 for particle_type in self.get_type_map().values(): 3510 check_list=[] 3511 for key in compulsory_parameters_in_df: 3512 value_in_df=self.find_value_from_es_type(es_type=particle_type, 3513 column_name=key) 3514 check_list.append(pd.isna(value_in_df)) 3515 if any(check_list): 3516 non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 3517 column_name='label')) 3518 else: 3519 particles_types_with_LJ_parameters.append(particle_type) 3520 # Set up LJ interactions between all particle types 3521 for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 3522 particle_name1 = self.find_value_from_es_type(es_type=type_pair[0], 3523 column_name="name") 3524 particle_name2 = self.find_value_from_es_type(es_type=type_pair[1], 3525 column_name="name") 3526 lj_parameters= self.get_lj_parameters(particle_name1 = particle_name1, 3527 particle_name2 = particle_name2, 3528 combining_rule = combining_rule) 3529 3530 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 3531 if not lj_parameters: 3532 continue 3533 espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = lj_parameters["epsilon"].to('reduced_energy').magnitude, 3534 sigma = lj_parameters["sigma"].to('reduced_length').magnitude, 3535 cutoff = lj_parameters["cutoff"].to('reduced_length').magnitude, 3536 offset = lj_parameters["offset"].to("reduced_length").magnitude, 3537 shift = shift) 3538 index = len(self.df) 3539 label1 = self.find_value_from_es_type(es_type=type_pair[0], column_name="label") 3540 label2 = self.find_value_from_es_type(es_type=type_pair[1], column_name="label") 3541 self.df.at [index, 'name'] = f'LJ: {label1}-{label2}' 3542 lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params() 3543 3544 self.add_value_to_df(index=index, 3545 key=('pmb_type',''), 3546 new_value='LennardJones') 3547 3548 self.add_value_to_df(index=index, 3549 key=('parameters_of_the_potential',''), 3550 new_value=lj_params, 3551 non_standard_value=True) 3552 if non_parametrized_labels: 3553 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.') 3554 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 calculatesigma
andepsilon
for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'. - warning(
bool
, optional): switch to activate/deactivate warning messages. Defaults to True.
Note:
- LJ interactions will only be set up between particles with defined values of
sigma
andepsilon
in the pmb.df.- Currently, the only
combining_rule
supported is Lorentz-Berthelot.- Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html
3556 def write_pmb_df (self, filename): 3557 ''' 3558 Writes the pyMBE dataframe into a csv file 3559 3560 Args: 3561 filename(`str`): Path to the csv file 3562 ''' 3563 3564 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence'] 3565 df = self.df.copy(deep=True) 3566 for column_name in columns_with_list_or_dict: 3567 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) 3568 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) 3569 df.fillna(pd.NA, inplace=True) 3570 df.to_csv(filename) 3571 return
Writes the pyMBE dataframe into a csv file
Arguments:
- filename(
str
): Path to the csv file
43 class NumpyEncoder(json.JSONEncoder): 44 """ 45 Custom JSON encoder that converts NumPy arrays to Python lists 46 and NumPy scalars to Python scalars. 47 """ 48 def default(self, obj): 49 if isinstance(obj, np.ndarray): 50 return obj.tolist() 51 if isinstance(obj, np.generic): 52 return obj.item() 53 return super().default(obj)
Custom JSON encoder that converts NumPy arrays to Python lists and NumPy scalars to Python scalars.
48 def default(self, obj): 49 if isinstance(obj, np.ndarray): 50 return obj.tolist() 51 if isinstance(obj, np.generic): 52 return obj.item() 53 return super().default(obj)
Implement this method in a subclass such that it returns
a serializable object for o
, or calls the base implementation
(to raise a TypeError
).
For example, to support arbitrary iterators, you could implement default like this::
def default(self, o):
try:
iterable = iter(o)
except TypeError:
pass
else:
return list(iterable)
# Let the base class default method raise the TypeError
return JSONEncoder.default(self, o)