def test_compare_with_cp2k():
    """
    Test overlap matrix transformation from cartesian to spherical
    """
    # Overlap matrix in cartesian coordinates
    basisname = "DZVP-MOLOPT-SR-GTH"
    # Molecular geometry in a.u.
    atoms = change_mol_units(readXYZ('test/test_files/ethylene.xyz'))

    dictCGFs = create_dict_CGFs(path_hdf5, basisname, atoms)

    # Compute the overlap matrix using the general multipole expression
    rs = calcMtxOverlapP(atoms, dictCGFs)
    number_of_contracted = sum(len(dictCGFs[at.symbol]) for at in atoms)
    mtx_overlap = triang2mtx(rs, number_of_contracted)

    dict_global_norms = compute_normalization_sphericals(dictCGFs)

    with h5py.File(path_hdf5, 'r') as f5:
        transf_mtx = calc_transf_matrix(f5, atoms, basisname,
                                        dict_global_norms, 'cp2k')

    # Use a sparse representation of the transformation matrix
    transf_mtx = sparse.csr_matrix(transf_mtx)
    transpose = transf_mtx.transpose()

    # Compare the results with CP2K overlap matrix
    test = transf_mtx.dot(sparse.csr_matrix.dot(mtx_overlap, transpose))

    sum_diag = np.sum(np.diag(test))
    number_of_sphericals = transf_mtx.shape[0]

    assert abs(sum_diag - number_of_sphericals) < 1e-16
Ejemplo n.º 2
0
def test_compare_with_cp2k():
    """
    Test overlap matrix transformation from cartesian to spherical
    """
    # Overlap matrix in cartesian coordinates
    basisname = "DZVP-MOLOPT-SR-GTH"
    # Molecular geometry in a.u.
    atoms = change_mol_units(
        readXYZ('test/test_files/Cd33Se33_fivePoints.xyz'))

    dictCGFs = create_dict_CGFs(path_hdf5, basisname, atoms)

    # Compute the overlap matrix using the general multipole expression
    rs = calcMtxOverlapP(atoms, dictCGFs)
    mtx_overlap = triang2mtx(rs, 1452)  # there are 1452 Cartesian basis CGFs

    dict_global_norms = compute_normalization_sphericals(dictCGFs)

    with h5py.File(path_hdf5, 'r') as f5:
        transf_mtx = calc_transf_matrix(f5, atoms, basisname,
                                        dict_global_norms, 'cp2k')

    # Use a sparse representation of the transformation matrix
    transf_mtx = sparse.csr_matrix(transf_mtx)
    transpose = transf_mtx.transpose()

    # Compare the results with CP2K overlap matrix
    test = transf_mtx.dot(sparse.csr_matrix.dot(mtx_overlap, transpose))
    expected = np.load('test/test_files/overlap_Cd33Se33_cp2k.npy')

    assert np.allclose(test, expected, atol=1e-5)
Ejemplo n.º 3
0
def test_cube():
    """
    Test the density compute to create a cube file.
    """
    if not os.path.exists(scratch_path):
        os.makedirs(scratch_path)

    # Overlap matrix in cartesian coordinates
    basisname = "DZVP-MOLOPT-SR-GTH"

    # Read coordinates
    molecule = change_mol_units(readXYZ(path_xyz))

    # String representation of the molecule
    geometries = split_file_geometries(path_xyz)

    try:
        shutil.copy(path_original_hdf5, path_test_hdf5)
        # Contracted Gauss functions
        dictCGFs = create_dict_CGFs(path_test_hdf5, basisname, molecule)

        dict_global_norms = compute_normalization_sphericals(dictCGFs)

        # Compute the transformation matrix and store it in the HDF5
        with h5py.File(path_test_hdf5, 'r') as f5:
            transf_mtx = calc_transf_matrix(
                f5, molecule, basisname, dict_global_norms, 'cp2k')

        path_transf_mtx = join(project_name, 'trans_mtx')
        store_arrays_in_hdf5(
            path_test_hdf5, path_transf_mtx, transf_mtx)

        # voxel size and Number of steps
        grid_data = GridCube(0.300939, 80)

        rs = workflow_compute_cubes(
            'cp2k', project_name, package_args, path_time_coeffs=None,
            grid_data=grid_data, guess_args=None, geometries=geometries,
            dictCGFs=dictCGFs, calc_new_wf_guess_on_points=[],
            path_hdf5=path_test_hdf5, enumerate_from=0, package_config=None,
            traj_folders=None, work_dir=scratch_path, basisname=basisname,
            hdf5_trans_mtx=path_transf_mtx, nHOMO=6, orbitals_range=(6, 6),
            ignore_warnings=False)

        xs = retrieve_hdf5_data(path_test_hdf5, rs)[0]
        np.save("grid_density", xs)

    finally:
        # Remove intermediate results
        shutil.rmtree(scratch_path)
Ejemplo n.º 4
0
def main(parser):
    """
    These calculation is based on the paper:
    `Theoretical analysis of electronic band structure of
    2- to 3-nm Si nanocrystals`
    PHYSICAL REVIEW B 87, 195420 (2013)
    """
    # Parse Command line
    project_name, path_hdf5, path_xyz, lattice_cte, basis_name, lower, \
        upper = read_cmd_line(parser)
    # Coordinates transformation
    atoms = readXYZ(path_xyz)
    symbols = np.array([at.symbol for at in atoms])
    coords_angstrom = np.concatenate([at.xyz for at in atoms])
    angstroms_to_au = 1.889725989
    coords = angstroms_to_au * coords_angstrom
    lattice_cte = lattice_cte * angstroms_to_au

    # Dictionary containing as key the atomic symbols and as values the set of CGFs
    home = os.path.expanduser('~')
    basiscp2k = join(home, "Cp2k/cp2k_basis/BASIS_MOLOPT")
    potcp2k = join(home, "Cp2k/cp2k_basis/GTH_POTENTIALS")
    cp2k_config = {"basis": basiscp2k, "potential": potcp2k}
    dictCGFs = create_dict_CGFs(path_hdf5,
                                basis_name,
                                atoms,
                                package_config=cp2k_config)
    count_cgfs = np.vectorize(lambda s: len(dictCGFs[s]))
    number_of_basis = np.sum(np.apply_along_axis(count_cgfs, 0, symbols))

    # K-space grid to calculate the fuzzy band
    nPoints = 20
    # grid_k_vectors = grid_kspace(initial, final, nPoints)
    map_grid_kspace = lambda ps: [grid_kspace(i, f, nPoints) for i, f in ps]

    # Grid
    grids_alpha = map_grid_kspace(create_alpha_paths(lattice_cte))
    grids_beta = map_grid_kspace(create_beta_paths(lattice_cte))

    # Apply the whole fourier transform to the subset of the grid
    # correspoding to each process
    momentum_density = partial(compute_momentum_density, project_name, symbols,
                               coords, dictCGFs, number_of_basis, path_hdf5)

    orbitals = list(range(lower, upper + 1))
    dim_x = len(orbitals)
    result = np.empty((dim_x, nPoints))
    with Pool() as p:
        for i, orb in enumerate(orbitals):
            # compute density
            density = momentum_density(orb)
            alphas = [
                normalize(p.map(density, grid_k)) for grid_k in grids_alpha
            ]
            betas = [
                normalize(p.map(density, grid_k)) for grid_k in grids_beta
            ]
            rs_alphas = normalize(np.stack(alphas).sum(axis=0))
            rs_betas = normalize(np.stack(betas).sum(axis=0))
            rss = normalize(rs_alphas + rs_betas)
            result[i] = rss
            np.save('alphas', rs_alphas)
            np.save('betas', rs_betas)
            print("Orb: ", orb)
            print(rss)
Ejemplo n.º 5
0
def initialize(project_name: str,
               path_traj_xyz: str,
               basisname: str,
               enumerate_from: int = 0,
               calculate_guesses: int = 'first',
               path_hdf5: str = None,
               scratch_path: str = None,
               path_basis: str = None,
               path_potential: str = None,
               geometry_units: str = 'angstrom') -> Dict:
    """
    Initialize all the data required to schedule the workflows associated with
    the nonadaibatic coupling
    """
    # User variables
    username = getpass.getuser()

    # Scratch
    if scratch_path is None:
        scratch_path = join('/tmp', username, project_name)
        logger.warning(
            "path to scratch was not defined, using: {}".format(scratch_path))

    # If the directory does not exist create it
    if not os.path.exists(scratch_path):
        os.makedirs(scratch_path)

    # Cp2k configuration files
    cp2k_config = {"basis": path_basis, "potential": path_potential}

    # HDF5 path
    if path_hdf5 is None:
        path_hdf5 = join(scratch_path, 'quantum.hdf5')
        logger.warning(
            "path to the HDF5 was not defined, using: {}".format(path_hdf5))

    # all_geometries type :: [String]
    geometries = split_file_geometries(path_traj_xyz)

    # Create a folder for each point the the dynamics
    traj_folders = create_point_folder(scratch_path, len(geometries),
                                       enumerate_from)
    if calculate_guesses is None:
        points_guess = []
    elif calculate_guesses.lower() in 'first':
        # Calculate new Guess in the first geometry
        points_guess = [enumerate_from]
        msg = "An initial Calculation will be computed as guess for the wave function"
        logger.info(msg)
    else:
        # Calculate new Guess in each geometry
        points_guess = [enumerate_from + i for i in range(len(geometries))]
        msg = "A guess calculation will be done for each geometry"
        logger.info(msg)

    # Generate a list of tuples containing the atomic label
    # and the coordinates to generate
    # the primitive CGFs
    atoms = parse_string_xyz(geometries[0])
    if 'angstrom' in geometry_units.lower():
        atoms = change_mol_units(atoms)

    # CGFs per element
    dictCGFs = create_dict_CGFs(path_hdf5,
                                basisname,
                                atoms,
                                package_config=cp2k_config)

    # Calculcate the matrix to transform from cartesian to spherical
    # representation of the overlap matrix
    hdf5_trans_mtx = store_transf_matrix(path_hdf5, atoms, dictCGFs, basisname,
                                         project_name)

    d = {
        'package_config': cp2k_config,
        'path_hdf5': path_hdf5,
        'calc_new_wf_guess_on_points': points_guess,
        'geometries': geometries,
        'enumerate_from': enumerate_from,
        'dictCGFs': dictCGFs,
        'work_dir': scratch_path,
        'traj_folders': traj_folders,
        'basisname': basisname,
        'hdf5_trans_mtx': hdf5_trans_mtx
    }

    return d