Functions

Index

Data_manager

PeriLab.Data_manager.check_propertyFunction
check_property(block_id::Int64, property::String)

Checks if the specified property exists for the given block_id.

Arguments

  • block_id::Int64: The ID of the block.
  • property::String: The name of the property to check.

Returns

  • Bool: true if the property exists, false otherwise.
source
PeriLab.Data_manager.create_bond_fieldFunction
create_bond_field(name::String, type::Type, dof::Int64)

Creates a bond field with the given name, data type, and degree of freedom.

Arguments

  • name::String: The name of the bond field.
  • vartype::Type: The data type of the bond field.
  • dof::Int64: The degrees of freedom per bond.
  • VectorOrArray::String (optional) - Vector or Materix; Default is vector

Returns

  • bond_field::Field: The created bond field for the current time step.
  • bond_field_np1::Field: The created bond field for the next time step.

Example:

create_bond_field("stress", Float64, 6)  # creates a stress bond field with 6 degrees of freedom
source
PeriLab.Data_manager.create_constant_bond_fieldFunction
create_constant_bond_field(name::String, type::Type, dof::Int64, default_value::Union{Int64,Float64,Bool}=0))

Creates a constant bond field with the given name, data type, and degree of freedom.

Arguments

  • name::String: The name of the bond field.
  • vartype::Type: The data type of the bond field.
  • dof::Int64: The degrees of freedom per bond.
  • default_value::Union{Int64,Float64,Bool}=0) (optional) - filled with zero or false

Returns

  • constant_bond_field::Field: The created constant bond field.

Example:

create_constant_bond_field("density", Float64, 1)  # creates a density constant bond field
source
PeriLab.Data_manager.create_constant_node_fieldFunction
create_constant_node_field(name::String, type::Type, dof::Int64)

Creates a constant node field with the given name, data type, and degree of freedom.

Arguments

  • name::String: The name of the node field.
  • vartype::Type: The data type of the node field.
  • dof::Int64: The degrees of freedom per node.
  • VectorOrArray::String (optional) - Vector or Materix; Default is vector

Returns

  • constant_node_field::Field: The created constant node field.

Example:

create_constant_node_field("temperature", Float64, 1)  # creates a temperature constant node field
source
PeriLab.Data_manager.create_fieldFunction
create_field(name::String, vartype::Type, bondNode::String, dof::Int64, default_value::Any=0)

Create a field with the given name for the specified vartype. If the field already exists, return the existing field. If the field does not exist, create a new field with the specified characteristics.

Arguments

  • name::String: The name of the field.
  • vartype::Type: The data type of the field.
  • dof::Int64: The degrees of freedom per node.
  • default_value::Any: The default value of the field.

Returns

The field with the given name and specified characteristics.

source
PeriLab.Data_manager.create_node_fieldFunction
create_node_field(name::String, type::Type, dof::Int64)

Creates a node field with the given name, data type, and degree of freedom.

Arguments

  • name::String: The name of the node field.
  • type::Type: The data type of the node field.
  • dof::Int64: The degree of freedom of each node.
  • VectorOrArray::String (optional) - Vector or Materix; Default is vector

Returns

  • node_field::Field: The created node field for the current time step.
  • node_field_np1::Field: The created node field for the next time step.

Example:

create_node_field("displacement", Float64, 3)  # creates a displacement node field with 3 degrees of freedom
source
PeriLab.Data_manager.get_dofFunction
get_dof()

Retrieves the degree of freedom (dof) value.

Returns

  • dof (integer): The current degree of freedom value.

Example:

get_dof()  # returns the current degree of freedom
source
PeriLab.Data_manager.get_fieldFunction
get_field(name::String, time::String)

Returns the field with the given name and time.

Arguments

  • name::String: The name of the field.
  • time::String: The time of the field.

Returns

  • field::Field: The field with the given name and time.
source
get_field(name::String, throw_error::Bool=true)

Returns the field with the given name.

Arguments

  • name::String: The name of the field.
  • throw_error::Bool=true: Whether to throw an error if the field does not exist.

Returns

  • field::Field: The field with the given name.
source
PeriLab.Data_manager.get_local_nodesFunction
get_local_nodes()

Determines the local node numbering.

Returns

  • get_local_nodes (array): returns local nodes.

Example:

get_local_nodes()  # returns local nodes or if they do not exist at the core an empty array
source
PeriLab.Data_manager.get_nnodesFunction
get_nnodes()

Retrieves the number of nodes.

Returns

  • num_controller::Int64 : The current number of nodes.

Example:

get_nnodes()  # returns the current number of controler nodes. The neighbors are not included
source
PeriLab.Data_manager.get_propertiesFunction
get_properties(block_id::Int64, property::String)

This function retrieves the value of a specified property for a given block_id if it exists in the properties dictionary.

Arguments

  • block_id::Int64: The identifier of the block for which to retrieve the property.
  • property::String: The dictionary entrycontaining the properties for the blocks.

Returns

  • property_value::Any: The value associated with the specified property for the given block_id.
  • Dict(): An empty dictionary if the specified property does not exist for the given block_id.

Example

```julia block_properties = Dict( 1 => Dict("color" => "red", "size" => 10), 2 => Dict("color" => "blue", "height" => 20) )

Retrieve the 'color' property for block 1

colorvalue = getproperties(1, "color") # Returns "red"

Try to retrieve a non-existent property for block 2

nonexistentvalue = get_properties(2, "width") # Returns an empty dictionary

source
PeriLab.Data_manager.get_propertyFunction
get_property(block_id::Int64, property::String, value_name::String)

This function retrieves a specific value_name associated with a specified property for a given block_id if it exists in the properties dictionary.

Arguments

  • block_id::Int64: The identifier of the block for which to retrieve the property.
  • property::String: The String property type (e.g. Material model) for the blocks.
  • value_name::String: The name of the value within the specified property.

Returns

  • value::Any: The value associated with the specified value_name within the property for the given block_id.
  • nothing: If the specified block_id, property, or value_name does not exist in the dictionary.

Example

```julia

source
PeriLab.Data_manager.get_rankFunction
get_rank()

This function returns the rank of the core.

Returns

  • rank::Any: The value of the rank variable.

Example

```julia currentrank = getrank()

source
PeriLab.Data_manager.get_max_rankFunction
get_max_rank()

This function returns the maximal rank of MPI the max_rank.

Returns

  • max_rank::Number: The value of the max_rank variable.

Example

```julia rank = getmaxrank()

source
PeriLab.Data_manager.loc_to_globFunction
loc_to_glob(range::UnitRange{Int64})

Converts the local index to the global index.

Arguments

  • range::UnitRange{Int64}: The range of the local index.

Example:

loc_to_glob(1:10)  # converts the local index to the global index
source
PeriLab.Data_manager.init_propertiesFunction
init_properties()

This function initializes the properties dictionary. Order of dictionary defines, in which order the models are called later on.

Returns

  • keys(properties[1]): The keys of the properties dictionary in defined order for the Model_Factory.jl.
source
PeriLab.Data_manager.set_dofFunction
set_dof(n::Int64)

Sets the degree of freedom (dof) value globally.

Arguments

  • n::Int64: The value to set as the degree of freedom.

Example:

set_dof(3)  # sets the degree of freedom to 3
source
PeriLab.Data_manager.set_glob_to_locFunction
set_glob_to_loc(dict)

Sets the global-to-local mapping dict globally.

Arguments

  • dict (array): The dict representing the global-to-local mapping.

Example:

set_glob_to_loc([1, 3, 5])  # sets the global-to-local mapping dict
source
PeriLab.Data_manager.set_num_controllerFunction
set_num_controller(n::Int64)

Sets the number of controller nodes globally. For one core the number of nodes is equal to the number of controller nodes.

Arguments

  • n::Int64: The value to set as the number of nodes.

Example:

set_num_controller(10)  # sets the number of nodes to 10
source
PeriLab.Data_manager.set_nsetFunction
set_nset(name, nodes)

Set the nodes associated with a named node set.

Arguments

  • name::String: The name of the node set.
  • nodes::Vector{Int}: The node indices associated with the node set.
source
PeriLab.Data_manager.set_num_responderFunction
set_num_responder(n::Int64)

Sets the number of responder nodes globally. For one core the number of responder is zero. responder hold the information of the neighbors, of one node, but are not evaluated.

Arguments

  • n::Int64: The value to set as the number of nodes.

Example:

set_num_responder(10)  # sets the number of responder nodes to 10
source
PeriLab.Data_manager.set_propertyFunction
set_property(block_id, property, value_name, value)

Sets the value of a specified property for a given block_id.

Arguments

  • block_id::Int64: The identifier of the block for which to set the property.
  • property::String: The name of the property.
  • value_name::String: The name of the value within the specified property.
  • value::Any: The value to set for the specified value_name.
source
PeriLab.Data_manager.set_propertiesFunction
set_properties(block_id, property, values)

Sets the values of a specified property for a given block_id.

Arguments

  • block_id::Int64: The identifier of the block for which to set the property.
  • property::String: The name of the property.
  • values::Any: The values to set for the specified property.
source
set_properties(property, values)

Sets the values of a specified property for a all blocks. E.g. for FEM, because it corresponds not to a block yet,

Arguments

  • property::String: The name of the property.
  • values::Any: The values to set for the specified property.
source
PeriLab.Data_manager.set_synchFunction
set_synch(name, download_from_cores, upload_to_cores)

Sets the synchronization dictionary globally.

Arguments

  • name::String: The name of the field.
  • download_from_cores::Bool: Whether to download the field from the cores.
  • upload_to_cores::Bool: Whether to upload the field to the cores.
source
PeriLab.Data_manager.synch_managerFunction
synch_manager(synchronise_field, direction::String)

Synchronises the fields.

Arguments

  • synchronise_field: The function to synchronise the field.
  • direction::String: The direction of the synchronisation.
source

Solver

PeriLab.Solver.initFunction
init(params::Dict, datamanager::Module)

Initialize the solver

Arguments

  • params::Dict: The parameters
  • datamanager::Module: Datamanager
  • to::TimerOutputs.TimerOutput: A timer output

Returns

  • block_nodes::Dict{Int64,Vector{Int64}}: A dictionary mapping block IDs to collections of nodes.
  • bcs::Dict{Any,Any}: A dictionary containing boundary conditions.
  • datamanager::Module: The data manager module that provides access to data fields and properties.
  • solver_options::Dict{String,Any}: A dictionary containing solver options.
source
PeriLab.Solver.get_block_nodesFunction
get_block_nodes(block_ids, nnodes)

Returns a dictionary mapping block IDs to collections of nodes.

Arguments

  • block_ids::Vector{Int64}: A vector of block IDs
  • nnodes::Int64: The number of nodes

Returns

  • block_nodes::Dict{Int64,Vector{Int64}}: A dictionary mapping block IDs to collections of nodes
source
PeriLab.Solver.set_densityFunction
set_density(params::Dict, block_nodes::Dict, density::SubArray)

Sets the density of the nodes in the dictionary.

Arguments

  • params::Dict: The parameters
  • block_nodes::Dict: A dictionary mapping block IDs to collections of nodes
  • density::SubArray: The density

Returns

  • density::SubArray: The density
source
PeriLab.Solver.set_horizonFunction
set_horizon(params::Dict, block_nodes::Dict, horizon::SubArray)

Sets the horizon of the nodes in the dictionary.

Arguments

  • params::Dict: The parameters
  • block_nodes::Dict: A dictionary mapping block IDs to collections of nodes
  • horizon::SubArray: The horizon

Returns

  • horizon::SubArray: The horizon
source
PeriLab.Solver.solverFunction
solver(solver_options::Dict{String,Any}, block_nodes::Dict{Int64,Vector{Int64}}, bcs::Dict{Any,Any}, datamanager::Module, outputs::Dict{Int64,Dict{}}, result_files::Vector{Any}, write_results, to, silent::Bool)

Runs the solver.

Arguments

  • solver_options::Dict{String,Any}: The solver options
  • block_nodes::Dict{Int64,Vector{Int64}}: A dictionary mapping block IDs to collections of nodes
  • bcs::Dict{Any,Any}: The boundary conditions
  • datamanager::Module: The data manager module
  • outputs::Dict{Int64,Dict{}}: A dictionary for output settings
  • result_files::Vector{Any}: A vector of result files
  • write_results: A function to write simulation results
  • to::TimerOutputs.TimerOutput: A timer output
  • silent::Bool: A boolean flag to suppress progress bars

Returns

  • result_files: A vector of updated result files
source
PeriLab.Solver.synchronise_fieldFunction
synchronise_field(comm, synch_fields::Dict, overlap_map, get_field, synch_field::String, direction::String)

Synchronises field.

Arguments

  • comm: The MPI communicator
  • synch_fields::Dict: A dictionary of fields
  • overlap_map: The overlap map
  • get_field: The function to get the field
  • synch_field::String: The field
  • direction::String: The direction

Returns

  • nothing
source

Verlet

PeriLab.Solver.Verlet.compute_thermodynamic_critical_time_stepFunction
compute_thermodynamic_critical_time_step(nodes::Union{SubArray,Vector{Int64}}, datamanager::Module, lambda::Float64, Cv::Float64)

Calculate the critical time step for a thermodynamic simulation based on [10].

This function iterates over a collection of nodes and computes the critical time step for each node using provided input data and parameters.

Arguments

  • nodes::Union{SubArray, Vector{Int64}}: The collection of nodes to calculate the critical time step for.
  • datamanager::Module: The data manager module that provides access to required data fields.
  • lambda::Float64: The material parameter used in the calculations.
  • Cv::Float64: The heat capacity at constant volume used in the calculations.

Returns

  • Float64: The calculated critical time step for the thermodynamic simulation.

Dependencies

This function depends on the following data fields from the datamanager module:

  • get_nlist(): Returns the neighbor list.
  • get_field("Density"): Returns the density field.
  • get_field("Bond Length"): Returns the bond distance field.
  • get_field("Volume"): Returns the volume field.
  • get_field("Number of Neighbors"): Returns the number of neighbors field.
source
PeriLab.Solver.Verlet.get_cs_denominatorFunction
get_cs_denominator(volume::Union{SubArray,Vector{Float64},Vector{Int64}}, undeformed_bond::Union{SubArray,Vector{Float64},Vector{Int64}})

Calculate the denominator for the critical time step calculation.

Arguments

  • volume::Union{SubArray,Vector{Float64},Vector{Int64}}: The volume field.
  • undeformed_bond::Union{SubArray,Vector{Float64},Vector{Int64}}: The undeformed bond field.

Returns

  • Float64: The denominator for the critical time step calculation.
source
PeriLab.Solver.Verlet.compute_mechanical_critical_time_stepFunction
compute_mechanical_critical_time_step(nodes::Union{SubArray,Vector{Int64}}, datamanager::Module, bulk_modulus::Float64)

Calculate the critical time step for a mechanical simulation using a bond-based approximation [11].

This function iterates over a collection of nodes and computes the critical time step for each node based on the given input data and parameters.

Arguments

  • nodes::Union{SubArray, Vector{Int64}}: The collection of nodes to calculate the critical time step for.
  • datamanager::Module: The data manager module that provides access to required data fields.
  • bulk_modulus::Float64: The bulk modulus used in the calculations.

Returns

  • Float64: The calculated critical time step for the mechanical simulation.

Dependencies

This function depends on the following data fields from the datamanager module:

  • get_nlist(): Returns the neighbor list.
  • get_field("Density"): Returns the density field.
  • get_field("Bond Length"): Returns the bond distance field.
  • get_field("Volume"): Returns the volume field.
  • get_field("Horizon"): Returns the horizon field.
source
PeriLab.Solver.Verlet.test_timestepFunction
test_timestep(t::Float64, critical_time_step::Float64)

Compare a time step t with a critical time step critical_time_step and update critical_time_step if t is smaller.

Arguments

  • t::Float64: The time step to compare with critical_time_step.
  • critical_time_step::Float64: The current critical time step.

Returns

  • critical_time_step::Float64: The updated critical time step, which is either the original critical_time_step or t, whichever is smaller.
source
PeriLab.Solver.Verlet.compute_crititical_time_stepFunction
compute_crititical_time_step(datamanager::Module, block_nodes::Dict{Int64,Vector{Int64}}, mechanical::Bool, thermo::Bool)

Calculate the critical time step for a simulation considering both mechanical and thermodynamic aspects.

This function computes the critical time step by considering mechanical and thermodynamic properties of different blocks. The resulting critical time step is based on the smallest critical time step found among the blocks.

Arguments

  • datamanager::Module: The data manager module that provides access to required data fields and properties.
  • block_nodes::Dict{Int64, Vector{Int64}}: A dictionary mapping block IDs to collections of nodes.
  • mechanical::Bool: If true, mechanical properties are considered in the calculation.
  • thermo::Bool: If true, thermodynamic properties are considered in the calculation.

Returns

  • Float64: The calculated critical time step based on the smallest critical time step found among the blocks.

Dependencies

This function may depend on the following functions:

  • compute_thermodynamic_critical_time_step: Used if thermo is true to calculate thermodynamic critical time steps.
  • compute_mechanical_critical_time_step: Used if mechanical is true to calculate mechanical critical time steps.
  • The availability of specific properties from the data manager module.

Errors

  • If required properties are not available in the data manager, it may raise an error message.
source
PeriLab.Solver.Verlet.init_solverFunction
init_solver(params::Dict, datamanager::Module, block_nodes::Dict{Int64,Vector{Int64}}, mechanical::Bool, thermo::Bool)

Initialize the Verlet solver for a simulation.

This function sets up the Verlet solver for a simulation by initializing various parameters and calculating the time step based on provided parameters or critical time step calculations.

Arguments

  • params::Dict: A dictionary containing simulation parameters.
  • datamanager::Module: The data manager module that provides access to required data fields and properties.
  • block_nodes::Dict{Int64,Vector{Int64}}: A dictionary mapping block IDs to collections of nodes.
  • mechanical::Bool: If true, mechanical properties are considered in the calculation.
  • thermo::Bool: If true, thermodynamic properties are considered in the calculation.

Returns

A tuple (initial_time, dt, nsteps, numerical_damping) where:

  • initial_time::Float64: The initial time for the simulation.
  • dt::Float64: The time step for the simulation.
  • nsteps::Int64: The number of time integration steps.
  • numerical_damping::Float64: The numerical damping factor.
  • max_damage::Float64: The maximum damage in the simulation.

Dependencies

This function may depend on the following functions:

  • get_initial_time, get_final_time, get_safety_factor, get_fixed_dt: Used to retrieve simulation parameters.
  • compute_crititical_time_step: Used to calculate the critical time step if dt is not fixed.
  • get_integration_steps: Used to determine the number of integration steps and adjust the time step.
  • find_and_set_core_value_min and find_and_set_core_value_max: Used to set core values in a distributed computing environment.
source
PeriLab.Solver.Verlet.get_integration_stepsFunction
get_integration_steps(initial_time::Float64, end_time::Float64, dt::Float64)

Calculate the number of integration steps and the adjusted time step for a numerical integration process.

Arguments

  • initial_time::Float64: The initial time for the integration.
  • end_time::Float64: The final time for the integration.
  • dt::Float64: The time step size.

Returns

A tuple (nsteps, dt) where:

  • nsteps::Int64: The number of integration steps required to cover the specified time range.
  • dt::Float64: The adjusted time step size to evenly divide the time range.

Errors

  • Throws an error if the dt is less than or equal to zero.
source
PeriLab.Solver.Verlet.run_solverFunction
run_solver(
    solver_options::Dict{Any,Any},
    block_nodes::Dict{Int64,Vector{Int64}},
    bcs::Dict{Any,Any},
    datamanager::Module,
    outputs::Dict{Int64,Dict{}},
    result_files::Vector{Any},
    synchronise_field,
    write_results,
    to::TimerOutputs.TimerOutput,
    silent::Bool
)

Run the Verlet solver for a simulation based on the strategy provided in [1] and [2].

This function performs the Verlet solver simulation, updating various data fields and properties over a specified number of time steps.

Arguments

  • solver_options::Dict{String,Any}: A dictionary containing solver options and parameters.
  • block_nodes::Dict{Int64,Vector{Int64}}: A dictionary mapping block IDs to collections of nodes.
  • bcs::Dict{Any,Any}: A dictionary containing boundary conditions.
  • datamanager::Module: The data manager module that provides access to data fields and properties.
  • outputs::Dict{Int64,Dict{}}: A dictionary for output settings.
  • result_files::Vector{Any}: A vector of result files.
  • synchronise_field: A function for synchronization.
  • write_results: A function to write simulation results.
  • to::TimerOutputs.TimerOutput: A timer output.
  • silent::Bool: A boolean flag to suppress progress bars.

Returns

  • result_files: A vector of updated result files.

Dependencies

This function depends on various data fields and properties from the datamanager module, as well as several helper functions. It also relies on solver options and boundary conditions provided in the input parameters.

Function Workflow

  1. Initialize simulation parameters and data fields.
  2. Perform Verlet integration over a specified number of time steps.
  3. Update data fields and properties based on the solver options.
  4. Write simulation results using the write_results function.
source

Boundary_conditions

PeriLab.Solver.Boundary_conditions.check_valid_bcsFunction
check_valid_bcs(bcs::Dict{String,Any}, datamanager::Module

Check if the boundary conditions are valid

Arguments

  • bcs::Dict{String,Any}: The boundary conditions
  • datamanager::Module: The data manager module

Returns

  • working_bcs::Dict{String,Any}: The valid boundary conditions
source
PeriLab.Solver.Boundary_conditions.init_BCsFunction
init_BCs(params::Dict, datamanager)

Initialize the boundary conditions

Arguments

  • params::Dict: The parameters
  • datamanager::Module: Datamanager

Returns

  • bcs::Dict{Any,Any}: The boundary conditions
source
PeriLab.Solver.Boundary_conditions.apply_bc_dirichletFunction
apply_bc_dirichlet(bcs::Dict, datamanager::Module, time::Float64)

Apply the boundary conditions

Arguments

  • bcs::Dict{Any,Any}: The boundary conditions
  • datamanager::Module: Datamanager
  • time::Float64: The current time

Returns

  • datamanager::Module: Datamanager
source
PeriLab.Solver.Boundary_conditions.apply_bc_neumannFunction
apply_bc_neumann(bcs::Dict, datamanager::Module, time::Float64)

Apply the boundary conditions

Arguments

  • bcs::Dict{Any,Any}: The boundary conditions
  • datamanager::Module: Datamanager
  • time::Float64: The current time

Returns

  • datamanager::Module: Datamanager
source
PeriLab.Solver.Boundary_conditions.eval_bcFunction
eval_bc(field_values::Union{SubArray,Vector{Float64},Vector{Int64}}, bc::Union{Float64,Float64,Int64,String}, coordinates::Matrix{Float64}, time::Float64, dof::Int64)

Working with if-statements "if t>2 0 else 20 end" works for scalars. If you want to evaluate a vector, please use the Julia notation as input "ifelse.(x .> y, 10, 20)"

source