Executing Simulations

Overview

Questions

  • How can I execute a series of workflow steps on many simulations?

Objectives

  • Introduce workflows.

  • Demonstrate how to use signac-flow to define workflow steps and their associated pre-conditions and post-conditions.

  • Execute the workflow to randomize and compress all state points in the data space.

Boilerplate code

[1]:
import math

import hoomd
import signac

Workflow steps

The Introducing HOOMD-blue tutorial employs distinct workflow steps that must be performed in sequence to complete the self-assembly study at a given state point. These are initialization, randomization, compression, equilibration, and analysis. The previous section in this tutorial initialized every state point in the data space. This section will randomize and compress them and the next section will equilibrate them. Analysis can also be implemented as a workflow step, but this is left as an exercise for the reader.

Use signac-flow to define these workflow steps as Python functions and execute them in the proper sequence on all state points in the data space.

[2]:
import flow

Define a function that creates a Simulation object based on the signac job:

[3]:
def create_simulation(job):
    cpu = hoomd.device.CPU()

    # Set the simulation seed from the statepoint.
    sim = hoomd.Simulation(device=cpu, seed=job.statepoint.seed)
    mc = hoomd.hpmc.integrate.ConvexPolyhedron()
    mc.shape['octahedron'] = dict(vertices=[
        (-0.5, 0, 0),
        (0.5, 0, 0),
        (0, -0.5, 0),
        (0, 0.5, 0),
        (0, 0, -0.5),
        (0, 0, 0.5),
    ])
    sim.operations.integrator = mc

    return sim

This is the code from Introducing HOOMD-blue, adapted to use the seed for the current signac job with job.statepoint.seed.

Subclass FlowProject by convention.

[4]:
class Project(flow.FlowProject):
    pass

Define a function that executes the randomization step:

[5]:
@Project.operation  # Workflow step.
@Project.pre.true('initialized')  # Pre-condition in job document.
@Project.post.true('randomized')  # Post-condition in job document.
def randomize(job):
    sim = create_simulation(job)

    # Read `lattice.gsd` from the signac job's directory.
    sim.create_state_from_gsd(filename=job.fn('lattice.gsd'))
    sim.run(10e3)

    # Write `random.gsd` to the signac job's directory.
    hoomd.write.GSD.write(state=sim.state,
                          mode='xb',
                          filename=job.fn('random.gsd'))

    # Set the 'randomized' to satisfy the post-condition.
    job.document['randomized'] = True

In signac-flow’s terminology, a project operation is a complete workflow step that modifies the signac job’s state. Recall that the similarly named HOOMD-blue Operation is a class that acts on the state of the simulation at defined time steps. This tutorial uses the term workflow step wherever possible to avoid ambiguity.

@Project.operation is a decorator that declares the randomize function is a workflow step. @Project.pre and @Project.post define pre-conditions and post-conditions for this step, which return a boolean to indicate whether a workflow step is ready to start (pre-condition) or complete (post-condition). In this block, both conditions are evaluated by pre.true and post.true which examine the job document and check whether the item with the given key evaluates to True. Use pre- and post-conditions to define the sequence in which the workflow steps will execute. Here, the pre-condition checking 'initialized' is satisfied for those signac jobs that were initialized in the previous section which set job.document['initialize'] = True.

The body of the function creates the Simulation object using the create_simulation method from above, completes the randomization as in the Introducing HOOMD-blue tutorial, and sets the 'randomized' item in the job document to True. The function writes randomized.gsd to the signac job’s assigned directory using the job.fn method.

Similarly define a function that executes the compression step:

[6]:
@Project.operation
@Project.pre.after(randomize)  # Execute after randomize completes.
@Project.post.true('compressed_step')
def compress(job):
    sim = create_simulation(job)

    # Read `random.gsd` from the signac job directory.
    sim.create_state_from_gsd(filename=job.fn('random.gsd'))

    a = math.sqrt(2) / 2
    V_particle = 1 / 3 * math.sqrt(2) * a**3

    initial_box = sim.state.box
    final_box = hoomd.Box.from_box(initial_box)

    # Set the final box volume to the volume fraction for this signac job.
    final_box.volume = (sim.state.N_particles * V_particle
                        / job.statepoint.volume_fraction)
    compress = hoomd.hpmc.update.QuickCompress(
        trigger=hoomd.trigger.Periodic(10), target_box=final_box)
    sim.operations.updaters.append(compress)

    periodic = hoomd.trigger.Periodic(10)
    tune = hoomd.hpmc.tune.MoveSize.scale_solver(moves=['a', 'd'],
                                                 target=0.2,
                                                 trigger=periodic,
                                                 max_translation_move=0.2,
                                                 max_rotation_move=0.2)
    sim.operations.tuners.append(tune)

    while not compress.complete and sim.timestep < 1e6:
        sim.run(1000)

    if not compress.complete:
        raise RuntimeError("Compression failed to complete")

    # Write `compressed.gsd` to the job document.
    hoomd.write.GSD.write(state=sim.state,
                          mode='xb',
                          filename=job.fn('compressed.gsd'))

    # Set 'compressed step' in the signac job document.
    job.document['compressed_step'] = sim.timestep

This workflow step executes after the randomize step completes using the pre.after(randomize) pre-condition. The body of the function contains the code from the Introducing HOOMD-blue tutorial, changed to use the volume fraction for the current signac job with job.statepoint.volume_fraction and to read and write files from the signac job’s directory with job.fn.

The compress operation sets the compressed_step item in the job document and uses that to evaluate the post-condition. The next section of the tutorial will use the value of compressed_step.

Run the workflow

Now that you have defined the workflow steps, check the status of the workflow:

[7]:
project = Project()
project.print_status(overview=False,
                     detailed=True,
                     parameters=['volume_fraction'])

Detailed View:

job id                            operation        volume_fraction  labels
--------------------------------  -------------  -----------------  --------
59363805e6f46a715bc154b38dffc4e4  randomize [U]                0.6
972b10bd6b308f65f0bc3a06db58cf9d  randomize [U]                0.4
c1a59a95a0e8b4526b28cf12aa0a689e  randomize [U]                0.5

[U]:unknown [R]:registered [I]:inactive [S]:submitted [H]:held [Q]:queued [A]:active [E]:error


Each signac job is ready to execute randomize, the first step in the workflow. Run it:

[8]:
project.run(names=['randomize'])

Every signac job directory in the data space now has a random.gsd file produced by randomize:

[9]:
!ls workspace/*
workspace/59363805e6f46a715bc154b38dffc4e4:
lattice.gsd              signac_job_document.json
random.gsd               signac_statepoint.json

workspace/972b10bd6b308f65f0bc3a06db58cf9d:
lattice.gsd              signac_job_document.json
random.gsd               signac_statepoint.json

workspace/c1a59a95a0e8b4526b28cf12aa0a689e:
lattice.gsd              signac_job_document.json
random.gsd               signac_statepoint.json

Now, the status shows that the compress step is ready:

[10]:
project.print_status(overview=False,
                     detailed=True,
                     parameters=['volume_fraction'])

Detailed View:

job id                            operation       volume_fraction  labels
--------------------------------  ------------  -----------------  --------
59363805e6f46a715bc154b38dffc4e4  compress [U]                0.6
972b10bd6b308f65f0bc3a06db58cf9d  compress [U]                0.4
c1a59a95a0e8b4526b28cf12aa0a689e  compress [U]                0.5

[U]:unknown [R]:registered [I]:inactive [S]:submitted [H]:held [Q]:queued [A]:active [E]:error


Execute it:

[11]:
project.run(names=['compress'])

Every signac job directory in the data space now has a compressed.gsd file produced by compress:

[12]:
!ls workspace/*
workspace/59363805e6f46a715bc154b38dffc4e4:
compressed.gsd           random.gsd               signac_statepoint.json
lattice.gsd              signac_job_document.json

workspace/972b10bd6b308f65f0bc3a06db58cf9d:
compressed.gsd           random.gsd               signac_statepoint.json
lattice.gsd              signac_job_document.json

workspace/c1a59a95a0e8b4526b28cf12aa0a689e:
compressed.gsd           random.gsd               signac_statepoint.json
lattice.gsd              signac_job_document.json

Summary

In this section of the tutorial, you defined the workflow steps to randomize and compress the initial configuration using signac-flow, along with the pre- and post-conditions needed to sequence the steps. Then you executed the workflow steps on all state points in the dataset. The directory for each simulation now contains compressed.gsd and is ready for equilibration at the target volume fraction.

The next section in this tutorial teaches you how to write a workflow step that can continue itself and complete over several submissions.

This tutorial only teaches the basics of signac-flow. Read the signac-flow documentation to learn more.