WARPXM v1.10.0
Loading...
Searching...
No Matches
Warpy: the WARPXM Python interface

Warpy is a Python interface for setting up and running WARPXM simulations.

Features:

  • Compatibility with Python 2 and 3 (although Python 3 is strongly recommended).
  • Input file generation for FVM and DG simulations using pure python code.
  • Significantly reduces specifying redundant information and annoying boilerplate.
  • Automates run directory setup, and copying of important files (input files, meshes, Git version information).
  • Handles job submission on HPC machines.
  • Can handle basic post processing (using matplotlib package directly or by generating a Paraview xdmf file.)
  • Allows generating multiple cases for performing parametric case studies/convergence tests.
  • Functionality for running WARPXM with or without MPI

Getting started

To get started using warpy, import the warpy package

import warpy

Note: The path to warpy must be in your PYTHONPATH (see install instructions).

TODO: Make warpy an installable package with pip.

Configuring warpy locally

The warpy_user_config.py file contains common user-adjustable settings for running WARPXM simulations with warpy, such as number of MPI process, where to place simulation run directories, and queue names on HPC clusters, to name a few.

Warpy will look for the user config file in $HOME/.config. If it is not already there, warpy will copy a default version to that location.

Configuration options are regularly added. Please see warpy_user_config.py for current options and description of use.

Extending warpy

See warpy developer information.

How warpy relates to the WARPXM input File

WARPXM simulations are setup based on instructions provided in an xml-like formatted input file with the extension .inp (see input file documentation for more details). These input files become long and difficult to read, even for modestly complex simulations. Warpy was developed to set up WARPXM simulations in a more intuitive, object-oriented way, as well as automate directory set up, copying of useful files, and performing some basic post processing. Warpy can also automate the submission of jobs on HPC machines.

The next section discusses the creation of a WARPXM simulation input file using the warpy library.

Constructing a simulation input file with warpy

WARPy is a regular Python package; that means you can do anything Python lets you do. The general process for writing a warpy script is:

Define a mesh

WARPXM/WARPy allows defining two different types of meshes: rectangular block generated meshes, or Abaqus meshes. Periodic boundaries are associated topologically with meshes, so these are specified at this stage rather than when you specify other boundary conditions. This is done by specifying pairs of nodesets are periodic.

sorder = 'secondOrder'
bmesh = warpy.mesh.block(Bounds=[-1,1,0,1],
NumCells=[128,128],
NodeSets=['Left','Right','Bottom','Top'],
NumLayers=1,
PeriodicBoundaries=['Left','Right','Bottom','Top'],
basis_array_set=sorder)
amesh = warpy.mesh.arbitrary('mesh_filename.inp',
NumLayers=1,
basis_array_set=sorder,
PeriodicBoundaries=['Left','Right'])
general unstructured mesh
Definition: mesh.py:56
block mesh generator
Definition: mesh.py:87

Define variables:

Variables represent a collection of scalar components expressed on some basis. They can be restricted to particular subdomains, as well as defined to exist virtually on a subdomain. A variable can also be defined to be explodable or not. Explodable variables are expanded depending on the number of stages required to implement a particular RK method; non-explodable variables only exist once. An example usage of non-explodable variables is implementing a static electric field while evolving other variables' response. By default if no subdomains are specified the variable exists globally in the domain.

sorder = 'firstOrder'
fluida = warpy.variable(name='a',
components=['rho', 'px', 'py', 'pz', 'e'],
basis_array_set=sorder,
subdomains=['regiona'])
efield = warpy.variable(name='E',
components=['Ex','Ey','Ez'],
basis_array_set=sorder,
subdomains=['regionb'],
virtual_subdomains=['regionc'],
explodable=False)
Definition: variable.py:4

Define initial conditions

TODO: This section is out-of-date. Needs updating.

Initial conditions are defined in two stages: first by creating a functional operation, then creating an applicator to apply the function to a variable. You can re-use previously defined function operations multiple times to initialize different variables. When creating applicators, you can specify a few components, or use the default which passes in all components.

# define some functions
null_func = warpy.functions.set_to(name='null', value=0)
idens_ic = warpy.functions.heaviside(name='idens_ic',
direction=[1,0,0],
center=[0,0,0],
neg=Ai*n_0*(1-jump),
pos=Ai*n_0*(1+jump))
# apply functions to variables
applicators = [
warpy.applicator(spatial_order=sorder,
fun=null_func,
var=field,
components=['Ex','Ey','Ez'],
spatial_scheme='Nodal'),
warpy.applicator(spatial_order=sorder,
fun=null_func,
var=field2,
spatial_scheme='Nodal'),
warpy.applicator(spatial_order=sorder,
fun=idens_ic,
var=ions,
components=['rho'])]
Heaviside step function.
Definition: heaviside.py:4
Function which sets field to a given uniform value.
Definition: set_to.py:4

Define equation set applications and associated spatial solvers

Applications in warpy.apps are used to implement the physics of a particular model. See the individual apps for what parameters to pass. You can define a spatial solver to be restricted to certain subdomains.

ssi_apps = [warpy.apps.five_moment.euler(name='ion_euler', gamma=5/3, fluid=ions),
warpy.apps.five_moment.field_source(name='ion_fsource', charge=1, mass=1,
skin_depth_norm=1, fluid=ions, field=em_fields)]
# ... similar for electrons and E&M fields
# create the spatial solver for all apps
ss = warpy.spatial_solvers.dg(name='dg', spatial_order=sorder, applications=ssi_apps+sse_apps+ssm_apps, on_subdomains=['plasma'])
Euler fluid flux.
Definition: five_moment.py:4
Lorentz force operator for 5-moment model.
Definition: five_moment.py:234

Define a writer

The writer defines what data to export and how. There are two primary types of writers:

Define boundary condition applications and variable adjusters

Non-periodic boundary conditions are defined using applications in warpy.apps. Once the apps are defined, they can be added to a boundary condition variable adjuster to restrict them to a given subdomain.

bc_apps = [warpy.apps.simple.bc_dirichlet(name='left boundary', values=[1], on_boundaries=['Left'], variable=fluid, components=['rho'])]
vas = [warpy.variable_adjusters.boundary_conditions(name='bcs', priority=1,spatial_order=sorder, applications=bc_apps, on_subdomains=['FluidSubdomain'])]
Dirichlet boundary condition.
Definition: simple.py:151
spatial solver - currently hardcoded for rk methods
Definition: boundary_conditions.py:4

Create a temporal solver

There are two type of temporal solvers currently supported by WARPXM: explicit and implicit Runge-Kutta schemes. These are both host actions.

temporal_solver = warpy.host_actions.erk(name='rk', scheme='SSPRK3', spatial_solvers=ss_list,
variable_adjusters=va_list)
Explicit Runge-Kutta temporal solver Note: Dormand45 currently will not work correctly with limiters ...
Definition: erk.py:5

Define a timestep controller

The timestep controller defines how the timestep can be updated (or not), as well as determining if a timestep needs to be re-tried.

dt = 1e-3
dt_controller = warpy.dt_calc.fixed_dt(dt)
Definition: fixed_dt.py:4

Create the DG simulation

The DG simulation object ties everything together, and represents the entire input file. There are a lot of extra optional parameters for configuring a simulation.

sim = warpy.dg_sim(name='advection',
meshes=[mesh],
initial_conditions=applicators,
temporal_solvers=[temporal_solver],
writers=[writer],
time=[0,1],
dt_controller=dt_controller,
flexible_writeout=False,
write_steps=100,
verbosity='info')
Discontinuous finite element RK simulation.
Definition: dg_sim.py:11

Run the simulation

The simulation object can be used to generate a directory where WARPXM will be run for a given test case. If gen_xdmf=True, then an XDMF file will be generated which can be used to post-process in ParaView or VisIt. Note that if VisIt is being used, then you need detect_nonscalar to be set to false.

sim.run('example', gen_xdmf=True, detect_nonscalar=True)

Debugging

When making modifications to warpy, it is a good idea to review that the generated WARPXM input file is being populated as expected.

For more targeted debugging, warpy objects can also be printed out individually as a string using the generate member function to get an idea of what the generated input file section would look like. Note that this may not work as intended for applications, since these depend on a variable index map which is generated from the associated spatial solver/variable adjuster. If you print out an application directly it will generate a dummy components map to use for indices. Once a spatial solver/variable adjuster has been added to a temporal solver it should have the correct variable index mapping.

# generate returns a list of blocks which will be added to the input file. Can pass a start/stop stage for items which may expand into multiple blocks.
print(fluid.generate(start=0,stop=3))
[<fluid>
Type = variable
Kind = distributed
ComponentNames = [rho]
BasisArraySet = firstOrder
</fluid>, <fluid_1>
Type = variable
Kind = distributed
ComponentNames = [rho]
BasisArraySet = firstOrder
</fluid_1>, <fluid_2>
Type = variable
Kind = distributed
ComponentNames = [rho]
BasisArraySet = firstOrder
</fluid_2>]

Examples

Setting up a generic MHD and/or Euler multiphysics simulation

See this user guide on things to consider when setting up a multiphysics plasma simulation with WARPXM.

More example problems

A collection of example warpy input files can be found in the examples folder.

Recommended examples include: