def test_get_data_quantities_AL_NI_VA_interaction():
    """Test that an interaction with a VA produces the correct data quantities

    We just have a template database that has the phase defined. We then hot
    patch the Model object to have the GM from the fixed model we printed out
    and the data we printed out. The hot patch is needed because this is
    formation enthalpy data and the model needs to have the lower order terms
    in composition.

    One possible issue is that the new GM in the fixed model does not have any
    individual contributions, so it cannot be used to test excluded model
    contributions. The only excluded model contributions in this data are
    idmix, but the property we are testing is HM_FORM, so the feature transform
    of the idmix property should be zero.

    """
    # Hack the namespace to make the copy-pasted Gibbs energy function work
    from sympy import log, Piecewise
    T = v.T

    data = [{'components': ['AL', 'NI', 'VA'], 'phases': ['BCC_B2'], 'solver': {'mode': 'manual', 'sublattice_occupancies': [[1.0, [0.5, 0.5], 1.0], [1.0, [0.75, 0.25], 1.0]], 'sublattice_site_ratios': [0.5, 0.5, 1.0], 'sublattice_configurations': (('AL', ('NI', 'VA'), 'VA'), ('AL', ('NI', 'VA'), 'VA')), 'comment': 'BCC_B2 sublattice configuration (2SL)'}, 'conditions': {'P': 101325.0, 'T': np.array([300.])}, 'reference_state': 'SGTE91', 'output': 'HM_FORM', 'values': np.array([[[-40316.61077, -56361.58554]]]), 'reference': 'C. Jiang 2009 (constrained SQS)', 'excluded_model_contributions': ['idmix']}, {'components': ['AL', 'NI', 'VA'], 'phases': ['BCC_B2'], 'solver': {'mode': 'manual', 'sublattice_occupancies': [[1.0, [0.5, 0.5], 1.0], [1.0, [0.75, 0.25], 1.0]], 'sublattice_site_ratios': [0.5, 0.5, 1.0], 'sublattice_configurations': (('AL', ('NI', 'VA'), 'VA'), ('AL', ('NI', 'VA'), 'VA')), 'comment': 'BCC_B2 sublattice configuration (2SL)'}, 'conditions': {'P': 101325.0, 'T': np.array([300.])}, 'reference_state': 'SGTE91', 'output': 'HM_FORM', 'values': np.array([[[-41921.43363, -57769.49473]]]), 'reference': 'C. Jiang 2009 (relaxed SQS)', 'excluded_model_contributions': ['idmix']}]
    NEW_GM = 8.3145*T*(0.5*Piecewise((v.SiteFraction("BCC_B2", 0, "AL")*log(v.SiteFraction("BCC_B2", 0, "AL")), v.SiteFraction("BCC_B2", 0, "AL") > 1.0e-16), (0, True))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI")) + 0.5*Piecewise((v.SiteFraction("BCC_B2", 0, "NI")*log(v.SiteFraction("BCC_B2", 0, "NI")), v.SiteFraction("BCC_B2", 0, "NI") > 1.0e-16), (0, True))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI")) + 0.5*Piecewise((v.SiteFraction("BCC_B2", 0, "VA")*log(v.SiteFraction("BCC_B2", 0, "VA")), v.SiteFraction("BCC_B2", 0, "VA") > 1.0e-16), (0, True))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI")) + 0.5*Piecewise((v.SiteFraction("BCC_B2", 1, "AL")*log(v.SiteFraction("BCC_B2", 1, "AL")), v.SiteFraction("BCC_B2", 1, "AL") > 1.0e-16), (0, True))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI")) + 0.5*Piecewise((v.SiteFraction("BCC_B2", 1, "NI")*log(v.SiteFraction("BCC_B2", 1, "NI")), v.SiteFraction("BCC_B2", 1, "NI") > 1.0e-16), (0, True))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI")) + 0.5*Piecewise((v.SiteFraction("BCC_B2", 1, "VA")*log(v.SiteFraction("BCC_B2", 1, "VA")), v.SiteFraction("BCC_B2", 1, "VA") > 1.0e-16), (0, True))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI")) + Piecewise((v.SiteFraction("BCC_B2", 2, "VA")*log(v.SiteFraction("BCC_B2", 2, "VA")), v.SiteFraction("BCC_B2", 2, "VA") > 1.0e-16), (0, True))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI"))) + (45262.9*v.SiteFraction("BCC_B2", 0, "AL")*v.SiteFraction("BCC_B2", 0, "NI")*v.SiteFraction("BCC_B2", 1, "AL")*v.SiteFraction("BCC_B2", 2, "VA") + 45262.9*v.SiteFraction("BCC_B2", 0, "AL")*v.SiteFraction("BCC_B2", 1, "AL")*v.SiteFraction("BCC_B2", 1, "NI")*v.SiteFraction("BCC_B2", 2, "VA"))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI")) + (1.0*v.SiteFraction("BCC_B2", 0, "AL")*v.SiteFraction("BCC_B2", 1, "AL")*v.SiteFraction("BCC_B2", 2, "VA")*Piecewise((10083 - 4.813*T, (T >= 298.15) & (T < 2900.0)), (0, True)) + v.SiteFraction("BCC_B2", 0, "AL")*v.SiteFraction("BCC_B2", 1, "NI")*v.SiteFraction("BCC_B2", 2, "VA")*(9.52839e-8*T**3 + 0.00123463*T**2 + 0.000871898*T*log(T) + 1.31471*T - 64435.3 + 23095.2/T) + v.SiteFraction("BCC_B2", 0, "AL")*v.SiteFraction("BCC_B2", 1, "VA")*v.SiteFraction("BCC_B2", 2, "VA")*(10.0*T + 16432.5) + v.SiteFraction("BCC_B2", 0, "NI")*v.SiteFraction("BCC_B2", 1, "AL")*v.SiteFraction("BCC_B2", 2, "VA")*(9.52839e-8*T**3 + 0.00123463*T**2 + 0.000871898*T*log(T) + 1.31471*T - 64435.3 + 23095.2/T) + 1.0*v.SiteFraction("BCC_B2", 0, "NI")*v.SiteFraction("BCC_B2", 1, "NI")*v.SiteFraction("BCC_B2", 2, "VA")*Piecewise((8715.084 - 3.556*T, (T >= 298.15) & (T < 3000.0)), (0, True)) + 32790.6*v.SiteFraction("BCC_B2", 0, "NI")*v.SiteFraction("BCC_B2", 1, "VA")*v.SiteFraction("BCC_B2", 2, "VA") + v.SiteFraction("BCC_B2", 0, "VA")*v.SiteFraction("BCC_B2", 1, "AL")*v.SiteFraction("BCC_B2", 2, "VA")*(10.0*T + 16432.5) + 32790.6*v.SiteFraction("BCC_B2", 0, "VA")*v.SiteFraction("BCC_B2", 1, "NI")*v.SiteFraction("BCC_B2", 2, "VA"))/(0.5*v.SiteFraction("BCC_B2", 0, "AL") + 0.5*v.SiteFraction("BCC_B2", 0, "NI") + 0.5*v.SiteFraction("BCC_B2", 1, "AL") + 0.5*v.SiteFraction("BCC_B2", 1, "NI"))

    dbf = Database("""$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
    $ Date: 2019-12-08 18:05
    $ Components: AL, NI, VA
    $ Phases: BCC_B2
    $ Generated by brandon (pycalphad 0.8.1.post1)
    $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$

    ELEMENT AL FCC_A1 26.982 4577.3 28.322 !
    ELEMENT NI FCC_A1 58.69 4787.0 29.796 !
    ELEMENT VA VACUUM 0.0 0.0 0.0 !

    TYPE_DEFINITION % SEQ * !
    DEFINE_SYSTEM_DEFAULT ELEMENT 2 !
    DEFAULT_COMMAND DEFINE_SYSTEM_ELEMENT VA !

    PHASE BCC_B2 %  3 0.5 0.5 1 !
    CONSTITUENT BCC_B2 :AL,NI,VA:AL,NI,VA:VA: !

    """)
    mod = Model(dbf, ['AL', 'NI', 'VA'], 'BCC_B2')
    dd = {ky: 0.0 for ky in mod.models.keys()}
    dd['GM'] = NEW_GM
    mod.models = dd
    print(mod.HM)
    config_tup = (('AL',), ('NI', 'VA'), ('VA',))
    calculate_dict = get_prop_samples(data, config_tup)
    sample_condition_dicts = _get_sample_condition_dicts(calculate_dict, list(map(len, config_tup)))
    qty = get_data_quantities('HM_FORM', mod, [0], data, sample_condition_dicts)
    print(qty)
    assert np.all(np.isclose([-6254.7802775, -5126.1206475, -7458.3974225, -6358.04118875], qty))
Exemple #2
0
def _compare_data_to_parameters(dbf,
                                comps,
                                phase_name,
                                desired_data,
                                mod,
                                configuration,
                                x,
                                y,
                                ax=None):
    """
    Return one set of plotted Axes with data compared to calculated parameters

    Parameters
    ----------
    dbf : Database
        pycalphad thermodynamic database containing the relevant parameters.
    comps : list
        Names of components to consider in the calculation.
    phase_name : str
        Name of the considered phase phase
    desired_data :
    mod : Model
        A pycalphad Model. The Model may or may not have the reference state zeroed out for formation properties.
    configuration :
    x : str
        Model property to plot on the x-axis e.g. 'T', 'HM_MIX', 'SM_FORM'
    y : str
        Model property to plot on the y-axis e.g. 'T', 'HM_MIX', 'SM_FORM'
    ax : matplotlib.Axes
        Default axes used if not specified.

    Returns
    -------
    matplotlib.Axes

    """
    all_samples = np.array(get_samples(desired_data), dtype=np.object)
    endpoints = endmembers_from_interaction(configuration)
    interacting_subls = [
        c for c in list_to_tuple(configuration) if isinstance(c, tuple)
    ]
    disordered_config = False
    if (len(set(interacting_subls)) == 1) and (len(interacting_subls[0]) == 2):
        # This configuration describes all sublattices with the same two elements interacting
        # In general this is a high-dimensional space; just plot the diagonal to see the disordered mixing
        endpoints = [endpoints[0], endpoints[-1]]
        disordered_config = True
    if not ax:
        fig = plt.figure(figsize=plt.figaspect(1))
        ax = fig.gca()
    bar_chart = False
    bar_labels = []
    bar_data = []
    if y.endswith('_FORM'):
        # We were passed a Model object with zeroed out reference states
        yattr = y[:-5]
    else:
        yattr = y
    if len(endpoints) == 1:
        # This is an endmember so we can just compute T-dependent stuff
        temperatures = np.array([i[0] for i in all_samples], dtype=np.float)
        if temperatures.min() != temperatures.max():
            temperatures = np.linspace(temperatures.min(),
                                       temperatures.max(),
                                       num=100)
        else:
            # We only have one temperature: let's do a bar chart instead
            bar_chart = True
            temperatures = temperatures.min()
        endmember = _translate_endmember_to_array(
            endpoints[0], mod.ast.atoms(v.SiteFraction))[None, None]
        predicted_quantities = calculate(dbf,
                                         comps, [phase_name],
                                         output=yattr,
                                         T=temperatures,
                                         P=101325,
                                         points=endmember,
                                         model=mod,
                                         mode='numpy')
        if y == 'HM' and x == 'T':
            # Shift enthalpy data so that value at minimum T is zero
            predicted_quantities[yattr] -= predicted_quantities[yattr].sel(
                T=temperatures[0]).values.flatten()
        response_data = predicted_quantities[yattr].values.flatten()
        if not bar_chart:
            extra_kwargs = {}
            if len(response_data) < 10:
                extra_kwargs['markersize'] = 20
                extra_kwargs['marker'] = '.'
                extra_kwargs['linestyle'] = 'none'
                extra_kwargs['clip_on'] = False
            ax.plot(temperatures,
                    response_data,
                    label='This work',
                    color='k',
                    **extra_kwargs)
            ax.set_xlabel(plot_mapping.get(x, x))
            ax.set_ylabel(plot_mapping.get(y, y))
        else:
            bar_labels.append('This work')
            bar_data.append(response_data[0])
    elif len(endpoints) == 2:
        # Binary interaction parameter
        first_endpoint = _translate_endmember_to_array(
            endpoints[0], mod.ast.atoms(v.SiteFraction))
        second_endpoint = _translate_endmember_to_array(
            endpoints[1], mod.ast.atoms(v.SiteFraction))
        point_matrix = np.linspace(0, 1, num=100)[None].T * second_endpoint + \
            (1 - np.linspace(0, 1, num=100))[None].T * first_endpoint
        # TODO: Real temperature support
        point_matrix = point_matrix[None, None]
        predicted_quantities = calculate(dbf,
                                         comps, [phase_name],
                                         output=yattr,
                                         T=300,
                                         P=101325,
                                         points=point_matrix,
                                         model=mod,
                                         mode='numpy')
        response_data = predicted_quantities[yattr].values.flatten()
        if not bar_chart:
            extra_kwargs = {}
            if len(response_data) < 10:
                extra_kwargs['markersize'] = 20
                extra_kwargs['marker'] = '.'
                extra_kwargs['linestyle'] = 'none'
                extra_kwargs['clip_on'] = False
            ax.plot(np.linspace(0, 1, num=100),
                    response_data,
                    label='This work',
                    color='k',
                    **extra_kwargs)
            ax.set_xlim((0, 1))
            ax.set_xlabel(
                str(':'.join(endpoints[0])) + ' to ' +
                str(':'.join(endpoints[1])))
            ax.set_ylabel(plot_mapping.get(y, y))
        else:
            bar_labels.append('This work')
            bar_data.append(response_data[0])
    else:
        raise NotImplementedError(
            'No support for plotting configuration {}'.format(configuration))

    bib_reference_keys = sorted(
        list({entry['reference']
              for entry in desired_data}))
    symbol_map = bib_marker_map(bib_reference_keys)

    for data in desired_data:
        indep_var_data = None
        response_data = np.zeros_like(data['values'], dtype=np.float)
        if x == 'T' or x == 'P':
            indep_var_data = np.array(data['conditions'][x],
                                      dtype=np.float).flatten()
        elif x == 'Z':
            if disordered_config:
                # Take the second element of the first interacting sublattice as the coordinate
                # Because it's disordered all sublattices should be equivalent
                # TODO: Fix this to filter because we need to guarantee the plot points are disordered
                occ = data['solver']['sublattice_occupancies']
                subl_idx = np.nonzero(
                    [isinstance(c, (list, tuple)) for c in occ[0]])[0]
                if len(subl_idx) > 1:
                    subl_idx = int(subl_idx[0])
                else:
                    subl_idx = int(subl_idx)
                indep_var_data = [c[subl_idx][1] for c in occ]
            else:
                interactions = np.array([i[1][1] for i in get_samples([data])],
                                        dtype=np.float)
                indep_var_data = 1 - (interactions + 1) / 2
            if y.endswith('_MIX') and data['output'].endswith('_FORM'):
                # All the _FORM data we have still has the lattice stability contribution
                # Need to zero it out to shift formation data to mixing
                mod_latticeonly = Model(
                    dbf,
                    comps,
                    phase_name,
                    parameters={'GHSER' + c.upper(): 0
                                for c in comps})
                mod_latticeonly.models = {
                    key: value
                    for key, value in mod_latticeonly.models.items()
                    if key == 'ref'
                }
                temps = data['conditions'].get('T', 300)
                pressures = data['conditions'].get('P', 101325)
                points = build_sitefractions(
                    phase_name, data['solver']['sublattice_configurations'],
                    data['solver']['sublattice_occupancies'])
                for point_idx in range(len(points)):
                    missing_variables = mod_latticeonly.ast.atoms(
                        v.SiteFraction) - set(points[point_idx].keys())
                    # Set unoccupied values to zero
                    points[point_idx].update(
                        {key: 0
                         for key in missing_variables})
                    # Change entry to a sorted array of site fractions
                    points[point_idx] = list(
                        OrderedDict(sorted(points[point_idx].items(),
                                           key=str)).values())
                points = np.array(points, dtype=np.float)
                # TODO: Real temperature support
                points = points[None, None]
                stability = calculate(dbf,
                                      comps, [phase_name],
                                      output=data['output'][:-5],
                                      T=temps,
                                      P=pressures,
                                      points=points,
                                      model=mod_latticeonly,
                                      mode='numpy')
                response_data -= stability[data['output'][:-5]].values

        response_data += np.array(data['values'], dtype=np.float)
        response_data = response_data.flatten()
        if not bar_chart:
            extra_kwargs = {}
            if len(response_data) < 10:
                extra_kwargs['markersize'] = 8
                extra_kwargs['linestyle'] = 'none'
                extra_kwargs['clip_on'] = False
            ref = data.get('reference', '')
            mark = symbol_map[ref]['markers']
            ax.plot(indep_var_data,
                    response_data,
                    label=symbol_map[ref]['formatted'],
                    marker=mark['marker'],
                    fillstyle=mark['fillstyle'],
                    **extra_kwargs)
        else:
            bar_labels.append(data.get('reference', None))
            bar_data.append(response_data[0])
    if bar_chart:
        ax.barh(0.02 * np.arange(len(bar_data)),
                bar_data,
                color='k',
                height=0.01)
        endmember_title = ' to '.join([':'.join(i) for i in endpoints])
        ax.get_figure().suptitle('{} (T = {} K)'.format(
            endmember_title, temperatures),
                                 fontsize=20)
        ax.set_yticks(0.02 * np.arange(len(bar_data)))
        ax.set_yticklabels(bar_labels, fontsize=20)
        # This bar chart is rotated 90 degrees, so "y" is now x
        ax.set_xlabel(plot_mapping.get(y, y))
    else:
        ax.set_frame_on(False)
        leg = ax.legend(loc='best')
        leg.get_frame().set_edgecolor('black')
    return ax
Exemple #3
0
def plot_interaction(dbf, comps, phase_name, configuration, output, datasets=None, symmetry=None, ax=None, plot_kwargs=None, dataplot_kwargs=None) -> plt.Axes:
    """
    Return one set of plotted Axes with data compared to calculated parameters

    Parameters
    ----------
    dbf : Database
        pycalphad thermodynamic database containing the relevant parameters.
    comps : Sequence[str]
        Names of components to consider in the calculation.
    phase_name : str
        Name of the considered phase phase
    configuration : Tuple[Tuple[str]]
        ESPEI-style configuration
    output : str
        Model property to plot on the y-axis e.g. ``'HM_MIX'``, or ``'SM_MIX'``.
        Must be a ``'_MIX'`` property.
    datasets : tinydb.TinyDB
    symmetry : list
        List of lists containing indices of symmetric sublattices e.g. [[0, 1], [2, 3]]
    ax : plt.Axes
        Default axes used if not specified.
    plot_kwargs : Optional[Dict[str, Any]]
        Keyword arguments to ``ax.plot`` for the predicted data.
    dataplot_kwargs : Optional[Dict[str, Any]]
        Keyword arguments to ``ax.plot`` the observed data.

    Returns
    -------
    plt.Axes

    """
    if not output.endswith('_MIX'):
        raise ValueError("`plot_interaction` only supports HM_MIX, SM_MIX, or CPM_MIX outputs.")
    if not plot_kwargs:
        plot_kwargs = {}
    if not dataplot_kwargs:
        dataplot_kwargs = {}

    if not ax:
        ax = plt.subplot()

    # Plot predicted values from the database
    grid, predicted_values = _get_interaction_predicted_values(dbf, comps, phase_name, configuration, output)
    plot_kwargs.setdefault('label', 'This work')
    plot_kwargs.setdefault('color', 'k')
    ax.plot(grid, predicted_values, **plot_kwargs)

    # Plot the observed values from the datasets
    # TODO: model exclusions handling
    # TODO: better reference state handling
    mod_srf = Model(dbf, comps, phase_name, parameters={'GHSER'+c.upper(): 0 for c in comps})
    mod_srf.models = {'ref': mod_srf.models['ref']}

    # _MIX assumption
    prop = output.split('_MIX')[0]
    desired_props = (f"{prop}_MIX", f"{prop}_FORM")
    if datasets is not None:
        solver_qry = (tinydb.where('solver').test(symmetry_filter, configuration, recursive_tuplify(symmetry) if symmetry else symmetry))
        desired_data = get_prop_data(comps, phase_name, desired_props, datasets, additional_query=solver_qry)
        desired_data = filter_configurations(desired_data, configuration, symmetry)
        desired_data = filter_temperatures(desired_data)
    else:
        desired_data = []

    species = unpack_components(dbf, comps)
    # phase constituents are Species objects, so we need to be doing intersections with those
    phase_constituents = dbf.phases[phase_name].constituents
    # phase constituents must be filtered to only active
    constituents = [[sp.name for sp in sorted(subl_constituents.intersection(species))] for subl_constituents in phase_constituents]
    subl_dof = list(map(len, constituents))
    calculate_dict = get_prop_samples(desired_data, constituents)
    sample_condition_dicts = _get_sample_condition_dicts(calculate_dict, subl_dof)
    interacting_subls = [c for c in recursive_tuplify(configuration) if isinstance(c, tuple)]
    if (len(set(interacting_subls)) == 1) and (len(interacting_subls[0]) == 2):
        # This configuration describes all sublattices with the same two elements interacting
        # In general this is a high-dimensional space; just plot the diagonal to see the disordered mixing
        endpoints = endmembers_from_interaction(configuration)
        endpoints = [endpoints[0], endpoints[-1]]
        disordered_config = True
    else:
        disordered_config = False
    bib_reference_keys = sorted({entry.get('reference', '') for entry in desired_data})
    symbol_map = bib_marker_map(bib_reference_keys)
    for data in desired_data:
        indep_var_data = None
        response_data = np.zeros_like(data['values'], dtype=np.float_)
        if disordered_config:
            # Take the second element of the first interacting sublattice as the coordinate
            # Because it's disordered all sublattices should be equivalent
            # TODO: Fix this to filter because we need to guarantee the plot points are disordered
            occ = data['solver']['sublattice_occupancies']
            subl_idx = np.nonzero([isinstance(c, (list, tuple)) for c in occ[0]])[0]
            if len(subl_idx) > 1:
                subl_idx = int(subl_idx[0])
            else:
                subl_idx = int(subl_idx)
            indep_var_data = [c[subl_idx][1] for c in occ]
        else:
            interactions = np.array([cond_dict[Symbol('YS')] for cond_dict in sample_condition_dicts])
            indep_var_data = 1 - (interactions+1)/2
        if data['output'].endswith('_FORM'):
            # All the _FORM data we have still has the lattice stability contribution
            # Need to zero it out to shift formation data to mixing
            temps = data['conditions'].get('T', 298.15)
            pressures = data['conditions'].get('P', 101325)
            points = build_sitefractions(phase_name, data['solver']['sublattice_configurations'],
                                            data['solver']['sublattice_occupancies'])
            for point_idx in range(len(points)):
                missing_variables = mod_srf.ast.atoms(v.SiteFraction) - set(points[point_idx].keys())
                # Set unoccupied values to zero
                points[point_idx].update({key: 0 for key in missing_variables})
                # Change entry to a sorted array of site fractions
                points[point_idx] = list(OrderedDict(sorted(points[point_idx].items(), key=str)).values())
            points = np.array(points, dtype=np.float_)
            # TODO: Real temperature support
            points = points[None, None]
            stability = calculate(dbf, comps, [phase_name], output=data['output'][:-5],
                                    T=temps, P=pressures, points=points,
                                    model=mod_srf)
            response_data -= stability[data['output'][:-5]].values.squeeze()

        response_data += np.array(data['values'], dtype=np.float_)
        response_data = response_data.flatten()
        ref = data.get('reference', '')
        dataplot_kwargs.setdefault('markersize', 8)
        dataplot_kwargs.setdefault('linestyle', 'none')
        dataplot_kwargs.setdefault('clip_on', False)
        # Cannot use setdefault because it won't overwrite previous iterations
        dataplot_kwargs['label'] = symbol_map[ref]['formatted']
        dataplot_kwargs['marker'] = symbol_map[ref]['markers']['marker']
        dataplot_kwargs['fillstyle'] = symbol_map[ref]['markers']['fillstyle']
        ax.plot(indep_var_data, response_data, **dataplot_kwargs)
    ax.set_xlim((0, 1))
    ax.set_xlabel(str(':'.join(endpoints[0])) + ' to ' + str(':'.join(endpoints[1])))
    ax.set_ylabel(plot_mapping.get(output, output))
    leg = ax.legend(loc='center left', bbox_to_anchor=(1, 0.5))  # legend outside
    leg.get_frame().set_edgecolor('black')
    return ax