diff --git a/pymatgen/analysis/gb/grain.py b/pymatgen/analysis/gb/grain.py index 71613ea8ea0..22692b82e5e 100644 --- a/pymatgen/analysis/gb/grain.py +++ b/pymatgen/analysis/gb/grain.py @@ -67,6 +67,7 @@ def __init__( oriented_unit_cell: Structure, validate_proximity: bool = False, coords_are_cartesian: bool = False, + properties: dict | None = None, ) -> None: """ Makes a GB structure, a structure object with additional information @@ -109,6 +110,8 @@ def __init__( that are less than 0.01 Ang apart. Defaults to False. coords_are_cartesian (bool): Set to True if you are providing coordinates in Cartesian coordinates. Defaults to False. + properties (dict): dictionary containing properties associated + with the whole GrainBoundary. """ self.oriented_unit_cell = oriented_unit_cell self.rotation_axis = rotation_axis @@ -125,6 +128,7 @@ def __init__( validate_proximity=validate_proximity, coords_are_cartesian=coords_are_cartesian, site_properties=site_properties, + properties=properties, ) def copy(self): diff --git a/pymatgen/analysis/magnetism/analyzer.py b/pymatgen/analysis/magnetism/analyzer.py index 8ab84d75256..30171a07b3b 100644 --- a/pymatgen/analysis/magnetism/analyzer.py +++ b/pymatgen/analysis/magnetism/analyzer.py @@ -553,16 +553,16 @@ def matches_ordering(self, other: Structure) -> bool: b_positive = CollinearMagneticStructureAnalyzer(other, overwrite_magmom_mode="normalize", make_primitive=False) b_negative = b_positive.structure.copy() - b_negative.add_site_property("magmom", np.multiply(-1, b_negative.site_properties["magmom"])) + b_negative.add_site_property("magmom", -np.array(b_negative.site_properties["magmom"])) - b_negative = CollinearMagneticStructureAnalyzer( + analyzer = CollinearMagneticStructureAnalyzer( b_negative, overwrite_magmom_mode="normalize", make_primitive=False ) b_positive = b_positive.get_structure_with_spin() - b_negative = b_negative.get_structure_with_spin() + analyzer = analyzer.get_structure_with_spin() - return a.matches(b_positive) or a.matches(b_negative) + return a.matches(b_positive) or a.matches(analyzer) def __str__(self): """ diff --git a/pymatgen/core/structure.py b/pymatgen/core/structure.py index 983305296b9..9e9875b736a 100644 --- a/pymatgen/core/structure.py +++ b/pymatgen/core/structure.py @@ -455,7 +455,7 @@ def from_file(cls, filename: str): """Reads in SiteCollection from a filename.""" raise NotImplementedError - def add_site_property(self, property_name: str, values: list): + def add_site_property(self, property_name: str, values: Sequence | ArrayLike): """Adds a property to a site. Note: This is the preferred method for adding magnetic moments, selective dynamics, and related site-specific properties to a structure/molecule object. @@ -593,12 +593,11 @@ def add_spin_by_element(self, spins: dict[str, float]) -> None: new_species[species] = occu site.species = Composition(new_species) - def add_spin_by_site(self, spins: list[float]) -> None: + def add_spin_by_site(self, spins: Sequence[float]) -> None: """Add spin states to structure by site. Args: - spins (list): List of spins - E.g., [+5, -5, 0, 0] + spins (list): e.g. [+5, -5, 0, 0] """ if len(spins) != len(self): raise ValueError(f"Spins for all sites must be specified, expected {len(self)} spins, got {len(spins)}") @@ -848,6 +847,7 @@ def __init__( coords_are_cartesian: bool = False, site_properties: dict | None = None, labels: Sequence[str | None] | None = None, + properties: dict | None = None, ) -> None: """Create a periodic structure. @@ -886,6 +886,9 @@ def __init__( list of strings, e.g. ['Li1', 'Li2']. Must have the same length as the species and fractional coords. Defaults to None for no labels. + properties (dict): Properties associated with the whole structure. + Will be serialized when writing the structure to JSON or YAML but is + lost when converting to other formats. """ if len(species) != len(coords): raise StructureError("atomic species and fractional coordinates must have same length") @@ -914,6 +917,7 @@ def __init__( if validate_proximity and not self.is_valid(): raise StructureError("Structure contains sites that are less than 0.01 Angstrom apart!") self._charge = charge + self.properties = properties or {} @classmethod def from_sites( @@ -922,6 +926,7 @@ def from_sites( charge: float | None = None, validate_proximity: bool = False, to_unit_cell: bool = False, + properties: dict | None = None, ) -> IStructure: """Convenience constructor to make a Structure from a list of sites. @@ -933,6 +938,9 @@ def from_sites( that are less than 0.01 Ang apart. Defaults to False. to_unit_cell (bool): Whether to translate sites into the unit cell. + properties (dict): Properties associated with the whole structure. + Will be serialized when writing the structure to JSON or YAML but is + lost when converting to other formats. Returns: (Structure) Note that missing properties are set as None. @@ -963,6 +971,7 @@ def from_sites( validate_proximity=validate_proximity, to_unit_cell=to_unit_cell, labels=labels, + properties=properties, ) @classmethod @@ -1982,7 +1991,7 @@ def get_sorted_structure(self, key: Callable | None = None, reverse: bool = Fals as if each comparison were reversed. """ sites = sorted(self, key=key, reverse=reverse) - return type(self).from_sites(sites, charge=self._charge) + return type(self).from_sites(sites, charge=self._charge, properties=self.properties) def get_reduced_structure(self, reduction_algo: Literal["niggli", "LLL"] = "niggli") -> IStructure | Structure: """Get a reduced structure. @@ -2011,7 +2020,7 @@ def get_reduced_structure(self, reduction_algo: Literal["niggli", "LLL"] = "nigg ) return self.copy() - def copy(self, site_properties=None, sanitize=False) -> Structure: + def copy(self, site_properties=None, sanitize=False, properties=None) -> Structure: """Convenience method to get a copy of the structure, with options to add site properties. @@ -2027,29 +2036,34 @@ def copy(self, site_properties=None, sanitize=False) -> Structure: carried out to obtain a relatively orthogonalized cell, (iii) all fractional coords for sites are mapped into the unit cell. + properties (dict): General properties to add or override. Returns: A copy of the Structure, with optionally new site_properties and optionally sanitized. """ - props = self.site_properties + new_site_props = self.site_properties if site_properties: - props.update(site_properties) + new_site_props.update(site_properties) + props = self.properties + if properties: + props.update(properties) if not sanitize: return self.__class__( self._lattice, self.species_and_occu, self.frac_coords, charge=self._charge, - site_properties=props, + site_properties=new_site_props, labels=self.labels, + properties=props, ) reduced_latt = self._lattice.get_lll_reduced_lattice() new_sites = [] for idx, site in enumerate(self): frac_coords = reduced_latt.get_fractional_coords(site.coords) site_props = {} - for prop, val in props.items(): + for prop, val in new_site_props.items(): site_props[prop] = val[idx] new_sites.append( PeriodicSite( @@ -2063,7 +2077,7 @@ def copy(self, site_properties=None, sanitize=False) -> Structure: ) ) new_sites = sorted(new_sites) - return type(self).from_sites(new_sites, charge=self._charge) + return type(self).from_sites(new_sites, charge=self._charge, properties=props) def interpolate( self, @@ -2519,6 +2533,7 @@ def as_dict(self, verbosity=1, fmt=None, **kwargs) -> dict[str, Any]: "@class": type(self).__name__, "charge": self.charge, "lattice": latt_dict, + "properties": self.properties, } for site in self: site_dict = site.as_dict(verbosity=verbosity) # type: ignore[call-arg] @@ -2572,7 +2587,7 @@ def from_dict(cls, d: dict[str, Any], fmt: Literal["abivars"] | None = None) -> lattice = Lattice.from_dict(d["lattice"]) sites = [PeriodicSite.from_dict(sd, lattice) for sd in d["sites"]] charge = d.get("charge") - return cls.from_sites(sites, charge=charge) + return cls.from_sites(sites, charge=charge, properties=d.get("properties")) def to(self, filename: str = "", fmt: str = "", **kwargs) -> str: """Outputs the structure to a file or string. @@ -2613,11 +2628,11 @@ def to(self, filename: str = "", fmt: str = "", **kwargs) -> str: writer = Cssr(self) # type: ignore elif fmt == "json" or fnmatch(filename.lower(), "*.json*"): - dct = json.dumps(self.as_dict()) + json_str = json.dumps(self.as_dict()) if filename: with zopen(filename, "wt") as file: - file.write(dct) - return dct + file.write(json_str) + return json_str elif fmt == "xsf" or fnmatch(filename.lower(), "*.xsf*"): from pymatgen.io.xcrysden import XSF @@ -2752,7 +2767,7 @@ def from_str( struct = struct.get_sorted_structure() if merge_tol: struct.merge_sites(merge_tol) - return cls.from_sites(struct) + return cls.from_sites(struct, properties=struct.properties) @classmethod def from_file(cls, filename, primitive=False, sort=False, merge_tol=0.0, **kwargs) -> Structure | IStructure: @@ -2853,6 +2868,7 @@ def __init__( site_properties: dict | None = None, labels: Sequence[str | None] | None = None, charge_spin_check: bool = True, + properties: dict | None = None, ) -> None: """Create a Molecule. @@ -2880,6 +2896,8 @@ def __init__( charge_spin_check (bool): Whether to check that the charge and spin multiplicity are compatible with each other. Defaults to True. + properties (dict): dictionary containing properties associated + with the whole molecule. """ if len(species) != len(coords): raise StructureError( @@ -2912,6 +2930,7 @@ def __init__( self._spin_multiplicity = spin_multiplicity else: self._spin_multiplicity = 1 if n_electrons % 2 == 0 else 2 + self.properties = properties or {} @property def charge(self) -> float: @@ -2951,7 +2970,7 @@ def copy(self) -> IMolecule | Molecule: Returns: IMolecule | Molecule """ - return type(self).from_sites(self) + return type(self).from_sites(self, properties=self.properties) @classmethod def from_sites( @@ -2961,6 +2980,7 @@ def from_sites( spin_multiplicity: int | None = None, validate_proximity: bool = False, charge_spin_check: bool = True, + properties: dict | None = None, ) -> IMolecule | Molecule: """Convenience constructor to make a Molecule from a list of sites. @@ -2974,6 +2994,8 @@ def from_sites( charge_spin_check (bool): Whether to check that the charge and spin multiplicity are compatible with each other. Defaults to True. + properties (dict): dictionary containing properties associated + with the whole molecule. """ props = collections.defaultdict(list) for site in sites: @@ -2989,6 +3011,7 @@ def from_sites( site_properties=props, labels=labels, charge_spin_check=charge_spin_check, + properties=properties, ) def break_bond(self, ind1: int, ind2: int, tol: float = 0.2) -> tuple[IMolecule | Molecule, ...]: @@ -3126,6 +3149,7 @@ def as_dict(self): "charge": self.charge, "spin_multiplicity": self.spin_multiplicity, "sites": [], + "properties": self.properties, } for site in self: site_dict = site.as_dict() @@ -3148,7 +3172,8 @@ def from_dict(cls, d) -> IMolecule | Molecule: sites = [Site.from_dict(sd) for sd in d["sites"]] charge = d.get("charge", 0) spin_multiplicity = d.get("spin_multiplicity") - return cls.from_sites(sites, charge=charge, spin_multiplicity=spin_multiplicity) + properties = d.get("properties") + return cls.from_sites(sites, charge=charge, spin_multiplicity=spin_multiplicity, properties=properties) def get_distance(self, i: int, j: int) -> float: """Get distance between site i and j. @@ -3343,6 +3368,7 @@ def get_centered_molecule(self) -> IMolecule | Molecule: site_properties=self.site_properties, charge_spin_check=self._charge_spin_check, labels=self.labels, + properties=self.properties, ) def to(self, filename: str = "", fmt: str = "") -> str | None: @@ -3378,7 +3404,7 @@ def to(self, filename: str = "", fmt: str = "") -> str | None: with zopen(filename, "wt", encoding="utf8") as file: file.write(json_str) return json_str - elif fmt == "yaml" or fnmatch(filename, "*.yaml*"): + elif fmt in ("yaml", "yml") or fnmatch(filename, "*.yaml*") or fnmatch(filename, "*.yml*"): yaml = YAML() str_io = StringIO() yaml.dump(self.as_dict(), str_io) @@ -3418,14 +3444,15 @@ def from_str( from pymatgen.io.gaussian import GaussianInput from pymatgen.io.xyz import XYZ - if fmt.lower() == "xyz": + fmt = fmt.lower() # type: ignore[assignment] + if fmt == "xyz": mol = XYZ.from_str(input_string).molecule elif fmt in ["gjf", "g03", "g09", "com", "inp"]: mol = GaussianInput.from_str(input_string).molecule elif fmt == "json": dct = json.loads(input_string) return cls.from_dict(dct) - elif fmt == "yaml": + elif fmt in ("yaml", "yml"): yaml = YAML() dct = yaml.load(input_string) return cls.from_dict(dct) @@ -3433,7 +3460,7 @@ def from_str( from pymatgen.io.babel import BabelMolAdaptor mol = BabelMolAdaptor.from_str(input_string, file_format=fmt).pymatgen_mol - return cls.from_sites(mol) + return cls.from_sites(mol, properties=mol.properties) @classmethod def from_file(cls, filename): @@ -3463,7 +3490,7 @@ def from_file(cls, filename): return GaussianOutput(filename).final_structure if fnmatch(fname, "*.json*") or fnmatch(fname, "*.mson*"): return cls.from_str(contents, fmt="json") - if fnmatch(fname, "*.yaml*"): + if fnmatch(fname, "*.yaml*") or fnmatch(filename, "*.yml*"): return cls.from_str(contents, fmt="yaml") from pymatgen.io.babel import BabelMolAdaptor @@ -3491,6 +3518,7 @@ def __init__( coords_are_cartesian: bool = False, site_properties: dict | None = None, labels: Sequence[str | None] | None = None, + properties: dict | None = None, ): """Create a periodic structure. @@ -3528,6 +3556,9 @@ def __init__( list of strings, e.g. ['Li1', 'Li2']. Must have the same length as the species and fractional coords. Defaults to None for no labels. + properties (dict): Properties associated with the whole structure. + Will be serialized when writing the structure to JSON or YAML but is + lost when converting to other formats. """ super().__init__( lattice, @@ -3539,6 +3570,7 @@ def __init__( coords_are_cartesian=coords_are_cartesian, site_properties=site_properties, labels=labels, + properties=properties, ) self._sites: list[PeriodicSite] = list(self._sites) # type: ignore @@ -4261,6 +4293,7 @@ def __init__( site_properties: dict | None = None, labels: Sequence[str | None] | None = None, charge_spin_check: bool = True, + properties: dict | None = None, ) -> None: """Creates a MutableMolecule. @@ -4288,6 +4321,8 @@ def __init__( charge_spin_check (bool): Whether to check that the charge and spin multiplicity are compatible with each other. Defaults to True. + properties (dict): dictionary containing properties associated + with the whole molecule. """ super().__init__( species, @@ -4298,6 +4333,7 @@ def __init__( site_properties=site_properties, labels=labels, charge_spin_check=charge_spin_check, + properties=properties, ) self._sites: list[Site] = list(self._sites) # type: ignore diff --git a/pymatgen/core/surface.py b/pymatgen/core/surface.py index aae750e9679..752949b0b0b 100644 --- a/pymatgen/core/surface.py +++ b/pymatgen/core/surface.py @@ -95,6 +95,7 @@ def __init__( coords_are_cartesian=False, site_properties=None, energy=None, + properties=None, ): """Makes a Slab structure, a structure object with additional information and methods pertaining to slabs. @@ -141,6 +142,8 @@ def __init__( have to be the same length as the atomic species and fractional_coords. Defaults to None for no properties. energy (float): A value for the energy. + properties (dict): dictionary containing properties associated + with the whole slab. """ self.oriented_unit_cell = oriented_unit_cell self.miller_index = tuple(miller_index) @@ -514,6 +517,7 @@ def from_dict(cls, d): scale_factor=d["scale_factor"], site_properties=struct.site_properties, energy=d["energy"], + properties=d.get("properties"), ) def get_surface_sites(self, tag=False): diff --git a/pymatgen/io/babel.py b/pymatgen/io/babel.py index 1b41111ec89..aa6270f4003 100644 --- a/pymatgen/io/babel.py +++ b/pymatgen/io/babel.py @@ -29,18 +29,21 @@ __date__ = "Apr 28, 2012" +needs_openbabel = requires( + openbabel, + "BabelMolAdaptor requires openbabel to be installed with Python bindings. " + "Please get it at http://openbabel.org (version >=3.0.0).", +) + + class BabelMolAdaptor: """ Adaptor serves as a bridge between OpenBabel's Molecule and pymatgen's Molecule. """ - @requires( - openbabel, - "BabelMolAdaptor requires openbabel to be installed with Python bindings. " - "Please get it at http://openbabel.org (version >=3.0.0).", - ) - def __init__(self, mol): + @needs_openbabel + def __init__(self, mol: Molecule | openbabel.OBMol | pybel.Molecule) -> None: """ Initializes with pymatgen Molecule or OpenBabel's OBMol. @@ -342,6 +345,7 @@ def from_molecule_graph(mol): def from_string(cls, *args, **kwargs): return cls.from_str(*args, **kwargs) + @needs_openbabel @staticmethod def from_str(string_data, file_format="xyz"): """ diff --git a/pymatgen/io/vasp/inputs.py b/pymatgen/io/vasp/inputs.py index d8ec7ecebab..04794df19f6 100644 --- a/pymatgen/io/vasp/inputs.py +++ b/pymatgen/io/vasp/inputs.py @@ -127,7 +127,8 @@ def __init__( self.structure = self.structure.get_sorted_structure() self.true_names = true_names self.comment = structure.formula if comment is None else comment - self.predictor_corrector_preamble = predictor_corrector_preamble + if predictor_corrector_preamble: + self.structure.properties["predictor_corrector_preamble"] = predictor_corrector_preamble else: raise ValueError("Structure with partial occupancies cannot be converted into POSCAR!") @@ -148,6 +149,11 @@ def predictor_corrector(self): """Predictor corrector in Poscar.""" return self.structure.site_properties.get("predictor_corrector") + @property + def predictor_corrector_preamble(self): + """Predictor corrector preamble in Poscar.""" + return self.structure.properties.get("predictor_corrector_preamble") + @velocities.setter # type: ignore def velocities(self, velocities): """Setter for Poscar.velocities.""" @@ -163,6 +169,11 @@ def predictor_corrector(self, predictor_corrector): """Setter for Poscar.predictor_corrector.""" self.structure.add_site_property("predictor_corrector", predictor_corrector) + @predictor_corrector_preamble.setter # type: ignore + def predictor_corrector_preamble(self, predictor_corrector_preamble): + """Setter for Poscar.predictor_corrector.""" + self.structure.properties["predictor_corrector"] = predictor_corrector_preamble + @property def site_symbols(self): """ diff --git a/pymatgen/symmetry/structure.py b/pymatgen/symmetry/structure.py index c84fceee593..62bd413695c 100644 --- a/pymatgen/symmetry/structure.py +++ b/pymatgen/symmetry/structure.py @@ -47,6 +47,7 @@ def __init__( [site.species for site in structure], structure.frac_coords, site_properties=structure.site_properties, + properties=structure.properties, ) equivalent_indices: list[list[int]] = [[] for _ in range(len(uniq))] diff --git a/tests/core/test_structure.py b/tests/core/test_structure.py index c55f00666e9..e27232c6cfb 100644 --- a/tests/core/test_structure.py +++ b/tests/core/test_structure.py @@ -74,7 +74,9 @@ def setUp(self): coords = [[0, 0, 0], [0.0, 0, 0.0000001]] with pytest.raises(StructureError, match="Structure contains sites that are less than 0.01 Angstrom apart"): IStructure(self.lattice, ["Si"] * 2, coords, validate_proximity=True) - self.propertied_structure = IStructure(self.lattice, ["Si"] * 2, coords, site_properties={"magmom": [5, -5]}) + self.propertied_structure = IStructure( + self.lattice, ["Si"] * 2, coords, site_properties={"magmom": [5, -5]}, properties={"test_property": "test"} + ) self.labeled_structure = IStructure(self.lattice, ["Si"] * 2, coords, labels=["Si1", "Si2"]) self.lattice_pbc = Lattice( @@ -216,10 +218,12 @@ def test_as_dict(self): ], coords, site_properties={"magmom": [5, -5]}, + properties={"general_property": "test"}, ) d = struct.as_dict() assert d["sites"][0]["properties"]["magmom"] == 5 assert d["sites"][0]["species"][0]["spin"] == 3 + assert d["properties"]["general_property"] == "test" d = struct.as_dict(0) assert "volume" not in d["lattice"] @@ -280,10 +284,12 @@ def test_from_dict(self): "xyz": [3.8401979336749994, 1.2247250003039056e-06, 2.351631795225], }, ], + "properties": {"test_property": "test"}, } struct = IStructure.from_dict(d) assert struct[0].magmom == 5 assert struct[0].specie.spin == 3 + assert struct.properties["test_property"] == "test" assert isinstance(struct, IStructure) def test_site_properties(self): @@ -292,25 +298,36 @@ def test_site_properties(self): assert self.propertied_structure[0].magmom == 5 assert self.propertied_structure[1].magmom == -5 + def test_properties_dict(self): + assert self.propertied_structure.properties == {"test_property": "test"} + def test_copy(self): - new_struct = self.propertied_structure.copy(site_properties={"charge": [2, 3]}) + new_struct = self.propertied_structure.copy( + site_properties={"charge": [2, 3]}, properties={"another_prop": "test"} + ) assert new_struct[0].magmom == 5 assert new_struct[1].magmom == -5 assert new_struct[0].charge == 2 assert new_struct[1].charge == 3 + assert new_struct.properties["another_prop"] == "test" + assert new_struct.properties["test_property"] == "test" coords = [] coords.append([0, 0, 0]) coords.append([0.0, 0, 0.0000001]) - struct = IStructure(self.lattice, ["O", "Si"], coords, site_properties={"magmom": [5, -5]}) + struct = IStructure( + self.lattice, ["O", "Si"], coords, site_properties={"magmom": [5, -5]}, properties={"test_property": "test"} + ) - new_struct = struct.copy(site_properties={"charge": [2, 3]}, sanitize=True) + new_struct = struct.copy(site_properties={"charge": [2, 3]}, sanitize=True, properties={"another_prop": "test"}) assert new_struct[0].magmom == -5 assert new_struct[1].magmom == 5 assert new_struct[0].charge == 3 assert new_struct[1].charge == 2 assert new_struct.volume == approx(struct.volume) + assert new_struct.properties["another_prop"] == "test" + assert new_struct.properties["test_property"] == "test" def test_interpolate(self): coords = [] @@ -925,6 +942,26 @@ def test_propertied_structure(self): self.struct.append("Li", [0.3, 0.3, 0.3]) assert len(self.struct.site_properties["charge"]) == 3 + props = {"test_property": 42} + struct_with_props = Structure( + lattice=Lattice.cubic(3), + species=("Fe", "Fe"), + coords=((0, 0, 0), (0.5, 0.5, 0.5)), + site_properties={"magmom": [5, -5]}, + properties=props, + ) + + dct = struct_with_props.as_dict() + struct = Structure.from_dict(dct) + assert struct.properties == props + assert dct == struct.as_dict() + + json_str = struct_with_props.to(fmt="json") + assert '"test_property": 42' in json_str + struct = Structure.from_str(json_str, fmt="json") + assert struct.properties == props + assert dct == struct.as_dict() + def test_perturb(self): dist = 0.1 pre_perturbation_sites = self.struct.copy() @@ -1657,6 +1694,11 @@ def test_properties(self): assert self.mol.is_ordered assert self.mol.formula == "H4 C1" + def test_properties_dict(self): + properties = {"test_property": "test"} + self.mol.properties = properties + assert self.mol.properties == properties + def test_repr_str(self): expected = """Full Formula (H4 C1) Reduced Formula: H4C @@ -1821,6 +1863,7 @@ def test_as_from_dict(self): self.coords, charge=1, site_properties={"magmom": [0.5, -0.5, 1, 2, 3]}, + properties={"test_properties": "test"}, ) dct = propertied_mol.as_dict() assert dct["sites"][0]["properties"]["magmom"] == 0.5 @@ -1829,6 +1872,7 @@ def test_as_from_dict(self): assert mol[0].magmom == 0.5 assert mol.formula == "H4 C1" assert mol.charge == 1 + assert mol.properties == {"test_properties": "test"} def test_default_dict_attrs(self): d = self.mol.as_dict() @@ -1836,12 +1880,15 @@ def test_default_dict_attrs(self): assert d["spin_multiplicity"] == 1 def test_to_from_file_string(self): - for fmt in ["xyz", "json", "g03", "yaml"]: + self.mol.properties["test_prop"] = 42 + for fmt in ("xyz", "json", "g03", "yaml", "yml"): mol = self.mol.to(fmt=fmt) assert isinstance(mol, str) mol = IMolecule.from_str(mol, fmt=fmt) assert mol == self.mol assert isinstance(mol, IMolecule) + if fmt in ("json", "yaml", "yml"): + assert mol.properties.get("test_prop") == 42 ch4_xyz_str = self.mol.to(filename=f"{self.tmp_path}/CH4_testing.xyz") with open("CH4_testing.xyz") as xyz_file: diff --git a/tests/io/vasp/test_sets.py b/tests/io/vasp/test_sets.py index a228f2c74c2..b1012b15469 100644 --- a/tests/io/vasp/test_sets.py +++ b/tests/io/vasp/test_sets.py @@ -589,7 +589,7 @@ def test_valid_magmom_struct(self): struct.insert(0, "Li", [0, 0, 0]) vis = MPRelaxSet(struct, user_potcar_settings={"Fe": "Fe"}, validate_magmom=False) - with pytest.raises(TypeError, match=r"float\(\) argument must be a string or a real number"): + with pytest.raises(TypeError, match="argument must be a string or a real number, not 'NoneType'"): print(vis.get_vasp_input()) vis = MPRelaxSet(struct, user_potcar_settings={"Fe": "Fe"}, validate_magmom=True)