pyMBE
1import re 2import sys 3import ast 4import json 5import pint 6import numpy as np 7import pandas as pd 8import scipy.optimize 9 10 11class pymbe_library(): 12 13 """ 14 The library for the Molecular Builder for ESPResSo (pyMBE) 15 16 Attributes: 17 N_A(`obj`): Avogadro number using the `pmb.units` UnitRegistry. 18 Kb(`obj`): Boltzmann constant using the `pmb.units` UnitRegistry. 19 e(`obj`): Elemental charge using the `pmb.units` UnitRegistry. 20 df(`obj`): PandasDataframe used to bookkeep all the information stored in pyMBE. Typically refered as `pmb.df`. 21 kT(`obj`): Thermal energy using the `pmb.units` UnitRegistry. It is used as the unit of reduced energy. 22 Kw(`obj`): Ionic product of water using the `pmb.units` UnitRegistry. Used in the setup of the G-RxMC method. 23 """ 24 units = pint.UnitRegistry() 25 N_A=6.02214076e23 / units.mol 26 Kb=1.38064852e-23 * units.J / units.K 27 e=1.60217662e-19 *units.C 28 df=None 29 kT=None 30 Kw=None 31 SEED=None 32 rng=None 33 34 35 class NumpyEncoder(json.JSONEncoder): 36 37 """ 38 Custom JSON encoder that converts NumPy arrays to Python lists 39 and NumPy scalars to Python scalars. 40 """ 41 42 def default(self, obj): 43 if isinstance(obj, np.ndarray): 44 return obj.tolist() 45 if isinstance(obj, np.generic): 46 return obj.item() 47 return super().default(obj) 48 49 50 def __init__(self, SEED, temperature=None, unit_length=None, unit_charge=None, Kw=None): 51 """ 52 Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 53 and sets up the `pmb.df` for bookkeeping. 54 55 Args: 56 temperature(`obj`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. 57 unit_length(`obj`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. 58 unit_charge (`obj`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 59 Kw (`obj`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 60 61 Note: 62 - If no `temperature` is given, a value of 298.15 K is assumed by default. 63 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 64 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 65 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 66 """ 67 # Seed and RNG 68 self.SEED=SEED 69 self.rng = np.random.default_rng(SEED) 70 self.set_reduced_units(unit_length=unit_length, unit_charge=unit_charge, 71 temperature=temperature, Kw=Kw, verbose=False) 72 self.setup_df() 73 return 74 75 def enable_motion_of_rigid_object(self, name, espresso_system): 76 ''' 77 Enables the motion of the rigid object `name` in the `espresso_system`. 78 79 Args: 80 name(`str`): Label of the object. 81 espresso_system(`obj`): Instance of a system object from the espressomd library. 82 83 Note: 84 - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]. 85 ''' 86 print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]') 87 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 88 if pmb_type != 'protein': 89 raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein') 90 molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list() 91 for molecule_id in molecule_ids_list: 92 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 93 center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system) 94 rigid_object_center = espresso_system.part.add(pos=center_of_mass, 95 rotation=[True,True,True], 96 type=self.propose_unused_type()) 97 for particle_id in particle_ids_list: 98 pid = espresso_system.part.by_id(particle_id) 99 pid.vs_auto_relate_to(rigid_object_center.id) 100 return 101 102 def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False): 103 """ 104 Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles. 105 106 Args: 107 particle_id1 (`int`): particle_id of the type of the first particle type of the bonded particles 108 particle_id2 (`int`): particle_id of the type of the second particle type of the bonded particles 109 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. 110 """ 111 particle_name1 = self.df.loc[self.df['particle_id']==particle_id1].name.values[0] 112 particle_name2 = self.df.loc[self.df['particle_id']==particle_id2].name.values[0] 113 bond_key = self.find_bond_key(particle_name1=particle_name1, 114 particle_name2=particle_name2, 115 use_default_bond=use_default_bond) 116 if not bond_key: 117 return 118 self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1) 119 indexs = np.where(self.df['name']==bond_key) 120 index_list = list (indexs[0]) 121 used_bond_df = self.df.loc[self.df['particle_id2'].notnull()] 122 #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict 123 used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 ) 124 used_bond_index = used_bond_df.index.to_list() 125 126 for index in index_list: 127 if index not in used_bond_index: 128 self.clean_df_row(index=int(index)) 129 self.df.at[index,'particle_id'] = particle_id1 130 self.df.at[index,'particle_id2'] = particle_id2 131 break 132 return 133 134 def add_bonds_to_espresso(self, espresso_system) : 135 """ 136 Adds all bonds defined in `pmb.df` to `espresso_system`. 137 138 Args: 139 espresso_system (str): system object of espressomd library 140 """ 141 142 if 'bond' in self.df.values: 143 bond_df = self.df.loc[self.df ['pmb_type'] == 'bond'] 144 bond_list = bond_df.bond_object.values.tolist() 145 for bond in bond_list: 146 espresso_system.bonded_inter.add(bond) 147 else: 148 print ('WARNING: There are no bonds defined in pymbe.df') 149 150 return 151 152 def add_value_to_df(self,index,key,new_value, verbose=True, non_standard_value=False): 153 """ 154 Adds a value to a cell in the `pmb.df` DataFrame. 155 156 Args: 157 index(`int`): index of the row to add the value to. 158 key(`str`): the column label to add the value to. 159 verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 160 non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False. 161 """ 162 163 token = "#protected:" 164 165 def protect(obj): 166 if non_standard_value: 167 return token + json.dumps(obj, cls=self.NumpyEncoder) 168 return obj 169 170 def deprotect(obj): 171 if non_standard_value and isinstance(obj, str) and obj.startswith(token): 172 return json.loads(obj.removeprefix(token)) 173 return obj 174 175 # Make sure index is a scalar integer value 176 index = int (index) 177 assert isinstance(index, int), '`index` should be a scalar integer value.' 178 idx = pd.IndexSlice 179 if self.check_if_df_cell_has_a_value(index=index,key=key): 180 old_value= self.df.loc[index,idx[key]] 181 if verbose: 182 if protect(old_value) != protect(new_value): 183 name=self.df.loc[index,('name','')] 184 pmb_type=self.df.loc[index,('pmb_type','')] 185 print(f"WARNING: you are redefining the properties of {name} of pmb_type {pmb_type}") 186 print(f'WARNING: overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}') 187 self.df.loc[index,idx[key]] = protect(new_value) 188 if non_standard_value: 189 self.df[key] = self.df[key].apply(deprotect) 190 return 191 192 def assign_molecule_id(self, name, molecule_index, pmb_type, used_molecules_id): 193 """ 194 Assigns the `molecule_id` of the pmb object given by `pmb_type` 195 196 Args: 197 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 198 pmb_type(`str`): pmb_object_type to assign the `molecule_id` 199 molecule_index (`int`): index of the current `pmb_object_type` to assign the `molecule_id` 200 used_molecules_id (`lst`): list with the `molecule_id` values already used. 201 202 Returns: 203 molecule_id(`int`): Id of the molecule 204 """ 205 206 self.clean_df_row(index=int(molecule_index)) 207 208 if self.df['molecule_id'].isnull().values.all(): 209 molecule_id = 0 210 else: 211 # check if a residue is part of another molecule 212 check_residue_name = self.df[self.df['residue_list'].astype(str).str.contains(name)] 213 mol_pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 214 if not check_residue_name.empty and mol_pmb_type == pmb_type: 215 for value in check_residue_name.index.to_list(): 216 if value not in used_molecules_id: 217 molecule_id = self.df.loc[value].molecule_id.values[0] 218 break 219 else: 220 molecule_id = self.df['molecule_id'].max() +1 221 222 self.add_value_to_df (key=('molecule_id',''), 223 index=int(molecule_index), 224 new_value=molecule_id, 225 verbose=False) 226 227 return molecule_id 228 229 def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system): 230 """ 231 Calculates the center of the molecule with a given molecule_id. 232 233 Args: 234 molecule_id(`int`): Id of the molecule whose center of mass is to be calculated. 235 espresso_system(`obj`): Instance of a system object from the espressomd library. 236 237 Returns: 238 center_of_mass(`lst`): Coordinates of the center of mass. 239 """ 240 total_beads = 0 241 center_of_mass = np.zeros(3) 242 axis_list = [0,1,2] 243 particle_id_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 244 for pid in particle_id_list: 245 total_beads +=1 246 for axis in axis_list: 247 center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis] 248 center_of_mass = center_of_mass /total_beads 249 return center_of_mass 250 251 def calculate_HH(self, molecule_name, pH_list=None, pka_set=None): 252 """ 253 Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 254 for molecules with the name `molecule_name`. 255 256 Args: 257 molecule_name (`str`): name of the molecule to calculate the ideal charge for 258 pH_list(`lst`): pH-values to calculate. 259 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 260 261 Returns: 262 Z_HH (`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list` 263 264 Note: 265 - This function supports objects with pmb types: "molecule", "peptide" and "protein". 266 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 267 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 268 - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan` should be used. 269 """ 270 if pH_list is None: 271 pH_list=np.linspace(2,12,50) 272 if pka_set is None: 273 pka_set=self.get_pka_set() 274 self.check_pka_set(pka_set=pka_set) 275 charge_map = self.get_charge_map() 276 Z_HH=[] 277 for pH_value in pH_list: 278 Z=0 279 index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 280 residue_list = self.df.at [index,('residue_list','')] 281 sequence = self.df.at [index,('sequence','')] 282 if np.any(pd.isnull(sequence)): 283 # Molecule has no sequence 284 for residue in residue_list: 285 list_of_particles_in_residue = self.search_particles_in_residue(residue) 286 for particle in list_of_particles_in_residue: 287 if particle in pka_set.keys(): 288 if pka_set[particle]['acidity'] == 'acidic': 289 psi=-1 290 elif pka_set[particle]['acidity']== 'basic': 291 psi=+1 292 else: 293 psi=0 294 Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value']))) 295 Z_HH.append(Z) 296 else: 297 # Molecule has a sequence 298 if not isinstance(sequence, list): 299 # If the df has been read by file, the sequence needs to be parsed. 300 sequence = self.parse_sequence_from_file(sequence=sequence) 301 for name in sequence: 302 if name in pka_set.keys(): 303 if pka_set[name]['acidity'] == 'acidic': 304 psi=-1 305 elif pka_set[name]['acidity']== 'basic': 306 psi=+1 307 else: 308 psi=0 309 Z+=psi/(1+10**(psi*(pH_value-pka_set[name]['pka_value']))) 310 else: 311 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 312 Z+=charge_map[state_one_type] 313 Z_HH.append(Z) 314 315 return Z_HH 316 317 def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None): 318 """ 319 Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning. 320 321 Args: 322 c_macro ('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 323 c_salt ('float'): Salt concentration in the reservoir. 324 pH_list ('lst'): List of pH-values in the reservoir. 325 pka_set ('dict'): {"name": {"pka_value": pka, "acidity": acidity}}. 326 327 Returns: 328 {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 329 pH_system_list ('lst'): List of pH_values in the system. 330 partition_coefficients_list ('lst'): List of partition coefficients of cations. 331 332 Note: 333 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 334 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 335 - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result. 336 """ 337 if pH_list is None: 338 pH_list=np.linspace(2,12,50) 339 if pka_set is None: 340 pka_set=self.get_pka_set() 341 self.check_pka_set(pka_set=pka_set) 342 343 partition_coefficients_list = [] 344 pH_system_list = [] 345 Z_HH_Donnan={} 346 for key in c_macro: 347 Z_HH_Donnan[key] = [] 348 349 def calc_charges(c_macro, pH): 350 """ 351 Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation. 352 353 Args: 354 c_macro ('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 355 pH ('float'): pH-value that is used in the HH equation. 356 357 Returns: 358 charge ('dict'): {"molecule_name": charge} 359 """ 360 charge = {} 361 for name in c_macro: 362 charge[name] = self.calculate_HH(name, [pH], pka_set)[0] 363 return charge 364 365 def calc_partition_coefficient(charge, c_macro): 366 """ 367 Calculates the partition coefficients of positive ions according to the ideal Donnan theory. 368 369 Args: 370 charge ('dict'): {"molecule_name": charge} 371 c_macro ('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 372 """ 373 nonlocal ionic_strength_res 374 charge_density = 0.0 375 for key in charge: 376 charge_density += charge[key] * c_macro[key] 377 return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude 378 379 for pH_value in pH_list: 380 # calculate the ionic strength of the reservoir 381 if pH_value <= 7.0: 382 ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 383 elif pH_value > 7.0: 384 ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt 385 386 #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations 387 #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation. 388 #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects. 389 equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro)) 390 logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq") 391 partition_coefficient = 10**logxi.root 392 393 charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient)) 394 for key in c_macro: 395 Z_HH_Donnan[key].append(charges_temp[key]) 396 397 pH_system_list.append(pH_value - np.log10(partition_coefficient)) 398 partition_coefficients_list.append(partition_coefficient) 399 400 return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 401 402 def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff): 403 """ 404 Calculates the initial bond length that is used when setting up molecules, 405 based on the minimum of the sum of bonded and short-range (LJ) interactions. 406 407 Args: 408 bond_object (`cls`): instance of a bond object from espressomd library 409 bond_type (`str`): label identifying the used bonded potential 410 epsilon (`float`): LJ epsilon of the interaction between the particles 411 sigma (`float`): LJ sigma of the interaction between the particles 412 cutoff (`float`): cutoff-radius of the LJ interaction 413 """ 414 def truncated_lj_potential(x, epsilon, sigma, cutoff): 415 if x>cutoff: 416 return 0.0 417 else: 418 return 4*epsilon*((sigma/x)**12-(sigma/x)**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6) 419 420 epsilon_red=epsilon.to('reduced_energy').magnitude 421 sigma_red=sigma.to('reduced_length').magnitude 422 cutoff_red=cutoff.to('reduced_length').magnitude 423 424 if bond_type == "harmonic": 425 r_0 = bond_object.params.get('r_0') 426 k = bond_object.params.get('k') 427 l0 = scipy.optimize.minimize(lambda x: 0.5*k*(x-r_0)**2 + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red), x0=r_0).x 428 elif bond_type == "FENE": 429 r_0 = bond_object.params.get('r_0') 430 k = bond_object.params.get('k') 431 d_r_max = bond_object.params.get('d_r_max') 432 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), x0=1.0).x 433 return l0 434 435 def calculate_net_charge (self, espresso_system, molecule_name): 436 ''' 437 Calculates the net charge per molecule of molecules with `name` = molecule_name. 438 Returns the net charge per molecule and a maps with the net charge per residue and molecule. 439 440 Args: 441 espresso_system: system information 442 molecule_name (str): name of the molecule to calculate the net charge 443 444 Returns: 445 {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }} 446 447 Note: 448 - The net charge of the molecule is averaged over all molecules of type `name` 449 - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name` 450 ''' 451 valid_pmb_types = ["molecule", "protein"] 452 pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0] 453 if pmb_type not in valid_pmb_types: 454 raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}") 455 456 id_map = self.get_particle_id_map(object_name=molecule_name) 457 def create_charge_map(espresso_system,id_map,label): 458 charge_map={} 459 for super_id in id_map[label].keys(): 460 net_charge=0 461 for pid in id_map[label][super_id]: 462 net_charge+=espresso_system.part.by_id(pid).q 463 charge_map[super_id]=net_charge 464 return charge_map 465 net_charge_molecules=create_charge_map(label="molecule_map", 466 espresso_system=espresso_system, 467 id_map=id_map) 468 net_charge_residues=create_charge_map(label="residue_map", 469 espresso_system=espresso_system, 470 id_map=id_map) 471 mean_charge=np.mean(np.array(list(net_charge_molecules.values()))) 472 return {"mean": mean_charge, "molecules": net_charge_molecules, "residues": net_charge_residues} 473 474 def center_molecule_in_simulation_box(self, molecule_id, espresso_system): 475 """ 476 Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`. 477 478 Args: 479 molecule_id(`int`): Id of the molecule to be centered. 480 espresso_system(`obj`): Instance of a system object from the espressomd library. 481 """ 482 if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0: 483 raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.") 484 center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system) 485 box_center = [espresso_system.box_l[0]/2.0, 486 espresso_system.box_l[1]/2.0, 487 espresso_system.box_l[2]/2.0] 488 particle_id_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 489 for pid in particle_id_list: 490 es_pos = espresso_system.part.by_id(pid).pos 491 espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center 492 return 493 494 def check_dimensionality(self, variable, expected_dimensionality): 495 """ 496 Checks if the dimensionality of `variable` matches `expected_dimensionality`. 497 498 Args: 499 `variable`(`pint.Quantity`): Quantity to be checked. 500 `expected_dimensionality(`str`): Expected dimension of the variable. 501 502 Returns: 503 `bool`: `True` if the variable if of the expected dimensionality, `False` otherwise. 504 505 Note: 506 - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality). 507 - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"` 508 """ 509 correct_dimensionality=variable.check(f"{expected_dimensionality}") 510 if not correct_dimensionality: 511 raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}") 512 return correct_dimensionality 513 514 def check_if_df_cell_has_a_value(self, index,key): 515 """ 516 Checks if a cell in the `pmb.df` at the specified index and column has a value. 517 518 Args: 519 index(`int`): Index of the row to check. 520 key(`str`): Column label to check. 521 522 Returns: 523 `bool`: `True` if the cell has a value, `False` otherwise. 524 """ 525 idx = pd.IndexSlice 526 import warnings 527 with warnings.catch_warnings(): 528 warnings.simplefilter("ignore") 529 return not pd.isna(self.df.loc[index, idx[key]]) 530 531 def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined): 532 """ 533 Checks if `name` is defined in `pmb.df`. 534 535 Args: 536 name(`str`): label to check if defined in `pmb.df`. 537 pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`. 538 539 Returns: 540 `bool`: `True` for success, `False` otherwise. 541 """ 542 if name in self.df['name'].unique(): 543 current_object_type = self.df[self.df['name']==name].pmb_type.values[0] 544 if current_object_type != pmb_type_to_be_defined: 545 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") 546 return True 547 else: 548 return False 549 550 551 def check_pka_set(self, pka_set): 552 """ 553 Checks that `pka_set` has the formatting expected by the pyMBE library. 554 555 Args: 556 pka_set (dict): {"name" : {"pka_value": pka, "acidity": acidity}} 557 """ 558 required_keys=['pka_value','acidity'] 559 for required_key in required_keys: 560 for pka_entry in pka_set.values(): 561 if required_key not in pka_entry.keys(): 562 raise ValueError(f'missing a required key "{required_key}" in the following entry of pka_set "{pka_entry}"') 563 return 564 565 def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")): 566 """ 567 Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a np.nan value. 568 569 Args: 570 index(`int`): Index of the row to clean. 571 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`]. 572 """ 573 for column_key in columns_keys_to_clean: 574 self.add_value_to_df(key=(column_key,''),index=index,new_value=np.nan, verbose=False) 575 return 576 577 def convert_columns_to_original_format (self, df): 578 """ 579 Converts the columns of the Dataframe to the original format in pyMBE. 580 581 Args: 582 df(`DataFrame`): dataframe with pyMBE information as a string 583 584 """ 585 586 columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', 'model',('state_one','es_type'),('state_two','es_type'),('state_one','charge'),('state_two','charge') ] 587 588 columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset'] 589 590 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence'] 591 592 for column_name in columns_dtype_int: 593 df[column_name] = df[column_name].astype(object) 594 595 for column_name in columns_with_list_or_dict: 596 if df[column_name].isnull().all(): 597 df[column_name] = df[column_name].astype(object) 598 else: 599 df[column_name] = df[column_name].apply(lambda x: ast.literal_eval(str(x)) if pd.notnull(x) else x) 600 601 for column_name in columns_with_units: 602 df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x) 603 604 df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x) 605 606 return df 607 608 def convert_str_to_bond_object (self, bond_str): 609 610 """ 611 Convert a row read as a `str` to the corresponding bond object. There are two supported bonds: HarmonicBond and FeneBond 612 613 Args: 614 bond_str (`str`): string with the information of a bond object 615 616 Returns: 617 bond_object(`obj`): EsPRESSo bond object 618 """ 619 620 from espressomd.interactions import HarmonicBond 621 from espressomd.interactions import FeneBond 622 623 supported_bonds = ['HarmonicBond', 'FeneBond'] 624 625 for bond in supported_bonds: 626 variable = re.subn(f'{bond}', '', bond_str) 627 if variable[1] == 1: 628 params = ast.literal_eval(variable[0]) 629 if bond == 'HarmonicBond': 630 bond_object = HarmonicBond(r_cut =params['r_cut'], k = params['k'], r_0=params['r_0']) 631 elif bond == 'FeneBond': 632 bond_object = FeneBond(k = params['k'], d_r_max =params['d_r_max'], r_0=params['r_0']) 633 634 return bond_object 635 636 def copy_df_entry(self, name, column_name, number_of_copies): 637 ''' 638 Creates 'number_of_copies' of a given 'name' in `pymbe.df`. 639 640 Args: 641 name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df` 642 column_name(`str`): Column name to use as a filter. 643 number_of_copies(`int`): number of copies of `name` to be created. 644 645 Note: 646 - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 647 ''' 648 649 valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ] 650 if column_name not in valid_column_names: 651 raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}") 652 df_by_name = self.df.loc[self.df.name == name] 653 if number_of_copies != 1: 654 if df_by_name[column_name].isnull().values.any(): 655 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True) 656 else: 657 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 658 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 659 df_by_name_repeated[column_name] =np.NaN 660 # Concatenate the new particle rows to `df` 661 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 662 else: 663 if not df_by_name[column_name].isnull().values.any(): 664 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 665 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 666 df_by_name_repeated[column_name] =np.NaN 667 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 668 return 669 670 def create_added_salt (self, espresso_system, cation_name, anion_name, c_salt, verbose=True): 671 """ 672 Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`. 673 674 Args: 675 espresso_system (`obj`): instance of an espresso system object. 676 cation_name(`str`): `name` of a particle with a positive charge. 677 anion_name(`str`): `name` of a particle with a negative charge. 678 c_salt(`float`): Salt concentration. 679 verbose(`bool`): switch to activate/deactivate verbose. Defaults to True. 680 681 Returns: 682 c_salt_calculated (float): Calculated salt concentration added to `espresso_system`. 683 """ 684 cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.values[0] 685 anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.values[0] 686 if cation_name_charge <= 0: 687 raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge) 688 if anion_name_charge >= 0: 689 raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge) 690 # Calculate the number of ions in the simulation box 691 volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3') 692 if c_salt.check('[substance] [length]**-3'): 693 N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude) 694 c_salt_calculated=N_ions/(volume*self.N_A) 695 elif c_salt.check('[length]**-3'): 696 N_ions= int((volume*c_salt.to('reduced_length**-3')*self.N_A).magnitude) 697 c_salt_calculated=N_ions/volume 698 else: 699 raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt) 700 N_cation = N_ions*abs(anion_name_charge) 701 N_anion = N_ions*abs(cation_name_charge) 702 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation) 703 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion) 704 if verbose: 705 print(f"\n Added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions") 706 return c_salt_calculated 707 708 def create_bond_in_espresso(self, bond_type, bond_parameters): 709 ''' 710 Creates either a harmonic or a FENE bond in ESPREesSo 711 712 Args: 713 bond_type (`str`) : label to identify the potential to model the bond. 714 bond_parameters (`dict`) : parameters of the potential of the bond. 715 716 Note: 717 Currently, only HARMONIC and FENE bonds are supported. 718 719 For a HARMONIC bond the dictionary must contain: 720 721 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 722 using the `pmb.units` UnitRegistry. 723 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 724 the `pmb.units` UnitRegistry. 725 726 For a FENE bond the dictionay must additionally contain: 727 728 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 729 units of length using the `pmb.units` UnitRegistry. Default 'None'. 730 731 Returns: 732 bond_object (`obj`): a harmonic or a FENE bond object in ESPREesSo 733 ''' 734 from espressomd import interactions 735 736 valid_bond_types = ["harmonic", "FENE"] 737 738 if 'k' in bond_parameters: 739 bond_magnitude = bond_parameters['k'].to('reduced_energy / reduced_length**2') 740 else: 741 raise ValueError("Magnitud of the potential (k) is missing") 742 743 if bond_type == 'harmonic': 744 if 'r_0' in bond_parameters: 745 bond_length = bond_parameters['r_0'].to('reduced_length') 746 else: 747 raise ValueError("Equilibrium bond length (r_0) is missing") 748 bond_object = interactions.HarmonicBond(k = bond_magnitude.magnitude, 749 r_0 = bond_length.magnitude) 750 elif bond_type == 'FENE': 751 if 'r_0' in bond_parameters: 752 bond_length = bond_parameters['r_0'].to('reduced_length').magnitude 753 else: 754 print("WARNING: No value provided for r_0. Defaulting to r_0 = 0") 755 bond_length=0 756 if 'd_r_max' in bond_parameters: 757 max_bond_stret = bond_parameters['d_r_max'].to('reduced_length') 758 else: 759 raise ValueError("Maximal stretching length (d_r_max) is missing") 760 bond_object = interactions.FeneBond(r_0 = bond_length, 761 k = bond_magnitude.magnitude, 762 d_r_max = max_bond_stret.magnitude) 763 else: 764 raise ValueError(f"Bond type {bond_type} currently not implemented in pyMBE, valid types are {valid_bond_types}") 765 766 return bond_object 767 768 769 def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True): 770 """ 771 Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`. 772 773 Args: 774 object_name(`str`): `name` of a pymbe object. 775 espresso_system(`obj`): Instance of a system object from the espressomd library. 776 cation_name(`str`): `name` of a particle with a positive charge. 777 anion_name(`str`): `name` of a particle with a negative charge. 778 verbose(`bool`): switch to activate/deactivate verbose. Defaults to True. 779 780 Returns: 781 counterion_number (dict): {"name": number} 782 """ 783 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.iloc[0] 784 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.iloc[0] 785 object_ids = self.get_particle_id_map(object_name=object_name)["all"] 786 counterion_number={} 787 object_charge={} 788 for name in ['positive', 'negative']: 789 object_charge[name]=0 790 for id in object_ids: 791 if espresso_system.part.by_id(id).q > 0: 792 object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q )) 793 elif espresso_system.part.by_id(id).q < 0: 794 object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q )) 795 if object_charge['positive'] % abs(anion_charge) == 0: 796 counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge)) 797 else: 798 raise ValueError('The number of positive charges in the pmb_object must be divisible by the charge of the anion') 799 if object_charge['negative'] % abs(cation_charge) == 0: 800 counterion_number[cation_name]=int(object_charge['negative']/cation_charge) 801 else: 802 raise ValueError('The number of negative charges in the pmb_object must be divisible by the charge of the cation') 803 if counterion_number[cation_name] > 0: 804 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name]) 805 else: 806 counterion_number[cation_name]=0 807 if counterion_number[anion_name] > 0: 808 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name]) 809 else: 810 counterion_number[anion_name] = 0 811 if verbose: 812 print('The following counter-ions have been created: ') 813 for name in counterion_number.keys(): 814 print(f'Ion type: {name} created number: {counterion_number[name]}') 815 return counterion_number 816 817 def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, use_default_bond=False): 818 """ 819 Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 820 821 Args: 822 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 823 espresso_system(`obj`): Instance of a system object from espressomd library. 824 number_of_molecules(`int`): Number of molecules of type `name` to be created. 825 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 826 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df` 827 Returns: 828 829 molecules_info (`dict`): {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 830 """ 831 if list_of_first_residue_positions is not None: 832 for item in list_of_first_residue_positions: 833 if not isinstance(item, list): 834 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.") 835 elif len(item) != 3: 836 raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.") 837 838 if len(list_of_first_residue_positions) != number_of_molecules: 839 raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}") 840 if number_of_molecules <= 0: 841 return 842 if not self.check_if_name_is_defined_in_df(name=name, 843 pmb_type_to_be_defined='molecule'): 844 raise ValueError(f"{name} must correspond to a label of a pmb_type='molecule' defined on df") 845 first_residue = True 846 molecules_info = {} 847 residue_list = self.df[self.df['name']==name].residue_list.values [0] 848 849 self.copy_df_entry(name=name, 850 column_name='molecule_id', 851 number_of_copies=number_of_molecules) 852 853 molecules_index = np.where(self.df['name']==name) 854 molecule_index_list =list(molecules_index[0])[-number_of_molecules:] 855 used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist() 856 pos_index = 0 857 for molecule_index in molecule_index_list: 858 molecule_id = self.assign_molecule_id(name=name,pmb_type='molecule',used_molecules_id=used_molecules_id,molecule_index=molecule_index) 859 molecules_info[molecule_id] = {} 860 for residue in residue_list: 861 if first_residue: 862 if list_of_first_residue_positions is None: 863 residue_position = None 864 else: 865 for item in list_of_first_residue_positions: 866 residue_position = [np.array(list_of_first_residue_positions[pos_index])] 867 # Generate an arbitrary random unit vector 868 backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0], 869 radius=1, 870 n_samples=1, 871 on_surface=True)[0] 872 residues_info = self.create_residue(name=residue, 873 number_of_residues=1, 874 espresso_system=espresso_system, 875 central_bead_position=residue_position, 876 use_default_bond= use_default_bond, 877 backbone_vector=backbone_vector) 878 residue_id = next(iter(residues_info)) 879 for index in self.df[self.df['residue_id']==residue_id].index: 880 self.add_value_to_df(key=('molecule_id',''), 881 index=int (index), 882 new_value=molecule_id) 883 central_bead_id = residues_info[residue_id]['central_bead_id'] 884 previous_residue = residue 885 residue_position = espresso_system.part.by_id(central_bead_id).pos 886 previous_residue_id = central_bead_id 887 first_residue = False 888 else: 889 previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0] 890 new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0] 891 bond = self.search_bond(particle_name1=previous_central_bead_name, 892 particle_name2=new_central_bead_name, 893 hard_check=True, 894 use_default_bond=use_default_bond) 895 l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 896 particle_name2=new_central_bead_name, 897 hard_check=True, 898 use_default_bond=use_default_bond) 899 residue_position = residue_position+backbone_vector*l0 900 residues_info = self.create_residue(name=residue, 901 number_of_residues=1, 902 espresso_system=espresso_system, 903 central_bead_position=[residue_position], 904 use_default_bond= use_default_bond, 905 backbone_vector=backbone_vector) 906 residue_id = next(iter(residues_info)) 907 for index in self.df[self.df['residue_id']==residue_id].index: 908 if not self.check_if_df_cell_has_a_value(index=index,key=('molecule_id','')): 909 self.df.at[index,'molecule_id'] = molecule_id 910 self.add_value_to_df(key=('molecule_id',''), 911 index=int (index), 912 new_value=molecule_id, 913 verbose=False) 914 central_bead_id = residues_info[residue_id]['central_bead_id'] 915 espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id)) 916 self.add_bond_in_df(particle_id1=central_bead_id, 917 particle_id2=previous_residue_id, 918 use_default_bond=use_default_bond) 919 previous_residue_id = central_bead_id 920 previous_residue = residue 921 molecules_info[molecule_id][residue_id] = residues_info[residue_id] 922 first_residue = True 923 pos_index+=1 924 925 return molecules_info 926 927 def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False): 928 """ 929 Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`. 930 931 Args: 932 name (`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 933 espresso_system (`cls`): Instance of a system object from the espressomd library. 934 number_of_particles (`int`): Number of particles to be created. 935 position (list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None. 936 fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False. 937 Returns: 938 created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`. 939 """ 940 if number_of_particles <=0: 941 return 942 if not self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'): 943 raise ValueError(f"{name} must correspond to a label of a pmb_type='particle' defined on df") 944 # Copy the data of the particle `number_of_particles` times in the `df` 945 self.copy_df_entry(name=name,column_name='particle_id',number_of_copies=number_of_particles) 946 # Get information from the particle type `name` from the df 947 q = self.df.loc[self.df['name']==name].state_one.charge.values[0] 948 es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 949 # Get a list of the index in `df` corresponding to the new particles to be created 950 index = np.where(self.df['name']==name) 951 index_list =list(index[0])[-number_of_particles:] 952 # Create the new particles into `espresso_system` 953 created_pid_list=[] 954 for index in range (number_of_particles): 955 df_index=int (index_list[index]) 956 self.clean_df_row(index=df_index) 957 if position is None: 958 particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l) 959 else: 960 particle_position = position[index] 961 if len(espresso_system.part.all()) == 0: 962 bead_id = 0 963 else: 964 bead_id = max (espresso_system.part.all().id) + 1 965 created_pid_list.append(bead_id) 966 967 if fix: 968 espresso_system.part.add (id=bead_id, pos = particle_position, type = es_type, q = q,fix =[fix,fix,fix]) 969 else: 970 espresso_system.part.add (id=bead_id, pos = particle_position, type = es_type, q = q) 971 self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id, verbose=False) 972 return created_pid_list 973 974 def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False): 975 """ 976 Creates all `particle`s associated to `pmb object` into `espresso` a number of times equal to `number_of_objects`. 977 978 Args: 979 name(`str`): Unique label of the `pmb object` to be created. 980 number_of_objects(`int`): Number of `pmb object`s to be created. 981 espresso_system(`obj`): Instance of an espresso system object from espressomd library. 982 position(`list`): Coordinates where the particles should be created. 983 use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`. 984 985 Note: 986 - 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. 987 """ 988 allowed_objects=['particle','residue','molecule'] 989 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 990 if pmb_type not in allowed_objects: 991 raise ValueError('Object type not supported, supported types are ', allowed_objects) 992 if pmb_type == 'particle': 993 self.create_particle(name=name, number_of_particles=number_of_objects, espresso_system=espresso_system, position=position) 994 elif pmb_type == 'residue': 995 self.create_residue(name=name, number_of_residues=number_of_objects, espresso_system=espresso_system, central_bead_position=position,use_default_bond=use_default_bond) 996 elif pmb_type == 'molecule': 997 self.create_molecule(name=name, number_of_molecules=number_of_objects, espresso_system=espresso_system, use_default_bond=use_default_bond, list_of_first_residue_positions=position) 998 return 999 1000 def create_protein(self, name, number_of_proteins, espresso_system, topology_dict): 1001 """ 1002 Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions` 1003 1004 Args: 1005 name(`str`): Label of the protein to be created. 1006 espresso_system(`obj`): Instance of a system object from the espressomd library. 1007 number_of_proteins(`int`): Number of proteins to be created. 1008 positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}} 1009 """ 1010 1011 if number_of_proteins <=0: 1012 return 1013 if not self.check_if_name_is_defined_in_df(name=name, 1014 pmb_type_to_be_defined='protein'): 1015 raise ValueError(f"{name} must correspond to a name of a pmb_type='protein' defined on df") 1016 1017 self.copy_df_entry(name=name,column_name='molecule_id',number_of_copies=number_of_proteins) 1018 1019 protein_index = np.where(self.df['name']==name) 1020 protein_index_list =list(protein_index[0])[-number_of_proteins:] 1021 used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist() 1022 1023 box_half=espresso_system.box_l[0]/2.0 1024 1025 for molecule_index in protein_index_list: 1026 1027 molecule_id = self.assign_molecule_id (name=name,pmb_type='protein',used_molecules_id=used_molecules_id,molecule_index=molecule_index) 1028 1029 protein_center = self.generate_coordinates_outside_sphere(radius = 1, 1030 max_dist=box_half, 1031 n_samples=1, 1032 center=[box_half]*3)[0] 1033 1034 for residue in topology_dict.keys(): 1035 1036 residue_name = re.split(r'\d+', residue)[0] 1037 residue_number = re.split(r'(\d+)', residue)[1] 1038 residue_position = topology_dict[residue]['initial_pos'] 1039 position = residue_position + protein_center 1040 1041 particle_id = self.create_particle(name=residue_name, 1042 espresso_system=espresso_system, 1043 number_of_particles=1, 1044 position=[position], 1045 fix = True) 1046 1047 index = self.df[self.df['particle_id']==particle_id[0]].index.values[0] 1048 1049 self.add_value_to_df(key=('residue_id',''), 1050 index=int (index), 1051 new_value=residue_number) 1052 1053 self.add_value_to_df(key=('molecule_id',''), 1054 index=int (index), 1055 new_value=molecule_id) 1056 1057 return 1058 1059 def create_residue(self, name, espresso_system, number_of_residues, central_bead_position=None,use_default_bond=False, backbone_vector=None): 1060 """ 1061 Creates `number_of_residues` residues of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1062 1063 Args: 1064 name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df` 1065 espresso_system(`obj`): Instance of a system object from espressomd library. 1066 number_of_residue(`int`): Number of residues of type `name` to be created. 1067 central_bead_position(`list` of `float`): Position of the central bead. 1068 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` 1069 backbone_vector (`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`. 1070 1071 Returns: 1072 residues_info (`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1073 """ 1074 if number_of_residues <= 0: 1075 return 1076 if not self.check_if_name_is_defined_in_df(name=name, 1077 pmb_type_to_be_defined='residue'): 1078 raise ValueError(f"{name} must correspond to a label of a pmb_type='residue' defined on df") 1079 # Copy the data of a residue a `number_of_residues` times in the `df 1080 self.copy_df_entry(name=name, 1081 column_name='residue_id', 1082 number_of_copies=number_of_residues) 1083 residues_index = np.where(self.df['name']==name) 1084 residue_index_list =list(residues_index[0])[-number_of_residues:] 1085 # Internal bookkepping of the residue info (important for side-chain residues) 1086 # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1087 residues_info={} 1088 for residue_index in residue_index_list: 1089 self.clean_df_row(index=int(residue_index)) 1090 # Assign a residue_id 1091 if self.df['residue_id'].isnull().all(): 1092 residue_id=0 1093 else: 1094 residue_id = self.df['residue_id'].max() + 1 1095 self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id, verbose=False) 1096 # create the principal bead 1097 if self.df.loc[self.df['name']==name].central_bead.values[0] is np.NaN: 1098 raise ValueError("central_bead must contain a particle name") 1099 central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0] 1100 central_bead_id = self.create_particle(name=central_bead_name, 1101 espresso_system=espresso_system, 1102 position=central_bead_position, 1103 number_of_particles = 1)[0] 1104 central_bead_position=espresso_system.part.by_id(central_bead_id).pos 1105 #assigns same residue_id to the central_bead particle created. 1106 index = self.df[self.df['particle_id']==central_bead_id].index.values[0] 1107 self.df.at [index,'residue_id'] = residue_id 1108 # Internal bookkeeping of the central bead id 1109 residues_info[residue_id]={} 1110 residues_info[residue_id]['central_bead_id']=central_bead_id 1111 # create the lateral beads 1112 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1113 side_chain_beads_ids = [] 1114 for side_chain_element in side_chain_list: 1115 if side_chain_element not in self.df.values: 1116 raise ValueError (side_chain_element +'is not defined') 1117 pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 1118 if pmb_type == 'particle': 1119 bond = self.search_bond(particle_name1=central_bead_name, 1120 particle_name2=side_chain_element, 1121 hard_check=True, 1122 use_default_bond=use_default_bond) 1123 l0 = self.get_bond_length(particle_name1=central_bead_name, 1124 particle_name2=side_chain_element, 1125 hard_check=True, 1126 use_default_bond=use_default_bond) 1127 1128 if backbone_vector is None: 1129 bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1130 radius=l0, 1131 n_samples=1, 1132 on_surface=True)[0] 1133 else: 1134 bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1135 magnitude=l0) 1136 1137 side_bead_id = self.create_particle(name=side_chain_element, 1138 espresso_system=espresso_system, 1139 position=[bead_position], 1140 number_of_particles=1)[0] 1141 index = self.df[self.df['particle_id']==side_bead_id].index.values[0] 1142 self.add_value_to_df(key=('residue_id',''), 1143 index=int (index), 1144 new_value=residue_id, 1145 verbose=False) 1146 side_chain_beads_ids.append(side_bead_id) 1147 espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id)) 1148 self.add_bond_in_df(particle_id1=central_bead_id, 1149 particle_id2=side_bead_id, 1150 use_default_bond=use_default_bond) 1151 elif pmb_type == 'residue': 1152 central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0] 1153 bond = self.search_bond(particle_name1=central_bead_name, 1154 particle_name2=central_bead_side_chain, 1155 hard_check=True, 1156 use_default_bond=use_default_bond) 1157 l0 = self.get_bond_length(particle_name1=central_bead_name, 1158 particle_name2=central_bead_side_chain, 1159 hard_check=True, 1160 use_default_bond=use_default_bond) 1161 if backbone_vector is None: 1162 residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1163 radius=l0, 1164 n_samples=1, 1165 on_surface=True)[0] 1166 else: 1167 residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1168 magnitude=l0) 1169 lateral_residue_info = self.create_residue(name=side_chain_element, espresso_system=espresso_system, 1170 number_of_residues=1, central_bead_position=[residue_position],use_default_bond=use_default_bond) 1171 lateral_residue_dict=list(lateral_residue_info.values())[0] 1172 central_bead_side_chain_id=lateral_residue_dict['central_bead_id'] 1173 lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids'] 1174 residue_id_side_chain=list(lateral_residue_info.keys())[0] 1175 # Change the residue_id of the residue in the side chain to the one of the biger residue 1176 index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0] 1177 self.add_value_to_df(key=('residue_id',''), 1178 index=int(index), 1179 new_value=residue_id, 1180 verbose=False) 1181 # Change the residue_id of the particles in the residue in the side chain 1182 side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids 1183 for particle_id in side_chain_beads_ids: 1184 index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0] 1185 self.add_value_to_df(key=('residue_id',''), 1186 index=int (index), 1187 new_value=residue_id, 1188 verbose=False) 1189 espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id)) 1190 self.add_bond_in_df(particle_id1=central_bead_id, 1191 particle_id2=central_bead_side_chain_id, 1192 use_default_bond=use_default_bond) 1193 # Internal bookkeeping of the side chain beads ids 1194 residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids 1195 return residues_info 1196 1197 def create_variable_with_units(self, variable): 1198 """ 1199 Returns a pint object with the value and units defined in `variable`. 1200 1201 Args: 1202 variable(`dict` or `str`): {'value': value, 'units': units} 1203 Returns: 1204 variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry. 1205 """ 1206 1207 if isinstance(variable, dict): 1208 1209 value=variable.pop('value') 1210 units=variable.pop('units') 1211 1212 elif isinstance(variable, str): 1213 1214 value = float(re.split(r'\s+', variable)[0]) 1215 units = re.split(r'\s+', variable)[1] 1216 1217 variable_with_units=value*self.units(units) 1218 1219 return variable_with_units 1220 1221 def define_AA_particles_in_sequence(self, sequence, sigma_dict=None): 1222 ''' 1223 Defines in `pmb.df` all the different particles in `sequence`. 1224 1225 Args: 1226 sequence(`lst`): Sequence of the peptide or protein. 1227 1228 Note: 1229 - It assumes that the names in `sequence` correspond to amino acid names using the standard one letter code. 1230 ''' 1231 1232 already_defined_AA=[] 1233 1234 for residue_name in sequence: 1235 if residue_name in already_defined_AA: 1236 continue 1237 self.define_particle (name=residue_name, q=0) 1238 1239 if sigma_dict: 1240 self.define_particles_parameter_from_dict(param_dict = sigma_dict, 1241 param_name = 'sigma') 1242 return 1243 1244 def define_AA_residues(self, sequence, model): 1245 """ 1246 Defines in `pmb.df` all the different residues in `sequence`. 1247 1248 Args: 1249 sequence(`lst`): Sequence of the peptide or protein. 1250 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1251 1252 Returns: 1253 residue_list (`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1254 """ 1255 residue_list = [] 1256 for residue_name in sequence: 1257 if model == '1beadAA': 1258 central_bead = residue_name 1259 side_chains = [] 1260 elif model == '2beadAA': 1261 if residue_name in ['c','n', 'G']: 1262 central_bead = residue_name 1263 side_chains = [] 1264 else: 1265 central_bead = 'CA' 1266 side_chains = [residue_name] 1267 if residue_name not in residue_list: 1268 self.define_residue(name = 'AA-'+residue_name, 1269 central_bead = central_bead, 1270 side_chains = side_chains) 1271 residue_list.append('AA-'+residue_name) 1272 return residue_list 1273 1274 def define_bond(self, bond_type, bond_parameters, particle_pairs): 1275 1276 ''' 1277 Defines a pmb object of type `bond` in `pymbe.df`. 1278 1279 Args: 1280 bond_type (`str`) : label to identify the potential to model the bond. 1281 bond_parameters (`dict`) : parameters of the potential of the bond. 1282 particle_pairs (`lst`) : list of the `names` of the `particles` to be bonded. 1283 1284 Note: 1285 Currently, only HARMONIC and FENE bonds are supported. 1286 1287 For a HARMONIC bond the dictionary must contain the following parameters: 1288 1289 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 1290 using the `pmb.units` UnitRegistry. 1291 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 1292 the `pmb.units` UnitRegistry. 1293 1294 For a FENE bond the dictionay must contain the same parameters as for a HARMONIC bond and: 1295 1296 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 1297 units of length using the `pmb.units` UnitRegistry. Default 'None'. 1298 ''' 1299 1300 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1301 1302 for particle_name1, particle_name2 in particle_pairs: 1303 1304 lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1, 1305 particle_name2 = particle_name2, 1306 combining_rule = 'Lorentz-Berthelot') 1307 1308 cutoff = 2**(1.0/6.0) * lj_parameters["sigma"] 1309 1310 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1311 bond_type = bond_type, 1312 epsilon = lj_parameters["epsilon"], 1313 sigma = lj_parameters["sigma"], 1314 cutoff = cutoff) 1315 index = len(self.df) 1316 for label in [particle_name1+'-'+particle_name2,particle_name2+'-'+particle_name1]: 1317 self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond") 1318 self.df.at [index,'name']= particle_name1+'-'+particle_name2 1319 self.df.at [index,'bond_object'] = bond_object 1320 self.df.at [index,'l0'] = l0 1321 self.add_value_to_df(index=index, 1322 key=('pmb_type',''), 1323 new_value='bond') 1324 self.add_value_to_df(index=index, 1325 key=('parameters_of_the_potential',''), 1326 new_value=bond_object.get_params(), 1327 non_standard_value=True) 1328 1329 return 1330 1331 1332 def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None): 1333 """ 1334 Asigns `bond` in `pmb.df` as the default bond. 1335 The LJ parameters can be optionally provided to calculate the initial bond length 1336 1337 Args: 1338 bond_type (`str`) : label to identify the potential to model the bond. 1339 bond_parameters (`dict`) : parameters of the potential of the bond. 1340 sigma (`float`, optional) : LJ sigma of the interaction between the particles 1341 epsilon (`float`, optional): LJ epsilon for the interaction between the particles 1342 cutoff (`float`, optional) : cutoff-radius of the LJ interaction 1343 1344 Note: 1345 - Currently, only harmonic and FENE bonds are supported. 1346 """ 1347 1348 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1349 1350 if epsilon is None: 1351 epsilon=1*self.units('reduced_energy') 1352 if sigma is None: 1353 sigma=1*self.units('reduced_length') 1354 if cutoff is None: 1355 cutoff=2**(1.0/6.0)*self.units('reduced_length') 1356 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1357 bond_type = bond_type, 1358 epsilon = epsilon, 1359 sigma = sigma, 1360 cutoff = cutoff) 1361 1362 if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'): 1363 return 1364 if len(self.df.index) != 0: 1365 index = max(self.df.index)+1 1366 else: 1367 index = 0 1368 self.df.at [index,'name'] = 'default' 1369 self.df.at [index,'bond_object'] = bond_object 1370 self.df.at [index,'l0'] = l0 1371 self.add_value_to_df(index = index, 1372 key = ('pmb_type',''), 1373 new_value = 'bond') 1374 self.add_value_to_df(index = index, 1375 key = ('parameters_of_the_potential',''), 1376 new_value = bond_object.get_params(), 1377 non_standard_value=True) 1378 return 1379 1380 def define_epsilon_value_of_particles(self, eps_dict): 1381 ''' 1382 Defines the epsilon value of the particles in `eps_dict` into the `pmb.df`. 1383 1384 Args: 1385 eps_dict(`dict`): {'name': epsilon} 1386 ''' 1387 for residue in eps_dict.keys(): 1388 label_list = self.df[self.df['name'] == residue].index.tolist() 1389 for index in label_list: 1390 epsilon = eps_dict[residue] 1391 self.add_value_to_df(key= ('epsilon',''), 1392 index=int (index), 1393 new_value=epsilon) 1394 return 1395 1396 def define_molecule(self, name, residue_list): 1397 """ 1398 Defines a pyMBE object of type `molecule` in `pymbe.df`. 1399 1400 Args: 1401 name (`str`): Unique label that identifies the `molecule`. 1402 residue_list (`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1403 """ 1404 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'): 1405 return 1406 index = len(self.df) 1407 self.df.at [index,'name'] = name 1408 self.df.at [index,'pmb_type'] = 'molecule' 1409 self.df.at [index,('residue_list','')] = residue_list 1410 return 1411 1412 def define_particle(self, name, q=0, acidity='inert', pka=None, sigma=None, epsilon=None, cutoff=None, offset=None,verbose=True): 1413 """ 1414 Defines the properties of a particle object. 1415 1416 Args: 1417 name (`str`): Unique label that identifies this particle type. 1418 q (`int`, optional): Permanent charge of this particle type. Defaults to 0. 1419 acidity (`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to `inert`. 1420 pka (`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to None. 1421 sigma (`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1422 cutoff (`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1423 offset (`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1424 epsilon (`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to None. 1425 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 1426 1427 Note: 1428 - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units. 1429 - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units. 1430 - `cutoff` defaults to `2**(1./6.) reduced_length`. 1431 - `offset` defaults to 0. 1432 - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions. 1433 - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`. 1434 """ 1435 1436 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'): 1437 index = self.df[self.df['name']==name].index[0] 1438 else: 1439 index = len(self.df) 1440 self.df.at [index, 'name'] = name 1441 self.df.at [index,'pmb_type'] = 'particle' 1442 1443 # If `cutoff` and `offset` are not defined, default them to 0 1444 1445 if cutoff is None: 1446 cutoff=self.units.Quantity(2**(1./6.), "reduced_length") 1447 if offset is None: 1448 offset=self.units.Quantity(0, "reduced_length") 1449 1450 # Define LJ parameters 1451 parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"}, 1452 "cutoff":{"value": cutoff, "dimensionality": "[length]"}, 1453 "offset":{"value": offset, "dimensionality": "[length]"}, 1454 "epsilon":{"value": epsilon, "dimensionality": "[energy]"},} 1455 1456 for parameter_key in parameters_with_dimensionality.keys(): 1457 if parameters_with_dimensionality[parameter_key]["value"] is not None: 1458 self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 1459 expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"]) 1460 self.add_value_to_df(key=(parameter_key,''), 1461 index=index, 1462 new_value=parameters_with_dimensionality[parameter_key]["value"], 1463 verbose=verbose) 1464 1465 # Define particle acid/base properties 1466 self.set_particle_acidity (name=name, 1467 acidity = acidity , 1468 default_charge=q, 1469 pka=pka, 1470 verbose=verbose) 1471 return 1472 1473 def define_particles_parameter_from_dict(self, param_dict, param_name): 1474 ''' 1475 Defines the `param_name` value of the particles in `param_dict` into the `pmb.df`. 1476 1477 Args: 1478 param_dict(`dict`): {'name': `param_name` value} 1479 ''' 1480 for residue in param_dict.keys(): 1481 label_list = self.df[self.df['name'] == residue].index.tolist() 1482 for index in label_list: 1483 value = param_dict[residue] 1484 self.add_value_to_df(key= (param_name,''), 1485 index=int (index), 1486 new_value=value) 1487 return 1488 1489 def define_peptide(self, name, sequence, model): 1490 """ 1491 Defines a pyMBE object of type `peptide` in the `pymbe.df`. 1492 1493 Args: 1494 name (`str`): Unique label that identifies the `peptide`. 1495 sequence (`string`): Sequence of the `peptide`. 1496 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1497 """ 1498 if not self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'): 1499 valid_keys = ['1beadAA','2beadAA'] 1500 if model not in valid_keys: 1501 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1502 1503 clean_sequence = self.protein_sequence_parser(sequence=sequence) 1504 residue_list = self.define_AA_residues(sequence=clean_sequence,model=model) 1505 1506 self.define_molecule(name = name, residue_list=residue_list) 1507 index = self.df.loc[self.df['name'] == name].index.item() 1508 self.df.at [index,'model'] = model 1509 self.df.at [index,('sequence','')] = clean_sequence 1510 return 1511 1512 def define_protein(self, name,model, topology_dict): 1513 """ 1514 Defines a pyMBE object of type `protein` in `pymbe.df`. 1515 1516 Args: 1517 name (`str`): Unique label that identifies the `protein`. 1518 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1519 topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 1520 """ 1521 1522 protein_seq_list = [] 1523 1524 valid_keys = ['1beadAA','2beadAA'] 1525 1526 if model not in valid_keys: 1527 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1528 1529 if model == '2beadAA': 1530 self.define_particle(name='CA') 1531 1532 sigma_dict = {} 1533 1534 for residue in topology_dict.keys(): 1535 residue_name = re.split(r'\d+', residue)[0] 1536 1537 if residue_name not in sigma_dict.keys(): 1538 sigma_dict [residue_name] = topology_dict[residue]['sigma'] 1539 if residue_name not in ('CA', 'Ca'): 1540 protein_seq_list.append(residue_name) 1541 1542 protein_sequence = ''.join(protein_seq_list) 1543 clean_sequence = self.protein_sequence_parser(sequence=protein_sequence) 1544 1545 self.define_AA_particles_in_sequence (sequence=clean_sequence, sigma_dict = sigma_dict) 1546 residue_list = self.define_AA_residues(sequence=clean_sequence, model=model) 1547 1548 index = len(self.df) 1549 self.df.at [index,'name'] = name 1550 self.df.at [index,'pmb_type'] = 'protein' 1551 self.df.at [index,'model'] = model 1552 self.df.at [index,('sequence','')] = clean_sequence 1553 self.df.at [index,('residue_list','')] = residue_list 1554 1555 return 1556 1557 def define_residue(self, name, central_bead, side_chains): 1558 """ 1559 Defines a pyMBE object of type `residue` in `pymbe.df`. 1560 1561 Args: 1562 name (`str`): Unique label that identifies the `residue`. 1563 central_bead (`str`): `name` of the `particle` to be placed as central_bead of the `residue`. 1564 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. 1565 """ 1566 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'): 1567 return 1568 index = len(self.df) 1569 self.df.at [index, 'name'] = name 1570 self.df.at [index,'pmb_type'] = 'residue' 1571 self.df.at [index,'central_bead'] = central_bead 1572 self.df.at [index,('side_chains','')] = side_chains 1573 return 1574 1575 def destroy_pmb_object_in_system(self, name, espresso_system): 1576 """ 1577 Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 1578 1579 Args: 1580 name (str): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`. 1581 espresso_system (cls): Instance of a system class from espressomd library. 1582 1583 Note: 1584 - 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. 1585 """ 1586 allowed_objects = ['particle','residue','molecule'] 1587 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 1588 if pmb_type not in allowed_objects: 1589 raise ValueError('Object type not supported, supported types are ', allowed_objects) 1590 if pmb_type == 'particle': 1591 particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1592 & (self.df['molecule_id'].isna())] 1593 particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1594 & (self.df['molecule_id'].isna())].particle_id.tolist() 1595 for particle_id in particle_ids_list: 1596 espresso_system.part.by_id(particle_id).remove() 1597 self.df = self.df.drop(index=particles_index) 1598 if pmb_type == 'residue': 1599 residues_id = self.df.loc[self.df['name']== name].residue_id.to_list() 1600 for residue_id in residues_id: 1601 particle_ids_list = self.df.loc[self.df['residue_id']==residue_id].particle_id.dropna().to_list() 1602 self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index) 1603 for particle_id in particle_ids_list: 1604 espresso_system.part.by_id(particle_id).remove() 1605 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1606 if pmb_type == 'molecule': 1607 molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list() 1608 for molecule_id in molecules_id: 1609 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 1610 1611 self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index) 1612 for particle_id in particle_ids_list: 1613 espresso_system.part.by_id(particle_id).remove() 1614 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1615 1616 self.df.reset_index(drop=True,inplace=True) 1617 1618 return 1619 1620 def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200): 1621 """ 1622 Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration. 1623 To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an 1624 ion pair and the concentrations of the various species is solved numerically using a self-consistent approach. 1625 More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831). 1626 This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237). 1627 1628 Args: 1629 pH_res ('float'): pH-value in the reservoir. 1630 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 1631 activity_coefficient_monovalent_pair ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 1632 1633 Returns: 1634 cH_res ('float'): Concentration of H+ ions. 1635 cOH_res ('float'): Concentration of OH- ions. 1636 cNa_res ('float'): Concentration of Na+ ions. 1637 cCl_res ('float'): Concentration of Cl- ions. 1638 """ 1639 1640 self_consistent_run = 0 1641 cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+ 1642 1643 def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res): 1644 #Calculate and initial guess for the concentrations of various species based on ideal gas estimate 1645 cOH_res = self.Kw / cH_res 1646 cNa_res = None 1647 cCl_res = None 1648 if cOH_res>=cH_res: 1649 #adjust the concentration of sodium if there is excess OH- in the reservoir: 1650 cNa_res = c_salt_res + (cOH_res-cH_res) 1651 cCl_res = c_salt_res 1652 else: 1653 # adjust the concentration of chloride if there is excess H+ in the reservoir 1654 cCl_res = c_salt_res + (cH_res-cOH_res) 1655 cNa_res = c_salt_res 1656 1657 def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res): 1658 nonlocal max_number_sc_runs, self_consistent_run 1659 if self_consistent_run<max_number_sc_runs: 1660 self_consistent_run+=1 1661 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1662 cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res)) 1663 if cOH_res>=cH_res: 1664 #adjust the concentration of sodium if there is excess OH- in the reservoir: 1665 cNa_res = c_salt_res + (cOH_res-cH_res) 1666 cCl_res = c_salt_res 1667 else: 1668 # adjust the concentration of chloride if there is excess H+ in the reservoir 1669 cCl_res = c_salt_res + (cH_res-cOH_res) 1670 cNa_res = c_salt_res 1671 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 1672 else: 1673 return cH_res, cOH_res, cNa_res, cCl_res 1674 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 1675 1676 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 1677 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1678 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 1679 1680 while abs(determined_pH-pH_res)>1e-6: 1681 if determined_pH > pH_res: 1682 cH_res *= 1.005 1683 else: 1684 cH_res /= 1.003 1685 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 1686 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1687 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 1688 self_consistent_run=0 1689 1690 return cH_res, cOH_res, cNa_res, cCl_res 1691 1692 def filter_df(self, pmb_type): 1693 """ 1694 Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns. 1695 1696 Args: 1697 pmb_type(`str`): pmb_object_type to filter in `pmb.df`. 1698 1699 Returns: 1700 pmb_type_df(`obj`): filtered `pmb.df`. 1701 """ 1702 pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type] 1703 pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1) 1704 return pmb_type_df 1705 1706 def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False): 1707 """ 1708 Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 1709 1710 Args: 1711 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 1712 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 1713 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'. 1714 1715 Returns: 1716 bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` 1717 1718 Note: 1719 - If `use_default_bond`=`True`, it returns "default" if no key is found. 1720 """ 1721 bond_keys = [particle_name1 +'-'+ particle_name2, particle_name2 +'-'+ particle_name1 ] 1722 bond_defined=False 1723 for bond_key in bond_keys: 1724 if bond_key in self.df.values: 1725 bond_defined=True 1726 correct_key=bond_key 1727 break 1728 if bond_defined: 1729 return correct_key 1730 elif use_default_bond: 1731 return 'default' 1732 else: 1733 return 1734 1735 def find_value_from_es_type(self, es_type, column_name): 1736 """ 1737 Finds a value in `pmb.df` for a `column_name` and `es_type` pair. 1738 1739 Args: 1740 es_type(`int`): value of the espresso type 1741 column_name(`str`): name of the column in `pymbe.df` 1742 1743 Returns: 1744 Value in `pymbe.df` matching `column_name` and `es_type` 1745 """ 1746 idx = pd.IndexSlice 1747 for state in ['state_one', 'state_two']: 1748 index = np.where(self.df[(state, 'es_type')] == es_type)[0] 1749 if len(index) > 0: 1750 if column_name == 'label': 1751 label = self.df.loc[idx[index[0]], idx[(state,column_name)]] 1752 return label 1753 else: 1754 column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]] 1755 return column_name_value 1756 return None 1757 1758 def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples): 1759 """ 1760 Generates coordinates outside a sphere centered at `center`. 1761 1762 Args: 1763 center(`lst`): Coordinates of the center of the sphere. 1764 radius(`float`): Radius of the sphere. 1765 max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates. 1766 n_samples(`int`): Number of sample points. 1767 1768 Returns: 1769 coord_list(`lst`): Coordinates of the sample points. 1770 """ 1771 if not radius > 0: 1772 raise ValueError (f'The value of {radius} must be a positive value') 1773 if not radius < max_dist: 1774 raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))') 1775 coord_list = [] 1776 counter = 0 1777 while counter<n_samples: 1778 coord = self.generate_random_points_in_a_sphere(center=center, 1779 radius=max_dist, 1780 n_samples=1)[0] 1781 if np.linalg.norm(coord-np.asarray(center))>=radius: 1782 coord_list.append (coord) 1783 counter += 1 1784 return coord_list 1785 1786 def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False): 1787 """ 1788 Uniformly samples points from a hypersphere. If on_surface is set to True, the points are 1789 uniformly sampled from the surface of the hypersphere. 1790 1791 Args: 1792 center(`lst`): Array with the coordinates of the center of the spheres. 1793 radius(`float`): Radius of the sphere. 1794 n_samples(`int`): Number of sample points to generate inside the sphere. 1795 on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere. 1796 1797 Returns: 1798 samples(`list`): Coordinates of the sample points inside the hypersphere. 1799 1800 Note: 1801 - Algorithm from: https://baezortega.github.io/2018/10/14/hypersphere-sampling/ 1802 """ 1803 # initial values 1804 center=np.array(center) 1805 d = center.shape[0] 1806 # sample n_samples points in d dimensions from a standard normal distribution 1807 samples = self.rng.normal(size=(n_samples, d)) 1808 # make the samples lie on the surface of the unit hypersphere 1809 normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis] 1810 samples /= normalize_radii 1811 if not on_surface: 1812 # make the samples lie inside the hypersphere with the correct density 1813 uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis] 1814 new_radii = np.power(uniform_points, 1/d) 1815 samples *= new_radii 1816 # scale the points to have the correct radius and center 1817 samples = samples * radius + center 1818 return samples 1819 1820 def generate_trial_perpendicular_vector(self,vector,magnitude): 1821 """ 1822 Generates an orthogonal vector to the input `vector`. 1823 1824 Args: 1825 vector(`lst`): arbitrary vector. 1826 magnitude(`float`): magnitude of the orthogonal vector. 1827 1828 Returns: 1829 (`lst`): Orthogonal vector with the same magnitude as the input vector. 1830 """ 1831 np_vec = np.array(vector) 1832 np_vec /= np.linalg.norm(np_vec) 1833 if np.all(np_vec == 0): 1834 raise ValueError('Zero vector') 1835 # Generate a random vector 1836 random_vector = self.generate_random_points_in_a_sphere(radius=1, 1837 center=[0,0,0], 1838 n_samples=1, 1839 on_surface=True)[0] 1840 # Project the random vector onto the input vector and subtract the projection 1841 projection = np.dot(random_vector, np_vec) * np_vec 1842 perpendicular_vector = random_vector - projection 1843 # Normalize the perpendicular vector to have the same magnitude as the input vector 1844 perpendicular_vector /= np.linalg.norm(perpendicular_vector) 1845 return perpendicular_vector*magnitude 1846 1847 def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 1848 """ 1849 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length. 1850 If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead. 1851 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. 1852 1853 Args: 1854 particle_name1 (str): label of the type of the first particle type of the bonded particles. 1855 particle_name2 (str): label of the type of the second particle type of the bonded particles. 1856 hard_check (bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 1857 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. 1858 1859 Returns: 1860 l0 (`float`): bond length 1861 1862 Note: 1863 - 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`. 1864 - If `hard_check`=`True` stops the code when no bond is found. 1865 """ 1866 bond_key = self.find_bond_key(particle_name1=particle_name1, 1867 particle_name2=particle_name2, 1868 use_default_bond=use_default_bond) 1869 if bond_key: 1870 return self.df[self.df['name']==bond_key].l0.values[0] 1871 else: 1872 print("Bond not defined between particles ", particle_name1, " and ", particle_name2) 1873 if hard_check: 1874 sys.exit(1) 1875 else: 1876 return 1877 1878 def get_charge_map(self): 1879 ''' 1880 Gets the charge of each `espresso_type` in `pymbe.df`. 1881 1882 Returns: 1883 charge_map(`dict`): {espresso_type: charge}. 1884 ''' 1885 if self.df.state_one['es_type'].isnull().values.any(): 1886 df_state_one = self.df.state_one.dropna() 1887 df_state_two = self.df.state_two.dropna() 1888 else: 1889 df_state_one = self.df.state_one 1890 if self.df.state_two['es_type'].isnull().values.any(): 1891 df_state_two = self.df.state_two.dropna() 1892 else: 1893 df_state_two = self.df.state_two 1894 state_one = pd.Series (df_state_one.charge.values,index=df_state_one.es_type.values) 1895 state_two = pd.Series (df_state_two.charge.values,index=df_state_two.es_type.values) 1896 charge_map = pd.concat([state_one,state_two],axis=0).to_dict() 1897 return charge_map 1898 1899 def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'): 1900 """ 1901 Returns the Lennard-Jones parameters for the interaction between the particle types given by 1902 `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule. 1903 1904 Args: 1905 particle_name1 (str): label of the type of the first particle type 1906 particle_name2 (str): label of the type of the second particle type 1907 combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'. 1908 1909 Returns: 1910 {"epsilon": LJ epsilon, "sigma": LJ sigma} 1911 1912 Note: 1913 Currently, the only `combining_rule` supported is Lorentz-Berthelot. 1914 """ 1915 supported_combining_rules=["Lorentz-Berthelot"] 1916 1917 if combining_rule not in supported_combining_rules: 1918 raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}") 1919 1920 eps1 = self.df.loc[self.df["name"]==particle_name1].epsilon.iloc[0] 1921 sigma1 = self.df.loc[self.df["name"]==particle_name1].sigma.iloc[0] 1922 eps2 = self.df.loc[self.df["name"]==particle_name2].epsilon.iloc[0] 1923 sigma2 = self.df.loc[self.df["name"]==particle_name2].sigma.iloc[0] 1924 1925 return {"epsilon": np.sqrt(eps1*eps2), "sigma": (sigma1+sigma2)/2.0} 1926 1927 def get_particle_id_map(self, object_name): 1928 ''' 1929 Gets all the ids associated with the object with name `object_name` in `pmb.df` 1930 1931 Args: 1932 object_name(`str`): name of the object 1933 1934 Returns: 1935 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]}, } 1936 ''' 1937 object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0] 1938 valid_types = ["particle", "molecule", "residue", "protein"] 1939 if object_type not in valid_types: 1940 raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}") 1941 id_list = [] 1942 mol_map = {} 1943 res_map = {} 1944 def do_res_map(res_ids): 1945 for res_id in res_ids: 1946 res_list=self.df.loc[self.df['residue_id']== res_id].particle_id.dropna().tolist() 1947 res_map[res_id]=res_list 1948 return res_map 1949 if object_type in ['molecule', 'protein']: 1950 mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist() 1951 for mol_id in mol_ids: 1952 res_ids = set(self.df.loc[self.df['molecule_id']== mol_id].residue_id.dropna().tolist()) 1953 res_map=do_res_map(res_ids=res_ids) 1954 mol_list=self.df.loc[self.df['molecule_id']== mol_id].particle_id.dropna().tolist() 1955 id_list+=mol_list 1956 mol_map[mol_id]=mol_list 1957 elif object_type == 'residue': 1958 res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist() 1959 res_map=do_res_map(res_ids=res_ids) 1960 id_list=[] 1961 for res_id_list in res_map.values(): 1962 id_list+=res_id_list 1963 elif object_type == 'particle': 1964 id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist() 1965 return {"all": id_list, "molecule_map": mol_map, "residue_map": res_map} 1966 1967 def get_pka_set(self): 1968 ''' 1969 Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df` 1970 1971 Returns: 1972 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 1973 ''' 1974 titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna() 1975 pka_set = {} 1976 for index in titratables_AA_df.name.keys(): 1977 name = titratables_AA_df.name[index] 1978 pka_value = titratables_AA_df.pka[index] 1979 acidity = titratables_AA_df.acidity[index] 1980 pka_set[name] = {'pka_value':pka_value,'acidity':acidity} 1981 return pka_set 1982 1983 def get_radius_map(self): 1984 ''' 1985 Gets the effective radius of each `espresso type` in `pmb.df`. 1986 1987 Returns: 1988 radius_map(`dict`): {espresso_type: radius}. 1989 ''' 1990 1991 df_state_one = self.df[[('sigma',''),('state_one','es_type')]].dropna().drop_duplicates() 1992 df_state_two = self.df[[('sigma',''),('state_two','es_type')]].dropna().drop_duplicates() 1993 1994 state_one = pd.Series(df_state_one.sigma.values/2.0,index=df_state_one.state_one.es_type.values) 1995 state_two = pd.Series(df_state_two.sigma.values/2.0,index=df_state_two.state_two.es_type.values) 1996 radius_map = pd.concat([state_one,state_two],axis=0).to_dict() 1997 1998 return radius_map 1999 2000 def get_resource(self, path): 2001 ''' 2002 Locate a file resource of the pyMBE package. 2003 2004 Args: 2005 path(`str`): Relative path to the resource 2006 2007 Returns: 2008 path(`str`): Absolute path to the resource 2009 2010 ''' 2011 import os 2012 return os.path.join(os.path.dirname(__file__), path) 2013 2014 def get_type_map(self): 2015 """ 2016 Gets all different espresso types assigned to particles in `pmb.df`. 2017 2018 Returns: 2019 type_map(`dict`): {"name": espresso_type}. 2020 """ 2021 if self.df.state_one['es_type'].isnull().values.any(): 2022 df_state_one = self.df.state_one.dropna(how='all') 2023 df_state_two = self.df.state_two.dropna(how='all') 2024 else: 2025 df_state_one = self.df.state_one 2026 if self.df.state_two['es_type'].isnull().values.any(): 2027 df_state_two = self.df.state_two.dropna(how='all') 2028 else: 2029 df_state_two = self.df.state_two 2030 state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label) 2031 state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label) 2032 type_map = pd.concat([state_one,state_two],axis=0).to_dict() 2033 return type_map 2034 2035 def load_interaction_parameters(self, filename, verbose=True): 2036 """ 2037 Loads the interaction parameters stored in `filename` into `pmb.df` 2038 2039 Args: 2040 filename(`str`): name of the file to be read 2041 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2042 """ 2043 without_units = ['q','es_type','acidity'] 2044 with_units = ['sigma','epsilon'] 2045 2046 with open(filename, 'r') as f: 2047 interaction_data = json.load(f) 2048 interaction_parameter_set = interaction_data["data"] 2049 2050 for key in interaction_parameter_set: 2051 param_dict=interaction_parameter_set[key] 2052 object_type=param_dict['object_type'] 2053 if object_type == 'particle': 2054 not_requiered_attributes={} 2055 for not_requiered_key in without_units+with_units: 2056 if not_requiered_key in param_dict.keys(): 2057 if not_requiered_key in with_units: 2058 not_requiered_attributes[not_requiered_key]=self.create_variable_with_units(variable=param_dict.pop(not_requiered_key)) 2059 elif not_requiered_key in without_units: 2060 not_requiered_attributes[not_requiered_key]=param_dict.pop(not_requiered_key) 2061 else: 2062 if not_requiered_key == 'acidity': 2063 not_requiered_attributes[not_requiered_key] = 'inert' 2064 else: 2065 not_requiered_attributes[not_requiered_key]=None 2066 self.define_particle(name=param_dict.pop('name'), 2067 q=not_requiered_attributes.pop('q'), 2068 sigma=not_requiered_attributes.pop('sigma'), 2069 acidity=not_requiered_attributes.pop('acidity'), 2070 epsilon=not_requiered_attributes.pop('epsilon'), 2071 verbose=verbose) 2072 elif object_type == 'residue': 2073 self.define_residue (name = param_dict.pop('name'), 2074 central_bead = param_dict.pop('central_bead_name'), 2075 side_chains = param_dict.pop('side_chains_names')) 2076 elif object_type == 'molecule': 2077 self.define_molecule(name=param_dict.pop('name'), 2078 residue_list=param_dict.pop('residue_name_list')) 2079 elif object_type == 'peptide': 2080 self.define_peptide(name=param_dict.pop('name'), 2081 sequence=param_dict.pop('sequence'), 2082 model=param_dict.pop('model')) 2083 elif object_type == 'bond': 2084 particle_pairs = param_dict.pop('particle_pairs') 2085 bond_parameters = param_dict.pop('bond_parameters') 2086 bond_type = param_dict.pop('bond_type') 2087 if bond_type == 'harmonic': 2088 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2089 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2090 bond = {'r_0' : r_0, 2091 'k' : k, 2092 } 2093 2094 elif bond_type == 'FENE': 2095 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2096 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2097 d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max')) 2098 bond = {'r_0' : r_0, 2099 'k' : k, 2100 'd_r_max': d_r_max, 2101 } 2102 else: 2103 raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds") 2104 self.define_bond(bond_type=bond_type, bond_parameters=bond, particle_pairs=particle_pairs) 2105 else: 2106 raise ValueError(object_type+' is not a known pmb object type') 2107 2108 return 2109 2110 def load_pka_set(self, filename, verbose=True): 2111 """ 2112 Loads the pka_set stored in `filename` into `pmb.df`. 2113 2114 Args: 2115 filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}. 2116 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2117 2118 Note: 2119 - If `name` is already defined in the `pymbe.df`, it prints a warning. 2120 """ 2121 with open(filename, 'r') as f: 2122 pka_data = json.load(f) 2123 pka_set = pka_data["data"] 2124 2125 self.check_pka_set(pka_set=pka_set) 2126 2127 for key in pka_set: 2128 acidity = pka_set[key]['acidity'] 2129 pka_value = pka_set[key]['pka_value'] 2130 self.define_particle( 2131 name=key, 2132 acidity=acidity, 2133 pka=pka_value, 2134 verbose=verbose) 2135 return 2136 2137 def parse_sequence_from_file(self,sequence): 2138 """ 2139 Parses the given sequence such that it can be used in pyMBE. This function has to be used if the df was read from a file. 2140 2141 Args: 2142 sequence(`str`): sequence to be parsed 2143 2144 Returns: 2145 sequence(`lst`): parsed sequence 2146 """ 2147 sequence = sequence.replace(' ', '') 2148 sequence = sequence.replace("'", '') 2149 sequence = sequence.split(",")[1:-1] 2150 return sequence 2151 2152 def print_reduced_units(self): 2153 """ 2154 Prints the current set of reduced units defined in pyMBE.units. 2155 """ 2156 print("\nCurrent set of reduced units:") 2157 unit_length=self.units.Quantity(1,'reduced_length') 2158 unit_energy=self.units.Quantity(1,'reduced_energy') 2159 unit_charge=self.units.Quantity(1,'reduced_charge') 2160 print(unit_length.to('nm'), "=", unit_length) 2161 print(unit_energy.to('J'), "=", unit_energy) 2162 print('Temperature:', (self.kT/self.Kb).to("K")) 2163 print(unit_charge.to('C'), "=", unit_charge) 2164 print() 2165 2166 def propose_unused_type(self): 2167 """ 2168 Searches in `pmb.df` all the different particle types defined and returns a new one. 2169 2170 Returns: 2171 unused_type(`int`): unused particle type 2172 """ 2173 type_map=self.get_type_map() 2174 if type_map == {}: 2175 unused_type = 0 2176 else: 2177 unused_type=max(type_map.values())+1 2178 return unused_type 2179 2180 def protein_sequence_parser(self, sequence): 2181 ''' 2182 Parses `sequence` to the one letter code for amino acids. 2183 2184 Args: 2185 sequence(`str` or `lst`): Sequence of the amino acid. 2186 2187 Returns: 2188 clean_sequence(`lst`): `sequence` using the one letter code. 2189 2190 Note: 2191 - Accepted formats for `sequence` are: 2192 - `lst` with one letter or three letter code of each aminoacid in each element 2193 - `str` with the sequence using the one letter code 2194 - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-" 2195 2196 ''' 2197 # Aminoacid key 2198 keys={"ALA": "A", 2199 "ARG": "R", 2200 "ASN": "N", 2201 "ASP": "D", 2202 "CYS": "C", 2203 "GLU": "E", 2204 "GLN": "Q", 2205 "GLY": "G", 2206 "HIS": "H", 2207 "ILE": "I", 2208 "LEU": "L", 2209 "LYS": "K", 2210 "MET": "M", 2211 "PHE": "F", 2212 "PRO": "P", 2213 "SER": "S", 2214 "THR": "T", 2215 "TRP": "W", 2216 "TYR": "Y", 2217 "VAL": "V", 2218 "PSER": "J", 2219 "PTHR": "U", 2220 "PTyr": "Z", 2221 "NH2": "n", 2222 "COOH": "c"} 2223 clean_sequence=[] 2224 if isinstance(sequence, str): 2225 if sequence.find("-") != -1: 2226 splited_sequence=sequence.split("-") 2227 for residue in splited_sequence: 2228 if len(residue) == 1: 2229 if residue in keys.values(): 2230 residue_ok=residue 2231 else: 2232 if residue.upper() in keys.values(): 2233 residue_ok=residue.upper() 2234 else: 2235 raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence") 2236 clean_sequence.append(residue_ok) 2237 else: 2238 if residue in keys.keys(): 2239 clean_sequence.append(keys[residue]) 2240 else: 2241 if residue.upper() in keys.keys(): 2242 clean_sequence.append(keys[residue.upper()]) 2243 else: 2244 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2245 else: 2246 for residue in sequence: 2247 if residue in keys.values(): 2248 residue_ok=residue 2249 else: 2250 if residue.upper() in keys.values(): 2251 residue_ok=residue.upper() 2252 else: 2253 raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence") 2254 clean_sequence.append(residue_ok) 2255 if isinstance(sequence, list): 2256 for residue in sequence: 2257 if residue in keys.values(): 2258 residue_ok=residue 2259 else: 2260 if residue.upper() in keys.values(): 2261 residue_ok=residue.upper() 2262 elif (residue.upper() in keys.keys()): 2263 clean_sequence.append(keys[residue.upper()]) 2264 else: 2265 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2266 clean_sequence.append(residue_ok) 2267 return clean_sequence 2268 2269 def read_pmb_df (self,filename): 2270 """ 2271 Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE. 2272 2273 Args: 2274 filename(str): path to file. 2275 2276 Note: 2277 This function only accepts files with CSV format. 2278 """ 2279 2280 if filename.rsplit(".", 1)[1] != "csv": 2281 raise ValueError("Only files with CSV format are supported") 2282 df = pd.read_csv (filename,header=[0, 1], index_col=0) 2283 columns_names = self.setup_df() 2284 2285 multi_index = pd.MultiIndex.from_tuples(columns_names) 2286 df.columns = multi_index 2287 2288 self.df = self.convert_columns_to_original_format(df) 2289 2290 return self.df 2291 2292 def read_protein_vtf_in_df (self,filename,unit_length=None): 2293 """ 2294 Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`. 2295 2296 Args: 2297 filename(`str`): Path to the vtf file with the coarse-grained model. 2298 unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None. 2299 2300 Returns: 2301 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 2302 2303 Note: 2304 - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom. 2305 """ 2306 2307 print (f'Loading protein coarse grain model file: {filename}') 2308 2309 coord_list = [] 2310 particles_dict = {} 2311 2312 if unit_length is None: 2313 unit_length = 1 * self.units.angstrom 2314 2315 with open (filename,'r') as protein_model: 2316 for line in protein_model : 2317 line_split = line.split() 2318 if line_split: 2319 line_header = line_split[0] 2320 if line_header == 'atom': 2321 atom_id = line_split[1] 2322 atom_name = line_split[3] 2323 atom_resname = line_split[5] 2324 chain_id = line_split[9] 2325 radius = float(line_split [11])*unit_length 2326 sigma = 2*radius 2327 particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, sigma] 2328 elif line_header.isnumeric(): 2329 atom_coord = line_split[1:] 2330 atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord] 2331 coord_list.append (atom_coord) 2332 2333 numbered_label = [] 2334 i = 0 2335 2336 for atom_id in particles_dict.keys(): 2337 2338 if atom_id == 1: 2339 atom_name = particles_dict[atom_id][0] 2340 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2341 numbered_label.append(numbered_name) 2342 2343 elif atom_id != 1: 2344 2345 if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]: 2346 i += 1 2347 count = 1 2348 atom_name = particles_dict[atom_id][0] 2349 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2350 numbered_label.append(numbered_name) 2351 2352 elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]: 2353 if count == 2 or particles_dict[atom_id][1] == 'GLY': 2354 i +=1 2355 count = 0 2356 atom_name = particles_dict[atom_id][0] 2357 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2358 numbered_label.append(numbered_name) 2359 count +=1 2360 2361 topology_dict = {} 2362 2363 for i in range (0, len(numbered_label)): 2364 topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] , 2365 'chain_id':numbered_label[i][1], 2366 'sigma':numbered_label[i][2] } 2367 2368 return topology_dict 2369 2370 def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2371 """ 2372 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2373 If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead. 2374 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. 2375 2376 Args: 2377 particle_name1 (str): label of the type of the first particle type of the bonded particles. 2378 particle_name2 (str): label of the type of the second particle type of the bonded particles. 2379 hard_check (bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2380 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. 2381 2382 Returns: 2383 bond (cls): bond object from the espressomd library. 2384 2385 Note: 2386 - 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`. 2387 - If `hard_check`=`True` stops the code when no bond is found. 2388 """ 2389 2390 bond_key = self.find_bond_key(particle_name1=particle_name1, 2391 particle_name2=particle_name2, 2392 use_default_bond=use_default_bond) 2393 if use_default_bond: 2394 if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'): 2395 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") 2396 if bond_key: 2397 return self.df[self.df['name']==bond_key].bond_object.values[0] 2398 else: 2399 print("Bond not defined between particles ", particle_name1, " and ", particle_name2) 2400 if hard_check: 2401 sys.exit(1) 2402 else: 2403 return 2404 2405 def search_particles_in_residue(self, residue_name): 2406 ''' 2407 Searches for all particles in a given residue of name `residue_name`. 2408 2409 Args: 2410 residue_name (`str`): name of the residue to be searched 2411 2412 Returns: 2413 list_of_particles_in_residue (`lst`): list of the names of all particles in the residue 2414 2415 Note: 2416 - 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. 2417 2418 ''' 2419 index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 2420 central_bead = self.df.at [index_residue, ('central_bead', '')] 2421 list_of_side_chains = self.df.at [index_residue, ('side_chains', '')] 2422 2423 list_of_particles_in_residue = [] 2424 list_of_particles_in_residue.append(central_bead) 2425 2426 for side_chain in list_of_side_chains: 2427 object_type = self.df[self.df['name']==side_chain].pmb_type.values[0] 2428 2429 if object_type == "residue": 2430 list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain) 2431 list_of_particles_in_residue += list_of_particles_in_side_chain_residue 2432 elif object_type == "particle": 2433 list_of_particles_in_residue.append(side_chain) 2434 2435 return list_of_particles_in_residue 2436 2437 def set_particle_acidity(self, name, acidity='inert', default_charge=0, pka=None, verbose=True): 2438 """ 2439 Sets the particle acidity if it is acidic or basic, creates `state_one` and `state_two` with the protonated and 2440 deprotonated states. In each state is set: `label`,`charge` and `es_type`. If it is inert, it will define only `state_one`. 2441 2442 Args: 2443 name (`str`): Unique label that identifies the `particle`. 2444 acidity (`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to `inert`. 2445 default_charge (`int`): Charge of the particle. Defaults to 0. 2446 pka (`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to None. 2447 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2448 """ 2449 acidity_valid_keys = ['inert','acidic', 'basic'] 2450 if acidity not in acidity_valid_keys: 2451 raise ValueError(f"Acidity {acidity} provided for particle name {name} is not supproted. Valid keys are: {acidity_valid_keys}") 2452 if acidity in ['acidic', 'basic'] and pka is None: 2453 raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.") 2454 for index in self.df[self.df['name']==name].index: 2455 if pka: 2456 self.add_value_to_df(key=('pka',''), 2457 index=index, 2458 new_value=pka, 2459 verbose=verbose) 2460 self.add_value_to_df(key=('acidity',''), 2461 index=index, 2462 new_value=acidity, 2463 verbose=verbose) 2464 if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')): 2465 self.add_value_to_df(key=('state_one','es_type'), 2466 index=index, 2467 new_value=self.propose_unused_type(), 2468 verbose=verbose) 2469 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'inert': 2470 self.add_value_to_df(key=('state_one','charge'), 2471 index=index, 2472 new_value=default_charge, 2473 verbose=verbose) 2474 self.add_value_to_df(key=('state_one','label'), 2475 index=index, 2476 new_value=name, 2477 verbose=verbose) 2478 else: 2479 protonated_label = f'{name}H' 2480 self.add_value_to_df(key=('state_one','label'), 2481 index=index, 2482 new_value=protonated_label, 2483 verbose=verbose) 2484 self.add_value_to_df(key=('state_two','label'), 2485 index=index, 2486 new_value=name, 2487 verbose=verbose) 2488 if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')): 2489 self.add_value_to_df(key=('state_two','es_type'), 2490 index=index, 2491 new_value=self.propose_unused_type(), 2492 verbose=verbose) 2493 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'acidic': 2494 self.add_value_to_df(key=('state_one','charge'), 2495 index=index,new_value=0, 2496 verbose=verbose) 2497 self.add_value_to_df(key=('state_two','charge'), 2498 index=index, 2499 new_value=-1, 2500 verbose=verbose) 2501 elif self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'basic': 2502 self.add_value_to_df(key=('state_one','charge'), 2503 index=index,new_value=+1, 2504 verbose=verbose) 2505 self.add_value_to_df(key=('state_two','charge'), 2506 index=index, 2507 new_value=0, 2508 verbose=verbose) 2509 return 2510 2511 def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True): 2512 """ 2513 Sets the set of reduced units used by pyMBE.units and it prints it. 2514 2515 Args: 2516 unit_length (`obj`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 2517 unit_charge (`obj`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 2518 temperature (`obj`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 2519 Kw (`obj`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 2520 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2521 2522 Note: 2523 - If no `temperature` is given, a value of 298.15 K is assumed by default. 2524 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 2525 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 2526 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 2527 """ 2528 self.units=pint.UnitRegistry() 2529 if unit_length is None: 2530 unit_length=0.355*self.units.nm 2531 if temperature is None: 2532 temperature=298.15 * self.units.K 2533 if unit_charge is None: 2534 unit_charge=self.units.e 2535 if Kw is None: 2536 Kw = 1e-14 2537 self.N_A=6.02214076e23 / self.units.mol 2538 self.Kb=1.38064852e-23 * self.units.J / self.units.K 2539 self.e=1.60217662e-19 *self.units.C 2540 self.kT=temperature*self.Kb 2541 self.Kw=Kw*self.units.mol**2 / (self.units.l**2) 2542 self.units.define(f'reduced_energy = {self.kT} ') 2543 self.units.define(f'reduced_length = {unit_length}') 2544 self.units.define(f'reduced_charge = {unit_charge}') 2545 if verbose: 2546 self.print_reduced_units() 2547 return 2548 2549 def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2550 """ 2551 Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 2552 2553 Args: 2554 counter_ion(`str`): `name` of the counter_ion `particle`. 2555 constant_pH(`float`): pH-value. 2556 exclusion_range(`float`, optional): Below this value, no particles will be inserted. 2557 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2558 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2559 2560 Returns: 2561 RE (`obj`): Instance of a reaction_ensemble.ConstantpHEnsemble object from the espressomd library. 2562 sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2563 """ 2564 from espressomd import reaction_methods 2565 if exclusion_range is None: 2566 exclusion_range = max(self.get_radius_map().values())*2.0 2567 if pka_set is None: 2568 pka_set=self.get_pka_set() 2569 self.check_pka_set(pka_set=pka_set) 2570 if use_exclusion_radius_per_type: 2571 exclusion_radius_per_type = self.get_radius_map() 2572 else: 2573 exclusion_radius_per_type = {} 2574 2575 RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2576 exclusion_range=exclusion_range.magnitude, 2577 seed=self.SEED, 2578 constant_pH=constant_pH, 2579 exclusion_radius_per_type = exclusion_radius_per_type 2580 ) 2581 sucessfull_reactions_labels=[] 2582 charge_map = self.get_charge_map() 2583 for name in pka_set.keys(): 2584 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2585 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in sg.type_map') 2586 continue 2587 gamma=10**-pka_set[name]['pka_value'] 2588 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2589 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2590 counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0] 2591 RE.add_reaction(gamma=gamma, 2592 reactant_types=[state_one_type], 2593 product_types=[state_two_type, counter_ion_type], 2594 default_charges={state_one_type: charge_map[state_one_type], 2595 state_two_type: charge_map[state_two_type], 2596 counter_ion_type: charge_map[counter_ion_type]}) 2597 sucessfull_reactions_labels.append(name) 2598 return RE, sucessfull_reactions_labels 2599 2600 def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient=None, exclusion_range=None, use_exclusion_radius_per_type = False): 2601 """ 2602 Sets up grand-canonical coupling to a reservoir of salt. 2603 For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead. 2604 2605 Args: 2606 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2607 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 2608 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 2609 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2610 exclusion_range(`float`, optional): For distances shorter than this value, no particles will be inserted. 2611 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2612 2613 Returns: 2614 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2615 """ 2616 from espressomd import reaction_methods 2617 if exclusion_range is None: 2618 exclusion_range = max(self.get_radius_map().values())*2.0 2619 if use_exclusion_radius_per_type: 2620 exclusion_radius_per_type = self.get_radius_map() 2621 else: 2622 exclusion_radius_per_type = {} 2623 2624 #If no function for the activity coefficient is provided, the ideal case is assumed. 2625 if activity_coefficient is None: 2626 activity_coefficient = lambda x: 1.0 2627 2628 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2629 exclusion_range=exclusion_range.magnitude, 2630 seed=self.SEED, 2631 exclusion_radius_per_type = exclusion_radius_per_type 2632 ) 2633 2634 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2635 determined_activity_coefficient = activity_coefficient(c_salt_res) 2636 K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient 2637 2638 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 2639 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 2640 2641 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.charge.values[0] 2642 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.charge.values[0] 2643 2644 if salt_cation_charge <= 0: 2645 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 2646 if salt_anion_charge >= 0: 2647 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 2648 2649 # Grand-canonical coupling to the reservoir 2650 RE.add_reaction( 2651 gamma = K_salt.magnitude, 2652 reactant_types = [], 2653 reactant_coefficients = [], 2654 product_types = [ salt_cation_es_type, salt_anion_es_type ], 2655 product_coefficients = [ 1, 1 ], 2656 default_charges = { 2657 salt_cation_es_type: salt_cation_charge, 2658 salt_anion_es_type: salt_anion_charge, 2659 } 2660 ) 2661 2662 return RE 2663 2664 def setup_grxmc_reactions(self, pH_res, c_salt_res, proton_name, hydroxide_name, salt_cation_name, salt_anion_name, activity_coefficient=None, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2665 """ 2666 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 2667 reservoir of small ions. 2668 This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1]. 2669 2670 [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. 2671 2672 Args: 2673 pH_res ('float): pH-value in the reservoir. 2674 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2675 proton_name ('str'): Name of the proton (H+) particle. 2676 hydroxide_name ('str'): Name of the hydroxide (OH-) particle. 2677 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 2678 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 2679 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2680 exclusion_range(`float`, optional): For distances shorter than this value, no particles will be inserted. 2681 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2682 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2683 2684 Returns: 2685 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2686 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2687 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 2688 """ 2689 from espressomd import reaction_methods 2690 if exclusion_range is None: 2691 exclusion_range = max(self.get_radius_map().values())*2.0 2692 if pka_set is None: 2693 pka_set=self.get_pka_set() 2694 self.check_pka_set(pka_set=pka_set) 2695 if use_exclusion_radius_per_type: 2696 exclusion_radius_per_type = self.get_radius_map() 2697 else: 2698 exclusion_radius_per_type = {} 2699 2700 #If no function for the activity coefficient is provided, the ideal case is assumed. 2701 if activity_coefficient is None: 2702 activity_coefficient = lambda x: 1.0 2703 2704 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2705 exclusion_range=exclusion_range.magnitude, 2706 seed=self.SEED, 2707 exclusion_radius_per_type = exclusion_radius_per_type 2708 ) 2709 2710 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2711 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 2712 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2713 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 2714 K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2715 K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2716 K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2717 2718 proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0] 2719 hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0] 2720 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 2721 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 2722 2723 proton_charge = self.df.loc[self.df['name']==proton_name].state_one.charge.values[0] 2724 hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.charge.values[0] 2725 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.charge.values[0] 2726 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.charge.values[0] 2727 2728 if proton_charge <= 0: 2729 raise ValueError('ERROR proton charge must be positive, charge ', proton_charge) 2730 if salt_cation_charge <= 0: 2731 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 2732 if hydroxide_charge >= 0: 2733 raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge) 2734 if salt_anion_charge >= 0: 2735 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 2736 2737 # Grand-canonical coupling to the reservoir 2738 # 0 = H+ + OH- 2739 RE.add_reaction( 2740 gamma = K_W.magnitude, 2741 reactant_types = [], 2742 reactant_coefficients = [], 2743 product_types = [ proton_es_type, hydroxide_es_type ], 2744 product_coefficients = [ 1, 1 ], 2745 default_charges = { 2746 proton_es_type: proton_charge, 2747 hydroxide_es_type: hydroxide_charge, 2748 } 2749 ) 2750 2751 # 0 = Na+ + Cl- 2752 RE.add_reaction( 2753 gamma = K_NACL.magnitude, 2754 reactant_types = [], 2755 reactant_coefficients = [], 2756 product_types = [ salt_cation_es_type, salt_anion_es_type ], 2757 product_coefficients = [ 1, 1 ], 2758 default_charges = { 2759 salt_cation_es_type: salt_cation_charge, 2760 salt_anion_es_type: salt_anion_charge, 2761 } 2762 ) 2763 2764 # 0 = Na+ + OH- 2765 RE.add_reaction( 2766 gamma = (K_NACL * K_W / K_HCL).magnitude, 2767 reactant_types = [], 2768 reactant_coefficients = [], 2769 product_types = [ salt_cation_es_type, hydroxide_es_type ], 2770 product_coefficients = [ 1, 1 ], 2771 default_charges = { 2772 salt_cation_es_type: salt_cation_charge, 2773 hydroxide_es_type: hydroxide_charge, 2774 } 2775 ) 2776 2777 # 0 = H+ + Cl- 2778 RE.add_reaction( 2779 gamma = K_HCL.magnitude, 2780 reactant_types = [], 2781 reactant_coefficients = [], 2782 product_types = [ proton_es_type, salt_anion_es_type ], 2783 product_coefficients = [ 1, 1 ], 2784 default_charges = { 2785 proton_es_type: proton_charge, 2786 salt_anion_es_type: salt_anion_charge, 2787 } 2788 ) 2789 2790 # Annealing moves to ensure sufficient sampling 2791 # Cation annealing H+ = Na+ 2792 RE.add_reaction( 2793 gamma = (K_NACL / K_HCL).magnitude, 2794 reactant_types = [proton_es_type], 2795 reactant_coefficients = [ 1 ], 2796 product_types = [ salt_cation_es_type ], 2797 product_coefficients = [ 1 ], 2798 default_charges = { 2799 proton_es_type: proton_charge, 2800 salt_cation_es_type: salt_cation_charge, 2801 } 2802 ) 2803 2804 # Anion annealing OH- = Cl- 2805 RE.add_reaction( 2806 gamma = (K_HCL / K_W).magnitude, 2807 reactant_types = [hydroxide_es_type], 2808 reactant_coefficients = [ 1 ], 2809 product_types = [ salt_anion_es_type ], 2810 product_coefficients = [ 1 ], 2811 default_charges = { 2812 hydroxide_es_type: hydroxide_charge, 2813 salt_anion_es_type: salt_anion_charge, 2814 } 2815 ) 2816 2817 sucessful_reactions_labels=[] 2818 charge_map = self.get_charge_map() 2819 for name in pka_set.keys(): 2820 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2821 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.') 2822 continue 2823 2824 Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 2825 2826 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2827 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2828 2829 # Reaction in terms of proton: HA = A + H+ 2830 RE.add_reaction(gamma=Ka.magnitude, 2831 reactant_types=[state_one_type], 2832 reactant_coefficients=[1], 2833 product_types=[state_two_type, proton_es_type], 2834 product_coefficients=[1, 1], 2835 default_charges={state_one_type: charge_map[state_one_type], 2836 state_two_type: charge_map[state_two_type], 2837 proton_es_type: proton_charge}) 2838 2839 # Reaction in terms of salt cation: HA = A + Na+ 2840 RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude, 2841 reactant_types=[state_one_type], 2842 reactant_coefficients=[1], 2843 product_types=[state_two_type, salt_cation_es_type], 2844 product_coefficients=[1, 1], 2845 default_charges={state_one_type: charge_map[state_one_type], 2846 state_two_type: charge_map[state_two_type], 2847 salt_cation_es_type: salt_cation_charge}) 2848 2849 # Reaction in terms of hydroxide: OH- + HA = A 2850 RE.add_reaction(gamma=(Ka / K_W).magnitude, 2851 reactant_types=[state_one_type, hydroxide_es_type], 2852 reactant_coefficients=[1, 1], 2853 product_types=[state_two_type], 2854 product_coefficients=[1], 2855 default_charges={state_one_type: charge_map[state_one_type], 2856 state_two_type: charge_map[state_two_type], 2857 hydroxide_es_type: hydroxide_charge}) 2858 2859 # Reaction in terms of salt anion: Cl- + HA = A 2860 RE.add_reaction(gamma=(Ka / K_HCL).magnitude, 2861 reactant_types=[state_one_type, salt_anion_es_type], 2862 reactant_coefficients=[1, 1], 2863 product_types=[state_two_type], 2864 product_coefficients=[1], 2865 default_charges={state_one_type: charge_map[state_one_type], 2866 state_two_type: charge_map[state_two_type], 2867 salt_anion_es_type: salt_anion_charge}) 2868 2869 sucessful_reactions_labels.append(name) 2870 return RE, sucessful_reactions_labels, ionic_strength_res 2871 2872 def setup_grxmc_unified (self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient=None, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2873 """ 2874 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 2875 reservoir of small ions. 2876 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-}. 2877 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'. 2878 2879 [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4). 2880 [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. 2881 2882 Args: 2883 pH_res ('float'): pH-value in the reservoir. 2884 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2885 cation_name ('str'): Name of the cationic particle. 2886 anion_name ('str'): Name of the anionic particle. 2887 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2888 exclusion_range(`float`, optional): Below this value, no particles will be inserted. 2889 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2890 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`. 2891 2892 Returns: 2893 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2894 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2895 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 2896 """ 2897 from espressomd import reaction_methods 2898 if exclusion_range is None: 2899 exclusion_range = max(self.get_radius_map().values())*2.0 2900 if pka_set is None: 2901 pka_set=self.get_pka_set() 2902 self.check_pka_set(pka_set=pka_set) 2903 if use_exclusion_radius_per_type: 2904 exclusion_radius_per_type = self.get_radius_map() 2905 else: 2906 exclusion_radius_per_type = {} 2907 2908 #If no function for the activity coefficient is provided, the ideal case is assumed. 2909 if activity_coefficient is None: 2910 activity_coefficient = lambda x: 1.0 2911 2912 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2913 exclusion_range=exclusion_range.magnitude, 2914 seed=self.SEED, 2915 exclusion_radius_per_type = exclusion_radius_per_type 2916 ) 2917 2918 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2919 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 2920 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2921 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 2922 a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 2923 a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 2924 a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 2925 K_XX = a_cation * a_anion 2926 2927 cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0] 2928 anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0] 2929 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.values[0] 2930 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.values[0] 2931 if cation_charge <= 0: 2932 raise ValueError('ERROR cation charge must be positive, charge ', cation_charge) 2933 if anion_charge >= 0: 2934 raise ValueError('ERROR anion charge must be negative, charge ', anion_charge) 2935 2936 # Coupling to the reservoir: 0 = X+ + X- 2937 RE.add_reaction( 2938 gamma = K_XX.magnitude, 2939 reactant_types = [], 2940 reactant_coefficients = [], 2941 product_types = [ cation_es_type, anion_es_type ], 2942 product_coefficients = [ 1, 1 ], 2943 default_charges = { 2944 cation_es_type: cation_charge, 2945 anion_es_type: anion_charge, 2946 } 2947 ) 2948 2949 sucessful_reactions_labels=[] 2950 charge_map = self.get_charge_map() 2951 for name in pka_set.keys(): 2952 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2953 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in sg.type_map') 2954 continue 2955 2956 Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l 2957 gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen 2958 2959 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2960 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2961 2962 # Reaction in terms of small cation: HA = A + X+ 2963 RE.add_reaction(gamma=gamma_K_AX.magnitude, 2964 reactant_types=[state_one_type], 2965 reactant_coefficients=[1], 2966 product_types=[state_two_type, cation_es_type], 2967 product_coefficients=[1, 1], 2968 default_charges={state_one_type: charge_map[state_one_type], 2969 state_two_type: charge_map[state_two_type], 2970 cation_es_type: cation_charge}) 2971 2972 # Reaction in terms of small anion: X- + HA = A 2973 RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude, 2974 reactant_types=[state_one_type, anion_es_type], 2975 reactant_coefficients=[1, 1], 2976 product_types=[state_two_type], 2977 product_coefficients=[1], 2978 default_charges={state_one_type: charge_map[state_one_type], 2979 state_two_type: charge_map[state_two_type], 2980 anion_es_type: anion_charge}) 2981 2982 sucessful_reactions_labels.append(name) 2983 return RE, sucessful_reactions_labels, ionic_strength_res 2984 2985 def setup_df (self): 2986 """ 2987 Sets up the pyMBE's dataframe `pymbe.df`. 2988 2989 Returns: 2990 columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe 2991 """ 2992 2993 columns_dtypes = { 2994 'name': { 2995 '': str}, 2996 'pmb_type': { 2997 '': str}, 2998 'particle_id': { 2999 '': object}, 3000 'particle_id2': { 3001 '': object}, 3002 'residue_id': { 3003 '': object}, 3004 'molecule_id': { 3005 '': object}, 3006 'acidity': { 3007 '': str}, 3008 'pka': { 3009 '': float}, 3010 'central_bead': { 3011 '': object}, 3012 'side_chains': { 3013 '': object}, 3014 'residue_list': { 3015 '': object}, 3016 'model': { 3017 '': str}, 3018 'sigma': { 3019 '': object}, 3020 'cutoff': { 3021 '': object}, 3022 'offset': { 3023 '': object}, 3024 'epsilon': { 3025 '': object}, 3026 'state_one': { 3027 'label': str, 3028 'es_type': object, 3029 'charge': object }, 3030 'state_two': { 3031 'label': str, 3032 'es_type': object, 3033 'charge': object }, 3034 'sequence': { 3035 '': object}, 3036 'bond_object': { 3037 '': object}, 3038 'parameters_of_the_potential':{ 3039 '': object}, 3040 'l0': { 3041 '': float}, 3042 } 3043 3044 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()])) 3045 3046 for level1, sub_dtypes in columns_dtypes.items(): 3047 for level2, dtype in sub_dtypes.items(): 3048 self.df[level1, level2] = self.df[level1, level2].astype(dtype) 3049 3050 columns_names = pd.MultiIndex.from_frame(self.df) 3051 columns_names = columns_names.names 3052 3053 return columns_names 3054 3055 def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True): 3056 """ 3057 Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`. 3058 3059 Args: 3060 espresso_system(`obj`): Instance of a system object from the espressomd library. 3061 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. 3062 combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'. 3063 warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True. 3064 3065 Note: 3066 - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 3067 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 3068 - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html 3069 3070 """ 3071 from itertools import combinations_with_replacement 3072 implemented_combining_rules = ['Lorentz-Berthelot'] 3073 compulsory_parameters_in_df = ['sigma','epsilon'] 3074 # Sanity check 3075 if combining_rule not in implemented_combining_rules: 3076 raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules) 3077 if shift_potential: 3078 shift="auto" 3079 else: 3080 shift=0 3081 # List which particles have sigma and epsilon values defined in pmb.df and which ones don't 3082 particles_types_with_LJ_parameters = [] 3083 non_parametrized_labels= [] 3084 for particle_type in self.get_type_map().values(): 3085 check_list=[] 3086 for key in compulsory_parameters_in_df: 3087 value_in_df=self.find_value_from_es_type(es_type=particle_type, 3088 column_name=key) 3089 check_list.append(np.isnan(value_in_df)) 3090 if any(check_list): 3091 non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 3092 column_name='label')) 3093 else: 3094 particles_types_with_LJ_parameters.append(particle_type) 3095 lj_parameters_keys=["label","sigma","epsilon","offset","cutoff"] 3096 # Set up LJ interactions between all particle types 3097 for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 3098 lj_parameters={} 3099 for key in lj_parameters_keys: 3100 lj_parameters[key]=[] 3101 # Search the LJ parameters of the type pair 3102 for ptype in type_pair: 3103 for key in lj_parameters_keys: 3104 lj_parameters[key].append(self.find_value_from_es_type(es_type=ptype, column_name=key)) 3105 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 3106 if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]): 3107 continue 3108 # Apply combining rule 3109 if combining_rule == 'Lorentz-Berthelot': 3110 sigma=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2 3111 cutoff=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2 3112 offset=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2 3113 epsilon=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1]) 3114 espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = epsilon.to('reduced_energy').magnitude, 3115 sigma = sigma.to('reduced_length').magnitude, 3116 cutoff = cutoff.to('reduced_length').magnitude, 3117 offset = offset.to("reduced_length").magnitude, 3118 shift = shift) 3119 index = len(self.df) 3120 self.df.at [index, 'name'] = f'LJ: {lj_parameters["label"][0]}-{lj_parameters["label"][1]}' 3121 lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params() 3122 3123 self.add_value_to_df(index=index, 3124 key=('pmb_type',''), 3125 new_value='LennardJones') 3126 3127 self.add_value_to_df(index=index, 3128 key=('parameters_of_the_potential',''), 3129 new_value=lj_params, 3130 non_standard_value=True) 3131 if non_parametrized_labels and warnings: 3132 print(f'WARNING: The following particles do not have a defined value of sigma or epsilon in pmb.df: {non_parametrized_labels}. No LJ interaction has been added in ESPResSo for those particles.') 3133 3134 return 3135 3136 def setup_particle_sigma(self, topology_dict): 3137 ''' 3138 Sets up sigma of the particles in `topology_dict`. 3139 3140 Args: 3141 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 3142 ''' 3143 for residue in topology_dict.keys(): 3144 residue_name = re.split(r'\d+', residue)[0] 3145 residue_number = re.split(r'(\d+)', residue)[1] 3146 residue_sigma = topology_dict[residue]['sigma'] 3147 sigma = residue_sigma*self.units.nm 3148 index = self.df[(self.df['residue_id']==residue_number) & (self.df['name']==residue_name) ].index.values[0] 3149 self.add_value_to_df(key= ('sigma',''), 3150 index=int (index), 3151 new_value=sigma) 3152 return 3153 3154 def write_pmb_df (self, filename): 3155 ''' 3156 Writes the pyMBE dataframe into a csv file 3157 3158 Args: 3159 filename (`str`): Path to the csv file 3160 ''' 3161 3162 self.df.to_csv(filename) 3163 3164 return 3165 3166 def write_output_vtf_file(self, espresso_system, filename): 3167 ''' 3168 Writes a snapshot of `espresso_system` on the vtf file `filename`. 3169 3170 Args: 3171 espresso_system(`obj`): Instance of a system object from the espressomd library. 3172 filename(`str`): Path to the vtf file. 3173 3174 ''' 3175 box = espresso_system.box_l[0] 3176 with open(filename, mode='w+t') as coordinates: 3177 coordinates.write (f'unitcell {box} {box} {box} \n') 3178 for particle in espresso_system.part: 3179 type_label = self.find_value_from_es_type(es_type=particle.type, column_name='label') 3180 coordinates.write (f'atom {particle.id} radius 1 name {type_label} type {type_label}\n' ) 3181 coordinates.write ('timestep indexed\n') 3182 for particle in espresso_system.part: 3183 coordinates.write (f'{particle.id} \t {particle.pos[0]} \t {particle.pos[1]} \t {particle.pos[2]}\n') 3184 return
12class pymbe_library(): 13 14 """ 15 The library for the Molecular Builder for ESPResSo (pyMBE) 16 17 Attributes: 18 N_A(`obj`): Avogadro number using the `pmb.units` UnitRegistry. 19 Kb(`obj`): Boltzmann constant using the `pmb.units` UnitRegistry. 20 e(`obj`): Elemental charge using the `pmb.units` UnitRegistry. 21 df(`obj`): PandasDataframe used to bookkeep all the information stored in pyMBE. Typically refered as `pmb.df`. 22 kT(`obj`): Thermal energy using the `pmb.units` UnitRegistry. It is used as the unit of reduced energy. 23 Kw(`obj`): Ionic product of water using the `pmb.units` UnitRegistry. Used in the setup of the G-RxMC method. 24 """ 25 units = pint.UnitRegistry() 26 N_A=6.02214076e23 / units.mol 27 Kb=1.38064852e-23 * units.J / units.K 28 e=1.60217662e-19 *units.C 29 df=None 30 kT=None 31 Kw=None 32 SEED=None 33 rng=None 34 35 36 class NumpyEncoder(json.JSONEncoder): 37 38 """ 39 Custom JSON encoder that converts NumPy arrays to Python lists 40 and NumPy scalars to Python scalars. 41 """ 42 43 def default(self, obj): 44 if isinstance(obj, np.ndarray): 45 return obj.tolist() 46 if isinstance(obj, np.generic): 47 return obj.item() 48 return super().default(obj) 49 50 51 def __init__(self, SEED, temperature=None, unit_length=None, unit_charge=None, Kw=None): 52 """ 53 Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 54 and sets up the `pmb.df` for bookkeeping. 55 56 Args: 57 temperature(`obj`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. 58 unit_length(`obj`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. 59 unit_charge (`obj`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 60 Kw (`obj`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 61 62 Note: 63 - If no `temperature` is given, a value of 298.15 K is assumed by default. 64 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 65 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 66 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 67 """ 68 # Seed and RNG 69 self.SEED=SEED 70 self.rng = np.random.default_rng(SEED) 71 self.set_reduced_units(unit_length=unit_length, unit_charge=unit_charge, 72 temperature=temperature, Kw=Kw, verbose=False) 73 self.setup_df() 74 return 75 76 def enable_motion_of_rigid_object(self, name, espresso_system): 77 ''' 78 Enables the motion of the rigid object `name` in the `espresso_system`. 79 80 Args: 81 name(`str`): Label of the object. 82 espresso_system(`obj`): Instance of a system object from the espressomd library. 83 84 Note: 85 - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]. 86 ''' 87 print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]') 88 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 89 if pmb_type != 'protein': 90 raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein') 91 molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list() 92 for molecule_id in molecule_ids_list: 93 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 94 center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system) 95 rigid_object_center = espresso_system.part.add(pos=center_of_mass, 96 rotation=[True,True,True], 97 type=self.propose_unused_type()) 98 for particle_id in particle_ids_list: 99 pid = espresso_system.part.by_id(particle_id) 100 pid.vs_auto_relate_to(rigid_object_center.id) 101 return 102 103 def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False): 104 """ 105 Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles. 106 107 Args: 108 particle_id1 (`int`): particle_id of the type of the first particle type of the bonded particles 109 particle_id2 (`int`): particle_id of the type of the second particle type of the bonded particles 110 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. 111 """ 112 particle_name1 = self.df.loc[self.df['particle_id']==particle_id1].name.values[0] 113 particle_name2 = self.df.loc[self.df['particle_id']==particle_id2].name.values[0] 114 bond_key = self.find_bond_key(particle_name1=particle_name1, 115 particle_name2=particle_name2, 116 use_default_bond=use_default_bond) 117 if not bond_key: 118 return 119 self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1) 120 indexs = np.where(self.df['name']==bond_key) 121 index_list = list (indexs[0]) 122 used_bond_df = self.df.loc[self.df['particle_id2'].notnull()] 123 #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict 124 used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 ) 125 used_bond_index = used_bond_df.index.to_list() 126 127 for index in index_list: 128 if index not in used_bond_index: 129 self.clean_df_row(index=int(index)) 130 self.df.at[index,'particle_id'] = particle_id1 131 self.df.at[index,'particle_id2'] = particle_id2 132 break 133 return 134 135 def add_bonds_to_espresso(self, espresso_system) : 136 """ 137 Adds all bonds defined in `pmb.df` to `espresso_system`. 138 139 Args: 140 espresso_system (str): system object of espressomd library 141 """ 142 143 if 'bond' in self.df.values: 144 bond_df = self.df.loc[self.df ['pmb_type'] == 'bond'] 145 bond_list = bond_df.bond_object.values.tolist() 146 for bond in bond_list: 147 espresso_system.bonded_inter.add(bond) 148 else: 149 print ('WARNING: There are no bonds defined in pymbe.df') 150 151 return 152 153 def add_value_to_df(self,index,key,new_value, verbose=True, non_standard_value=False): 154 """ 155 Adds a value to a cell in the `pmb.df` DataFrame. 156 157 Args: 158 index(`int`): index of the row to add the value to. 159 key(`str`): the column label to add the value to. 160 verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 161 non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False. 162 """ 163 164 token = "#protected:" 165 166 def protect(obj): 167 if non_standard_value: 168 return token + json.dumps(obj, cls=self.NumpyEncoder) 169 return obj 170 171 def deprotect(obj): 172 if non_standard_value and isinstance(obj, str) and obj.startswith(token): 173 return json.loads(obj.removeprefix(token)) 174 return obj 175 176 # Make sure index is a scalar integer value 177 index = int (index) 178 assert isinstance(index, int), '`index` should be a scalar integer value.' 179 idx = pd.IndexSlice 180 if self.check_if_df_cell_has_a_value(index=index,key=key): 181 old_value= self.df.loc[index,idx[key]] 182 if verbose: 183 if protect(old_value) != protect(new_value): 184 name=self.df.loc[index,('name','')] 185 pmb_type=self.df.loc[index,('pmb_type','')] 186 print(f"WARNING: you are redefining the properties of {name} of pmb_type {pmb_type}") 187 print(f'WARNING: overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}') 188 self.df.loc[index,idx[key]] = protect(new_value) 189 if non_standard_value: 190 self.df[key] = self.df[key].apply(deprotect) 191 return 192 193 def assign_molecule_id(self, name, molecule_index, pmb_type, used_molecules_id): 194 """ 195 Assigns the `molecule_id` of the pmb object given by `pmb_type` 196 197 Args: 198 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 199 pmb_type(`str`): pmb_object_type to assign the `molecule_id` 200 molecule_index (`int`): index of the current `pmb_object_type` to assign the `molecule_id` 201 used_molecules_id (`lst`): list with the `molecule_id` values already used. 202 203 Returns: 204 molecule_id(`int`): Id of the molecule 205 """ 206 207 self.clean_df_row(index=int(molecule_index)) 208 209 if self.df['molecule_id'].isnull().values.all(): 210 molecule_id = 0 211 else: 212 # check if a residue is part of another molecule 213 check_residue_name = self.df[self.df['residue_list'].astype(str).str.contains(name)] 214 mol_pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 215 if not check_residue_name.empty and mol_pmb_type == pmb_type: 216 for value in check_residue_name.index.to_list(): 217 if value not in used_molecules_id: 218 molecule_id = self.df.loc[value].molecule_id.values[0] 219 break 220 else: 221 molecule_id = self.df['molecule_id'].max() +1 222 223 self.add_value_to_df (key=('molecule_id',''), 224 index=int(molecule_index), 225 new_value=molecule_id, 226 verbose=False) 227 228 return molecule_id 229 230 def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system): 231 """ 232 Calculates the center of the molecule with a given molecule_id. 233 234 Args: 235 molecule_id(`int`): Id of the molecule whose center of mass is to be calculated. 236 espresso_system(`obj`): Instance of a system object from the espressomd library. 237 238 Returns: 239 center_of_mass(`lst`): Coordinates of the center of mass. 240 """ 241 total_beads = 0 242 center_of_mass = np.zeros(3) 243 axis_list = [0,1,2] 244 particle_id_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 245 for pid in particle_id_list: 246 total_beads +=1 247 for axis in axis_list: 248 center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis] 249 center_of_mass = center_of_mass /total_beads 250 return center_of_mass 251 252 def calculate_HH(self, molecule_name, pH_list=None, pka_set=None): 253 """ 254 Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 255 for molecules with the name `molecule_name`. 256 257 Args: 258 molecule_name (`str`): name of the molecule to calculate the ideal charge for 259 pH_list(`lst`): pH-values to calculate. 260 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 261 262 Returns: 263 Z_HH (`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list` 264 265 Note: 266 - This function supports objects with pmb types: "molecule", "peptide" and "protein". 267 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 268 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 269 - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan` should be used. 270 """ 271 if pH_list is None: 272 pH_list=np.linspace(2,12,50) 273 if pka_set is None: 274 pka_set=self.get_pka_set() 275 self.check_pka_set(pka_set=pka_set) 276 charge_map = self.get_charge_map() 277 Z_HH=[] 278 for pH_value in pH_list: 279 Z=0 280 index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 281 residue_list = self.df.at [index,('residue_list','')] 282 sequence = self.df.at [index,('sequence','')] 283 if np.any(pd.isnull(sequence)): 284 # Molecule has no sequence 285 for residue in residue_list: 286 list_of_particles_in_residue = self.search_particles_in_residue(residue) 287 for particle in list_of_particles_in_residue: 288 if particle in pka_set.keys(): 289 if pka_set[particle]['acidity'] == 'acidic': 290 psi=-1 291 elif pka_set[particle]['acidity']== 'basic': 292 psi=+1 293 else: 294 psi=0 295 Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value']))) 296 Z_HH.append(Z) 297 else: 298 # Molecule has a sequence 299 if not isinstance(sequence, list): 300 # If the df has been read by file, the sequence needs to be parsed. 301 sequence = self.parse_sequence_from_file(sequence=sequence) 302 for name in sequence: 303 if name in pka_set.keys(): 304 if pka_set[name]['acidity'] == 'acidic': 305 psi=-1 306 elif pka_set[name]['acidity']== 'basic': 307 psi=+1 308 else: 309 psi=0 310 Z+=psi/(1+10**(psi*(pH_value-pka_set[name]['pka_value']))) 311 else: 312 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 313 Z+=charge_map[state_one_type] 314 Z_HH.append(Z) 315 316 return Z_HH 317 318 def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None): 319 """ 320 Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning. 321 322 Args: 323 c_macro ('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 324 c_salt ('float'): Salt concentration in the reservoir. 325 pH_list ('lst'): List of pH-values in the reservoir. 326 pka_set ('dict'): {"name": {"pka_value": pka, "acidity": acidity}}. 327 328 Returns: 329 {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 330 pH_system_list ('lst'): List of pH_values in the system. 331 partition_coefficients_list ('lst'): List of partition coefficients of cations. 332 333 Note: 334 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 335 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 336 - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result. 337 """ 338 if pH_list is None: 339 pH_list=np.linspace(2,12,50) 340 if pka_set is None: 341 pka_set=self.get_pka_set() 342 self.check_pka_set(pka_set=pka_set) 343 344 partition_coefficients_list = [] 345 pH_system_list = [] 346 Z_HH_Donnan={} 347 for key in c_macro: 348 Z_HH_Donnan[key] = [] 349 350 def calc_charges(c_macro, pH): 351 """ 352 Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation. 353 354 Args: 355 c_macro ('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 356 pH ('float'): pH-value that is used in the HH equation. 357 358 Returns: 359 charge ('dict'): {"molecule_name": charge} 360 """ 361 charge = {} 362 for name in c_macro: 363 charge[name] = self.calculate_HH(name, [pH], pka_set)[0] 364 return charge 365 366 def calc_partition_coefficient(charge, c_macro): 367 """ 368 Calculates the partition coefficients of positive ions according to the ideal Donnan theory. 369 370 Args: 371 charge ('dict'): {"molecule_name": charge} 372 c_macro ('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 373 """ 374 nonlocal ionic_strength_res 375 charge_density = 0.0 376 for key in charge: 377 charge_density += charge[key] * c_macro[key] 378 return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude 379 380 for pH_value in pH_list: 381 # calculate the ionic strength of the reservoir 382 if pH_value <= 7.0: 383 ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 384 elif pH_value > 7.0: 385 ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt 386 387 #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations 388 #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation. 389 #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects. 390 equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro)) 391 logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq") 392 partition_coefficient = 10**logxi.root 393 394 charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient)) 395 for key in c_macro: 396 Z_HH_Donnan[key].append(charges_temp[key]) 397 398 pH_system_list.append(pH_value - np.log10(partition_coefficient)) 399 partition_coefficients_list.append(partition_coefficient) 400 401 return {"charges_dict": Z_HH_Donnan, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 402 403 def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff): 404 """ 405 Calculates the initial bond length that is used when setting up molecules, 406 based on the minimum of the sum of bonded and short-range (LJ) interactions. 407 408 Args: 409 bond_object (`cls`): instance of a bond object from espressomd library 410 bond_type (`str`): label identifying the used bonded potential 411 epsilon (`float`): LJ epsilon of the interaction between the particles 412 sigma (`float`): LJ sigma of the interaction between the particles 413 cutoff (`float`): cutoff-radius of the LJ interaction 414 """ 415 def truncated_lj_potential(x, epsilon, sigma, cutoff): 416 if x>cutoff: 417 return 0.0 418 else: 419 return 4*epsilon*((sigma/x)**12-(sigma/x)**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6) 420 421 epsilon_red=epsilon.to('reduced_energy').magnitude 422 sigma_red=sigma.to('reduced_length').magnitude 423 cutoff_red=cutoff.to('reduced_length').magnitude 424 425 if bond_type == "harmonic": 426 r_0 = bond_object.params.get('r_0') 427 k = bond_object.params.get('k') 428 l0 = scipy.optimize.minimize(lambda x: 0.5*k*(x-r_0)**2 + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red), x0=r_0).x 429 elif bond_type == "FENE": 430 r_0 = bond_object.params.get('r_0') 431 k = bond_object.params.get('k') 432 d_r_max = bond_object.params.get('d_r_max') 433 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), x0=1.0).x 434 return l0 435 436 def calculate_net_charge (self, espresso_system, molecule_name): 437 ''' 438 Calculates the net charge per molecule of molecules with `name` = molecule_name. 439 Returns the net charge per molecule and a maps with the net charge per residue and molecule. 440 441 Args: 442 espresso_system: system information 443 molecule_name (str): name of the molecule to calculate the net charge 444 445 Returns: 446 {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }} 447 448 Note: 449 - The net charge of the molecule is averaged over all molecules of type `name` 450 - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name` 451 ''' 452 valid_pmb_types = ["molecule", "protein"] 453 pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0] 454 if pmb_type not in valid_pmb_types: 455 raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}") 456 457 id_map = self.get_particle_id_map(object_name=molecule_name) 458 def create_charge_map(espresso_system,id_map,label): 459 charge_map={} 460 for super_id in id_map[label].keys(): 461 net_charge=0 462 for pid in id_map[label][super_id]: 463 net_charge+=espresso_system.part.by_id(pid).q 464 charge_map[super_id]=net_charge 465 return charge_map 466 net_charge_molecules=create_charge_map(label="molecule_map", 467 espresso_system=espresso_system, 468 id_map=id_map) 469 net_charge_residues=create_charge_map(label="residue_map", 470 espresso_system=espresso_system, 471 id_map=id_map) 472 mean_charge=np.mean(np.array(list(net_charge_molecules.values()))) 473 return {"mean": mean_charge, "molecules": net_charge_molecules, "residues": net_charge_residues} 474 475 def center_molecule_in_simulation_box(self, molecule_id, espresso_system): 476 """ 477 Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`. 478 479 Args: 480 molecule_id(`int`): Id of the molecule to be centered. 481 espresso_system(`obj`): Instance of a system object from the espressomd library. 482 """ 483 if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0: 484 raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.") 485 center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system) 486 box_center = [espresso_system.box_l[0]/2.0, 487 espresso_system.box_l[1]/2.0, 488 espresso_system.box_l[2]/2.0] 489 particle_id_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 490 for pid in particle_id_list: 491 es_pos = espresso_system.part.by_id(pid).pos 492 espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center 493 return 494 495 def check_dimensionality(self, variable, expected_dimensionality): 496 """ 497 Checks if the dimensionality of `variable` matches `expected_dimensionality`. 498 499 Args: 500 `variable`(`pint.Quantity`): Quantity to be checked. 501 `expected_dimensionality(`str`): Expected dimension of the variable. 502 503 Returns: 504 `bool`: `True` if the variable if of the expected dimensionality, `False` otherwise. 505 506 Note: 507 - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality). 508 - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"` 509 """ 510 correct_dimensionality=variable.check(f"{expected_dimensionality}") 511 if not correct_dimensionality: 512 raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}") 513 return correct_dimensionality 514 515 def check_if_df_cell_has_a_value(self, index,key): 516 """ 517 Checks if a cell in the `pmb.df` at the specified index and column has a value. 518 519 Args: 520 index(`int`): Index of the row to check. 521 key(`str`): Column label to check. 522 523 Returns: 524 `bool`: `True` if the cell has a value, `False` otherwise. 525 """ 526 idx = pd.IndexSlice 527 import warnings 528 with warnings.catch_warnings(): 529 warnings.simplefilter("ignore") 530 return not pd.isna(self.df.loc[index, idx[key]]) 531 532 def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined): 533 """ 534 Checks if `name` is defined in `pmb.df`. 535 536 Args: 537 name(`str`): label to check if defined in `pmb.df`. 538 pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`. 539 540 Returns: 541 `bool`: `True` for success, `False` otherwise. 542 """ 543 if name in self.df['name'].unique(): 544 current_object_type = self.df[self.df['name']==name].pmb_type.values[0] 545 if current_object_type != pmb_type_to_be_defined: 546 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") 547 return True 548 else: 549 return False 550 551 552 def check_pka_set(self, pka_set): 553 """ 554 Checks that `pka_set` has the formatting expected by the pyMBE library. 555 556 Args: 557 pka_set (dict): {"name" : {"pka_value": pka, "acidity": acidity}} 558 """ 559 required_keys=['pka_value','acidity'] 560 for required_key in required_keys: 561 for pka_entry in pka_set.values(): 562 if required_key not in pka_entry.keys(): 563 raise ValueError(f'missing a required key "{required_key}" in the following entry of pka_set "{pka_entry}"') 564 return 565 566 def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")): 567 """ 568 Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a np.nan value. 569 570 Args: 571 index(`int`): Index of the row to clean. 572 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`]. 573 """ 574 for column_key in columns_keys_to_clean: 575 self.add_value_to_df(key=(column_key,''),index=index,new_value=np.nan, verbose=False) 576 return 577 578 def convert_columns_to_original_format (self, df): 579 """ 580 Converts the columns of the Dataframe to the original format in pyMBE. 581 582 Args: 583 df(`DataFrame`): dataframe with pyMBE information as a string 584 585 """ 586 587 columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', 'model',('state_one','es_type'),('state_two','es_type'),('state_one','charge'),('state_two','charge') ] 588 589 columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset'] 590 591 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence'] 592 593 for column_name in columns_dtype_int: 594 df[column_name] = df[column_name].astype(object) 595 596 for column_name in columns_with_list_or_dict: 597 if df[column_name].isnull().all(): 598 df[column_name] = df[column_name].astype(object) 599 else: 600 df[column_name] = df[column_name].apply(lambda x: ast.literal_eval(str(x)) if pd.notnull(x) else x) 601 602 for column_name in columns_with_units: 603 df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x) 604 605 df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x) 606 607 return df 608 609 def convert_str_to_bond_object (self, bond_str): 610 611 """ 612 Convert a row read as a `str` to the corresponding bond object. There are two supported bonds: HarmonicBond and FeneBond 613 614 Args: 615 bond_str (`str`): string with the information of a bond object 616 617 Returns: 618 bond_object(`obj`): EsPRESSo bond object 619 """ 620 621 from espressomd.interactions import HarmonicBond 622 from espressomd.interactions import FeneBond 623 624 supported_bonds = ['HarmonicBond', 'FeneBond'] 625 626 for bond in supported_bonds: 627 variable = re.subn(f'{bond}', '', bond_str) 628 if variable[1] == 1: 629 params = ast.literal_eval(variable[0]) 630 if bond == 'HarmonicBond': 631 bond_object = HarmonicBond(r_cut =params['r_cut'], k = params['k'], r_0=params['r_0']) 632 elif bond == 'FeneBond': 633 bond_object = FeneBond(k = params['k'], d_r_max =params['d_r_max'], r_0=params['r_0']) 634 635 return bond_object 636 637 def copy_df_entry(self, name, column_name, number_of_copies): 638 ''' 639 Creates 'number_of_copies' of a given 'name' in `pymbe.df`. 640 641 Args: 642 name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df` 643 column_name(`str`): Column name to use as a filter. 644 number_of_copies(`int`): number of copies of `name` to be created. 645 646 Note: 647 - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 648 ''' 649 650 valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ] 651 if column_name not in valid_column_names: 652 raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}") 653 df_by_name = self.df.loc[self.df.name == name] 654 if number_of_copies != 1: 655 if df_by_name[column_name].isnull().values.any(): 656 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True) 657 else: 658 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 659 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 660 df_by_name_repeated[column_name] =np.NaN 661 # Concatenate the new particle rows to `df` 662 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 663 else: 664 if not df_by_name[column_name].isnull().values.any(): 665 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 666 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 667 df_by_name_repeated[column_name] =np.NaN 668 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 669 return 670 671 def create_added_salt (self, espresso_system, cation_name, anion_name, c_salt, verbose=True): 672 """ 673 Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`. 674 675 Args: 676 espresso_system (`obj`): instance of an espresso system object. 677 cation_name(`str`): `name` of a particle with a positive charge. 678 anion_name(`str`): `name` of a particle with a negative charge. 679 c_salt(`float`): Salt concentration. 680 verbose(`bool`): switch to activate/deactivate verbose. Defaults to True. 681 682 Returns: 683 c_salt_calculated (float): Calculated salt concentration added to `espresso_system`. 684 """ 685 cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.values[0] 686 anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.values[0] 687 if cation_name_charge <= 0: 688 raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge) 689 if anion_name_charge >= 0: 690 raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge) 691 # Calculate the number of ions in the simulation box 692 volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3') 693 if c_salt.check('[substance] [length]**-3'): 694 N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude) 695 c_salt_calculated=N_ions/(volume*self.N_A) 696 elif c_salt.check('[length]**-3'): 697 N_ions= int((volume*c_salt.to('reduced_length**-3')*self.N_A).magnitude) 698 c_salt_calculated=N_ions/volume 699 else: 700 raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt) 701 N_cation = N_ions*abs(anion_name_charge) 702 N_anion = N_ions*abs(cation_name_charge) 703 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation) 704 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion) 705 if verbose: 706 print(f"\n Added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions") 707 return c_salt_calculated 708 709 def create_bond_in_espresso(self, bond_type, bond_parameters): 710 ''' 711 Creates either a harmonic or a FENE bond in ESPREesSo 712 713 Args: 714 bond_type (`str`) : label to identify the potential to model the bond. 715 bond_parameters (`dict`) : parameters of the potential of the bond. 716 717 Note: 718 Currently, only HARMONIC and FENE bonds are supported. 719 720 For a HARMONIC bond the dictionary must contain: 721 722 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 723 using the `pmb.units` UnitRegistry. 724 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 725 the `pmb.units` UnitRegistry. 726 727 For a FENE bond the dictionay must additionally contain: 728 729 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 730 units of length using the `pmb.units` UnitRegistry. Default 'None'. 731 732 Returns: 733 bond_object (`obj`): a harmonic or a FENE bond object in ESPREesSo 734 ''' 735 from espressomd import interactions 736 737 valid_bond_types = ["harmonic", "FENE"] 738 739 if 'k' in bond_parameters: 740 bond_magnitude = bond_parameters['k'].to('reduced_energy / reduced_length**2') 741 else: 742 raise ValueError("Magnitud of the potential (k) is missing") 743 744 if bond_type == 'harmonic': 745 if 'r_0' in bond_parameters: 746 bond_length = bond_parameters['r_0'].to('reduced_length') 747 else: 748 raise ValueError("Equilibrium bond length (r_0) is missing") 749 bond_object = interactions.HarmonicBond(k = bond_magnitude.magnitude, 750 r_0 = bond_length.magnitude) 751 elif bond_type == 'FENE': 752 if 'r_0' in bond_parameters: 753 bond_length = bond_parameters['r_0'].to('reduced_length').magnitude 754 else: 755 print("WARNING: No value provided for r_0. Defaulting to r_0 = 0") 756 bond_length=0 757 if 'd_r_max' in bond_parameters: 758 max_bond_stret = bond_parameters['d_r_max'].to('reduced_length') 759 else: 760 raise ValueError("Maximal stretching length (d_r_max) is missing") 761 bond_object = interactions.FeneBond(r_0 = bond_length, 762 k = bond_magnitude.magnitude, 763 d_r_max = max_bond_stret.magnitude) 764 else: 765 raise ValueError(f"Bond type {bond_type} currently not implemented in pyMBE, valid types are {valid_bond_types}") 766 767 return bond_object 768 769 770 def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True): 771 """ 772 Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`. 773 774 Args: 775 object_name(`str`): `name` of a pymbe object. 776 espresso_system(`obj`): Instance of a system object from the espressomd library. 777 cation_name(`str`): `name` of a particle with a positive charge. 778 anion_name(`str`): `name` of a particle with a negative charge. 779 verbose(`bool`): switch to activate/deactivate verbose. Defaults to True. 780 781 Returns: 782 counterion_number (dict): {"name": number} 783 """ 784 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.iloc[0] 785 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.iloc[0] 786 object_ids = self.get_particle_id_map(object_name=object_name)["all"] 787 counterion_number={} 788 object_charge={} 789 for name in ['positive', 'negative']: 790 object_charge[name]=0 791 for id in object_ids: 792 if espresso_system.part.by_id(id).q > 0: 793 object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q )) 794 elif espresso_system.part.by_id(id).q < 0: 795 object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q )) 796 if object_charge['positive'] % abs(anion_charge) == 0: 797 counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge)) 798 else: 799 raise ValueError('The number of positive charges in the pmb_object must be divisible by the charge of the anion') 800 if object_charge['negative'] % abs(cation_charge) == 0: 801 counterion_number[cation_name]=int(object_charge['negative']/cation_charge) 802 else: 803 raise ValueError('The number of negative charges in the pmb_object must be divisible by the charge of the cation') 804 if counterion_number[cation_name] > 0: 805 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name]) 806 else: 807 counterion_number[cation_name]=0 808 if counterion_number[anion_name] > 0: 809 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name]) 810 else: 811 counterion_number[anion_name] = 0 812 if verbose: 813 print('The following counter-ions have been created: ') 814 for name in counterion_number.keys(): 815 print(f'Ion type: {name} created number: {counterion_number[name]}') 816 return counterion_number 817 818 def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, use_default_bond=False): 819 """ 820 Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 821 822 Args: 823 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 824 espresso_system(`obj`): Instance of a system object from espressomd library. 825 number_of_molecules(`int`): Number of molecules of type `name` to be created. 826 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 827 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df` 828 Returns: 829 830 molecules_info (`dict`): {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 831 """ 832 if list_of_first_residue_positions is not None: 833 for item in list_of_first_residue_positions: 834 if not isinstance(item, list): 835 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.") 836 elif len(item) != 3: 837 raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.") 838 839 if len(list_of_first_residue_positions) != number_of_molecules: 840 raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}") 841 if number_of_molecules <= 0: 842 return 843 if not self.check_if_name_is_defined_in_df(name=name, 844 pmb_type_to_be_defined='molecule'): 845 raise ValueError(f"{name} must correspond to a label of a pmb_type='molecule' defined on df") 846 first_residue = True 847 molecules_info = {} 848 residue_list = self.df[self.df['name']==name].residue_list.values [0] 849 850 self.copy_df_entry(name=name, 851 column_name='molecule_id', 852 number_of_copies=number_of_molecules) 853 854 molecules_index = np.where(self.df['name']==name) 855 molecule_index_list =list(molecules_index[0])[-number_of_molecules:] 856 used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist() 857 pos_index = 0 858 for molecule_index in molecule_index_list: 859 molecule_id = self.assign_molecule_id(name=name,pmb_type='molecule',used_molecules_id=used_molecules_id,molecule_index=molecule_index) 860 molecules_info[molecule_id] = {} 861 for residue in residue_list: 862 if first_residue: 863 if list_of_first_residue_positions is None: 864 residue_position = None 865 else: 866 for item in list_of_first_residue_positions: 867 residue_position = [np.array(list_of_first_residue_positions[pos_index])] 868 # Generate an arbitrary random unit vector 869 backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0], 870 radius=1, 871 n_samples=1, 872 on_surface=True)[0] 873 residues_info = self.create_residue(name=residue, 874 number_of_residues=1, 875 espresso_system=espresso_system, 876 central_bead_position=residue_position, 877 use_default_bond= use_default_bond, 878 backbone_vector=backbone_vector) 879 residue_id = next(iter(residues_info)) 880 for index in self.df[self.df['residue_id']==residue_id].index: 881 self.add_value_to_df(key=('molecule_id',''), 882 index=int (index), 883 new_value=molecule_id) 884 central_bead_id = residues_info[residue_id]['central_bead_id'] 885 previous_residue = residue 886 residue_position = espresso_system.part.by_id(central_bead_id).pos 887 previous_residue_id = central_bead_id 888 first_residue = False 889 else: 890 previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0] 891 new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0] 892 bond = self.search_bond(particle_name1=previous_central_bead_name, 893 particle_name2=new_central_bead_name, 894 hard_check=True, 895 use_default_bond=use_default_bond) 896 l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 897 particle_name2=new_central_bead_name, 898 hard_check=True, 899 use_default_bond=use_default_bond) 900 residue_position = residue_position+backbone_vector*l0 901 residues_info = self.create_residue(name=residue, 902 number_of_residues=1, 903 espresso_system=espresso_system, 904 central_bead_position=[residue_position], 905 use_default_bond= use_default_bond, 906 backbone_vector=backbone_vector) 907 residue_id = next(iter(residues_info)) 908 for index in self.df[self.df['residue_id']==residue_id].index: 909 if not self.check_if_df_cell_has_a_value(index=index,key=('molecule_id','')): 910 self.df.at[index,'molecule_id'] = molecule_id 911 self.add_value_to_df(key=('molecule_id',''), 912 index=int (index), 913 new_value=molecule_id, 914 verbose=False) 915 central_bead_id = residues_info[residue_id]['central_bead_id'] 916 espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id)) 917 self.add_bond_in_df(particle_id1=central_bead_id, 918 particle_id2=previous_residue_id, 919 use_default_bond=use_default_bond) 920 previous_residue_id = central_bead_id 921 previous_residue = residue 922 molecules_info[molecule_id][residue_id] = residues_info[residue_id] 923 first_residue = True 924 pos_index+=1 925 926 return molecules_info 927 928 def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False): 929 """ 930 Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`. 931 932 Args: 933 name (`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 934 espresso_system (`cls`): Instance of a system object from the espressomd library. 935 number_of_particles (`int`): Number of particles to be created. 936 position (list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None. 937 fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False. 938 Returns: 939 created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`. 940 """ 941 if number_of_particles <=0: 942 return 943 if not self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'): 944 raise ValueError(f"{name} must correspond to a label of a pmb_type='particle' defined on df") 945 # Copy the data of the particle `number_of_particles` times in the `df` 946 self.copy_df_entry(name=name,column_name='particle_id',number_of_copies=number_of_particles) 947 # Get information from the particle type `name` from the df 948 q = self.df.loc[self.df['name']==name].state_one.charge.values[0] 949 es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 950 # Get a list of the index in `df` corresponding to the new particles to be created 951 index = np.where(self.df['name']==name) 952 index_list =list(index[0])[-number_of_particles:] 953 # Create the new particles into `espresso_system` 954 created_pid_list=[] 955 for index in range (number_of_particles): 956 df_index=int (index_list[index]) 957 self.clean_df_row(index=df_index) 958 if position is None: 959 particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l) 960 else: 961 particle_position = position[index] 962 if len(espresso_system.part.all()) == 0: 963 bead_id = 0 964 else: 965 bead_id = max (espresso_system.part.all().id) + 1 966 created_pid_list.append(bead_id) 967 968 if fix: 969 espresso_system.part.add (id=bead_id, pos = particle_position, type = es_type, q = q,fix =[fix,fix,fix]) 970 else: 971 espresso_system.part.add (id=bead_id, pos = particle_position, type = es_type, q = q) 972 self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id, verbose=False) 973 return created_pid_list 974 975 def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False): 976 """ 977 Creates all `particle`s associated to `pmb object` into `espresso` a number of times equal to `number_of_objects`. 978 979 Args: 980 name(`str`): Unique label of the `pmb object` to be created. 981 number_of_objects(`int`): Number of `pmb object`s to be created. 982 espresso_system(`obj`): Instance of an espresso system object from espressomd library. 983 position(`list`): Coordinates where the particles should be created. 984 use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`. 985 986 Note: 987 - 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. 988 """ 989 allowed_objects=['particle','residue','molecule'] 990 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 991 if pmb_type not in allowed_objects: 992 raise ValueError('Object type not supported, supported types are ', allowed_objects) 993 if pmb_type == 'particle': 994 self.create_particle(name=name, number_of_particles=number_of_objects, espresso_system=espresso_system, position=position) 995 elif pmb_type == 'residue': 996 self.create_residue(name=name, number_of_residues=number_of_objects, espresso_system=espresso_system, central_bead_position=position,use_default_bond=use_default_bond) 997 elif pmb_type == 'molecule': 998 self.create_molecule(name=name, number_of_molecules=number_of_objects, espresso_system=espresso_system, use_default_bond=use_default_bond, list_of_first_residue_positions=position) 999 return 1000 1001 def create_protein(self, name, number_of_proteins, espresso_system, topology_dict): 1002 """ 1003 Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions` 1004 1005 Args: 1006 name(`str`): Label of the protein to be created. 1007 espresso_system(`obj`): Instance of a system object from the espressomd library. 1008 number_of_proteins(`int`): Number of proteins to be created. 1009 positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}} 1010 """ 1011 1012 if number_of_proteins <=0: 1013 return 1014 if not self.check_if_name_is_defined_in_df(name=name, 1015 pmb_type_to_be_defined='protein'): 1016 raise ValueError(f"{name} must correspond to a name of a pmb_type='protein' defined on df") 1017 1018 self.copy_df_entry(name=name,column_name='molecule_id',number_of_copies=number_of_proteins) 1019 1020 protein_index = np.where(self.df['name']==name) 1021 protein_index_list =list(protein_index[0])[-number_of_proteins:] 1022 used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist() 1023 1024 box_half=espresso_system.box_l[0]/2.0 1025 1026 for molecule_index in protein_index_list: 1027 1028 molecule_id = self.assign_molecule_id (name=name,pmb_type='protein',used_molecules_id=used_molecules_id,molecule_index=molecule_index) 1029 1030 protein_center = self.generate_coordinates_outside_sphere(radius = 1, 1031 max_dist=box_half, 1032 n_samples=1, 1033 center=[box_half]*3)[0] 1034 1035 for residue in topology_dict.keys(): 1036 1037 residue_name = re.split(r'\d+', residue)[0] 1038 residue_number = re.split(r'(\d+)', residue)[1] 1039 residue_position = topology_dict[residue]['initial_pos'] 1040 position = residue_position + protein_center 1041 1042 particle_id = self.create_particle(name=residue_name, 1043 espresso_system=espresso_system, 1044 number_of_particles=1, 1045 position=[position], 1046 fix = True) 1047 1048 index = self.df[self.df['particle_id']==particle_id[0]].index.values[0] 1049 1050 self.add_value_to_df(key=('residue_id',''), 1051 index=int (index), 1052 new_value=residue_number) 1053 1054 self.add_value_to_df(key=('molecule_id',''), 1055 index=int (index), 1056 new_value=molecule_id) 1057 1058 return 1059 1060 def create_residue(self, name, espresso_system, number_of_residues, central_bead_position=None,use_default_bond=False, backbone_vector=None): 1061 """ 1062 Creates `number_of_residues` residues of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1063 1064 Args: 1065 name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df` 1066 espresso_system(`obj`): Instance of a system object from espressomd library. 1067 number_of_residue(`int`): Number of residues of type `name` to be created. 1068 central_bead_position(`list` of `float`): Position of the central bead. 1069 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` 1070 backbone_vector (`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`. 1071 1072 Returns: 1073 residues_info (`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1074 """ 1075 if number_of_residues <= 0: 1076 return 1077 if not self.check_if_name_is_defined_in_df(name=name, 1078 pmb_type_to_be_defined='residue'): 1079 raise ValueError(f"{name} must correspond to a label of a pmb_type='residue' defined on df") 1080 # Copy the data of a residue a `number_of_residues` times in the `df 1081 self.copy_df_entry(name=name, 1082 column_name='residue_id', 1083 number_of_copies=number_of_residues) 1084 residues_index = np.where(self.df['name']==name) 1085 residue_index_list =list(residues_index[0])[-number_of_residues:] 1086 # Internal bookkepping of the residue info (important for side-chain residues) 1087 # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1088 residues_info={} 1089 for residue_index in residue_index_list: 1090 self.clean_df_row(index=int(residue_index)) 1091 # Assign a residue_id 1092 if self.df['residue_id'].isnull().all(): 1093 residue_id=0 1094 else: 1095 residue_id = self.df['residue_id'].max() + 1 1096 self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id, verbose=False) 1097 # create the principal bead 1098 if self.df.loc[self.df['name']==name].central_bead.values[0] is np.NaN: 1099 raise ValueError("central_bead must contain a particle name") 1100 central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0] 1101 central_bead_id = self.create_particle(name=central_bead_name, 1102 espresso_system=espresso_system, 1103 position=central_bead_position, 1104 number_of_particles = 1)[0] 1105 central_bead_position=espresso_system.part.by_id(central_bead_id).pos 1106 #assigns same residue_id to the central_bead particle created. 1107 index = self.df[self.df['particle_id']==central_bead_id].index.values[0] 1108 self.df.at [index,'residue_id'] = residue_id 1109 # Internal bookkeeping of the central bead id 1110 residues_info[residue_id]={} 1111 residues_info[residue_id]['central_bead_id']=central_bead_id 1112 # create the lateral beads 1113 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1114 side_chain_beads_ids = [] 1115 for side_chain_element in side_chain_list: 1116 if side_chain_element not in self.df.values: 1117 raise ValueError (side_chain_element +'is not defined') 1118 pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 1119 if pmb_type == 'particle': 1120 bond = self.search_bond(particle_name1=central_bead_name, 1121 particle_name2=side_chain_element, 1122 hard_check=True, 1123 use_default_bond=use_default_bond) 1124 l0 = self.get_bond_length(particle_name1=central_bead_name, 1125 particle_name2=side_chain_element, 1126 hard_check=True, 1127 use_default_bond=use_default_bond) 1128 1129 if backbone_vector is None: 1130 bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1131 radius=l0, 1132 n_samples=1, 1133 on_surface=True)[0] 1134 else: 1135 bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1136 magnitude=l0) 1137 1138 side_bead_id = self.create_particle(name=side_chain_element, 1139 espresso_system=espresso_system, 1140 position=[bead_position], 1141 number_of_particles=1)[0] 1142 index = self.df[self.df['particle_id']==side_bead_id].index.values[0] 1143 self.add_value_to_df(key=('residue_id',''), 1144 index=int (index), 1145 new_value=residue_id, 1146 verbose=False) 1147 side_chain_beads_ids.append(side_bead_id) 1148 espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id)) 1149 self.add_bond_in_df(particle_id1=central_bead_id, 1150 particle_id2=side_bead_id, 1151 use_default_bond=use_default_bond) 1152 elif pmb_type == 'residue': 1153 central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0] 1154 bond = self.search_bond(particle_name1=central_bead_name, 1155 particle_name2=central_bead_side_chain, 1156 hard_check=True, 1157 use_default_bond=use_default_bond) 1158 l0 = self.get_bond_length(particle_name1=central_bead_name, 1159 particle_name2=central_bead_side_chain, 1160 hard_check=True, 1161 use_default_bond=use_default_bond) 1162 if backbone_vector is None: 1163 residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1164 radius=l0, 1165 n_samples=1, 1166 on_surface=True)[0] 1167 else: 1168 residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1169 magnitude=l0) 1170 lateral_residue_info = self.create_residue(name=side_chain_element, espresso_system=espresso_system, 1171 number_of_residues=1, central_bead_position=[residue_position],use_default_bond=use_default_bond) 1172 lateral_residue_dict=list(lateral_residue_info.values())[0] 1173 central_bead_side_chain_id=lateral_residue_dict['central_bead_id'] 1174 lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids'] 1175 residue_id_side_chain=list(lateral_residue_info.keys())[0] 1176 # Change the residue_id of the residue in the side chain to the one of the biger residue 1177 index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0] 1178 self.add_value_to_df(key=('residue_id',''), 1179 index=int(index), 1180 new_value=residue_id, 1181 verbose=False) 1182 # Change the residue_id of the particles in the residue in the side chain 1183 side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids 1184 for particle_id in side_chain_beads_ids: 1185 index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0] 1186 self.add_value_to_df(key=('residue_id',''), 1187 index=int (index), 1188 new_value=residue_id, 1189 verbose=False) 1190 espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id)) 1191 self.add_bond_in_df(particle_id1=central_bead_id, 1192 particle_id2=central_bead_side_chain_id, 1193 use_default_bond=use_default_bond) 1194 # Internal bookkeeping of the side chain beads ids 1195 residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids 1196 return residues_info 1197 1198 def create_variable_with_units(self, variable): 1199 """ 1200 Returns a pint object with the value and units defined in `variable`. 1201 1202 Args: 1203 variable(`dict` or `str`): {'value': value, 'units': units} 1204 Returns: 1205 variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry. 1206 """ 1207 1208 if isinstance(variable, dict): 1209 1210 value=variable.pop('value') 1211 units=variable.pop('units') 1212 1213 elif isinstance(variable, str): 1214 1215 value = float(re.split(r'\s+', variable)[0]) 1216 units = re.split(r'\s+', variable)[1] 1217 1218 variable_with_units=value*self.units(units) 1219 1220 return variable_with_units 1221 1222 def define_AA_particles_in_sequence(self, sequence, sigma_dict=None): 1223 ''' 1224 Defines in `pmb.df` all the different particles in `sequence`. 1225 1226 Args: 1227 sequence(`lst`): Sequence of the peptide or protein. 1228 1229 Note: 1230 - It assumes that the names in `sequence` correspond to amino acid names using the standard one letter code. 1231 ''' 1232 1233 already_defined_AA=[] 1234 1235 for residue_name in sequence: 1236 if residue_name in already_defined_AA: 1237 continue 1238 self.define_particle (name=residue_name, q=0) 1239 1240 if sigma_dict: 1241 self.define_particles_parameter_from_dict(param_dict = sigma_dict, 1242 param_name = 'sigma') 1243 return 1244 1245 def define_AA_residues(self, sequence, model): 1246 """ 1247 Defines in `pmb.df` all the different residues in `sequence`. 1248 1249 Args: 1250 sequence(`lst`): Sequence of the peptide or protein. 1251 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1252 1253 Returns: 1254 residue_list (`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1255 """ 1256 residue_list = [] 1257 for residue_name in sequence: 1258 if model == '1beadAA': 1259 central_bead = residue_name 1260 side_chains = [] 1261 elif model == '2beadAA': 1262 if residue_name in ['c','n', 'G']: 1263 central_bead = residue_name 1264 side_chains = [] 1265 else: 1266 central_bead = 'CA' 1267 side_chains = [residue_name] 1268 if residue_name not in residue_list: 1269 self.define_residue(name = 'AA-'+residue_name, 1270 central_bead = central_bead, 1271 side_chains = side_chains) 1272 residue_list.append('AA-'+residue_name) 1273 return residue_list 1274 1275 def define_bond(self, bond_type, bond_parameters, particle_pairs): 1276 1277 ''' 1278 Defines a pmb object of type `bond` in `pymbe.df`. 1279 1280 Args: 1281 bond_type (`str`) : label to identify the potential to model the bond. 1282 bond_parameters (`dict`) : parameters of the potential of the bond. 1283 particle_pairs (`lst`) : list of the `names` of the `particles` to be bonded. 1284 1285 Note: 1286 Currently, only HARMONIC and FENE bonds are supported. 1287 1288 For a HARMONIC bond the dictionary must contain the following parameters: 1289 1290 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 1291 using the `pmb.units` UnitRegistry. 1292 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 1293 the `pmb.units` UnitRegistry. 1294 1295 For a FENE bond the dictionay must contain the same parameters as for a HARMONIC bond and: 1296 1297 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 1298 units of length using the `pmb.units` UnitRegistry. Default 'None'. 1299 ''' 1300 1301 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1302 1303 for particle_name1, particle_name2 in particle_pairs: 1304 1305 lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1, 1306 particle_name2 = particle_name2, 1307 combining_rule = 'Lorentz-Berthelot') 1308 1309 cutoff = 2**(1.0/6.0) * lj_parameters["sigma"] 1310 1311 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1312 bond_type = bond_type, 1313 epsilon = lj_parameters["epsilon"], 1314 sigma = lj_parameters["sigma"], 1315 cutoff = cutoff) 1316 index = len(self.df) 1317 for label in [particle_name1+'-'+particle_name2,particle_name2+'-'+particle_name1]: 1318 self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond") 1319 self.df.at [index,'name']= particle_name1+'-'+particle_name2 1320 self.df.at [index,'bond_object'] = bond_object 1321 self.df.at [index,'l0'] = l0 1322 self.add_value_to_df(index=index, 1323 key=('pmb_type',''), 1324 new_value='bond') 1325 self.add_value_to_df(index=index, 1326 key=('parameters_of_the_potential',''), 1327 new_value=bond_object.get_params(), 1328 non_standard_value=True) 1329 1330 return 1331 1332 1333 def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None): 1334 """ 1335 Asigns `bond` in `pmb.df` as the default bond. 1336 The LJ parameters can be optionally provided to calculate the initial bond length 1337 1338 Args: 1339 bond_type (`str`) : label to identify the potential to model the bond. 1340 bond_parameters (`dict`) : parameters of the potential of the bond. 1341 sigma (`float`, optional) : LJ sigma of the interaction between the particles 1342 epsilon (`float`, optional): LJ epsilon for the interaction between the particles 1343 cutoff (`float`, optional) : cutoff-radius of the LJ interaction 1344 1345 Note: 1346 - Currently, only harmonic and FENE bonds are supported. 1347 """ 1348 1349 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1350 1351 if epsilon is None: 1352 epsilon=1*self.units('reduced_energy') 1353 if sigma is None: 1354 sigma=1*self.units('reduced_length') 1355 if cutoff is None: 1356 cutoff=2**(1.0/6.0)*self.units('reduced_length') 1357 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1358 bond_type = bond_type, 1359 epsilon = epsilon, 1360 sigma = sigma, 1361 cutoff = cutoff) 1362 1363 if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'): 1364 return 1365 if len(self.df.index) != 0: 1366 index = max(self.df.index)+1 1367 else: 1368 index = 0 1369 self.df.at [index,'name'] = 'default' 1370 self.df.at [index,'bond_object'] = bond_object 1371 self.df.at [index,'l0'] = l0 1372 self.add_value_to_df(index = index, 1373 key = ('pmb_type',''), 1374 new_value = 'bond') 1375 self.add_value_to_df(index = index, 1376 key = ('parameters_of_the_potential',''), 1377 new_value = bond_object.get_params(), 1378 non_standard_value=True) 1379 return 1380 1381 def define_epsilon_value_of_particles(self, eps_dict): 1382 ''' 1383 Defines the epsilon value of the particles in `eps_dict` into the `pmb.df`. 1384 1385 Args: 1386 eps_dict(`dict`): {'name': epsilon} 1387 ''' 1388 for residue in eps_dict.keys(): 1389 label_list = self.df[self.df['name'] == residue].index.tolist() 1390 for index in label_list: 1391 epsilon = eps_dict[residue] 1392 self.add_value_to_df(key= ('epsilon',''), 1393 index=int (index), 1394 new_value=epsilon) 1395 return 1396 1397 def define_molecule(self, name, residue_list): 1398 """ 1399 Defines a pyMBE object of type `molecule` in `pymbe.df`. 1400 1401 Args: 1402 name (`str`): Unique label that identifies the `molecule`. 1403 residue_list (`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1404 """ 1405 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'): 1406 return 1407 index = len(self.df) 1408 self.df.at [index,'name'] = name 1409 self.df.at [index,'pmb_type'] = 'molecule' 1410 self.df.at [index,('residue_list','')] = residue_list 1411 return 1412 1413 def define_particle(self, name, q=0, acidity='inert', pka=None, sigma=None, epsilon=None, cutoff=None, offset=None,verbose=True): 1414 """ 1415 Defines the properties of a particle object. 1416 1417 Args: 1418 name (`str`): Unique label that identifies this particle type. 1419 q (`int`, optional): Permanent charge of this particle type. Defaults to 0. 1420 acidity (`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to `inert`. 1421 pka (`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to None. 1422 sigma (`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1423 cutoff (`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1424 offset (`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1425 epsilon (`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to None. 1426 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 1427 1428 Note: 1429 - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units. 1430 - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units. 1431 - `cutoff` defaults to `2**(1./6.) reduced_length`. 1432 - `offset` defaults to 0. 1433 - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions. 1434 - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`. 1435 """ 1436 1437 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'): 1438 index = self.df[self.df['name']==name].index[0] 1439 else: 1440 index = len(self.df) 1441 self.df.at [index, 'name'] = name 1442 self.df.at [index,'pmb_type'] = 'particle' 1443 1444 # If `cutoff` and `offset` are not defined, default them to 0 1445 1446 if cutoff is None: 1447 cutoff=self.units.Quantity(2**(1./6.), "reduced_length") 1448 if offset is None: 1449 offset=self.units.Quantity(0, "reduced_length") 1450 1451 # Define LJ parameters 1452 parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"}, 1453 "cutoff":{"value": cutoff, "dimensionality": "[length]"}, 1454 "offset":{"value": offset, "dimensionality": "[length]"}, 1455 "epsilon":{"value": epsilon, "dimensionality": "[energy]"},} 1456 1457 for parameter_key in parameters_with_dimensionality.keys(): 1458 if parameters_with_dimensionality[parameter_key]["value"] is not None: 1459 self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 1460 expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"]) 1461 self.add_value_to_df(key=(parameter_key,''), 1462 index=index, 1463 new_value=parameters_with_dimensionality[parameter_key]["value"], 1464 verbose=verbose) 1465 1466 # Define particle acid/base properties 1467 self.set_particle_acidity (name=name, 1468 acidity = acidity , 1469 default_charge=q, 1470 pka=pka, 1471 verbose=verbose) 1472 return 1473 1474 def define_particles_parameter_from_dict(self, param_dict, param_name): 1475 ''' 1476 Defines the `param_name` value of the particles in `param_dict` into the `pmb.df`. 1477 1478 Args: 1479 param_dict(`dict`): {'name': `param_name` value} 1480 ''' 1481 for residue in param_dict.keys(): 1482 label_list = self.df[self.df['name'] == residue].index.tolist() 1483 for index in label_list: 1484 value = param_dict[residue] 1485 self.add_value_to_df(key= (param_name,''), 1486 index=int (index), 1487 new_value=value) 1488 return 1489 1490 def define_peptide(self, name, sequence, model): 1491 """ 1492 Defines a pyMBE object of type `peptide` in the `pymbe.df`. 1493 1494 Args: 1495 name (`str`): Unique label that identifies the `peptide`. 1496 sequence (`string`): Sequence of the `peptide`. 1497 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1498 """ 1499 if not self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'): 1500 valid_keys = ['1beadAA','2beadAA'] 1501 if model not in valid_keys: 1502 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1503 1504 clean_sequence = self.protein_sequence_parser(sequence=sequence) 1505 residue_list = self.define_AA_residues(sequence=clean_sequence,model=model) 1506 1507 self.define_molecule(name = name, residue_list=residue_list) 1508 index = self.df.loc[self.df['name'] == name].index.item() 1509 self.df.at [index,'model'] = model 1510 self.df.at [index,('sequence','')] = clean_sequence 1511 return 1512 1513 def define_protein(self, name,model, topology_dict): 1514 """ 1515 Defines a pyMBE object of type `protein` in `pymbe.df`. 1516 1517 Args: 1518 name (`str`): Unique label that identifies the `protein`. 1519 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1520 topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 1521 """ 1522 1523 protein_seq_list = [] 1524 1525 valid_keys = ['1beadAA','2beadAA'] 1526 1527 if model not in valid_keys: 1528 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1529 1530 if model == '2beadAA': 1531 self.define_particle(name='CA') 1532 1533 sigma_dict = {} 1534 1535 for residue in topology_dict.keys(): 1536 residue_name = re.split(r'\d+', residue)[0] 1537 1538 if residue_name not in sigma_dict.keys(): 1539 sigma_dict [residue_name] = topology_dict[residue]['sigma'] 1540 if residue_name not in ('CA', 'Ca'): 1541 protein_seq_list.append(residue_name) 1542 1543 protein_sequence = ''.join(protein_seq_list) 1544 clean_sequence = self.protein_sequence_parser(sequence=protein_sequence) 1545 1546 self.define_AA_particles_in_sequence (sequence=clean_sequence, sigma_dict = sigma_dict) 1547 residue_list = self.define_AA_residues(sequence=clean_sequence, model=model) 1548 1549 index = len(self.df) 1550 self.df.at [index,'name'] = name 1551 self.df.at [index,'pmb_type'] = 'protein' 1552 self.df.at [index,'model'] = model 1553 self.df.at [index,('sequence','')] = clean_sequence 1554 self.df.at [index,('residue_list','')] = residue_list 1555 1556 return 1557 1558 def define_residue(self, name, central_bead, side_chains): 1559 """ 1560 Defines a pyMBE object of type `residue` in `pymbe.df`. 1561 1562 Args: 1563 name (`str`): Unique label that identifies the `residue`. 1564 central_bead (`str`): `name` of the `particle` to be placed as central_bead of the `residue`. 1565 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. 1566 """ 1567 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'): 1568 return 1569 index = len(self.df) 1570 self.df.at [index, 'name'] = name 1571 self.df.at [index,'pmb_type'] = 'residue' 1572 self.df.at [index,'central_bead'] = central_bead 1573 self.df.at [index,('side_chains','')] = side_chains 1574 return 1575 1576 def destroy_pmb_object_in_system(self, name, espresso_system): 1577 """ 1578 Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 1579 1580 Args: 1581 name (str): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`. 1582 espresso_system (cls): Instance of a system class from espressomd library. 1583 1584 Note: 1585 - 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. 1586 """ 1587 allowed_objects = ['particle','residue','molecule'] 1588 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 1589 if pmb_type not in allowed_objects: 1590 raise ValueError('Object type not supported, supported types are ', allowed_objects) 1591 if pmb_type == 'particle': 1592 particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1593 & (self.df['molecule_id'].isna())] 1594 particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1595 & (self.df['molecule_id'].isna())].particle_id.tolist() 1596 for particle_id in particle_ids_list: 1597 espresso_system.part.by_id(particle_id).remove() 1598 self.df = self.df.drop(index=particles_index) 1599 if pmb_type == 'residue': 1600 residues_id = self.df.loc[self.df['name']== name].residue_id.to_list() 1601 for residue_id in residues_id: 1602 particle_ids_list = self.df.loc[self.df['residue_id']==residue_id].particle_id.dropna().to_list() 1603 self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index) 1604 for particle_id in particle_ids_list: 1605 espresso_system.part.by_id(particle_id).remove() 1606 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1607 if pmb_type == 'molecule': 1608 molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list() 1609 for molecule_id in molecules_id: 1610 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 1611 1612 self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index) 1613 for particle_id in particle_ids_list: 1614 espresso_system.part.by_id(particle_id).remove() 1615 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1616 1617 self.df.reset_index(drop=True,inplace=True) 1618 1619 return 1620 1621 def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200): 1622 """ 1623 Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration. 1624 To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an 1625 ion pair and the concentrations of the various species is solved numerically using a self-consistent approach. 1626 More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831). 1627 This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237). 1628 1629 Args: 1630 pH_res ('float'): pH-value in the reservoir. 1631 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 1632 activity_coefficient_monovalent_pair ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 1633 1634 Returns: 1635 cH_res ('float'): Concentration of H+ ions. 1636 cOH_res ('float'): Concentration of OH- ions. 1637 cNa_res ('float'): Concentration of Na+ ions. 1638 cCl_res ('float'): Concentration of Cl- ions. 1639 """ 1640 1641 self_consistent_run = 0 1642 cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+ 1643 1644 def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res): 1645 #Calculate and initial guess for the concentrations of various species based on ideal gas estimate 1646 cOH_res = self.Kw / cH_res 1647 cNa_res = None 1648 cCl_res = None 1649 if cOH_res>=cH_res: 1650 #adjust the concentration of sodium if there is excess OH- in the reservoir: 1651 cNa_res = c_salt_res + (cOH_res-cH_res) 1652 cCl_res = c_salt_res 1653 else: 1654 # adjust the concentration of chloride if there is excess H+ in the reservoir 1655 cCl_res = c_salt_res + (cH_res-cOH_res) 1656 cNa_res = c_salt_res 1657 1658 def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res): 1659 nonlocal max_number_sc_runs, self_consistent_run 1660 if self_consistent_run<max_number_sc_runs: 1661 self_consistent_run+=1 1662 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1663 cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res)) 1664 if cOH_res>=cH_res: 1665 #adjust the concentration of sodium if there is excess OH- in the reservoir: 1666 cNa_res = c_salt_res + (cOH_res-cH_res) 1667 cCl_res = c_salt_res 1668 else: 1669 # adjust the concentration of chloride if there is excess H+ in the reservoir 1670 cCl_res = c_salt_res + (cH_res-cOH_res) 1671 cNa_res = c_salt_res 1672 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 1673 else: 1674 return cH_res, cOH_res, cNa_res, cCl_res 1675 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 1676 1677 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 1678 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1679 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 1680 1681 while abs(determined_pH-pH_res)>1e-6: 1682 if determined_pH > pH_res: 1683 cH_res *= 1.005 1684 else: 1685 cH_res /= 1.003 1686 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 1687 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1688 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 1689 self_consistent_run=0 1690 1691 return cH_res, cOH_res, cNa_res, cCl_res 1692 1693 def filter_df(self, pmb_type): 1694 """ 1695 Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns. 1696 1697 Args: 1698 pmb_type(`str`): pmb_object_type to filter in `pmb.df`. 1699 1700 Returns: 1701 pmb_type_df(`obj`): filtered `pmb.df`. 1702 """ 1703 pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type] 1704 pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1) 1705 return pmb_type_df 1706 1707 def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False): 1708 """ 1709 Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 1710 1711 Args: 1712 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 1713 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 1714 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'. 1715 1716 Returns: 1717 bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` 1718 1719 Note: 1720 - If `use_default_bond`=`True`, it returns "default" if no key is found. 1721 """ 1722 bond_keys = [particle_name1 +'-'+ particle_name2, particle_name2 +'-'+ particle_name1 ] 1723 bond_defined=False 1724 for bond_key in bond_keys: 1725 if bond_key in self.df.values: 1726 bond_defined=True 1727 correct_key=bond_key 1728 break 1729 if bond_defined: 1730 return correct_key 1731 elif use_default_bond: 1732 return 'default' 1733 else: 1734 return 1735 1736 def find_value_from_es_type(self, es_type, column_name): 1737 """ 1738 Finds a value in `pmb.df` for a `column_name` and `es_type` pair. 1739 1740 Args: 1741 es_type(`int`): value of the espresso type 1742 column_name(`str`): name of the column in `pymbe.df` 1743 1744 Returns: 1745 Value in `pymbe.df` matching `column_name` and `es_type` 1746 """ 1747 idx = pd.IndexSlice 1748 for state in ['state_one', 'state_two']: 1749 index = np.where(self.df[(state, 'es_type')] == es_type)[0] 1750 if len(index) > 0: 1751 if column_name == 'label': 1752 label = self.df.loc[idx[index[0]], idx[(state,column_name)]] 1753 return label 1754 else: 1755 column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]] 1756 return column_name_value 1757 return None 1758 1759 def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples): 1760 """ 1761 Generates coordinates outside a sphere centered at `center`. 1762 1763 Args: 1764 center(`lst`): Coordinates of the center of the sphere. 1765 radius(`float`): Radius of the sphere. 1766 max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates. 1767 n_samples(`int`): Number of sample points. 1768 1769 Returns: 1770 coord_list(`lst`): Coordinates of the sample points. 1771 """ 1772 if not radius > 0: 1773 raise ValueError (f'The value of {radius} must be a positive value') 1774 if not radius < max_dist: 1775 raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))') 1776 coord_list = [] 1777 counter = 0 1778 while counter<n_samples: 1779 coord = self.generate_random_points_in_a_sphere(center=center, 1780 radius=max_dist, 1781 n_samples=1)[0] 1782 if np.linalg.norm(coord-np.asarray(center))>=radius: 1783 coord_list.append (coord) 1784 counter += 1 1785 return coord_list 1786 1787 def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False): 1788 """ 1789 Uniformly samples points from a hypersphere. If on_surface is set to True, the points are 1790 uniformly sampled from the surface of the hypersphere. 1791 1792 Args: 1793 center(`lst`): Array with the coordinates of the center of the spheres. 1794 radius(`float`): Radius of the sphere. 1795 n_samples(`int`): Number of sample points to generate inside the sphere. 1796 on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere. 1797 1798 Returns: 1799 samples(`list`): Coordinates of the sample points inside the hypersphere. 1800 1801 Note: 1802 - Algorithm from: https://baezortega.github.io/2018/10/14/hypersphere-sampling/ 1803 """ 1804 # initial values 1805 center=np.array(center) 1806 d = center.shape[0] 1807 # sample n_samples points in d dimensions from a standard normal distribution 1808 samples = self.rng.normal(size=(n_samples, d)) 1809 # make the samples lie on the surface of the unit hypersphere 1810 normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis] 1811 samples /= normalize_radii 1812 if not on_surface: 1813 # make the samples lie inside the hypersphere with the correct density 1814 uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis] 1815 new_radii = np.power(uniform_points, 1/d) 1816 samples *= new_radii 1817 # scale the points to have the correct radius and center 1818 samples = samples * radius + center 1819 return samples 1820 1821 def generate_trial_perpendicular_vector(self,vector,magnitude): 1822 """ 1823 Generates an orthogonal vector to the input `vector`. 1824 1825 Args: 1826 vector(`lst`): arbitrary vector. 1827 magnitude(`float`): magnitude of the orthogonal vector. 1828 1829 Returns: 1830 (`lst`): Orthogonal vector with the same magnitude as the input vector. 1831 """ 1832 np_vec = np.array(vector) 1833 np_vec /= np.linalg.norm(np_vec) 1834 if np.all(np_vec == 0): 1835 raise ValueError('Zero vector') 1836 # Generate a random vector 1837 random_vector = self.generate_random_points_in_a_sphere(radius=1, 1838 center=[0,0,0], 1839 n_samples=1, 1840 on_surface=True)[0] 1841 # Project the random vector onto the input vector and subtract the projection 1842 projection = np.dot(random_vector, np_vec) * np_vec 1843 perpendicular_vector = random_vector - projection 1844 # Normalize the perpendicular vector to have the same magnitude as the input vector 1845 perpendicular_vector /= np.linalg.norm(perpendicular_vector) 1846 return perpendicular_vector*magnitude 1847 1848 def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 1849 """ 1850 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length. 1851 If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead. 1852 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. 1853 1854 Args: 1855 particle_name1 (str): label of the type of the first particle type of the bonded particles. 1856 particle_name2 (str): label of the type of the second particle type of the bonded particles. 1857 hard_check (bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 1858 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. 1859 1860 Returns: 1861 l0 (`float`): bond length 1862 1863 Note: 1864 - 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`. 1865 - If `hard_check`=`True` stops the code when no bond is found. 1866 """ 1867 bond_key = self.find_bond_key(particle_name1=particle_name1, 1868 particle_name2=particle_name2, 1869 use_default_bond=use_default_bond) 1870 if bond_key: 1871 return self.df[self.df['name']==bond_key].l0.values[0] 1872 else: 1873 print("Bond not defined between particles ", particle_name1, " and ", particle_name2) 1874 if hard_check: 1875 sys.exit(1) 1876 else: 1877 return 1878 1879 def get_charge_map(self): 1880 ''' 1881 Gets the charge of each `espresso_type` in `pymbe.df`. 1882 1883 Returns: 1884 charge_map(`dict`): {espresso_type: charge}. 1885 ''' 1886 if self.df.state_one['es_type'].isnull().values.any(): 1887 df_state_one = self.df.state_one.dropna() 1888 df_state_two = self.df.state_two.dropna() 1889 else: 1890 df_state_one = self.df.state_one 1891 if self.df.state_two['es_type'].isnull().values.any(): 1892 df_state_two = self.df.state_two.dropna() 1893 else: 1894 df_state_two = self.df.state_two 1895 state_one = pd.Series (df_state_one.charge.values,index=df_state_one.es_type.values) 1896 state_two = pd.Series (df_state_two.charge.values,index=df_state_two.es_type.values) 1897 charge_map = pd.concat([state_one,state_two],axis=0).to_dict() 1898 return charge_map 1899 1900 def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'): 1901 """ 1902 Returns the Lennard-Jones parameters for the interaction between the particle types given by 1903 `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule. 1904 1905 Args: 1906 particle_name1 (str): label of the type of the first particle type 1907 particle_name2 (str): label of the type of the second particle type 1908 combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'. 1909 1910 Returns: 1911 {"epsilon": LJ epsilon, "sigma": LJ sigma} 1912 1913 Note: 1914 Currently, the only `combining_rule` supported is Lorentz-Berthelot. 1915 """ 1916 supported_combining_rules=["Lorentz-Berthelot"] 1917 1918 if combining_rule not in supported_combining_rules: 1919 raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}") 1920 1921 eps1 = self.df.loc[self.df["name"]==particle_name1].epsilon.iloc[0] 1922 sigma1 = self.df.loc[self.df["name"]==particle_name1].sigma.iloc[0] 1923 eps2 = self.df.loc[self.df["name"]==particle_name2].epsilon.iloc[0] 1924 sigma2 = self.df.loc[self.df["name"]==particle_name2].sigma.iloc[0] 1925 1926 return {"epsilon": np.sqrt(eps1*eps2), "sigma": (sigma1+sigma2)/2.0} 1927 1928 def get_particle_id_map(self, object_name): 1929 ''' 1930 Gets all the ids associated with the object with name `object_name` in `pmb.df` 1931 1932 Args: 1933 object_name(`str`): name of the object 1934 1935 Returns: 1936 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]}, } 1937 ''' 1938 object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0] 1939 valid_types = ["particle", "molecule", "residue", "protein"] 1940 if object_type not in valid_types: 1941 raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}") 1942 id_list = [] 1943 mol_map = {} 1944 res_map = {} 1945 def do_res_map(res_ids): 1946 for res_id in res_ids: 1947 res_list=self.df.loc[self.df['residue_id']== res_id].particle_id.dropna().tolist() 1948 res_map[res_id]=res_list 1949 return res_map 1950 if object_type in ['molecule', 'protein']: 1951 mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist() 1952 for mol_id in mol_ids: 1953 res_ids = set(self.df.loc[self.df['molecule_id']== mol_id].residue_id.dropna().tolist()) 1954 res_map=do_res_map(res_ids=res_ids) 1955 mol_list=self.df.loc[self.df['molecule_id']== mol_id].particle_id.dropna().tolist() 1956 id_list+=mol_list 1957 mol_map[mol_id]=mol_list 1958 elif object_type == 'residue': 1959 res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist() 1960 res_map=do_res_map(res_ids=res_ids) 1961 id_list=[] 1962 for res_id_list in res_map.values(): 1963 id_list+=res_id_list 1964 elif object_type == 'particle': 1965 id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist() 1966 return {"all": id_list, "molecule_map": mol_map, "residue_map": res_map} 1967 1968 def get_pka_set(self): 1969 ''' 1970 Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df` 1971 1972 Returns: 1973 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 1974 ''' 1975 titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna() 1976 pka_set = {} 1977 for index in titratables_AA_df.name.keys(): 1978 name = titratables_AA_df.name[index] 1979 pka_value = titratables_AA_df.pka[index] 1980 acidity = titratables_AA_df.acidity[index] 1981 pka_set[name] = {'pka_value':pka_value,'acidity':acidity} 1982 return pka_set 1983 1984 def get_radius_map(self): 1985 ''' 1986 Gets the effective radius of each `espresso type` in `pmb.df`. 1987 1988 Returns: 1989 radius_map(`dict`): {espresso_type: radius}. 1990 ''' 1991 1992 df_state_one = self.df[[('sigma',''),('state_one','es_type')]].dropna().drop_duplicates() 1993 df_state_two = self.df[[('sigma',''),('state_two','es_type')]].dropna().drop_duplicates() 1994 1995 state_one = pd.Series(df_state_one.sigma.values/2.0,index=df_state_one.state_one.es_type.values) 1996 state_two = pd.Series(df_state_two.sigma.values/2.0,index=df_state_two.state_two.es_type.values) 1997 radius_map = pd.concat([state_one,state_two],axis=0).to_dict() 1998 1999 return radius_map 2000 2001 def get_resource(self, path): 2002 ''' 2003 Locate a file resource of the pyMBE package. 2004 2005 Args: 2006 path(`str`): Relative path to the resource 2007 2008 Returns: 2009 path(`str`): Absolute path to the resource 2010 2011 ''' 2012 import os 2013 return os.path.join(os.path.dirname(__file__), path) 2014 2015 def get_type_map(self): 2016 """ 2017 Gets all different espresso types assigned to particles in `pmb.df`. 2018 2019 Returns: 2020 type_map(`dict`): {"name": espresso_type}. 2021 """ 2022 if self.df.state_one['es_type'].isnull().values.any(): 2023 df_state_one = self.df.state_one.dropna(how='all') 2024 df_state_two = self.df.state_two.dropna(how='all') 2025 else: 2026 df_state_one = self.df.state_one 2027 if self.df.state_two['es_type'].isnull().values.any(): 2028 df_state_two = self.df.state_two.dropna(how='all') 2029 else: 2030 df_state_two = self.df.state_two 2031 state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label) 2032 state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label) 2033 type_map = pd.concat([state_one,state_two],axis=0).to_dict() 2034 return type_map 2035 2036 def load_interaction_parameters(self, filename, verbose=True): 2037 """ 2038 Loads the interaction parameters stored in `filename` into `pmb.df` 2039 2040 Args: 2041 filename(`str`): name of the file to be read 2042 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2043 """ 2044 without_units = ['q','es_type','acidity'] 2045 with_units = ['sigma','epsilon'] 2046 2047 with open(filename, 'r') as f: 2048 interaction_data = json.load(f) 2049 interaction_parameter_set = interaction_data["data"] 2050 2051 for key in interaction_parameter_set: 2052 param_dict=interaction_parameter_set[key] 2053 object_type=param_dict['object_type'] 2054 if object_type == 'particle': 2055 not_requiered_attributes={} 2056 for not_requiered_key in without_units+with_units: 2057 if not_requiered_key in param_dict.keys(): 2058 if not_requiered_key in with_units: 2059 not_requiered_attributes[not_requiered_key]=self.create_variable_with_units(variable=param_dict.pop(not_requiered_key)) 2060 elif not_requiered_key in without_units: 2061 not_requiered_attributes[not_requiered_key]=param_dict.pop(not_requiered_key) 2062 else: 2063 if not_requiered_key == 'acidity': 2064 not_requiered_attributes[not_requiered_key] = 'inert' 2065 else: 2066 not_requiered_attributes[not_requiered_key]=None 2067 self.define_particle(name=param_dict.pop('name'), 2068 q=not_requiered_attributes.pop('q'), 2069 sigma=not_requiered_attributes.pop('sigma'), 2070 acidity=not_requiered_attributes.pop('acidity'), 2071 epsilon=not_requiered_attributes.pop('epsilon'), 2072 verbose=verbose) 2073 elif object_type == 'residue': 2074 self.define_residue (name = param_dict.pop('name'), 2075 central_bead = param_dict.pop('central_bead_name'), 2076 side_chains = param_dict.pop('side_chains_names')) 2077 elif object_type == 'molecule': 2078 self.define_molecule(name=param_dict.pop('name'), 2079 residue_list=param_dict.pop('residue_name_list')) 2080 elif object_type == 'peptide': 2081 self.define_peptide(name=param_dict.pop('name'), 2082 sequence=param_dict.pop('sequence'), 2083 model=param_dict.pop('model')) 2084 elif object_type == 'bond': 2085 particle_pairs = param_dict.pop('particle_pairs') 2086 bond_parameters = param_dict.pop('bond_parameters') 2087 bond_type = param_dict.pop('bond_type') 2088 if bond_type == 'harmonic': 2089 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2090 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2091 bond = {'r_0' : r_0, 2092 'k' : k, 2093 } 2094 2095 elif bond_type == 'FENE': 2096 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2097 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2098 d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max')) 2099 bond = {'r_0' : r_0, 2100 'k' : k, 2101 'd_r_max': d_r_max, 2102 } 2103 else: 2104 raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds") 2105 self.define_bond(bond_type=bond_type, bond_parameters=bond, particle_pairs=particle_pairs) 2106 else: 2107 raise ValueError(object_type+' is not a known pmb object type') 2108 2109 return 2110 2111 def load_pka_set(self, filename, verbose=True): 2112 """ 2113 Loads the pka_set stored in `filename` into `pmb.df`. 2114 2115 Args: 2116 filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}. 2117 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2118 2119 Note: 2120 - If `name` is already defined in the `pymbe.df`, it prints a warning. 2121 """ 2122 with open(filename, 'r') as f: 2123 pka_data = json.load(f) 2124 pka_set = pka_data["data"] 2125 2126 self.check_pka_set(pka_set=pka_set) 2127 2128 for key in pka_set: 2129 acidity = pka_set[key]['acidity'] 2130 pka_value = pka_set[key]['pka_value'] 2131 self.define_particle( 2132 name=key, 2133 acidity=acidity, 2134 pka=pka_value, 2135 verbose=verbose) 2136 return 2137 2138 def parse_sequence_from_file(self,sequence): 2139 """ 2140 Parses the given sequence such that it can be used in pyMBE. This function has to be used if the df was read from a file. 2141 2142 Args: 2143 sequence(`str`): sequence to be parsed 2144 2145 Returns: 2146 sequence(`lst`): parsed sequence 2147 """ 2148 sequence = sequence.replace(' ', '') 2149 sequence = sequence.replace("'", '') 2150 sequence = sequence.split(",")[1:-1] 2151 return sequence 2152 2153 def print_reduced_units(self): 2154 """ 2155 Prints the current set of reduced units defined in pyMBE.units. 2156 """ 2157 print("\nCurrent set of reduced units:") 2158 unit_length=self.units.Quantity(1,'reduced_length') 2159 unit_energy=self.units.Quantity(1,'reduced_energy') 2160 unit_charge=self.units.Quantity(1,'reduced_charge') 2161 print(unit_length.to('nm'), "=", unit_length) 2162 print(unit_energy.to('J'), "=", unit_energy) 2163 print('Temperature:', (self.kT/self.Kb).to("K")) 2164 print(unit_charge.to('C'), "=", unit_charge) 2165 print() 2166 2167 def propose_unused_type(self): 2168 """ 2169 Searches in `pmb.df` all the different particle types defined and returns a new one. 2170 2171 Returns: 2172 unused_type(`int`): unused particle type 2173 """ 2174 type_map=self.get_type_map() 2175 if type_map == {}: 2176 unused_type = 0 2177 else: 2178 unused_type=max(type_map.values())+1 2179 return unused_type 2180 2181 def protein_sequence_parser(self, sequence): 2182 ''' 2183 Parses `sequence` to the one letter code for amino acids. 2184 2185 Args: 2186 sequence(`str` or `lst`): Sequence of the amino acid. 2187 2188 Returns: 2189 clean_sequence(`lst`): `sequence` using the one letter code. 2190 2191 Note: 2192 - Accepted formats for `sequence` are: 2193 - `lst` with one letter or three letter code of each aminoacid in each element 2194 - `str` with the sequence using the one letter code 2195 - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-" 2196 2197 ''' 2198 # Aminoacid key 2199 keys={"ALA": "A", 2200 "ARG": "R", 2201 "ASN": "N", 2202 "ASP": "D", 2203 "CYS": "C", 2204 "GLU": "E", 2205 "GLN": "Q", 2206 "GLY": "G", 2207 "HIS": "H", 2208 "ILE": "I", 2209 "LEU": "L", 2210 "LYS": "K", 2211 "MET": "M", 2212 "PHE": "F", 2213 "PRO": "P", 2214 "SER": "S", 2215 "THR": "T", 2216 "TRP": "W", 2217 "TYR": "Y", 2218 "VAL": "V", 2219 "PSER": "J", 2220 "PTHR": "U", 2221 "PTyr": "Z", 2222 "NH2": "n", 2223 "COOH": "c"} 2224 clean_sequence=[] 2225 if isinstance(sequence, str): 2226 if sequence.find("-") != -1: 2227 splited_sequence=sequence.split("-") 2228 for residue in splited_sequence: 2229 if len(residue) == 1: 2230 if residue in keys.values(): 2231 residue_ok=residue 2232 else: 2233 if residue.upper() in keys.values(): 2234 residue_ok=residue.upper() 2235 else: 2236 raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence") 2237 clean_sequence.append(residue_ok) 2238 else: 2239 if residue in keys.keys(): 2240 clean_sequence.append(keys[residue]) 2241 else: 2242 if residue.upper() in keys.keys(): 2243 clean_sequence.append(keys[residue.upper()]) 2244 else: 2245 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2246 else: 2247 for residue in sequence: 2248 if residue in keys.values(): 2249 residue_ok=residue 2250 else: 2251 if residue.upper() in keys.values(): 2252 residue_ok=residue.upper() 2253 else: 2254 raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence") 2255 clean_sequence.append(residue_ok) 2256 if isinstance(sequence, list): 2257 for residue in sequence: 2258 if residue in keys.values(): 2259 residue_ok=residue 2260 else: 2261 if residue.upper() in keys.values(): 2262 residue_ok=residue.upper() 2263 elif (residue.upper() in keys.keys()): 2264 clean_sequence.append(keys[residue.upper()]) 2265 else: 2266 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2267 clean_sequence.append(residue_ok) 2268 return clean_sequence 2269 2270 def read_pmb_df (self,filename): 2271 """ 2272 Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE. 2273 2274 Args: 2275 filename(str): path to file. 2276 2277 Note: 2278 This function only accepts files with CSV format. 2279 """ 2280 2281 if filename.rsplit(".", 1)[1] != "csv": 2282 raise ValueError("Only files with CSV format are supported") 2283 df = pd.read_csv (filename,header=[0, 1], index_col=0) 2284 columns_names = self.setup_df() 2285 2286 multi_index = pd.MultiIndex.from_tuples(columns_names) 2287 df.columns = multi_index 2288 2289 self.df = self.convert_columns_to_original_format(df) 2290 2291 return self.df 2292 2293 def read_protein_vtf_in_df (self,filename,unit_length=None): 2294 """ 2295 Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`. 2296 2297 Args: 2298 filename(`str`): Path to the vtf file with the coarse-grained model. 2299 unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None. 2300 2301 Returns: 2302 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 2303 2304 Note: 2305 - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom. 2306 """ 2307 2308 print (f'Loading protein coarse grain model file: {filename}') 2309 2310 coord_list = [] 2311 particles_dict = {} 2312 2313 if unit_length is None: 2314 unit_length = 1 * self.units.angstrom 2315 2316 with open (filename,'r') as protein_model: 2317 for line in protein_model : 2318 line_split = line.split() 2319 if line_split: 2320 line_header = line_split[0] 2321 if line_header == 'atom': 2322 atom_id = line_split[1] 2323 atom_name = line_split[3] 2324 atom_resname = line_split[5] 2325 chain_id = line_split[9] 2326 radius = float(line_split [11])*unit_length 2327 sigma = 2*radius 2328 particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, sigma] 2329 elif line_header.isnumeric(): 2330 atom_coord = line_split[1:] 2331 atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord] 2332 coord_list.append (atom_coord) 2333 2334 numbered_label = [] 2335 i = 0 2336 2337 for atom_id in particles_dict.keys(): 2338 2339 if atom_id == 1: 2340 atom_name = particles_dict[atom_id][0] 2341 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2342 numbered_label.append(numbered_name) 2343 2344 elif atom_id != 1: 2345 2346 if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]: 2347 i += 1 2348 count = 1 2349 atom_name = particles_dict[atom_id][0] 2350 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2351 numbered_label.append(numbered_name) 2352 2353 elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]: 2354 if count == 2 or particles_dict[atom_id][1] == 'GLY': 2355 i +=1 2356 count = 0 2357 atom_name = particles_dict[atom_id][0] 2358 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2359 numbered_label.append(numbered_name) 2360 count +=1 2361 2362 topology_dict = {} 2363 2364 for i in range (0, len(numbered_label)): 2365 topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] , 2366 'chain_id':numbered_label[i][1], 2367 'sigma':numbered_label[i][2] } 2368 2369 return topology_dict 2370 2371 def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2372 """ 2373 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2374 If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead. 2375 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. 2376 2377 Args: 2378 particle_name1 (str): label of the type of the first particle type of the bonded particles. 2379 particle_name2 (str): label of the type of the second particle type of the bonded particles. 2380 hard_check (bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2381 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. 2382 2383 Returns: 2384 bond (cls): bond object from the espressomd library. 2385 2386 Note: 2387 - 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`. 2388 - If `hard_check`=`True` stops the code when no bond is found. 2389 """ 2390 2391 bond_key = self.find_bond_key(particle_name1=particle_name1, 2392 particle_name2=particle_name2, 2393 use_default_bond=use_default_bond) 2394 if use_default_bond: 2395 if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'): 2396 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") 2397 if bond_key: 2398 return self.df[self.df['name']==bond_key].bond_object.values[0] 2399 else: 2400 print("Bond not defined between particles ", particle_name1, " and ", particle_name2) 2401 if hard_check: 2402 sys.exit(1) 2403 else: 2404 return 2405 2406 def search_particles_in_residue(self, residue_name): 2407 ''' 2408 Searches for all particles in a given residue of name `residue_name`. 2409 2410 Args: 2411 residue_name (`str`): name of the residue to be searched 2412 2413 Returns: 2414 list_of_particles_in_residue (`lst`): list of the names of all particles in the residue 2415 2416 Note: 2417 - 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. 2418 2419 ''' 2420 index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 2421 central_bead = self.df.at [index_residue, ('central_bead', '')] 2422 list_of_side_chains = self.df.at [index_residue, ('side_chains', '')] 2423 2424 list_of_particles_in_residue = [] 2425 list_of_particles_in_residue.append(central_bead) 2426 2427 for side_chain in list_of_side_chains: 2428 object_type = self.df[self.df['name']==side_chain].pmb_type.values[0] 2429 2430 if object_type == "residue": 2431 list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain) 2432 list_of_particles_in_residue += list_of_particles_in_side_chain_residue 2433 elif object_type == "particle": 2434 list_of_particles_in_residue.append(side_chain) 2435 2436 return list_of_particles_in_residue 2437 2438 def set_particle_acidity(self, name, acidity='inert', default_charge=0, pka=None, verbose=True): 2439 """ 2440 Sets the particle acidity if it is acidic or basic, creates `state_one` and `state_two` with the protonated and 2441 deprotonated states. In each state is set: `label`,`charge` and `es_type`. If it is inert, it will define only `state_one`. 2442 2443 Args: 2444 name (`str`): Unique label that identifies the `particle`. 2445 acidity (`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to `inert`. 2446 default_charge (`int`): Charge of the particle. Defaults to 0. 2447 pka (`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to None. 2448 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2449 """ 2450 acidity_valid_keys = ['inert','acidic', 'basic'] 2451 if acidity not in acidity_valid_keys: 2452 raise ValueError(f"Acidity {acidity} provided for particle name {name} is not supproted. Valid keys are: {acidity_valid_keys}") 2453 if acidity in ['acidic', 'basic'] and pka is None: 2454 raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.") 2455 for index in self.df[self.df['name']==name].index: 2456 if pka: 2457 self.add_value_to_df(key=('pka',''), 2458 index=index, 2459 new_value=pka, 2460 verbose=verbose) 2461 self.add_value_to_df(key=('acidity',''), 2462 index=index, 2463 new_value=acidity, 2464 verbose=verbose) 2465 if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')): 2466 self.add_value_to_df(key=('state_one','es_type'), 2467 index=index, 2468 new_value=self.propose_unused_type(), 2469 verbose=verbose) 2470 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'inert': 2471 self.add_value_to_df(key=('state_one','charge'), 2472 index=index, 2473 new_value=default_charge, 2474 verbose=verbose) 2475 self.add_value_to_df(key=('state_one','label'), 2476 index=index, 2477 new_value=name, 2478 verbose=verbose) 2479 else: 2480 protonated_label = f'{name}H' 2481 self.add_value_to_df(key=('state_one','label'), 2482 index=index, 2483 new_value=protonated_label, 2484 verbose=verbose) 2485 self.add_value_to_df(key=('state_two','label'), 2486 index=index, 2487 new_value=name, 2488 verbose=verbose) 2489 if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')): 2490 self.add_value_to_df(key=('state_two','es_type'), 2491 index=index, 2492 new_value=self.propose_unused_type(), 2493 verbose=verbose) 2494 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'acidic': 2495 self.add_value_to_df(key=('state_one','charge'), 2496 index=index,new_value=0, 2497 verbose=verbose) 2498 self.add_value_to_df(key=('state_two','charge'), 2499 index=index, 2500 new_value=-1, 2501 verbose=verbose) 2502 elif self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'basic': 2503 self.add_value_to_df(key=('state_one','charge'), 2504 index=index,new_value=+1, 2505 verbose=verbose) 2506 self.add_value_to_df(key=('state_two','charge'), 2507 index=index, 2508 new_value=0, 2509 verbose=verbose) 2510 return 2511 2512 def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True): 2513 """ 2514 Sets the set of reduced units used by pyMBE.units and it prints it. 2515 2516 Args: 2517 unit_length (`obj`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 2518 unit_charge (`obj`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 2519 temperature (`obj`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 2520 Kw (`obj`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 2521 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2522 2523 Note: 2524 - If no `temperature` is given, a value of 298.15 K is assumed by default. 2525 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 2526 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 2527 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 2528 """ 2529 self.units=pint.UnitRegistry() 2530 if unit_length is None: 2531 unit_length=0.355*self.units.nm 2532 if temperature is None: 2533 temperature=298.15 * self.units.K 2534 if unit_charge is None: 2535 unit_charge=self.units.e 2536 if Kw is None: 2537 Kw = 1e-14 2538 self.N_A=6.02214076e23 / self.units.mol 2539 self.Kb=1.38064852e-23 * self.units.J / self.units.K 2540 self.e=1.60217662e-19 *self.units.C 2541 self.kT=temperature*self.Kb 2542 self.Kw=Kw*self.units.mol**2 / (self.units.l**2) 2543 self.units.define(f'reduced_energy = {self.kT} ') 2544 self.units.define(f'reduced_length = {unit_length}') 2545 self.units.define(f'reduced_charge = {unit_charge}') 2546 if verbose: 2547 self.print_reduced_units() 2548 return 2549 2550 def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2551 """ 2552 Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 2553 2554 Args: 2555 counter_ion(`str`): `name` of the counter_ion `particle`. 2556 constant_pH(`float`): pH-value. 2557 exclusion_range(`float`, optional): Below this value, no particles will be inserted. 2558 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2559 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2560 2561 Returns: 2562 RE (`obj`): Instance of a reaction_ensemble.ConstantpHEnsemble object from the espressomd library. 2563 sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2564 """ 2565 from espressomd import reaction_methods 2566 if exclusion_range is None: 2567 exclusion_range = max(self.get_radius_map().values())*2.0 2568 if pka_set is None: 2569 pka_set=self.get_pka_set() 2570 self.check_pka_set(pka_set=pka_set) 2571 if use_exclusion_radius_per_type: 2572 exclusion_radius_per_type = self.get_radius_map() 2573 else: 2574 exclusion_radius_per_type = {} 2575 2576 RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2577 exclusion_range=exclusion_range.magnitude, 2578 seed=self.SEED, 2579 constant_pH=constant_pH, 2580 exclusion_radius_per_type = exclusion_radius_per_type 2581 ) 2582 sucessfull_reactions_labels=[] 2583 charge_map = self.get_charge_map() 2584 for name in pka_set.keys(): 2585 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2586 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in sg.type_map') 2587 continue 2588 gamma=10**-pka_set[name]['pka_value'] 2589 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2590 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2591 counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0] 2592 RE.add_reaction(gamma=gamma, 2593 reactant_types=[state_one_type], 2594 product_types=[state_two_type, counter_ion_type], 2595 default_charges={state_one_type: charge_map[state_one_type], 2596 state_two_type: charge_map[state_two_type], 2597 counter_ion_type: charge_map[counter_ion_type]}) 2598 sucessfull_reactions_labels.append(name) 2599 return RE, sucessfull_reactions_labels 2600 2601 def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient=None, exclusion_range=None, use_exclusion_radius_per_type = False): 2602 """ 2603 Sets up grand-canonical coupling to a reservoir of salt. 2604 For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead. 2605 2606 Args: 2607 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2608 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 2609 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 2610 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2611 exclusion_range(`float`, optional): For distances shorter than this value, no particles will be inserted. 2612 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2613 2614 Returns: 2615 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2616 """ 2617 from espressomd import reaction_methods 2618 if exclusion_range is None: 2619 exclusion_range = max(self.get_radius_map().values())*2.0 2620 if use_exclusion_radius_per_type: 2621 exclusion_radius_per_type = self.get_radius_map() 2622 else: 2623 exclusion_radius_per_type = {} 2624 2625 #If no function for the activity coefficient is provided, the ideal case is assumed. 2626 if activity_coefficient is None: 2627 activity_coefficient = lambda x: 1.0 2628 2629 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2630 exclusion_range=exclusion_range.magnitude, 2631 seed=self.SEED, 2632 exclusion_radius_per_type = exclusion_radius_per_type 2633 ) 2634 2635 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2636 determined_activity_coefficient = activity_coefficient(c_salt_res) 2637 K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient 2638 2639 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 2640 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 2641 2642 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.charge.values[0] 2643 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.charge.values[0] 2644 2645 if salt_cation_charge <= 0: 2646 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 2647 if salt_anion_charge >= 0: 2648 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 2649 2650 # Grand-canonical coupling to the reservoir 2651 RE.add_reaction( 2652 gamma = K_salt.magnitude, 2653 reactant_types = [], 2654 reactant_coefficients = [], 2655 product_types = [ salt_cation_es_type, salt_anion_es_type ], 2656 product_coefficients = [ 1, 1 ], 2657 default_charges = { 2658 salt_cation_es_type: salt_cation_charge, 2659 salt_anion_es_type: salt_anion_charge, 2660 } 2661 ) 2662 2663 return RE 2664 2665 def setup_grxmc_reactions(self, pH_res, c_salt_res, proton_name, hydroxide_name, salt_cation_name, salt_anion_name, activity_coefficient=None, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2666 """ 2667 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 2668 reservoir of small ions. 2669 This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1]. 2670 2671 [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. 2672 2673 Args: 2674 pH_res ('float): pH-value in the reservoir. 2675 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2676 proton_name ('str'): Name of the proton (H+) particle. 2677 hydroxide_name ('str'): Name of the hydroxide (OH-) particle. 2678 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 2679 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 2680 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2681 exclusion_range(`float`, optional): For distances shorter than this value, no particles will be inserted. 2682 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2683 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2684 2685 Returns: 2686 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2687 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2688 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 2689 """ 2690 from espressomd import reaction_methods 2691 if exclusion_range is None: 2692 exclusion_range = max(self.get_radius_map().values())*2.0 2693 if pka_set is None: 2694 pka_set=self.get_pka_set() 2695 self.check_pka_set(pka_set=pka_set) 2696 if use_exclusion_radius_per_type: 2697 exclusion_radius_per_type = self.get_radius_map() 2698 else: 2699 exclusion_radius_per_type = {} 2700 2701 #If no function for the activity coefficient is provided, the ideal case is assumed. 2702 if activity_coefficient is None: 2703 activity_coefficient = lambda x: 1.0 2704 2705 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2706 exclusion_range=exclusion_range.magnitude, 2707 seed=self.SEED, 2708 exclusion_radius_per_type = exclusion_radius_per_type 2709 ) 2710 2711 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2712 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 2713 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2714 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 2715 K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2716 K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2717 K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2718 2719 proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0] 2720 hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0] 2721 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 2722 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 2723 2724 proton_charge = self.df.loc[self.df['name']==proton_name].state_one.charge.values[0] 2725 hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.charge.values[0] 2726 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.charge.values[0] 2727 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.charge.values[0] 2728 2729 if proton_charge <= 0: 2730 raise ValueError('ERROR proton charge must be positive, charge ', proton_charge) 2731 if salt_cation_charge <= 0: 2732 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 2733 if hydroxide_charge >= 0: 2734 raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge) 2735 if salt_anion_charge >= 0: 2736 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 2737 2738 # Grand-canonical coupling to the reservoir 2739 # 0 = H+ + OH- 2740 RE.add_reaction( 2741 gamma = K_W.magnitude, 2742 reactant_types = [], 2743 reactant_coefficients = [], 2744 product_types = [ proton_es_type, hydroxide_es_type ], 2745 product_coefficients = [ 1, 1 ], 2746 default_charges = { 2747 proton_es_type: proton_charge, 2748 hydroxide_es_type: hydroxide_charge, 2749 } 2750 ) 2751 2752 # 0 = Na+ + Cl- 2753 RE.add_reaction( 2754 gamma = K_NACL.magnitude, 2755 reactant_types = [], 2756 reactant_coefficients = [], 2757 product_types = [ salt_cation_es_type, salt_anion_es_type ], 2758 product_coefficients = [ 1, 1 ], 2759 default_charges = { 2760 salt_cation_es_type: salt_cation_charge, 2761 salt_anion_es_type: salt_anion_charge, 2762 } 2763 ) 2764 2765 # 0 = Na+ + OH- 2766 RE.add_reaction( 2767 gamma = (K_NACL * K_W / K_HCL).magnitude, 2768 reactant_types = [], 2769 reactant_coefficients = [], 2770 product_types = [ salt_cation_es_type, hydroxide_es_type ], 2771 product_coefficients = [ 1, 1 ], 2772 default_charges = { 2773 salt_cation_es_type: salt_cation_charge, 2774 hydroxide_es_type: hydroxide_charge, 2775 } 2776 ) 2777 2778 # 0 = H+ + Cl- 2779 RE.add_reaction( 2780 gamma = K_HCL.magnitude, 2781 reactant_types = [], 2782 reactant_coefficients = [], 2783 product_types = [ proton_es_type, salt_anion_es_type ], 2784 product_coefficients = [ 1, 1 ], 2785 default_charges = { 2786 proton_es_type: proton_charge, 2787 salt_anion_es_type: salt_anion_charge, 2788 } 2789 ) 2790 2791 # Annealing moves to ensure sufficient sampling 2792 # Cation annealing H+ = Na+ 2793 RE.add_reaction( 2794 gamma = (K_NACL / K_HCL).magnitude, 2795 reactant_types = [proton_es_type], 2796 reactant_coefficients = [ 1 ], 2797 product_types = [ salt_cation_es_type ], 2798 product_coefficients = [ 1 ], 2799 default_charges = { 2800 proton_es_type: proton_charge, 2801 salt_cation_es_type: salt_cation_charge, 2802 } 2803 ) 2804 2805 # Anion annealing OH- = Cl- 2806 RE.add_reaction( 2807 gamma = (K_HCL / K_W).magnitude, 2808 reactant_types = [hydroxide_es_type], 2809 reactant_coefficients = [ 1 ], 2810 product_types = [ salt_anion_es_type ], 2811 product_coefficients = [ 1 ], 2812 default_charges = { 2813 hydroxide_es_type: hydroxide_charge, 2814 salt_anion_es_type: salt_anion_charge, 2815 } 2816 ) 2817 2818 sucessful_reactions_labels=[] 2819 charge_map = self.get_charge_map() 2820 for name in pka_set.keys(): 2821 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2822 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.') 2823 continue 2824 2825 Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 2826 2827 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2828 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2829 2830 # Reaction in terms of proton: HA = A + H+ 2831 RE.add_reaction(gamma=Ka.magnitude, 2832 reactant_types=[state_one_type], 2833 reactant_coefficients=[1], 2834 product_types=[state_two_type, proton_es_type], 2835 product_coefficients=[1, 1], 2836 default_charges={state_one_type: charge_map[state_one_type], 2837 state_two_type: charge_map[state_two_type], 2838 proton_es_type: proton_charge}) 2839 2840 # Reaction in terms of salt cation: HA = A + Na+ 2841 RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude, 2842 reactant_types=[state_one_type], 2843 reactant_coefficients=[1], 2844 product_types=[state_two_type, salt_cation_es_type], 2845 product_coefficients=[1, 1], 2846 default_charges={state_one_type: charge_map[state_one_type], 2847 state_two_type: charge_map[state_two_type], 2848 salt_cation_es_type: salt_cation_charge}) 2849 2850 # Reaction in terms of hydroxide: OH- + HA = A 2851 RE.add_reaction(gamma=(Ka / K_W).magnitude, 2852 reactant_types=[state_one_type, hydroxide_es_type], 2853 reactant_coefficients=[1, 1], 2854 product_types=[state_two_type], 2855 product_coefficients=[1], 2856 default_charges={state_one_type: charge_map[state_one_type], 2857 state_two_type: charge_map[state_two_type], 2858 hydroxide_es_type: hydroxide_charge}) 2859 2860 # Reaction in terms of salt anion: Cl- + HA = A 2861 RE.add_reaction(gamma=(Ka / K_HCL).magnitude, 2862 reactant_types=[state_one_type, salt_anion_es_type], 2863 reactant_coefficients=[1, 1], 2864 product_types=[state_two_type], 2865 product_coefficients=[1], 2866 default_charges={state_one_type: charge_map[state_one_type], 2867 state_two_type: charge_map[state_two_type], 2868 salt_anion_es_type: salt_anion_charge}) 2869 2870 sucessful_reactions_labels.append(name) 2871 return RE, sucessful_reactions_labels, ionic_strength_res 2872 2873 def setup_grxmc_unified (self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient=None, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2874 """ 2875 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 2876 reservoir of small ions. 2877 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-}. 2878 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'. 2879 2880 [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4). 2881 [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. 2882 2883 Args: 2884 pH_res ('float'): pH-value in the reservoir. 2885 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2886 cation_name ('str'): Name of the cationic particle. 2887 anion_name ('str'): Name of the anionic particle. 2888 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2889 exclusion_range(`float`, optional): Below this value, no particles will be inserted. 2890 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2891 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`. 2892 2893 Returns: 2894 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2895 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2896 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 2897 """ 2898 from espressomd import reaction_methods 2899 if exclusion_range is None: 2900 exclusion_range = max(self.get_radius_map().values())*2.0 2901 if pka_set is None: 2902 pka_set=self.get_pka_set() 2903 self.check_pka_set(pka_set=pka_set) 2904 if use_exclusion_radius_per_type: 2905 exclusion_radius_per_type = self.get_radius_map() 2906 else: 2907 exclusion_radius_per_type = {} 2908 2909 #If no function for the activity coefficient is provided, the ideal case is assumed. 2910 if activity_coefficient is None: 2911 activity_coefficient = lambda x: 1.0 2912 2913 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2914 exclusion_range=exclusion_range.magnitude, 2915 seed=self.SEED, 2916 exclusion_radius_per_type = exclusion_radius_per_type 2917 ) 2918 2919 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2920 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 2921 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2922 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 2923 a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 2924 a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 2925 a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 2926 K_XX = a_cation * a_anion 2927 2928 cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0] 2929 anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0] 2930 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.values[0] 2931 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.values[0] 2932 if cation_charge <= 0: 2933 raise ValueError('ERROR cation charge must be positive, charge ', cation_charge) 2934 if anion_charge >= 0: 2935 raise ValueError('ERROR anion charge must be negative, charge ', anion_charge) 2936 2937 # Coupling to the reservoir: 0 = X+ + X- 2938 RE.add_reaction( 2939 gamma = K_XX.magnitude, 2940 reactant_types = [], 2941 reactant_coefficients = [], 2942 product_types = [ cation_es_type, anion_es_type ], 2943 product_coefficients = [ 1, 1 ], 2944 default_charges = { 2945 cation_es_type: cation_charge, 2946 anion_es_type: anion_charge, 2947 } 2948 ) 2949 2950 sucessful_reactions_labels=[] 2951 charge_map = self.get_charge_map() 2952 for name in pka_set.keys(): 2953 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2954 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in sg.type_map') 2955 continue 2956 2957 Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l 2958 gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen 2959 2960 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2961 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2962 2963 # Reaction in terms of small cation: HA = A + X+ 2964 RE.add_reaction(gamma=gamma_K_AX.magnitude, 2965 reactant_types=[state_one_type], 2966 reactant_coefficients=[1], 2967 product_types=[state_two_type, cation_es_type], 2968 product_coefficients=[1, 1], 2969 default_charges={state_one_type: charge_map[state_one_type], 2970 state_two_type: charge_map[state_two_type], 2971 cation_es_type: cation_charge}) 2972 2973 # Reaction in terms of small anion: X- + HA = A 2974 RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude, 2975 reactant_types=[state_one_type, anion_es_type], 2976 reactant_coefficients=[1, 1], 2977 product_types=[state_two_type], 2978 product_coefficients=[1], 2979 default_charges={state_one_type: charge_map[state_one_type], 2980 state_two_type: charge_map[state_two_type], 2981 anion_es_type: anion_charge}) 2982 2983 sucessful_reactions_labels.append(name) 2984 return RE, sucessful_reactions_labels, ionic_strength_res 2985 2986 def setup_df (self): 2987 """ 2988 Sets up the pyMBE's dataframe `pymbe.df`. 2989 2990 Returns: 2991 columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe 2992 """ 2993 2994 columns_dtypes = { 2995 'name': { 2996 '': str}, 2997 'pmb_type': { 2998 '': str}, 2999 'particle_id': { 3000 '': object}, 3001 'particle_id2': { 3002 '': object}, 3003 'residue_id': { 3004 '': object}, 3005 'molecule_id': { 3006 '': object}, 3007 'acidity': { 3008 '': str}, 3009 'pka': { 3010 '': float}, 3011 'central_bead': { 3012 '': object}, 3013 'side_chains': { 3014 '': object}, 3015 'residue_list': { 3016 '': object}, 3017 'model': { 3018 '': str}, 3019 'sigma': { 3020 '': object}, 3021 'cutoff': { 3022 '': object}, 3023 'offset': { 3024 '': object}, 3025 'epsilon': { 3026 '': object}, 3027 'state_one': { 3028 'label': str, 3029 'es_type': object, 3030 'charge': object }, 3031 'state_two': { 3032 'label': str, 3033 'es_type': object, 3034 'charge': object }, 3035 'sequence': { 3036 '': object}, 3037 'bond_object': { 3038 '': object}, 3039 'parameters_of_the_potential':{ 3040 '': object}, 3041 'l0': { 3042 '': float}, 3043 } 3044 3045 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()])) 3046 3047 for level1, sub_dtypes in columns_dtypes.items(): 3048 for level2, dtype in sub_dtypes.items(): 3049 self.df[level1, level2] = self.df[level1, level2].astype(dtype) 3050 3051 columns_names = pd.MultiIndex.from_frame(self.df) 3052 columns_names = columns_names.names 3053 3054 return columns_names 3055 3056 def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True): 3057 """ 3058 Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`. 3059 3060 Args: 3061 espresso_system(`obj`): Instance of a system object from the espressomd library. 3062 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. 3063 combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'. 3064 warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True. 3065 3066 Note: 3067 - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 3068 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 3069 - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html 3070 3071 """ 3072 from itertools import combinations_with_replacement 3073 implemented_combining_rules = ['Lorentz-Berthelot'] 3074 compulsory_parameters_in_df = ['sigma','epsilon'] 3075 # Sanity check 3076 if combining_rule not in implemented_combining_rules: 3077 raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules) 3078 if shift_potential: 3079 shift="auto" 3080 else: 3081 shift=0 3082 # List which particles have sigma and epsilon values defined in pmb.df and which ones don't 3083 particles_types_with_LJ_parameters = [] 3084 non_parametrized_labels= [] 3085 for particle_type in self.get_type_map().values(): 3086 check_list=[] 3087 for key in compulsory_parameters_in_df: 3088 value_in_df=self.find_value_from_es_type(es_type=particle_type, 3089 column_name=key) 3090 check_list.append(np.isnan(value_in_df)) 3091 if any(check_list): 3092 non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 3093 column_name='label')) 3094 else: 3095 particles_types_with_LJ_parameters.append(particle_type) 3096 lj_parameters_keys=["label","sigma","epsilon","offset","cutoff"] 3097 # Set up LJ interactions between all particle types 3098 for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 3099 lj_parameters={} 3100 for key in lj_parameters_keys: 3101 lj_parameters[key]=[] 3102 # Search the LJ parameters of the type pair 3103 for ptype in type_pair: 3104 for key in lj_parameters_keys: 3105 lj_parameters[key].append(self.find_value_from_es_type(es_type=ptype, column_name=key)) 3106 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 3107 if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]): 3108 continue 3109 # Apply combining rule 3110 if combining_rule == 'Lorentz-Berthelot': 3111 sigma=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2 3112 cutoff=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2 3113 offset=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2 3114 epsilon=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1]) 3115 espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = epsilon.to('reduced_energy').magnitude, 3116 sigma = sigma.to('reduced_length').magnitude, 3117 cutoff = cutoff.to('reduced_length').magnitude, 3118 offset = offset.to("reduced_length").magnitude, 3119 shift = shift) 3120 index = len(self.df) 3121 self.df.at [index, 'name'] = f'LJ: {lj_parameters["label"][0]}-{lj_parameters["label"][1]}' 3122 lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params() 3123 3124 self.add_value_to_df(index=index, 3125 key=('pmb_type',''), 3126 new_value='LennardJones') 3127 3128 self.add_value_to_df(index=index, 3129 key=('parameters_of_the_potential',''), 3130 new_value=lj_params, 3131 non_standard_value=True) 3132 if non_parametrized_labels and warnings: 3133 print(f'WARNING: The following particles do not have a defined value of sigma or epsilon in pmb.df: {non_parametrized_labels}. No LJ interaction has been added in ESPResSo for those particles.') 3134 3135 return 3136 3137 def setup_particle_sigma(self, topology_dict): 3138 ''' 3139 Sets up sigma of the particles in `topology_dict`. 3140 3141 Args: 3142 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 3143 ''' 3144 for residue in topology_dict.keys(): 3145 residue_name = re.split(r'\d+', residue)[0] 3146 residue_number = re.split(r'(\d+)', residue)[1] 3147 residue_sigma = topology_dict[residue]['sigma'] 3148 sigma = residue_sigma*self.units.nm 3149 index = self.df[(self.df['residue_id']==residue_number) & (self.df['name']==residue_name) ].index.values[0] 3150 self.add_value_to_df(key= ('sigma',''), 3151 index=int (index), 3152 new_value=sigma) 3153 return 3154 3155 def write_pmb_df (self, filename): 3156 ''' 3157 Writes the pyMBE dataframe into a csv file 3158 3159 Args: 3160 filename (`str`): Path to the csv file 3161 ''' 3162 3163 self.df.to_csv(filename) 3164 3165 return 3166 3167 def write_output_vtf_file(self, espresso_system, filename): 3168 ''' 3169 Writes a snapshot of `espresso_system` on the vtf file `filename`. 3170 3171 Args: 3172 espresso_system(`obj`): Instance of a system object from the espressomd library. 3173 filename(`str`): Path to the vtf file. 3174 3175 ''' 3176 box = espresso_system.box_l[0] 3177 with open(filename, mode='w+t') as coordinates: 3178 coordinates.write (f'unitcell {box} {box} {box} \n') 3179 for particle in espresso_system.part: 3180 type_label = self.find_value_from_es_type(es_type=particle.type, column_name='label') 3181 coordinates.write (f'atom {particle.id} radius 1 name {type_label} type {type_label}\n' ) 3182 coordinates.write ('timestep indexed\n') 3183 for particle in espresso_system.part: 3184 coordinates.write (f'{particle.id} \t {particle.pos[0]} \t {particle.pos[1]} \t {particle.pos[2]}\n') 3185 return
The library for the Molecular Builder for ESPResSo (pyMBE)
Attributes:
- N_A(
obj
): Avogadro number using thepmb.units
UnitRegistry. - Kb(
obj
): Boltzmann constant using thepmb.units
UnitRegistry. - e(
obj
): Elemental charge using thepmb.units
UnitRegistry. - df(
obj
): PandasDataframe used to bookkeep all the information stored in pyMBE. Typically refered aspmb.df
. - kT(
obj
): Thermal energy using thepmb.units
UnitRegistry. It is used as the unit of reduced energy. - Kw(
obj
): Ionic product of water using thepmb.units
UnitRegistry. Used in the setup of the G-RxMC method.
51 def __init__(self, SEED, temperature=None, unit_length=None, unit_charge=None, Kw=None): 52 """ 53 Initializes the pymbe_library by setting up the reduced unit system with `temperature` and `reduced_length` 54 and sets up the `pmb.df` for bookkeeping. 55 56 Args: 57 temperature(`obj`,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. 58 unit_length(`obj`, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. 59 unit_charge (`obj`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 60 Kw (`obj`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 61 62 Note: 63 - If no `temperature` is given, a value of 298.15 K is assumed by default. 64 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 65 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 66 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 67 """ 68 # Seed and RNG 69 self.SEED=SEED 70 self.rng = np.random.default_rng(SEED) 71 self.set_reduced_units(unit_length=unit_length, unit_charge=unit_charge, 72 temperature=temperature, Kw=Kw, verbose=False) 73 self.setup_df() 74 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(
obj
,optional): Value of the temperature in the pyMBE UnitRegistry. Defaults to None. - unit_length(
obj
, optional): Value of the unit of length in the pyMBE UnitRegistry. Defaults to None. - unit_charge (
obj
,optional): Reduced unit of charge defined using thepmb.units
UnitRegistry. Defaults to None. - Kw (
obj
,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.
76 def enable_motion_of_rigid_object(self, name, espresso_system): 77 ''' 78 Enables the motion of the rigid object `name` in the `espresso_system`. 79 80 Args: 81 name(`str`): Label of the object. 82 espresso_system(`obj`): Instance of a system object from the espressomd library. 83 84 Note: 85 - It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]. 86 ''' 87 print ('enable_motion_of_rigid_object requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"]') 88 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 89 if pmb_type != 'protein': 90 raise ValueError (f'The pmb_type: {pmb_type} is not currently supported. The supported pmb_type is: protein') 91 molecule_ids_list = self.df.loc[self.df['name']==name].molecule_id.to_list() 92 for molecule_id in molecule_ids_list: 93 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 94 center_of_mass = self.calculate_center_of_mass_of_molecule ( molecule_id=molecule_id,espresso_system=espresso_system) 95 rigid_object_center = espresso_system.part.add(pos=center_of_mass, 96 rotation=[True,True,True], 97 type=self.propose_unused_type()) 98 for particle_id in particle_ids_list: 99 pid = espresso_system.part.by_id(particle_id) 100 pid.vs_auto_relate_to(rigid_object_center.id) 101 return
Enables the motion of the rigid object name
in the espresso_system
.
Arguments:
- name(
str
): Label of the object. - espresso_system(
obj
): Instance of a system object from the espressomd library.
Note:
- It requires that espressomd has the following features activated: ["VIRTUAL_SITES_RELATIVE", "MASS"].
103 def add_bond_in_df(self, particle_id1, particle_id2, use_default_bond=False): 104 """ 105 Adds a bond entry on the `pymbe.df` storing the particle_ids of the two bonded particles. 106 107 Args: 108 particle_id1 (`int`): particle_id of the type of the first particle type of the bonded particles 109 particle_id2 (`int`): particle_id of the type of the second particle type of the bonded particles 110 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. 111 """ 112 particle_name1 = self.df.loc[self.df['particle_id']==particle_id1].name.values[0] 113 particle_name2 = self.df.loc[self.df['particle_id']==particle_id2].name.values[0] 114 bond_key = self.find_bond_key(particle_name1=particle_name1, 115 particle_name2=particle_name2, 116 use_default_bond=use_default_bond) 117 if not bond_key: 118 return 119 self.copy_df_entry(name=bond_key,column_name='particle_id2',number_of_copies=1) 120 indexs = np.where(self.df['name']==bond_key) 121 index_list = list (indexs[0]) 122 used_bond_df = self.df.loc[self.df['particle_id2'].notnull()] 123 #without this drop the program crashes when dropping duplicates because the 'bond' column is a dict 124 used_bond_df = used_bond_df.drop([('bond_object','')],axis =1 ) 125 used_bond_index = used_bond_df.index.to_list() 126 127 for index in index_list: 128 if index not in used_bond_index: 129 self.clean_df_row(index=int(index)) 130 self.df.at[index,'particle_id'] = particle_id1 131 self.df.at[index,'particle_id2'] = particle_id2 132 break 133 return
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.
135 def add_bonds_to_espresso(self, espresso_system) : 136 """ 137 Adds all bonds defined in `pmb.df` to `espresso_system`. 138 139 Args: 140 espresso_system (str): system object of espressomd library 141 """ 142 143 if 'bond' in self.df.values: 144 bond_df = self.df.loc[self.df ['pmb_type'] == 'bond'] 145 bond_list = bond_df.bond_object.values.tolist() 146 for bond in bond_list: 147 espresso_system.bonded_inter.add(bond) 148 else: 149 print ('WARNING: There are no bonds defined in pymbe.df') 150 151 return
Adds all bonds defined in pmb.df
to espresso_system
.
Arguments:
- espresso_system (str): system object of espressomd library
153 def add_value_to_df(self,index,key,new_value, verbose=True, non_standard_value=False): 154 """ 155 Adds a value to a cell in the `pmb.df` DataFrame. 156 157 Args: 158 index(`int`): index of the row to add the value to. 159 key(`str`): the column label to add the value to. 160 verbose(`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 161 non_standard_value(`bool`, optional): Switch to enable insertion of non-standard values, such as `dict` objects. Defaults to False. 162 """ 163 164 token = "#protected:" 165 166 def protect(obj): 167 if non_standard_value: 168 return token + json.dumps(obj, cls=self.NumpyEncoder) 169 return obj 170 171 def deprotect(obj): 172 if non_standard_value and isinstance(obj, str) and obj.startswith(token): 173 return json.loads(obj.removeprefix(token)) 174 return obj 175 176 # Make sure index is a scalar integer value 177 index = int (index) 178 assert isinstance(index, int), '`index` should be a scalar integer value.' 179 idx = pd.IndexSlice 180 if self.check_if_df_cell_has_a_value(index=index,key=key): 181 old_value= self.df.loc[index,idx[key]] 182 if verbose: 183 if protect(old_value) != protect(new_value): 184 name=self.df.loc[index,('name','')] 185 pmb_type=self.df.loc[index,('pmb_type','')] 186 print(f"WARNING: you are redefining the properties of {name} of pmb_type {pmb_type}") 187 print(f'WARNING: overwritting the value of the entry `{key}`: old_value = {old_value} new_value = {new_value}') 188 self.df.loc[index,idx[key]] = protect(new_value) 189 if non_standard_value: 190 self.df[key] = self.df[key].apply(deprotect) 191 return
Adds a value to a cell in the pmb.df
DataFrame.
Arguments:
- index(
int
): index of the row to add the value to. - key(
str
): the column label to add the value to. - verbose(
bool
, optional): Switch to activate/deactivate verbose. Defaults to True. - non_standard_value(
bool
, optional): Switch to enable insertion of non-standard values, such asdict
objects. Defaults to False.
193 def assign_molecule_id(self, name, molecule_index, pmb_type, used_molecules_id): 194 """ 195 Assigns the `molecule_id` of the pmb object given by `pmb_type` 196 197 Args: 198 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 199 pmb_type(`str`): pmb_object_type to assign the `molecule_id` 200 molecule_index (`int`): index of the current `pmb_object_type` to assign the `molecule_id` 201 used_molecules_id (`lst`): list with the `molecule_id` values already used. 202 203 Returns: 204 molecule_id(`int`): Id of the molecule 205 """ 206 207 self.clean_df_row(index=int(molecule_index)) 208 209 if self.df['molecule_id'].isnull().values.all(): 210 molecule_id = 0 211 else: 212 # check if a residue is part of another molecule 213 check_residue_name = self.df[self.df['residue_list'].astype(str).str.contains(name)] 214 mol_pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 215 if not check_residue_name.empty and mol_pmb_type == pmb_type: 216 for value in check_residue_name.index.to_list(): 217 if value not in used_molecules_id: 218 molecule_id = self.df.loc[value].molecule_id.values[0] 219 break 220 else: 221 molecule_id = self.df['molecule_id'].max() +1 222 223 self.add_value_to_df (key=('molecule_id',''), 224 index=int(molecule_index), 225 new_value=molecule_id, 226 verbose=False) 227 228 return molecule_id
Assigns the molecule_id
of the pmb object given by pmb_type
Arguments:
- name(
str
): Label of the molecule type to be created.name
must be defined inpmb.df
- pmb_type(
str
): pmb_object_type to assign themolecule_id
- molecule_index (
int
): index of the currentpmb_object_type
to assign themolecule_id
- used_molecules_id (
lst
): list with themolecule_id
values already used.
Returns:
molecule_id(
int
): Id of the molecule
230 def calculate_center_of_mass_of_molecule(self, molecule_id, espresso_system): 231 """ 232 Calculates the center of the molecule with a given molecule_id. 233 234 Args: 235 molecule_id(`int`): Id of the molecule whose center of mass is to be calculated. 236 espresso_system(`obj`): Instance of a system object from the espressomd library. 237 238 Returns: 239 center_of_mass(`lst`): Coordinates of the center of mass. 240 """ 241 total_beads = 0 242 center_of_mass = np.zeros(3) 243 axis_list = [0,1,2] 244 particle_id_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 245 for pid in particle_id_list: 246 total_beads +=1 247 for axis in axis_list: 248 center_of_mass [axis] += espresso_system.part.by_id(pid).pos[axis] 249 center_of_mass = center_of_mass /total_beads 250 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(
obj
): Instance of a system object from the espressomd library.
Returns:
center_of_mass(
lst
): Coordinates of the center of mass.
252 def calculate_HH(self, molecule_name, pH_list=None, pka_set=None): 253 """ 254 Calculates the charge per molecule according to the ideal Henderson-Hasselbalch titration curve 255 for molecules with the name `molecule_name`. 256 257 Args: 258 molecule_name (`str`): name of the molecule to calculate the ideal charge for 259 pH_list(`lst`): pH-values to calculate. 260 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 261 262 Returns: 263 Z_HH (`lst`): Henderson-Hasselbalch prediction of the charge of `sequence` in `pH_list` 264 265 Note: 266 - This function supports objects with pmb types: "molecule", "peptide" and "protein". 267 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 268 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 269 - This function should only be used for single-phase systems. For two-phase systems `pmb.calculate_HH_Donnan` should be used. 270 """ 271 if pH_list is None: 272 pH_list=np.linspace(2,12,50) 273 if pka_set is None: 274 pka_set=self.get_pka_set() 275 self.check_pka_set(pka_set=pka_set) 276 charge_map = self.get_charge_map() 277 Z_HH=[] 278 for pH_value in pH_list: 279 Z=0 280 index = self.df.loc[self.df['name'] == molecule_name].index[0].item() 281 residue_list = self.df.at [index,('residue_list','')] 282 sequence = self.df.at [index,('sequence','')] 283 if np.any(pd.isnull(sequence)): 284 # Molecule has no sequence 285 for residue in residue_list: 286 list_of_particles_in_residue = self.search_particles_in_residue(residue) 287 for particle in list_of_particles_in_residue: 288 if particle in pka_set.keys(): 289 if pka_set[particle]['acidity'] == 'acidic': 290 psi=-1 291 elif pka_set[particle]['acidity']== 'basic': 292 psi=+1 293 else: 294 psi=0 295 Z+=psi/(1+10**(psi*(pH_value-pka_set[particle]['pka_value']))) 296 Z_HH.append(Z) 297 else: 298 # Molecule has a sequence 299 if not isinstance(sequence, list): 300 # If the df has been read by file, the sequence needs to be parsed. 301 sequence = self.parse_sequence_from_file(sequence=sequence) 302 for name in sequence: 303 if name in pka_set.keys(): 304 if pka_set[name]['acidity'] == 'acidic': 305 psi=-1 306 elif pka_set[name]['acidity']== 'basic': 307 psi=+1 308 else: 309 psi=0 310 Z+=psi/(1+10**(psi*(pH_value-pka_set[name]['pka_value']))) 311 else: 312 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 313 Z+=charge_map[state_one_type] 314 Z_HH.append(Z) 315 316 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.
318 def calculate_HH_Donnan(self, c_macro, c_salt, pH_list=None, pka_set=None): 319 """ 320 Calculates the charge on the different molecules according to the Henderson-Hasselbalch equation coupled to the Donnan partitioning. 321 322 Args: 323 c_macro ('dict'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 324 c_salt ('float'): Salt concentration in the reservoir. 325 pH_list ('lst'): List of pH-values in the reservoir. 326 pka_set ('dict'): {"name": {"pka_value": pka, "acidity": acidity}}. 327 328 Returns: 329 {"charges_dict": {"name": charges}, "pH_system_list": pH_system_list, "partition_coefficients": partition_coefficients_list} 330 pH_system_list ('lst'): List of pH_values in the system. 331 partition_coefficients_list ('lst'): List of partition coefficients of cations. 332 333 Note: 334 - If no `pH_list` is given, 50 equispaced pH-values ranging from 2 to 12 are calculated 335 - If no `pka_set` is given, the pKa values are taken from `pmb.df` 336 - If `c_macro` does not contain all charged molecules in the system, this function is likely to provide the wrong result. 337 """ 338 if pH_list is None: 339 pH_list=np.linspace(2,12,50) 340 if pka_set is None: 341 pka_set=self.get_pka_set() 342 self.check_pka_set(pka_set=pka_set) 343 344 partition_coefficients_list = [] 345 pH_system_list = [] 346 Z_HH_Donnan={} 347 for key in c_macro: 348 Z_HH_Donnan[key] = [] 349 350 def calc_charges(c_macro, pH): 351 """ 352 Calculates the charges of the different kinds of molecules according to the Henderson-Hasselbalch equation. 353 354 Args: 355 c_macro ('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 356 pH ('float'): pH-value that is used in the HH equation. 357 358 Returns: 359 charge ('dict'): {"molecule_name": charge} 360 """ 361 charge = {} 362 for name in c_macro: 363 charge[name] = self.calculate_HH(name, [pH], pka_set)[0] 364 return charge 365 366 def calc_partition_coefficient(charge, c_macro): 367 """ 368 Calculates the partition coefficients of positive ions according to the ideal Donnan theory. 369 370 Args: 371 charge ('dict'): {"molecule_name": charge} 372 c_macro ('dic'): {"name": concentration} - A dict containing the concentrations of all charged macromolecular species in the system. 373 """ 374 nonlocal ionic_strength_res 375 charge_density = 0.0 376 for key in charge: 377 charge_density += charge[key] * c_macro[key] 378 return (-charge_density / (2 * ionic_strength_res) + np.sqrt((charge_density / (2 * ionic_strength_res))**2 + 1)).magnitude 379 380 for pH_value in pH_list: 381 # calculate the ionic strength of the reservoir 382 if pH_value <= 7.0: 383 ionic_strength_res = 10 ** (-pH_value) * self.units.mol/self.units.l + c_salt 384 elif pH_value > 7.0: 385 ionic_strength_res = 10 ** (-(14-pH_value)) * self.units.mol/self.units.l + c_salt 386 387 #Determine the partition coefficient of positive ions by solving the system of nonlinear, coupled equations 388 #consisting of the partition coefficient given by the ideal Donnan theory and the Henderson-Hasselbalch equation. 389 #The nonlinear equation is formulated for log(xi) since log-operations are not supported for RootResult objects. 390 equation = lambda logxi: logxi - np.log10(calc_partition_coefficient(calc_charges(c_macro, pH_value - logxi), c_macro)) 391 logxi = scipy.optimize.root_scalar(equation, bracket=[-1e2, 1e2], method="brentq") 392 partition_coefficient = 10**logxi.root 393 394 charges_temp = calc_charges(c_macro, pH_value-np.log10(partition_coefficient)) 395 for key in c_macro: 396 Z_HH_Donnan[key].append(charges_temp[key]) 397 398 pH_system_list.append(pH_value - np.log10(partition_coefficient)) 399 partition_coefficients_list.append(partition_coefficient) 400 401 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.
403 def calculate_initial_bond_length(self, bond_object, bond_type, epsilon, sigma, cutoff): 404 """ 405 Calculates the initial bond length that is used when setting up molecules, 406 based on the minimum of the sum of bonded and short-range (LJ) interactions. 407 408 Args: 409 bond_object (`cls`): instance of a bond object from espressomd library 410 bond_type (`str`): label identifying the used bonded potential 411 epsilon (`float`): LJ epsilon of the interaction between the particles 412 sigma (`float`): LJ sigma of the interaction between the particles 413 cutoff (`float`): cutoff-radius of the LJ interaction 414 """ 415 def truncated_lj_potential(x, epsilon, sigma, cutoff): 416 if x>cutoff: 417 return 0.0 418 else: 419 return 4*epsilon*((sigma/x)**12-(sigma/x)**6) - 4*epsilon*((sigma/cutoff)**12-(sigma/cutoff)**6) 420 421 epsilon_red=epsilon.to('reduced_energy').magnitude 422 sigma_red=sigma.to('reduced_length').magnitude 423 cutoff_red=cutoff.to('reduced_length').magnitude 424 425 if bond_type == "harmonic": 426 r_0 = bond_object.params.get('r_0') 427 k = bond_object.params.get('k') 428 l0 = scipy.optimize.minimize(lambda x: 0.5*k*(x-r_0)**2 + truncated_lj_potential(x, epsilon_red, sigma_red, cutoff_red), x0=r_0).x 429 elif bond_type == "FENE": 430 r_0 = bond_object.params.get('r_0') 431 k = bond_object.params.get('k') 432 d_r_max = bond_object.params.get('d_r_max') 433 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), x0=1.0).x 434 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 (
cls
): instance of a bond object from espressomd library - bond_type (
str
): label identifying the used bonded potential - epsilon (
float
): LJ epsilon of the interaction between the particles - sigma (
float
): LJ sigma of the interaction between the particles - cutoff (
float
): cutoff-radius of the LJ interaction
436 def calculate_net_charge (self, espresso_system, molecule_name): 437 ''' 438 Calculates the net charge per molecule of molecules with `name` = molecule_name. 439 Returns the net charge per molecule and a maps with the net charge per residue and molecule. 440 441 Args: 442 espresso_system: system information 443 molecule_name (str): name of the molecule to calculate the net charge 444 445 Returns: 446 {"mean": mean_net_charge, "molecules": {mol_id: net_charge_of_mol_id, }, "residues": {res_id: net_charge_of_res_id, }} 447 448 Note: 449 - The net charge of the molecule is averaged over all molecules of type `name` 450 - The net charge of each particle type is averaged over all particle of the same type in all molecules of type `name` 451 ''' 452 valid_pmb_types = ["molecule", "protein"] 453 pmb_type=self.df.loc[self.df['name']==molecule_name].pmb_type.values[0] 454 if pmb_type not in valid_pmb_types: 455 raise ValueError("The pyMBE object with name {molecule_name} has a pmb_type {pmb_type}. This function only supports pyMBE types {valid_pmb_types}") 456 457 id_map = self.get_particle_id_map(object_name=molecule_name) 458 def create_charge_map(espresso_system,id_map,label): 459 charge_map={} 460 for super_id in id_map[label].keys(): 461 net_charge=0 462 for pid in id_map[label][super_id]: 463 net_charge+=espresso_system.part.by_id(pid).q 464 charge_map[super_id]=net_charge 465 return charge_map 466 net_charge_molecules=create_charge_map(label="molecule_map", 467 espresso_system=espresso_system, 468 id_map=id_map) 469 net_charge_residues=create_charge_map(label="residue_map", 470 espresso_system=espresso_system, 471 id_map=id_map) 472 mean_charge=np.mean(np.array(list(net_charge_molecules.values()))) 473 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: system information
- molecule_name (str): name of the molecule to calculate the net charge
Returns:
{"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
475 def center_molecule_in_simulation_box(self, molecule_id, espresso_system): 476 """ 477 Centers the pmb object matching `molecule_id` in the center of the simulation box in `espresso_md`. 478 479 Args: 480 molecule_id(`int`): Id of the molecule to be centered. 481 espresso_system(`obj`): Instance of a system object from the espressomd library. 482 """ 483 if len(self.df.loc[self.df['molecule_id']==molecule_id].pmb_type) == 0: 484 raise ValueError("The provided molecule_id is not present in the pyMBE dataframe.") 485 center_of_mass = self.calculate_center_of_mass_of_molecule(molecule_id=molecule_id,espresso_system=espresso_system) 486 box_center = [espresso_system.box_l[0]/2.0, 487 espresso_system.box_l[1]/2.0, 488 espresso_system.box_l[2]/2.0] 489 particle_id_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 490 for pid in particle_id_list: 491 es_pos = espresso_system.part.by_id(pid).pos 492 espresso_system.part.by_id(pid).pos = es_pos - center_of_mass + box_center 493 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(
obj
): Instance of a system object from the espressomd library.
495 def check_dimensionality(self, variable, expected_dimensionality): 496 """ 497 Checks if the dimensionality of `variable` matches `expected_dimensionality`. 498 499 Args: 500 `variable`(`pint.Quantity`): Quantity to be checked. 501 `expected_dimensionality(`str`): Expected dimension of the variable. 502 503 Returns: 504 `bool`: `True` if the variable if of the expected dimensionality, `False` otherwise. 505 506 Note: 507 - `expected_dimensionality` takes dimensionality following the Pint standards [docs](https://pint.readthedocs.io/en/0.10.1/wrapping.html?highlight=dimensionality#checking-dimensionality). 508 - For example, to check for a variable corresponding to a velocity `expected_dimensionality = "[length]/[time]"` 509 """ 510 correct_dimensionality=variable.check(f"{expected_dimensionality}") 511 if not correct_dimensionality: 512 raise ValueError(f"The variable {variable} should have a dimensionality of {expected_dimensionality}, instead the variable has a dimensionality of {variable.dimensionality}") 513 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]"
515 def check_if_df_cell_has_a_value(self, index,key): 516 """ 517 Checks if a cell in the `pmb.df` at the specified index and column has a value. 518 519 Args: 520 index(`int`): Index of the row to check. 521 key(`str`): Column label to check. 522 523 Returns: 524 `bool`: `True` if the cell has a value, `False` otherwise. 525 """ 526 idx = pd.IndexSlice 527 import warnings 528 with warnings.catch_warnings(): 529 warnings.simplefilter("ignore") 530 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.
532 def check_if_name_is_defined_in_df(self, name, pmb_type_to_be_defined): 533 """ 534 Checks if `name` is defined in `pmb.df`. 535 536 Args: 537 name(`str`): label to check if defined in `pmb.df`. 538 pmb_type_to_be_defined(`str`): pmb object type corresponding to `name`. 539 540 Returns: 541 `bool`: `True` for success, `False` otherwise. 542 """ 543 if name in self.df['name'].unique(): 544 current_object_type = self.df[self.df['name']==name].pmb_type.values[0] 545 if current_object_type != pmb_type_to_be_defined: 546 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") 547 return True 548 else: 549 return False
Checks if name
is defined in pmb.df
.
Arguments:
- name(
str
): label to check if defined inpmb.df
. - pmb_type_to_be_defined(
str
): pmb object type corresponding toname
.
Returns:
bool
:True
for success,False
otherwise.
552 def check_pka_set(self, pka_set): 553 """ 554 Checks that `pka_set` has the formatting expected by the pyMBE library. 555 556 Args: 557 pka_set (dict): {"name" : {"pka_value": pka, "acidity": acidity}} 558 """ 559 required_keys=['pka_value','acidity'] 560 for required_key in required_keys: 561 for pka_entry in pka_set.values(): 562 if required_key not in pka_entry.keys(): 563 raise ValueError(f'missing a required key "{required_key}" in the following entry of pka_set "{pka_entry}"') 564 return
Checks that pka_set
has the formatting expected by the pyMBE library.
Arguments:
- pka_set (dict): {"name" : {"pka_value": pka, "acidity": acidity}}
566 def clean_df_row(self, index, columns_keys_to_clean=("particle_id", "particle_id2", "residue_id", "molecule_id")): 567 """ 568 Cleans the columns of `pmb.df` in `columns_keys_to_clean` of the row with index `index` by assigning them a np.nan value. 569 570 Args: 571 index(`int`): Index of the row to clean. 572 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`]. 573 """ 574 for column_key in columns_keys_to_clean: 575 self.add_value_to_df(key=(column_key,''),index=index,new_value=np.nan, verbose=False) 576 return
Cleans the columns of pmb.df
in columns_keys_to_clean
of the row with index index
by assigning them a np.nan value.
Arguments:
- index(
int
): Index of the row to clean. - columns_keys_to_clean(
list
ofstr
, optional): List with the column keys to be cleaned. Defaults to [particle_id
,particle_id2
,residue_id
,molecule_id
].
578 def convert_columns_to_original_format (self, df): 579 """ 580 Converts the columns of the Dataframe to the original format in pyMBE. 581 582 Args: 583 df(`DataFrame`): dataframe with pyMBE information as a string 584 585 """ 586 587 columns_dtype_int = ['particle_id','particle_id2', 'residue_id','molecule_id', 'model',('state_one','es_type'),('state_two','es_type'),('state_one','charge'),('state_two','charge') ] 588 589 columns_with_units = ['sigma', 'epsilon', 'cutoff', 'offset'] 590 591 columns_with_list_or_dict = ['residue_list','side_chains', 'parameters_of_the_potential','sequence'] 592 593 for column_name in columns_dtype_int: 594 df[column_name] = df[column_name].astype(object) 595 596 for column_name in columns_with_list_or_dict: 597 if df[column_name].isnull().all(): 598 df[column_name] = df[column_name].astype(object) 599 else: 600 df[column_name] = df[column_name].apply(lambda x: ast.literal_eval(str(x)) if pd.notnull(x) else x) 601 602 for column_name in columns_with_units: 603 df[column_name] = df[column_name].apply(lambda x: self.create_variable_with_units(x) if pd.notnull(x) else x) 604 605 df['bond_object'] = df['bond_object'].apply(lambda x: self.convert_str_to_bond_object(x) if pd.notnull(x) else x) 606 607 return df
Converts the columns of the Dataframe to the original format in pyMBE.
Arguments:
- df(
DataFrame
): dataframe with pyMBE information as a string
609 def convert_str_to_bond_object (self, bond_str): 610 611 """ 612 Convert a row read as a `str` to the corresponding bond object. There are two supported bonds: HarmonicBond and FeneBond 613 614 Args: 615 bond_str (`str`): string with the information of a bond object 616 617 Returns: 618 bond_object(`obj`): EsPRESSo bond object 619 """ 620 621 from espressomd.interactions import HarmonicBond 622 from espressomd.interactions import FeneBond 623 624 supported_bonds = ['HarmonicBond', 'FeneBond'] 625 626 for bond in supported_bonds: 627 variable = re.subn(f'{bond}', '', bond_str) 628 if variable[1] == 1: 629 params = ast.literal_eval(variable[0]) 630 if bond == 'HarmonicBond': 631 bond_object = HarmonicBond(r_cut =params['r_cut'], k = params['k'], r_0=params['r_0']) 632 elif bond == 'FeneBond': 633 bond_object = FeneBond(k = params['k'], d_r_max =params['d_r_max'], r_0=params['r_0']) 634 635 return bond_object
Convert a row read as a str
to the corresponding bond object. There are two supported bonds: HarmonicBond and FeneBond
Arguments:
- bond_str (
str
): string with the information of a bond object
Returns:
bond_object(
obj
): EsPRESSo bond object
637 def copy_df_entry(self, name, column_name, number_of_copies): 638 ''' 639 Creates 'number_of_copies' of a given 'name' in `pymbe.df`. 640 641 Args: 642 name(`str`): Label of the particle/residue/molecule type to be created. `name` must be defined in `pmb.df` 643 column_name(`str`): Column name to use as a filter. 644 number_of_copies(`int`): number of copies of `name` to be created. 645 646 Note: 647 - Currently, column_name only supports "particle_id", "particle_id2", "residue_id" and "molecule_id" 648 ''' 649 650 valid_column_names=["particle_id", "residue_id", "molecule_id", "particle_id2" ] 651 if column_name not in valid_column_names: 652 raise ValueError(f"{column_name} is not a valid column_name, currently only the following are supported: {valid_column_names}") 653 df_by_name = self.df.loc[self.df.name == name] 654 if number_of_copies != 1: 655 if df_by_name[column_name].isnull().values.any(): 656 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies-1), ignore_index=True) 657 else: 658 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 659 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 660 df_by_name_repeated[column_name] =np.NaN 661 # Concatenate the new particle rows to `df` 662 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 663 else: 664 if not df_by_name[column_name].isnull().values.any(): 665 df_by_name = df_by_name[df_by_name.index == df_by_name.index.min()] 666 df_by_name_repeated = pd.concat ([df_by_name]*(number_of_copies), ignore_index=True) 667 df_by_name_repeated[column_name] =np.NaN 668 self.df = pd.concat ([self.df,df_by_name_repeated], ignore_index=True) 669 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"
671 def create_added_salt (self, espresso_system, cation_name, anion_name, c_salt, verbose=True): 672 """ 673 Creates a `c_salt` concentration of `cation_name` and `anion_name` ions into the `espresso_system`. 674 675 Args: 676 espresso_system (`obj`): instance of an espresso system object. 677 cation_name(`str`): `name` of a particle with a positive charge. 678 anion_name(`str`): `name` of a particle with a negative charge. 679 c_salt(`float`): Salt concentration. 680 verbose(`bool`): switch to activate/deactivate verbose. Defaults to True. 681 682 Returns: 683 c_salt_calculated (float): Calculated salt concentration added to `espresso_system`. 684 """ 685 cation_name_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.values[0] 686 anion_name_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.values[0] 687 if cation_name_charge <= 0: 688 raise ValueError('ERROR cation charge must be positive, charge ',cation_name_charge) 689 if anion_name_charge >= 0: 690 raise ValueError('ERROR anion charge must be negative, charge ', anion_name_charge) 691 # Calculate the number of ions in the simulation box 692 volume=self.units.Quantity(espresso_system.volume(), 'reduced_length**3') 693 if c_salt.check('[substance] [length]**-3'): 694 N_ions= int((volume*c_salt.to('mol/reduced_length**3')*self.N_A).magnitude) 695 c_salt_calculated=N_ions/(volume*self.N_A) 696 elif c_salt.check('[length]**-3'): 697 N_ions= int((volume*c_salt.to('reduced_length**-3')*self.N_A).magnitude) 698 c_salt_calculated=N_ions/volume 699 else: 700 raise ValueError('Unknown units for c_salt, please provided it in [mol / volume] or [particle / volume]', c_salt) 701 N_cation = N_ions*abs(anion_name_charge) 702 N_anion = N_ions*abs(cation_name_charge) 703 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=N_cation) 704 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=N_anion) 705 if verbose: 706 print(f"\n Added salt concentration of {c_salt_calculated.to('mol/L')} given by {N_cation} cations and {N_anion} anions") 707 return c_salt_calculated
Creates a c_salt
concentration of cation_name
and anion_name
ions into the espresso_system
.
Arguments:
- espresso_system (
obj
): instance of an espresso system object. - cation_name(
str
):name
of a particle with a positive charge. - anion_name(
str
):name
of a particle with a negative charge. - c_salt(
float
): Salt concentration. - verbose(
bool
): switch to activate/deactivate verbose. Defaults to True.
Returns:
c_salt_calculated (float): Calculated salt concentration added to
espresso_system
.
709 def create_bond_in_espresso(self, bond_type, bond_parameters): 710 ''' 711 Creates either a harmonic or a FENE bond in ESPREesSo 712 713 Args: 714 bond_type (`str`) : label to identify the potential to model the bond. 715 bond_parameters (`dict`) : parameters of the potential of the bond. 716 717 Note: 718 Currently, only HARMONIC and FENE bonds are supported. 719 720 For a HARMONIC bond the dictionary must contain: 721 722 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 723 using the `pmb.units` UnitRegistry. 724 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 725 the `pmb.units` UnitRegistry. 726 727 For a FENE bond the dictionay must additionally contain: 728 729 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 730 units of length using the `pmb.units` UnitRegistry. Default 'None'. 731 732 Returns: 733 bond_object (`obj`): a harmonic or a FENE bond object in ESPREesSo 734 ''' 735 from espressomd import interactions 736 737 valid_bond_types = ["harmonic", "FENE"] 738 739 if 'k' in bond_parameters: 740 bond_magnitude = bond_parameters['k'].to('reduced_energy / reduced_length**2') 741 else: 742 raise ValueError("Magnitud of the potential (k) is missing") 743 744 if bond_type == 'harmonic': 745 if 'r_0' in bond_parameters: 746 bond_length = bond_parameters['r_0'].to('reduced_length') 747 else: 748 raise ValueError("Equilibrium bond length (r_0) is missing") 749 bond_object = interactions.HarmonicBond(k = bond_magnitude.magnitude, 750 r_0 = bond_length.magnitude) 751 elif bond_type == 'FENE': 752 if 'r_0' in bond_parameters: 753 bond_length = bond_parameters['r_0'].to('reduced_length').magnitude 754 else: 755 print("WARNING: No value provided for r_0. Defaulting to r_0 = 0") 756 bond_length=0 757 if 'd_r_max' in bond_parameters: 758 max_bond_stret = bond_parameters['d_r_max'].to('reduced_length') 759 else: 760 raise ValueError("Maximal stretching length (d_r_max) is missing") 761 bond_object = interactions.FeneBond(r_0 = bond_length, 762 k = bond_magnitude.magnitude, 763 d_r_max = max_bond_stret.magnitude) 764 else: 765 raise ValueError(f"Bond type {bond_type} currently not implemented in pyMBE, valid types are {valid_bond_types}") 766 767 return bond_object
Creates either a harmonic or a FENE bond in ESPREesSo
Arguments:
- bond_type (
str
) : label to identify the potential to model the bond. - bond_parameters (
dict
) : parameters of the potential of the bond.
Note:
Currently, only HARMONIC and FENE bonds are supported.
For a HARMONIC bond the dictionary must contain:
- k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 using the `pmb.units` UnitRegistry. - r_0 (`obj`) : Equilibrium bond length. It should have units of length using the `pmb.units` UnitRegistry.
For a FENE bond the dictionay must additionally contain:
- d_r_max (`obj`): Maximal stretching length for FENE. It should have units of length using the `pmb.units` UnitRegistry. Default 'None'.
Returns:
bond_object (
obj
): a harmonic or a FENE bond object in ESPREesSo
770 def create_counterions(self, object_name, cation_name, anion_name, espresso_system,verbose=True): 771 """ 772 Creates particles of `cation_name` and `anion_name` in `espresso_system` to counter the net charge of `pmb_object`. 773 774 Args: 775 object_name(`str`): `name` of a pymbe object. 776 espresso_system(`obj`): Instance of a system object from the espressomd library. 777 cation_name(`str`): `name` of a particle with a positive charge. 778 anion_name(`str`): `name` of a particle with a negative charge. 779 verbose(`bool`): switch to activate/deactivate verbose. Defaults to True. 780 781 Returns: 782 counterion_number (dict): {"name": number} 783 """ 784 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.iloc[0] 785 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.iloc[0] 786 object_ids = self.get_particle_id_map(object_name=object_name)["all"] 787 counterion_number={} 788 object_charge={} 789 for name in ['positive', 'negative']: 790 object_charge[name]=0 791 for id in object_ids: 792 if espresso_system.part.by_id(id).q > 0: 793 object_charge['positive']+=1*(np.abs(espresso_system.part.by_id(id).q )) 794 elif espresso_system.part.by_id(id).q < 0: 795 object_charge['negative']+=1*(np.abs(espresso_system.part.by_id(id).q )) 796 if object_charge['positive'] % abs(anion_charge) == 0: 797 counterion_number[anion_name]=int(object_charge['positive']/abs(anion_charge)) 798 else: 799 raise ValueError('The number of positive charges in the pmb_object must be divisible by the charge of the anion') 800 if object_charge['negative'] % abs(cation_charge) == 0: 801 counterion_number[cation_name]=int(object_charge['negative']/cation_charge) 802 else: 803 raise ValueError('The number of negative charges in the pmb_object must be divisible by the charge of the cation') 804 if counterion_number[cation_name] > 0: 805 self.create_particle(espresso_system=espresso_system, name=cation_name, number_of_particles=counterion_number[cation_name]) 806 else: 807 counterion_number[cation_name]=0 808 if counterion_number[anion_name] > 0: 809 self.create_particle(espresso_system=espresso_system, name=anion_name, number_of_particles=counterion_number[anion_name]) 810 else: 811 counterion_number[anion_name] = 0 812 if verbose: 813 print('The following counter-ions have been created: ') 814 for name in counterion_number.keys(): 815 print(f'Ion type: {name} created number: {counterion_number[name]}') 816 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(
obj
): Instance of a system object from the espressomd library. - cation_name(
str
):name
of a particle with a positive charge. - anion_name(
str
):name
of a particle with a negative charge. - verbose(
bool
): switch to activate/deactivate verbose. Defaults to True.
Returns: counterion_number (dict): {"name": number}
818 def create_molecule(self, name, number_of_molecules, espresso_system, list_of_first_residue_positions=None, use_default_bond=False): 819 """ 820 Creates `number_of_molecules` molecule of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 821 822 Args: 823 name(`str`): Label of the molecule type to be created. `name` must be defined in `pmb.df` 824 espresso_system(`obj`): Instance of a system object from espressomd library. 825 number_of_molecules(`int`): Number of molecules of type `name` to be created. 826 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 827 use_default_bond(`bool`, optional): Controls if a bond of type `default` is used to bond particle with undefined bonds in `pymbe.df` 828 Returns: 829 830 molecules_info (`dict`): {molecule_id: {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids": [particle_id1, ...]}}} 831 """ 832 if list_of_first_residue_positions is not None: 833 for item in list_of_first_residue_positions: 834 if not isinstance(item, list): 835 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.") 836 elif len(item) != 3: 837 raise ValueError("The provided input position is formatted wrong. The elements in the provided list does not have 3 coordinates, corresponding to xyz coord.") 838 839 if len(list_of_first_residue_positions) != number_of_molecules: 840 raise ValueError(f"Number of positions provided in {list_of_first_residue_positions} does not match number of molecules desired, {number_of_molecules}") 841 if number_of_molecules <= 0: 842 return 843 if not self.check_if_name_is_defined_in_df(name=name, 844 pmb_type_to_be_defined='molecule'): 845 raise ValueError(f"{name} must correspond to a label of a pmb_type='molecule' defined on df") 846 first_residue = True 847 molecules_info = {} 848 residue_list = self.df[self.df['name']==name].residue_list.values [0] 849 850 self.copy_df_entry(name=name, 851 column_name='molecule_id', 852 number_of_copies=number_of_molecules) 853 854 molecules_index = np.where(self.df['name']==name) 855 molecule_index_list =list(molecules_index[0])[-number_of_molecules:] 856 used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist() 857 pos_index = 0 858 for molecule_index in molecule_index_list: 859 molecule_id = self.assign_molecule_id(name=name,pmb_type='molecule',used_molecules_id=used_molecules_id,molecule_index=molecule_index) 860 molecules_info[molecule_id] = {} 861 for residue in residue_list: 862 if first_residue: 863 if list_of_first_residue_positions is None: 864 residue_position = None 865 else: 866 for item in list_of_first_residue_positions: 867 residue_position = [np.array(list_of_first_residue_positions[pos_index])] 868 # Generate an arbitrary random unit vector 869 backbone_vector = self.generate_random_points_in_a_sphere(center=[0,0,0], 870 radius=1, 871 n_samples=1, 872 on_surface=True)[0] 873 residues_info = self.create_residue(name=residue, 874 number_of_residues=1, 875 espresso_system=espresso_system, 876 central_bead_position=residue_position, 877 use_default_bond= use_default_bond, 878 backbone_vector=backbone_vector) 879 residue_id = next(iter(residues_info)) 880 for index in self.df[self.df['residue_id']==residue_id].index: 881 self.add_value_to_df(key=('molecule_id',''), 882 index=int (index), 883 new_value=molecule_id) 884 central_bead_id = residues_info[residue_id]['central_bead_id'] 885 previous_residue = residue 886 residue_position = espresso_system.part.by_id(central_bead_id).pos 887 previous_residue_id = central_bead_id 888 first_residue = False 889 else: 890 previous_central_bead_name=self.df[self.df['name']==previous_residue].central_bead.values[0] 891 new_central_bead_name=self.df[self.df['name']==residue].central_bead.values[0] 892 bond = self.search_bond(particle_name1=previous_central_bead_name, 893 particle_name2=new_central_bead_name, 894 hard_check=True, 895 use_default_bond=use_default_bond) 896 l0 = self.get_bond_length(particle_name1=previous_central_bead_name, 897 particle_name2=new_central_bead_name, 898 hard_check=True, 899 use_default_bond=use_default_bond) 900 residue_position = residue_position+backbone_vector*l0 901 residues_info = self.create_residue(name=residue, 902 number_of_residues=1, 903 espresso_system=espresso_system, 904 central_bead_position=[residue_position], 905 use_default_bond= use_default_bond, 906 backbone_vector=backbone_vector) 907 residue_id = next(iter(residues_info)) 908 for index in self.df[self.df['residue_id']==residue_id].index: 909 if not self.check_if_df_cell_has_a_value(index=index,key=('molecule_id','')): 910 self.df.at[index,'molecule_id'] = molecule_id 911 self.add_value_to_df(key=('molecule_id',''), 912 index=int (index), 913 new_value=molecule_id, 914 verbose=False) 915 central_bead_id = residues_info[residue_id]['central_bead_id'] 916 espresso_system.part.by_id(central_bead_id).add_bond((bond, previous_residue_id)) 917 self.add_bond_in_df(particle_id1=central_bead_id, 918 particle_id2=previous_residue_id, 919 use_default_bond=use_default_bond) 920 previous_residue_id = central_bead_id 921 previous_residue = residue 922 molecules_info[molecule_id][residue_id] = residues_info[residue_id] 923 first_residue = True 924 pos_index+=1 925 926 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(
obj
): 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 - 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, ...]}}}
928 def create_particle(self, name, espresso_system, number_of_particles, position=None, fix=False): 929 """ 930 Creates `number_of_particles` particles of type `name` into `espresso_system` and bookkeeps them into `pymbe.df`. 931 932 Args: 933 name (`str`): Label of the particle type to be created. `name` must be a `particle` defined in `pmb_df`. 934 espresso_system (`cls`): Instance of a system object from the espressomd library. 935 number_of_particles (`int`): Number of particles to be created. 936 position (list of [`float`,`float`,`float`], optional): Initial positions of the particles. If not given, particles are created in random positions. Defaults to None. 937 fix(`bool`, optional): Controls if the particle motion is frozen in the integrator, it is used to create rigid objects. Defaults to False. 938 Returns: 939 created_pid_list(`list` of `float`): List with the ids of the particles created into `espresso_system`. 940 """ 941 if number_of_particles <=0: 942 return 943 if not self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'): 944 raise ValueError(f"{name} must correspond to a label of a pmb_type='particle' defined on df") 945 # Copy the data of the particle `number_of_particles` times in the `df` 946 self.copy_df_entry(name=name,column_name='particle_id',number_of_copies=number_of_particles) 947 # Get information from the particle type `name` from the df 948 q = self.df.loc[self.df['name']==name].state_one.charge.values[0] 949 es_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 950 # Get a list of the index in `df` corresponding to the new particles to be created 951 index = np.where(self.df['name']==name) 952 index_list =list(index[0])[-number_of_particles:] 953 # Create the new particles into `espresso_system` 954 created_pid_list=[] 955 for index in range (number_of_particles): 956 df_index=int (index_list[index]) 957 self.clean_df_row(index=df_index) 958 if position is None: 959 particle_position = self.rng.random((1, 3))[0] *np.copy(espresso_system.box_l) 960 else: 961 particle_position = position[index] 962 if len(espresso_system.part.all()) == 0: 963 bead_id = 0 964 else: 965 bead_id = max (espresso_system.part.all().id) + 1 966 created_pid_list.append(bead_id) 967 968 if fix: 969 espresso_system.part.add (id=bead_id, pos = particle_position, type = es_type, q = q,fix =[fix,fix,fix]) 970 else: 971 espresso_system.part.add (id=bead_id, pos = particle_position, type = es_type, q = q) 972 self.add_value_to_df(key=('particle_id',''),index=df_index,new_value=bead_id, verbose=False) 973 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 (
cls
): 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
.
975 def create_pmb_object(self, name, number_of_objects, espresso_system, position=None, use_default_bond=False): 976 """ 977 Creates all `particle`s associated to `pmb object` into `espresso` a number of times equal to `number_of_objects`. 978 979 Args: 980 name(`str`): Unique label of the `pmb object` to be created. 981 number_of_objects(`int`): Number of `pmb object`s to be created. 982 espresso_system(`obj`): Instance of an espresso system object from espressomd library. 983 position(`list`): Coordinates where the particles should be created. 984 use_default_bond(`bool`,optional): Controls if a `default` bond is used to bond particles with undefined bonds in `pmb.df`. Defaults to `False`. 985 986 Note: 987 - 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. 988 """ 989 allowed_objects=['particle','residue','molecule'] 990 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 991 if pmb_type not in allowed_objects: 992 raise ValueError('Object type not supported, supported types are ', allowed_objects) 993 if pmb_type == 'particle': 994 self.create_particle(name=name, number_of_particles=number_of_objects, espresso_system=espresso_system, position=position) 995 elif pmb_type == 'residue': 996 self.create_residue(name=name, number_of_residues=number_of_objects, espresso_system=espresso_system, central_bead_position=position,use_default_bond=use_default_bond) 997 elif pmb_type == 'molecule': 998 self.create_molecule(name=name, number_of_molecules=number_of_objects, espresso_system=espresso_system, use_default_bond=use_default_bond, list_of_first_residue_positions=position) 999 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(
obj
): 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
.
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.
1001 def create_protein(self, name, number_of_proteins, espresso_system, topology_dict): 1002 """ 1003 Creates `number_of_proteins` molecules of type `name` into `espresso_system` at the coordinates in `positions` 1004 1005 Args: 1006 name(`str`): Label of the protein to be created. 1007 espresso_system(`obj`): Instance of a system object from the espressomd library. 1008 number_of_proteins(`int`): Number of proteins to be created. 1009 positions(`dict`): {'ResidueNumber': {'initial_pos': [], 'chain_id': ''}} 1010 """ 1011 1012 if number_of_proteins <=0: 1013 return 1014 if not self.check_if_name_is_defined_in_df(name=name, 1015 pmb_type_to_be_defined='protein'): 1016 raise ValueError(f"{name} must correspond to a name of a pmb_type='protein' defined on df") 1017 1018 self.copy_df_entry(name=name,column_name='molecule_id',number_of_copies=number_of_proteins) 1019 1020 protein_index = np.where(self.df['name']==name) 1021 protein_index_list =list(protein_index[0])[-number_of_proteins:] 1022 used_molecules_id = self.df.molecule_id.dropna().drop_duplicates().tolist() 1023 1024 box_half=espresso_system.box_l[0]/2.0 1025 1026 for molecule_index in protein_index_list: 1027 1028 molecule_id = self.assign_molecule_id (name=name,pmb_type='protein',used_molecules_id=used_molecules_id,molecule_index=molecule_index) 1029 1030 protein_center = self.generate_coordinates_outside_sphere(radius = 1, 1031 max_dist=box_half, 1032 n_samples=1, 1033 center=[box_half]*3)[0] 1034 1035 for residue in topology_dict.keys(): 1036 1037 residue_name = re.split(r'\d+', residue)[0] 1038 residue_number = re.split(r'(\d+)', residue)[1] 1039 residue_position = topology_dict[residue]['initial_pos'] 1040 position = residue_position + protein_center 1041 1042 particle_id = self.create_particle(name=residue_name, 1043 espresso_system=espresso_system, 1044 number_of_particles=1, 1045 position=[position], 1046 fix = True) 1047 1048 index = self.df[self.df['particle_id']==particle_id[0]].index.values[0] 1049 1050 self.add_value_to_df(key=('residue_id',''), 1051 index=int (index), 1052 new_value=residue_number) 1053 1054 self.add_value_to_df(key=('molecule_id',''), 1055 index=int (index), 1056 new_value=molecule_id) 1057 1058 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(
obj
): 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': ''}}
1060 def create_residue(self, name, espresso_system, number_of_residues, central_bead_position=None,use_default_bond=False, backbone_vector=None): 1061 """ 1062 Creates `number_of_residues` residues of type `name` into `espresso_system` and bookkeeps them into `pmb.df`. 1063 1064 Args: 1065 name(`str`): Label of the residue type to be created. `name` must be defined in `pmb.df` 1066 espresso_system(`obj`): Instance of a system object from espressomd library. 1067 number_of_residue(`int`): Number of residues of type `name` to be created. 1068 central_bead_position(`list` of `float`): Position of the central bead. 1069 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` 1070 backbone_vector (`list` of `float`): Backbone vector of the molecule. All side chains are created perpendicularly to `backbone_vector`. 1071 1072 Returns: 1073 residues_info (`dict`): {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1074 """ 1075 if number_of_residues <= 0: 1076 return 1077 if not self.check_if_name_is_defined_in_df(name=name, 1078 pmb_type_to_be_defined='residue'): 1079 raise ValueError(f"{name} must correspond to a label of a pmb_type='residue' defined on df") 1080 # Copy the data of a residue a `number_of_residues` times in the `df 1081 self.copy_df_entry(name=name, 1082 column_name='residue_id', 1083 number_of_copies=number_of_residues) 1084 residues_index = np.where(self.df['name']==name) 1085 residue_index_list =list(residues_index[0])[-number_of_residues:] 1086 # Internal bookkepping of the residue info (important for side-chain residues) 1087 # Dict structure {residue_id:{"central_bead_id":central_bead_id, "side_chain_ids":[particle_id1, ...]}} 1088 residues_info={} 1089 for residue_index in residue_index_list: 1090 self.clean_df_row(index=int(residue_index)) 1091 # Assign a residue_id 1092 if self.df['residue_id'].isnull().all(): 1093 residue_id=0 1094 else: 1095 residue_id = self.df['residue_id'].max() + 1 1096 self.add_value_to_df(key=('residue_id',''),index=int (residue_index),new_value=residue_id, verbose=False) 1097 # create the principal bead 1098 if self.df.loc[self.df['name']==name].central_bead.values[0] is np.NaN: 1099 raise ValueError("central_bead must contain a particle name") 1100 central_bead_name = self.df.loc[self.df['name']==name].central_bead.values[0] 1101 central_bead_id = self.create_particle(name=central_bead_name, 1102 espresso_system=espresso_system, 1103 position=central_bead_position, 1104 number_of_particles = 1)[0] 1105 central_bead_position=espresso_system.part.by_id(central_bead_id).pos 1106 #assigns same residue_id to the central_bead particle created. 1107 index = self.df[self.df['particle_id']==central_bead_id].index.values[0] 1108 self.df.at [index,'residue_id'] = residue_id 1109 # Internal bookkeeping of the central bead id 1110 residues_info[residue_id]={} 1111 residues_info[residue_id]['central_bead_id']=central_bead_id 1112 # create the lateral beads 1113 side_chain_list = self.df.loc[self.df.index[residue_index]].side_chains.values[0] 1114 side_chain_beads_ids = [] 1115 for side_chain_element in side_chain_list: 1116 if side_chain_element not in self.df.values: 1117 raise ValueError (side_chain_element +'is not defined') 1118 pmb_type = self.df[self.df['name']==side_chain_element].pmb_type.values[0] 1119 if pmb_type == 'particle': 1120 bond = self.search_bond(particle_name1=central_bead_name, 1121 particle_name2=side_chain_element, 1122 hard_check=True, 1123 use_default_bond=use_default_bond) 1124 l0 = self.get_bond_length(particle_name1=central_bead_name, 1125 particle_name2=side_chain_element, 1126 hard_check=True, 1127 use_default_bond=use_default_bond) 1128 1129 if backbone_vector is None: 1130 bead_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1131 radius=l0, 1132 n_samples=1, 1133 on_surface=True)[0] 1134 else: 1135 bead_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1136 magnitude=l0) 1137 1138 side_bead_id = self.create_particle(name=side_chain_element, 1139 espresso_system=espresso_system, 1140 position=[bead_position], 1141 number_of_particles=1)[0] 1142 index = self.df[self.df['particle_id']==side_bead_id].index.values[0] 1143 self.add_value_to_df(key=('residue_id',''), 1144 index=int (index), 1145 new_value=residue_id, 1146 verbose=False) 1147 side_chain_beads_ids.append(side_bead_id) 1148 espresso_system.part.by_id(central_bead_id).add_bond((bond, side_bead_id)) 1149 self.add_bond_in_df(particle_id1=central_bead_id, 1150 particle_id2=side_bead_id, 1151 use_default_bond=use_default_bond) 1152 elif pmb_type == 'residue': 1153 central_bead_side_chain = self.df[self.df['name']==side_chain_element].central_bead.values[0] 1154 bond = self.search_bond(particle_name1=central_bead_name, 1155 particle_name2=central_bead_side_chain, 1156 hard_check=True, 1157 use_default_bond=use_default_bond) 1158 l0 = self.get_bond_length(particle_name1=central_bead_name, 1159 particle_name2=central_bead_side_chain, 1160 hard_check=True, 1161 use_default_bond=use_default_bond) 1162 if backbone_vector is None: 1163 residue_position=self.generate_random_points_in_a_sphere(center=central_bead_position, 1164 radius=l0, 1165 n_samples=1, 1166 on_surface=True)[0] 1167 else: 1168 residue_position=central_bead_position+self.generate_trial_perpendicular_vector(vector=backbone_vector, 1169 magnitude=l0) 1170 lateral_residue_info = self.create_residue(name=side_chain_element, espresso_system=espresso_system, 1171 number_of_residues=1, central_bead_position=[residue_position],use_default_bond=use_default_bond) 1172 lateral_residue_dict=list(lateral_residue_info.values())[0] 1173 central_bead_side_chain_id=lateral_residue_dict['central_bead_id'] 1174 lateral_beads_side_chain_ids=lateral_residue_dict['side_chain_ids'] 1175 residue_id_side_chain=list(lateral_residue_info.keys())[0] 1176 # Change the residue_id of the residue in the side chain to the one of the biger residue 1177 index = self.df[(self.df['residue_id']==residue_id_side_chain) & (self.df['pmb_type']=='residue') ].index.values[0] 1178 self.add_value_to_df(key=('residue_id',''), 1179 index=int(index), 1180 new_value=residue_id, 1181 verbose=False) 1182 # Change the residue_id of the particles in the residue in the side chain 1183 side_chain_beads_ids+=[central_bead_side_chain_id]+lateral_beads_side_chain_ids 1184 for particle_id in side_chain_beads_ids: 1185 index = self.df[(self.df['particle_id']==particle_id) & (self.df['pmb_type']=='particle')].index.values[0] 1186 self.add_value_to_df(key=('residue_id',''), 1187 index=int (index), 1188 new_value=residue_id, 1189 verbose=False) 1190 espresso_system.part.by_id(central_bead_id).add_bond((bond, central_bead_side_chain_id)) 1191 self.add_bond_in_df(particle_id1=central_bead_id, 1192 particle_id2=central_bead_side_chain_id, 1193 use_default_bond=use_default_bond) 1194 # Internal bookkeeping of the side chain beads ids 1195 residues_info[residue_id]['side_chain_ids']=side_chain_beads_ids 1196 return residues_info
Creates number_of_residues
residues 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(
obj
): Instance of a system object from espressomd library. - number_of_residue(
int
): Number of residues of typename
to be created. - 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, ...]}}
1198 def create_variable_with_units(self, variable): 1199 """ 1200 Returns a pint object with the value and units defined in `variable`. 1201 1202 Args: 1203 variable(`dict` or `str`): {'value': value, 'units': units} 1204 Returns: 1205 variable_with_units(`obj`): variable with units using the pyMBE UnitRegistry. 1206 """ 1207 1208 if isinstance(variable, dict): 1209 1210 value=variable.pop('value') 1211 units=variable.pop('units') 1212 1213 elif isinstance(variable, str): 1214 1215 value = float(re.split(r'\s+', variable)[0]) 1216 units = re.split(r'\s+', variable)[1] 1217 1218 variable_with_units=value*self.units(units) 1219 1220 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.
1222 def define_AA_particles_in_sequence(self, sequence, sigma_dict=None): 1223 ''' 1224 Defines in `pmb.df` all the different particles in `sequence`. 1225 1226 Args: 1227 sequence(`lst`): Sequence of the peptide or protein. 1228 1229 Note: 1230 - It assumes that the names in `sequence` correspond to amino acid names using the standard one letter code. 1231 ''' 1232 1233 already_defined_AA=[] 1234 1235 for residue_name in sequence: 1236 if residue_name in already_defined_AA: 1237 continue 1238 self.define_particle (name=residue_name, q=0) 1239 1240 if sigma_dict: 1241 self.define_particles_parameter_from_dict(param_dict = sigma_dict, 1242 param_name = 'sigma') 1243 return
Defines in pmb.df
all the different particles in sequence
.
Arguments:
- sequence(
lst
): Sequence of the peptide or protein.
Note:
- It assumes that the names in
sequence
correspond to amino acid names using the standard one letter code.
1245 def define_AA_residues(self, sequence, model): 1246 """ 1247 Defines in `pmb.df` all the different residues in `sequence`. 1248 1249 Args: 1250 sequence(`lst`): Sequence of the peptide or protein. 1251 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1252 1253 Returns: 1254 residue_list (`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1255 """ 1256 residue_list = [] 1257 for residue_name in sequence: 1258 if model == '1beadAA': 1259 central_bead = residue_name 1260 side_chains = [] 1261 elif model == '2beadAA': 1262 if residue_name in ['c','n', 'G']: 1263 central_bead = residue_name 1264 side_chains = [] 1265 else: 1266 central_bead = 'CA' 1267 side_chains = [residue_name] 1268 if residue_name not in residue_list: 1269 self.define_residue(name = 'AA-'+residue_name, 1270 central_bead = central_bead, 1271 side_chains = side_chains) 1272 residue_list.append('AA-'+residue_name) 1273 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
.
1275 def define_bond(self, bond_type, bond_parameters, particle_pairs): 1276 1277 ''' 1278 Defines a pmb object of type `bond` in `pymbe.df`. 1279 1280 Args: 1281 bond_type (`str`) : label to identify the potential to model the bond. 1282 bond_parameters (`dict`) : parameters of the potential of the bond. 1283 particle_pairs (`lst`) : list of the `names` of the `particles` to be bonded. 1284 1285 Note: 1286 Currently, only HARMONIC and FENE bonds are supported. 1287 1288 For a HARMONIC bond the dictionary must contain the following parameters: 1289 1290 - k (`obj`) : Magnitude of the bond. It should have units of energy/length**2 1291 using the `pmb.units` UnitRegistry. 1292 - r_0 (`obj`) : Equilibrium bond length. It should have units of length using 1293 the `pmb.units` UnitRegistry. 1294 1295 For a FENE bond the dictionay must contain the same parameters as for a HARMONIC bond and: 1296 1297 - d_r_max (`obj`): Maximal stretching length for FENE. It should have 1298 units of length using the `pmb.units` UnitRegistry. Default 'None'. 1299 ''' 1300 1301 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1302 1303 for particle_name1, particle_name2 in particle_pairs: 1304 1305 lj_parameters=self.get_lj_parameters(particle_name1 = particle_name1, 1306 particle_name2 = particle_name2, 1307 combining_rule = 'Lorentz-Berthelot') 1308 1309 cutoff = 2**(1.0/6.0) * lj_parameters["sigma"] 1310 1311 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1312 bond_type = bond_type, 1313 epsilon = lj_parameters["epsilon"], 1314 sigma = lj_parameters["sigma"], 1315 cutoff = cutoff) 1316 index = len(self.df) 1317 for label in [particle_name1+'-'+particle_name2,particle_name2+'-'+particle_name1]: 1318 self.check_if_name_is_defined_in_df(name=label, pmb_type_to_be_defined="bond") 1319 self.df.at [index,'name']= particle_name1+'-'+particle_name2 1320 self.df.at [index,'bond_object'] = bond_object 1321 self.df.at [index,'l0'] = l0 1322 self.add_value_to_df(index=index, 1323 key=('pmb_type',''), 1324 new_value='bond') 1325 self.add_value_to_df(index=index, 1326 key=('parameters_of_the_potential',''), 1327 new_value=bond_object.get_params(), 1328 non_standard_value=True) 1329 1330 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 dictionay must contain the same parameters as for a HARMONIC bond and:
- d_r_max (`obj`): Maximal stretching length for FENE. It should have units of length using the `pmb.units` UnitRegistry. Default 'None'.
1333 def define_default_bond(self, bond_type, bond_parameters, epsilon=None, sigma=None, cutoff=None): 1334 """ 1335 Asigns `bond` in `pmb.df` as the default bond. 1336 The LJ parameters can be optionally provided to calculate the initial bond length 1337 1338 Args: 1339 bond_type (`str`) : label to identify the potential to model the bond. 1340 bond_parameters (`dict`) : parameters of the potential of the bond. 1341 sigma (`float`, optional) : LJ sigma of the interaction between the particles 1342 epsilon (`float`, optional): LJ epsilon for the interaction between the particles 1343 cutoff (`float`, optional) : cutoff-radius of the LJ interaction 1344 1345 Note: 1346 - Currently, only harmonic and FENE bonds are supported. 1347 """ 1348 1349 bond_object=self.create_bond_in_espresso(bond_type, bond_parameters) 1350 1351 if epsilon is None: 1352 epsilon=1*self.units('reduced_energy') 1353 if sigma is None: 1354 sigma=1*self.units('reduced_length') 1355 if cutoff is None: 1356 cutoff=2**(1.0/6.0)*self.units('reduced_length') 1357 l0 = self.calculate_initial_bond_length(bond_object = bond_object, 1358 bond_type = bond_type, 1359 epsilon = epsilon, 1360 sigma = sigma, 1361 cutoff = cutoff) 1362 1363 if self.check_if_name_is_defined_in_df(name='default',pmb_type_to_be_defined='bond'): 1364 return 1365 if len(self.df.index) != 0: 1366 index = max(self.df.index)+1 1367 else: 1368 index = 0 1369 self.df.at [index,'name'] = 'default' 1370 self.df.at [index,'bond_object'] = bond_object 1371 self.df.at [index,'l0'] = l0 1372 self.add_value_to_df(index = index, 1373 key = ('pmb_type',''), 1374 new_value = 'bond') 1375 self.add_value_to_df(index = index, 1376 key = ('parameters_of_the_potential',''), 1377 new_value = bond_object.get_params(), 1378 non_standard_value=True) 1379 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
Note:
- Currently, only harmonic and FENE bonds are supported.
1381 def define_epsilon_value_of_particles(self, eps_dict): 1382 ''' 1383 Defines the epsilon value of the particles in `eps_dict` into the `pmb.df`. 1384 1385 Args: 1386 eps_dict(`dict`): {'name': epsilon} 1387 ''' 1388 for residue in eps_dict.keys(): 1389 label_list = self.df[self.df['name'] == residue].index.tolist() 1390 for index in label_list: 1391 epsilon = eps_dict[residue] 1392 self.add_value_to_df(key= ('epsilon',''), 1393 index=int (index), 1394 new_value=epsilon) 1395 return
Defines the epsilon value of the particles in eps_dict
into the pmb.df
.
Arguments:
- eps_dict(
dict
): {'name': epsilon}
1397 def define_molecule(self, name, residue_list): 1398 """ 1399 Defines a pyMBE object of type `molecule` in `pymbe.df`. 1400 1401 Args: 1402 name (`str`): Unique label that identifies the `molecule`. 1403 residue_list (`list` of `str`): List of the `name`s of the `residue`s in the sequence of the `molecule`. 1404 """ 1405 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='molecule'): 1406 return 1407 index = len(self.df) 1408 self.df.at [index,'name'] = name 1409 self.df.at [index,'pmb_type'] = 'molecule' 1410 self.df.at [index,('residue_list','')] = residue_list 1411 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
.
1413 def define_particle(self, name, q=0, acidity='inert', pka=None, sigma=None, epsilon=None, cutoff=None, offset=None,verbose=True): 1414 """ 1415 Defines the properties of a particle object. 1416 1417 Args: 1418 name (`str`): Unique label that identifies this particle type. 1419 q (`int`, optional): Permanent charge of this particle type. Defaults to 0. 1420 acidity (`str`, optional): Identifies whether if the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to `inert`. 1421 pka (`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to None. 1422 sigma (`pint.Quantity`, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1423 cutoff (`pint.Quantity`, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1424 offset (`pint.Quantity`, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. 1425 epsilon (`pint.Quantity`, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to None. 1426 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 1427 1428 Note: 1429 - `sigma`, `cutoff` and `offset` must have a dimensitonality of `[length]` and should be defined using pmb.units. 1430 - `epsilon` must have a dimensitonality of `[energy]` and should be defined using pmb.units. 1431 - `cutoff` defaults to `2**(1./6.) reduced_length`. 1432 - `offset` defaults to 0. 1433 - The default setup corresponds to the Weeks−Chandler−Andersen (WCA) model, corresponding to purely steric interactions. 1434 - For more information on `sigma`, `epsilon`, `cutoff` and `offset` check `pmb.setup_lj_interactions()`. 1435 """ 1436 1437 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='particle'): 1438 index = self.df[self.df['name']==name].index[0] 1439 else: 1440 index = len(self.df) 1441 self.df.at [index, 'name'] = name 1442 self.df.at [index,'pmb_type'] = 'particle' 1443 1444 # If `cutoff` and `offset` are not defined, default them to 0 1445 1446 if cutoff is None: 1447 cutoff=self.units.Quantity(2**(1./6.), "reduced_length") 1448 if offset is None: 1449 offset=self.units.Quantity(0, "reduced_length") 1450 1451 # Define LJ parameters 1452 parameters_with_dimensionality={"sigma":{"value": sigma, "dimensionality": "[length]"}, 1453 "cutoff":{"value": cutoff, "dimensionality": "[length]"}, 1454 "offset":{"value": offset, "dimensionality": "[length]"}, 1455 "epsilon":{"value": epsilon, "dimensionality": "[energy]"},} 1456 1457 for parameter_key in parameters_with_dimensionality.keys(): 1458 if parameters_with_dimensionality[parameter_key]["value"] is not None: 1459 self.check_dimensionality(variable=parameters_with_dimensionality[parameter_key]["value"], 1460 expected_dimensionality=parameters_with_dimensionality[parameter_key]["dimensionality"]) 1461 self.add_value_to_df(key=(parameter_key,''), 1462 index=index, 1463 new_value=parameters_with_dimensionality[parameter_key]["value"], 1464 verbose=verbose) 1465 1466 # Define particle acid/base properties 1467 self.set_particle_acidity (name=name, 1468 acidity = acidity , 1469 default_charge=q, 1470 pka=pka, 1471 verbose=verbose) 1472 return
Defines the properties of a particle object.
Arguments:
- name (
str
): Unique label that identifies this particle type. - q (
int
, optional): Permanent charge of this particle type. Defaults to 0. - acidity (
str
, optional): Identifies whether if the particle isacidic
orbasic
, used to setup constant pH simulations. Defaults toinert
. - pka (
float
, optional): Ifparticle
is an acid or a base, it defines its pka-value. Defaults to None. - sigma (
pint.Quantity
, optional): Sigma parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. - cutoff (
pint.Quantity
, optional): Cutoff parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. - offset (
pint.Quantity
, optional): Offset parameter used to set up Lennard-Jones interactions for this particle type. Defaults to None. - epsilon (
pint.Quantity
, optional): Epsilon parameter used to setup Lennard-Jones interactions for this particle tipe. Defaults to None. - verbose (
bool
, optional): Switch to activate/deactivate verbose. Defaults to True.
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()
.
1474 def define_particles_parameter_from_dict(self, param_dict, param_name): 1475 ''' 1476 Defines the `param_name` value of the particles in `param_dict` into the `pmb.df`. 1477 1478 Args: 1479 param_dict(`dict`): {'name': `param_name` value} 1480 ''' 1481 for residue in param_dict.keys(): 1482 label_list = self.df[self.df['name'] == residue].index.tolist() 1483 for index in label_list: 1484 value = param_dict[residue] 1485 self.add_value_to_df(key= (param_name,''), 1486 index=int (index), 1487 new_value=value) 1488 return
Defines the param_name
value of the particles in param_dict
into the pmb.df
.
Arguments:
- param_dict(
dict
): {'name':param_name
value}
1490 def define_peptide(self, name, sequence, model): 1491 """ 1492 Defines a pyMBE object of type `peptide` in the `pymbe.df`. 1493 1494 Args: 1495 name (`str`): Unique label that identifies the `peptide`. 1496 sequence (`string`): Sequence of the `peptide`. 1497 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1498 """ 1499 if not self.check_if_name_is_defined_in_df(name = name, pmb_type_to_be_defined='peptide'): 1500 valid_keys = ['1beadAA','2beadAA'] 1501 if model not in valid_keys: 1502 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1503 1504 clean_sequence = self.protein_sequence_parser(sequence=sequence) 1505 residue_list = self.define_AA_residues(sequence=clean_sequence,model=model) 1506 1507 self.define_molecule(name = name, residue_list=residue_list) 1508 index = self.df.loc[self.df['name'] == name].index.item() 1509 self.df.at [index,'model'] = model 1510 self.df.at [index,('sequence','')] = clean_sequence 1511 return
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.
1513 def define_protein(self, name,model, topology_dict): 1514 """ 1515 Defines a pyMBE object of type `protein` in `pymbe.df`. 1516 1517 Args: 1518 name (`str`): Unique label that identifies the `protein`. 1519 model (`string`): Model name. Currently only models with 1 bead '1beadAA' or with 2 beads '2beadAA' per amino acid are supported. 1520 topology_dict (`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 1521 """ 1522 1523 protein_seq_list = [] 1524 1525 valid_keys = ['1beadAA','2beadAA'] 1526 1527 if model not in valid_keys: 1528 raise ValueError('Invalid label for the peptide model, please choose between 1beadAA or 2beadAA') 1529 1530 if model == '2beadAA': 1531 self.define_particle(name='CA') 1532 1533 sigma_dict = {} 1534 1535 for residue in topology_dict.keys(): 1536 residue_name = re.split(r'\d+', residue)[0] 1537 1538 if residue_name not in sigma_dict.keys(): 1539 sigma_dict [residue_name] = topology_dict[residue]['sigma'] 1540 if residue_name not in ('CA', 'Ca'): 1541 protein_seq_list.append(residue_name) 1542 1543 protein_sequence = ''.join(protein_seq_list) 1544 clean_sequence = self.protein_sequence_parser(sequence=protein_sequence) 1545 1546 self.define_AA_particles_in_sequence (sequence=clean_sequence, sigma_dict = sigma_dict) 1547 residue_list = self.define_AA_residues(sequence=clean_sequence, model=model) 1548 1549 index = len(self.df) 1550 self.df.at [index,'name'] = name 1551 self.df.at [index,'pmb_type'] = 'protein' 1552 self.df.at [index,'model'] = model 1553 self.df.at [index,('sequence','')] = clean_sequence 1554 self.df.at [index,('residue_list','')] = residue_list 1555 1556 return
Defines a pyMBE object of type protein
in pymbe.df
.
Arguments:
- name (
str
): Unique label that identifies theprotein
. - 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, 'sigma': sigma_value}
1558 def define_residue(self, name, central_bead, side_chains): 1559 """ 1560 Defines a pyMBE object of type `residue` in `pymbe.df`. 1561 1562 Args: 1563 name (`str`): Unique label that identifies the `residue`. 1564 central_bead (`str`): `name` of the `particle` to be placed as central_bead of the `residue`. 1565 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. 1566 """ 1567 if self.check_if_name_is_defined_in_df(name=name,pmb_type_to_be_defined='residue'): 1568 return 1569 index = len(self.df) 1570 self.df.at [index, 'name'] = name 1571 self.df.at [index,'pmb_type'] = 'residue' 1572 self.df.at [index,'central_bead'] = central_bead 1573 self.df.at [index,('side_chains','')] = side_chains 1574 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.
1576 def destroy_pmb_object_in_system(self, name, espresso_system): 1577 """ 1578 Destroys all particles associated with `name` in `espresso_system` amd removes the destroyed pmb_objects from `pmb.df` 1579 1580 Args: 1581 name (str): Label of the pmb object to be destroyed. The pmb object must be defined in `pymbe.df`. 1582 espresso_system (cls): Instance of a system class from espressomd library. 1583 1584 Note: 1585 - 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. 1586 """ 1587 allowed_objects = ['particle','residue','molecule'] 1588 pmb_type = self.df.loc[self.df['name']==name].pmb_type.values[0] 1589 if pmb_type not in allowed_objects: 1590 raise ValueError('Object type not supported, supported types are ', allowed_objects) 1591 if pmb_type == 'particle': 1592 particles_index = self.df.index[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1593 & (self.df['molecule_id'].isna())] 1594 particle_ids_list= self.df.loc[(self.df['name'] == name) & (self.df['residue_id'].isna()) 1595 & (self.df['molecule_id'].isna())].particle_id.tolist() 1596 for particle_id in particle_ids_list: 1597 espresso_system.part.by_id(particle_id).remove() 1598 self.df = self.df.drop(index=particles_index) 1599 if pmb_type == 'residue': 1600 residues_id = self.df.loc[self.df['name']== name].residue_id.to_list() 1601 for residue_id in residues_id: 1602 particle_ids_list = self.df.loc[self.df['residue_id']==residue_id].particle_id.dropna().to_list() 1603 self.df = self.df.drop(self.df[self.df['residue_id'] == residue_id].index) 1604 for particle_id in particle_ids_list: 1605 espresso_system.part.by_id(particle_id).remove() 1606 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1607 if pmb_type == 'molecule': 1608 molecules_id = self.df.loc[self.df['name']== name].molecule_id.to_list() 1609 for molecule_id in molecules_id: 1610 particle_ids_list = self.df.loc[self.df['molecule_id']==molecule_id].particle_id.dropna().to_list() 1611 1612 self.df = self.df.drop(self.df[self.df['molecule_id'] == molecule_id].index) 1613 for particle_id in particle_ids_list: 1614 espresso_system.part.by_id(particle_id).remove() 1615 self.df= self.df.drop(self.df[self.df['particle_id']==particle_id].index) 1616 1617 self.df.reset_index(drop=True,inplace=True) 1618 1619 return
Destroys all particles associated with name
in espresso_system
amd removes the destroyed pmb_objects from pmb.df
Arguments:
- name (str): Label of the pmb object to be destroyed. The pmb object must be defined in
pymbe.df
. - espresso_system (cls): 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.
1621 def determine_reservoir_concentrations(self, pH_res, c_salt_res, activity_coefficient_monovalent_pair, max_number_sc_runs=200): 1622 """ 1623 Determines the concentrations of the various species in the reservoir for given values of the pH and salt concentration. 1624 To do this, a system of nonlinear equations involving the pH, the ionic product of water, the activity coefficient of an 1625 ion pair and the concentrations of the various species is solved numerically using a self-consistent approach. 1626 More details can be found in chapter 5.3 of Landsgesell (doi.org/10.18419/opus-10831). 1627 This is a modified version of the code by Landsgesell et al. (doi.org/10.18419/darus-2237). 1628 1629 Args: 1630 pH_res ('float'): pH-value in the reservoir. 1631 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 1632 activity_coefficient_monovalent_pair ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 1633 1634 Returns: 1635 cH_res ('float'): Concentration of H+ ions. 1636 cOH_res ('float'): Concentration of OH- ions. 1637 cNa_res ('float'): Concentration of Na+ ions. 1638 cCl_res ('float'): Concentration of Cl- ions. 1639 """ 1640 1641 self_consistent_run = 0 1642 cH_res = 10**(-pH_res) * self.units.mol/self.units.l #initial guess for the concentration of H+ 1643 1644 def determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res): 1645 #Calculate and initial guess for the concentrations of various species based on ideal gas estimate 1646 cOH_res = self.Kw / cH_res 1647 cNa_res = None 1648 cCl_res = None 1649 if cOH_res>=cH_res: 1650 #adjust the concentration of sodium if there is excess OH- in the reservoir: 1651 cNa_res = c_salt_res + (cOH_res-cH_res) 1652 cCl_res = c_salt_res 1653 else: 1654 # adjust the concentration of chloride if there is excess H+ in the reservoir 1655 cCl_res = c_salt_res + (cH_res-cOH_res) 1656 cNa_res = c_salt_res 1657 1658 def calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res): 1659 nonlocal max_number_sc_runs, self_consistent_run 1660 if self_consistent_run<max_number_sc_runs: 1661 self_consistent_run+=1 1662 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1663 cOH_res = self.Kw / (cH_res * activity_coefficient_monovalent_pair(ionic_strength_res)) 1664 if cOH_res>=cH_res: 1665 #adjust the concentration of sodium if there is excess OH- in the reservoir: 1666 cNa_res = c_salt_res + (cOH_res-cH_res) 1667 cCl_res = c_salt_res 1668 else: 1669 # adjust the concentration of chloride if there is excess H+ in the reservoir 1670 cCl_res = c_salt_res + (cH_res-cOH_res) 1671 cNa_res = c_salt_res 1672 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 1673 else: 1674 return cH_res, cOH_res, cNa_res, cCl_res 1675 return calculate_concentrations_self_consistently(cH_res, cOH_res, cNa_res, cCl_res) 1676 1677 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 1678 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1679 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 1680 1681 while abs(determined_pH-pH_res)>1e-6: 1682 if determined_pH > pH_res: 1683 cH_res *= 1.005 1684 else: 1685 cH_res /= 1.003 1686 cH_res, cOH_res, cNa_res, cCl_res = determine_reservoir_concentrations_selfconsistently(cH_res, c_salt_res) 1687 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 1688 determined_pH = -np.log10(cH_res.to('mol/L').magnitude * np.sqrt(activity_coefficient_monovalent_pair(ionic_strength_res))) 1689 self_consistent_run=0 1690 1691 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 ('float'): 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 ('float'): Concentration of H+ ions. cOH_res ('float'): Concentration of OH- ions. cNa_res ('float'): Concentration of Na+ ions. cCl_res ('float'): Concentration of Cl- ions.
1693 def filter_df(self, pmb_type): 1694 """ 1695 Filters `pmb.df` and returns a sub-set of it containing only rows with pmb_object_type=`pmb_type` and non-NaN columns. 1696 1697 Args: 1698 pmb_type(`str`): pmb_object_type to filter in `pmb.df`. 1699 1700 Returns: 1701 pmb_type_df(`obj`): filtered `pmb.df`. 1702 """ 1703 pmb_type_df = self.df.loc[self.df['pmb_type']== pmb_type] 1704 pmb_type_df = pmb_type_df.dropna( axis=1, thresh=1) 1705 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(
obj
): filteredpmb.df
.
1707 def find_bond_key(self, particle_name1, particle_name2, use_default_bond=False): 1708 """ 1709 Searches for the `name` of the bond between `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 1710 1711 Args: 1712 particle_name1(`str`): label of the type of the first particle type of the bonded particles. 1713 particle_name2(`str`): label of the type of the second particle type of the bonded particles. 1714 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'. 1715 1716 Returns: 1717 bond_key (str): `name` of the bond between `particle_name1` and `particle_name2` 1718 1719 Note: 1720 - If `use_default_bond`=`True`, it returns "default" if no key is found. 1721 """ 1722 bond_keys = [particle_name1 +'-'+ particle_name2, particle_name2 +'-'+ particle_name1 ] 1723 bond_defined=False 1724 for bond_key in bond_keys: 1725 if bond_key in self.df.values: 1726 bond_defined=True 1727 correct_key=bond_key 1728 break 1729 if bond_defined: 1730 return correct_key 1731 elif use_default_bond: 1732 return 'default' 1733 else: 1734 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
Note:
- If
use_default_bond
=True
, it returns "default" if no key is found.
1736 def find_value_from_es_type(self, es_type, column_name): 1737 """ 1738 Finds a value in `pmb.df` for a `column_name` and `es_type` pair. 1739 1740 Args: 1741 es_type(`int`): value of the espresso type 1742 column_name(`str`): name of the column in `pymbe.df` 1743 1744 Returns: 1745 Value in `pymbe.df` matching `column_name` and `es_type` 1746 """ 1747 idx = pd.IndexSlice 1748 for state in ['state_one', 'state_two']: 1749 index = np.where(self.df[(state, 'es_type')] == es_type)[0] 1750 if len(index) > 0: 1751 if column_name == 'label': 1752 label = self.df.loc[idx[index[0]], idx[(state,column_name)]] 1753 return label 1754 else: 1755 column_name_value = self.df.loc[idx[index[0]], idx[(column_name,'')]] 1756 return column_name_value 1757 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
1759 def generate_coordinates_outside_sphere(self, center, radius, max_dist, n_samples): 1760 """ 1761 Generates coordinates outside a sphere centered at `center`. 1762 1763 Args: 1764 center(`lst`): Coordinates of the center of the sphere. 1765 radius(`float`): Radius of the sphere. 1766 max_dist(`float`): Maximum distance from the center of the spahre to generate coordinates. 1767 n_samples(`int`): Number of sample points. 1768 1769 Returns: 1770 coord_list(`lst`): Coordinates of the sample points. 1771 """ 1772 if not radius > 0: 1773 raise ValueError (f'The value of {radius} must be a positive value') 1774 if not radius < max_dist: 1775 raise ValueError(f'The min_dist ({radius} must be lower than the max_dist ({max_dist}))') 1776 coord_list = [] 1777 counter = 0 1778 while counter<n_samples: 1779 coord = self.generate_random_points_in_a_sphere(center=center, 1780 radius=max_dist, 1781 n_samples=1)[0] 1782 if np.linalg.norm(coord-np.asarray(center))>=radius: 1783 coord_list.append (coord) 1784 counter += 1 1785 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.
1787 def generate_random_points_in_a_sphere(self, center, radius, n_samples, on_surface=False): 1788 """ 1789 Uniformly samples points from a hypersphere. If on_surface is set to True, the points are 1790 uniformly sampled from the surface of the hypersphere. 1791 1792 Args: 1793 center(`lst`): Array with the coordinates of the center of the spheres. 1794 radius(`float`): Radius of the sphere. 1795 n_samples(`int`): Number of sample points to generate inside the sphere. 1796 on_surface (`bool`, optional): If set to True, points will be uniformly sampled on the surface of the hypersphere. 1797 1798 Returns: 1799 samples(`list`): Coordinates of the sample points inside the hypersphere. 1800 1801 Note: 1802 - Algorithm from: https://baezortega.github.io/2018/10/14/hypersphere-sampling/ 1803 """ 1804 # initial values 1805 center=np.array(center) 1806 d = center.shape[0] 1807 # sample n_samples points in d dimensions from a standard normal distribution 1808 samples = self.rng.normal(size=(n_samples, d)) 1809 # make the samples lie on the surface of the unit hypersphere 1810 normalize_radii = np.linalg.norm(samples, axis=1)[:, np.newaxis] 1811 samples /= normalize_radii 1812 if not on_surface: 1813 # make the samples lie inside the hypersphere with the correct density 1814 uniform_points = self.rng.uniform(size=n_samples)[:, np.newaxis] 1815 new_radii = np.power(uniform_points, 1/d) 1816 samples *= new_radii 1817 # scale the points to have the correct radius and center 1818 samples = samples * radius + center 1819 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.
Note:
- Algorithm from: https://baezortega.github.io/2018/10/14/hypersphere-sampling/
1821 def generate_trial_perpendicular_vector(self,vector,magnitude): 1822 """ 1823 Generates an orthogonal vector to the input `vector`. 1824 1825 Args: 1826 vector(`lst`): arbitrary vector. 1827 magnitude(`float`): magnitude of the orthogonal vector. 1828 1829 Returns: 1830 (`lst`): Orthogonal vector with the same magnitude as the input vector. 1831 """ 1832 np_vec = np.array(vector) 1833 np_vec /= np.linalg.norm(np_vec) 1834 if np.all(np_vec == 0): 1835 raise ValueError('Zero vector') 1836 # Generate a random vector 1837 random_vector = self.generate_random_points_in_a_sphere(radius=1, 1838 center=[0,0,0], 1839 n_samples=1, 1840 on_surface=True)[0] 1841 # Project the random vector onto the input vector and subtract the projection 1842 projection = np.dot(random_vector, np_vec) * np_vec 1843 perpendicular_vector = random_vector - projection 1844 # Normalize the perpendicular vector to have the same magnitude as the input vector 1845 perpendicular_vector /= np.linalg.norm(perpendicular_vector) 1846 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.
1848 def get_bond_length(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 1849 """ 1850 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns the initial bond length. 1851 If `use_default_bond` is activated and a "default" bond is defined, returns the length of that default bond instead. 1852 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. 1853 1854 Args: 1855 particle_name1 (str): label of the type of the first particle type of the bonded particles. 1856 particle_name2 (str): label of the type of the second particle type of the bonded particles. 1857 hard_check (bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 1858 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. 1859 1860 Returns: 1861 l0 (`float`): bond length 1862 1863 Note: 1864 - 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`. 1865 - If `hard_check`=`True` stops the code when no bond is found. 1866 """ 1867 bond_key = self.find_bond_key(particle_name1=particle_name1, 1868 particle_name2=particle_name2, 1869 use_default_bond=use_default_bond) 1870 if bond_key: 1871 return self.df[self.df['name']==bond_key].l0.values[0] 1872 else: 1873 print("Bond not defined between particles ", particle_name1, " and ", particle_name2) 1874 if hard_check: 1875 sys.exit(1) 1876 else: 1877 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 (
float
): 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.
1879 def get_charge_map(self): 1880 ''' 1881 Gets the charge of each `espresso_type` in `pymbe.df`. 1882 1883 Returns: 1884 charge_map(`dict`): {espresso_type: charge}. 1885 ''' 1886 if self.df.state_one['es_type'].isnull().values.any(): 1887 df_state_one = self.df.state_one.dropna() 1888 df_state_two = self.df.state_two.dropna() 1889 else: 1890 df_state_one = self.df.state_one 1891 if self.df.state_two['es_type'].isnull().values.any(): 1892 df_state_two = self.df.state_two.dropna() 1893 else: 1894 df_state_two = self.df.state_two 1895 state_one = pd.Series (df_state_one.charge.values,index=df_state_one.es_type.values) 1896 state_two = pd.Series (df_state_two.charge.values,index=df_state_two.es_type.values) 1897 charge_map = pd.concat([state_one,state_two],axis=0).to_dict() 1898 return charge_map
Gets the charge of each espresso_type
in pymbe.df
.
Returns:
charge_map(
dict
): {espresso_type: charge}.
1900 def get_lj_parameters(self, particle_name1, particle_name2, combining_rule='Lorentz-Berthelot'): 1901 """ 1902 Returns the Lennard-Jones parameters for the interaction between the particle types given by 1903 `particle_name1` and `particle_name2` in `pymbe.df`, calculated according to the provided combining rule. 1904 1905 Args: 1906 particle_name1 (str): label of the type of the first particle type 1907 particle_name2 (str): label of the type of the second particle type 1908 combining_rule (`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential betwen a pair of particles. Defaults to 'Lorentz-Berthelot'. 1909 1910 Returns: 1911 {"epsilon": LJ epsilon, "sigma": LJ sigma} 1912 1913 Note: 1914 Currently, the only `combining_rule` supported is Lorentz-Berthelot. 1915 """ 1916 supported_combining_rules=["Lorentz-Berthelot"] 1917 1918 if combining_rule not in supported_combining_rules: 1919 raise ValueError(f"Combining_rule {combining_rule} currently not implemented in pyMBE, valid keys are {supported_combining_rules}") 1920 1921 eps1 = self.df.loc[self.df["name"]==particle_name1].epsilon.iloc[0] 1922 sigma1 = self.df.loc[self.df["name"]==particle_name1].sigma.iloc[0] 1923 eps2 = self.df.loc[self.df["name"]==particle_name2].epsilon.iloc[0] 1924 sigma2 = self.df.loc[self.df["name"]==particle_name2].sigma.iloc[0] 1925 1926 return {"epsilon": np.sqrt(eps1*eps2), "sigma": (sigma1+sigma2)/2.0}
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": LJ epsilon, "sigma": LJ sigma}
Note:
Currently, the only
combining_rule
supported is Lorentz-Berthelot.
1928 def get_particle_id_map(self, object_name): 1929 ''' 1930 Gets all the ids associated with the object with name `object_name` in `pmb.df` 1931 1932 Args: 1933 object_name(`str`): name of the object 1934 1935 Returns: 1936 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]}, } 1937 ''' 1938 object_type = self.df.loc[self.df['name']== object_name].pmb_type.values[0] 1939 valid_types = ["particle", "molecule", "residue", "protein"] 1940 if object_type not in valid_types: 1941 raise ValueError(f"{object_name} is of pmb_type {object_type}, which is not supported by this function. Supported types are {valid_types}") 1942 id_list = [] 1943 mol_map = {} 1944 res_map = {} 1945 def do_res_map(res_ids): 1946 for res_id in res_ids: 1947 res_list=self.df.loc[self.df['residue_id']== res_id].particle_id.dropna().tolist() 1948 res_map[res_id]=res_list 1949 return res_map 1950 if object_type in ['molecule', 'protein']: 1951 mol_ids = self.df.loc[self.df['name']== object_name].molecule_id.dropna().tolist() 1952 for mol_id in mol_ids: 1953 res_ids = set(self.df.loc[self.df['molecule_id']== mol_id].residue_id.dropna().tolist()) 1954 res_map=do_res_map(res_ids=res_ids) 1955 mol_list=self.df.loc[self.df['molecule_id']== mol_id].particle_id.dropna().tolist() 1956 id_list+=mol_list 1957 mol_map[mol_id]=mol_list 1958 elif object_type == 'residue': 1959 res_ids = self.df.loc[self.df['name']== object_name].residue_id.dropna().tolist() 1960 res_map=do_res_map(res_ids=res_ids) 1961 id_list=[] 1962 for res_id_list in res_map.values(): 1963 id_list+=res_id_list 1964 elif object_type == 'particle': 1965 id_list = self.df.loc[self.df['name']== object_name].particle_id.dropna().tolist() 1966 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]}, }
1968 def get_pka_set(self): 1969 ''' 1970 Gets the pka-values and acidities of the particles with acid/base properties in `pmb.df` 1971 1972 Returns: 1973 pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}} 1974 ''' 1975 titratables_AA_df = self.df[[('name',''),('pka',''),('acidity','')]].drop_duplicates().dropna() 1976 pka_set = {} 1977 for index in titratables_AA_df.name.keys(): 1978 name = titratables_AA_df.name[index] 1979 pka_value = titratables_AA_df.pka[index] 1980 acidity = titratables_AA_df.acidity[index] 1981 pka_set[name] = {'pka_value':pka_value,'acidity':acidity} 1982 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}}
1984 def get_radius_map(self): 1985 ''' 1986 Gets the effective radius of each `espresso type` in `pmb.df`. 1987 1988 Returns: 1989 radius_map(`dict`): {espresso_type: radius}. 1990 ''' 1991 1992 df_state_one = self.df[[('sigma',''),('state_one','es_type')]].dropna().drop_duplicates() 1993 df_state_two = self.df[[('sigma',''),('state_two','es_type')]].dropna().drop_duplicates() 1994 1995 state_one = pd.Series(df_state_one.sigma.values/2.0,index=df_state_one.state_one.es_type.values) 1996 state_two = pd.Series(df_state_two.sigma.values/2.0,index=df_state_two.state_two.es_type.values) 1997 radius_map = pd.concat([state_one,state_two],axis=0).to_dict() 1998 1999 return radius_map
Gets the effective radius of each espresso type
in pmb.df
.
Returns:
radius_map(
dict
): {espresso_type: radius}.
2001 def get_resource(self, path): 2002 ''' 2003 Locate a file resource of the pyMBE package. 2004 2005 Args: 2006 path(`str`): Relative path to the resource 2007 2008 Returns: 2009 path(`str`): Absolute path to the resource 2010 2011 ''' 2012 import os 2013 return os.path.join(os.path.dirname(__file__), path)
Locate a file resource of the pyMBE package.
Arguments:
- path(
str
): Relative path to the resource
Returns:
path(
str
): Absolute path to the resource
2015 def get_type_map(self): 2016 """ 2017 Gets all different espresso types assigned to particles in `pmb.df`. 2018 2019 Returns: 2020 type_map(`dict`): {"name": espresso_type}. 2021 """ 2022 if self.df.state_one['es_type'].isnull().values.any(): 2023 df_state_one = self.df.state_one.dropna(how='all') 2024 df_state_two = self.df.state_two.dropna(how='all') 2025 else: 2026 df_state_one = self.df.state_one 2027 if self.df.state_two['es_type'].isnull().values.any(): 2028 df_state_two = self.df.state_two.dropna(how='all') 2029 else: 2030 df_state_two = self.df.state_two 2031 state_one = pd.Series (df_state_one.es_type.values,index=df_state_one.label) 2032 state_two = pd.Series (df_state_two.es_type.values,index=df_state_two.label) 2033 type_map = pd.concat([state_one,state_two],axis=0).to_dict() 2034 return type_map
Gets all different espresso types assigned to particles in pmb.df
.
Returns:
type_map(
dict
): {"name": espresso_type}.
2036 def load_interaction_parameters(self, filename, verbose=True): 2037 """ 2038 Loads the interaction parameters stored in `filename` into `pmb.df` 2039 2040 Args: 2041 filename(`str`): name of the file to be read 2042 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2043 """ 2044 without_units = ['q','es_type','acidity'] 2045 with_units = ['sigma','epsilon'] 2046 2047 with open(filename, 'r') as f: 2048 interaction_data = json.load(f) 2049 interaction_parameter_set = interaction_data["data"] 2050 2051 for key in interaction_parameter_set: 2052 param_dict=interaction_parameter_set[key] 2053 object_type=param_dict['object_type'] 2054 if object_type == 'particle': 2055 not_requiered_attributes={} 2056 for not_requiered_key in without_units+with_units: 2057 if not_requiered_key in param_dict.keys(): 2058 if not_requiered_key in with_units: 2059 not_requiered_attributes[not_requiered_key]=self.create_variable_with_units(variable=param_dict.pop(not_requiered_key)) 2060 elif not_requiered_key in without_units: 2061 not_requiered_attributes[not_requiered_key]=param_dict.pop(not_requiered_key) 2062 else: 2063 if not_requiered_key == 'acidity': 2064 not_requiered_attributes[not_requiered_key] = 'inert' 2065 else: 2066 not_requiered_attributes[not_requiered_key]=None 2067 self.define_particle(name=param_dict.pop('name'), 2068 q=not_requiered_attributes.pop('q'), 2069 sigma=not_requiered_attributes.pop('sigma'), 2070 acidity=not_requiered_attributes.pop('acidity'), 2071 epsilon=not_requiered_attributes.pop('epsilon'), 2072 verbose=verbose) 2073 elif object_type == 'residue': 2074 self.define_residue (name = param_dict.pop('name'), 2075 central_bead = param_dict.pop('central_bead_name'), 2076 side_chains = param_dict.pop('side_chains_names')) 2077 elif object_type == 'molecule': 2078 self.define_molecule(name=param_dict.pop('name'), 2079 residue_list=param_dict.pop('residue_name_list')) 2080 elif object_type == 'peptide': 2081 self.define_peptide(name=param_dict.pop('name'), 2082 sequence=param_dict.pop('sequence'), 2083 model=param_dict.pop('model')) 2084 elif object_type == 'bond': 2085 particle_pairs = param_dict.pop('particle_pairs') 2086 bond_parameters = param_dict.pop('bond_parameters') 2087 bond_type = param_dict.pop('bond_type') 2088 if bond_type == 'harmonic': 2089 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2090 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2091 bond = {'r_0' : r_0, 2092 'k' : k, 2093 } 2094 2095 elif bond_type == 'FENE': 2096 k = self.create_variable_with_units(variable=bond_parameters.pop('k')) 2097 r_0 = self.create_variable_with_units(variable=bond_parameters.pop('r_0')) 2098 d_r_max = self.create_variable_with_units(variable=bond_parameters.pop('d_r_max')) 2099 bond = {'r_0' : r_0, 2100 'k' : k, 2101 'd_r_max': d_r_max, 2102 } 2103 else: 2104 raise ValueError("Current implementation of pyMBE only supports harmonic and FENE bonds") 2105 self.define_bond(bond_type=bond_type, bond_parameters=bond, particle_pairs=particle_pairs) 2106 else: 2107 raise ValueError(object_type+' is not a known pmb object type') 2108 2109 return
Loads the interaction parameters stored in filename
into pmb.df
Arguments:
- filename(
str
): name of the file to be read - verbose (
bool
, optional): Switch to activate/deactivate verbose. Defaults to True.
2111 def load_pka_set(self, filename, verbose=True): 2112 """ 2113 Loads the pka_set stored in `filename` into `pmb.df`. 2114 2115 Args: 2116 filename(`str`): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}. 2117 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2118 2119 Note: 2120 - If `name` is already defined in the `pymbe.df`, it prints a warning. 2121 """ 2122 with open(filename, 'r') as f: 2123 pka_data = json.load(f) 2124 pka_set = pka_data["data"] 2125 2126 self.check_pka_set(pka_set=pka_set) 2127 2128 for key in pka_set: 2129 acidity = pka_set[key]['acidity'] 2130 pka_value = pka_set[key]['pka_value'] 2131 self.define_particle( 2132 name=key, 2133 acidity=acidity, 2134 pka=pka_value, 2135 verbose=verbose) 2136 return
Loads the pka_set stored in filename
into pmb.df
.
Arguments:
- filename(
str
): name of the file with the pka set to be loaded. Expected format is {name:{"acidity": acidity, "pka_value":pka_value}}. - verbose (
bool
, optional): Switch to activate/deactivate verbose. Defaults to True.
Note:
- If
name
is already defined in thepymbe.df
, it prints a warning.
2138 def parse_sequence_from_file(self,sequence): 2139 """ 2140 Parses the given sequence such that it can be used in pyMBE. This function has to be used if the df was read from a file. 2141 2142 Args: 2143 sequence(`str`): sequence to be parsed 2144 2145 Returns: 2146 sequence(`lst`): parsed sequence 2147 """ 2148 sequence = sequence.replace(' ', '') 2149 sequence = sequence.replace("'", '') 2150 sequence = sequence.split(",")[1:-1] 2151 return sequence
Parses the given sequence such that it can be used in pyMBE. This function has to be used if the df was read from a file.
Arguments:
- sequence(
str
): sequence to be parsed
Returns:
sequence(
lst
): parsed sequence
2153 def print_reduced_units(self): 2154 """ 2155 Prints the current set of reduced units defined in pyMBE.units. 2156 """ 2157 print("\nCurrent set of reduced units:") 2158 unit_length=self.units.Quantity(1,'reduced_length') 2159 unit_energy=self.units.Quantity(1,'reduced_energy') 2160 unit_charge=self.units.Quantity(1,'reduced_charge') 2161 print(unit_length.to('nm'), "=", unit_length) 2162 print(unit_energy.to('J'), "=", unit_energy) 2163 print('Temperature:', (self.kT/self.Kb).to("K")) 2164 print(unit_charge.to('C'), "=", unit_charge) 2165 print()
Prints the current set of reduced units defined in pyMBE.units.
2167 def propose_unused_type(self): 2168 """ 2169 Searches in `pmb.df` all the different particle types defined and returns a new one. 2170 2171 Returns: 2172 unused_type(`int`): unused particle type 2173 """ 2174 type_map=self.get_type_map() 2175 if type_map == {}: 2176 unused_type = 0 2177 else: 2178 unused_type=max(type_map.values())+1 2179 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
2181 def protein_sequence_parser(self, sequence): 2182 ''' 2183 Parses `sequence` to the one letter code for amino acids. 2184 2185 Args: 2186 sequence(`str` or `lst`): Sequence of the amino acid. 2187 2188 Returns: 2189 clean_sequence(`lst`): `sequence` using the one letter code. 2190 2191 Note: 2192 - Accepted formats for `sequence` are: 2193 - `lst` with one letter or three letter code of each aminoacid in each element 2194 - `str` with the sequence using the one letter code 2195 - `str` with the squence using the three letter code, each aminoacid must be separated by a hyphen "-" 2196 2197 ''' 2198 # Aminoacid key 2199 keys={"ALA": "A", 2200 "ARG": "R", 2201 "ASN": "N", 2202 "ASP": "D", 2203 "CYS": "C", 2204 "GLU": "E", 2205 "GLN": "Q", 2206 "GLY": "G", 2207 "HIS": "H", 2208 "ILE": "I", 2209 "LEU": "L", 2210 "LYS": "K", 2211 "MET": "M", 2212 "PHE": "F", 2213 "PRO": "P", 2214 "SER": "S", 2215 "THR": "T", 2216 "TRP": "W", 2217 "TYR": "Y", 2218 "VAL": "V", 2219 "PSER": "J", 2220 "PTHR": "U", 2221 "PTyr": "Z", 2222 "NH2": "n", 2223 "COOH": "c"} 2224 clean_sequence=[] 2225 if isinstance(sequence, str): 2226 if sequence.find("-") != -1: 2227 splited_sequence=sequence.split("-") 2228 for residue in splited_sequence: 2229 if len(residue) == 1: 2230 if residue in keys.values(): 2231 residue_ok=residue 2232 else: 2233 if residue.upper() in keys.values(): 2234 residue_ok=residue.upper() 2235 else: 2236 raise ValueError("Unknown one letter code for a residue given: ", residue, " please review the input sequence") 2237 clean_sequence.append(residue_ok) 2238 else: 2239 if residue in keys.keys(): 2240 clean_sequence.append(keys[residue]) 2241 else: 2242 if residue.upper() in keys.keys(): 2243 clean_sequence.append(keys[residue.upper()]) 2244 else: 2245 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2246 else: 2247 for residue in sequence: 2248 if residue in keys.values(): 2249 residue_ok=residue 2250 else: 2251 if residue.upper() in keys.values(): 2252 residue_ok=residue.upper() 2253 else: 2254 raise ValueError("Unknown one letter code for a residue: ", residue, " please review the input sequence") 2255 clean_sequence.append(residue_ok) 2256 if isinstance(sequence, list): 2257 for residue in sequence: 2258 if residue in keys.values(): 2259 residue_ok=residue 2260 else: 2261 if residue.upper() in keys.values(): 2262 residue_ok=residue.upper() 2263 elif (residue.upper() in keys.keys()): 2264 clean_sequence.append(keys[residue.upper()]) 2265 else: 2266 raise ValueError("Unknown code for a residue: ", residue, " please review the input sequence") 2267 clean_sequence.append(residue_ok) 2268 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 "-"
2270 def read_pmb_df (self,filename): 2271 """ 2272 Reads a pyMBE's Dataframe stored in `filename` and stores the information into pyMBE. 2273 2274 Args: 2275 filename(str): path to file. 2276 2277 Note: 2278 This function only accepts files with CSV format. 2279 """ 2280 2281 if filename.rsplit(".", 1)[1] != "csv": 2282 raise ValueError("Only files with CSV format are supported") 2283 df = pd.read_csv (filename,header=[0, 1], index_col=0) 2284 columns_names = self.setup_df() 2285 2286 multi_index = pd.MultiIndex.from_tuples(columns_names) 2287 df.columns = multi_index 2288 2289 self.df = self.convert_columns_to_original_format(df) 2290 2291 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.
2293 def read_protein_vtf_in_df (self,filename,unit_length=None): 2294 """ 2295 Loads a coarse-grained protein model in a vtf file `filename` into the `pmb.df` and it labels it with `name`. 2296 2297 Args: 2298 filename(`str`): Path to the vtf file with the coarse-grained model. 2299 unit_length(`obj`): unit of length of the the coordinates in `filename` using the pyMBE UnitRegistry. Defaults to None. 2300 2301 Returns: 2302 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 2303 2304 Note: 2305 - If no `unit_length` is provided, it is assumed that the coordinates are in Angstrom. 2306 """ 2307 2308 print (f'Loading protein coarse grain model file: {filename}') 2309 2310 coord_list = [] 2311 particles_dict = {} 2312 2313 if unit_length is None: 2314 unit_length = 1 * self.units.angstrom 2315 2316 with open (filename,'r') as protein_model: 2317 for line in protein_model : 2318 line_split = line.split() 2319 if line_split: 2320 line_header = line_split[0] 2321 if line_header == 'atom': 2322 atom_id = line_split[1] 2323 atom_name = line_split[3] 2324 atom_resname = line_split[5] 2325 chain_id = line_split[9] 2326 radius = float(line_split [11])*unit_length 2327 sigma = 2*radius 2328 particles_dict [int(atom_id)] = [atom_name , atom_resname, chain_id, sigma] 2329 elif line_header.isnumeric(): 2330 atom_coord = line_split[1:] 2331 atom_coord = [(float(i)*unit_length).to('reduced_length').magnitude for i in atom_coord] 2332 coord_list.append (atom_coord) 2333 2334 numbered_label = [] 2335 i = 0 2336 2337 for atom_id in particles_dict.keys(): 2338 2339 if atom_id == 1: 2340 atom_name = particles_dict[atom_id][0] 2341 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2342 numbered_label.append(numbered_name) 2343 2344 elif atom_id != 1: 2345 2346 if particles_dict[atom_id-1][1] != particles_dict[atom_id][1]: 2347 i += 1 2348 count = 1 2349 atom_name = particles_dict[atom_id][0] 2350 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2351 numbered_label.append(numbered_name) 2352 2353 elif particles_dict[atom_id-1][1] == particles_dict[atom_id][1]: 2354 if count == 2 or particles_dict[atom_id][1] == 'GLY': 2355 i +=1 2356 count = 0 2357 atom_name = particles_dict[atom_id][0] 2358 numbered_name = [f'{atom_name}{i}',particles_dict[atom_id][2],particles_dict[atom_id][3]] 2359 numbered_label.append(numbered_name) 2360 count +=1 2361 2362 topology_dict = {} 2363 2364 for i in range (0, len(numbered_label)): 2365 topology_dict [numbered_label[i][0]] = {'initial_pos': coord_list[i] , 2366 'chain_id':numbered_label[i][1], 2367 'sigma':numbered_label[i][2] } 2368 2369 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.
2371 def search_bond(self, particle_name1, particle_name2, hard_check=False, use_default_bond=False) : 2372 """ 2373 Searches for bonds between the particle types given by `particle_name1` and `particle_name2` in `pymbe.df` and returns it. 2374 If `use_default_bond` is activated and a "default" bond is defined, returns that default bond instead. 2375 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. 2376 2377 Args: 2378 particle_name1 (str): label of the type of the first particle type of the bonded particles. 2379 particle_name2 (str): label of the type of the second particle type of the bonded particles. 2380 hard_check (bool, optional): If it is activated, the code stops if no bond is found. Defaults to False. 2381 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. 2382 2383 Returns: 2384 bond (cls): bond object from the espressomd library. 2385 2386 Note: 2387 - 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`. 2388 - If `hard_check`=`True` stops the code when no bond is found. 2389 """ 2390 2391 bond_key = self.find_bond_key(particle_name1=particle_name1, 2392 particle_name2=particle_name2, 2393 use_default_bond=use_default_bond) 2394 if use_default_bond: 2395 if not self.check_if_name_is_defined_in_df(name="default",pmb_type_to_be_defined='bond'): 2396 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") 2397 if bond_key: 2398 return self.df[self.df['name']==bond_key].bond_object.values[0] 2399 else: 2400 print("Bond not defined between particles ", particle_name1, " and ", particle_name2) 2401 if hard_check: 2402 sys.exit(1) 2403 else: 2404 return
Searches for bonds between the particle types given by particle_name1
and particle_name2
in pymbe.df
and returns it.
If use_default_bond
is activated and a "default" bond is defined, returns that default bond instead.
If no bond is found, it prints a message and it does not return anything. If hard_check
is activated, the code stops if no bond is found.
Arguments:
- particle_name1 (str): label of the type of the first particle type of the bonded particles.
- particle_name2 (str): label of the type of the second particle type of the bonded particles.
- hard_check (bool, optional): If it is activated, the code stops if no bond is found. Defaults to False.
- use_default_bond (bool, optional): If it is activated, the "default" bond is returned if no bond is found between
particle_name1
andparticle_name2
. Defaults to False.
Returns:
bond (cls): 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.
2406 def search_particles_in_residue(self, residue_name): 2407 ''' 2408 Searches for all particles in a given residue of name `residue_name`. 2409 2410 Args: 2411 residue_name (`str`): name of the residue to be searched 2412 2413 Returns: 2414 list_of_particles_in_residue (`lst`): list of the names of all particles in the residue 2415 2416 Note: 2417 - 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. 2418 2419 ''' 2420 index_residue = self.df.loc[self.df['name'] == residue_name].index[0].item() 2421 central_bead = self.df.at [index_residue, ('central_bead', '')] 2422 list_of_side_chains = self.df.at [index_residue, ('side_chains', '')] 2423 2424 list_of_particles_in_residue = [] 2425 list_of_particles_in_residue.append(central_bead) 2426 2427 for side_chain in list_of_side_chains: 2428 object_type = self.df[self.df['name']==side_chain].pmb_type.values[0] 2429 2430 if object_type == "residue": 2431 list_of_particles_in_side_chain_residue = self.search_particles_in_residue(side_chain) 2432 list_of_particles_in_residue += list_of_particles_in_side_chain_residue 2433 elif object_type == "particle": 2434 list_of_particles_in_residue.append(side_chain) 2435 2436 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.
2438 def set_particle_acidity(self, name, acidity='inert', default_charge=0, pka=None, verbose=True): 2439 """ 2440 Sets the particle acidity if it is acidic or basic, creates `state_one` and `state_two` with the protonated and 2441 deprotonated states. In each state is set: `label`,`charge` and `es_type`. If it is inert, it will define only `state_one`. 2442 2443 Args: 2444 name (`str`): Unique label that identifies the `particle`. 2445 acidity (`str`): Identifies whether the particle is `acidic` or `basic`, used to setup constant pH simulations. Defaults to `inert`. 2446 default_charge (`int`): Charge of the particle. Defaults to 0. 2447 pka (`float`, optional): If `particle` is an acid or a base, it defines its pka-value. Defaults to None. 2448 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2449 """ 2450 acidity_valid_keys = ['inert','acidic', 'basic'] 2451 if acidity not in acidity_valid_keys: 2452 raise ValueError(f"Acidity {acidity} provided for particle name {name} is not supproted. Valid keys are: {acidity_valid_keys}") 2453 if acidity in ['acidic', 'basic'] and pka is None: 2454 raise ValueError(f"pKa not provided for particle with name {name} with acidity {acidity}. pKa must be provided for acidic or basic particles.") 2455 for index in self.df[self.df['name']==name].index: 2456 if pka: 2457 self.add_value_to_df(key=('pka',''), 2458 index=index, 2459 new_value=pka, 2460 verbose=verbose) 2461 self.add_value_to_df(key=('acidity',''), 2462 index=index, 2463 new_value=acidity, 2464 verbose=verbose) 2465 if not self.check_if_df_cell_has_a_value(index=index,key=('state_one','es_type')): 2466 self.add_value_to_df(key=('state_one','es_type'), 2467 index=index, 2468 new_value=self.propose_unused_type(), 2469 verbose=verbose) 2470 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'inert': 2471 self.add_value_to_df(key=('state_one','charge'), 2472 index=index, 2473 new_value=default_charge, 2474 verbose=verbose) 2475 self.add_value_to_df(key=('state_one','label'), 2476 index=index, 2477 new_value=name, 2478 verbose=verbose) 2479 else: 2480 protonated_label = f'{name}H' 2481 self.add_value_to_df(key=('state_one','label'), 2482 index=index, 2483 new_value=protonated_label, 2484 verbose=verbose) 2485 self.add_value_to_df(key=('state_two','label'), 2486 index=index, 2487 new_value=name, 2488 verbose=verbose) 2489 if not self.check_if_df_cell_has_a_value(index=index,key=('state_two','es_type')): 2490 self.add_value_to_df(key=('state_two','es_type'), 2491 index=index, 2492 new_value=self.propose_unused_type(), 2493 verbose=verbose) 2494 if self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'acidic': 2495 self.add_value_to_df(key=('state_one','charge'), 2496 index=index,new_value=0, 2497 verbose=verbose) 2498 self.add_value_to_df(key=('state_two','charge'), 2499 index=index, 2500 new_value=-1, 2501 verbose=verbose) 2502 elif self.df.loc [self.df['name'] == name].acidity.iloc[0] == 'basic': 2503 self.add_value_to_df(key=('state_one','charge'), 2504 index=index,new_value=+1, 2505 verbose=verbose) 2506 self.add_value_to_df(key=('state_two','charge'), 2507 index=index, 2508 new_value=0, 2509 verbose=verbose) 2510 return
Sets the particle acidity if it is acidic or basic, creates state_one
and state_two
with the protonated and
deprotonated states. In each state is set: label
,charge
and es_type
. If it is inert, it will define only state_one
.
Arguments:
- name (
str
): Unique label that identifies theparticle
. - acidity (
str
): Identifies whether the particle isacidic
orbasic
, used to setup constant pH simulations. Defaults toinert
. - default_charge (
int
): Charge of the particle. Defaults to 0. - pka (
float
, optional): Ifparticle
is an acid or a base, it defines its pka-value. Defaults to None. - verbose (
bool
, optional): Switch to activate/deactivate verbose. Defaults to True.
2512 def set_reduced_units(self, unit_length=None, unit_charge=None, temperature=None, Kw=None, verbose=True): 2513 """ 2514 Sets the set of reduced units used by pyMBE.units and it prints it. 2515 2516 Args: 2517 unit_length (`obj`,optional): Reduced unit of length defined using the `pmb.units` UnitRegistry. Defaults to None. 2518 unit_charge (`obj`,optional): Reduced unit of charge defined using the `pmb.units` UnitRegistry. Defaults to None. 2519 temperature (`obj`,optional): Temperature of the system, defined using the `pmb.units` UnitRegistry. Defaults to None. 2520 Kw (`obj`,optional): Ionic product of water in mol^2/l^2. Defaults to None. 2521 verbose (`bool`, optional): Switch to activate/deactivate verbose. Defaults to True. 2522 2523 Note: 2524 - If no `temperature` is given, a value of 298.15 K is assumed by default. 2525 - If no `unit_length` is given, a value of 0.355 nm is assumed by default. 2526 - If no `unit_charge` is given, a value of 1 elementary charge is assumed by default. 2527 - If no `Kw` is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default. 2528 """ 2529 self.units=pint.UnitRegistry() 2530 if unit_length is None: 2531 unit_length=0.355*self.units.nm 2532 if temperature is None: 2533 temperature=298.15 * self.units.K 2534 if unit_charge is None: 2535 unit_charge=self.units.e 2536 if Kw is None: 2537 Kw = 1e-14 2538 self.N_A=6.02214076e23 / self.units.mol 2539 self.Kb=1.38064852e-23 * self.units.J / self.units.K 2540 self.e=1.60217662e-19 *self.units.C 2541 self.kT=temperature*self.Kb 2542 self.Kw=Kw*self.units.mol**2 / (self.units.l**2) 2543 self.units.define(f'reduced_energy = {self.kT} ') 2544 self.units.define(f'reduced_length = {unit_length}') 2545 self.units.define(f'reduced_charge = {unit_charge}') 2546 if verbose: 2547 self.print_reduced_units() 2548 return
Sets the set of reduced units used by pyMBE.units and it prints it.
Arguments:
- unit_length (
obj
,optional): Reduced unit of length defined using thepmb.units
UnitRegistry. Defaults to None. - unit_charge (
obj
,optional): Reduced unit of charge defined using thepmb.units
UnitRegistry. Defaults to None. - temperature (
obj
,optional): Temperature of the system, defined using thepmb.units
UnitRegistry. Defaults to None. - Kw (
obj
,optional): Ionic product of water in mol^2/l^2. Defaults to None. - verbose (
bool
, optional): Switch to activate/deactivate verbose. Defaults to True.
Note:
- If no
temperature
is given, a value of 298.15 K is assumed by default.- If no
unit_length
is given, a value of 0.355 nm is assumed by default.- If no
unit_charge
is given, a value of 1 elementary charge is assumed by default.- If no
Kw
is given, a value of 10^(-14) * mol^2 / l^2 is assumed by default.
2550 def setup_cpH (self, counter_ion, constant_pH, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2551 """ 2552 Sets up the Acid/Base reactions for acidic/basic `particles` defined in `pmb.df` to be sampled in the constant pH ensemble. 2553 2554 Args: 2555 counter_ion(`str`): `name` of the counter_ion `particle`. 2556 constant_pH(`float`): pH-value. 2557 exclusion_range(`float`, optional): Below this value, no particles will be inserted. 2558 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2559 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2560 2561 Returns: 2562 RE (`obj`): Instance of a reaction_ensemble.ConstantpHEnsemble object from the espressomd library. 2563 sucessfull_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2564 """ 2565 from espressomd import reaction_methods 2566 if exclusion_range is None: 2567 exclusion_range = max(self.get_radius_map().values())*2.0 2568 if pka_set is None: 2569 pka_set=self.get_pka_set() 2570 self.check_pka_set(pka_set=pka_set) 2571 if use_exclusion_radius_per_type: 2572 exclusion_radius_per_type = self.get_radius_map() 2573 else: 2574 exclusion_radius_per_type = {} 2575 2576 RE = reaction_methods.ConstantpHEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2577 exclusion_range=exclusion_range.magnitude, 2578 seed=self.SEED, 2579 constant_pH=constant_pH, 2580 exclusion_radius_per_type = exclusion_radius_per_type 2581 ) 2582 sucessfull_reactions_labels=[] 2583 charge_map = self.get_charge_map() 2584 for name in pka_set.keys(): 2585 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2586 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in sg.type_map') 2587 continue 2588 gamma=10**-pka_set[name]['pka_value'] 2589 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2590 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2591 counter_ion_type = self.df.loc[self.df['name']==counter_ion].state_one.es_type.values[0] 2592 RE.add_reaction(gamma=gamma, 2593 reactant_types=[state_one_type], 2594 product_types=[state_two_type, counter_ion_type], 2595 default_charges={state_one_type: charge_map[state_one_type], 2596 state_two_type: charge_map[state_two_type], 2597 counter_ion_type: charge_map[counter_ion_type]}) 2598 sucessfull_reactions_labels.append(name) 2599 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(
float
, 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 (
obj
): Instance of a reaction_ensemble.ConstantpHEnsemble object from the espressomd library. sucessfull_reactions_labels(lst
): Labels of the reactions set up by pyMBE.
2601 def setup_gcmc(self, c_salt_res, salt_cation_name, salt_anion_name, activity_coefficient=None, exclusion_range=None, use_exclusion_radius_per_type = False): 2602 """ 2603 Sets up grand-canonical coupling to a reservoir of salt. 2604 For reactive systems coupled to a reservoir, the grand-reaction method has to be used instead. 2605 2606 Args: 2607 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2608 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 2609 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 2610 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2611 exclusion_range(`float`, optional): For distances shorter than this value, no particles will be inserted. 2612 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2613 2614 Returns: 2615 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2616 """ 2617 from espressomd import reaction_methods 2618 if exclusion_range is None: 2619 exclusion_range = max(self.get_radius_map().values())*2.0 2620 if use_exclusion_radius_per_type: 2621 exclusion_radius_per_type = self.get_radius_map() 2622 else: 2623 exclusion_radius_per_type = {} 2624 2625 #If no function for the activity coefficient is provided, the ideal case is assumed. 2626 if activity_coefficient is None: 2627 activity_coefficient = lambda x: 1.0 2628 2629 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2630 exclusion_range=exclusion_range.magnitude, 2631 seed=self.SEED, 2632 exclusion_radius_per_type = exclusion_radius_per_type 2633 ) 2634 2635 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2636 determined_activity_coefficient = activity_coefficient(c_salt_res) 2637 K_salt = (c_salt_res.to('1/(N_A * reduced_length**3)')**2) * determined_activity_coefficient 2638 2639 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 2640 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 2641 2642 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.charge.values[0] 2643 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.charge.values[0] 2644 2645 if salt_cation_charge <= 0: 2646 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 2647 if salt_anion_charge >= 0: 2648 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 2649 2650 # Grand-canonical coupling to the reservoir 2651 RE.add_reaction( 2652 gamma = K_salt.magnitude, 2653 reactant_types = [], 2654 reactant_coefficients = [], 2655 product_types = [ salt_cation_es_type, salt_anion_es_type ], 2656 product_coefficients = [ 1, 1 ], 2657 default_charges = { 2658 salt_cation_es_type: salt_cation_charge, 2659 salt_anion_es_type: salt_anion_charge, 2660 } 2661 ) 2662 2663 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 ('float'): 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', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
- exclusion_range(
float
, 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 (
obj
): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library.
2665 def setup_grxmc_reactions(self, pH_res, c_salt_res, proton_name, hydroxide_name, salt_cation_name, salt_anion_name, activity_coefficient=None, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2666 """ 2667 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 2668 reservoir of small ions. 2669 This implementation uses the original formulation of the grand-reaction method by Landsgesell et al. [1]. 2670 2671 [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. 2672 2673 Args: 2674 pH_res ('float): pH-value in the reservoir. 2675 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2676 proton_name ('str'): Name of the proton (H+) particle. 2677 hydroxide_name ('str'): Name of the hydroxide (OH-) particle. 2678 salt_cation_name ('str'): Name of the salt cation (e.g. Na+) particle. 2679 salt_anion_name ('str'): Name of the salt anion (e.g. Cl-) particle. 2680 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2681 exclusion_range(`float`, optional): For distances shorter than this value, no particles will be inserted. 2682 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2683 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius for each espresso_type is used. Defaults to `False`. 2684 2685 Returns: 2686 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2687 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2688 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 2689 """ 2690 from espressomd import reaction_methods 2691 if exclusion_range is None: 2692 exclusion_range = max(self.get_radius_map().values())*2.0 2693 if pka_set is None: 2694 pka_set=self.get_pka_set() 2695 self.check_pka_set(pka_set=pka_set) 2696 if use_exclusion_radius_per_type: 2697 exclusion_radius_per_type = self.get_radius_map() 2698 else: 2699 exclusion_radius_per_type = {} 2700 2701 #If no function for the activity coefficient is provided, the ideal case is assumed. 2702 if activity_coefficient is None: 2703 activity_coefficient = lambda x: 1.0 2704 2705 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2706 exclusion_range=exclusion_range.magnitude, 2707 seed=self.SEED, 2708 exclusion_radius_per_type = exclusion_radius_per_type 2709 ) 2710 2711 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2712 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 2713 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2714 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 2715 K_W = cH_res.to('1/(N_A * reduced_length**3)') * cOH_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2716 K_NACL = cNa_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2717 K_HCL = cH_res.to('1/(N_A * reduced_length**3)') * cCl_res.to('1/(N_A * reduced_length**3)') * determined_activity_coefficient 2718 2719 proton_es_type = self.df.loc[self.df['name']==proton_name].state_one.es_type.values[0] 2720 hydroxide_es_type = self.df.loc[self.df['name']==hydroxide_name].state_one.es_type.values[0] 2721 salt_cation_es_type = self.df.loc[self.df['name']==salt_cation_name].state_one.es_type.values[0] 2722 salt_anion_es_type = self.df.loc[self.df['name']==salt_anion_name].state_one.es_type.values[0] 2723 2724 proton_charge = self.df.loc[self.df['name']==proton_name].state_one.charge.values[0] 2725 hydroxide_charge = self.df.loc[self.df['name']==hydroxide_name].state_one.charge.values[0] 2726 salt_cation_charge = self.df.loc[self.df['name']==salt_cation_name].state_one.charge.values[0] 2727 salt_anion_charge = self.df.loc[self.df['name']==salt_anion_name].state_one.charge.values[0] 2728 2729 if proton_charge <= 0: 2730 raise ValueError('ERROR proton charge must be positive, charge ', proton_charge) 2731 if salt_cation_charge <= 0: 2732 raise ValueError('ERROR salt cation charge must be positive, charge ', salt_cation_charge) 2733 if hydroxide_charge >= 0: 2734 raise ValueError('ERROR hydroxide charge must be negative, charge ', hydroxide_charge) 2735 if salt_anion_charge >= 0: 2736 raise ValueError('ERROR salt anion charge must be negative, charge ', salt_anion_charge) 2737 2738 # Grand-canonical coupling to the reservoir 2739 # 0 = H+ + OH- 2740 RE.add_reaction( 2741 gamma = K_W.magnitude, 2742 reactant_types = [], 2743 reactant_coefficients = [], 2744 product_types = [ proton_es_type, hydroxide_es_type ], 2745 product_coefficients = [ 1, 1 ], 2746 default_charges = { 2747 proton_es_type: proton_charge, 2748 hydroxide_es_type: hydroxide_charge, 2749 } 2750 ) 2751 2752 # 0 = Na+ + Cl- 2753 RE.add_reaction( 2754 gamma = K_NACL.magnitude, 2755 reactant_types = [], 2756 reactant_coefficients = [], 2757 product_types = [ salt_cation_es_type, salt_anion_es_type ], 2758 product_coefficients = [ 1, 1 ], 2759 default_charges = { 2760 salt_cation_es_type: salt_cation_charge, 2761 salt_anion_es_type: salt_anion_charge, 2762 } 2763 ) 2764 2765 # 0 = Na+ + OH- 2766 RE.add_reaction( 2767 gamma = (K_NACL * K_W / K_HCL).magnitude, 2768 reactant_types = [], 2769 reactant_coefficients = [], 2770 product_types = [ salt_cation_es_type, hydroxide_es_type ], 2771 product_coefficients = [ 1, 1 ], 2772 default_charges = { 2773 salt_cation_es_type: salt_cation_charge, 2774 hydroxide_es_type: hydroxide_charge, 2775 } 2776 ) 2777 2778 # 0 = H+ + Cl- 2779 RE.add_reaction( 2780 gamma = K_HCL.magnitude, 2781 reactant_types = [], 2782 reactant_coefficients = [], 2783 product_types = [ proton_es_type, salt_anion_es_type ], 2784 product_coefficients = [ 1, 1 ], 2785 default_charges = { 2786 proton_es_type: proton_charge, 2787 salt_anion_es_type: salt_anion_charge, 2788 } 2789 ) 2790 2791 # Annealing moves to ensure sufficient sampling 2792 # Cation annealing H+ = Na+ 2793 RE.add_reaction( 2794 gamma = (K_NACL / K_HCL).magnitude, 2795 reactant_types = [proton_es_type], 2796 reactant_coefficients = [ 1 ], 2797 product_types = [ salt_cation_es_type ], 2798 product_coefficients = [ 1 ], 2799 default_charges = { 2800 proton_es_type: proton_charge, 2801 salt_cation_es_type: salt_cation_charge, 2802 } 2803 ) 2804 2805 # Anion annealing OH- = Cl- 2806 RE.add_reaction( 2807 gamma = (K_HCL / K_W).magnitude, 2808 reactant_types = [hydroxide_es_type], 2809 reactant_coefficients = [ 1 ], 2810 product_types = [ salt_anion_es_type ], 2811 product_coefficients = [ 1 ], 2812 default_charges = { 2813 hydroxide_es_type: hydroxide_charge, 2814 salt_anion_es_type: salt_anion_charge, 2815 } 2816 ) 2817 2818 sucessful_reactions_labels=[] 2819 charge_map = self.get_charge_map() 2820 for name in pka_set.keys(): 2821 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2822 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in the type map.') 2823 continue 2824 2825 Ka = (10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 2826 2827 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2828 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2829 2830 # Reaction in terms of proton: HA = A + H+ 2831 RE.add_reaction(gamma=Ka.magnitude, 2832 reactant_types=[state_one_type], 2833 reactant_coefficients=[1], 2834 product_types=[state_two_type, proton_es_type], 2835 product_coefficients=[1, 1], 2836 default_charges={state_one_type: charge_map[state_one_type], 2837 state_two_type: charge_map[state_two_type], 2838 proton_es_type: proton_charge}) 2839 2840 # Reaction in terms of salt cation: HA = A + Na+ 2841 RE.add_reaction(gamma=(Ka * K_NACL / K_HCL).magnitude, 2842 reactant_types=[state_one_type], 2843 reactant_coefficients=[1], 2844 product_types=[state_two_type, salt_cation_es_type], 2845 product_coefficients=[1, 1], 2846 default_charges={state_one_type: charge_map[state_one_type], 2847 state_two_type: charge_map[state_two_type], 2848 salt_cation_es_type: salt_cation_charge}) 2849 2850 # Reaction in terms of hydroxide: OH- + HA = A 2851 RE.add_reaction(gamma=(Ka / K_W).magnitude, 2852 reactant_types=[state_one_type, hydroxide_es_type], 2853 reactant_coefficients=[1, 1], 2854 product_types=[state_two_type], 2855 product_coefficients=[1], 2856 default_charges={state_one_type: charge_map[state_one_type], 2857 state_two_type: charge_map[state_two_type], 2858 hydroxide_es_type: hydroxide_charge}) 2859 2860 # Reaction in terms of salt anion: Cl- + HA = A 2861 RE.add_reaction(gamma=(Ka / K_HCL).magnitude, 2862 reactant_types=[state_one_type, salt_anion_es_type], 2863 reactant_coefficients=[1, 1], 2864 product_types=[state_two_type], 2865 product_coefficients=[1], 2866 default_charges={state_one_type: charge_map[state_one_type], 2867 state_two_type: charge_map[state_two_type], 2868 salt_anion_es_type: salt_anion_charge}) 2869 2870 sucessful_reactions_labels.append(name) 2871 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 ('float'): 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', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
- exclusion_range(
float
, 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 ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).
2873 def setup_grxmc_unified (self, pH_res, c_salt_res, cation_name, anion_name, activity_coefficient=None, exclusion_range=None, pka_set=None, use_exclusion_radius_per_type = False): 2874 """ 2875 Sets up Acid/Base reactions for acidic/basic 'particles' defined in 'pmb.df', as well as a grand-canonical coupling to a 2876 reservoir of small ions. 2877 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-}. 2878 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'. 2879 2880 [1] Curk, T., Yuan, J., & Luijten, E. (2022). Accelerated simulation method for charge regulation effects. The Journal of Chemical Physics, 156(4). 2881 [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. 2882 2883 Args: 2884 pH_res ('float'): pH-value in the reservoir. 2885 c_salt_res ('float'): Concentration of monovalent salt (e.g. NaCl) in the reservoir. 2886 cation_name ('str'): Name of the cationic particle. 2887 anion_name ('str'): Name of the anionic particle. 2888 activity_coefficient ('callable', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength. 2889 exclusion_range(`float`, optional): Below this value, no particles will be inserted. 2890 pka_set(`dict`,optional): Desired pka_set, pka_set(`dict`): {"name" : {"pka_value": pka, "acidity": acidity}}. Defaults to None. 2891 use_exclusion_radius_per_type(`bool`,optional): Controls if one exclusion_radius per each espresso_type. Defaults to `False`. 2892 2893 Returns: 2894 RE (`obj`): Instance of a reaction_ensemble.ReactionEnsemble object from the espressomd library. 2895 sucessful_reactions_labels(`lst`): Labels of the reactions set up by pyMBE. 2896 ionic_strength_res ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients). 2897 """ 2898 from espressomd import reaction_methods 2899 if exclusion_range is None: 2900 exclusion_range = max(self.get_radius_map().values())*2.0 2901 if pka_set is None: 2902 pka_set=self.get_pka_set() 2903 self.check_pka_set(pka_set=pka_set) 2904 if use_exclusion_radius_per_type: 2905 exclusion_radius_per_type = self.get_radius_map() 2906 else: 2907 exclusion_radius_per_type = {} 2908 2909 #If no function for the activity coefficient is provided, the ideal case is assumed. 2910 if activity_coefficient is None: 2911 activity_coefficient = lambda x: 1.0 2912 2913 RE = reaction_methods.ReactionEnsemble(kT=self.kT.to('reduced_energy').magnitude, 2914 exclusion_range=exclusion_range.magnitude, 2915 seed=self.SEED, 2916 exclusion_radius_per_type = exclusion_radius_per_type 2917 ) 2918 2919 # Determine the concentrations of the various species in the reservoir and the equilibrium constants 2920 cH_res, cOH_res, cNa_res, cCl_res = self.determine_reservoir_concentrations(pH_res, c_salt_res, activity_coefficient) 2921 ionic_strength_res = 0.5*(cNa_res+cCl_res+cOH_res+cH_res) 2922 determined_activity_coefficient = activity_coefficient(ionic_strength_res) 2923 a_hydrogen = (10 ** (-pH_res) * self.units.mol/self.units.l).to('1/(N_A * reduced_length**3)') 2924 a_cation = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 2925 a_anion = (cH_res+cNa_res).to('1/(N_A * reduced_length**3)') * np.sqrt(determined_activity_coefficient) 2926 K_XX = a_cation * a_anion 2927 2928 cation_es_type = self.df.loc[self.df['name']==cation_name].state_one.es_type.values[0] 2929 anion_es_type = self.df.loc[self.df['name']==anion_name].state_one.es_type.values[0] 2930 cation_charge = self.df.loc[self.df['name']==cation_name].state_one.charge.values[0] 2931 anion_charge = self.df.loc[self.df['name']==anion_name].state_one.charge.values[0] 2932 if cation_charge <= 0: 2933 raise ValueError('ERROR cation charge must be positive, charge ', cation_charge) 2934 if anion_charge >= 0: 2935 raise ValueError('ERROR anion charge must be negative, charge ', anion_charge) 2936 2937 # Coupling to the reservoir: 0 = X+ + X- 2938 RE.add_reaction( 2939 gamma = K_XX.magnitude, 2940 reactant_types = [], 2941 reactant_coefficients = [], 2942 product_types = [ cation_es_type, anion_es_type ], 2943 product_coefficients = [ 1, 1 ], 2944 default_charges = { 2945 cation_es_type: cation_charge, 2946 anion_es_type: anion_charge, 2947 } 2948 ) 2949 2950 sucessful_reactions_labels=[] 2951 charge_map = self.get_charge_map() 2952 for name in pka_set.keys(): 2953 if self.check_if_name_is_defined_in_df(name,pmb_type_to_be_defined='particle') is False : 2954 print('WARNING: the acid-base reaction of ' + name +' has not been set up because its espresso type is not defined in sg.type_map') 2955 continue 2956 2957 Ka = 10**-pka_set[name]['pka_value'] * self.units.mol/self.units.l 2958 gamma_K_AX = Ka.to('1/(N_A * reduced_length**3)').magnitude * a_cation / a_hydrogen 2959 2960 state_one_type = self.df.loc[self.df['name']==name].state_one.es_type.values[0] 2961 state_two_type = self.df.loc[self.df['name']==name].state_two.es_type.values[0] 2962 2963 # Reaction in terms of small cation: HA = A + X+ 2964 RE.add_reaction(gamma=gamma_K_AX.magnitude, 2965 reactant_types=[state_one_type], 2966 reactant_coefficients=[1], 2967 product_types=[state_two_type, cation_es_type], 2968 product_coefficients=[1, 1], 2969 default_charges={state_one_type: charge_map[state_one_type], 2970 state_two_type: charge_map[state_two_type], 2971 cation_es_type: cation_charge}) 2972 2973 # Reaction in terms of small anion: X- + HA = A 2974 RE.add_reaction(gamma=gamma_K_AX.magnitude / K_XX.magnitude, 2975 reactant_types=[state_one_type, anion_es_type], 2976 reactant_coefficients=[1, 1], 2977 product_types=[state_two_type], 2978 product_coefficients=[1], 2979 default_charges={state_one_type: charge_map[state_one_type], 2980 state_two_type: charge_map[state_two_type], 2981 anion_es_type: anion_charge}) 2982 2983 sucessful_reactions_labels.append(name) 2984 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 ('float'): 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', optional): A function that calculates the activity coefficient of an ion pair as a function of the ionic strength.
- exclusion_range(
float
, 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 (
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 ('float'): Ionic strength of the reservoir (useful for calculating partition coefficients).
2986 def setup_df (self): 2987 """ 2988 Sets up the pyMBE's dataframe `pymbe.df`. 2989 2990 Returns: 2991 columns_names(`obj`): pandas multiindex object with the column names of the pyMBE's dataframe 2992 """ 2993 2994 columns_dtypes = { 2995 'name': { 2996 '': str}, 2997 'pmb_type': { 2998 '': str}, 2999 'particle_id': { 3000 '': object}, 3001 'particle_id2': { 3002 '': object}, 3003 'residue_id': { 3004 '': object}, 3005 'molecule_id': { 3006 '': object}, 3007 'acidity': { 3008 '': str}, 3009 'pka': { 3010 '': float}, 3011 'central_bead': { 3012 '': object}, 3013 'side_chains': { 3014 '': object}, 3015 'residue_list': { 3016 '': object}, 3017 'model': { 3018 '': str}, 3019 'sigma': { 3020 '': object}, 3021 'cutoff': { 3022 '': object}, 3023 'offset': { 3024 '': object}, 3025 'epsilon': { 3026 '': object}, 3027 'state_one': { 3028 'label': str, 3029 'es_type': object, 3030 'charge': object }, 3031 'state_two': { 3032 'label': str, 3033 'es_type': object, 3034 'charge': object }, 3035 'sequence': { 3036 '': object}, 3037 'bond_object': { 3038 '': object}, 3039 'parameters_of_the_potential':{ 3040 '': object}, 3041 'l0': { 3042 '': float}, 3043 } 3044 3045 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()])) 3046 3047 for level1, sub_dtypes in columns_dtypes.items(): 3048 for level2, dtype in sub_dtypes.items(): 3049 self.df[level1, level2] = self.df[level1, level2].astype(dtype) 3050 3051 columns_names = pd.MultiIndex.from_frame(self.df) 3052 columns_names = columns_names.names 3053 3054 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
3056 def setup_lj_interactions(self, espresso_system, shift_potential=True, combining_rule='Lorentz-Berthelot', warnings=True): 3057 """ 3058 Sets up the Lennard-Jones (LJ) potential between all pairs of particle types with values for `sigma`, `offset`, and `epsilon` stored in `pymbe.df`. 3059 3060 Args: 3061 espresso_system(`obj`): Instance of a system object from the espressomd library. 3062 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. 3063 combining_rule(`string`, optional): combining rule used to calculate `sigma` and `epsilon` for the potential between a pair of particles. Defaults to 'Lorentz-Berthelot'. 3064 warning(`bool`, optional): switch to activate/deactivate warning messages. Defaults to True. 3065 3066 Note: 3067 - LJ interactions will only be set up between particles with defined values of `sigma` and `epsilon` in the pmb.df. 3068 - Currently, the only `combining_rule` supported is Lorentz-Berthelot. 3069 - Check the documentation of ESPResSo for more info about the potential https://espressomd.github.io/doc4.2.0/inter_non-bonded.html 3070 3071 """ 3072 from itertools import combinations_with_replacement 3073 implemented_combining_rules = ['Lorentz-Berthelot'] 3074 compulsory_parameters_in_df = ['sigma','epsilon'] 3075 # Sanity check 3076 if combining_rule not in implemented_combining_rules: 3077 raise ValueError('In the current version of pyMBE, the only combinatorial rules implemented are ', implemented_combining_rules) 3078 if shift_potential: 3079 shift="auto" 3080 else: 3081 shift=0 3082 # List which particles have sigma and epsilon values defined in pmb.df and which ones don't 3083 particles_types_with_LJ_parameters = [] 3084 non_parametrized_labels= [] 3085 for particle_type in self.get_type_map().values(): 3086 check_list=[] 3087 for key in compulsory_parameters_in_df: 3088 value_in_df=self.find_value_from_es_type(es_type=particle_type, 3089 column_name=key) 3090 check_list.append(np.isnan(value_in_df)) 3091 if any(check_list): 3092 non_parametrized_labels.append(self.find_value_from_es_type(es_type=particle_type, 3093 column_name='label')) 3094 else: 3095 particles_types_with_LJ_parameters.append(particle_type) 3096 lj_parameters_keys=["label","sigma","epsilon","offset","cutoff"] 3097 # Set up LJ interactions between all particle types 3098 for type_pair in combinations_with_replacement(particles_types_with_LJ_parameters, 2): 3099 lj_parameters={} 3100 for key in lj_parameters_keys: 3101 lj_parameters[key]=[] 3102 # Search the LJ parameters of the type pair 3103 for ptype in type_pair: 3104 for key in lj_parameters_keys: 3105 lj_parameters[key].append(self.find_value_from_es_type(es_type=ptype, column_name=key)) 3106 # If one of the particle has sigma=0, no LJ interations are set up between that particle type and the others 3107 if not all(sigma_value.magnitude for sigma_value in lj_parameters["sigma"]): 3108 continue 3109 # Apply combining rule 3110 if combining_rule == 'Lorentz-Berthelot': 3111 sigma=(lj_parameters["sigma"][0]+lj_parameters["sigma"][1])/2 3112 cutoff=(lj_parameters["cutoff"][0]+lj_parameters["cutoff"][1])/2 3113 offset=(lj_parameters["offset"][0]+lj_parameters["offset"][1])/2 3114 epsilon=np.sqrt(lj_parameters["epsilon"][0]*lj_parameters["epsilon"][1]) 3115 espresso_system.non_bonded_inter[type_pair[0],type_pair[1]].lennard_jones.set_params(epsilon = epsilon.to('reduced_energy').magnitude, 3116 sigma = sigma.to('reduced_length').magnitude, 3117 cutoff = cutoff.to('reduced_length').magnitude, 3118 offset = offset.to("reduced_length").magnitude, 3119 shift = shift) 3120 index = len(self.df) 3121 self.df.at [index, 'name'] = f'LJ: {lj_parameters["label"][0]}-{lj_parameters["label"][1]}' 3122 lj_params=espresso_system.non_bonded_inter[type_pair[0], type_pair[1]].lennard_jones.get_params() 3123 3124 self.add_value_to_df(index=index, 3125 key=('pmb_type',''), 3126 new_value='LennardJones') 3127 3128 self.add_value_to_df(index=index, 3129 key=('parameters_of_the_potential',''), 3130 new_value=lj_params, 3131 non_standard_value=True) 3132 if non_parametrized_labels and warnings: 3133 print(f'WARNING: The following particles do not have a defined value of sigma or epsilon in pmb.df: {non_parametrized_labels}. No LJ interaction has been added in ESPResSo for those particles.') 3134 3135 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(
obj
): 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
3137 def setup_particle_sigma(self, topology_dict): 3138 ''' 3139 Sets up sigma of the particles in `topology_dict`. 3140 3141 Args: 3142 topology_dict(`dict`): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value} 3143 ''' 3144 for residue in topology_dict.keys(): 3145 residue_name = re.split(r'\d+', residue)[0] 3146 residue_number = re.split(r'(\d+)', residue)[1] 3147 residue_sigma = topology_dict[residue]['sigma'] 3148 sigma = residue_sigma*self.units.nm 3149 index = self.df[(self.df['residue_id']==residue_number) & (self.df['name']==residue_name) ].index.values[0] 3150 self.add_value_to_df(key= ('sigma',''), 3151 index=int (index), 3152 new_value=sigma) 3153 return
Sets up sigma of the particles in topology_dict
.
Arguments:
- topology_dict(
dict
): {'initial_pos': coords_list, 'chain_id': id, 'sigma': sigma_value}
3155 def write_pmb_df (self, filename): 3156 ''' 3157 Writes the pyMBE dataframe into a csv file 3158 3159 Args: 3160 filename (`str`): Path to the csv file 3161 ''' 3162 3163 self.df.to_csv(filename) 3164 3165 return
Writes the pyMBE dataframe into a csv file
Arguments:
- filename (
str
): Path to the csv file
3167 def write_output_vtf_file(self, espresso_system, filename): 3168 ''' 3169 Writes a snapshot of `espresso_system` on the vtf file `filename`. 3170 3171 Args: 3172 espresso_system(`obj`): Instance of a system object from the espressomd library. 3173 filename(`str`): Path to the vtf file. 3174 3175 ''' 3176 box = espresso_system.box_l[0] 3177 with open(filename, mode='w+t') as coordinates: 3178 coordinates.write (f'unitcell {box} {box} {box} \n') 3179 for particle in espresso_system.part: 3180 type_label = self.find_value_from_es_type(es_type=particle.type, column_name='label') 3181 coordinates.write (f'atom {particle.id} radius 1 name {type_label} type {type_label}\n' ) 3182 coordinates.write ('timestep indexed\n') 3183 for particle in espresso_system.part: 3184 coordinates.write (f'{particle.id} \t {particle.pos[0]} \t {particle.pos[1]} \t {particle.pos[2]}\n') 3185 return
Writes a snapshot of espresso_system
on the vtf file filename
.
Arguments:
- espresso_system(
obj
): Instance of a system object from the espressomd library. - filename(
str
): Path to the vtf file.
36 class NumpyEncoder(json.JSONEncoder): 37 38 """ 39 Custom JSON encoder that converts NumPy arrays to Python lists 40 and NumPy scalars to Python scalars. 41 """ 42 43 def default(self, obj): 44 if isinstance(obj, np.ndarray): 45 return obj.tolist() 46 if isinstance(obj, np.generic): 47 return obj.item() 48 return super().default(obj)
Custom JSON encoder that converts NumPy arrays to Python lists and NumPy scalars to Python scalars.
43 def default(self, obj): 44 if isinstance(obj, np.ndarray): 45 return obj.tolist() 46 if isinstance(obj, np.generic): 47 return obj.item() 48 return super().default(obj)
Implement this method in a subclass such that it returns
a serializable object for o
, or calls the base implementation
(to raise a TypeError
).
For example, to support arbitrary iterators, you could implement default like this::
def default(self, o):
try:
iterable = iter(o)
except TypeError:
pass
else:
return list(iterable)
# Let the base class default method raise the TypeError
return JSONEncoder.default(self, o)
Inherited Members
- json.encoder.JSONEncoder
- JSONEncoder
- item_separator
- key_separator
- skipkeys
- ensure_ascii
- check_circular
- allow_nan
- sort_keys
- indent
- encode
- iterencode