def Return_U_from_Aniso_Expand(inputs, new_crystal_matrix, coordinate_file,
                               Parameter_file, Program, output_file_name,
                               molecules_in_coord):

    # Converting the crystal matrix parameter to lattice parameters
    new_lattice_parameters = Ex.crystal_matrix_to_lattice_parameters(
        Ex.array_to_triangle_crystal_matrix(new_crystal_matrix))

    # Determining the file ending of the coordinate file
    file_ending = psf.assign_coordinate_file_ending(Program)

    # Determine the input coordinate files lattice parameters
    old_lattice_parameters = psf.Lattice_parameters(Program, coordinate_file)

    # Expand input coordinate file using the new lattice parameters
    Ex.Expand_Structure(inputs,
                        coordinate_file,
                        'lattice_parameters',
                        output_file_name,
                        dlattice_parameters=new_lattice_parameters[:6] -
                        old_lattice_parameters)
    # Computing the potential energy of the new expanded structure
    U = psf.Potential_energy(
        output_file_name + file_ending, Program,
        Parameter_file=Parameter_file) / molecules_in_coord
    return U
Beispiel #2
0
def anisotropic_gradient_settings_1D(inputs, dC_dLambda):
    # Determining the file ending based on the program
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    # Setting the energy cutoff
    cutoff = program_cutoff(inputs.program)

    # Lambda step sizes to take
    steps = np.array([5e-02, 1e-01, 5e-01, 1e01, 5e01, 1e02, 5e02, 1e03, 5e03])

    # Number of total step sizes
    n_steps = len(steps)

    # Potential energy of input file and a place to store the expanded structures potential energy
    U_0 = psf.Potential_energy(inputs.coordinate_file, inputs.program, Parameter_file=inputs.tinker_parameter_file) / \
          inputs.number_of_molecules
    U = np.zeros(n_steps)

    for i in range(n_steps):
        dlattice_matrix = Ex.array_to_triangle_crystal_matrix(steps[i] *
                                                              dC_dLambda)
        Ex.Expand_Structure(inputs,
                            inputs.coordinate_file,
                            'crystal_matrix',
                            inputs.output,
                            dcrystal_matrix=dlattice_matrix)
        U[i] = psf.Potential_energy(inputs.output + file_ending,
                                    inputs.program,
                                    Parameter_file=inputs.tinker_parameter_file
                                    ) / inputs.number_of_molecules
        subprocess.call(['rm', inputs.output + file_ending])
        if (U[i] - U_0) > cutoff:
            LocGrd_dLambda = steps[i]
            end_plot = i
            break

    # Plotting the results
    plt.plot(np.log10(steps[:end_plot + 1]),
             U[:end_plot + 1] - U_0,
             linestyle='--',
             marker='o')
    plt.xlabel('$\log({d\lambda})$', fontsize=22)
    plt.ylabel('$\Delta U$ [kcal/mol]', fontsize=22)
    plt.ylim((0., 2 * cutoff))
    plt.axhline(y=cutoff, c='grey', linestyle='--')
    plt.tight_layout()
    plt.savefig(inputs.output + '_LocGrd_Lambda_FracStep.pdf')
    plt.close()
    print('dLambda used: ', LocGrd_dLambda)
    # returning the value of dV
    return LocGrd_dLambda
def constrained_minimization(inputs,
                             Coordinate_file,
                             Program,
                             Parameter_file=''):
    # Determining the file ending of the coordinate file
    file_ending = psf.assign_coordinate_file_ending(Program)

    # Determining the lattice parameters and volume of the input coordinate file
    lp_0 = psf.Lattice_parameters(Program, Coordinate_file)
    V0 = Pr.Volume(lattice_parameters=lp_0)
    cm_0 = Ex.triangle_crystal_matrix_to_array(
        Ex.Lattice_parameters_to_Crystal_matrix(lp_0))

    bnds = np.matrix([[cm_0[0] - cm_0[0] * 0.2, cm_0[0] + cm_0[0] * 0.2],
                      [cm_0[1] - cm_0[1] * 0.2, cm_0[1] + cm_0[1] * 0.2],
                      [cm_0[2] - cm_0[2] * 0.2, cm_0[2] + cm_0[2] * 0.2],
                      [cm_0[3] - cm_0[0] * 0.2, cm_0[3] + cm_0[0] * 0.2],
                      [cm_0[4] - cm_0[0] * 0.2, cm_0[4] + cm_0[0] * 0.2],
                      [cm_0[5] - cm_0[0] * 0.2, cm_0[5] + cm_0[0] * 0.2]])

    output = scipy.optimize.minimize(
        Return_U_from_Aniso_Expand,
        cm_0, (inputs, Coordinate_file, Parameter_file, Program,
               'constV_minimize', 4),
        method='SLSQP',
        constraints=({
            'type':
            'eq',
            'fun':
            lambda cm: np.linalg.det(Ex.array_to_triangle_crystal_matrix(cm)) -
            V0
        }),
        bounds=bnds,
        tol=1e-08,
        options={
            'ftol': float(1e-6),
            'disp': True,
            'eps': float(1e-4)
        })

    dlattice_parameters = Ex.crystal_matrix_to_lattice_parameters(
        Ex.array_to_triangle_crystal_matrix(output.x)) - lp_0
    Ex.Expand_Structure(inputs,
                        Coordinate_file,
                        'lattice_parameters',
                        'constV_minimize',
                        dlattice_parameters=dlattice_parameters)
    subprocess.call(['mv', 'constV_minimize' + file_ending, Coordinate_file])
Beispiel #4
0
def Anisotropic_Local_Gradient_1D(inputs, coordinate_file, temperature,
                                  LocGrd_dLambda, dC_dLambda,
                                  **keyword_parameters):
    # Determining the file ending of the coordinate files
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    if temperature == 0.:
        # If temperature is zero, we assume that the local gradient is the same at 0.1K
        temperature = 1e-5

    # Retrieving the wavenumbers of the initial structure
    wavenumber_keywords = {
        'Gruneisen': keyword_parameters['Gruneisen'],
        'Wavenumber_Reference': keyword_parameters['Wavenumber_Reference'],
        'ref_crystal_matrix': keyword_parameters['ref_crystal_matrix']
    }
    wavenumbers = Wvn.Call_Wavenumbers(inputs,
                                       Coordinate_file=coordinate_file,
                                       **wavenumber_keywords)

    G = np.zeros(3)
    U = np.zeros(3)
    Av = np.zeros(3)
    S = np.zeros(3)

    G[1], U[1], Av[1] = Pr.Gibbs_Free_Energy(
        temperature,
        inputs.pressure,
        inputs.program,
        wavenumbers,
        coordinate_file,
        inputs.statistical_mechanics,
        inputs.number_of_molecules,
        Parameter_file=inputs.tinker_parameter_file)

    crystal_matrix_array = LocGrd_dLambda * dC_dLambda
    # Straining the input structure by the current diagonal
    delta = ('p', 'm')
    for d in delta:
        if d == 'm':
            cm_factor = -1.0
            out_factor = 0
        else:
            cm_factor = 1.0
            out_factor = 2

        # Expand the crystal structure
        Expand_Structure(inputs,
                         coordinate_file,
                         'crystal_matrix',
                         d,
                         dcrystal_matrix=array_to_triangle_crystal_matrix(
                             cm_factor * crystal_matrix_array))

        # Compute the wavenumbers
        wavenumbers_hold = Wvn.Call_Wavenumbers(inputs,
                                                Coordinate_file=d +
                                                file_ending,
                                                **wavenumber_keywords)

        # Compute the energy
        G[out_factor], U[out_factor], Av[out_factor] = \
            Pr.Gibbs_Free_Energy(temperature, inputs.pressure, inputs.program, wavenumbers_hold, d + file_ending,
                                 inputs.statistical_mechanics, inputs.number_of_molecules,
                                 Parameter_file=inputs.tinker_parameter_file)

        # Compute the entropy
        S[out_factor] = Pr.Vibrational_Entropy(temperature, wavenumbers_hold, inputs.statistical_mechanics) / \
                        inputs.number_of_molecules

        # Remove excess files
        subprocess.call(['rm', d + file_ending])

    # Computing the numerical gradient of entropy wrt lambda
    dS_dLambda = (S[2] - S[0]) / (2 * LocGrd_dLambda)

    # Calculating the finite difference of dG/deta for forward, central, and backwards
    dG_dLambda = np.zeros(3)
    dG_dLambda[0] = (G[1] - G[0]) / (LocGrd_dLambda)
    dG_dLambda[1] = (G[2] - G[0]) / (2 * LocGrd_dLambda)
    dG_dLambda[2] = (G[2] - G[1]) / (LocGrd_dLambda)

    # Computing the second derivative Gibbs free energy as a function of strain using a finite difference approach
    ddG_ddLambda = (G[2] - 2 * G[1] + G[0]) / (LocGrd_dLambda**2)

    # Computing the expansion of lambda
    dLambda = dS_dLambda / ddG_ddLambda

    # Determining if the system is still at a free energy minima wrt lambda
    left_minimum = NO.aniso_gradient_1D(dG_dLambda, ddG_ddLambda, dS_dLambda,
                                        dLambda)
    return dLambda, wavenumbers, left_minimum
Beispiel #5
0
def Anisotropic_Local_Gradient(inputs,
                               coordinate_file,
                               temperature,
                               LocGrd_dC,
                               zeta=-1.,
                               **keyword_parameters):
    """
    This function calculates the local gradient of anisotropic expansion for a given coordinate file

    :param Coordinate_file: file containing lattice parameters (and coordinates)
    :param Program: 'Tinker' Tinker molecular modeling
                    'Test' Test case
    :param Temperature: in Kelvin
    :param Pressure: in atm
    :param LocGrd_LatParam_FracStep: 
    :param molecules_in_coord: number of molecules in Coordinate_file
    :param Statistical_mechanics: 'Classical' Classical mechanics
                                  'Quantum' Quantum mechanics
    :param Method: 'GaQ' Gradient anisotropic QHA
                   'GaQg' Gradient anisotropic QHA with Gruneisen Parameter
    :param Hessian_number: 73 Hessians to calculate the complete anistropic gradient
                           25 for d**2G_dUdU only calculating the diagonals and off-diags. of the upper left 3x3 matrix
                           19 for d**2G_dUdU only calculating the upper left 3x3 matrix
                           13 for d**2G_dUdU only calculating the diagonals
                           7  for d**2G_dUdU only calculating the upper left 3x3 matrix daigonals
    :param keyword_parameters: Parameter_file, Gruneisen, Wavenumber_reference, Volume_reference

    Optional Parameters
    Parameter_file: program specific file containing force field parameters
    Gruneisen: isotropic Gruneisen parameter
    Wavenumber_reference: reference wavenumbers for the Gruneisen parameter
    Crystal_matrix_Reference
    """
    min_numerical_crystal_matrix = 1.0e-7

    # Determining the file ending of the coordinate files
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    # Preparing the matrix with each entry as d**2G/(dC*dC)
    ddG_ddC = np.zeros((6, 6))

    # Preparing the vector with each entry as dS/dC and ddG/dCd(zeta or T)
    # dS/dC will be more accurate than ddG/dCdT
    # ddG/dCd(zeta) must be used for the EZP expansion
    dS_dC = np.zeros(6)
    ddG_dCdzeta = np.zeros(6)

    # A place to save potential and helmholtz vibrational energies to output
    U = np.zeros((6, 2))
    Av = np.zeros((6, 2))
    PV = np.zeros((6, 2))

    # dictionaries for saved intermediate data
    dG_dict = dict()
    wavenumbers_dict = dict()

    # Making array for dG/dC
    dG_dC = np.zeros((6, 3))

    if temperature == 0. and zeta == -1. and inputs.statistical_mechanics == 'Classical':
        # If temperature is zero, we assume that the local gradient is the same at 0.1K
        temperature = 0.0001
    elif temperature == 0. and zeta == -1. and inputs.statistical_mechanics == 'Quantum':
        # If temperature is zero, we assume that the local gradient is the same at 0.1K
        temperature = 0.1
    elif zeta == 0.:
        zeta = 0.0001

    # Modified anisotropic Local Gradient
    if (inputs.anisotropic_type
            == '6D') or ((inputs.anisotropic_type == '1D')
                         and not os.path.isfile(inputs.output + '_dC_' +
                                                inputs.method + '.npy')):
        diag_limit = 6
        off_diag_limit = 6
    elif inputs.anisotropic_type == '3D':
        diag_limit = 3
        off_diag_limit = 3
    else:
        print("Aniso_LocGrad_Type = ", inputs.anisotropic_type,
              " is not a valid option.")
        sys.exit()

    # Retrieving the wavenumbers of the initial structure
    wavenumbers = Wvn.Call_Wavenumbers(
        inputs,
        Coordinate_file=coordinate_file,
        ref_crystal_matrix=keyword_parameters['ref_crystal_matrix'],
        Gruneisen=keyword_parameters['Gruneisen'],
        Wavenumber_Reference=keyword_parameters['Wavenumber_Reference'])

    G_hold, U_0, Av_0 = Pr.Gibbs_Free_Energy(
        temperature,
        inputs.pressure,
        inputs.program,
        wavenumbers,
        coordinate_file,
        inputs.statistical_mechanics,
        inputs.number_of_molecules,
        Parameter_file=inputs.tinker_parameter_file)

    PV_0 = G_hold - U_0 - Av_0
    dG_dict['0'] = U_0 + Av_0 * np.absolute(zeta) + PV_0

    Av_0 *= zeta
    for i in range(diag_limit):
        # Calculating the diagonals of ddG_ddeta and the vector dS_deta

        # Setting the additional strain to the input structure for the current diagonal element
        if LocGrd_dC[i] < min_numerical_crystal_matrix:
            continue  # they remain zero, we don't bother to calculate them

        cm_array = np.zeros(6)
        cm_array[i] = LocGrd_dC[i]

        # Straining the input structure by the current diagonal
        delta1 = ('p', 'm')
        for d in delta1:
            if d == 'm':
                cm_factor = -1.0
                out_factor = 0
            else:
                cm_factor = 1.0
                out_factor = 1
            Expand_Structure(inputs,
                             coordinate_file,
                             'crystal_matrix',
                             d,
                             dcrystal_matrix=array_to_triangle_crystal_matrix(
                                 cm_factor * cm_array))

            wavenumbers_dict[d] = Wvn.Call_Wavenumbers(
                inputs,
                Coordinate_file=d + file_ending,
                ref_crystal_matrix=keyword_parameters['ref_crystal_matrix'],
                Gruneisen=keyword_parameters['Gruneisen'],
                Wavenumber_Reference=keyword_parameters['Wavenumber_Reference']
            )
            G_hold, U[i, out_factor], Av[i, out_factor] = \
                Pr.Gibbs_Free_Energy(temperature, inputs.pressure, inputs.program, wavenumbers_dict[d], d + file_ending,
                                     inputs.statistical_mechanics, inputs.number_of_molecules,
                                     Parameter_file=inputs.tinker_parameter_file)
            PV[i, out_factor] = G_hold - U[i, out_factor] - Av[i, out_factor]
            dG_dict[d] = U[i, out_factor] + Av[i, out_factor] * np.absolute(
                zeta) + PV[i, out_factor]


#            dG_dict[d] = G_hold + Av[i, out_factor] * (np.absolute(zeta) - 1)

# Computing the derivative of entropy as a funciton of strain using a finite difference approach
        if zeta == -1.:
            if temperature < 1.0 and inputs.statistical_mechanics == 'Quantum':
                # The limit of quantum entropy blows up at low T so computing -d^2 G / d C^2
                Gpp, _, _ = Pr.Gibbs_Free_Energy(
                    temperature + 0.05,
                    inputs.pressure,
                    inputs.program,
                    wavenumbers_dict['p'],
                    'p' + file_ending,
                    inputs.statistical_mechanics,
                    inputs.number_of_molecules,
                    Parameter_file=inputs.tinker_parameter_file)
                Gpm, _, _ = Pr.Gibbs_Free_Energy(
                    temperature - 0.05,
                    inputs.pressure,
                    inputs.program,
                    wavenumbers_dict['p'],
                    'p' + file_ending,
                    inputs.statistical_mechanics,
                    inputs.number_of_molecules,
                    Parameter_file=inputs.tinker_parameter_file)
                Gmp, _, _ = Pr.Gibbs_Free_Energy(
                    temperature + 0.05,
                    inputs.pressure,
                    inputs.program,
                    wavenumbers_dict['m'],
                    'm' + file_ending,
                    inputs.statistical_mechanics,
                    inputs.number_of_molecules,
                    Parameter_file=inputs.tinker_parameter_file)
                Gmm, _, _ = Pr.Gibbs_Free_Energy(
                    temperature - 0.05,
                    inputs.pressure,
                    inputs.program,
                    wavenumbers_dict['m'],
                    'm' + file_ending,
                    inputs.statistical_mechanics,
                    inputs.number_of_molecules,
                    Parameter_file=inputs.tinker_parameter_file)
                dS_dC[i] = -(Gpp - Gpm - Gpm + Gmm) / (4 *
                                                       (temperature * 0.01) *
                                                       LocGrd_dC[i])
            else:
                Sp = Pr.Vibrational_Entropy(temperature, wavenumbers_dict['p'], inputs.statistical_mechanics) / \
                     inputs.number_of_molecules
                Sm = Pr.Vibrational_Entropy(temperature, wavenumbers_dict['m'], inputs.statistical_mechanics) / \
                     inputs.number_of_molecules
                dS_dC[i] = (Sp - Sm) / (2 * LocGrd_dC[i])
        else:
            dzeta = inputs.zeta_numerical_step * 0.05
            ddG_dCdzeta[i] = ((U[i, 1] + PV[i, 1] + Av[i, 1] *
                               (zeta + dzeta)) -
                              (U[i, 1] + PV[i, 1] + Av[i, 1] *
                               (zeta - dzeta)) -
                              (U[i, 0] + PV[i, 0] + Av[i, 0] *
                               (zeta + dzeta)) +
                              (U[i, 0] + PV[i, 0] + Av[i, 0] *
                               (zeta - dzeta))) / (4 * dzeta * LocGrd_dC[i])

        # Calculating the finite difference of dG/deta for forward, central, and backwards
        dG_dC[i, 0] = (dG_dict['0'] - dG_dict['m']) / (LocGrd_dC[i])
        dG_dC[i, 1] = (dG_dict['p'] - dG_dict['m']) / (2 * LocGrd_dC[i])
        dG_dC[i, 2] = (dG_dict['p'] - dG_dict['0']) / (LocGrd_dC[i])

        # Computing the second derivative Gibbs free energy as a function of strain using a finite difference approach
        ddG_ddC[i, i] = (dG_dict['p'] - 2 * dG_dict['0'] +
                         dG_dict['m']) / (LocGrd_dC[i]**2)

        # otherwise, these stay as zero
        if i < off_diag_limit:
            # Computing the off diagonals of d**2 G/ (deta*deta)
            for j in np.arange(i + 1, off_diag_limit):
                if LocGrd_dC[j] < min_numerical_crystal_matrix:
                    continue  # don't bother to calculate them distance changed is too small.
                # Setting the strain of the second element for the off diagonal
                cm_array_2 = np.zeros(6)
                cm_array_2[j] = LocGrd_dC[j]

                # Expanding the structure for a 4 different strains
                delta2 = list()
                for di in delta1:
                    for dj in delta1:
                        d2 = di + dj
                        delta2.append(
                            d2
                        )  # create the list of 2D changes as length 2 strings
                        if dj == 'm':  # if the second dimension is a minus, dcrystal_matrix is negative
                            cm_factor = -1.0
                        else:
                            cm_factor = 1.0

                        Expand_Structure(
                            inputs,
                            di + file_ending,
                            'crystal_matrix',
                            d2,
                            dcrystal_matrix=array_to_triangle_crystal_matrix(
                                cm_factor * cm_array_2))

                for d in delta2:
                    wavenumbers_dict[d] = \
                        Wvn.Call_Wavenumbers(inputs, Coordinate_file=d + file_ending,
                                             ref_crystal_matrix=keyword_parameters['ref_crystal_matrix'],
                                             Gruneisen=keyword_parameters['Gruneisen'],
                                             Wavenumber_Reference=keyword_parameters['Wavenumber_Reference'])

                    G_hold, ignore, Av_hold = \
                        Pr.Gibbs_Free_Energy(temperature, inputs.pressure, inputs.program, wavenumbers_dict[d],
                                             d + file_ending, inputs.statistical_mechanics, inputs.number_of_molecules,
                                             Parameter_file=inputs.tinker_parameter_file)
                    dG_dict[d] = G_hold + Av_hold * (np.absolute(zeta) - 1)

                # Calculating the diagonal elements of d**2 G/(deta*deta)
                ddG_ddC[i, j] = (dG_dict['pp'] - dG_dict['pm'] - dG_dict['mp'] + dG_dict['mm']) / \
                                (4 * LocGrd_dC[i] * LocGrd_dC[j])

                # Making d**2 G/(deta*deta) symetric
                ddG_ddC[j, i] = ddG_ddC[i, j]

                # Removing excess files
                for d in delta2:
                    subprocess.call(['rm', d + file_ending])
        for d in delta1:
            subprocess.call(['rm', d + file_ending])

    # Calculating deta/dT for all strains
    if zeta == -1.:
        dC_dT = np.linalg.lstsq(ddG_ddC, dS_dC, rcond=None)[0]
        # Saving numerical outputs
        NO.raw_energies(np.array([U_0]), np.array([Av_0]), U, Av)
        left_minimum = NO.aniso_gradient(dG_dC, ddG_ddC, dS_dC, dC_dT)
    else:
        dC_dT = np.linalg.lstsq(ddG_ddC, -ddG_dCdzeta, rcond=None)[0]
        print(dC_dT)
        left_minimum = False
    return dC_dT, wavenumbers, left_minimum
Beispiel #6
0
def Isotropic_Local_Gradient(inputs, coordinate_file, temperature, LocGrd_dV,
                             **keyword_parameters):
    """
    This function calculates the local gradient of isotropic expansion for a given coordinate file
    
    :param Coordinate_file: file containing lattice parameters (and coordinates)
    :param Program: 'Tinker' Tinker molecular modeling
                    'Test' Test case
    :param Temperature: in Kelvin
    :param Pressure: in atm
    :param LocGrd_Vol_FracStep: fractional volumetric step size for numerical gradient 
    :param molecules_in_coord: number of molecules in Coordinate_file
    :param Statistical_mechanics: 'Classical' Classical mechanics
                                  'Quantum' Quantum mechanics
    :param Method: 'GiQ' Gradient isotropic QHA
                   'GiQg' Gradient isotropic QHA with Gruneisen Parameter
    :param keyword_parameters: Parameter_file, Gruneisen, Wavenumber_reference, Volume_reference
    
    Optional Parameters
    Parameter_file: program specific file containing force field parameters
    Gruneisen: isotropic Gruneisen parameter
    Wavenumber_reference: reference wavenumbers for the Gruneisen parameter
    Volume_reference: reference volume for the Gruneisen parameter
    """
    # Assigning general names for expanded and compressed structures
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    coordinate_plus = 'plus' + file_ending
    coordinate_minus = 'minus' + file_ending

    # Determining the volume of Coordinate_file
    volume = Pr.Volume(Program=inputs.program, Coordinate_file=coordinate_file)

    # Determining the change in lattice parameter for isotropic expansion
    dlattice_parameters_p = Isotropic_Change_Lattice_Parameters(
        (volume + LocGrd_dV) / volume, inputs.program, coordinate_file)
    dlattice_parameters_m = Isotropic_Change_Lattice_Parameters(
        (volume - LocGrd_dV) / volume, inputs.program, coordinate_file)

    # Building the isotropically expanded and compressed strucutres
    Expand_Structure(inputs,
                     coordinate_file,
                     'lattice_parameters',
                     'plus',
                     dlattice_parameters=dlattice_parameters_p)
    Expand_Structure(inputs,
                     coordinate_file,
                     'lattice_parameters',
                     'minus',
                     dlattice_parameters=dlattice_parameters_m)

    # Calculating wavenumbers coordinate_file, plus.*, and minus.*
    if inputs.method == 'GiQ':
        wavenumbers = Wvn.Call_Wavenumbers(inputs,
                                           Coordinate_file=coordinate_file)
        wavenumbers_plus = Wvn.Call_Wavenumbers(
            inputs, Coordinate_file=coordinate_plus)
        wavenumbers_minus = Wvn.Call_Wavenumbers(
            inputs, Coordinate_file=coordinate_minus)
    else:
        wavenumbers = Wvn.Call_Wavenumbers(
            inputs,
            Gruneisen=keyword_parameters['Gruneisen'],
            Wavenumber_Reference=keyword_parameters['Wavenumber_Reference'],
            Volume_Reference=keyword_parameters['Volume_Reference'],
            New_Volume=volume)
        wavenumbers_plus = Wvn.Call_Wavenumbers(
            inputs,
            Gruneisen=keyword_parameters['Gruneisen'],
            Wavenumber_Reference=keyword_parameters['Wavenumber_Reference'],
            Volume_Reference=keyword_parameters['Volume_Reference'],
            New_Volume=volume + LocGrd_dV)
        wavenumbers_minus = Wvn.Call_Wavenumbers(
            inputs,
            Gruneisen=keyword_parameters['Gruneisen'],
            Wavenumber_Reference=keyword_parameters['Wavenumber_Reference'],
            Volume_Reference=keyword_parameters['Volume_Reference'],
            New_Volume=volume - LocGrd_dV)

    # If temperature is zero, we assume that the local gradient is the same at 0.001K
    if temperature == 0.:
        temperature = 1e-03

    # Calculating the numerator of the local gradient -dS/dV
    dS = (Pr.Vibrational_Entropy(temperature, wavenumbers_plus,
                                 inputs.statistical_mechanics) /
          inputs.number_of_molecules - Pr.Vibrational_Entropy(
              temperature, wavenumbers_minus, inputs.statistical_mechanics) /
          inputs.number_of_molecules) / (2 * LocGrd_dV)

    # Calculating the denominator of the local gradient d**2G/dV**2
    ddG = (
        Pr.Gibbs_Free_Energy(temperature,
                             inputs.pressure,
                             inputs.program,
                             wavenumbers_plus,
                             coordinate_plus,
                             inputs.statistical_mechanics,
                             inputs.number_of_molecules,
                             Parameter_file=inputs.tinker_parameter_file)[0] -
        2 *
        Pr.Gibbs_Free_Energy(temperature,
                             inputs.pressure,
                             inputs.program,
                             wavenumbers,
                             coordinate_file,
                             inputs.statistical_mechanics,
                             inputs.number_of_molecules,
                             Parameter_file=inputs.tinker_parameter_file)[0] +
        Pr.Gibbs_Free_Energy(
            temperature,
            inputs.pressure,
            inputs.program,
            wavenumbers_minus,
            coordinate_minus,
            inputs.statistical_mechanics,
            inputs.number_of_molecules,
            Parameter_file=inputs.tinker_parameter_file)[0]) / (LocGrd_dV**2)

    # Computing the backward, central, and forward finite difference of dG/dV
    dG = np.zeros(3)
    dG[0] = (Pr.Gibbs_Free_Energy(
        temperature,
        inputs.pressure,
        inputs.program,
        wavenumbers,
        coordinate_file,
        inputs.statistical_mechanics,
        inputs.number_of_molecules,
        Parameter_file=inputs.tinker_parameter_file)[0] - Pr.Gibbs_Free_Energy(
            temperature,
            inputs.pressure,
            inputs.program,
            wavenumbers_minus,
            coordinate_minus,
            inputs.statistical_mechanics,
            inputs.number_of_molecules,
            Parameter_file=inputs.tinker_parameter_file)[0]) / (LocGrd_dV)

    dG[1] = (Pr.Gibbs_Free_Energy(
        temperature,
        inputs.pressure,
        inputs.program,
        wavenumbers_plus,
        coordinate_plus,
        inputs.statistical_mechanics,
        inputs.number_of_molecules,
        Parameter_file=inputs.tinker_parameter_file)[0] - Pr.Gibbs_Free_Energy(
            temperature,
            inputs.pressure,
            inputs.program,
            wavenumbers_minus,
            coordinate_minus,
            inputs.statistical_mechanics,
            inputs.number_of_molecules,
            Parameter_file=inputs.tinker_parameter_file)[0]) / (2 * LocGrd_dV)

    dG[2] = (Pr.Gibbs_Free_Energy(
        temperature,
        inputs.pressure,
        inputs.program,
        wavenumbers_plus,
        coordinate_plus,
        inputs.statistical_mechanics,
        inputs.number_of_molecules,
        Parameter_file=inputs.tinker_parameter_file)[0] - Pr.Gibbs_Free_Energy(
            temperature,
            inputs.pressure,
            inputs.program,
            wavenumbers,
            coordinate_file,
            inputs.statistical_mechanics,
            inputs.number_of_molecules,
            Parameter_file=inputs.tinker_parameter_file)[0]) / (LocGrd_dV)

    # Saving numerical outputs
    left_minimum = NO.iso_gradient(dG, ddG, dS, dS / ddG)

    # Removing excess files
    subprocess.call(['rm', coordinate_plus, coordinate_minus])
    return dS / ddG, wavenumbers, volume, left_minimum
Beispiel #7
0
def pressure_setup(Temperature=[0.0, 25.0, 50.0, 75.0, 100.0],
                   Pressure=1.,
                   Method='HA',
                   Program='Test',
                   Output='out',
                   Coordinate_file='molecule.xyz',
                   Parameter_file='keyfile.key',
                   molecules_in_coord=1,
                   properties_to_save=['G', 'T'],
                   NumAnalysis_method='RK4',
                   NumAnalysis_step=25.0,
                   LocGrd_Vol_FracStep='',
                   LocGrd_CMatrix_FracStep='',
                   StepWise_Vol_StepFrac=1.5e-3,
                   StepWise_Vol_LowerFrac=0.97,
                   StepWise_Vol_UpperFrac=1.16,
                   Statistical_mechanics='Classical',
                   Gruneisen_Vol_FracStep=1.5e-3,
                   Gruneisen_Lat_FracStep=1.0e-3,
                   Wavenum_Tol=-1.,
                   Gradient_MaxTemp=300.0,
                   Aniso_LocGrad_Type='6D',
                   min_RMS_gradient=0.01,
                   eq_of_state='Murnaghan',
                   gru_from_0T_0P=True,
                   cp2kroot='BNZ_NMA_p2'):

    if Method == 'HA':
        print(
            "Pressure vs. Temperature methods have not been implimented for the Harmonic Approximation"
        )
        sys.exit()

    elif (Method == 'GaQ') or (Method == 'GaQg'):
        print(
            "Pressure vs. Temperature methods have not been implimented for anisotropic expansion"
        )
        sys.exit()

    file_ending = psf.assign_coordinate_file_ending(Program)

    # Making an array of volume fractions
    V_frac = np.arange(StepWise_Vol_LowerFrac, 1.0, StepWise_Vol_StepFrac)
    V_frac = np.arange(StepWise_Vol_LowerFrac, StepWise_Vol_UpperFrac,
                       StepWise_Vol_StepFrac)

    # Volume of lattice minimum strucutre
    V0 = Pr.Volume(Program=Program, Coordinate_file=Coordinate_file)

    # Making an array to store the potential energy and volume of each structure
    U = np.zeros(len(V_frac))
    V = np.zeros(len(V_frac))

    for i in range(len(V_frac)):
        # Expanding the structures and saving the required data
        Ex.Call_Expansion(Method,
                          'expand',
                          Program,
                          Coordinate_file,
                          molecules_in_coord,
                          min_RMS_gradient,
                          volume_fraction_change=V_frac[i],
                          Output='temporary',
                          Parameter_file=Parameter_file)
        U[i] = psf.Potential_energy('temporary' + file_ending,
                                    Program,
                                    Parameter_file=Parameter_file)
        V[i] = Pr.Volume(Program=Program,
                         Coordinate_file='temporary' + file_ending)
        subprocess.call(['rm', 'temporary' + file_ending])

    V0 = Pr.Volume(Program=Program, Coordinate_file=Coordinate_file)
    E0 = psf.Potential_energy(Coordinate_file,
                              Program,
                              Parameter_file=Parameter_file)

    if eq_of_state != 'None':
        [B, dB], _ = scipy.optimize.curve_fit(
            lambda V, B, dB: eos.EV_EOS(V, V0, B, dB, E0, eq_of_state),
            V,
            U,
            p0=[2., 2.])
        np.save(Output + '_EOS', [V0, B, dB, E0])

        plt.plot(V, (eos.EV_EOS(V, V0, B, dB, E0, eq_of_state) - U) /
                 molecules_in_coord)
        plt.xlabel('Volume [Ang.$^{3}$]', fontsize=18)
        plt.ylabel('$\delta(U_{EOS})$ [kcal/mol]', fontsize=18)
        plt.tight_layout()
        plt.savefig(Output + '_EOS_dU.pdf')
        plt.close()

        plt.plot(
            V,
            eos.EV_EOS(V, V0, B, dB, E0, eq_of_state) / molecules_in_coord)
        plt.scatter(V, U / molecules_in_coord)
        plt.xlabel('Volume [Ang.$^{3}$]', fontsize=18)
        plt.ylabel('$U$ [kcal/mol]', fontsize=18)
        plt.tight_layout()
        plt.savefig(Output + '_EOS_U.pdf')
        plt.close()

        if gru_from_0T_0P == True:
            EOS_TvP_Gru_0T_0P(Method, Temperature, Pressure, Program, Output,
                              Coordinate_file, Parameter_file,
                              molecules_in_coord, properties_to_save,
                              Statistical_mechanics, Gruneisen_Vol_FracStep,
                              Wavenum_Tol, min_RMS_gradient, eq_of_state,
                              cp2kroot, V0, B, dB, E0)
            sys.exit()

        else:
            lattice_parameter_name = Coordinate_file.split('.')[0]
            for i in range(len(Pressure)):
                v_p = scipy.optimize.minimize(pressure_minimization_at_0T,
                                              V0,
                                              args=(Pressure[i], V0, B, dB,
                                                    eq_of_state),
                                              method='Nelder-Mead',
                                              tol=1.e-15).x

                print(Pressure[i], v_p)
                V = np.arange(100, 5000, 10)
                plt.plot(
                    V,
                    eos.EV_EOS(V, V0, B, dB, E0, eq_of_state) +
                    Pressure[i] * V)
                plt.scatter(
                    v_p,
                    eos.EV_EOS(v_p, V0, B, dB, E0, eq_of_state) +
                    Pressure[i] * v_p)
                plt.show()
    else:
        print("WARNING: Option of None for EOS is not set up yet, exiting.")
        sys.exit()
def Setup_Anisotropic_Gruneisen(inputs):
    r"""
    Computes the anisotropic Gruneisen parameters of the lattice minimum structure

    Parameters
    ----------
    inputs: Class
        Contains all user defined values and filled in with default program values

    Returns
    -------
    gruneisen: List[float]
        Anisotropic Gruneisen parameters.
    wavenumbers_ref: List[float]
        Wavenumbers [cm$^{-1}$] of the lattice minimum structure.
    """
    # Determining the file ending for the program used
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    # Determining if the anisotropic Gruneisen parameters have already been calculated
    # Typically, this take a long time. It may be advantageous to parallelize this in the future.
    re_run = False
    if os.path.isfile('GRU_wvn.npy') and os.path.isfile('GRU_eigen.npy'):
        # Loading in previously computed Gruneisen parameters
        wavenumbers = np.load('GRU_wvn.npy')
        eigenvectors = np.load('GRU_eigen.npy')
        re_run = True
    else:
        # Computing the reference wavenumbers and eigenvectors (the lattice minimum structure fed in)
        wavenumbers_ref, eigenvectors_ref = psf.Wavenumber_and_Vectors(
            inputs.program, inputs.coordinate_file,
            inputs.tinker_parameter_file)

        # Setting the number of vibrational modes
        number_of_modes = int(len(wavenumbers_ref))
        number_of_atoms = psf.atoms_count(inputs.program,
                                          inputs.coordinate_file)

        # Setting a place to store all the wavenumbers and eigenvalues for the Gruenisen parameters
        wavenumbers = np.zeros((7, number_of_modes))
        eigenvectors = np.zeros((7, number_of_modes, 3 * number_of_atoms))
        wavenumbers[0] = wavenumbers_ref
        eigenvectors[0] = eigenvectors_ref

        # Out putting information for user output
        NO.start_anisoGru()

        # Saving that data computed thusfar in case the system crashes or times out
        np.save('GRU_eigen', eigenvectors)
        np.save('GRU_wvn', wavenumbers)

    # Cycling through the six principle six principal strains
    for i in range(6):
        if re_run and (wavenumbers[i + 1, 3] != 0.):
            # Skipping a given strain if the wavenumbers have previously been computed and loaded in
            pass
        else:
            # Expanding the lattice minimum structure in the direction of the i-th principal strain
            applied_strain = np.zeros(6)
            applied_strain[i] = inputs.gruneisen_matrix_strain_stepsize
            Ex.Expand_Structure(
                inputs,
                inputs.coordinate_file,
                'strain',
                'temp',
                strain=Ex.strain_matrix(applied_strain),
                crystal_matrix=Ex.Lattice_parameters_to_Crystal_matrix(
                    psf.Lattice_parameters(inputs.program,
                                           inputs.coordinate_file)))

            # Computing the strained wavenumbers and eigenvectors
            wavenumbers_unorganized, eigenvectors_unorganized = psf.Wavenumber_and_Vectors(
                inputs.program, 'temp' + file_ending,
                inputs.tinker_parameter_file)

            # Determining how the strained eigenvectors match up with the reference structure
            z, weight = matching_eigenvectors_of_modes(
                number_of_modes, eigenvectors[0], eigenvectors_unorganized)
            NO.GRU_weight(
                weight
            )  # Writing out the precision of matching the modes with one another

            # Re-organizing the expanded wavenumbers
            wavenumbers[i + 1], eigenvectors[i + 1] = reorder_modes(
                z, wavenumbers_unorganized, eigenvectors_unorganized)

            # Saving the eigenvectors and wavenumbers
            np.save('GRU_eigen', eigenvectors)
            np.save('GRU_wvn', wavenumbers)

            # Removing the strained coordinate file
            subprocess.call(['rm', 'temp' + file_ending])

    # Calculating the Gruneisen parameters due to the six principal strains
    gruneisen = np.zeros((number_of_modes, 6))
    for i in range(6):
        # Calculating the Gruneisen parameters
        gruneisen[3:, i] = -(np.log(wavenumbers[i + 1, 3:]) - np.log(wavenumbers[0, 3:])) \
                           / inputs.gruneisen_matrix_strain_stepsize
    return gruneisen, wavenumbers[0]
def Setup_Isotropic_Gruneisen(inputs):
    r"""
    Computes the isotropic Gruneisen parameters of the lattice minimum structure

    Parameters
    ----------
    inputs: Class
        Contains all user defined values and filled in with default program values

    Returns
    -------
    gruneisen: List[float]
        Isotropic Gruneisen parameters.
    wavenumbers_ref: List[float]
        Wavenumbers [cm$^{-1}$] of the lattice minimum structure.
    volume_ref: float
        Volume [$\AA^{3}$] of the lattice minimum structure.
    """

    # Change in lattice parameters for expanded structure
    dLattice_Parameters = Ex.Isotropic_Change_Lattice_Parameters(
        (1 + inputs.gruneisen_volume_fraction_stepsize), inputs.program,
        inputs.coordinate_file)

    # Determining wavenumbers of lattice strucutre and expanded structure
    lattice_parameters = psf.Lattice_parameters(inputs.program,
                                                inputs.coordinate_file)

    # Expanding structure
    Ex.Expand_Structure(inputs,
                        inputs.coordinate_file,
                        'lattice_parameters',
                        'temp',
                        dlattice_parameters=dLattice_Parameters)

    # File ending of the coordinate file
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    # Computing the reference wavenumbers and eigenvectors
    wavenumbers_ref, eigenvectors_ref = psf.Wavenumber_and_Vectors(
        inputs.program, inputs.coordinate_file, inputs.tinker_parameter_file)
    number_of_modes = len(wavenumbers_ref)  # Number of vibrational modes

    # Computing the strained wavenumbers and eigenvectors
    wavenumbers_unorganized, eigenvectors_unorganized = psf.Wavenumber_and_Vectors(
        inputs.program, 'temp' + file_ending, inputs.tinker_parameter_file)

    # Determining how the strained eigenvectors match up with the reference structure
    z, weight = matching_eigenvectors_of_modes(number_of_modes,
                                               eigenvectors_ref,
                                               eigenvectors_unorganized)

    # Outputing information for user output
    NO.start_isoGru()
    NO.GRU_weight(weight)

    # Re-organizing the expanded wavenumbers
    wavenumbers_organized, eigenvectors_organized = reorder_modes(
        z, wavenumbers_unorganized, eigenvectors_unorganized)

    # Calculating the volume of the lattice minimum and expanded structure
    volume_ref = Pr.Volume(lattice_parameters=lattice_parameters)
    volume_expand = volume_ref + inputs.gruneisen_volume_fraction_stepsize * volume_ref

    # Computing the Gruneisen parameters
    gruneisen = np.zeros(len(wavenumbers_ref))
    gruneisen[3:] = -(np.log(wavenumbers_ref[3:]) - np.log(wavenumbers_organized[3:])) / \
                    (np.log(volume_ref) - np.log(volume_expand))
    for x in range(0, len(gruneisen)):
        if wavenumbers_ref[x] == 0:
            gruneisen[x] = 0.0

    # Removing extra files created in process
    subprocess.call(['rm', 'temp' + file_ending])
    return gruneisen, wavenumbers_ref, volume_ref
Beispiel #10
0
def temperature_lattice_dynamics(inputs, data, input_file='input.yaml'):
    # Geometry and lattice optimizing the input structure
    if inputs.tinker_xtalmin and (inputs.program == 'Tinker'):
        psf.tinker_xtalmin(inputs)
        inputs.tinker_xtalmin = False

    # Running the harmonic approximation
    if inputs.method == 'HA':
        print("Performing Harmonic Approximation")

        # Determining if the wavenumbers have already been calculated
        if os.path.isfile(inputs.output + '_' + inputs.method + '_WVN.npy'):
            # Loading in te wavenumbers
            wavenumbers = np.load(inputs.output + '_' + inputs.method + '_WVN.npy')
            print("   Importing wavenumbers from:" + inputs.output + '_' + inputs.method + '_WVN.npy')
        else:
            # Computing the wavenumbers and saving them
            print("   Computing wavenumbers of coordinate file")
            wavenumbers = Wvn.Call_Wavenumbers(inputs, Coordinate_file=inputs.coordinate_file,
                                               Parameter_file=inputs.tinker_parameter_file)
            np.save(inputs.output + '_' + inputs.method + '_WVN', wavenumbers)

        # Making sure all wavenumbers are within the user specified tolerance
        if all(wavenumbers[:3] < inputs.wavenumber_tolerance) and \
                all(wavenumbers[:3] > -1. * inputs.wavenumber_tolerance):
            print("   All wavenumbers are greater than tolerance of: " + str(inputs.wavenumber_tolerance) + " cm^-1")
            properties = Pr.Properties_with_Temperature_and_Pressure(inputs, inputs.coordinate_file, wavenumbers)
            print("   All properties have been saved in " + inputs.output + "_raw.npy")
            np.save(inputs.output + '_raw', properties)
            print("   Saving user specified properties in indipendent files:")
            Pr.Save_Properties_2D(inputs, properties)
        exit()
    else:
        if not os.path.isdir('Cords'):
            print("Creating directory 'Cords/' to store structures along Gibbs free energy path")
            subprocess.call(['mkdir', 'Cords'])

    # Expanding the crystal with the zero point energy
    if (inputs.statistical_mechanics == 'Quantum') and (inputs.coordinate_file != 'ezp_minimum' + psf.assign_coordinate_file_ending(inputs.program)):
        if any(inputs.gradient_matrix_fractions != 0.):
            crystal_matrix_array = Ex.triangle_crystal_matrix_to_array(Ex.Lattice_parameters_to_Crystal_matrix(
                psf.Lattice_parameters(inputs.program, inputs.coordinate_file)))
            LocGrd_dC = np.absolute(inputs.gradient_matrix_fractions * crystal_matrix_array)
            for i in range(len(LocGrd_dC)):
                if LocGrd_dC[i] == 0.:
                    LocGrd_dC[i] = inputs.gradient_matrix_fractions[i]
        else:
            LocGrd_dC = Ss.anisotropic_gradient_settings(inputs, data, input_file)

        TNA.anisotropic_gradient_expansion_ezp(inputs, LocGrd_dC)
        inputs.coordinate_file = 'ezp_minimum' + psf.assign_coordinate_file_ending(inputs.program)

    # Running through QHA
    if (inputs.method == 'SiQ') or (inputs.method == 'SiQg'):
        # Stepwise Isotropic QHA
        print("Performing Stepwise Isotropic Quasi-Harmonic Approximation")
        properties = TNA.stepwise_expansion(inputs)

        print("   Saving user specified properties in indipendent files:")
        Pr.Save_Properties_1D(inputs, properties)

    elif (inputs.method == 'GiQ') or (inputs.method == 'GiQg'):
        if inputs.gradient_vol_fraction == (0. or None):
            LocGrd_dV = Ss.isotropic_gradient_settings(inputs)
        else:
            V_0 = Pr.Volume(Program=inputs.program, Coordinate_file=inputs.coordinate_file)
            LocGrd_dV = inputs.gradient_vol_fraction * V_0
        # Gradient Isotropic QHA
        print("Performing Gradient Isotropic Quasi-Harmonic Approximation")
        properties = TNA.Isotropic_Gradient_Expansion(inputs, LocGrd_dV)

        print("   Saving user specified properties in indipendent files:")
        Pr.Save_Properties_1D(inputs, properties)
    elif (inputs.method == 'GaQ') or (inputs.method == 'GaQg'):
        if inputs.statistical_mechanics == 'Classical':
            if any(inputs.gradient_matrix_fractions != 0.):
                crystal_matrix_array = Ex.triangle_crystal_matrix_to_array(Ex.Lattice_parameters_to_Crystal_matrix(
                    psf.Lattice_parameters(inputs.program, inputs.coordinate_file)))
                LocGrd_dC = np.absolute(inputs.gradient_matrix_fractions * crystal_matrix_array)
                for i in range(len(LocGrd_dC)):
                    if LocGrd_dC[i] == 0.:
                        LocGrd_dC[i] = inputs.gradient_matrix_fractions[i]
            else:
                LocGrd_dC = Ss.anisotropic_gradient_settings(inputs, data, input_file)

        if inputs.anisotropic_type != '1D':
            print("Performing Gradient Anisotropic Quasi-Harmonic Approximation")
            properties = TNA.Anisotropic_Gradient_Expansion(inputs, LocGrd_dC)
            print("   Saving user specified properties in independent files:")
            Pr.Save_Properties_1D(inputs, properties)
        else:
            print("Performing 1D-Gradient Anisotropic Quasi-Harmonic Approximation")
            properties = TNA.Anisotropic_Gradient_Expansion_1D(inputs, LocGrd_dC)
            print("   Saving user specified properties in independent files:")
            Pr.Save_Properties_1D(inputs, properties)
        
    elif inputs.method == 'SaQply':
        print("Performing Quasi-Anisotropic Quasi-Harmonic Approximation")
        properties = TNA.stepwise_expansion(inputs)
        print("   Saving user specified properties in independent files:")
        Pr.Save_Properties_1D(inputs, properties)
    print("Lattice dynamic calculation is complete!")
Beispiel #11
0
def isotropic_gradient_settings(inputs):
    #Coordinate_file, Program, Parameter_file, molecules_in_coord, min_RMS_gradient, Output,
    #                       Pressure):
    # Determining the file ending based on the program
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    # Setting the energy cutoff
    cutoff = program_cutoff(inputs.program)

    # Fractional step sizes to take
    steps = np.array(
        [5e-05, 1e-04, 5e-04, 1e-03, 5e-03, 1e-02, 5e-02, 1e-01, 5e-01])

    # Number of total step sizes
    n_steps = len(steps)

    # Potential energy of input file and a place to store the expanded structures potential energy
    U_0 = (psf.Potential_energy(inputs.coordinate_file, inputs.program, Parameter_file=inputs.tinker_parameter_file) \
           + Pr.PV_energy(inputs.pressure, Pr.Volume(Program=inputs.program, Coordinate_file=inputs.coordinate_file))) / \
          inputs.number_of_molecules
    U = np.zeros((n_steps))

    for i in range(n_steps):
        # Setting how much the lattice parameters must be changed
        dlattice_parameters = Ex.Isotropic_Change_Lattice_Parameters(
            1. + steps[i], inputs.program, inputs.coordinate_file)

        # Expanding the strucutre
        Ex.Expand_Structure(inputs,
                            inputs.coordinate_file,
                            'lattice_parameters',
                            inputs.output,
                            dlattice_parameters=dlattice_parameters)

        # Computing the potential energy
        U[i] = (psf.Potential_energy(inputs.output + file_ending, inputs.program,
                                    Parameter_file=inputs.tinker_parameter_file) \
                + Pr.PV_energy(inputs.pressure, Pr.Volume(Program=inputs.program, Coordinate_file=inputs.output + file_ending))) / \
                inputs.number_of_molecules
        subprocess.call(['rm', inputs.output + file_ending])

        if (U[i] - U_0) > cutoff:
            # Ending the run if we've exceeded the energy cut-off
            LocGrd_Vol_FracStep = steps[i]
            steps = steps[:i + 1]
            U = U[:i + 1]
            break

    # Plotting the results
    plt.plot(np.log10(steps), U - U_0, linestyle='--', marker='o')
    plt.xlabel('$\log({dV/V_{0}})$', fontsize=22)
    plt.ylabel('$\Delta U$ [kcal/mol]', fontsize=22)
    plt.ylim((0., 2 * cutoff))
    plt.axhline(y=cutoff, c='grey', linestyle='--')
    plt.tight_layout()
    plt.savefig(inputs.output + '_LocGrd_Vol_FracStep.pdf')
    plt.close()

    # Printing step size
    print("After analysis, LocGrd_Vol_FracStep = ", LocGrd_Vol_FracStep)

    # initial volume
    V_0 = Pr.Volume(Program=inputs.program,
                    Coordinate_file=inputs.coordinate_file)

    # returning the value of dV
    return LocGrd_Vol_FracStep * V_0
Beispiel #12
0
def anisotropic_gradient_settings(inputs, data, input_file):
    # Determining the file ending based on the program
    file_ending = psf.assign_coordinate_file_ending(inputs.program)

    # Setting the energy cutoff
    cutoff = program_cutoff(inputs.program)

    # Fractional step sizes to take
    steps = np.array([
        5e-05, 1e-04, 5e-04, 1e-03, 5e-03, 1e-02, 5e-02, 1e-01, 5e-01, 1., 5.,
        1e01, 5e01, 1e02, 5e02, 1e03, 5e03, 1e04, 5e04
    ])

    # Number of total step sizes
    n_steps = len(steps)

    # Potential energy of input file and a place to store the expanded structures potential energy
    U_0 = psf.Potential_energy(inputs.coordinate_file, inputs.program, Parameter_file=inputs.tinker_parameter_file) / \
          inputs.number_of_molecules
    U = np.zeros((6, n_steps))

    # Determining the tensor parameters of the input file
    crystal_matrix = Ex.Lattice_parameters_to_Crystal_matrix(
        psf.Lattice_parameters(inputs.program, inputs.coordinate_file))
    crystal_matrix_array = Ex.triangle_crystal_matrix_to_array(crystal_matrix)

    LocGrd_CMatrix_FracStep = np.zeros(6)
    LocGrd_CMatrix_Step = np.zeros(6)
    for j in range(6):
        for i in range(n_steps):
            dlattice_matrix_array = np.zeros(6)
            dlattice_matrix_array[j] = np.absolute(crystal_matrix_array[j] *
                                                   steps[i])
            if np.absolute(crystal_matrix_array[j]) < 1e-4:
                dlattice_matrix_array[j] = steps[i]
            dlattice_matrix = Ex.array_to_triangle_crystal_matrix(
                dlattice_matrix_array)
            Ex.Expand_Structure(inputs,
                                inputs.coordinate_file,
                                'crystal_matrix',
                                inputs.output,
                                dcrystal_matrix=dlattice_matrix)
            U[j, i] = psf.Potential_energy(
                inputs.output + file_ending,
                inputs.program,
                Parameter_file=inputs.tinker_parameter_file
            ) / inputs.number_of_molecules
            subprocess.call(['rm', inputs.output + file_ending])
            if (U[j, i] - U_0) > cutoff:
                LocGrd_CMatrix_FracStep[j] = 10**np.interp(
                    cutoff, U[j, i - 1:i + 1] - U_0,
                    np.log10(steps[i - 1:i + 1]))
                LocGrd_CMatrix_Step[j] = np.absolute(
                    crystal_matrix_array[j] * LocGrd_CMatrix_FracStep[j])
                if np.absolute(crystal_matrix_array[j]) < 1e-4:
                    LocGrd_CMatrix_Step[j] = LocGrd_CMatrix_FracStep[j]
                break
        plt.scatter(np.log10(LocGrd_CMatrix_FracStep[j]),
                    cutoff,
                    marker='x',
                    color='r')
        plt.plot(np.log10(steps[:i + 1]),
                 U[j, :i + 1] - U_0,
                 linestyle='--',
                 marker='o',
                 label='C' + str(j + 1))

    # Plotting the results
    plt.xlabel('$\log({dC/C_{0}})$', fontsize=22)
    plt.ylabel('$\Delta U$ [kcal/mol]', fontsize=22)
    plt.ylim((0., 2 * cutoff))
    plt.axhline(y=cutoff, c='grey', linestyle='--')
    plt.legend(loc='upper right', ncol=2, fontsize=18)
    plt.tight_layout()
    plt.savefig(inputs.output + '_LocGrd_CMatrix_FracStep.pdf')
    plt.close()

    # Printing step size
    print("After analysis, LocGrd_CMatrix_FracStep = ",
          LocGrd_CMatrix_FracStep)
    data['gradient']['matrix_fractions'] = LocGrd_CMatrix_FracStep.tolist()
    with open(input_file, 'w') as yaml_file:
        yaml.dump(data, yaml_file, default_flow_style=False)
    # returning the value of dV
    return LocGrd_CMatrix_Step