Python binding: pyexadis
¶
This section documents the raw ExaDiS classes and functions available through the python interface, pyexadis
.
These modules are direct bindings to the backend C++ modules implemented in ExaDiS.
For documentation of the ExaDiS python modules available in OpenDiS format, see Python modules: pyexadis_base.
For documentation about the backend C++ modules please see the Developer guide section of the documentation.
General Attributes¶
- pyexadis.BCC_CRYSTAL: int
Constant for BCC crystal type.
- pyexadis.FCC_CRYSTAL: int
Constant for FCC crystal type.
General Classes¶
- class pyexadis.Params¶
Simulation parameters.
- __init__(crystal="", burgmag, mu, nu, a, maxseg, minseg, rann=-1.0, rtol=-1.0, maxdt=1e-7, nextdt=1e-12, split3node=1)
Initialize parameters.
- Parameters:
crystal – Crystal type (string)
burgmag – Burgers vector magnitude
mu – Shear modulus
nu – Poisson’s ratio
a – Dislocation core radius
maxseg – Maximum line discretization length
minseg – Minimum line discretization length
rann – Annihilation distance
rtol – Error tolerance
maxdt – Maximum timestep size
nextdt – Starting timestep size
split3node – Enable splitting of 3-nodes
- crystalparams: CrystalParams
Crystal parameters.
- burgmag: float
Burgers vector magnitude (scaling length).
- mu: float
Shear modulus.
- nu: float
Poisson’s ratio.
- a: float
Dislocation core radius.
- maxseg: float
Maximum line discretization length.
- minseg: float
Minimum line discretization length.
- rann: float
Annihilation distance.
- rtol: float
Error tolerance.
- maxdt: float
Maximum timestep size.
- nextdt: float
Starting timestep size.
- split3node: int
Enable splitting of 3-nodes.
—
- class pyexadis.CrystalParams¶
Parameters for crystal type, orientation, and glide planes.
- set_crystal_type(str)
Set the crystal type.
- R: array
Crystal orientation matrix.
- use_glide_planes: bool
Use and maintain dislocation glide planes.
- enforce_glide_planes: bool
Enforce glide planes option.
- num_bcc_plane_families: int
Number of BCC plane families (1, 2, or 3).
—
- class pyexadis.Crystal¶
Crystal type and orientation.
- __init__(type)
Initialize with crystal type.
- __init__(type, R)
Initialize with crystal type and orientation matrix.
- __init__(crystalparams)
Initialize with crystal parameters.
- type: int
Index of the crystal type.
- R: array
Crystal orientation matrix.
- set_orientation(R)
Set crystal orientation matrix.
- set_orientation(euler_angles)
Set crystal orientation via Euler angles.
—
- class pyexadis.Cell¶
Simulation cell and periodic boundary conditions.
- __init__(Lbox, centered=False)
Initialize cubic cell.
- __init__(Lvecbox, centered=False)
Initialize cell with vector box.
- __init__(bmin, bmax)
Initialize cell with bounds.
- __init__(h, origin=Vec3(0.0), is_periodic=[PBC_BOUND, PBC_BOUND, PBC_BOUND])
Initialize cell with matrix, origin, and periodicity.
- __init__(cell)
Copy constructor.
- h: array
Cell matrix.
- origin: array
Origin of the cell.
- center()
Returns the center of the cell.
- closest_image(Rref, R)
Returns the closest image of an array of positions from another reference.
- pbc_fold(R)
Fold an array of positions to the primary cell.
- is_inside(R)
Checks if a position is inside the primary cell.
- are_inside(R)
Checks if an array of positions are inside the primary cell.
- is_triclinic()
Returns if the box is triclinic.
- is_periodic()
Get the cell periodic boundary condition flags along the 3 dimensions.
- get_bounds()
Get the (orthorhombic) bounds of the cell.
- volume()
Returns the volume of the cell.
—
- class pyexadis.NodeTag¶
Tag identifying a node by domain and local index.
- __init__(domain, index)
- Parameters:
domain – Domain index (int)
index – Local index (int)
- domain: int
Domain index.
- index: int
Local index.
—
- class pyexadis.DisNode¶
Dislocation node, including position, constraint, and tag.
- __init__(pos, constraint)
- Parameters:
pos – Node position (Vec3)
constraint – Node constraint flag (int)
- tag: NodeTag
Node tag (domain, index).
- constraint: int
Node constraint flag.
- pos: Vec3
Node position (x, y, z).
—
- class pyexadis.DisSeg¶
Dislocation segment between two nodes.
- __init__(n1, n2, burg, plane=Vec3(0.0))
- Parameters:
n1 – Segment start node index (int)
n2 – Segment end node index (int)
burg – Segment Burgers vector (Vec3)
plane – Segment plane normal (Vec3, default zero vector)
- n1: int
Segment start node index.
- n2: int
Segment end node index.
- burg: Vec3
Segment Burgers vector.
- plane: Vec3
Segment plane normal.
—
- class pyexadis.Conn¶
Connectivity information for a node.
- __init__()
Default constructor.
- num: int
Number of connections.
- node(i)
Get the node index of the i-th connection.
- Parameters:
i – Connection index (int)
- Return type:
int
- seg(i)
Get the segment index of the i-th connection.
- Parameters:
i – Connection index (int)
- Return type:
int
- order(i)
Get the order value of the i-th connection.
- Parameters:
i – Connection index (int)
- Return type:
int
- add_connection(node, seg, order)
Add a connection.
- Parameters:
node – Node index (int)
seg – Segment index (int)
order – Order value (int)
- Return type:
bool
- remove_connection(i)
Remove the i-th connection.
- Parameters:
i – Connection index (int)
—
- class pyexadis.SerialDisNet¶
Serial dislocation network object. This is the direct binding to the C++ implementation. In the python interface, a wrapper ExaDisNet object is usually used instead as a proxy for a SerialDisNet object.
- __init__(cell)
- Parameters:
cell – Simulation cell (Cell)
- number_of_nodes()
Get the number of nodes.
- Return type:
int
- number_of_segs()
Get the number of segments.
- Return type:
int
- dislocation_density(burgmag)
Get the dislocation density.
- Return type:
float
- cell: Cell
Simulation cell.
- nodes(i)
Get the i-th node by reference.
- Parameters:
i – Node index (int)
- Return type:
- segs(i)
Get the i-th segment by reference.
- Parameters:
i – Segment index (int)
- Return type:
- conn(i)
Get the connectivity object for the i-th node by reference.
- Parameters:
i – Node index (int)
- Return type:
- find_connection(n1, n2)
Find connection index c in the Conn array of node 1 that connects to node 2, i.e. where Conn(n1).node(c) == n2.
- Parameters:
n1 – Node 1 index (int)
n2 – Node 2 index (int)
- Return type:
int
- generate_connectivity()
Generate internal connectivity array for the network.
- _add_node(pos, constraint=UNCONSTRAINED)
Add node (x, y, z[, constraint]) to the network.
- Parameters:
pos – Position (Vec3)
constraint – Constraint flag (int, default UNCONSTRAINED)
- _add_seg(n1, n2, burg, plane=Vec3(0.0))
Add segment (n1, n2, burg[, plane]) to the network.
- Parameters:
n1 – Start node index (int)
n2 – End node index (int)
burg – Burgers vector (Vec3)
plane – Plane normal (Vec3, default zero vector)
- move_node(node_index, new_pos, dEp)
Move a node to a new position.
- Parameters:
node_index – Index of node (int)
new_pos – New position (Vec3)
dEp – Reference to the plastic strain variable to be incremented (Mat33)
- split_seg(seg_index, new_pos)
Split a segment by inserting a new node at new_pos.
- Parameters:
seg_index – Segment index (int)
new_pos – New node position (Vec3)
- Returns:
Index of the new node (int)
- split_node(node_index, arms)
Split a node given the list of arms to be transferred to the new node.
- Parameters:
node_index – Node index (int)
arms – Node arms (array)
- Returns:
Index of the new node (int)
- merge_nodes(node_index1, node_index2, dEp)
Merge two nodes into the first node.
- Parameters:
node_index1 – First node index (int)
node_index2 – Second node index (int)
dEp – Reference to the plastic strain variable to be incremented (Mat33)
- Returns:
Merge error (bool)
- merge_nodes_position(node_index1, node_index2, new_pos, dEp)
Merge two nodes into the first node at the new position.
- Parameters:
node_index1 – First node index (int)
node_index2 – Second node index (int)
new_pos – New position (Vec3)
dEp – Reference to the plastic strain variable to be incremented (Mat33)
- Returns:
Merge error (bool)
- remove_segs(seg_indices)
Remove segments by indices.
- Parameters:
seg_indices – List of segment indices
- remove_nodes(node_indices)
Remove nodes by indices.
- Parameters:
node_indices – List of node indices
- purge_network()
Purge the network (remove unused nodes/segs).
- update()
Update network memory after modifications.
—
- class pyexadis.ExaDisNet¶
Dislocation network object, supporting import/export, manipulation, and analysis of network data.
- __init__()
Default constructor.
- __init__(cell, nodes, segs)
Construct network from cell, nodes, and segments.
- Parameters:
cell – Cell object
nodes – List of node data (list of lists: [dom, id, x, y, z, constraint])
segs – List of segment data (list of lists: [n1, n2, burg, plane])
- import_data(cell, nodes, segs)
Set the network with (cell, nodes, segs) data.
- Parameters:
cell – Cell object
nodes – List of node data
segs – List of segment data
- number_of_nodes()
Returns the number of nodes in the network.
- Return type:
int
- number_of_segs()
Returns the number of segments in the network.
- Return type:
int
- is_sane()
Checks if the network is sane (valid connectivity, no corrupt data).
- Return type:
bool
- get_cell()
Get the cell containing the network by value.
- Return type:
- get_nodes_array()
Get the list of nodes in the network.
- Returns:
List of [dom, id, x, y, z, constraint] for each node
- get_segs_array()
Get the list of segments in the network.
- Returns:
List of [n1, n2, burg, plane] for each segment
- get_forces()
Get the list of node forces (fx, fy, fz) in the network.
- Returns:
List of [fx, fy, fz] for each node
- get_velocities()
Get the list of node velocities (vx, vy, vz) in the network.
- Returns:
List of [vx, vy, vz] for each node
- set_positions(positions)
Set the list of node positions (x, y, z) in the network.
- Parameters:
positions – List of [x, y, z] for each node
- set_forces(forces)
Set the list of node forces (fx, fy, fz) in the network.
- Parameters:
forces – List of [fx, fy, fz] for each node
- set_velocities(velocities)
Set the list of node velocities (vx, vy, vz) in the network.
- Parameters:
velocities – List of [vx, vy, vz] for each node
- write_data(filename)
Write network data in ParaDiS format.
- Parameters:
filename – Output file path (str)
- get_plastic_strain()
Returns the plastic strain, plastic spin, and dislocation density, as computed since the last integration step.
- Returns:
Plastic strain array, plastic spin array, dislocation density
- physical_links()
Returns the list of segments for each physical dislocation link.
- Returns:
List of segment indices for each link
- _get_crystal()
Get the Crystal object by reference.
- Return type:
- _get_serial_network()
Get the underlying SerialDisNet object by reference.
- Return type:
—
- class pyexadis.System(pyexadis.ExaDisNet)¶
ExaDiS simulation system, combining a network and simulation parameters.
- __init__(net, params)
Construct a system from a network and parameters.
- Parameters:
net – ExaDisNet object
params – Params object
- set_neighbor_cutoff(cutoff)
Set neighbor cutoff of the system.
- Parameters:
cutoff – Cutoff value (float)
- set_applied_stress(stress)
Set applied stress of the system.
- Parameters:
stress – List or array of [xx, yy, zz, yz, xz, xy]
- print_timers(dev=False)
Print simulation timers.
- Parameters:
dev – If True, print device timers (bool, default False)
—
General Functions¶
- pyexadis.initialize(num_threads=-1, device_id=0)¶
Initialize the python binding module.
- Parameters:
num_threads – Number of OMP threads (default: -1)
device_id – Device GPU index (default: 0)
- pyexadis.finalize()¶
Finalize the python binding module.
- pyexadis.read_paradis(filename)¶
Read ParaDiS data file.
- Parameters:
filename – Path to ParaDiS data file
- Return type:
- pyexadis.generate_prismatic_config(crystal, Lbox, numsources, radius, maxseg=-1, seed=1234, uniform=False)
Generate a configuration made of prismatic loops in a cubic cell.
- Parameters:
crystal – Crystal object
Lbox – Box size
numsources – Number of sources
radius – Loop radius
maxseg – Maximum segment length (default: -1)
seed – Random seed (default: 1234)
uniform – Uniform distribution (default: False)
- Return type:
- pyexadis.generate_prismatic_config(crystal, cell, numsources, radius, maxseg=-1, seed=1234, uniform=False)¶
Generate a configuration made of prismatic loops in a custom cell.
- Parameters:
crystal – Crystal object
cell – Cell object
numsources – Number of sources
radius – Loop radius
maxseg – Maximum segment length (default: -1)
seed – Random seed (default: 1234)
uniform – Uniform distribution (default: False)
- Return type:
Force¶
Force Modules and Factories¶
- class pyexadis.Force.Force¶
ExaDiS base force class.
- name()
Get force name.
—
- class pyexadis.Force.CORE_SELF_PKEXT(pyexadis.Force.Force)¶
Core, self, and external PK force model.
- class Params
Parameters for CORE_SELF_PKEXT model.
- __init__(Ecore=-1.0, Ecore_junc_fact=1.0)
- Parameters:
Ecore – Core energy (float, default -1.0)
Ecore_junc_fact – Junction energy factor (float, default 1.0)
—
- class pyexadis.Force.ForceFFT(pyexadis.Force.Force)¶
FFT-based force model.
- class Params
Parameters for ForceFFT model.
- __init__(Ngrid)
- Parameters:
Ngrid – Grid dimensions (list of 3 ints)
- make(params, fparams, cell)
Factory to create a ForceFFT model.
- Parameters:
params – General simulation parameters
fparams – ForceFFT.Params object
cell – Cell object
- Return type:
- get_neighbor_cutoff()
Returns the neighbor cutoff used in ForceFFT.
- get_rcgrid()
Returns the grid cutoff used in ForceFFT.
- interpolate_stress(R)
Interpolates stress at query positions from FFT grid values.
- Parameters:
R – Array of positions
- export_stress_gridval()
Exports grid stress values.
—
- class pyexadis.Force.LINE_TENSION_MODEL(pyexadis.Force.Force)¶
Line tension force model.
- make(params, coreparams)
Factory to create a LINE_TENSION_MODEL.
- Parameters:
params – Global simulation parameters
coreparams – CORE_SELF_PKEXT.Params object
- Return type:
—
- class pyexadis.Force.CUTOFF_MODEL(pyexadis.Force.Force)¶
Cutoff-based force model.
- class Params
Parameters for CUTOFF_MODEL.
- __init__(coreparams, cutoff)
- Parameters:
coreparams – CORE_SELF_PKEXT.Params object
cutoff – Cutoff distance (float)
- make(params, fparams)
Factory to create a CUTOFF_MODEL.
- Parameters:
params – Global simulation parameters
fparams – CUTOFF_MODEL.Params object
- Return type:
—
- class pyexadis.Force.DDD_FFT_MODEL(pyexadis.Force.Force)¶
DDD-FFT force model.
- class Params
Parameters for DDD_FFT_MODEL.
- __init__(coreparams, Ngrid)
- Parameters:
coreparams – CORE_SELF_PKEXT.Params object
Ngrid – Grid dimensions (list of 3 ints)
- make(params, fparams, cell)
Factory to create a DDD_FFT_MODEL.
- Parameters:
params – Global simulation parameters
fparams – DDD_FFT_MODEL.Params object
cell – Cell object
- Return type:
—
- class pyexadis.Force.SUBCYCLING_MODEL(pyexadis.Force.Force)¶
Subcycling force model.
- class Params
Parameters for SUBCYCLING_MODEL.
- __init__(coreparams, Ngrid, drift=False, flong_group0=True)
- Parameters:
coreparams – CORE_SELF_PKEXT.Params object
Ngrid – Grid dimensions (list of 3 ints)
drift – Enable drift subcycling scheme (bool, default False)
flong_group0 – Integrate long-range force in group 0 (bool, default True)
- make(params, fparams, cell)
Factory to create a SUBCYCLING_MODEL.
- Parameters:
params – Global simulation parameters
fparams – SUBCYCLING_MODEL.Params object
cell – Cell object
- Return type:
—
- class pyexadis.Force.ForcePython(pyexadis.Force.Force)¶
Python-based force model.
- make(params, force)
Factory to create a Python-based force model.
- Parameters:
params – Global simulation parameters
force – CalForce python object implementing force methods
- Return type:
—
Force Binder¶
- class pyexadis.ForceBind¶
Wrapper for force computation and related operations.
- force: Force
Internal force object.
- neighbor_cutoff: float
Neighbor cutoff distance.
- _pre_compute(system)
Pre-compute force of the system (internal use).
- _compute(system)
Compute force of the system (internal use).
- _get_force_fft()
Return internal ForceFFT object if it exists in the internal force model.
- Return type:
- compute_force(net, applied_stress=[0.0, 0.0, 0.0, 0.0, 0.0, 0.0], pre_compute=True)
Wrapper to compute nodal forces of the system.
- Parameters:
net – ExaDisNet object
applied_stress – Applied stress tensor (list of 6 floats)
pre_compute – Whether to run pre-computation (bool)
- pre_compute_force(net)
Wrapper to perform pre-computations before compute_node_force().
- Parameters:
net – ExaDisNet object
- compute_node_force(net, i, applied_stress)
Wrapper to compute the force on a single node.
- Parameters:
net – ExaDisNet object
i – Node index (int)
applied_stress – Applied stress tensor (list of 6 floats)
—
Force Calculation Functions¶
- pyexadis.compute_force_n2(net, mu, nu, a)¶
Compute elastic forces using brute-force N^2 calculation.
- Parameters:
net – ExaDisNet object
mu – Shear modulus (float)
nu – Poisson ratio (float)
a – Core cutoff (float)
- pyexadis.compute_force_cutoff(net, mu, nu, a, cutoff, maxseg=0.0)¶
Compute elastic forces using a segment pair cutoff.
- Parameters:
net – ExaDisNet object
mu – Shear modulus (float)
nu – Poisson ratio (float)
a – Core cutoff (float)
cutoff – Cutoff distance (float)
maxseg – Max segment length (float, default 0.0)
- pyexadis.compute_force_segseglist(net, mu, nu, a, segseglist)¶
Compute elastic forces given a list of segment pairs.
- Parameters:
net – ExaDisNet object
mu – Shear modulus (float)
nu – Poisson ratio (float)
a – Core cutoff (float)
segseglist – List of segment pairs
—
Mobility¶
Mobility Binder¶
- class pyexadis.Mobility¶
Mobility law wrapper.
- compute(system)
Compute mobility of the system.
—
Mobility Law Factories¶
- pyexadis.make_mobility_glide(params, mobparams)¶
Instantiate a GLIDE mobility law.
- pyexadis.make_mobility_bcc_0b(params, mobparams)¶
Instantiate a BCC_0B mobility law.
- pyexadis.make_mobility_bcc_nl(params, mobparams)¶
Instantiate a BCC_NL mobility law.
- pyexadis.make_mobility_fcc_0(params, mobparams)¶
Instantiate a FCC_0 mobility law.
- pyexadis.make_mobility_fcc_0_fric(params, mobparams)¶
Instantiate a FCC_0_FRIC mobility law.
- pyexadis.make_mobility_fcc_0b(params, mobparams)¶
Instantiate a FCC_0B mobility law.
- pyexadis.make_mobility_python(params, mobility)¶
Instantiate a python-based mobility model.
—
Mobility Calculation Functions¶
- pyexadis.compute_mobility(net, mobility, nodeforces, nodetags=[])¶
Wrapper to compute nodal velocities.
- Parameters:
net – ExaDisNet object
mobility – Mobility object
nodeforces – List of node forces
nodetags – List of NodeTag objects (optional)
- pyexadis.compute_node_mobility(net, i, mobility, fi)¶
Wrapper to compute the mobility of a single node.
- Parameters:
net – ExaDisNet object
i – Node index (int)
mobility – Mobility object
fi – Force on node
—
Integrator¶
Integrator Binder¶
- class pyexadis.Integrator¶
Integrator wrapper.
- integrate(system)
Integrate the system.
—
Integrator Factories¶
- pyexadis.make_integrator_trapezoid(params, intparams, force, mobility)¶
Instantiate a trapezoid integrator.
- pyexadis.make_integrator_trapezoid_multi(params, intparams, force, mobility)¶
Instantiate a multi-step trapezoid integrator.
- pyexadis.make_integrator_rkf(params, intparams, force, mobility)¶
Instantiate a RKF integrator.
- pyexadis.make_integrator_rkf_multi(params, intparams, force, mobility)¶
Instantiate a multi-step RKF integrator.
- pyexadis.make_integrator_subcycling(params, intparams, force, mobility)¶
Instantiate a subcycling integrator.
—
Integration Functions¶
- pyexadis.integrate(net, integrator, nodevels, applied_stress, nodetags=[])¶
Wrapper to perform a time-integration step.
- pyexadis.integrate_euler(net, params, dt, nodevels, nodetags=[])¶
Time-integrate positions using the Euler integrator.
—
Collision¶
- class pyexadis.Collision¶
Collision handler.
- handle(system)
Handle collision of the system.
—
- pyexadis.make_collision(collision_mode, params)¶
Instantiate a collision class.
- pyexadis.handle_collision(net, collision, xold=[], dt=0.0)¶
Wrapper to handle collisions.
—
Topology¶
- class pyexadis.Topology¶
Topology handler.
- handle(system)
Handle topology of the system.
—
- pyexadis.make_topology(topology_mode, params, topolparams, force, mobility)¶
Instantiate a topology class.
- pyexadis.handle_topology(net, topology, dt)¶
Wrapper to handle topological operations.
—
Remesh¶
- class pyexadis.Remesh¶
Remesh handler.
- remesh(system)
Remesh the system.
—
- pyexadis.make_remesh(remesh_rule, params, remeshparams)¶
Instantiate a remesh class.
- pyexadis.remesh(net, remesh)¶
Wrapper to remesh the network.
—
Cross-slip¶
- class pyexadis.CrossSlip¶
Cross-slip handler.
- handle(system)
Handle cross-slip operations of the system.
—
- pyexadis.make_cross_slip(cross_slip_mode, params, force)¶
Instantiate a cross-slip class.
- pyexadis.handle_cross_slip(net, cross_slip)¶
Wrapper to handle cross-slip operations.
—
Simulation Driver¶
Drivers¶
- class pyexadis.ExaDiSApp¶
Base application class for pyexadis simulation driver.
—
- class pyexadis.Driver¶
Main simulation driver, manages simulation state, modules, and execution.
- __init__()
Default constructor.
- __init__(system)
Construct driver with a system.
- Parameters:
system – SystemBind object
- outputdir: str
Output directory path for the simulation.
- update_state(state)
Update the state dictionary with simulation state.
- read_restart(state, restart_file)
Read restart file.
- set_system(system)
Set system for the simulation.
- Parameters:
system – SystemBind object
- set_modules(force, mobility, integrator, collision, topology, remesh, cross_slip=CrossSlipBind())
Set modules for the simulation.
- Parameters:
force – ForceBind object
mobility – MobilityBind object
integrator – IntegratorBind object
collision – CollisionBind object
topology – TopologyBind object
remesh – RemeshBind object
cross_slip – CrossSlipBind object (optional)
- set_simulation(restart='')
Set things up before running the simulation.
- Parameters:
restart – Restart file path (str, optional)
- initialize(ctrl, check_modules=True)
Initialize simulation.
- Parameters:
ctrl – Driver.Control object
check_modules – Check modules before initialization (bool, default True)
- step(ctrl)
Execute a simulation step.
- run(ctrl)
Run the simulation.
- oprec_replay(ctrl, oprec_files)
Replay the simulation from OpRec files.
—
Driver Control¶
- class pyexadis.Driver.Control¶
Control object for simulation parameters.
- nsteps: int or Stepper
Number of steps or stepper object.
- loading: Loadings
Loading type.
- erate: float
Loading rate.
- edir: list or array
Loading direction.
- appstress: list or array
Applied stress.
- rotation: bool
Enable crystal rotation.
- printfreq: int
Print frequency.
- propfreq: int
Properties output frequency.
- outfreq: int
Configuration and restart output frequency.
- outfreqdt: float
Configuration and restart output time frequency.
- oprecwritefreq: int
OpRec write frequency.
- oprecfilefreq: int
OpRec new file frequency.
- oprecposfreq: int
OpRec nodal positions save frequency.
- set_props(prop_fields)
Set property fields for the output.
—
Driver Loadings Enum¶
- class pyexadis.Driver.Loadings¶
Enum for loading types.
- STRESS_CONTROL
Stress control loading.
- STRAIN_RATE_CONTROL
Strain rate control loading.
—
Driver Stepper¶
- class pyexadis.Driver.Stepper¶
Stepper object for simulation steps.
- __init__(type, stopval)
Construct stepper with stopping type and value.
- Parameters:
type – Stopping type
stopval – Stopping value
- NUM_STEPS(nsteps)
Construct stepper that iterates to a number of steps.
- MAX_STEPS(maxsteps)
Construct stepper that iterates to a maximum number of steps.
- MAX_STRAIN(strain)
Construct stepper that iterates to a maximum strain.
- MAX_TIME(time)
Construct stepper that iterates to a maximum simulation time.
- MAX_WALLTIME(time)
Construct stepper that iterates to a maximum wall clock time.
- iterate(driver)
Iterate a simulation step.
—