Skip to content

Commit

Permalink
Add option to export structure in LAMMPS format in abistruct.py export
Browse files Browse the repository at this point in the history
  • Loading branch information
gmatteo committed May 26, 2024
1 parent d25347b commit 3e5f1d1
Show file tree
Hide file tree
Showing 5 changed files with 69 additions and 76 deletions.
7 changes: 7 additions & 0 deletions abipy/core/structure.py
Original file line number Diff line number Diff line change
Expand Up @@ -1806,6 +1806,13 @@ def convert(self, fmt: str = "cif", **kwargs) -> str:
return Poscar(self).get_str(significant_figures=12)
except AttributeError:
return Poscar(self).get_string(significant_figures=12)

elif fmt.lower() == "lammps":
from pymatgen.io.lammps.data import LammpsData
# Convert the structure to a LAMMPS data file
lammps_data = LammpsData.from_structure(self)
return lammps_data.get_str()

else:
return super().to(fmt=fmt, **kwargs)

Expand Down
78 changes: 38 additions & 40 deletions abipy/eph/gstore.py
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,7 @@
from monty.functools import lazy_property
#from monty.termcolor import cprint
from abipy.core.structure import Structure
from abipy.core.mixins import AbinitNcFile, Has_Structure, Has_ElectronBands, Has_Header, NotebookWriter
from abipy.core.mixins import AbinitNcFile, Has_Structure, Has_ElectronBands, Has_Header #, NotebookWriter
from abipy.tools.typing import PathLike
#from abipy.tools.plotting import (add_fig_kwargs, get_ax_fig_plt, get_axarray_fig_plt, set_axlims, set_visible,
# rotate_ticklabels, ax_append_title, set_ax_xylabels, linestyles)
Expand All @@ -26,7 +26,6 @@
from abipy.eph.common import BaseEphReader



class GstoreFile(AbinitNcFile, Has_Header, Has_Structure, Has_ElectronBands): # , NotebookWriter):
"""
This file stores the e-ph matrix elements and provides methods to analyze and plot results.
Expand Down Expand Up @@ -178,15 +177,15 @@ def from_gstore(cls, gstore: GstoreFile, spin: int):

vk_cart_ibz, vkmat_cart_ibz = None, None
if ncr.with_vk == 1:
# NCF_CHECK(nctk_def_arrays(spin_ncid, nctkarr_t("vk_cart_ibz", "dp", "three, nb, gstore_nkibz")))
vk_cart_ibz = ncr.read_value("vk_cart_ibz", path=path)
# NCF_CHECK(nctk_def_arrays(spin_ncid, nctkarr_t("vk_cart_ibz", "dp", "three, nb, gstore_nkibz")))
vk_cart_ibz = ncr.read_value("vk_cart_ibz", path=path)

if ncr.with_vk == 2:
# Full (nb x nb) matrix.
# Have to transpose (nb_kq, nb_k) submatrix written by Fortran.
# NCF_CHECK(nctk_def_arrays(spin_ncid, nctkarr_t("vkmat_cart_ibz", "dp", "two, three, nb, nb, gstore_nkibz")))
vkmat_cart_ibz = ncr.read_value("vkmat_cart_ibz", path=path).transpose(0, 1, 3, 2, 4).copy()
vkmat_cart_ibz = vkmat_cart_ibz[...,0] + 1j*vkmat_cart_ibz[...,1]
# Full (nb x nb) matrix.
# Have to transpose (nb_kq, nb_k) submatrix written by Fortran.
# NCF_CHECK(nctk_def_arrays(spin_ncid, nctkarr_t("vkmat_cart_ibz", "dp", "two, three, nb, nb, gstore_nkibz")))
vkmat_cart_ibz = ncr.read_value("vkmat_cart_ibz", path=path).transpose(0, 1, 3, 2, 4).copy()
vkmat_cart_ibz = vkmat_cart_ibz[...,0] + 1j*vkmat_cart_ibz[...,1]

# Note conversion between Fortran and python convention.
bstart = ncr.read_value("bstart", path=path) - 1
Expand Down Expand Up @@ -382,7 +381,6 @@ def path2group(self) -> dict:
return self.rootgrp.groups



class GstoreRobot(Robot, RobotWithEbands):
"""
This robot analyzes the results contained in multiple GSTORE.nc files.
Expand Down Expand Up @@ -413,38 +411,38 @@ def compare(self) -> None:

@staticmethod
def compare_two_gstores(gstore1: GstoreFile, gstore2: GstoreFile):
"""
Helper function to compare two GSTORE files.
"""
aname_list = ["structure", "nsppol", "cplex", "nkbz", "nkibz",
"nqbz", "nqibz", "completed", "kzone", "qzone", "kfilter", "gmode",
"brange_spin", "erange_spin", "glob_spin_nq", "glob_nk_spin",
]

for aname in aname_list:
# Get attributes in gstore first, then in gstore.r, else raise.
if hasattr(gstore1, aname):
val1, val2 = getattr(gstore1, aname), getattr(gstore2, aname)
elif hasattr(gstore1.r, aname):
val1, val2 = getattr(gstore1.r, aname), getattr(gstore2.r, aname)
else:
raise AttributeError(f"Cannot find attribute `{aname=}` neither in gstore not in gstore.r")
"""
Helper function to compare two GSTORE files.
"""
aname_list = ["structure", "nsppol", "cplex", "nkbz", "nkibz",
"nqbz", "nqibz", "completed", "kzone", "qzone", "kfilter", "gmode",
"brange_spin", "erange_spin", "glob_spin_nq", "glob_nk_spin",
]

for aname in aname_list:
# Get attributes in gstore first, then in gstore.r, else raise.
if hasattr(gstore1, aname):
val1, val2 = getattr(gstore1, aname), getattr(gstore2, aname)
elif hasattr(gstore1.r, aname):
val1, val2 = getattr(gstore1.r, aname), getattr(gstore2.r, aname)
else:
raise AttributeError(f"Cannot find attribute `{aname=}` neither in gstore not in gstore.r")

# Now compare val1 and val2 taking into account the type.
if isinstance(val1, (str, int, float, Structure)):
eq = val1 == val2
elif isinstance(val1, np.ndarray):
eq = np.allclose(val1, val2)
else:
raise TypeError(f"Don't know how to handle comparison for type: {type(val1)}")

if not eq:
raise RuntimeError(f"Different values of {aname=}, {val1=}, {val2=}")

# Compare gkq objects for each spin.
for spin in range(gstore1.nsppol):
gqk1, gqk2 = gstore1.gqk_spin[spin], gstore2.gqk_spin[spin]
gqk1.compare(gqk2)
if isinstance(val1, (str, int, float, Structure)):
eq = val1 == val2
elif isinstance(val1, np.ndarray):
eq = np.allclose(val1, val2)
else:
raise TypeError(f"Don't know how to handle comparison for type: {type(val1)}")

if not eq:
raise RuntimeError(f"Different values of {aname=}, {val1=}, {val2=}")

# Compare gkq objects for each spin.
for spin in range(gstore1.nsppol):
gqk1, gqk2 = gstore1.gqk_spin[spin], gstore2.gqk_spin[spin]
gqk1.compare(gqk2)

def yield_figs(self, **kwargs): # pragma: no cover
"""
Expand Down
14 changes: 8 additions & 6 deletions abipy/flowtk/qutils.py
Original file line number Diff line number Diff line change
Expand Up @@ -308,11 +308,12 @@ def get_slurm_template(body) -> str:
"""
return """\
#!/bin/bash
# Please customize this section using your settings.
#
# Please CUSTOMIZE this section according to your cluster and the type of calculation
#
#SBATCH --job-name=my_job
#SBATCH --output=%j_%x.out
#SBATCH --error=%j_%x.err
#SBATCH --output=%j_%x.slurm.out
#SBATCH --error=%j_%x.slurm.err
#SBATCH --partition=debug
#SBATCH --nodes=1
#SBATCH --ntasks-per-node=64
Expand All @@ -337,11 +338,10 @@ def get_slurm_template(body) -> str:
echo "number of mpi tasks: $SLURM_NTASKS tasks"
echo "OMP_NUM_THREADS : $OMP_NUM_THREADS"
echo "number of gpus : $SLURM_GPUS_ON_NODE"
echo "------------------- Node list ------------------"
echo $SLURM_JOB_NODELIST
echo "---------------- Checking limits ---------------"
echo "---------------- Printing limits ---------------"
ulimit -s unlimited
ulimit -a
Expand All @@ -350,6 +350,7 @@ def get_slurm_template(body) -> str:
# ------------------------------------------------------------------------------
echo "----------------- Environment ------------------"
source $HOME/vasp.6.2.1/modules.sh
module list
Expand All @@ -359,6 +360,7 @@ def get_slurm_template(body) -> str:
{body}
echo -n "This run completed on: "
date
"""

Expand Down
44 changes: 15 additions & 29 deletions abipy/ml/extxyz_generator.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,14 +3,14 @@
from __future__ import annotations

import os
import numpy as np
#import numpy as np
import abipy.core.abinit_units as abu
#try:
# import ase
#except ImportError as exc:
# raise ImportError("ase not installed. Try `pip install ase`.") from exc
from pathlib import Path
from ase.atoms import Atoms
#from ase.atoms import Atoms
from ase.calculators.singlepoint import SinglePointCalculator
from monty.string import list_strings # marquee,
#from monty.functools import lazy_property
Expand All @@ -22,7 +22,6 @@
#from ase.calculators.calculator import Calculator
#from ase.io.vasp import write_vasp_xdatcar, write_vasp
from ase.stress import full_3x3_to_voigt_6_stress
from ase.calculators.singlepoint import SinglePointCalculator
from ase.io import write
from pymatgen.io.vasp.outputs import Vasprun, Outcar
from abipy.core import Structure
Expand Down Expand Up @@ -64,9 +63,9 @@ class ExtxyzIOWriter:

@classmethod
def from_top(cls, top: PathLike, ext: str):
from monty.os.path import find_exts
filepaths = find_exts(str(top), ext)
return cls(filepaths)
from monty.os.path import find_exts
filepaths = find_exts(str(top), ext)
return cls(filepaths)

def __init__(self, filepaths: list[PathLike]):
self.filepaths = list_strings(filepaths)
Expand All @@ -90,7 +89,7 @@ def to_string(self, verbose=0) -> str:
def __str__(self) -> str:
return self.to_string()

def write(self, xyz_filepath: PathLike, overwrite: bool=False):
def write(self, xyz_filepath: PathLike, overwrite: bool = False):
"""
"""
if not overwrite and os.path.isfile(xyz_filepath):
Expand Down Expand Up @@ -134,7 +133,7 @@ def yield_atoms(self):
yield atoms


def check_vasp_success(vasprun, outcar, verbose: int = 1) -> bool:
def check_vasp_success(vasprun: Vasprun, outcar: Outcar, verbose: int = 1) -> bool:
"""
Check if a VASP calculation completed successfully.
Expand Down Expand Up @@ -223,22 +222,7 @@ def to_string(self, verbose=0) -> str:

return "\n".join(lines)

#def get_custom_incar_settings():
# # Define custom INCAR settings
# custom_incar_settings = {
# 'ENCUT': 600, # Override plane-wave cutoff energy
# 'ISMEAR': -5, # Use tetrahedron method with Blöchl corrections
# 'SIGMA': 0.01, # Smearing width
# 'EDIFF': 1E-8, # Electronic energy convergence criterion
# 'NSW': 0, # Number of ionic steps (static calculation)
# 'ISIF': 2, # Stress tensor calculation
# 'LREAL': 'Auto', # Projection operators (automatic)
# 'LWAVE': False, # Do not write WAVECAR
# 'LCHARG': True # Write CHGCAR
# }
# return custom_incar_settings

def sbatch(self, max_jobs: int=100) -> list[int]:
def sbatch(self, max_jobs: int = 100) -> list[int]:
"""
"""
if not self.topdir.exists(): self.topdir.mkdir()
Expand All @@ -247,7 +231,7 @@ def sbatch(self, max_jobs: int=100) -> list[int]:
for index in self.traj_range:
workdir = self.topdir / f"SINGLEPOINT_{index}"
if workdir.exists():
print(f"{workdir=} already exists. Ignoring it")
print(f"{str(workdir)} already exists. Ignoring it")
continue

try:
Expand All @@ -262,9 +246,12 @@ def sbatch(self, max_jobs: int=100) -> list[int]:
if self.code == "vasp":
# Generate VASP input files using the Materials Project settings for a single-point calculation
from pymatgen.io.vasp.sets import MPStaticSet
user_incar_settings = {"NCORE": 2}
user_incar_settings = {
"NCORE": 2
'LWAVE': False, # Do not write WAVECAR
'LCHARG': False, # Do not Write CHGCAR
}
vasp_input_set = MPStaticSet(structure, user_incar_settings=user_incar_settings)

vasp_input_set.write_input(workdir)
with open(workdir / "run_custodian.py", "wt") as fh:
fh.write(self.custodian_script_str)
Expand All @@ -282,7 +269,7 @@ def sbatch(self, max_jobs: int=100) -> list[int]:

return job_ids

def write_xyz(self, xyz_filepath: PathLike, dryrun=False) -> None:
def write_xyz(self, xyz_filepath: PathLike, dry_run=False) -> None:
"""
"""
ext = {
Expand All @@ -292,4 +279,3 @@ def write_xyz(self, xyz_filepath: PathLike, dryrun=False) -> None:

writer = ExtxyzIOWriter.from_top(self.topdir, ext)
writer.write(xyz_filepath)

2 changes: 1 addition & 1 deletion abipy/scripts/abistruct.py
Original file line number Diff line number Diff line change
Expand Up @@ -199,7 +199,7 @@ def add_primitive_options(parser):
group.add_argument('--primitive-standard', default=False, action='store_true',
help="Enforce primitive standard cell.")

supported_formats = "(abivars, cif, xsf, poscar, qe, siesta, wannier90, cssr, json, None)"
supported_formats = "(abivars, cif, xsf, poscar, qe, siesta, wannier90, cssr, json, lammps, fleur-inpgen, None)"

def add_format_arg(parser, default, option=True, formats=None):
"""
Expand Down

0 comments on commit 3e5f1d1

Please sign in to comment.