def interpolate(start_node, end_node, num_interp): ''' ''' nifty.printcool(" interpolate") num_nodes = num_interp + 2 nodes = [None] * (num_nodes) nodes[0] = start_node nodes[-1] = end_node sign = 1 nR = 1 nP = 1 nn = nR + nP for n in range(num_interp): if num_nodes - nn > 1: stepsize = 1. / float(num_nodes - nn) else: stepsize = 0.5 if sign == 1: iR = nR - 1 iP = num_nodes - nP iN = nR nodes[nR] = GSM.add_node(nodes[iR], nodes[iP], stepsize, iN) if nodes[nR] is None: raise RuntimeError # print(" Energy of node {} is {:5.4}".format(nR,nodes[nR].energy-E0)) nR += 1 nn += 1 else: n1 = num_nodes - nP n2 = n1 - 1 n3 = nR - 1 nodes[n2] = GSM.add_node(nodes[n1], nodes[n3], stepsize, n2) if nodes[n2] is None: raise RuntimeError # print(" Energy of node {} is {:5.4}".format(nR,nodes[nR].energy-E0)) nP += 1 nn += 1 sign *= -1 return nodes
def check_if_grown(self): ''' Check if the string is grown Returns True if grown ''' self.pastts = self.past_ts() isDone = False # TODO break planes condition1 = (abs(self.nodes[self.nR-1].bdist) <= (1-self.BDIST_RATIO)*abs(self.nodes[0].bdist)) print(" bdist %.3f" % self.nodes[self.nR-1].bdist) fp = self.find_peaks('growing') if self.pastts and self.current_nnodes > 3 and condition1: # TODO extra criterion here print(" pastts is ", self.pastts) if self.TSnode == self.nR-1: print(" The highest energy node is the last") print(" not continuing with TS optimization.") self.tscontinue = False nifty.printcool("Over the hill") isDone = True elif fp == -1 and self.energies[self.nR-1] > 200. and self.nodes[self.nR-1].gradrms > self.options['CONV_TOL']*5: print("growth_iters over: all uphill and high energy") self.end_early = 2 self.tscontinue = False self.nnodes = self.nR isDone = True elif fp == -2: print("growth_iters over: all uphill and flattening out") self.end_early = 2 self.tscontinue = False self.nnodes = self.nR isDone = True # ADD extra criteria here to check if TS is higher energy than product return isDone
if __name__ == "__main__": import psiw from utilities import nifty ##### => Job Data <= ##### states = [(1, 0), (1, 1)] charge = 0 nocc = 7 nactive = 2 basis = '6-31gs' filepath = '../../data/ethylene.xyz' #### => PSIW Obj <= ###### nifty.printcool("Build resources") resources = ls.ResourceList.build() nifty.printcool('{}'.format(resources)) molecule = ls.Molecule.from_xyz_file(filepath) geom = psiw.geometry.Geometry.build( resources=resources, molecule=molecule, basisname=basis, ) nifty.printcool('{}'.format(geom)) ref = psiw.RHF.from_options( geometry=geom, g_convergence=1.0E-6, fomo=True,
def __init__(self, options): super(OpenMM, self).__init__(options) # get simulation from options if it exists # self.options['job_data']['simulation'] = self.options['job_data'].get('simulation',None) self.simulation = self.options['job_data'].get('simulation', None) if self.lot_inp_file is not None and self.simulation is None: # Now go through the logic of determining which FILE options are activated. self.file_options.set_active('use_crystal', False, bool, "Use crystal unit parameters") self.file_options.set_active( 'use_pme', False, bool, '', "Use particle mesh ewald-- requires periodic boundary conditions" ) self.file_options.set_active('cutoff', 1.0, float, '', depend=(self.file_options.use_pme), msg="Requires PME") self.file_options.set_active('prmtopfile', None, str, "parameter file") self.file_options.set_active('inpcrdfile', None, str, "inpcrd file") self.file_options.set_active('restrain_bondfile', None, str, 'list of bonds to restrain') self.file_options.set_active('restrain_torfile', None, str, "list of torsions to restrain") self.file_options.set_active('restrain_tranfile', None, str, "list of translations to restrain") for line in self.file_options.record(): print(line) # set all active values to self for easy access for key in self.file_options.ActiveOptions: setattr(self, key, self.file_options.ActiveOptions[key]) nifty.printcool(" Options for OpenMM") for val in [self.prmtopfile, self.inpcrdfile]: assert val is not None, "Missing prmtop or inpcrdfile" # Integrator will never be used (Simulation requires one) integrator = openmm.VerletIntegrator(1.0) # create simulation object if self.use_crystal: crystal = load_file(self.prmtopfile, self.inpcrdfile) if self.use_pme: system = crystal.createSystem( nonbondedMethod=openmm_app.PME, nonbondedCutoff=self.cutoff * openmm_units.nanometer, ) else: system = crystal.createSystem( nonbondedMethod=openmm_app.NoCutoff, ) # Add restraints self.add_restraints(system) self.simulation = openmm_app.Simulation( crystal.topology, system, integrator) # set the box vectors inpcrd = openmm_app.AmberInpcrdFile(self.inpcrdfile) if inpcrd.boxVectors is not None: print(" setting box vectors") print(inpcrd.boxVectors) self.simulation.context.setPeriodicBoxVectors( *inpcrd.boxVectors) else: # Do not use crystal parameters prmtop = openmm_app.AmberPrmtopFile(self.prmtopfile) if self.use_pme: system = prmtop.createSystem( nonbondedMethod=openmm_app.PME, nonbondedCutoff=self.cutoff * openmm_units.nanometer, ) else: system = prmtop.createSystem( nonbondedMethod=openmm_app.NoCutoff, ) # add restraints self.add_restraints(system) self.simulation = openmm_app.Simulation( prmtop.topology, system, integrator, )
def add_restraints(self, system): # Bond Restraints if self.restrain_bondfile is not None: nifty.printcool(" Adding bonding restraints!") # Harmonic constraint flat_bottom_force = openmm.CustomBondForce( 'step(r-r0) * (k/2) * (r-r0)^2') flat_bottom_force.addPerBondParameter('r0') flat_bottom_force.addPerBondParameter('k') system.addForce(flat_bottom_force) with open(self.restrain_bondfile, 'r') as input_file: for line in input_file: print(line) columns = line.split() atom_index_i = int(columns[0]) atom_index_j = int(columns[1]) r0 = float(columns[2]) k = float(columns[3]) flat_bottom_force.addBond(atom_index_i, atom_index_j, [r0, k]) # Torsion restraint if self.restrain_torfile is not None: nifty.printcool(" Adding torsional restraints!") # Harmonic constraint tforce = openmm.CustomTorsionForce( "0.5*k*min(dtheta, 2*pi-dtheta)^2; dtheta = abs(theta-theta0); pi = 3.1415926535" ) tforce.addPerTorsionParameter("k") tforce.addPerTorsionParameter("theta0") system.addForce(tforce) xyz = manage_xyz.xyz_to_np(self.geom) with open(self.restrain_torfile, 'r') as input_file: for line in input_file: columns = line.split() a = int(columns[0]) b = int(columns[1]) c = int(columns[2]) d = int(columns[3]) k = float(columns[4]) dih = Dihedral(a, b, c, d) theta0 = dih.value(xyz) tforce.addTorsion(a, b, c, d, [k, theta0]) # Translation restraint if self.restrain_tranfile is not None: nifty.printcool(" Adding translational restraints!") trforce = openmm.CustomExternalForce( "k*periodicdistance(x, y, z, x0, y0, z0)^2") trforce.addPerParticleParameter("k") trforce.addPerParticleParameter("x0") trforce.addPerParticleParameter("y0") trforce.addPerParticleParameter("z0") system.addForce(trforce) xyz = manage_xyz.xyz_to_np(self.geom) with open(self.restrain_tranfile, 'r') as input_file: for line in input_file: columns = line.split() a = int(columns[0]) k = float(columns[1]) x0 = xyz[a, 0] * 0.1 # Units are in nm y0 = xyz[a, 1] * 0.1 # Units are in nm z0 = xyz[a, 2] * 0.1 # Units are in nm trforce.addParticle(a, [k, x0, y0, z0])
def __init__(self, options): super(pDynamo, self).__init__(options) #pDynamo self.system = self.options['job_data'].get('system', None) print(" making folder scratch/{}".format(self.node_id)) os.system('mkdir -p scratch/{}'.format(self.node_id)) # if simulation doesn't exist create it if self.lot_inp_file is not None and self.simulation is None: # Now go through the logic of determining which FILE options are activated. # DO NOT DUPLICATE OPTIONS WHICH ARE ALREADY PART OF LOT OPTIONS (e.g. charge) # ORCA options self.file_options.set_active( 'use_orca', False, bool, "Use ORCA to evaluate energies and gradients") self.file_options.set_active( 'orca_method', "B3LYP", str, "Method to use in ORCA e.g. HF, MP2, or density functional", depend=(self.file_options.use_orca), msg="Must use ORCA to specify density functional") self.file_options.set_active( 'basis', "6-31g", str, "Basis set for wavefunction or density functional", depend=(self.file_options.use_orca), msg="Must use ORCA to specify density functional") self.file_options.set_active('tole', 1e-4, float, "Energy tolerance for convergence") self.file_options.set_active('maxiter', 100, int, "Number of SCF cycles") self.file_options.set_active('slowconv', False, bool, "Convergence option for ORCA") self.file_options.set_active( 'scfconvtol', 'NormalSCF', str, "Convergence option for ORCA", allowed=['NormalSCF', 'TightSCF', 'ExtremeSCF']) self.file_options.set_active('d3', False, bool, "Use Grimme's D3 dispersion") # QM/MM CHARMM self.file_options.set_active('qmatom_file', None, str, '') self.file_options.set_active( 'use_charmm_qmmm', False, bool, 'Use CHARMM molecular mechanics parameters to perform QMMM', depend=(self.file_options.qmatom_file is not None), msg="Must define qm atoms") self.file_options.set_active( 'path_to_prm', None, str, 'path to folder containing Charmm parameter files') self.file_options.set_active( 'path_to_str', None, str, 'path to folder containing to Charmm str files') self.file_options.set_active('psf_file', None, str, 'Path to file containing CHARMM PSF') self.file_options.set_active('crd_file', None, str, 'Path to file containing CHARMM CRD') # DFTB options self.file_options.set_active( 'use_dftb', False, bool, "Use DFTB to evaluate energies and gradients", clash=(self.file_options.use_orca), msg="We're not using DFTB+") self.file_options.set_active( 'path_to_skf', None, str, 'path to folder containing skf files') self.file_options.set_active('use_scc', True, bool, "Use self-consistent charge") # General options self.file_options.set_active( 'command', None, str, 'pDynamo requires a path to an executable like ORCA or DFTB+') self.file_options.set_active('scratch', None, str, 'Folder to store temporary files') self.file_options.force_active('scratch', 'scratch/{}'.format(self.node_id), 'Setting scratch folder') nifty.printcool(" Options for pdynamo") for line in self.file_options.record(): print(line) # Build system self.build_system() # set all active values to self for easy access for key in self.file_options.ActiveOptions: setattr(self, key, self.file_options.ActiveOptions[key])
def ic_reparam(nodes, energies, climbing=False, ic_reparam_steps=8, print_level=1, NUM_CORE=1, MAXRE=0.25): ''' Reparameterizes the string using Delocalizedin internal coordinatesusing three-way tangents at the TS node Only pushes nodes outwards during reparameterization because otherwise too many things change. Be careful, however, if the path is allup or alldown then this can cause Parameters ---------- nodes : list of molecule objects energies : list of energies in kcal/mol ic_reparam_steps : int max number of reparameterization steps print_level : int verbosity ''' nifty.printcool("reparametrizing string nodes") nnodes = len(nodes) rpart = np.zeros(nnodes) for n in range(1, nnodes): rpart[n] = 1. / (nnodes - 1) deltadqs = np.zeros(nnodes) TSnode = np.argmax(energies) disprms = 100 if ((TSnode == nnodes - 1) or (TSnode == 0)) and climbing: raise RuntimeError(" TS node shouldn't be the first or last node") ideal_progress_gained = np.zeros(nnodes) if climbing: for n in range(1, TSnode): ideal_progress_gained[n] = 1. / (TSnode) for n in range(TSnode + 1, nnodes): ideal_progress_gained[n] = 1. / (nnodes - TSnode - 1) ideal_progress_gained[TSnode] = 0. else: for n in range(1, nnodes): ideal_progress_gained[n] = 1. / (nnodes - 1) for i in range(ic_reparam_steps): ictan, dqmaga = GSM.get_tangents(nodes) totaldqmag = np.sum(dqmaga) if climbing: progress = np.zeros(nnodes) progress_gained = np.zeros(nnodes) h1dqmag = np.sum(dqmaga[:TSnode + 1]) h2dqmag = np.sum(dqmaga[TSnode + 1:nnodes]) if print_level > 0: print(" h1dqmag, h2dqmag: %3.2f %3.2f" % (h1dqmag, h2dqmag)) progress_gained[:TSnode] = dqmaga[:TSnode] / h1dqmag progress_gained[TSnode + 1:] = dqmaga[TSnode + 1:] / h2dqmag progress[:TSnode] = np.cumsum(progress_gained[:TSnode]) progress[TSnode:] = np.cumsum(progress_gained[TSnode:]) else: progress = np.cumsum(dqmaga) / totaldqmag progress_gained = dqmaga / totaldqmag if i == 0: orig_dqmaga = copy(dqmaga) orig_progress_gained = copy(progress_gained) if climbing: difference = np.zeros(nnodes) for n in range(TSnode): difference[ n] = ideal_progress_gained[n] - progress_gained[n] deltadqs[n] = difference[n] * h1dqmag for n in range(TSnode + 1, nnodes): difference[ n] = ideal_progress_gained[n] - progress_gained[n] deltadqs[n] = difference[n] * h2dqmag else: difference = ideal_progress_gained - progress_gained deltadqs = difference * totaldqmag if print_level > 1: print(" ideal progress gained per step", end=' ') for n in range(nnodes): print(" step [{}]: {:1.3f}".format( n, ideal_progress_gained[n]), end=' ') print() print(" path progress ", end=' ') for n in range(nnodes): print(" step [{}]: {:1.3f}".format(n, progress_gained[n]), end=' ') print() print(" difference ", end=' ') for n in range(nnodes): print(" step [{}]: {:1.3f}".format(n, difference[n]), end=' ') print() print(" deltadqs ", end=' ') for n in range(nnodes): print(" step [{}]: {:1.3f}".format(n, deltadqs[n]), end=' ') print() # disprms = np.linalg.norm(deltadqs)/np.sqrt(nnodes-1) disprms = np.linalg.norm(deltadqs) / np.sqrt(nnodes - 1) print(" disprms: {:1.3}\n".format(disprms)) if disprms < 0.02: break # Move nodes if climbing: deltadqs[TSnode - 2] -= deltadqs[TSnode - 1] deltadqs[nnodes - 2] -= deltadqs[nnodes - 1] for n in range(1, nnodes - 1): if abs(deltadqs[n]) > MAXRE: deltadqs[n] = np.sign(deltadqs[n]) * MAXRE for n in range(TSnode - 1): deltadqs[n + 1] += deltadqs[n] for n in range(TSnode + 1, nnodes - 2): deltadqs[n + 1] += deltadqs[n] for n in range(nnodes): if abs(deltadqs[n]) > MAXRE: deltadqs[n] = np.sign(deltadqs[n]) * MAXRE if NUM_CORE > 1: # 5/14/2021 TS node f***s this up?! tans = [ ictan[n] if deltadqs[n] < 0 else ictan[n + 1] for n in chain(range(1, TSnode), range(TSnode + 1, nnodes - 1)) ] # + [ ictan[n] if deltadqs[n]<0 else ictan[n+1] for n in range(TSnode+1,nnodes-1)] pool = mp.Pool(NUM_CORE) Vecs = pool.map( worker, ((nodes[0].coord_obj, "build_dlc", node.xyz, tan) for node, tan in zip( nodes[1:TSnode] + nodes[TSnode + 1:nnodes - 1], tans))) pool.close() pool.join() for n, node in enumerate(nodes[1:TSnode] + nodes[TSnode + 1:nnodes - 1]): node.coord_basis = Vecs[n] # move the positions dqs = [ deltadqs[n] * nodes[n].constraints[:, 0] for n in chain(range(1, TSnode), range(TSnode + 1, nnodes - 1)) ] pool = mp.Pool(NUM_CORE) newXyzs = pool.map( worker, ((node.coord_obj, "newCartesian", node.xyz, dq) for node, dq in zip( nodes[1:TSnode] + nodes[TSnode + 1:nnodes - 1], dqs))) pool.close() pool.join() for n, node in enumerate(nodes[1:TSnode] + nodes[TSnode + 1:nnodes - 1]): node.xyz = newXyzs[n] else: for n in chain(range(1, TSnode), range(TSnode + 1, nnodes - 1)): if deltadqs[n] < 0: # print(f" Moving node {n} along tan[{n}] this much {deltadqs[n]}") print(" Moving node {} along tan[{}] this much {}". format(n, n, deltadqs[n])) nodes[n].update_coordinate_basis(ictan[n]) constraint = nodes[n].constraints[:, 0] dq = deltadqs[n] * constraint nodes[n].update_xyz(dq, verbose=(print_level > 1)) elif deltadqs[n] > 0: print(" Moving node {} along tan[{}] this much {}". format(n, n + 1, deltadqs[n])) nodes[n].update_coordinate_basis(ictan[n + 1]) constraint = nodes[n].constraints[:, 0] dq = deltadqs[n] * constraint nodes[n].update_xyz(dq, verbose=(print_level > 1)) else: # e.g 11-2 = 9, deltadq[9] -= deltadqs[10] deltadqs[nnodes - 2] -= deltadqs[nnodes - 1] for n in range(1, nnodes - 1): if abs(deltadqs[n]) > MAXRE: deltadqs[n] = np.sign(deltadqs[n]) * MAXRE for n in range(1, nnodes - 2): deltadqs[n + 1] += deltadqs[n] for n in range(1, nnodes - 1): if abs(deltadqs[n]) > MAXRE: deltadqs[n] = np.sign(deltadqs[n]) * MAXRE if NUM_CORE > 1: # Update the coordinate basis tans = [ ictan[n] if deltadqs[n] < 0 else ictan[n + 1] for n in range(1, nnodes - 1) ] pool = mp.Pool(NUM_CORE) Vecs = pool.map( worker, ((nodes[0].coord_obj, "build_dlc", node.xyz, tan) for node, tan in zip(nodes[1:nnodes - 1], tans))) pool.close() pool.join() for n, node in enumerate(nodes[1:nnodes - 1]): node.coord_basis = Vecs[n] # move the positions dqs = [ deltadqs[n] * nodes[n].constraints[:, 0] for n in range(1, nnodes - 1) ] pool = mp.Pool(NUM_CORE) newXyzs = pool.map( worker, ((node.coord_obj, "newCartesian", node.xyz, dq) for node, dq in zip(nodes[1:nnodes - 1], dqs))) pool.close() pool.join() for n, node in enumerate(nodes[1:nnodes - 1]): node.xyz = newXyzs[n] else: for n in range(1, nnodes - 1): if deltadqs[n] < 0: # print(f" Moving node {n} along tan[{n}] this much {deltadqs[n]}") print(" Moving node {} along tan[{}] this much {}". format(n, n, deltadqs[n])) nodes[n].update_coordinate_basis(ictan[n]) constraint = nodes[n].constraints[:, 0] dq = deltadqs[n] * constraint nodes[n].update_xyz(dq, verbose=(print_level > 1)) elif deltadqs[n] > 0: print(" Moving node {} along tan[{}] this much {}". format(n, n + 1, deltadqs[n])) nodes[n].update_coordinate_basis(ictan[n + 1]) constraint = nodes[n].constraints[:, 0] dq = deltadqs[n] * constraint nodes[n].update_xyz(dq, verbose=(print_level > 1)) if climbing: ictan, dqmaga = GSM.get_tangents(nodes) h1dqmag = np.sum(dqmaga[:TSnode + 1]) h2dqmag = np.sum(dqmaga[TSnode + 1:nnodes]) if print_level > 0: print(" h1dqmag, h2dqmag: %3.2f %3.2f" % (h1dqmag, h2dqmag)) progress_gained[:TSnode] = dqmaga[:TSnode] / h1dqmag progress_gained[TSnode + 1:] = dqmaga[TSnode + 1:] / h2dqmag progress[:TSnode] = np.cumsum(progress_gained[:TSnode]) progress[TSnode:] = np.cumsum(progress_gained[TSnode:]) else: ictan, dqmaga = GSM.get_tangents(nodes) totaldqmag = np.sum(dqmaga) progress = np.cumsum(dqmaga) / totaldqmag progress_gained = dqmaga / totaldqmag print() if print_level > 0: print(" ideal progress gained per step", end=' ') for n in range(nnodes): print(" step [{}]: {:1.3f}".format(n, ideal_progress_gained[n]), end=' ') print() print(" original path progress ", end=' ') for n in range(nnodes): print(" step [{}]: {:1.3f}".format(n, orig_progress_gained[n]), end=' ') print() print(" reparameterized path progress ", end=' ') for n in range(nnodes): print(" step [{}]: {:1.3f}".format(n, progress_gained[n]), end=' ') print() print(" spacings (begin ic_reparam, steps", end=' ') for n in range(nnodes): print(" {:1.2}".format(orig_dqmaga[n]), end=' ') print() print(" spacings (end ic_reparam, steps: {}/{}):".format( i + 1, ic_reparam_steps), end=' ') for n in range(nnodes): print(" {:1.2}".format(dqmaga[n]), end=' ') print("\n disprms: {:1.3}".format(disprms)) return
def go_gsm(self, max_iters=50, opt_steps=3, rtype=0): """rtype=0 MECI search rtype=1 MESX search """ assert rtype in [0, 1], "rtype not defined" if rtype == 0: nifty.printcool("Doing SE-MECI search") else: nifty.printcool("Doing SE-MESX search") self.nodes[0].gradrms = 0. self.nodes[0].V0 = self.nodes[0].energy print(' Initial energy is {:1.4f}'.format(self.nodes[0].energy)) sys.stdout.flush() # stash bdist for node 0 _, self.nodes[0].bdist = self.get_tangent( self.nodes[0], None, driving_coords=self.driving_coords) print(" Initial bdist is %1.3f" % self.nodes[0].bdist) # interpolate first node self.add_GSM_nodeR() # grow string self.grow_string(max_iters=max_iters, max_opt_steps=opt_steps) print(' SE_Cross growth phase over') print(' Warning last node still not fully optimized') if True: path = os.path.join( os.getcwd(), 'scratch/{:03d}/{}'.format(self.ID, self.nR - 1)) # doing extra constrained penalty optimization for MECI print(" extra constrained optimization for the nnR-1 = %d" % (self.nR - 1)) self.optimizer[self.nR - 1].conv_grms = self.options['CONV_TOL'] * 5 ictan = self.get_tangent_xyz( self.nodes[self.nR - 1].xyz, self.nodes[self.nR - 2].xyz, self.newic.primitive_internal_coordinates) self.nodes[self.nR - 1].PES.sigma = 1.5 self.optimizer[self.nR - 1].optimize( molecule=self.nodes[self.nR - 1], refE=self.nodes[0].V0, opt_type='ICTAN', opt_steps=5, ictan=ictan, path=path, ) ictan = self.get_tangent_xyz( self.nodes[self.nR - 1].xyz, self.nodes[self.nR - 2].xyz, self.newic.primitive_internal_coordinates) self.nodes[self.nR - 1].PES.sigma = 2.5 self.optimizer[self.nR - 1].optimize( molecule=self.nodes[self.nR - 1], refE=self.nodes[0].V0, opt_type='ICTAN', opt_steps=5, ictan=ictan, path=path, ) ictan = self.get_tangent_xyz( self.nodes[self.nR - 1].xyz, self.nodes[self.nR - 2].xyz, self.newic.primitive_internal_coordinates) self.nodes[self.nR - 1].PES.sigma = 3.5 self.optimizer[self.nR - 1].optimize( molecule=self.nodes[self.nR - 1], refE=self.nodes[0].V0, opt_type='ICTAN', opt_steps=5, ictan=ictan, path=path, ) self.xyz_writer('after_penalty_{:03}.xyz'.format(self.ID), self.geometries, self.energies, self.gradrmss, self.dEs) self.optimizer[self.nR].opt_cross = True self.nodes[0].V0 = self.nodes[0].PES.PES2.energy if rtype == 0: # MECI optimization self.nodes[self.nR] = Molecule.copy_from_options( self.nodes[self.nR - 1], new_node_id=self.nR) avg_pes = Avg_PES.create_pes_from(self.nodes[self.nR].PES) self.nodes[self.nR].PES = avg_pes path = os.path.join(os.getcwd(), 'scratch/{:03d}/{}'.format(self.ID, self.nR)) self.optimizer[self.nR].conv_grms = self.options['CONV_TOL'] self.optimizer[ self.nR].conv_gmax = 0.1 # self.options['CONV_gmax'] self.optimizer[self.nR].conv_Ediff = self.options['CONV_Ediff'] self.optimizer[self.nR].conv_dE = self.options['CONV_dE'] self.optimizer[self.nR].optimize( molecule=self.nodes[self.nR], refE=self.nodes[0].V0, opt_type='MECI', opt_steps=100, verbose=True, path=path, ) if not self.optimizer[self.nR].converged: print( "doing extra optimization in hopes that the MECI will converge." ) if self.nodes[self.nR].PES.PES2.energy - self.nodes[0].V0 < 20: self.optimizer[self.nR].optimize( molecule=self.nodes[self.nR], refE=self.nodes[0].V0, opt_type='MECI', opt_steps=100, verbose=True, path=path, ) else: # unconstrained penalty optimization # TODO make unctonstrained "CROSSING" which checks for dE convergence self.nodes[self.nR] = Molecule.copy_from_options( self.nodes[self.nR - 1], new_node_id=self.nR) self.nodes[self.nR].PES.sigma = 10.0 print(" sigma for node %d is %.3f" % (self.nR, self.nodes[self.nR].PES.sigma)) path = os.path.join(os.getcwd(), 'scratch/{:03d}/{}'.format(self.ID, self.nR)) self.optimizer[self.nR].opt_cross = True self.optimizer[self.nR].conv_grms = self.options['CONV_TOL'] # self.optimizer[self.nR].conv_gmax = self.options['CONV_gmax'] self.optimizer[self.nR].conv_Ediff = self.options['CONV_Ediff'] self.optimizer[self.nR].conv_dE = self.options['CONV_dE'] self.optimizer[self.nR].optimize( molecule=self.nodes[self.nR], refE=self.nodes[0].V0, opt_type='UNCONSTRAINED', opt_steps=200, verbose=True, path=path, ) self.xyz_writer('grown_string_{:03}.xyz'.format(self.ID), self.geometries, self.energies, self.gradrmss, self.dEs) if self.optimizer[self.nR].converged: self.nnodes = self.nR + 1 self.nodes = self.nodes[:self.nnodes] print("Setting all interior nodes to active") for n in range(1, self.nnodes - 1): self.active[n] = True self.active[self.nnodes - 1] = False self.active[0] = False # Convert all the PES to excited-states for n in range(self.nnodes): self.nodes[n].PES = PES.create_pes_from( self.nodes[n].PES.PES2, options={'gradient_states': [(1, 1)]}) print(" initial ic_reparam") self.reparameterize(ic_reparam_steps=25) print(" V_profile (after reparam): ", end=' ') energies = self.energies for n in range(self.nnodes): print(" {:7.3f}".format(float(energies[n])), end=' ') print() self.xyz_writer('grown_string1_{:03}.xyz'.format(self.ID), self.geometries, self.energies, self.gradrmss, self.dEs) deltaE = energies[-1] - energies[0] if deltaE > 20: print( " MECI energy is too high %5.4f. Don't try to optimize pathway" % deltaE) print("Exiting early") self.end_early = True else: print(" deltaE s1-minimum and MECI %5.4f" % deltaE) try: self.optimize_string(max_iter=max_iters, opt_steps=3, rtype=1) except Exception as error: if str(error) == "Ran out of iterations": print(error) self.end_early = True else: print(error) self.end_early = True else: print("Exiting early") self.end_early = True
def go_gsm(self, max_iters=50, opt_steps=3, rtype=2): """ rtype=2 Find and Climb TS, 1 Climb with no exact find, 0 turning of climbing image and TS search """ self.set_V0() if not self.isRestarted: if self.growth_direction == 0: self.add_GSM_nodes(2) elif self.growth_direction == 1: self.add_GSM_nodeR(1) elif self.growth_direction == 2: self.add_GSM_nodeP(1) # Grow String self.grow_string(max_iters=max_iters, max_opt_steps=opt_steps) nifty.printcool("Done Growing the String!!!") self.done_growing = True # nifty.printcool("initial ic_reparam") self.reparameterize() self.xyz_writer('grown_string_{:03}.xyz'.format(self.ID), self.geometries, self.energies, self.gradrmss, self.dEs) # Can check for intermediate at beginning but not doing that now. # else: # if self.has_intermediate(self.noise): # nifty.printcool(f" WARNING THIS REACTION HAS AN INTERMEDIATE within noise {self.noise}, opting out") # try: # self.optimize_string(max_iter=3,opt_steps=opt_steps,rtype=0) # except Exception as error: # print(" Done optimizing 3 times, checking if intermediate still exists") # if self.has_intermediate(self.noise): # self.tscontinue=False if self.tscontinue: try: self.optimize_string(max_iter=max_iters, opt_steps=opt_steps, rtype=rtype) except Exception as error: if str(error) == "Ran out of iterations": print(error) self.end_early = True else: print(error) self.end_early = True else: print("Exiting early") self.end_early = True filename = "opt_converged_{:03d}.xyz".format(self.ID) print(" Printing string to " + filename) self.xyz_writer(filename, self.geometries, self.energies, self.gradrmss, self.dEs) print("Finished GSM!") return
def build_topology(xyz, atoms, add_bond=None, hybrid_indices=None, bondlistfile=None, prim_idx_start_stop=None, **kwargs): """ Create topology and fragments; these are graph representations of the individual molecule fragments contained in the Molecule object. Parameters ---------- force_bonds : bool Build the bonds from interatomic distances. If the user calls build_topology from outside, assume this is the default behavior. If creating a Molecule object using __init__, do not force the building of bonds by default (only build bonds if not read from file.) topframe : int, optional Provide the frame number used for reading the bonds. If not provided, this will be taken from the top_settings field. If provided, this will take priority and write the value into top_settings. """ natoms = len(atoms) # can do an assert for xyz here CRA TODO if natoms > 100000: nifty.logger.warning( "Warning: Large number of atoms (%i), topology building may take a long time" % natoms) # Get hybrid indices hybrid_indices = hybrid_indices hybrid_idx_start_stop = [] if hybrid_indices is None: primitive_indices = range(len(atoms)) else: # specify Hybrid TRIC we need to specify which atoms to build topology for primitive_indices = [] for i in range(len(atoms)): if i not in hybrid_indices: primitive_indices.append(i) # print("non-cartesian indices") # print(primitive_indices) # get the hybrid start and stop indices new = True for i in range(natoms + 1): if i in hybrid_indices: if new: start = i new = False else: if not new: end = i - 1 new = True hybrid_idx_start_stop.append((start, end)) if not bondlistfile: nifty.printcool(" building bonds") print(prim_idx_start_stop) bonds = Topology.build_bonds(xyz, atoms, primitive_indices, prim_idx_start_stop) # print(" done") assert bondlistfile is None else: # bondlistfile: # prim_idx_start_stop = kwargs.get('prim_idx_start_stop',None) try: bonds = Topology.read_bonds_from_file( bondlistfile) # ,prim_idx_start_stop) except: raise RuntimeError if add_bond: print(" adding extra bonds") for bond in add_bond: bonds.append(bond) # print(bonds) # Create a NetworkX graph object to hold the bonds. G = MyG() for i in primitive_indices: element = atoms[i] a = element.symbol G.add_node(i) if parse_version(nx.__version__) >= parse_version('2.0'): nx.set_node_attributes(G, {i: a}, name='e') nx.set_node_attributes(G, {i: xyz[i]}, name='x') else: nx.set_node_attributes(G, 'e', {i: a}) nx.set_node_attributes(G, 'x', {i: xyz[i]}) for (i, j) in bonds: G.add_edge(i, j) # The Topology is simply the NetworkX graph object. # topology = G fragments = [G.subgraph(c).copy() for c in nx.connected_components(G)] for g in fragments: g.__class__ = MyG # print(len(fragments)) # for frag in fragments: # print(frag.L()) # print("nodes of Graph") # print(topology.L()) # Deprecated in networkx 2.2 # fragments = list(nx.connected_component_subgraphs(G)) # show topology # import matplotlib as mpl # mpl.use('Agg') # import matplotlib.pyplot as plt # plt.plot() # nx.draw(topology,with_labels=True,font_weight='bold') # plt.show() # plt.savefig('tmp.png') return G