Table Of Contents

Previous topic

CoreMirror class

Next topic

Pysic Fortran module

This Page

Pysic Utility modules

Modules for providing utility tools and parameters.

Plot

The plot utility defines a group of tools for exploring and plotting energy and force landscapes. It uses the matplotlib library.

pysic.utility.plot.plot_abs_force_on_line(index, system, direction=None, length=None, steps=100, start=None, end=None, lims=[-10000000000.0, 10000000000.0])[source]

Plots the absolute value of the force on a particle as a function of the position.

The method probes the system by moving a single particle on a line and recording the force. A plot is drawn. Also a tuple containing arrays of the distance traveled and the recorded forces is returned.

After the operation is complete, the initial structure is restored.

Parameters:

index: integer
index of the particle to be moved
system: ASE Atoms object
the structure to be explored
direction: double 3-vector
the direction where the atom is moved
length: double
the distance moved
steps: integer
number of points (taken uniformly on the movement path) for measuring the force
start: double 3-vector (array or list)
starting point for the trajectory - if not specified, the position of the particle in ‘system’ is used
end: double 3-vector (array or list)
end point for the trajectory - alternative for direction and length (will override them)
lims: double 2-vector (array or list)
lower and upper truncation limits - if a recorded value is smaller than the lower limit or larger than the upper, it is replaced by the corresponding truncation value
pysic.utility.plot.plot_abs_force_on_plane(index, system, directions, lengths, steps=100, start=None, lims=[-10000000000.0, 10000000000.0])[source]

Plots the absolute value of force on a particle as a function of the position.

The method probes the system by moving a single particle on a plane and recording the force. A contour plot is drawn. Also a tuple containing arrays of the distances traveled on the plane and the recorded forces is returned.

After the operation is complete, the initial structure is restored.

Parameters:

index: integer
index of the particle to be moved
system: ASE Atoms object
the structure to be explored
directions: double 2x3-matrix
The directions where the atom is moved - i.e., the vectors defining the plane. If the second vector is not perpendicular to the first, the normal component is automatically used instead.
lengths: double 2-vector
the distances moved in the given directions
steps: integer
number of points (taken uniformly on the movement plane) for measuring the force
start: double 3-vector (array or list)
starting point for the trajectory, i.e., a corner for the plane to be probed - if not specified, the position of the particle in ‘system’ is used
lims: double 2-vector (array or list)
lower and upper truncation limits - if a recorded value is smaller than the lower limit or larger than the upper, it is replaced by the corresponding truncation value
pysic.utility.plot.plot_energy_on_line(index, system, direction=None, length=None, steps=100, start=None, end=None, lims=[-10000000000.0, 10000000000.0])[source]

Plots the energy of the system as a function of the position of a single particle.

The method probes the system by moving a single particle on a line and recording the energy. A plot is drawn. Also a tuple containing arrays of the distance traveled and the recorded energies is returned.

After the operation is complete, the initial structure is restored.

Parameters:

index: integer
index of the particle to be moved
system: ASE Atoms object
the structure to be explored
direction: double 3-vector
the direction where the atom is moved
length: double
the distance moved
steps: integer
number of points (taken uniformly on the movement path) for measuring the energy
start: double 3-vector (array or list)
starting point for the trajectory - if not specified, the position of the particle in ‘system’ is used
end: double 3-vector (array or list)
end point for the trajectory - alternative for direction and length (will override them)
lims: double 2-vector (array or list)
lower and upper truncation limits - if a recorded value is smaller than the lower limit or larger than the upper, it is replaced by the corresponding truncation value
pysic.utility.plot.plot_energy_on_plane(index, system, directions, lengths, steps=100, start=None, lims=[-10000000000.0, 10000000000.0])[source]

Plots the energy of the system as a function of the position of a particle.

The method probes the system by moving a single particle on a plane and recording the energy. A contour plot is drawn. Also a tuple containing arrays of the distances traveled on the plane and the recorded energies is returned.

After the operation is complete, the initial structure is restored.

Parameters:

index: integer
index of the particle to be moved
system: ASE Atoms object
the structure to be explored
directions: double 2x3-matrix
The directions where the atom is moved - i.e., the vectors defining the plane. If the second vector is not perpendicular to the first, the normal component is automatically used instead.
lengths: double 2-vector
the distances moved in the given directions
steps: integer
number of points (taken uniformly on the movement plane) for measuring the energy
start: double 3-vector (array or list)
starting point for the trajectory, i.e., a corner for the plane to be probed - if not specified, the position of the particle in ‘system’ is used
lims: double 2-vector (array or list)
lower and upper truncation limits - if a recorded value is smaller than the lower limit or larger than the upper, it is replaced by the corresponding truncation value
pysic.utility.plot.plot_force_component_on_plane(index, system, directions, lengths, component, steps=100, start=None, lims=[-10000000000.0, 10000000000.0])[source]

Plots the projected component of force on a particle as a function of the position.

The method probes the system by moving a single particle on a plane and recording the force. The component of the force projected on a given vector is recorded. A contour plot is drawn. Also a tuple containing arrays of the distances traveled on the plane and the recorded forces is returned.

After the operation is complete, the initial structure is restored.

Parameters:

index: integer
index of the particle to be moved
system: ASE Atoms object
the structure to be explored
directions: double 2x3-matrix
The directions where the atom is moved - i.e., the vectors defining the plane. If the second vector is not perpendicular to the first, the normal component is automatically used instead.
lengths: double 2-vector
the distances moved in the given directions
component: double 3-vector
the direction on which the force is projected - e.g., if component is [1,0,0], the x-component is recorded
steps: integer
number of points (taken uniformly on the movement plane) for measuring the force
start: double 3-vector (array or list)
starting point for the trajectory, i.e., a corner for the plane to be probed - if not specified, the position of the particle in ‘system’ is used
lims: double 2-vector (array or list)
lower and upper truncation limits - if a recorded value is smaller than the lower limit or larger than the upper, it is replaced by the corresponding truncation value
pysic.utility.plot.plot_tangent_force_on_line(index, system, direction=None, length=None, steps=100, start=None, end=None, lims=[-10000000000.0, 10000000000.0])[source]

Plots the tangential force on a particle as a function of the position.

The method probes the system by moving a single particle on a line and recording the force tangent. A plot is drawn. Also a tuple containing arrays of the distance traveled and the recorded forces is returned.

After the operation is complete, the initial structure is restored.

Parameters:

index: integer
index of the particle to be moved
system: ASE Atoms object
the structure to be explored
direction: double 3-vector
the direction where the atom is moved
length: double
the distance moved
steps: integer
number of points (taken uniformly on the movement path) for measuring the energy
start: double 3-vector (array or list)
starting point for the trajectory - if not specified, the position of the particle in ‘system’ is used
end: double 3-vector (array or list)
end point for the trajectory - alternative for direction and length (will override them)
lims: double 2-vector (array or list)
lower and upper truncation limits - if a recorded value is smaller than the lower limit or larger than the upper, it is replaced by the corresponding truncation value
pysic.utility.plot.plot_tangent_force_on_plane(index, system, directions, lengths, steps=100, start=None, lims=[-10000000000.0, 10000000000.0])[source]

Plots the absolute value of the tangent component of force on a particle as a function of the position.

The method probes the system by moving a single particle on a plane and recording the force. The force is projected on the same plane, and the absolute value of the projection is calculated. A contour plot is drawn. Also a tuple containing arrays of the distances traveled on the plane and the recorded forces is returned.

After the operation is complete, the initial structure is restored.

Parameters:

index: integer
index of the particle to be moved
system: ASE Atoms object
the structure to be explored
directions: double 2x3-matrix
The directions where the atom is moved - i.e., the vectors defining the plane. If the second vector is not perpendicular to the first, the normal component is automatically used instead.
lengths: double 2-vector
the distances moved in the given directions
steps: integer
number of points (taken uniformly on the movement plane) for measuring the force
start: double 3-vector (array or list)
starting point for the trajectory, i.e., a corner for the plane to be probed - if not specified, the position of the particle in ‘system’ is used
lims: double 2-vector (array or list)
lower and upper truncation limits - if a recorded value is smaller than the lower limit or larger than the upper, it is replaced by the corresponding truncation value

Outliers

This module contains tools for statistical analysis of atomic structures. It can be used for calculating the bond length and angle distributions in a given structure according to bond rules given by the user. Log-likelihood distributions are generated based on these distributions, which can be further grouped or used to create false color representations of the structure.

This is an example of how to create the distributions:

import pysic
from pysic.utility.outliers import *
from ase.io import read

def do_everything(inputfile, outputfile, radii, periodic_directions, cell_lengths):
   """Stitch everything together.
   """

   system = read(inputfile)
   system.set_pbc(periodic_directions)
   system.set_cell(cell_lengths)
   structure = Structure(system)

   structure.add_bond(['Si','Si'], radii['Si']*2*0.6)
   structure.add_bond(['Si','O'], (radii['Si']+radii['O'])*0.6)
   structure.add_bond(['O','O'], radii['O']*2*0.6)
   structure.add_bond(['Si','H'], (radii['Si']+radii['H'])*0.6)
   structure.add_bond(['H','O'], (radii['H']+radii['O'])*0.6)
   structure.create_neighbor_lists()
   angles = structure.get_all_angles()
   distances = structure.get_all_distances()
   max_n = np.shape(system)[0]

   for a in angles:
       print "a", a.center_index, a.type1, a.type2, a.type3, a.value

   for d in distances:
       print "d", d.primary_index, d.type1, d.type2, d.value

   angle_distribs, dist_distribs = get_distributions(angles, distances, radii)
   a_logls, d_logls = get_log_likelihoods(angles, distances, angle_distribs,
                                          dist_distribs, max_n)
   write_to_file(outputfile, cell_lengths, system.get_chemical_symbols(),
                   system.get_positions(), a_logls, d_logls)


inputfile           = 'sio2.xyz'
outputfile          = 'sio2_analysis.txt'
periodic_directions = [True, True, True]
cell_lengths = [14.835, 14.835, 14.835]

radii = {'O' : 1.52,
         'Si': 2.10}


# All set? Run.
do_everything(inputfile, outputfile, radii, periodic_directions, cell_lengths)

The atomic system is first read from an xyz file using ASE read method. It is then translated to a Structure class, which the outliers tools use, and bonding rules are created based on element types. Once bonding has been defined, bond lengths, bond angles, and their distributions are calculated with the outliers tools. The results are printed in a file in this example, but it would of course be possible to access them directly in the script for further analysis.

Structure class

class pysic.utility.outliers.Structure(atoms)[source]
add_bond(elements, cutoff)[source]

Adds a bond between atoms of the given elements, up to the cutoff separation.

create_neighbor_lists()[source]

Creates neighbor lists for the structure.

get_all_angles()[source]

Returns a list of all 3-atom angles as Angle objects.

get_all_distances()[source]

Returns a list of all 2-atom distances as Distance objects.

get_bond_length(elements)[source]

Returns the stored bond length for the given pair of elements.

get_distances(index)[source]

Returns the distances between a given atom and its neighbors. (Sorted.)

get_neighbors(index)[source]

Returns the indices of neighboring atoms for the given atom. (Sorted according to distance.)

get_separations(index)[source]

Returns the separation vectors from a given atom to its neighbors. (Sorted according to distance.)

Distance class

class pysic.utility.outliers.Distance(primary_index, type1, type2, distance)[source]

Angle class

class pysic.utility.outliers.Angle(center_index, type1, type2, type3, angle)[source]

Outliers module

Detection of irregular atomic neighborhoods by bond distance and angle analysis.

pysic.utility.outliers.angle(A, O, B)[source]

Return the angle between vectors OA and OB

Parameters O, A, B: coordinates in 3d-space

pysic.utility.outliers.get_distributions(angles, distances, radii)[source]

Return observed log-distributions of angles and distances

Distributions are obtained for all combinations of observed atom types. They are histogram-based, using the breakpoints specified by parameters ‘angle.grid’ and ‘dist.grid’. The distributions are represented by vectors giving the log-probability of each “bin”. These are not normalized by the “bin” widths. This is not a problem as these would cancel out later.

pysic.utility.outliers.get_log_likelihoods(angles, distances, angle_distribs, dist_distribs, n_atoms)[source]

Calculate log-likelihoods of each atom

All observed angles and distances are considered independent. Parameters ‘angles’ and .distances’ supplie the observations, ‘angle_distribs’ and ‘dist_distribs’’ the histogram-based distributions and ‘angle.grid’ & ‘dist.grid’ the breakpoints of these histograms. Parameter ‘n_atoms’ gives the number of atoms in the original data.

pysic.utility.outliers.vec_angle(A, B)[source]

Return the angle between vectors A and B

Parameters A, B: coordinates in 3d-space

pysic.utility.outliers.write_to_file(filename, boxsize, atoms, coordinates, a_logls, d_logls)[source]

Write original data + results into a file

Creates/overwrites an xyz-like file with two additional columns. These additional columns will contain the log-likelihood contributions from (1) angles and (2) distances.

MPI

The mpi utility defines a group of tools for parallel computations.

Defines MPI safe routines for printing, file access etc.

pysic.utility.mpi.append_file(lines, filename)[source]

Appends the given lines to a text file so that only the master cpu writes.

Parameters:

lines: list of strings
lines to be written
filename: string
the name of the file
pysic.utility.mpi.cd(dir)[source]

Changes to the given directory on all cpus.

Parameters:

dir: string
the name of the directory
pysic.utility.mpi.cpu_id()[source]

Gets the cpu ID from the Fortran MPI.

Equivalent to get_cpu_id().

pysic.utility.mpi.finish_mpi()[source]

Terminates the MPI framework.

If the Fortran core is compiled in MPI mode, pysic will automatically initialize MPI upon being imported. This method terminates the MPI.

pysic.utility.mpi.get_cpu_id()[source]

Gets the cpu ID from the Fortran MPI.

pysic.utility.mpi.get_number_of_cpus()[source]

Gets the number of cpus from the Fortran MPI.

pysic.utility.mpi.is_master()[source]

Returns True for the cpu with ID 0.

pysic.utility.mpi.mkdir(dir)[source]

Creates a new directory only with the master cpu.

Parameters:

dir: string
the name of the directory
pysic.utility.mpi.mpi_barrier()[source]

Calls MPI barrier from the Fortran core MPI framework.

This is equivalent to sync_mpi()

pysic.utility.mpi.mprint(string)[source]

Prints the string to stdout only from the master cpu.

Parameters:

string: string
the string to be written
pysic.utility.mpi.sync_mpi()[source]

Calls MPI barrier from the Fortran core MPI framework.

pysic.utility.mpi.write_file(lines, filename)[source]

Writes the given lines to a text file so that only the master cpu writes.

Parameters:

lines: list of strings
lines to be written
filename: string
the name of the file

Archive

The archive utility defines a group of functions for archiving and retrieving simulation data in the hdf5 format. This requires hdf5 and the h5py Python library.

class pysic.utility.archive.Archive(filename)[source]

A class representing a data archive in the hdf5 format.

The archive is built on h5py. See its documentation to access further functionality.

Parameters:

filename: string
path to the file storing the data
access_dataset(name)[source]
add_dataset_metadata(attribute, metadata)[source]
add_group_metadata(attribute, metadata)[source]
create_dataset(name, data, **kwds)[source]
create_subgroup(name)[source]
data()[source]
delete_current_dataset()[source]
delete_dataset(name)[source]
delete_group()[source]
delete_subgroup(name)[source]
get_contents()[source]
get_dataset()[source]
get_dataset_metadata()[source]
get_dataset_name()[source]
get_group()[source]
get_group_metadata()[source]
get_group_name()[source]
get_system(name)[source]

Restore a stored system.

The function returns a tuple with the data:

system, energy, forces, stress, electronegativities

Parameters:

name: string
name of the group storing the system
group()[source]
list()[source]
move(source, destination)[source]
move_to_group(name)[source]
move_to_parent_group()[source]
move_to_root_group()[source]
remove_dataset_metadata(attribute)[source]
remove_group_metadata(attribute)[source]
state()[source]
store_system(name, atoms, calculate=True)[source]

Stores the information of an atomic system.

The routine creates a group with the given name and stores the atomic information there as named datasets:

'atomic numbers'
'positions'
'charges'
'magnetic moments'
'tags'
'cell'
'pbc'
'momenta'
'potential energy'
'forces'
'stress'
'electronegativites'

The data requiring calculation is stored only if a calculator is attached to the atoms, electronegativities only if the calculator is Pysic.

Parameters:

name: string
the name of the group for storing the data
atoms: Atoms object
the structure to be stored
calculate: boolean
if True, the data requiring calculations will not be stored

Convenience

The convenience utility defines functions to ease the handling of complicated data structures.

pysic.utility.convenience.expand_symbols_string(symbol_string)

Expands a string of chemical symbols to list.

The function parses a string of chemical symbols and turns it to a list such as those expected by Potential.

Examples:

>>> print expand_symbols_string("HH")
[['H', 'H']]
>>> print expand_symbols_string("SiSi,SiO,SiH")
[['Si', 'Si'], ['Si', 'O'], ['Si', 'H']]

Parameters:

symbol_string: string
the string to be expanded
pysic.utility.convenience.expand_symbols_table(symbol_list, type=None)[source]

Creates a table of symbols object.

The syntax for defining the targets is precise but somewhat cumbersome due to the large number of permutations one gets when the number of bodies increases. Oftentimes one does not need such fine control over all the parameters since many of them have the same numerical values. Therefore it is convenient to be able to define the targets in a more compact way.

This method generates the detailed target tables from compact syntax. By default, the method takes a list of list and multiplies each list with the others (note the call for a static method):

>>> pysic.utility.convenience.expand_symbols_table([  'Si',
...                                                  ['O', 'C'],
...                                                  ['H', 'O'] ])
[['Si', 'O', 'H'],
 ['Si', 'C', 'H'],
 ['Si', 'O', 'O'],
 ['Si', 'C', 'O']]

Other custom types of formatting can be defined with the type parameter.

For type ‘triplet’, the target list is created for triplets A-B-C from an input list of the form:

['A', 'B', 'C']

Remember that in the symbol table accepted by the BondOrderParameters, one needs to define the B-A and B-C bonds separately and so B appears as the first symbol in the output and the other two appear as second and third (both cases):

[['B', 'A', 'C'],
 ['B', 'C', 'A']]

However, for an A-B-A triplet, the A-B bond should only be defined once to prevent double counting. Like the default function, also here several triplets can be defined at once:

>>> pysic.BondOrderParameters.expand_symbols_table([ ['H', 'O'],
...                                                   'Si',
...                                                  ['O', 'C'] ],
...                                                type='triplet')
[['Si', 'H', 'O'],
 ['Si', 'O', 'H'],
 ['Si', 'H', 'C'],
 ['Si', 'C', 'H'],
 ['Si', 'O', 'O'],
 ['Si', 'O', 'C'],
 ['Si', 'C', 'O']]

Parameters:

symbol_list: list of strings list to be expanded to a table type: string specifies a custom way of generating the table

Geometry

The geometry utility defines tools for geometric operations.

class pysic.utility.geometry.Cell(vector1, vector2, vector3, pbc)[source]

Cell describing the simulation volume of a subvolume.

This class can be used by the user for coordinate manipulation. Note however, that ASE does not use on this class, as it is part of Pysic. The class is merely a tool for examining the geometry.

Parameters:

vector1: list of doubles
3-vector specifying the first vector spanning the cell \(\mathbf{v}_1\)
vector2: list of doubles
3-vector specifying the second vector spanning the cell \(\mathbf{v}_2\)
vector3: list of doubles
3-vector specifying the third vector spanning the cell \(\mathbf{v}_3\)
pbc: list of logicals
three logic switches for specifying periodic boundaries - True denotes periodicity
get_absolute_coordinates(fractional)[source]

For the given fractional coordinates, returns the absolute coordinates.

The absolute coordinates are the cell vectors multiplied by the fractional coordinates.

Parameters:

fractional: numpy double 3-vector
the fractional coordinates
get_distance(atom1, atom2, offsets=None)[source]

Calculates the distance between two atoms.

Offsets are multipliers for the cell vectors to be added to the plain separation vector r1-r2 between the atoms.

Parameters:

atom1: ASE Atoms object
first atom
atom2: ASE Atoms object
second atom
offsets: Numpy integer 3-vector
the periodic boundary offsets
get_relative_coordinates(coordinates)[source]

Returns the coordinates of the given atom in fractional coordinates.

The absolute position of the atom is given by multiplying the cell vectors by the fractional coordinates.

Parameters:

coordinates: numpy double 3-vector
the absolute coordinates
get_separation(atom1, atom2, offsets=None)[source]

Returns the separation vector between two atoms, r1-r2.

Offsets are multipliers for the cell vectors to be added to the plain separation vector r1-r2 between the atoms.

Parameters:

atom1: ASE Atoms object
first atom
atom2: ASE Atoms object
second atom
offsets: Numpy integer 3-vector
the periodic boundary offsets
get_wrapped_coordinates(coordinates)[source]

Wraps the coordinates of the given atom inside the simulation cell.

This method return the equivalent position (with respect to the periodic boundaries) of the atom inside the cell.

For instance, if the cell is spanned by vectors of length 10.0 in directions of \(x\), \(y\), and \(z\), an the coordinates [-1.0, 12.0, 3.0] wrap to [9.0, 2.0, 3.0].

Parameters:

coordinates: numpy double 3-vector
the absolute coordinates

Error

The error utility defines a group of intrinsic errors to describe situations where one tries to use or set up the calculator with errorneous or insufficient information.

The module also defines the Warning class and related routines for displaying warnings for the user upon suspicious but non-critical behavior.

exception pysic.utility.error.InvalidCoordinatorError(message='', coordinator=None)[source]

An error raised when an invalid coordination calculator is about to be created or used.

Parameters:

message: string
information describing why the error occurred
coordinator: Coordinator
the errorneous coordinator
exception pysic.utility.error.InvalidParametersError(message='', params=None)[source]

An error raised when an invalid set of parameters is about to be created.

Parameters:

message: string
information describing why the error occurred
params: BondOrderParameters
the errorneous parameters
exception pysic.utility.error.InvalidPotentialError(message='', potential=None)[source]

An error raised when an invalid potential is about to be created or used.

Parameters:

message: string
information describing why the error occurred
potential: Potential
the errorneous potential
exception pysic.utility.error.InvalidRelaxationError(message='', relaxation=None)[source]

An error raised when an invalid charge relaxation is about to be created.

Parameters:

message: string
information describing why the error occurred
params: ChargeRelaxation
the errorneous parameters
exception pysic.utility.error.InvalidSummationError(message='', summer=None)[source]

An error raised when an invalid coulomb summation is about to be created.

Parameters:

message: string
information describing why the error occurred
params: CoulombSummation
the errorneous summation
exception pysic.utility.error.LockedCoreError(message='')[source]

An error raised when a Pysic tries to access the core which is locked by another calculator.

Parameters:

message: string
information describing why the error occurred
exception pysic.utility.error.MissingAtomsError(message='')[source]

An error raised when the core is being updated with per atom information before updating the atoms.

Parameters:

message: string
information describing why the error occurred
exception pysic.utility.error.MissingNeighborsError(message='')[source]

An error raised when a calculation is initiated without initializing the neighbor lists.

In principle Pysic should always take care of handling the neighbors automatically. This error is an indication that there is loophole in the built-in preparations.

Parameters:

message: string
information describing why the error occurred
class pysic.utility.error.Warning(message, level)

A warning raised due to a potentially dangerous action.

The warning is not an exception, so it doesn’t by default interrupt execution. It will, however, display a message or perform an action depending on the warning settings.

The warning levels are:

1: a condition that leads to unwanted behaviour (e.g., creating redundant potential in core) 2: a condition that is likely unwanted, but possibly a hack 3: a condition that is likely unwanted, but is a featured hack (e.g., bond order mixing) 4: a condition that is harmless but may lead to errors later (e.g., defining a potential without targets - often you’ll specify the targets later) 5: a condition that may call for attention (notes to user)

Parameters:

message: string
information describing the cause for the warning
level: int
severity of the warning (1-5, 1 being most severe)
display()

Displays the message related to this warning, but only if the level of the warning is smaller (more severe) than the global warning_level.

exception pysic.utility.error.WarningInterruptException(message)

An error raised automatically for any warning when strict warnings are in use (warning_level = 6).

pysic.utility.error.error(message)

Raises an error and displays the reason. Independent of the warning level.

Parameters:

message: string
information describing the cause for the warning
pysic.utility.error.set_warning_level(level)

Set the warning level.

Parameters:

level: int
The level of warnings displayed. Should be an integer between 0 (no warnings) and 6 (warnings are treated as errors).
pysic.utility.error.style_message(header, message, width=80, x_pad=2, y_pad=1)

Styles a message to be printed into console.

The message can be a list of string, where each list item corresponds to a row. If a single string is provided, it is converted to a list. Rows are cut into pieces so that they will fit into the defined width.

pysic.utility.error.warn(message, level)

Raise and display a warning.

Parameters:

message: string
information describing the cause for the warning
level: int
severity of the warning (1-5, 1 being most severe)

Debug

The debug utility defines debugging tools.

F2py

The f2py utility defines tools for the Python-Fortran interfacing.

pysic.utility.f2py.char2int(char_in)[source]

Codes a single character to an integer.

pysic.utility.f2py.int2char(int_in)[source]

Decodes an integer to a single character.

pysic.utility.f2py.ints2str(ints_in)[source]

Decodes a list of integers to a string.

pysic.utility.f2py.str2ints(string_in, target_length=0)[source]

Codes a string to a list of integers.

Turns a string to a list of integers for f2py interfacing. If required, the length of the list can be specified and trailing spaces will be added to the end.