Esempio n. 1
0
    def add_last_node(self, rtype):
        assert rtype == 1 or rtype == 2, "rtype must be 1 or 2"
        samegeom = False
        noptsteps = 100
        if self.nodes[self.nR - 1].PES.lot.do_coupling:
            opt_type = 'MECI'
        else:
            opt_type = 'UNCONSTRAINED'

        if rtype == 1:
            print(" copying last node, opting")
            #self.nodes[self.nR] = DLC.copy_node(self.nodes[self.nR-1],self.nR)
            self.nodes[self.nR] = Molecule.copy_from_options(
                self.nodes[self.nR - 1], new_node_id=self.nR)
            print(" Optimizing node %i" % self.nR)
            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']
            path = os.path.join(os.getcwd(),
                                'scratch/{:03d}/{}'.format(self.ID, self.nR))
            self.optimizer[self.nR].optimize(
                molecule=self.nodes[self.nR],
                refE=self.nodes[0].V0,
                opt_steps=noptsteps,
                opt_type=opt_type,
                path=path,
            )
            self.active[self.nR] = True
            if (self.nodes[self.nR].xyz == self.nodes[self.nR - 1].xyz).all():
                print(" Opt did not produce new geometry")
            else:
                self.nR += 1
        elif rtype == 2:
            print(" already created node, opting")
            self.optimizer[self.nR - 1].conv_grms = self.options['CONV_TOL']
            self.optimizer[self.nR - 1].conv_gmax = self.options['CONV_gmax']
            self.optimizer[self.nR - 1].conv_Ediff = self.options['CONV_Ediff']
            self.optimizer[self.nR - 1].conv_dE = self.options['CONV_dE']
            path = os.path.join(
                os.getcwd(), 'scratch/{:03d}/{}'.format(self.ID, self.nR - 1))
            self.optimizer[self.nR - 1].optimize(
                molecule=self.nodes[self.nR - 1],
                refE=self.nodes[0].V0,
                opt_steps=noptsteps,
                opt_type=opt_type,
                path=path,
            )
        #print(" Aligning")
        #self.nodes[self.nR-1].xyz = self.com_rotate_move(self.nR-2,self.nR,self.nR-1)
        return
Esempio n. 2
0
    def check_if_grown(self):
        isDone = False
        if self.nn == self.nnodes:
            isDone = True

            # TODO should something be done for growthdirection 2?
            if self.growth_direction == 1:
                print("Setting LOT of last node")
                self.nodes[-1] = Molecule.copy_from_options(
                    MoleculeA=self.nodes[-2],
                    xyz=self.nodes[-1].xyz,
                    new_node_id=self.nnodes - 1)

        return isDone
Esempio n. 3
0
def main():
    parser = argparse.ArgumentParser(
        description="Reaction path transition state and photochemistry tool",
        formatter_class=argparse.RawDescriptionHelpFormatter,
        epilog=textwrap.dedent('''\
                Example of use:
                --------------------------------
                gsm -mode DE_GSM -xyzfile yourfile.xyz -package QChem -lot_inp_file qstart -ID 1
                ''')
        )
    parser.add_argument('-xyzfile', help='XYZ file containing reactant and, if DE-GSM, product.',  required=True)
    parser.add_argument('-isomers', help='driving coordinate file', type=str, required=False)
    parser.add_argument('-mode', default="DE_GSM",help='GSM Type (default: %(default)s)',choices=["DE_GSM","SE_GSM","SE_Cross"], type=str, required=True)
    parser.add_argument('-only_drive',action='store_true',help='')
    parser.add_argument('-package',default="QChem",type=str,help="Electronic structure theory package (default: %(default)s)",choices=["QChem","Orca","Molpro","PyTC","TeraChemCloud","OpenMM","DFTB","TeraChem","BAGEL","xTB_lot"])
    parser.add_argument('-lot_inp_file',type=str,default=None, help='external file to specify calculation e.g. qstart,gstart,etc. Highly package specific.',required=False)
    parser.add_argument('-ID',default=0, type=int,help='string identification number (default: %(default)s)',required=False)
    parser.add_argument('-num_nodes',type=int,default=11,help='number of nodes for string (defaults: 9 DE-GSM, 20 SE-GSM)',required=False)
    parser.add_argument('-pes_type',type=str,default='PES',help='Potential energy surface (default: %(default)s)',choices=['PES','Avg_PES','Penalty_PES'])
    parser.add_argument('-adiabatic_index',nargs="*",type=int,default=[0],help='Adiabatic index (default: %(default)s)',required=False)
    parser.add_argument('-multiplicity',nargs="*",type=int,default=[1],help='Multiplicity (default: %(default)s)')
    parser.add_argument('-FORCE_FILE',type=str,default=None,help='Constant force between atoms in AU,e.g. [(1,2,0.1214)]. Negative is tensile, positive is compresive')
    parser.add_argument('-RESTRAINT_FILE',type=str,default=None,help='Harmonic translational restraints')
    parser.add_argument('-optimizer',type=str,default='eigenvector_follow',help='The optimizer object. (default: %(default)s Recommend LBFGS for large molecules >1000 atoms)',required=False)
    parser.add_argument('-opt_print_level',type=int,default=1,help='Printout for optimization. 2 prints everything in opt.',required=False)
    parser.add_argument('-gsm_print_level',type=int,default=1,help='Printout for gsm. 1 prints ?',required=False)
    parser.add_argument('-linesearch',type=str,default='NoLineSearch',help='default: %(default)s',choices=['NoLineSearch','backtrack'])
    parser.add_argument('-coordinate_type',type=str,default='TRIC',help='Coordinate system (default %(default)s)',choices=['TRIC','DLC','HDLC'])
    parser.add_argument('-ADD_NODE_TOL',type=float,default=0.01,help='Convergence tolerance for adding new node (default: %(default)s)',required=False)
    parser.add_argument('-DQMAG_MAX',type=float,default=0.8,help='Maximum step size in single-ended mode (default: %(default)s)',required=False)
    parser.add_argument('-BDIST_RATIO',type=float,default=0.5,help='Reaction completion convergence in SE modes (default: %(default)s)')
    parser.add_argument('-CONV_TOL',type=float,default=0.0005,help='Convergence tolerance for optimizing nodes (default: %(default)s)',required=False)
    parser.add_argument('-growth_direction',type=int,default=0,help='Direction adding new nodes (default: %(default)s)',choices=[0,1,2])
    parser.add_argument('-reactant_geom_fixed',action='store_true',help='Fix reactant geometry i.e. do not pre-optimize')
    parser.add_argument('-product_geom_fixed',action='store_true',help='Fix product geometry i.e. do not pre-optimize')
    parser.add_argument('-nproc',type=int,default=1,help='Processors for calculation. Python will detect OMP_NUM_THREADS, only use this if you want to force the number of processors')
    parser.add_argument('-charge',type=int,default=0,help='Total system charge (default: %(default)s)')
    parser.add_argument('-max_gsm_iters',type=int,default=100,help='The maximum number of GSM cycles (default: %(default)s)')
    parser.add_argument('-max_opt_steps',type=int,help='The maximum number of node optimizations per GSM cycle (defaults: 3 DE-GSM, 20 SE-GSM)')
    parser.add_argument('-only_climb',action='store_true',help="Only use climbing image to optimize TS")
    parser.add_argument('-no_climb',action='store_true',help="Don't climb to the TS")
    parser.add_argument('-optimize_mesx',action='store_true',help='optimize to the MESX')
    parser.add_argument('-optimize_meci',action='store_true',help='optimize to the MECI')
    parser.add_argument('-restart_file',help='restart file',type=str)
    parser.add_argument('-mp_cores',type=int,default=1,help="Use python multiprocessing to parallelize jobs on a single compute node. Set OMP_NUM_THREADS, ncpus accordingly.")
    parser.add_argument('-dont_analyze_ICs',action='store_false',help="Don't post-print the internal coordinates primitives and values") #defaults to true
    parser.add_argument('-hybrid_coord_idx_file',type=str,default=None,help="A filename containing a list of  indices to use in hybrid coordinates. 0-Based indexed")
    parser.add_argument('-frozen_coord_idx_file',type=str,default=None,help="A filename containing a list of  indices to be frozen. 0-Based indexed")
    parser.add_argument('-conv_Ediff',default=100.,type=float,help='Energy difference convergence of optimization.')
    parser.add_argument('-conv_dE',default=1.,type=float,help='State difference energy convergence')
    parser.add_argument('-conv_gmax',default=100.,type=float,help='Max grad rms threshold')
    parser.add_argument('-DMAX',default=.1,type=float,help='')
    parser.add_argument('-sigma',default=1.,type=float,help='The strength of the difference energy penalty in Penalty_PES')
    parser.add_argument('-prim_idx_file',type=str,help="A filename containing a list of indices to define fragments. 0-Based indexed")
    parser.add_argument('-reparametrize',action='store_true',help='Reparametrize restart string equally along path')
    parser.add_argument('-interp_method',default='DLC',type=str,help='')
    parser.add_argument('-bonds_file',type=str,help="A file which contains the bond indices (0-based)")


    args = parser.parse_args()

    print_msg()

    if args.nproc>1:
        force_num_procs=True
        print("forcing number of processors to be {}!!!".format(args.nproc))
    else:
        force_num_procs=False
    if force_num_procs:
        nproc = args.nproc
    else:
        #nproc = get_nproc()
        try:
            nproc = int(os.environ['OMP_NUM_THREADS'])
        except: 
            nproc = 1
        print(" Using {} processors".format(nproc))

    inpfileq = {
               # LOT
              'lot_inp_file': args.lot_inp_file,
              'xyzfile' : args.xyzfile,
              'EST_Package': args.package,
              'reactant_geom_fixed' : args.reactant_geom_fixed,
              'nproc': args.nproc,
              'states': None,
              
              #PES
              'PES_type': args.pes_type,
              'adiabatic_index': args.adiabatic_index,
              'multiplicity': args.multiplicity,
              'FORCE_FILE': args.FORCE_FILE,
              'RESTRAINT_FILE': args.RESTRAINT_FILE,
              'FORCE': None,
              'RESTRAINTS': None,

              #optimizer
              'optimizer' : args.optimizer,
              'opt_print_level' : args.opt_print_level,
              'linesearch' : args.linesearch,
              'DMAX'    :   args.DMAX,

              #molecule
              'coordinate_type' : args.coordinate_type,
              'hybrid_coord_idx_file' : args.hybrid_coord_idx_file,
              'frozen_coord_idx_file' : args.frozen_coord_idx_file,
              'prim_idx_file' : args.prim_idx_file,

              # GSM
              'gsm_type': args.mode, # SE_GSM, SE_Cross
              'num_nodes' : args.num_nodes,
              'isomers_file': args.isomers,
              'ADD_NODE_TOL': args.ADD_NODE_TOL,
              'CONV_TOL': args.CONV_TOL,
              'conv_Ediff': args.conv_Ediff,
              'conv_dE': args.conv_dE,
              'conv_gmax': args.conv_gmax,
              'BDIST_RATIO':args.BDIST_RATIO,
              'DQMAG_MAX': args.DQMAG_MAX,
              'growth_direction': args.growth_direction,
              'ID':args.ID,
              'product_geom_fixed' : args.product_geom_fixed,
              'gsm_print_level' : args.gsm_print_level,
              'max_gsm_iters' : args.max_gsm_iters,
              'max_opt_steps' : args.max_opt_steps,
              #'use_multiprocessing': args.use_multiprocessing,
              'sigma'   :   args.sigma,
              }

    nifty.printcool_dictionary(inpfileq,title='Parsed GSM Keys : Values')


    #LOT
    nifty.printcool("Build the {} level of theory (LOT) object".format(inpfileq['EST_Package']))
    est_package=importlib.import_module("level_of_theories."+inpfileq['EST_Package'].lower())
    lot_class = getattr(est_package,inpfileq['EST_Package'])

    geoms = manage_xyz.read_xyzs(inpfileq['xyzfile'])
    if args.restart_file:
        geoms = manage_xyz.read_molden_geoms(args.restart_file)

    inpfileq['states'] = [ (int(m),int(s)) for m,s in zip(args.multiplicity,args.adiabatic_index)]
    if inpfileq['PES_type']!="PES":
        assert len(args.adiabatic_index)>1, "need more states"
        assert len(args.multiplicity)>1, "need more spins"
    if args.charge != 0:
        print("Warning: charge is not implemented for all level of theories. Make sure this is correct for your package.")
    if inpfileq['num_nodes'] is None:
        if inpfileq['gsm_type']=="DE_GSM":
            inpfileq['num_nodes']=9
        else:
            inpfileq['num_nodes']=20
    do_coupling = True if inpfileq['PES_type']=="Avg_PES" else False
    coupling_states = [ (int(m),int(s)) for m,s in zip(args.multiplicity,args.adiabatic_index)] if inpfileq['PES_type']=="Avg_PES" else []
    
    lot = lot_class.from_options(
            ID = inpfileq['ID'],
            lot_inp_file=inpfileq['lot_inp_file'],
            states=inpfileq['states'],
            gradient_states=inpfileq['states'],
            coupling_states=coupling_states,
            geom=geoms[0],
            nproc=nproc,
            charge=args.charge,
            do_coupling=do_coupling,
            )

    #PES
    if inpfileq['gsm_type'] == "SE_Cross":
        if inpfileq['PES_type']!="Penalty_PES":
            print(" setting PES type to Penalty")
            inpfileq['PES_type']="Penalty_PES"
    if args.optimize_mesx or args.optimize_meci  or inpfileq['gsm_type']=="SE_Cross":
        assert inpfileq['PES_type'] == "Penalty_PES", "Need penalty pes for optimizing MESX/MECI"
    if inpfileq['FORCE_FILE']:
        inpfileq['FORCE']=[]
        with open(inpfileq['FORCE_FILE'],'r') as f:
            tmp = filter(None, (line.rstrip() for line in f))
            lines=[]
            for line in tmp:
                lines.append(line)
        for line in lines:
            force=[]
            for i,elem in enumerate(line.split()):
                if i==0 or i==1:
                    force.append(int(elem))
                else:
                    force.append(float(elem))
            inpfileq['FORCE'].append(tuple(force))
        print(inpfileq['FORCE'])
    if inpfileq['RESTRAINT_FILE']:
        inpfileq['RESTRAINTS']=[]
        with open(inpfileq['RESTRAINT_FILE'],'r') as f:
            tmp = filter(None, (line.rstrip() for line in f))
            lines=[]
            for line in tmp:
                lines.append(line)
        for line in lines:
            restraint=[]
            for i,elem in enumerate(line.split()):
                if i==0:
                    restraint.append(int(elem))
                else:
                    restraint.append(float(elem))
            inpfileq['RESTRAINTS'].append(tuple(restraint))
        print(inpfileq['RESTRAINTS'])

    nifty.printcool("Building the {} objects".format(inpfileq['PES_type']))
    pes_class = getattr(sys.modules[__name__], inpfileq['PES_type'])
    if inpfileq['PES_type']=='PES':
        pes = pes_class.from_options(
                lot=lot,
                ad_idx=inpfileq['adiabatic_index'][0],
                multiplicity=inpfileq['multiplicity'][0],
                FORCE=inpfileq['FORCE'],
                RESTRAINTS=inpfileq['RESTRAINTS'],
                )
    else:
        pes1 = PES.from_options(
                lot=lot,multiplicity=inpfileq['states'][0][0],
                ad_idx=inpfileq['states'][0][1],
                FORCE=inpfileq['FORCE'],
                RESTRAINTS=inpfileq['RESTRAINTS'],
                )
        pes2 = PES.from_options(
                lot=lot,
                multiplicity=inpfileq['states'][1][0],
                ad_idx=inpfileq['states'][1][1],
                FORCE=inpfileq['FORCE'],
                RESTRAINTS=inpfileq['RESTRAINTS'],
                )
        if inpfileq['PES_type']=="Avg_PES":
            pes = pes_class(PES1=pes1,PES2=pes2,lot=lot)
        elif inpfileq['PES_type']=="Penalty_PES": 
            pes = pes_class(PES1=pes1,PES2=pes2,lot=lot,sigma=inpfileq['sigma'])
        else:
            raise NotImplementedError

    # Molecule
    nifty.printcool("Building the reactant object with {}".format(inpfileq['coordinate_type']))
    Form_Hessian = True if inpfileq['optimizer']=='eigenvector_follow' else False
    #form_primitives = True if inpfileq['gsm_type']!='DE_GSM' else False

    # hybrid coordinates
    if inpfileq['hybrid_coord_idx_file'] is not None:
        nifty.printcool(" Using Hybrid COORDINATES :)")
        assert inpfileq['coordinate_type']=="TRIC", "hybrid indices won't work (currently) with other coordinate systems"
        with open(inpfileq['hybrid_coord_idx_file']) as f:
            hybrid_indices = f.read().splitlines()
        hybrid_indices = [int(x) for x in hybrid_indices]
    else:
        hybrid_indices = None
    if inpfileq['frozen_coord_idx_file'] is not None:
        with open(inpfileq['frozen_coord_idx_file']) as f:
            frozen_indices = f.read().splitlines()
        frozen_indices = [int(x) for x in frozen_indices]
    else:
        frozen_indices = None


    # prim internal coordinates
    # The start and stop indexes of the primitive internal region, this defines the "fragments" so no large molecule is built        
    if inpfileq['prim_idx_file'] is not None:
        nifty.printcool(" Defining primitive internal region :)")
        assert inpfileq['coordinate_type']=="TRIC", "won't work (currently) with other coordinate systems"
        prim_indices = np.loadtxt(inpfileq['prim_idx_file'])
        if prim_indices.ndim==2:
            prim_indices = [ (int(prim_indices[i,0]), int(prim_indices[i,1])-1) for i in range(len(prim_indices))]
        elif prim_indices.ndim==1:
            prim_indices = [ (int(prim_indices[0]), int(prim_indices[1])-1) ]

        print(prim_indices)
        #with open(inpfileq['prim_idx_file']) as f:
        #    prim_indices = f.read().splitlines()
        #prim_indices = [int(x) for x in prim_indices]
    else:
        prim_indices = None



    # Build the topology
    nifty.printcool("Building the topology")
    atom_symbols  = manage_xyz.get_atoms(geoms[0])
    ELEMENT_TABLE = elements.ElementData()
    atoms = [ELEMENT_TABLE.from_symbol(atom) for atom in atom_symbols]
    xyz1 = manage_xyz.xyz_to_np(geoms[0])
    top1 = Topology.build_topology(
            xyz1,
            atoms,
            hybrid_indices=hybrid_indices,
            prim_idx_start_stop=prim_indices,
            bondlistfile=args.bonds_file,
            )

    if inpfileq['gsm_type'] == 'DE_GSM':
        # find union bonds
        xyz2 = manage_xyz.xyz_to_np(geoms[-1])
        top2 = Topology.build_topology(
                xyz2,
                atoms,
                hybrid_indices=hybrid_indices,
                prim_idx_start_stop=prim_indices,
                )

        # Add bonds to top1 that are present in top2
        # It's not clear if we should form the topology so the bonds
        # are the same since this might affect the Primitives of the xyz1 (slightly)
        # Later we stil need to form the union of bonds, angles and torsions
        # However, I think this is important, the way its formulated, for identifiyin 
        # the number of fragments and blocks, which is used in hybrid TRIC. 
        for bond in top2.edges():
            if bond in top1.edges:
                pass
            elif (bond[1],bond[0]) in top1.edges():
                pass
            else:
                print(" Adding bond {} to top1".format(bond))
                if bond[0]>bond[1]:
                    top1.add_edge(bond[0],bond[1])
                else:
                    top1.add_edge(bond[1],bond[0])
    elif inpfileq['gsm_type'] == 'SE_GSM' or inpfileq['gsm_type']=='SE_Cross':
        driving_coordinates = read_isomers_file(inpfileq['isomers_file'])

        driving_coord_prims=[]
        for dc in driving_coordinates:
            prim = get_driving_coord_prim(dc)
            if prim is not None:
                driving_coord_prims.append(prim)

        for prim in driving_coord_prims:
            if type(prim)==Distance:
                bond = (prim.atoms[0],prim.atoms[1])
                if bond in top1.edges:
                    pass
                elif (bond[1],bond[0]) in top1.edges():
                    pass
                else:
                    print(" Adding bond {} to top1".format(bond))
                    top1.add_edge(bond[0],bond[1])

    nifty.printcool("Building Primitive Internal Coordinates")
    connect=False
    addtr = False
    addcart=False
    if inpfileq['coordinate_type']=="DLC":
        connect=True
    elif inpfileq['coordinate_type']=="TRIC":
        addtr=True
    elif inpfileq['coordinate_type']=="HDLC":
        addcart=True
    p1 = PrimitiveInternalCoordinates.from_options(
            xyz=xyz1,
            atoms=atoms,
            connect=connect,
            addtr=addtr,
            addcart=addcart,
            topology=top1,
            )

    if inpfileq['gsm_type'] == 'DE_GSM':
        nifty.printcool("Building Primitive Internal Coordinates 2")
        p2 = PrimitiveInternalCoordinates.from_options(
                xyz=xyz2,
                atoms=atoms,
                addtr = addtr,
                addcart=addcart,
                connect=connect,
                topology=top1,  # Use the topology of 1 because we fixed it above
                ) 
        nifty.printcool("Forming Union of Primitives")
        # Form the union of primitives
        p1.add_union_primitives(p2)

        print("check {}".format(len(p1.Internals)))
    elif inpfileq['gsm_type'] == 'SE_GSM' or inpfileq['gsm_type']=='SE_Cross':
        for dc in driving_coord_prims:
            if type(dc)!=Distance: # Already handled in topology
                if dc not in p1.Internals:
                    print("Adding driving coord prim {} to Internals".format(dc))
                    p1.append_prim_to_block(dc)

    nifty.printcool("Building Delocalized Internal Coordinates")
    coord_obj1 = DelocalizedInternalCoordinates.from_options(
            xyz=xyz1,
            atoms=atoms,
            addtr = addtr,
            addcart=addcart,
            connect=connect,
            primitives=p1,
            ) 
    if inpfileq['gsm_type'] == 'DE_GSM':
        # TMP 
        pass


    nifty.printcool("Building the reactant")
    reactant = Molecule.from_options(
            geom=geoms[0],
            PES=pes,
            coord_obj = coord_obj1,
            Form_Hessian=Form_Hessian,
            frozen_atoms=frozen_indices,
            )

    if inpfileq['gsm_type']=='DE_GSM':
        nifty.printcool("Building the product object")
        xyz2 = manage_xyz.xyz_to_np(geoms[-1])
        product = Molecule.copy_from_options(
                reactant,
                xyz=xyz2,
                new_node_id = inpfileq['num_nodes']-1,
                copy_wavefunction=False,
                )
   
    # optimizer
    nifty.printcool("Building the Optimizer object")
    update_hess_in_bg = True
    if args.only_climb or inpfileq['optimizer']=="lbfgs":
        update_hess_in_bg = False
    opt_class = getattr(sys.modules[__name__], inpfileq['optimizer'])
    optimizer = opt_class.from_options(
            print_level=inpfileq['opt_print_level'],
            Linesearch=inpfileq['linesearch'],
            update_hess_in_bg = update_hess_in_bg,
            conv_Ediff = inpfileq['conv_Ediff'],
            conv_dE = inpfileq['conv_dE'],
            conv_gmax = inpfileq['conv_gmax'],
            DMAX = inpfileq['DMAX'],
            opt_climb = True if args.only_climb else False,
            )

    # GSM
    nifty.printcool("Building the GSM object")
    gsm_class = getattr(sys.modules[__name__], inpfileq['gsm_type'])
    if inpfileq['gsm_type']=="DE_GSM":
        gsm = gsm_class.from_options(
                reactant=reactant,
                product=product,
                nnodes=inpfileq['num_nodes'],
                CONV_TOL=inpfileq['CONV_TOL'],
                CONV_gmax=inpfileq['conv_gmax'],
                CONV_Ediff=inpfileq['conv_Ediff'],
                CONV_dE=inpfileq['conv_dE'],
                ADD_NODE_TOL=inpfileq['ADD_NODE_TOL'],
                growth_direction=inpfileq['growth_direction'],
                optimizer=optimizer,
                ID=inpfileq['ID'],
                print_level=inpfileq['gsm_print_level'],
                mp_cores=args.mp_cores,
                interp_method = args.interp_method,
                )
    else:
        gsm = gsm_class.from_options(
                reactant=reactant,
                nnodes=inpfileq['num_nodes'],
                DQMAG_MAX=inpfileq['DQMAG_MAX'],
                BDIST_RATIO=inpfileq['BDIST_RATIO'],
                CONV_TOL=inpfileq['CONV_TOL'],
                ADD_NODE_TOL=inpfileq['ADD_NODE_TOL'],
                optimizer=optimizer,
                print_level=inpfileq['gsm_print_level'],
                driving_coords=driving_coordinates,
                ID=inpfileq['ID'],
                mp_cores=args.mp_cores,
                interp_method = args.interp_method,
                )



    if args.only_drive:
        for i in range(gsm.nnodes-1):
            try:
                gsm.add_GSM_nodeR()
            except:
                break
        geoms=[]
        for node in gsm.nodes:
            if node is not None:
                geoms.append(node.geometry)
            else:
                break
        manage_xyz.write_xyzs('interpolated.xyz',geoms)
        return


    # For seam calculation
    if inpfileq['gsm_type']!='SE_Cross' and (inpfileq['PES_type'] =="Avg_PES" or inpfileq['PES_type']=="Penalty_PES"):
        optimizer.opt_cross = True

    if not inpfileq['reactant_geom_fixed'] and inpfileq['gsm_type']!='SE_Cross':
        path=os.path.join(os.getcwd(),'scratch/{:03}/{}/'.format(args.ID,0))
        nifty.printcool("REACTANT GEOMETRY NOT FIXED!!! OPTIMIZING")
        optimizer.optimize(
           molecule = reactant,
           refE = reactant.energy,
           opt_steps=100,
           path=path
           )

    if not inpfileq['product_geom_fixed'] and inpfileq['gsm_type']=='DE_GSM':
        path=os.path.join(os.getcwd(),'scratch/{:03}/{}/'.format(args.ID,args.num_nodes-1))
        nifty.printcool("PRODUCT GEOMETRY NOT FIXED!!! OPTIMIZING")
        optimizer.optimize(
           molecule = product,
           refE = reactant.energy,
           opt_steps=100,
           path=path
           )

    rtype=2
    if args.only_climb:
        rtype=1
    elif args.no_climb:
        rtype=0
    elif args.optimize_meci:
        rtype=0
    elif args.optimize_mesx:
        rtype=1
    elif inpfileq['gsm_type']=="SE_Cross":
        rtype=1

    if inpfileq['max_opt_steps'] is None:
        if inpfileq['gsm_type']=="DE_GSM":
            inpfileq['max_opt_steps']=3
        else:
            inpfileq['max_opt_steps']=20
   
    if args.restart_file is not None:
        gsm.setup_from_geometries(geoms,reparametrize=args.reparametrize)
    gsm.go_gsm(inpfileq['max_gsm_iters'],inpfileq['max_opt_steps'],rtype)
    if inpfileq['gsm_type']=='SE_Cross':
        post_processing(
                gsm,
                analyze_ICs=args.dont_analyze_ICs,
                have_TS=False,
                )
        manage_xyz.write_xyz(f'meci_{gsm.ID}.xyz',gsm.nodes[gsm.nR].geometry)

        if not gsm.end_early:
            manage_xyz.write_xyz(f'TSnode_{gsm.ID}.xyz',gsm.nodes[gsm.TSnode].geometry)
    else:
        post_processing(
                gsm,
                analyze_ICs=args.dont_analyze_ICs,
                have_TS=True,
                )
        manage_xyz.write_xyz(f'TSnode_{gsm.ID}.xyz',gsm.nodes[gsm.TSnode].geometry)

    cleanup_scratch(gsm.ID)

    return
Esempio n. 4
0
    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 = Base_Method.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.growth_iters(iters=max_iters, maxopt=opt_steps, nconstraints=1)
        print(' SE_Cross growth phase over')
        print(' Warning last node still not fully optimized')

        if True:
            # 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 = 0.01
            ictan, _ = Base_Method.tangent(self.nodes[self.nR - 1],
                                           self.nodes[self.nR - 2])
            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=50,
                ictan=ictan,
            )

        self.optimizer[self.nR].opt_cross = True
        if rtype == 0:
            # MECI optimization
            self.write_xyz_files(iters=1, base="after_penalty", nconstraints=1)
            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
            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='MECI',
                opt_steps=100,
            )
            self.write_xyz_files(iters=1, base="grown_string", nconstraints=1)
        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))
            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,
            )
            self.write_xyz_files(iters=1, base="grown_string", nconstraints=1)
Esempio n. 5
0
    def add_node(nodeR, nodeP, stepsize, node_id, **kwargs):
        '''
        Add a node between  nodeR and nodeP or if nodeP is none use driving coordinate to add new node
        '''

        # get driving coord
        driving_coords = kwargs.get('driving_coords', None)
        DQMAG_MAX = kwargs.get('DQMAG_MAX', 0.8)
        DQMAG_MIN = kwargs.get('DQMAG_MIN', 0.2)

        if nodeP is None:

            if driving_coords is None:
                raise RuntimeError(
                    "You didn't supply a driving coordinate and product node is None!"
                )

            BDISTMIN = 0.05
            ictan, bdist = GSM.get_tangent(nodeR,
                                           None,
                                           driving_coords=driving_coords)

            if bdist < BDISTMIN:
                print("bdist too small %.3f" % bdist)
                return None
            new_node = Molecule.copy_from_options(nodeR, new_node_id=node_id)
            new_node.update_coordinate_basis(constraints=ictan)
            constraint = new_node.constraints[:, 0]
            sign = -1.

            dqmag_scale = 1.5
            minmax = DQMAG_MAX - DQMAG_MIN
            a = bdist / dqmag_scale
            if a > 1.:
                a = 1.
            dqmag = sign * (DQMAG_MIN + minmax * a)
            if dqmag > DQMAG_MAX:
                dqmag = DQMAG_MAX
            print(" dqmag: %4.3f from bdist: %4.3f" % (dqmag, bdist))

            dq0 = dqmag * constraint
            print(" dq0[constraint]: %1.3f" % dqmag)

            new_node.update_xyz(dq0)
            new_node.bdist = bdist

        else:
            ictan, _ = GSM.get_tangent(nodeR, nodeP)
            nodeR.update_coordinate_basis(constraints=ictan)
            constraint = nodeR.constraints[:, 0]
            dqmag = np.linalg.norm(ictan)
            print(" dqmag: %1.3f" % dqmag)
            # sign=-1
            sign = 1.
            dqmag *= (sign * stepsize)
            print(" scaled dqmag: %1.3f" % dqmag)

            dq0 = dqmag * constraint
            old_xyz = nodeR.xyz.copy()
            new_xyz = nodeR.coord_obj.newCartesian(old_xyz, dq0)
            new_node = Molecule.copy_from_options(MoleculeA=nodeR,
                                                  xyz=new_xyz,
                                                  new_node_id=node_id)

        return new_node
Esempio n. 6
0
    def __init__(
        self,
        options,
    ):
        """ Constructor """
        self.options = options

        os.system('mkdir -p scratch')

        # Cache attributes
        self.nnodes = self.options['nnodes']
        self.nodes = [None] * self.nnodes
        self.nodes[0] = self.options['reactant']
        self.nodes[-1] = self.options['product']
        self.driving_coords = self.options['driving_coords']
        self.growth_direction = self.options['growth_direction']
        self.isRestarted = False
        self.DQMAG_MAX = self.options['DQMAG_MAX']
        self.DQMAG_MIN = self.options['DQMAG_MIN']
        self.BDIST_RATIO = self.options['BDIST_RATIO']
        self.ID = self.options['ID']
        self.optimizer = []
        self.interp_method = self.options['interp_method']
        self.CONV_TOL = self.options['CONV_TOL']
        self.noise = self.options['noise']
        self.mp_cores = self.options['mp_cores']
        self.xyz_writer = self.options['xyz_writer']

        optimizer = options['optimizer']
        for count in range(self.nnodes):
            self.optimizer.append(optimizer.__class__(
                optimizer.options.copy()))
        self.print_level = options['print_level']

        # Set initial values
        self.current_nnodes = 2
        self.nR = 1
        self.nP = 1
        self.climb = False
        self.find = False
        self.ts_exsteps = 3  # multiplier for ts node
        self.n0 = 1  # something to do with added nodes? "first node along current block"
        self.end_early = False
        self.tscontinue = True  # whether to continue with TS opt or not
        self.found_ts = False
        self.rn3m6 = np.sqrt(3. * self.nodes[0].natoms - 6.)
        self.gaddmax = self.options[
            'ADD_NODE_TOL']  # self.options['ADD_NODE_TOL']/self.rn3m6;
        print(" gaddmax:", self.gaddmax)
        self.ictan = [None] * self.nnodes
        self.active = [False] * self.nnodes
        self.climber = False  # is this string a climber?
        self.finder = False  # is this string a finder?
        self.done_growing = False
        self.nclimb = 0
        self.nhessreset = 10  # are these used??? TODO
        self.hessrcount = 0  # are these used?!  TODO
        self.hess_counter = 0  # it is probably good to reset the hessian
        self.newclimbscale = 2.
        self.TS_E_0 = None
        self.dE_iter = 100.  # change in max TS node

        self.nopt_intermediate = 0  # might be a duplicate of endearly_counter
        self.flag_intermediate = False
        self.endearly_counter = 0  # Find the intermediate x time
        self.pot_min = []
        self.ran_out = False  # if it ran out of iterations

        self.newic = Molecule.copy_from_options(
            self.nodes[0]
        )  # newic object is used for coordinate transformations
Esempio n. 7
0
if __name__ == '__main__':
    from .qchem import QChem
    from .pes import PES
    from .dlc_new import DelocalizedInternalCoordinates
    from .eigenvector_follow import eigenvector_follow
    from ._linesearch import backtrack, NoLineSearch
    from .molecule import Molecule

    basis = '6-31G'
    nproc = 8
    functional = 'B3LYP'
    filepath1 = "examples/tests/butadiene_ethene.xyz"
    lot1 = QChem.from_options(states=[(1, 0)],
                              charge=0,
                              basis=basis,
                              functional=functional,
                              nproc=nproc,
                              fnm=filepath1)
    pes1 = PES.from_options(lot=lot1, ad_idx=0, multiplicity=1)
    M1 = Molecule.from_options(fnm=filepath1, PES=pes1, coordinate_type="DLC")
    optimizer = eigenvector_follow.from_options(
        print_level=1
    )  #default parameters fine here/opt_type will get set by GSM

    gsm = SE_GSM.from_options(reactant=M1,
                              nnodes=20,
                              driving_coords=[("ADD", 6, 4), ("ADD", 5, 1)],
                              optimizer=optimizer,
                              print_level=1)
    gsm.go_gsm()
Esempio n. 8
0
    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,
            )

        write_molden_geoms('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,
            )
        write_molden_geoms('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()
            write_molden_geoms('grown_string1_{:03}.xyz'.format(self.ID),
                               self.geometries, self.energies, self.gradrmss,
                               self.dEs)

            deltaE = energies[-1] - energies[0]
            if deltaE > 20:
                print(
                    f" MECI energy is too high {deltaE}. Don't try to optimize pathway"
                )
                print("Exiting early")
                self.end_early = True
            else:
                print(f" deltaE s1-minimum and MECI {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
Esempio n. 9
0
    p1.add_union_primitives(p2)

    coord_obj1 = DelocalizedInternalCoordinates.from_options(
        xyz=xyz1,
        atoms=atoms,
        addtr=addtr,
        addcart=addcart,
        connect=connect,
        primitives=p1,
    )

    optimizer = eigenvector_follow.from_options()

    reactant = Molecule.from_options(
        geom=geoms[0],
        PES=pes,
        coord_obj=coord_obj1,
        Form_Hessian=True,
    )
    product = Molecule.copy_from_options(
        reactant,
        xyz=xyz2,
        new_node_id=len(geoms) - 1,
        copy_wavefunction=False,
    )

    gsm = DE_GSM.from_options(
        reactant=reactant,
        product=product,
        nnodes=len(geoms),
        optimizer=optimizer,
    )