Post-process¶
This page describes the usage of the python module happi
for extracting, viewing
and post-processing simulation data. First, you need to install happi.
The module can be imported directly in python:
import happi
Open a simulation¶
In a python command line (or script), call the following function to open your Smilei simulation. Note that several simulations can be opened at once, as long as they correspond to several restarts of the same simulation.
- happi.Open(results_path='.', reference_angular_frequency_SI=None, show=True, verbose=True, scan=True, pint=True)¶
results_path
: path or list of paths to the directory-ies where the results of the simulation-s are stored. It can also contain wildcards, such as*
and?
in order to include several simulations at once.reference_angular_frequency_SI
: overrides the value of the simulation parameterreference_angular_frequency_SI
, in order to re-scale units.show
: ifFalse
, figures will not plot on screen. Make sure that you have not loaded another simulation or the matplotlib package. You may need to restart python.verbose
: ifFalse
, less information is printed while post-processing.scan
: ifFalse
, HDF5 output files are not scanned initially, and the namelist is not read.pint
: ifTrue
, happi attempts to load the Pint package and to use it for managing units.
- Returns: An object containing various methods to extract and manipulate the simulation
outputs, as described below.
Example:
S = happi.Open("path/to/my/results")
Once a simulation is opened, several methods are available to find information on the namelist or open various diagnostics. Checkout the namelist documentation to find out which diagnostics are included in Smilei: scalars, fields, probes, particle binning, trajectories and performances.
Extract namelist information¶
Once a simulation is opened as shown above, you can access the content of the namelist
using the attribute namelist
:
S = happi.Open("path/to/my/results") # Open a simulation
print(S.namelist.Main.timestep) # print the timestep
print(S.namelist.Main.geometry) # print the simulation dimensions
All the variables defined in the original namelist are copied into this variable.
Concerning components like Species, External fields or Probe diagnostics, of which several instances may exist, you can directly iterate over them:
for species in S.namelist.Species:
print("species "+species.name+" has mass "+str(species.mass))
You can also access to a specific component by referencing its number:
F = S.namelist.ExternalField[0] # get the first external field
print("An external field "+F.field+" was applied")
In the case of the species, you can also obtain a given species by its name:
species = S.namelist.Species["electron1"]
print("species "+species.name+" has mass "+str(species.mass))
Obtain diagnostic information¶
Print available diagnostics
Commands S.Scalar
, S.Field
, S.Probe
(etc.) will display general information
about the corresponding diagnostics in the simulation.
List available diagnostics
- getDiags(diagType)¶
Returns a list of available diagnostics of the given type
diagType
: The diagnostic type ("Field"
,"Probe"
, etc.)
Information on specific diagnostics
- getScalars()¶
Returns a list of available scalars.
- getTrackSpecies()¶
Returns a list of available tracked species.
- fieldInfo(diag)¶
diag
: the number or name of a Field diagnostic
Returns a dictionnary containing:
"diagNumber"
: the diagnostic number"diagName"
: the diagnostic name"fields"
: list of the available fields in this diagnostic. In the case ofAMcylindrical
geometry, this is a dictionnary with a list of modes for each field.
- probeInfo(diag)¶
diag
: the number or name of a Probe diagnostic
Returns a dictionnary containing:
"probeNumber"
: the diagnostic number"probeName"
: the diagnostic name"fields"
: list of the available fields in this diagnostic
- performanceInfo()¶
Returns a dictionnary containing:
"quantities_uint"
: a list of the available integer quantities"quantities_double"
: a list of the available float quantities"patch_arrangement"
: the type of patch arrangement"timesteps"
: the list of timesteps
Open a Scalar diagnostic¶
- Scalar(scalar=None, timesteps=None, units=[''], data_log=False, data_transform=None, **kwargs)¶
scalar
: The name of the scalar, or an operation on scalars, such as"Uelm+Ukin"
.timesteps
ortimestep_indices
: The requested range of timesteps.If omitted, all timesteps are used.
If one number given, the nearest timestep available is used.
If two numbers given, all the timesteps in between are used.
When using
timesteps
, provide the timesteps themselves, but when usingtimestep_indices
, provide their indices in the list of the available timesteps.units
: A unit specification (see Specifying units)data_log
:- If
True
, then \(\log_{10}\) is applied to the output.
data_transform
:- If this is set to a function, the function is applied to the output before plotting.
See also Other arguments for diagnostics
Example:
S = happi.Open("path/to/my/results")
Diag = S.Scalar("Utot")
Open a Field diagnostic¶
- Field(diagNumber=None, field=None, timesteps=None, subset=None, average=None, units=[''], data_log=False, data_transform=None, moving=False, export_dir=None, **kwargs)¶
timesteps
(ortimestep_indices
),units
,data_log
,data_transform
: same as before.diagNumber
: number orname
of the fields diagnostic- If not given, then a list of available diagnostic numbers is printed.
field
: The name of a field ("Ex"
,"Ey"
, etc.)- If not given, then a list of available fields is printed.The string can also be an operation between several fields, such as
"Jx+Jy"
.
subset
: A selection of coordinates to be extracted.- Syntax 1:
subset = { axis : location, ... }
Syntax 2:subset = { axis : [start, stop] , ... }
Syntax 3:subset = { axis : [start, stop, step] , ... }
axis
must be"x"
,"y"
,"z"
or"r"
.Only the data within the chosen axes’ selections is extracted.WARNING: THE VALUE OFstep
IS A NUMBER OF CELLS.Example:subset = {"y":[10, 80, 4]}
average
: A selection of coordinates on which to average.- Syntax 1:
average = { axis : "all", ... }
Syntax 2:average = { axis : location, ... }
Syntax 3:average = { axis : [start, stop] , ... }
axis
must be"x"
,"y"
,"z"
or"r"
.The chosen axes will be removed:- With syntax 1, an average is performed over all the axis.- With syntax 2, only the bin closest tolocation
is kept.- With syntax 3, an average is performed fromstart
tostop
.Example:average = {"x":[4,5]}
will average for \(x\) within [4,5].
moving
: IfTrue
, plots will display the X coordinates evolving according to the moving windowexport_dir
: The directory where to export VTK files.See also Other arguments for diagnostics
In the case of an azimuthal mode cylindrical geometry (
AMcylindrical
), additional argument are available. You must choose one oftheta
orbuild3d
, defined below, in order to construct fields from their complex angular Fourier modes. In addition, themodes
argument is optional.theta
: An angle (in radians)- Calculates the field in a plane passing through the \(r=0\) axisand making an angle
theta
with the \(xy\) plane.
build3d
: A list of three ranges- Calculates the field interpolated in a 3D \(xyz\) grid.Each range is a list
[start, stop, step]
indicating the beginning,the end and the step of this grid.
modes
: An integer or a list of integers- Only these modes numbers will be used in the calculation. If omited, all modes are used.
Example:
S = happi.Open("path/to/my/results")
Diag = S.Field(0, "Ex", average = {"x":[4,5]}, theta=math.pi/4.)
Open a Probe diagnostic¶
- Probe(probeNumber=None, field=None, timesteps=None, subset=None, average=None, units=[''], data_log=False, data_transform=None, **kwargs)¶
timesteps
(ortimestep_indices
),units
,data_log
,data_transform
,export_dir
: same as before.probeNumber
: number orname
of the probe (the first one has number 0).- If not given, a list of available probes is printed.
field
: name of the field ("Bx"
,"By"
,"Bz"
,"Ex"
,"Ey"
,"Ez"
,"Jx"
,"Jy"
,"Jz"
or"Rho"
).- If not given, a list of available fields is printed.The string can also be an operation between several fields, such as
"Jx+Jy"
.
subset
andaverage
are very similar to those ofField()
, but they can only have the axes:"axis1"
,"axis2"
and"axis3"
. For instance,average={"axis1":"all"}
. Note that the axes are not necessarily \(x\), \(y\) or \(z\) because the probe mesh is arbitrary.See also Other arguments for diagnostics
Example:
S = happi.Open("path/to/my/results")
Diag = S.Probe(0, "Ex")
- Probe.changeField(field)¶
In cases where happi’s performance is an issue, it is possible to switch between different fields of an open
Probe
diagnostic using this method. Thefield
argument is the same as inProbe(...)
above.
Open a ParticleBinning diagnostic¶
- ParticleBinning(diagNumber=None, timesteps=None, subset=None, average=None, units=[''], data_log=False, data_transform=None, **kwargs)¶
timesteps
(ortimestep_indices
),units
,data_log
,data_transform
,export_dir
: same as before.diagNumber
: number orname
of the particle binning diagnostic (starts at 0).- If not given, a list of available diagnostics is printed.It can also be an operation between several diagnostics.For example,
"#0/#1"
computes the division by diagnostics 0 and 1.
subset
is similar to that ofField()
, although the axis must be one of"x"
,"y"
,"z"
,"px"
,"py"
,"pz"
,"p"
,"gamma"
,"ekin"
,"vx"
,"vy"
,"vz"
,"v"
or"charge"
.WARNING: With the syntax
subset={axis:[start, stop, step]}
, the value ofstep
is a number of bins.
average
: a selection of coordinates on which to average the data.- Syntax 1:
average = { axis : "all", ... }
Syntax 2:average = { axis : location, ... }
Syntax 3:average = { axis : [begin, end] , ... }
axis
must be"x"
,"y"
,"z"
,"px"
,"py"
,"pz"
,"p"
,"gamma"
,"ekin"
,"vx"
,"vy"
,"vz"
,"v"
or"charge"
.The chosen axes will be removed:- With syntax 1, an average is performed over all the axis.- With syntax 2, only the bin closest tolocation
is kept.- With syntax 3, an average is performed betweenbegin
andend
.Example:average={"x":[4,5]}
will average all the data for x within [4,5].
See also Other arguments for diagnostics
Example:
S = happi.Open("path/to/my/results")
Diag = S.ParticleBinning(1)
Units of the results:
Open a Screen diagnostic¶
- Screen(diagNumber=None, timesteps=None, subset=None, average=None, units=[''], data_log=False, data_transform=None, **kwargs)¶
timesteps
(ortimestep_indices
),units
,data_log
,data_transform
,export_dir
: same as before.diagNumber
,subset
andaverage
: identical to that of ParticleBinning diagnostics.See also Other arguments for diagnostics
Example:
S = happi.Open("path/to/my/results")
Diag = S.Screen(0)
Open a RadiationSpectrum diagnostic¶
- ParticleBinning(diagNumber=None, timesteps=None, subset=None, average=None, units=[''], data_log=False, data_transform=None, **kwargs)¶
timesteps
(ortimestep_indices
),units
,data_log
,data_transform
,export_dir
: same as before.diagNumber
,subset
andaverage
: identical to that of ParticleBinning diagnostics.See also Other arguments for diagnostics
Example:
S = happi.Open("path/to/my/results")
Diag = S.RadiationSpectrum(0)
Note
The resulting spectral power is in units of \(\omega_r\). If additional axes are used, the power spectrum is divided by the size of the bins of each axes.
Open a TrackParticles diagnostic¶
- TrackParticles(species=None, select='', axes=[], timesteps=None, sort=True, length=None, units=[''], **kwargs)¶
timesteps
(ortimestep_indices
),units
,export_dir
: same as before.species
: the name of a tracked-particle species. If omitted, a list of available tracked-particle species is printed.select
: Instructions for selecting particles among those available. A detailed explanation is provided belowaxes
: A list of axes for plotting the trajectories or obtaining particle data.Each axis is one of the
attributes
defined in the namelist. In addition, when there is a moving window, the axis"moving_x"
is automatically available.Example:axes = ["x"]
corresponds to \(x\) versus time.Example:axes = ["x","y"]
correspond to 2-D trajectories.Example:axes = ["x","px"]
correspond to phase-space trajectories.
sort
: may be eitherFalse
: the particles are not sorted by ID. This can save significant time, but prevents plotting, exporting to VTK, and theselect
argument. OnlygetData
anditerParticles
are available in this mode. Read this for more information on particle IDs.True
: the particles are sorted in a new file, unless this file already exists. If it does, sorted particles are directly read from the sorted file.A string for selecting particles (same syntax as
select
): only selected particles are sorted in a new file. The file name must be defined in the argumentsorted_as
. Iftimesteps
is used, only selected timesteps will be included in the created file.
sorted_as
: a keyword that defines the new sorted file name (whensort
is a selection) or refers to a previously user-defined sorted file name (whensort
is not given).length
: The length of each plotted trajectory, in number of timesteps.See also Other arguments for diagnostics
Example:
S = happi.Open("path/to/my/results")
Diag = S.TrackParticles("electrons", axes=["px","py"])
Detailed explanation of the select
parameter
times
is a condition on timesteps t
, for instance t>50
.condition
is a condition on particles properties (x
, y
, z
, px
, py
, pz
), for instance px>0
.- Syntax 1:
select="any(times, condition)"
Selects particles satisfyingcondition
for at least one of thetimes
.For example,select="any(t>0, px>1.)"
selects those reaching \(p_x>1\) at some point. - Syntax 2:
select="all(times, condition)"
Selects particles satisfyingcondition
at alltimes
.For example,select="all(t<40, px<1)"
selects those having \(p_x<1\) until timestep 40. - Syntax 3:
select=[ID1, ID2, ...]
Selects the provided particle IDs. - It is possible to make logical operations:
+
is OR;*
is AND;~
is NOT.For example,select="any((t>30)*(t<60), px>1) + all(t>0, (x>1)*(x<2))"
Open a NewParticles diagnostic¶
- NewParticles(species=None, select='', axes=[], units=[''], **kwargs)¶
units
: same as before.species
: same as forTrackParticles
axes
: same as forTrackParticles
, with the addition of another axist
that represents the time when the particle was born.select
: Instructions for selecting particles among those available. It must be a condition on particles propertiesaxes
, for instancepx>0
. It is possible to make logical operations:+
is OR;*
is AND;~
is NOT.Example:select="(x>1)*(x<2)"
It is also possible to select directly a list of IDs.
Example:select=[ID1, ID2, ...]
Open a Performances diagnostic¶
The post-processing of the performances diagnostic may be achieved in three different
modes: raw
, map
, or histogram
, described further below. You must choose one
and only one mode between those three.
- Performances(raw=None, map=None, histogram=None, timesteps=None, units=[''], data_log=False, data_transform=None, species=None, cumulative=True, **kwargs)¶
timesteps
,units
,data_log
,data_transform
,export_dir
: same as before.raw
: The name of a quantity, or an operation between them (see quantities below). The requested quantity is listed for each process.map
: The name of a quantity, or an operation between them (see quantities below). The requested quantity is mapped vs. space coordinates (1D and 2D only).histogram
: the list["quantity", min, max, nsteps]
. Makes a histogram of the requested quantity betweenmin
anmax
, withnsteps
bins. The"quantity"
may be an operation between the quantities listed further below.cumulative
: may beTrue
for timers accumulated for the duration of the simulation, orFalse
for timers reset to 0 at each output.See also Other arguments for diagnostics
Quantities at the MPI-process level (contain many patches):
hindex
: the starting index of each proc in the hilbert curve
number_of_cells
: the number of cells in each proc
number_of_particles
: the total number of non-frozen macro-particles in each proc (includes all species)
number_of_frozen_particles
: the number of frozen particles in each proc
total_load
: the load of each proc (number of macro-particles and cells weighted by cell_load coefficients)
timer_global
: global simulation time (only available for proc 0)
timer_particles
: time spent computing particles by each proc
timer_maxwell
: time spent solving maxwell by each proc
timer_envelope
: time spent solving the envelope propagation by each proc
timer_densities
: time spent projecting densities by each proc
timer_collisions
: time spent computing collisions by each proc
timer_movWindow
: time spent handling the moving window by each proc
timer_loadBal
: time spent balancing the load by each proc
timer_partMerging
: time spent merging particles by each proc
timer_syncPart
: time spent synchronzing particles by each proc
timer_syncField
: time spent synchronzing fields by each proc
timer_syncDens
: time spent synchronzing densities by each proc
timer_syncSusceptibility
: time spent synchronzing susceptibility by each proc
timer_diags
: time spent by each proc calculating and writing diagnostics
timer_total
: the sum of all timers above (except timer_global)
memory_total
: the total memory (RSS) used by the process in GB
memory_peak
: the peak memory (peak RSS) used by the process in GBWARNING: The timers
loadBal
anddiags
include global communications. This means they might contain time doing nothing, waiting for other processes. Thesync***
timers contain proc-to-proc communications, which also represents some waiting time.
Quantities at the patch level:
This requires
patch_information
in the namelist.
mpi_rank
: the MPI rank that contains the current patch
vecto
: the mode of the specified species in the current patch (vectorized of scalar) when the adaptive mode is activated. Here thespecies
argument has to be specified.WARNING: The patch quantities are only compatible with the
raw
mode and only in3Dcartesian
geometry
. The result is a patch matrix with the quantity on each patch.
Example: performance diagnostic at the MPI level:
S = happi.Open("path/to/my/results")
Diag = S.Performances(map="total_load")
Example: performance diagnostic at the patch level:
S = happi.Open("path/to/my/results")
Diag = S.Performances(raw="vecto", species="electron")
Specifying units¶
By default, all the diagnostics data is in code units (see Units).
To change the units, all the methods Scalar()
,
Field()
, Probe()
,
ParticleBinning()
and
TrackParticles()
support a units
argument.
It has three different syntaxes:
A list, for example
units = ["um/ns", "feet", "W/cm^2"]
In this case, any quantity found to be of the same dimension as one of these units will be converted.
A dictionary, for example
units = {"x":"um", "y":"um", "v":"Joule"}
In this case, we specify the units separately for axes
x
andy
, and for the data valuesv
.A
Units
object, for exampleunits = happi.Units("um/ns", "feet", x="um")
This version combines the two previous ones.
Requirements for changing units
The Pint module.
To obtain units in a non-normalized system (e.g. SI), the simulation must have the parameter
reference_angular_frequency_SI
set to a finite value. Otherwise, this parameter can be set during post-processing as an argument to thehappi.Open()
function.
Other arguments for diagnostics¶
All diagnostics above can use additional keyword arguments (kwargs
)
to manipulate the plotting options:
figure
: The figure number that is passed to matplotlib.vmin
,vmax
: data value limits.vsym
: makes data limits symmetric about 0 (vmin
andvmax
are ignored), and sets the colormap tosmileiD
.If
vsym = True
, autoscale symmetrically.If
vsym
is a number, limits are set to [-vsym
,vsym
].
xmin
,xmax
,ymin
,ymax
: axes limits.xfactor
,yfactor
: factors to rescale axes.xoffset
,yoffset
: numerical values to offset the coordinates. These values must be given in the original (normalized) units, i.e. not acounting for the factors above or for unit conversion.title
: a string that replaces the plot title (or the y-label in a 1D plot). The current simulation time can be included with the placeholders{time}
and{time_units}
, together with formatting instructions conforming to python’s string formatter. For instance:title = "Density @ $t = {time:.0f} {time_units}$"
.side
:"left"
(by default) or"right"
puts the y-axis on the left- or the right-hand-side.transparent
:None
(by default),"over"
,"under"
,"both"
, or a function. The colormap becomes transparent over, under, or outside both the boundaries set byvmin
andvmax
. This argument may be set instead to a function mapping the data value \(\in [0,1]\) to the transparency \(\in [0,1]\). For instancelambda x: 1-x
.Other Matplotlib arguments listed in Advanced plotting options.
Obtain the data¶
- Scalar.getData(timestep=None)¶
- Field.getData(timestep=None)¶
- Probe.getData(timestep=None)¶
- ParticleBinning.getData(timestep=None)¶
- Screen.getData(timestep=None)¶
- TrackParticles.getData(timestep=None)¶
Returns a list of the data arrays (one element for each timestep requested). In the case of
TrackParticles
, this method returns a dictionary containing one entry for each axis, and ifsort==False
, these entries are included inside an entry for each timestep.timestep
, if specified, is the only timestep number that is read and returned.
Example:
S = happi.Open("path/to/results") # Open the simulation Diag = S.Field(0, "Ex") # Open Ex in the first Field diag result = Diag.getData() # Get list of Ex arrays (one for each time)
- Scalar.getTimesteps()¶
- Field.getTimesteps()¶
- Probe.getTimesteps()¶
- ParticleBinning.getTimesteps()¶
- Screen.getTimesteps()¶
- TrackParticles.getTimesteps()¶
Returns a list of the timesteps requested.
- Scalar.getTimes()¶
- Field.getTimes()¶
- Probe.getTimes()¶
- ParticleBinning.getTimes()¶
- Screen.getTimes()¶
- TrackParticles.getTimes()¶
Returns the list of the times requested. By default, times are in the code’s units, but are converted to the diagnostic’s units defined by the
units
argument, if provided.
- Scalar.getAxis(axis)¶
- Field.getAxis(axis, timestep)¶
- Probe.getAxis(axis)¶
- ParticleBinning.getAxis(axis, timestep)¶
- Screen.getAxis(axis, timestep)¶
Returns the list of positions of the diagnostic data along the requested axis. If the axis is not available, returns an empty list. By default, axis positions are in the code’s units, but are converted to the diagnostic’s units defined by the
units
argument, if provided.axis
: the name of the requested axis.For
Field
: this is"x"
,"y"
or"z"
For
Probe
: this is"axis1"
,"axis2"
or"axis3"
For
ParticleBinning
andScreen
: this is thetype
of theaxes
defined in the namelist
timestep
: The timestep at which the axis is obtained. Only matters inParticleBinning
,Screen
andRadiationSpectrum
whenauto
axis limits are requested; or inField
whenmoving=True
.
- TrackParticles.iterParticles(timestep, chunksize=1)¶
This method, specific to the tracked particles, provides a fast iterator on chunks of particles for a given timestep. The argument
chunksize
is the number of particles in each chunk. Note that the data is not ordered by particle ID, meaning that particles are not ordered the same way from one timestep to another.The returned quantity for each iteration is a python dictionary containing key/value pairs
axis:array
, whereaxis
is the name of the particle characteristic ("x"
,"px"
, etc.) andarray
contains the corresponding particle values.Example:
S = happi.Open("path/to/my/results") # Open the simulation Diag = S.TrackParticles("my_particles") # Open the tracked particles npart = 0 sum_px = 0. # Loop particles of timestep 100 by chunks of 10000 for particle_chunk in Diag.iterParticles(100, chunksize=10000): npart += particle_chunk["px"].size sum_px += particle_chunk["px"].sum() # Calculate the average px mean_px = sum_px / npart
- Field.getXmoved(timestep)¶
Specific to Field diagnostics, this method returns the displacement of the moving window at the required
timestep
.
Export 2D or 3D data to VTK¶
- Field.toVTK(numberOfPieces=1)¶
- Probe.toVTK(numberOfPieces=1)¶
- ParticleBinning.toVTK(numberOfPieces=1)¶
- Performances.toVTK(numberOfPieces=1)¶
- Screen.toVTK(numberOfPieces=1)¶
- TrackParticles.toVTK(rendering='trajectory', data_format='xml')¶
Converts the data from a diagnostic object to the vtk format. Note the
export_dir
argument available for each diagnostic (see above).numberOfPieces
: the number of files into which the data will be split.rendering
: the type of output in the case ofTrackParticles()
:"trajectory"
: show particle trajectories. One file is generated for all trajectories. This rendering requires the particles to be sorted."cloud"
: show a cloud of particles. One file is generated for each iteration. This rendering can be used without sorting the particles.
data_format
: the data formatting in the case ofTrackParticles()
, either"vtk"
or"xml"
. The format"vtk"
results in ascii.
Example for tracked particles:
S = happi.Open("path/to/my/results") tracked_particles = S.TrackParticles("electron", axes=["x","y","z","px","py","pz","Id"], timesteps=[1,10]) # Create cloud of particles in separate files for each iteration tracked_particles.toVTK(rendering="cloud",data_format="xml"); # Create trajectory in a single file tracked_particles.toVTK(rendering="trajectory",data_format="xml");
Plot the data at one timestep¶
This is the first method to plot the data. It produces a static image of the data at one given timestep.
- Scalar.plot(...)¶
- Field.plot(...)¶
- Probe.plot(...)¶
- ParticleBinning.plot(...)¶
- TrackParticles.plot(...)¶
- Screen.plot(...)¶
All these methods have the same arguments described below.
- plot(timestep=None, saveAs=None, axes=None, dpi=200, **kwargs)¶
- If the data is 1D, it is plotted as a curve.If the data is 2D, it is plotted as a map.If the data is 0D, it is plotted as a curve as function of time.
timestep
: The iteration number at which to plot the data.saveAs
: name of a directory where to save each frame as figures. You can even specify a filename such asmydir/prefix.png
and it will automatically make successive files showing the timestep:mydir/prefix0.png
,mydir/prefix1.png
, etc.axes
: Matplotlib’s axes handle on which to plot. If None, make new axes.dpi
: the number of dots per inch forsaveAs
.
You may also have keyword-arguments (
kwargs
) described in Other arguments for diagnostics.
Example:
S = happi.Open("path/to/my/results")
S.ParticleBinning(1).plot(timestep=40, vmin=0, vmax=1e14)
Plot the data streaked over time¶
This second type of plot works only for 1D data. All available timesteps are streaked to produce a 2D image where the second axis is time.
- Scalar.streak(...)¶
- Field.streak(...)¶
- Probe.streak(...)¶
- ParticleBinning.streak(...)¶
- TrackParticles.streak(...)¶
- Screen.streak(...)¶
All these methods have the same arguments described below.
- streak(saveAs=None, axes=None, **kwargs)¶
All arguments are identical to those of
plot
, with the exception oftimestep
.
Example:
S = happi.Open("path/to/my/results")
S.ParticleBinning(1).streak()
Animated plot¶
This third plotting method animates the data over time.
- Scalar.animate(...)¶
- Field.animate(...)¶
- Probe.animate(...)¶
- ParticleBinning.animate(...)¶
- TrackParticles.animate(...)¶
- Screen.animate(...)¶
All these methods have the same arguments described below.
- animate(movie='', fps=15, dpi=200, saveAs=None, axes=None, **kwargs)¶
All arguments are identical to those of
streak
, with the addition of:movie
: name of a file to create a movie, such as"movie.avi"
or"movie.gif"
. Ifmovie=""
no movie is created.fps
: number of frames per second (only if movie requested).dpi
: number of dots per inch for bothmovie
andsaveAs
Example:
S = happi.Open("path/to/my/results")
S.ParticleBinning(1).animate()
Plot with a slider¶
This methods provides an interactive slider to change the time.
- Scalar.slide(...)¶
- Field.slide(...)¶
- Probe.slide(...)¶
- ParticleBinning.slide(...)¶
- TrackParticles.slide(...)¶
- Screen.slide(...)¶
All these methods have the same arguments described below.
- slide(axes=None, **kwargs)¶
See
plot
for the description of the arguments.
Example:
S = happi.Open("path/to/my/results")
S.ParticleBinning(1).slide(vmin=0)
Simultaneous plotting of multiple diagnostics¶
- happi.multiPlot(diag1, diag2, ..., **kwargs)¶
Makes an animated figure containing several plots (one for each diagnostic). If all diagnostics are of similar type, they may be overlayed on only one plot.
diag1
,diag2
, etc.- Diagnostics prepared by
Scalar()
,Field()
,Probe()
, etc.
Keyword-arguments
kwargs
are:figure
: The figure number that is passed to matplotlib (default is 1).shape
: The arrangement of plots inside the figure. For instance,[2, 1]
makes two plots stacked vertically, and[1, 2]
makes two plots stacked horizontally. If absent, stacks plots vertically.legend_font
: dictionnary to set the legend’s font properties, such as{'size':15, 'weight':'bold', 'family':'serif', 'color':'k'}
.movie
: filename to create a movie.fps
: frames per second for the movie.dpi
: resolution of themovie
orsaveAs
.saveAs
: name of a directory where to save each frame as figures. You can even specify a filename such asmydir/prefix.png
and it will automatically make successive files showing the timestep:mydir/prefix0.png
,mydir/prefix1.png
, etc.skipAnimation
: if True, plots only the last frame.timesteps
: same as thetimesteps
argument of theplot()
method.
- happi.multiSlide(diag1, diag2, ..., **kwargs)¶
Identical to
happi.multiPlot
but uses a time slider instead of an animation.diag1
,diag2
, etc.- Diagnostics prepared by
Scalar()
,Field()
,Probe()
, etc.
figure
,shape
, andlegend_font
: same as inhappi.multiPlot
.
Example:
S = happi.Open("path/to/my/results")
A = S.Probe(probeNumber=0, field="Ex")
B = S.ParticleBinning(diagNumber=1)
happi.multiPlot( A, B, figure=1 )
This plots a Probe and a ParticleBinning on the same figure, and makes an animation for all available timesteps.
Note
To plot several quantities on the same graph, you can try shape=[1,1]
.
One diagnostic may have the option side="right"
to use the right-hand-side axis.
Advanced plotting options¶
In addition to figure
, vmin
, vmax
, xmin
, xmax
, ymin
and ymax
,
there are many more optional arguments. They are directly passed to the matplotlib package.
For the figure: figsize
, dpi
, facecolor
, edgecolor
Please refer to matplotlib’s figure options.
For the axes frame: aspect
, axis_facecolor
, frame_on
, position
,
visible
, xlabel
, xscale
, xticklabels
, xticks
,
ylabel
, yscale
, yticklabels
, yticks
, zorder
Please refer to matplotlib’s axes options: the same as functions starting with
set_
listed here.
For the lines: color
, dashes
, drawstyle
, fillstyle
,
label
, linestyle
, linewidth
,
marker
, markeredgecolor
, markeredgewidth
,
markerfacecolor
, markerfacecoloralt
, markersize
, markevery
,
visible
, zorder
Please refer to matplotlib’s line options.
For the image: cmap
, aspect
, interpolation
, norm
Please refer to matplotlib’s image options.
For the colorbar: cbaspect
, orientation
, fraction
, pad
,
shrink
, anchor
, panchor
, extend
, extendfrac
, extendrect
,
spacing
, ticks
, format
, drawedges
, size
, clabel
Please refer to matplotlib’s colorbar options.
For the tick number format: style_x
, scilimits_x
, useOffset_x
,
style_y
, scilimits_y
, useOffset_y
Please refer to matplotlib’s tick label format.
For fonts: title_font
, xlabel_font
, xticklabels_font
,
ylabel_font
, yticklabels_font
, colorbar_font
These options are dictionnaries that may contain the entries available in matplotlib’s font properties, for instance:
title_font = {'size': 15, 'weight': 'bold', 'family':'serif', 'color': 'k'}
Example:
To choose a gray colormap of the image, use
cmap="gray"
:S = happi.Open("path/to/my/results") S.ParticleBinning(0, figure=1, cmap="gray") .plot()
Many colormaps are available from the matplotlib package. With
cmap=""
, you will get a list of available colormaps. Smilei’s default colormaps are:smilei
,smilei_r
,smileiD
andsmileiD_r
.
Update the plotting options¶
Other tools in happi
¶
- happi.openNamelist(namelist)¶
Reads a namelist and stores all its content in the returned object.
namelist
: the path to the namelist.
Example:
namelist = happi.openNamelist("path/no/my/namelist.py")
print namelist.Main.timestep