User questions and feedback can be directed to the pycalphad Google Group. Bugs can be reported to the GitHub repo.
0.5.2
. If this version is not available, please install the development version.Ag-Bi-Cu-Pb-Sb-Sn-nist-solders.tdb
thermodynamic database in this directory.
In [2]:
from pycalphad import Database
solder_dbf = Database('Ag-Bi-Cu-Pb-Sb-Sn-nist-solders.tdb')
phases
attribute that is a dictionary. Print this dictionary. The keys are the phase names and the values are Phase
objects. Each phase object has a constituents
attribute that returns a tuple of the sublattice model for that phase. Print the sublattice model of the BCC_A2
phase and confirm that the sublattice model is ((CU, SN), (VA))
. Note that the sublattices are wrapped as frozensets
. A frozenset
is an immutable version of set
.
In [3]:
print(solder_dbf.phases, end='\n\n')
print('BCC_A2 sublattice model:')
print(solder_dbf.phases['BCC_A2'].constituents)
SN-PB
system.binplot
from pycalphad and use the binplot
function to plot the binary phase diagram of the SN-PB
system described in the database you created. Use the phases ['LIQUID', 'BCT_A5', 'FCC_A1']
. Use conditions P=101325, T=(300, 700, 5), X(SN)=(0,1,0.005)
. Identify each single and two phase region.Remember to use the %matplotlib inline
magic if you are using a Jupyter Notebook!
In [4]:
%matplotlib inline
from pycalphad import binplot, variables as v
comps = ['PB', 'SN', 'VA']
phases = ['LIQUID', 'BCT_A5', 'FCC_A1']
conds = {v.P: 101325, v.T: (300, 700, 10), v.X('SN'): (0, 1, 0.01)}
binplot(solder_dbf, comps, phases, conds)
Out[4]:
calculate
from pycalphad to calculate Gibbs free energies for each of these phases. Use the xarray Dataset returned from the calculate
function and matplotlib to plot the free energy curves as a function of composition. Do this at 400 K, 450 K, 475 K and 600 K. Try to use the pycalphad.plot.utils.phase_legend
function to create a phase legendSee cell 6 in the pycalphad binary examples for an example.
In [5]:
from pycalphad import calculate
from pycalphad.plot.utils import phase_legend
import matplotlib.pyplot as plt
import numpy as np
legend_handles, colorlist = phase_legend(phases)
temperatures = 400, 450, 475, 600
for temp in temperatures:
fig = plt.figure()
ax = fig.gca()
ax.set_title('Pb-Sn {} K'.format(temp))
ax.set_ylabel('GM')
ax.set_xlabel('X(Sn)')
for phase in phases:
result = calculate(solder_dbf, comps, phase, T=temp, P=101325, output='GM')
ax.scatter(result.X.sel(component='SN'), result.GM,
marker='.', s=5, color=colorlist[phase.upper()])
ax.set_xlim((0, 1))
ax.legend(handles=legend_handles, loc='center left', bbox_to_anchor=(1, 0.6))
BI-PB
binary system.filter_active_phases
) to take a pycalphad Database and a list of compositions and return a list of only active phases. An active phase can be defined by a phase that has at least one of the desired components in each sublattice. E.g. a sublattice model of ((BI,SN), (PB,SN))
is valid for the components BI-PB-VA
, but the sublattice model ((BI,SN), (SN))
is not because the second sublattice does not contain any of the components. If you have done this right, running the function should return ['RHOMBO_A7', 'LIQUID', 'HCP_A3', 'FCC_A1', 'BCT_A5']
for BI-PB-VA
components. Hint: use set
methods
In [6]:
def filter_active_phases(dbf, comps):
active_phases = []
for phase_name, phase in dbf.phases.items():
sublattice_is_valid = [] # a list of booleans, e.g. [True, True, True] means that all 3 sublattices are valid
for sublattice in phase.constituents:
valid_sublattice = len(set(comps).intersection(sublattice)) > 0
sublattice_is_valid.append(valid_sublattice)
if all(sublattice_is_valid):
active_phases.append(phase_name)
return active_phases
phases = filter_active_phases(solder_dbf, ['BI', 'PB', 'VA'])
print('Active phase {}'.format(phases))
print('Inactive phases {}'.format(solder_dbf.phases.keys() - set(phases)))
In [7]:
from pycalphad import equilibrium
comps = ['BI', 'PB', 'VA']
phases = filter_active_phases(solder_dbf, comps)
conds = {v.P: 101325, v.T: (300, 700, 10), v.X('PB'): (0, 1, 0.02)}
eq_result = equilibrium(solder_dbf, comps, phases, conds)
print(eq_result)
In [8]:
from pycalphad.plot.eqplot import eqplot
eqplot(eq_result)
Out[8]:
calculate
. We can use calculate
to give us properties for individual phases, but often equilibrium properties are of interest. Use equilibrium
to calculate and plot the equilibrium free energy as a function of composition at 425 K. Bonus: plot the free energies of the phases on top of this plot.
In [9]:
from pycalphad import equilibrium
conds = {v.P: 101325, v.T: 425, v.X('PB'): (0, 1, 0.005)}
eq_result = equilibrium(solder_dbf, comps, phases, conds)
print(eq_result)
legend_handles, colorlist = phase_legend(phases)
fig = plt.figure()
ax = fig.gca()
ax.set_title('BI-PB 425 K')
ax.set_ylabel('GM')
ax.set_xlabel('X(Pb)')
for phase in phases:
result = calculate(solder_dbf, comps, phase, T=425, P=101325, output='GM')
ax.scatter(result.X.sel(component='PB'), result.GM,
marker='.', s=5, color=colorlist[phase.upper()])
# finally plot the equilibrium values
ax.scatter(eq_result.X_PB, eq_result.GM, c='r', lw=0)
ax.set_xlim((0, 1))
ax.legend(handles=legend_handles, loc='center left', bbox_to_anchor=(1, 0.6))
Out[9]:
DO3
phase.((CU, SN), (CU, SN))
with sublattice site ratios of (0.75, 0.25)
. That means there can be four different endmembers (an endmember has no mixing within sublattices) of this phase, pure Cu and pure Sn (((CU), (CU))
and ((SN), (SN))
), CU-25SN
(((CU), (SN))
), and SN-25CU
(((SN), (CU))
). Then with mixing, the entire free energy surface between these endmembers can be calculated. Use calculate
to plot the free energy surface of this phase at 500 K and find the endmembers.
In [11]:
fig = plt.figure()
ax = fig.gca()
ax.set_title('DO3 phase free energy 500 K')
ax.set_ylabel('GM')
ax.set_xlabel('X(Sn)')
result = calculate(solder_dbf, ['CU', 'SN'], 'DO3', T=500, P=101325, output='GM')
ax.scatter(result.X.sel(component='SN'), result.GM,
marker='.', s=5,)
ax.set_xlim((0, 1))
Out[11]:
In [12]:
import xarray as xr
eq_tern = xr.open_dataset('Ag-Bi-Cu-eq.nc')
In [13]:
fig = plt.figure()
ax = eqplot(eq_tern.sel(T=[500]))
ax.set_title('Ag-Bi-Cu 500 K')
fig = plt.figure()
ax = eqplot(eq_tern.sel(T=[600]))
ax.set_title('Ag-Bi-Cu 600 K')
fig = plt.figure()
ax = eqplot(eq_tern.sel(T=[700]))
ax.set_title('Ag-Bi-Cu 700 K')
fig = plt.figure()
ax = eqplot(eq_tern.sel(T=[800]))
ax.set_title('Ag-Bi-Cu 8000 K')
fig = plt.figure()
ax = eqplot(eq_tern.sel(T=[900]))
ax.set_title('Ag-Bi-Cu 900 K')
fig = plt.figure()
ax = eqplot(eq_tern.sel(T=[1000]))
ax.set_title('Ag-Bi-Cu 1000 K')
Out[13]: