38. truchas Python Module

Truchas provides a Python package for scripting. This is useful for generating custom restarts, performing automated optimization, and performing automated parameter studies, or computing custom metrics.

38.1. Examples

Custom Restart

Frequently it is necessary to generate custom mapped restarts. For example, if we need to map data to a new mesh, where there is a new block which needs to be initialized. Generally, we will need to initialize this block with a material and a temperature. The below script shows how to read in a Truchas h5 file, map it to a new mesh, and modify some field data before generating a restart file.

import numpy as np

import truchas

# Load truchas data from the preheat run, map it to a new mesh, and scale that mesh.
tdata = truchas.TruchasMappedData("preheat_output/preheat.h5", "mesh2.gen", 0.001)

# Select the last cycle in the preheat output. Further calls will read data from
# this cycle, overwrite data at this cycle, and generate a restart with our
# modifications.
sid = tdata.num_series()

# Set the VOF to [0, 1] (purely the final material in the input) in block 2
tdata.assign_value_block(sid, "Z_TEMP", 1, f)
tdata.assign_value_block(sid, "Z_TEMP", 3, f)

# Compute the function f = 1 + 2*x at cell centers across the mesh.
xc = tdata.centroids()
f = 1 + 2 * xc[0]

# Reassign the temperature to the function f, but only in blocks 1 and 3.
tdata.assign_value_block(sid, "Z_TEMP", 1, xc[:,0])

# Alternatively, we could read the whole field and modify it before reassigning.
temp = tdata.field(sid, "Z_TEMP")
graphite_region = tdata.region(1, 3)
temp[graphite_region] = f[graphite_region]
tdata.assign_field(sid, "Z_TEMP", temp)

# Alternatively, we could use Numpy functions
temp = tdata.field(sid, "Z_TEMP")
temp = np.where(tdata.blockid() == 1 or tdata.blockid() == 3, f, temp)
tdata.assign_field(sid, "Z_TEMP", temp)

# Write the custom restart file
tdata.write_restart("preheat.restart", sid)

Example Truchas Python script. This script reads a Truchas output file, generates a mapped restart onto a new mesh, and overwrites data for some variables on a specific block on that new mesh.

Optimization

Python scripting can be used to execute Truchas in a loop, for example to modify some heat flux boundary condition until the temperature at some point in the mesh fits an expected value. We can wrap this entire process into a Python function, then use an optimization function from the popular Scipy package to solve \(T(F) = 500\) for the temperature \(T\) at \(x = (0, 0, 0)\).

...

&PROBE
  coord = 0, 0, 0
  data = 'temperature'
  data_file = 'top-temperature.txt'
/

...

&THERMAL_BC
  name = 'heat-coil'
  face_set_ids = 14
  type = 'flux'
  flux = -{heat_flux:.4e} ! [W / m^2]
/

...

Part of a Truchas template input, here named template-heatup.inp. The text {heat_flux:.4e} will be replaced by the Python-generated input deck below.

import scipy.optimize as spopt

import truchas

def final_top_temperature(heat_flux, tenv, nprocs):
    # Generate input deck
    parameters = {"heat_flux": heat_flux}
    tenv.generate_input_deck(parameters, "template-heatup.inp", "heatup.inp")

    # Run truchas with the new input
    stdout, output = tenv.truchas(nprocs, "block-heatup.inp")

    # Read the temperature probe
    temperature = output.probe_data("top-temperature.txt")

    # It is useful to print both the input and output for every iteration,
    # since each iteration is a full simulation and can be slow.
    print(f"heat_flux = {heat_flux:.4e}, temperature = {temperature[-1,1]:.4e}\n")

    # Return the last element (final time) of the second column (temperature column)
    return temperature[-1,1]


if __name__ == "__main__":
    nprocs = 8

    # We need to allow overwriting the output directory, since each iteration of
    # the optimizer will re-run Truchas.
    tenv = truchas.TruchasEnvironment.default(overwrite_output=True)

    # We'll wrap our function final_top_temperature in a lambda, so we can pass
    # in the TruchasEnvironment and nprocs variables. We'll also offset by our
    # target value, since we're using a root finder.
    fx = lambda x: final_top_temperature(x, tenv, nprocs) - 500

    # The Scipy.Optimize.brentq root finder is used. Here we're passing our
    # lambda, an interval for our heat flux of [2.5e4, 5e4], a tolerance of 50
    # K, and absolute and relative tolerances on the heat flux.
    heat_flux = spopt.brentq(fx, 2.5e4, 5e4, xtol=50, rtol=5e-3)
    print(heat_flux)

Parameter Study

Suppose we wish to perform the following parameter study: Vary the heat flux on one boundary condition, and measure the maximum and average temperature in block 1. The following script will define the custom output metrics (max and avg temperature on block 1), the parameter dictionary, and perform the parameter study. This will produce a text file vary_heat_flux.txt which contains a table with the input and output values.

This script will also generate a database which caches Truchas outputs. This way, if we re-run the parameter study with more heat flux values later, we don’t need to re-run entire simulations which have already been computed.

import numpy as np

import truchas


# read the maximum temperature in block 1
def max_temperature_str(input_parameters, output):
    sid = output.num_series()
    temperature = output.field(sid, "Z_TEMP")
    block1 = output.region(1)
    max_temperature = max(temperature[block1])
    print(f"Max temperature = {max_temperature:.2f} K.")
    return f"{max_temperature:.2f}"


# compute the average temperature in block 1
def avg_temperature_str(input_parameters, output):
    sid = output.num_series()
    temperature = output.field(sid, "Z_TEMP")
    cell_volume = output.volumes()
    block1 = output.region(1)
    avg_temp = np.average(temperature[block1], weights=cell_volume[block1])
    print(f"Avg temperature = {avg_temp:.2f} K.")
    return f"{avg_tmp:.2f}"


if __name__ == "__main__":
    # Store the results from each run in a database. If we wish to run this
    # parameter study again with additional points, each run will be cached so
    # we don't need to re-run data we've already calculated.
    nproc = 8
    tenv = truchas.TruchasEnvironment.default(overwrite_output=True)
    tdb = truchas.TruchasDatabase("parameter_study")
    tstudy = truchas.TruchasStudy(tenv, tdb, nproc)

    # Suppose the template file has multiple values that need to be inserted,
    # but we're only doing a parameter study on "heat_flux2".
    input_parameters = {"heat_flux1": 4e4,
                        "heat_flux2": 1e4,
                        }
    heat_flux2_points = [1e4, 2e4, 4e4]

    # The output will be a text file with the following columns: "heat_flux2",
    # "max_temperature", and "average_temperature". The latter two are custom
    # output metrics, tied to functions defined above.
    output_metrics = {"max_temperature": max_temperature_str,
                      "average_temperature": avg_temperature_str,
                      }

    tstudy.do_1d_parameter_study("template-heatup.inp",
                                 input_parameters, "heat_flux2", heat_flux2_points,
                                 output_metrics, "vary_heat_flux.txt")

38.2. Components

TruchasEnvironment

class TruchasEnvironment(mpiexec, truchas_executable, input_dir, write_restart_executable, working_dir='.', overwrite_output=False)

This is the main interface for running Truchas via Python, grabbing outputs, and writing restarts. It is usually recommended to initialize your environment using the default initializer:

tenv = truchas.TruchasEnvironment.default()

See the default method for more parameters.

classmethod default(input_dir=None, overwrite_output=None)

Initialize using the CMake-generated default configuration. The working directory is initialized to the input_dir.

Parameters:
  • input_dir (str, optional) – The directory which contains the input file and where the output directory should be created. May be an absolute path, or relative to the current working directory. Defaults to the current working directory.

  • overwrite_output (bool, optional) – Choose whether or not to overwrite the output by default on subsequent calls to .truchas(…). Defaults to True if the TRUCHAS_OVERWRITE_OUTPUT environment variable is defined to a nonzero value, and False otherwise.

classmethod from_argv()

Initialize from command line arguments.

usage: test.py [-h] -c CONFIG [-o OUTPUT]

-h, --help

show this help message and exit

-c CONFIG, --config CONFIG

CMake-generated configuration file

-o OUTPUT, --output OUTPUT

Output directory

classmethod from_config_file(config_file, input_dir, working_dir='.')

Initialize from a json configuration file, input directory, and an optional working directory.

Parameters:
  • config_file (str) – Filename of the input JSON file containing the configuration parameters described below. May be relative to the current working directory, or an absolute path.

  • input_dir (str) – The directory which contains the input file. May be an absolute path, or relative to the current working directory.

  • working_dir (str, optional) – The directory where output directories should be generated. May be an absolute path, or relative to the current working directory. Defaults to the current working directory.

JSON configuration file parameters:

Parameters:
  • "mpiexec" (str) – Path to the mpiexec executable.

  • "truchas_executable" (str) – Path to the truchas executable.

  • "write_restart_executable" (str) – Path to the write-restart.py executable.

  • "overwrite_output" (bool) – Chooses whether to overwrite the output by default or not.

generate_input_deck(replacements, template_filenames, output_filename)

Generates a truchas input file from a template (or templates) and replacements. Expects a dictionary of replacements, e.g. {"sigma": 1.2e-3}. Will replace Python format-strings in the input template file, e.g. {sigma} is then replaced with .0012. Formatting can be specified in the input template in the typical Python way, such as {sigma:.2e}.

Parameters:
  • replacements (dict) – A dictionary of replacements to be inserted into the Truchas input deck according to the given template. Dictionary keys are the strings to replace, e.g. "sigma" or "temp", and the values are the items to insert.

  • template_filename (str or list of str) – Filename of the template for the Truchas input deck. Can also be a list of template filenames that should be concatenated. Should contain Python format-strings, such as {sigma} or {sigma:.2e} where replacements from the dictionary (e.g. {"sigma": 1.2e-3}) will be inserted.

  • output_filename (str) – Filename of the generated Truchas input deck.

input_directory()

Return the input directory

open_data(h5file)

Returns an output data object from an h5file, expected to be relative to the working directory.

Parameters:

h5file (str) – A Truchas output .h5 file. Expected to be relative to the working directory.

Return type:

TruchasData

output(output_file)

Returns an output data object from an h5file, expected to be relative to the input directory.

Parameters:

h5file (str) – A Truchas output .h5 file. Expected to be relative to the input directory.

Return type:

TruchasData

truchas(nprocs, input_file, restart_file=None, output_dir=None, overwrite_output=None)

Runs truchas with specified number of MPI ranks, and returns the terminal output and the output data.

stdout, output = tenv.truchas(4, "cast.inp")

Parameters:
  • nprocs (int) – Number of MPI ranks.

  • input_file (str) – Filename of a Truchas input deck.

  • restart_file (str, optional) – Filename of a Truchas restart file.

  • output_dir (str, optional) – User-specified Truchas output directory.

  • overwrite_output (bool, optional) – Overwrite the output directory if it already exists. Defaults to the value provided to the TruchasEnvironment initializer.

Returns:

stdout and Truchas output data.

Return type:

str, TruchasData

working_directory()

Return the working directory

write_restart(h5file, cycle_number, restart_file)

Write a restart file from the given H5 dump and cycle number to the given file. Files expected to be relative to working directory.

WARNING: This calls the write-restart.py program, and is primarily used for testing. Usually in Python scripts, it is better to use the TruchasData.write_restart() method.

Parameters:
  • h5file (str) – A Truchas output .h5 file. Expected to be relative to the working directory.

  • cycle_number (int) – The cycle used to generate the restart file.

  • restart_file (str) – The name of the restart file to be generated.

TruchasData

class TruchasData(filename)

A class for reading and interacting with Truchas output data. It provides access to fields, and allows users to manipulate them.

Field manipulation is in-memory, not saved to the H5 file. The Truchas output is read-only. Modifying data is intended for modified restarts, and is done through the assign_field and assign_value functions. The field routine will return a user-modified field, if one has been assigned via one of those two methods.

TruchasData objects may be directly initialized by its initializer. TruchasData objects are also provided by TruchasEnvironment.truchas() and TruchasEnvironment.output().

Parameters:

filename (str) – Filename of the Truchas h5 output file.

assign_field(series_id, field_name, field)

Reassign a field with a user-defined copy. Future calls to TruchasData.field() will return this copy, and this copy will be used in any restarts.

Note: VOF is a 2D array. The order of the elements is the same as listed in the PHYSICS namelist materials parameter and MATERIAL namelists phases parameter.

Parameters:
  • series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

  • field_name (str) – The name of a field in the Truchas output file. Valid values include "Z_TEMP", "Z_RHO" (density), "Z_ENTHALPY", "phi1" (species 1 concentration), "sigma" (stress tensor), "epsilon" (strain tensor), "epstherm" (thermal strain), "Rotation", "Displacement", "Gap Displacement", "Gap Normal Traction", "VOF", "Z_P" (pressure), and "Z_VC" (velocity). Valid values for the active h5 file can be found with the TruchasData.field_names() function.

  • field (numpy.ndarray or number) – When assigning to a scalar field (.e.g "Z_TEMP", "Z_RHO"), expected is either a 1D Numpy array with the same length as the field being overwritten (ncell or nnode), or a number. When assigning to a 2D array (e.g. "VOF" or "sigma"), expected is either a 2D Numpy array with the same shape as the field being written, or a 1D Numpy array with the same length as the second dimension of the field being overwritten (field.shape[1]).

assign_value_block(series_id, field_name, blockid, value)

Assign a value to a field within an element block.

Parameters:
  • series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

  • blockid (int) – The ID of cell block in the Truchas output file.

  • field_name (str) – The name of a field in the Truchas output file. Valid values include "Z_TEMP", "Z_RHO" (density), "Z_ENTHALPY", "phi1" (species 1 concentration), "sigma" (stress tensor), "epsilon" (strain tensor), "epstherm" (thermal strain), "Rotation", "Displacement", "Gap Displacement", "Gap Normal Traction", "VOF", "Z_P" (pressure), and "Z_VC" (velocity). Valid values for the active h5 file can be found with the TruchasData.field_names() function.

  • field (numpy.ndarray or number) – When assigning to a scalar field (.e.g "Z_TEMP", "Z_RHO"), expected is either a 1D Numpy array with the same length as the field being overwritten (ncell or nnode), or a number. When assigning to a 2D array (e.g. "VOF" or "sigma"), expected is either a 2D Numpy array with the same shape as the field being written, or a 1D Numpy array with the same length as the second dimension of the field being overwritten (field.shape[1]).

assign_value_cell(series_id, field_name, cellid, value)

Assign a value to a field at a given cell. This is equivalent to:

f = TruchasData.output(sid, "xxx")
f[cellid] = x
TruchasData.assign_field(sid, "xxx", f)
Parameters:
  • series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

  • blockid (int) – The ID of cell block in the Truchas output file.

  • field_name (str) – The name of a field in the Truchas output file. Valid values include "Z_TEMP", "Z_RHO" (density), "Z_ENTHALPY", "phi1" (species 1 concentration), "sigma" (stress tensor), "epsilon" (strain tensor), "epstherm" (thermal strain), "Rotation", "Displacement", "Gap Displacement", "Gap Normal Traction", "VOF", "Z_P" (pressure), and "Z_VC" (velocity). Valid values for the active h5 file can be found with the TruchasData.field_names() function.

  • field (numpy.ndarray or number) – When assigning to a scalar field (.e.g "Z_TEMP", "Z_RHO"), expected is a scalar number. When assigning to a 2D array (e.g. "VOF" or "sigma"), expected is a 1D Numpy array with the same length as the second dimension of the field being overwritten (field.shape[1]).

blockid()
Returns:

A ncell length field of integers, holding the Exodus block ID at for each cell.

Return type:

numpy.ndarray

centroids()
Returns:

A [ncell,3] field of cell centroids.

Return type:

numpy.ndarray

cycle(series_id)

Return the cycle number of a given series.

Parameters:

series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

Returns:

The cycle number at the given series_id

Return type:

int

field(series_id, field_name)

Return the requested field of the requested series as a ndarray. This will read the entire field from disk. If the field has been reassigned, returns the user-defined copy.

Parameters:
  • series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

  • field_name (str) – The name of a field in the Truchas output file. Valid values include "Z_TEMP", "Z_RHO" (density), "Z_ENTHALPY", "phi1" (species 1 concentration), "sigma" (stress tensor), "epsilon" (strain tensor), "epstherm" (thermal strain), "Rotation", "Displacement", "Gap Displacement", "Gap Normal Traction", "VOF", "Z_P" (pressure), and "Z_VC" (velocity). Valid values for the active h5 file can be found with the TruchasData.field_names() function.

Returns:

The requested field from the H5 file. If the field has been reassigned, returns the user-defined copy.

Return type:

numpy.ndarray

field_names(series_id=1)

Return a list of field names present in the first series.

Parameters:

series_id (int, optional) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()]. Defaults to 1. Usually the fields present does not change over the course of a simulation.

Result:

A list of field names present in the requested series.

Return type:

list of strings

filename

The filename of the HDF5 Truchas output.

node_coordinates()
Returns:

A [nnode,3] field of node coordinates.

Return type:

numpy.ndarray

num_series()
Returns:

The number of series in this file.

Return type:

int

num_species()
Returns:

The number of species in this file.

Return type:

int

probe_data(probe_filename)

Read data from a probe file.

Parameters:

probe_filename (str) – Name of a probe file. Expected to be relative to the directory holding the h5 file.

Returns:

2D array containing the text file data.

Return type:

numpy.ndarray

region(*region_blockids)

Return a list of bools indicating which cells are part of a given set of block ids. Input is an arbitrary number of block id integers.

Parameters:

region_blockids (int) –

An arbitrary number of block id integers. This function takes an arbitrary number of arguments, not a list. E.g.,

TruchasData.region(1)
TruchasData.region(1, 2, 3)

Returns:

A list of bools indicating which cells are part of a given set of block ids. E.g., [True, False, False, True, ...]. This is useful for masking off operations, for example:

mask = TruchasData.region(1, 2)
f1[mask] = f2[mask]
f1[mask] = 0

Return type:

numpy.ndarray of bools

region_node(*region_blockids)

Return a list of bools indicating which nodes are part of a given set of block ids. Input is an arbitrary number of block id integers. A node is considered part of the given block set if any one of the cells adjacent to it is part of one block in that set.

Parameters:

region_blockids (int) –

An arbitrary number of block id integers. This function takes an arbitrary number of arguments, not a list. E.g.,

TruchasData.region(1)
TruchasData.region(1, 2, 3)

Returns:

A list of bools indicating which nodes are part of a given set of block ids. E.g., [True, False, False, True, ...]. This is useful for masking off operations, for example:

mask = TruchasData.region(1, 2)
f1[mask] = f2[mask]
f1[mask] = 0

Return type:

numpy.ndarray of bools

series_id(cycle)

Return the series id of a given cycle. If not found, returns None.

Parameters:

cycle (int) – The number of a Truchas cycle.

Returns:

The series ID matching a given cycle. If there is no series exactly matching the cycle, returns None.

Return type:

int or None

stdout()
Returns:

Truchas stdout

Return type:

str

time(series_id)

Return the time at a given series number.

Parameters:

series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

Returns:

The time at the given series_id

Return type:

number

time_step(series_id)

Return the time step size at a given series number.

Parameters:

series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

Returns:

The time step size at the given series_id

Return type:

number

volumes()
Returns:

A length ncell field of cell volumes.

Return type:

numpy.ndarray

write_restart(outfile, series_id)

Write a restart file from given Truchas data object and series ID.

Parameters:
  • outfile (str) – The name of the restart file to be generated.

  • series_id (int) – The ID of a series in the Truchas output file. Valid values are in the range [1, TruchasData.num_series()].

TruchasMappedData

class TruchasMappedData(output_filename, exodus_filename, scale_factor=1, use_portage=False)

Bases: TruchasData

Interface for mapping Truchas Data to new exodus mesh. This class inherits all of the methods from the TruchasData class, and invisibly maps fields to the new mesh prior to returning them. Initialization takes extra parameters for mapping.

Parameters:
  • filename (str) – Filename of the Truchas h5 output file.

  • exodus_filename (str) – Filename of the exodus mesh file to which the data shall be mapped

  • scale_factor (number) – Scaling factor for the target exodus mesh.

  • use_portage (bool) – Flag to use the Portage backend rather than the Truchas built-in mapper. Defaults to False, using the Truchas built-in mapper.

finalize()

Free the memory used by the grid mapper. If mapping between multiple different mesh combinations, it’s best to call this directly rather than wait for Python garbage collection to clean it up.

TruchasStudy

class TruchasStudy(tenv, tdb, nprocs, njobs=1, restart_file=None)

A class for performing Truchas parameter studies which vary in one parameter (1D). This class expects a Truchas environment which will perform the tests (or generate input decks) and a Truchas Database which will cache the results.

One will generally want to execute either do_1d_parameter_study(), or generate_1d_parameter_study_inputs() and register_1d_parameter_study_outputs() and do_1d_parameter_study() in that order.

In the former case, the parameter study will be performed by sequentially running Truchas on the supplied metric. In the latter case, a collection of input files will be generated, which might be executed on a cluster via a SLURM array job. The results of these studies are then registered to the supplied database via the register_1d_parameter_study_outputs() method. At this point, running do_1d_parameter_study() will read outputs from the database of Truchas outputs instead of executing simulations.

To initialize:

Parameters:
  • tenv (TruchasEnvironment) – A TruchasEnvironment for running the parameter study.

  • tdb (TruchasDatabase) – A TruchasDatabase where results will be cached.

  • nproc (int) – Number of MPI ranks to use. When used with generate_1d_parameter_study_inputs(), this parameter is ignored.

  • njobs (int, optional (default 1)) – Number of Truchas simulations to run in simultaneously.

  • restart_file (str, optional (default None)) – Filename of Truchas restart file. This restart will be used for all runs.

do_1d_parameter_study(template_input_file, initial_parameters, variable, points, output_metrics, output_filename, extra_outputs=None)

Do a 1D parameter study of the given variable, across the given points. Output metrics are recorded to a text file.

Parameters:
  • initial_parameters (dict with str keys) – A dictionary of input parameters to supply to TruchasEnvironment.generate_input_deck(). Keys are the python format variables present in the template input deck. Values are either values to replace those format strings, or functions. If a function, it is evaluated with a single input: the dict of parameters in the current iteration. This allows values to be generated based on other values at this iteration. E.g., the “h2” key could hold a lambda which always returns double the value of the “h1” key.

  • variable (str) – The name of a variable in the initial_parameters dictionary which is to be varied for the parameter study.

  • points (list of values) – The values of the variable to iterate over for the parameter study.

  • output_metrics (dict with str keys and function values.) – A dictionary of output metrics to be printed to the output file. Functions are expected to accept two arguments: the input parameter dictionary and a TruchasData output set. These functions are to be written by the user to generate custom outputs. The result of each function must be a string.

  • output_filename (str) – A filename where the table of inputs and outputs will be written.

  • extra_outputs (function, optional) – A function of four arguments, the parameter dictionary, the Truchas output, the study filename (intended for generating new filenames), and the parameter study iteration index. This is for custom operations to be performed at each run. For instance, generating unique new files each containing the maximum temperature evolution for a different run.

Note

  • Toolpath filename should be specified in the template.inp file with the parameter {toolpath_file}. This parameter is skipped when generating the instance unique identifier, so both the Truchas input file and toolpath file are renamed to a unique identifier independent of the toolpath filename.

Warning

TODO Some current hardwired shortcomings…

  • Toolpath templates are hardcoded to be “template.json”.

  • Only a single json toolpath file is permitted, or none. The parameter key “toolpath_file” is disregarded when generating an input deck’s unique identifier for the TruchasDatabase. This is because the toolpath file ought to be named with the same hash as the input file.

generate_1d_parameter_study_inputs(template_input_file, initial_parameters, variable, points)

Generate input decks across a given span of data points. This is intended for generating inputs to then be simulated on a cluster via SLURM, not in a Python environment. These output files will be stored in the directory parameter_study_inputs. With this approach, it is expected to:

  1. Use generate_1d_parameter_study_inputs() to generate Truchas input files in the parameter_study_inputs directory.

  2. Run Truchas independent of Python, e.g. via a SLURM script.

  3. Use register_1d_parameter_study_inputs() to identify the Truchas outputs and move them into a TruchasDatabase.

  4. Use do_1d_parameter_study() to “run” the parameter study. This time, all of the outputs will be recognized as already-generated in the database, so only the postprocessing will be performed.

See do_1d_parameter_study() for argument documentation.

An example SLURM script is shown below for running each of these generated input files in a single array job.

#SBATCH -a 1-16  # Run 16 duplicates of this job in parallel

cd <path-to-parameter-study-inputs>

# exit if this run doesn't correspond to a file
nfiles=$(ls -1 *.inp | wc -l)
(( $SLURM_ARRAY_TASK_ID > $nfiles )) && exit

# select a unique input file from an alphabatized list and my job array ID
inputfile=$(ls -1 *.inp | sed -n "${SLURM_ARRAY_TASK_ID}p")
ml truchas
srun -c 1 truchas $inputfile
generate_input_deck(template_input_file, replacements, use_working_dir=False)

Generate a usable Truchas input deck from a given template file(s) and a dictionary of replacements.

generate_inputs(template_input_file, replacements, use_working_dir=False)

Generate Truchas input decks from a given template file(s) and a list of dictionaries of replacements.

Tip

If you need to delete the input decks generated by this input without deleting the templates, carefully use a command like the following:

ls -1 *inp | grep -v template | xargs rm
register_1d_parameter_study_outputs(initial_parameters, variable, points)

Register Truchas simulation outputs and input decks to a database. If generate_1d_parameter_study_inputs() is used to generate input files in the directory parameter_study_inputs, then Truchas is used externally to generate outputs for each of these input files, this function is used to register those outputs to a TruchasDatabase.

See do_1d_parameter_study() for argument documentation. See generate_1d_parameter_study_inputs() for the workflow description.

TruchasDatabase

class TruchasDatabase(directory)

Records Truchas outputs into a database for easy lookup. This class is used by TruchasStudy, and for most user code will not need to be used directly except for initialization.

Parameters:

directory (str) – A directory where the database will be stored. Should be either an empty directory, nonexistent directory, or the path of an already-existing database.

append_new_parameters(parameters)

Append a new set of parameters to the database. This requires rehashing all identifiers.

Parameters:

parameters (a dict with str keys and non-function values) – A dictionary of template replacements, including any new keys that were not present in previous database entries. For accuracy, any new keys must contain the value that was used in old runs (where presumably they were hardcoded in the input file).

backup_database()

Create a backup file database-backup.json in the database directory. This backs up the database structure, not the Truchas inputs nor output data.

delete(parameters)

Delete the data associated with the given parameter set from the database. Note this will delete all associated data from disk.

Parameters:

parameters (a dict with str keys and non-function values) – A dictionary of template replacements. Unlike TruchasStudy.do_1d_parameter_study(), this cannot contain functions (rather, it expects the functions have already been evaluated). This is so that a unique identifier can be produced.

donate(parameters, infile, dump_dir)

Add the input file and Truchas output directory associated with a set of parameters to the database. Note this will overwrite any existing item associated with the parameter set in the database.

Parameters:
  • parameters (a dict with str keys and non-function values) – A dictionary of template replacements. Unlike TruchasStudy.do_1d_parameter_study(), this cannot contain functions (rather, it expects the functions have already been evaluated). This is so that a unique identifier can be produced.

  • infile (str) – A Truchas input file

  • dump_dir (str) – A Truchas output dir associated with infile

exists(parameters)

Checks if the given parameter combination is present in the database.

Parameters:

parameters (a dict with str keys and non-function values) – A dictionary of template replacements. Unlike TruchasStudy.do_1d_parameter_study(), this cannot contain functions (rather, it expects the functions have already been evaluated). This is so that a unique identifier can be produced.

Returns:

True if the Truchas results associated with parameters is already present in the database, False otherwise.

Return type:

bool

static identifier(parameters)

Generates a unique identifier for the given parameter dictionary.

Parameters:

parameters (a dict with str keys and non-function values) – A dictionary of template replacements. Unlike TruchasStudy.do_1d_parameter_study(), this cannot contain functions (rather, it expects the functions have already been evaluated). This is so that a unique identifier can be produced.

Returns:

Unique identifier hash

Return type:

str

truchas_output(parameters)

Return the Truchas output corresponding to the given parameters. If such an output does not exist, return None.

Parameters:

parameters (a dict with str keys and non-function values) – A dictionary of template replacements. Unlike TruchasStudy.do_1d_parameter_study(), this cannot contain functions (rather, it expects the functions have already been evaluated). This is so that a unique identifier can be produced.

Returns:

Truchas output associated with parameters

Return type:

TruchasData

TruchasTest

compare_l2(field1, field2, tol, name, time)

Calculates the L2 norm between two fields, and compares against an input tolerance. Returns 0 if error is below tolerance and 1 otherwise. Prints result of comparison to terminal. Fields may be a scalar value.

compare_max(field1, field2, tol, name, time)

Calculates the maximum absolute difference between two fields, and compares against an input tolerance. Returns 0 if error is below tolerance and 1 otherwise. Prints result of comparison to terminal. Fields may be a scalar value.

compare_max_rel(field1, field2, tol, name, time)

Calculates the maximum relative absolute difference between two fields, and compares against an input tolerance. Returns 0 if error is below tolerance and 1 otherwise. Prints result of comparison to terminal. Fields may be a scalar value.