Skip to content

Commit 75b94f1

Browse files
authored
Merge pull request #226 from ICAMS/merge_pd_data
Merge pd data
2 parents 85abd6f + fc8bc24 commit 75b94f1

File tree

5 files changed

+67
-19
lines changed

5 files changed

+67
-19
lines changed

.bumpversion.cfg

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,5 @@
11
[bumpversion]
2-
current_version = 1.6.2
2+
current_version = 1.6.3
33
commit = True
44
tag = True
55

calphy/__init__.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -4,7 +4,7 @@
44
from calphy.alchemy import Alchemy
55
from calphy.routines import MeltingTemp
66

7-
__version__ = "1.6.2"
7+
__version__ = "1.6.3"
88

99
def addtest(a,b):
1010
return a+b

calphy/input.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -49,7 +49,7 @@
4949
from ase.io import read, write
5050
import shutil
5151

52-
__version__ = "1.6.2"
52+
__version__ = "1.6.3"
5353

5454

5555
def _check_equal(val):

calphy/phase_diagram.py

Lines changed: 63 additions & 15 deletions
Original file line numberDiff line numberDiff line change
@@ -1034,16 +1034,21 @@ def get_free_energy_mixing(dict_list, threshold=1E-3, boundary_trim=0.1):
10341034
right_partial = c_max < max_comp - threshold
10351035

10361036
if left_partial or right_partial:
1037-
trim = float(boundary_trim)
1037+
# Cap trim so it never exceeds 1/4 of the phase's own range,
1038+
# preventing narrow phases (e.g. ordered compounds) from being
1039+
# trimmed to an empty array.
1040+
phase_range = c_max - c_min
1041+
trim = min(float(boundary_trim), phase_range / 4.0)
10381042

10391043
mask = np.ones(len(comp), dtype=bool)
10401044
if left_partial:
10411045
mask &= comp >= c_min + trim
10421046
if right_partial:
10431047
mask &= comp <= c_max - trim
1044-
d["composition"] = comp[mask]
1045-
d["free_energy"] = d["free_energy"][mask]
1046-
d["free_energy_mix"] = d["free_energy_mix"][mask]
1048+
if mask.any():
1049+
d["composition"] = comp[mask]
1050+
d["free_energy"] = d["free_energy"][mask]
1051+
d["free_energy_mix"] = d["free_energy_mix"][mask]
10471052

10481053
return dict_list
10491054

@@ -1085,6 +1090,8 @@ def get_tangent_type(dict_list, tangent, energy):
10851090
right_phases = []
10861091

10871092
for d in dict_list:
1093+
if len(d["composition"]) == 0:
1094+
continue
10881095
diff = np.abs(left_c - d["composition"])
10891096
arg = np.argsort(diff)[0]
10901097
if diff[arg] < 1E-5:
@@ -1315,14 +1322,23 @@ class PhaseDiagram:
13151322
Parameters
13161323
----------
13171324
folders : dict
1318-
Mapping of phase name → folder path, e.g.
1319-
``{'cufcc': 'cufcc', 'agfcc': 'agfcc', 'lqd': 'lqd'}``.
1325+
Mapping of phase name → folder path (or list of folder paths).
1326+
Pass a list to merge multiple simulation folders into a single
1327+
phase — the raw data is combined *before* any fitting, so there
1328+
are no cross-folder inconsistencies::
1329+
1330+
# single folder per phase
1331+
{'cufcc': 'cufcc', 'agfcc': 'agfcc', 'lqd': 'lqd'}
1332+
1333+
# aufcc and cufcc merged into one 'fcc' phase
1334+
{'fcc': ['aufcc_folder', 'cufcc_folder'], 'lqd': 'lqd'}
1335+
13201336
reference_element : str
13211337
The element whose fraction is used as the composition axis
13221338
(e.g. ``'Ag'``).
13231339
composition_intervals : dict, optional
13241340
Per-phase composition bounds, e.g.
1325-
``{'cufcc': (0, 0.5), 'agfcc': (0.7, 1.0)}``.
1341+
``{'fcc': (0, 1), 'lqd': (0, 1)}``.
13261342
Phases not listed are auto-detected from the data: the
13271343
interval is set to the ``(min, max)`` of available
13281344
compositions for that phase.
@@ -1333,8 +1349,8 @@ class PhaseDiagram:
13331349
Examples
13341350
--------
13351351
>>> pd = PhaseDiagram(
1336-
... folders={'cufcc': 'cufcc', 'agfcc': 'agfcc', 'lqd': 'lqd'},
1337-
... reference_element='Ag',
1352+
... folders={'fcc': ['aufcc', 'cufcc'], 'lqd': 'lqd'},
1353+
... reference_element='Au',
13381354
... )
13391355
>>> pd.calculate(T_range=(400, 1400), T_step=5, fit_order=4,
13401356
... method='redlich-kister')
@@ -1351,11 +1367,33 @@ def __init__(self, folders, reference_element,
13511367
self.phases = list(folders.keys())
13521368

13531369
# ---- gather & clean ----
1370+
# Each value in `folders` can be a single path or a list of paths.
1371+
# IMPORTANT: fix_composition_scaling anchors each composition_scaling
1372+
# row to its phase's reference (is_reference=True) row. A single
1373+
# folder may contain calculations with *different* phase_name values
1374+
# (e.g. lqdAg + lqdCu in the same folder), each having its own
1375+
# reference point. We therefore assign a unique temp name per
1376+
# (folder_index, original_phase_name) group so that clean_df and
1377+
# fix_composition_scaling always see exactly one reference per group.
1378+
# Only after both steps are all groups merged under the user's key.
1379+
temp_to_user = {} # temporary_phase_name -> user key
13541380
dfs = []
1355-
for phase, folder in folders.items():
1356-
df = gather_results(folder, reduce_composition=True,
1357-
extract_phase_prefix=True)
1358-
dfs.append(df)
1381+
for phase, folder_spec in folders.items():
1382+
folder_list = ([folder_spec]
1383+
if isinstance(folder_spec, str)
1384+
else list(folder_spec))
1385+
for idx, folder in enumerate(folder_list):
1386+
df = gather_results(folder, reduce_composition=True,
1387+
extract_phase_prefix=True)
1388+
# Assign one temp name per distinct original phase_name so
1389+
# each group (with its own reference row) is processed
1390+
# independently through fix_composition_scaling.
1391+
def _make_temp(orig, phase=phase, idx=idx):
1392+
return f"__calphy_{phase}_{idx}_{orig}__"
1393+
for orig_name in df['phase_name'].unique():
1394+
temp_to_user[_make_temp(orig_name)] = phase
1395+
df['phase_name'] = df['phase_name'].apply(_make_temp)
1396+
dfs.append(df)
13591397

13601398
combined = pd.concat(dfs, ignore_index=True)
13611399
combined = combined.loc[combined.status == 'True']
@@ -1364,10 +1402,20 @@ def __init__(self, folders, reference_element,
13641402
combine_direct_calculations=True, smooth=smooth)
13651403
dfc = fix_composition_scaling(dfc, add_ideal_entropy=True)
13661404

1367-
for key, val in dfc.items():
1405+
# Merge temp-named groups back to user keys and label the phase column.
1406+
user_dfc = {}
1407+
for temp_name, phase_df in dfc.items():
1408+
user_key = temp_to_user[temp_name]
1409+
if user_key in user_dfc:
1410+
user_dfc[user_key] = pd.concat(
1411+
[user_dfc[user_key], phase_df], ignore_index=True)
1412+
else:
1413+
user_dfc[user_key] = phase_df
1414+
1415+
for key, val in user_dfc.items():
13681416
val['phase'] = key
13691417

1370-
self.df = pd.concat([val for val in dfc.values()])
1418+
self.df = pd.concat([val for val in user_dfc.values()])
13711419

13721420
# ---- auto-detect composition intervals from data ----
13731421
for phase in self.phases:

pyproject.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -4,7 +4,7 @@ build-backend = "setuptools.build_meta"
44

55
[project]
66
name = "calphy"
7-
version = "1.6.2"
7+
version = "1.6.3"
88
description = "free energy calculation for python"
99
readme = "README.md"
1010
license = "GPL-3.0-only"

0 commit comments

Comments
 (0)