Solver - Functions
Index
PeriLab.Solver.Boundary_conditions.apply_bc_dirichlet
PeriLab.Solver.Boundary_conditions.apply_bc_dirichlet_force
PeriLab.Solver.Boundary_conditions.apply_bc_neumann
PeriLab.Solver.Boundary_conditions.boundary_condition
PeriLab.Solver.Boundary_conditions.check_valid_bcs
PeriLab.Solver.Boundary_conditions.clean_up
PeriLab.Solver.Boundary_conditions.eval_bc!
PeriLab.Solver.Boundary_conditions.init_BCs
PeriLab.Solver.Verlet.calculate_strain
PeriLab.Solver.Verlet.calculate_stresses
PeriLab.Solver.Verlet.calculate_von_mises_stress
PeriLab.Solver.Verlet.compute_crititical_time_step
PeriLab.Solver.Verlet.compute_mechanical_critical_time_step
PeriLab.Solver.Verlet.compute_thermodynamic_critical_time_step
PeriLab.Solver.Verlet.find_and_set_core_value_avg
PeriLab.Solver.Verlet.find_and_set_core_value_max
PeriLab.Solver.Verlet.find_and_set_core_value_min
PeriLab.Solver.Verlet.find_and_set_core_value_sum
PeriLab.Solver.Verlet.gather_values
PeriLab.Solver.Verlet.get_cs_denominator
PeriLab.Solver.Verlet.get_forces_from_force_density
PeriLab.Solver.Verlet.get_integration_steps
PeriLab.Solver.Verlet.get_partial_stresses
PeriLab.Solver.Verlet.init_solver
PeriLab.Solver.Verlet.run_solver
PeriLab.Solver.Verlet.send_single_value_from_vector
PeriLab.Solver.Verlet.send_value
PeriLab.Solver.Verlet.send_vector_from_root_to_core_i
PeriLab.Solver.Verlet.split_vector
PeriLab.Solver.Verlet.synch_controller_bonds_to_responder
PeriLab.Solver.Verlet.synch_controller_bonds_to_responder_flattened
PeriLab.Solver.Verlet.synch_controller_to_responder
PeriLab.Solver.Verlet.synch_responder_to_controller
PeriLab.Solver.Verlet.test_timestep
PeriLab.Solver.find_and_set_core_value_avg
PeriLab.Solver.find_and_set_core_value_max
PeriLab.Solver.find_and_set_core_value_min
PeriLab.Solver.find_and_set_core_value_sum
PeriLab.Solver.gather_values
PeriLab.Solver.get_block_nodes
PeriLab.Solver.init
PeriLab.Solver.remove_models
PeriLab.Solver.send_single_value_from_vector
PeriLab.Solver.send_value
PeriLab.Solver.send_vector_from_root_to_core_i
PeriLab.Solver.set_density
PeriLab.Solver.set_fem_block
PeriLab.Solver.set_horizon
PeriLab.Solver.solver
PeriLab.Solver.split_vector
PeriLab.Solver.synch_controller_bonds_to_responder
PeriLab.Solver.synch_controller_bonds_to_responder_flattened
PeriLab.Solver.synch_controller_to_responder
PeriLab.Solver.synch_responder_to_controller
PeriLab.Solver.synchronise_field
Solver
PeriLab.Solver.find_and_set_core_value_avg
— Methodfind_and_set_core_value_avg(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value avg
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.find_and_set_core_value_max
— Methodfind_and_set_core_value_max(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value max
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.find_and_set_core_value_min
— Methodfind_and_set_core_value_min(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value min
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.find_and_set_core_value_sum
— Methodfind_and_set_core_value_sum(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value sum
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.gather_values
— Methodgather_values(comm::MPI.Comm, value::Any)
Gather values
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Any
: The value
Returns
recv_msg::Any
: The received message
PeriLab.Solver.get_block_nodes
— Methodget_block_nodes(block_ids, nnodes)
Returns a dictionary mapping block IDs to collections of nodes.
Arguments
block_ids::Vector{Int64}
: A vector of block IDsnnodes::Int64
: The number of nodes
Returns
block_nodes::Dict{Int64,Vector{Int64}}
: A dictionary mapping block IDs to collections of nodes
PeriLab.Solver.init
— Methodinit(params::Dict, datamanager::Module)
Initialize the solver
Arguments
params::Dict
: The parametersdatamanager::Module
: Datamanagerto::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.
PeriLab.Solver.remove_models
— Methodremove_models(datamanager::Module, solver_options::Vector{String})
Sets the active models to false if they are deactivated in the solver. They can be active, because they are defined as model and in the blocks.
Arguments
datamanager::Module
: The MPI communicatorsolver_options::Vector{String}
: A dictionary of fields
Returns
datamanager
PeriLab.Solver.send_single_value_from_vector
— Methodsend_single_value_from_vector(comm::MPI.Comm, controller::Int64, values::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}, type::Type)
Sends a single value from a vector to a controller
Arguments
comm::MPI.Comm
: The MPI communicatorcontroller::Int64
: The controllervalues::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The valuestype::Type
: The type
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.send_value
— Methodsend_value(comm::MPI.Comm, controller, send_msg)
Sends a value to a controller
Arguments
comm::MPI.Comm
: The MPI communicatorcontroller::Int64
: The controllersend_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The send message
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.send_vector_from_root_to_core_i
— Methodsend_vector_from_root_to_core_i(comm::MPI.Comm, send_msg, recv_msg, distribution)
Sends a vector from the root to the core i
Arguments
comm::MPI.Comm
: The MPI communicatorsend_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The send messagerecv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The receive messagedistribution::Vector{Int64}
: The distribution
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.set_density
— Methodset_density(params::Dict, block_nodes::Dict, density::Vector{Float64})
Sets the density of the nodes in the dictionary.
Arguments
params::Dict
: The parametersblock_nodes::Dict
: A dictionary mapping block IDs to collections of nodesdensity::Vector{Float64}
: The density
Returns
density::Vector{Float64}
: The density
PeriLab.Solver.set_fem_block
— Methodset_fem_block(params::Dict, block_nodes::Dict, fem_block::Vector{Bool})
Sets the fem_block of the nodes in the dictionary.
Arguments
params::Dict
: The parametersblock_nodes::Dict
: A dictionary mapping block IDs to collections of nodesfem_block::Vector{Bool}
: The fem_block
Returns
fem_block::Vector{Bool}
: The fem_block
PeriLab.Solver.set_horizon
— Methodset_horizon(params::Dict, block_nodes::Dict, horizon::Vector{Float64})
Sets the horizon of the nodes in the dictionary.
Arguments
params::Dict
: The parametersblock_nodes::Dict
: A dictionary mapping block IDs to collections of nodeshorizon::Vector{Float64}
: The horizon
Returns
horizon::Vector{Float64}
: The horizon
PeriLab.Solver.solver
— Methodsolver(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 optionsblock_nodes::Dict{Int64,Vector{Int64}}
: A dictionary mapping block IDs to collections of nodesbcs::Dict{Any,Any}
: The boundary conditionsdatamanager::Module
: The data manager moduleoutputs::Dict{Int64,Dict{}}
: A dictionary for output settingsresult_files::Vector{Any}
: A vector of result fileswrite_results
: A function to write simulation resultsto::TimerOutputs.TimerOutput
: A timer outputsilent::Bool
: A boolean flag to suppress progress bars
Returns
result_files
: A vector of updated result files
PeriLab.Solver.split_vector
— Methodsplit_vector(input, row_nums, dof)
Split a vector into a vector of matrices
Arguments
input::Vector
: The input vectorrow_nums::Vector
: The row numbersdof::Int
: The degree of freedom
Returns
result::Vector
: The result vector
PeriLab.Solver.synch_controller_bonds_to_responder
— Methodsynch_controller_bonds_to_responder(comm::MPI.Comm, overlapnodes, array, dof)
Synch the controller bonds to the responder
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesarray::Array
: The arraydof::Int
: The degree of freedom
Returns
array::Array
: The array
PeriLab.Solver.synch_controller_bonds_to_responder_flattened
— Methodsynch_controller_bonds_to_responder_flattened(comm::MPI.Comm, overlapnodes, array, dof)
Synch the controller bonds to the responder
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesarray::Array
: The arraydof::Int
: The degree of freedom
Returns
array::Array
: The array
PeriLab.Solver.synch_controller_to_responder
— Methodsynch_controller_to_responder(comm::MPI.Comm, overlapnodes, vector, dof)
Synch the controller to the responder
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesvector::Vector
: The vectordof::Int
: The degree of freedom
Returns
vector::Vector
: The vector
PeriLab.Solver.synch_responder_to_controller
— Methodsynch_responder_to_controller(comm::MPI.Comm, overlapnodes, vector, dof)
Synch the responder to the controller
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesvector::Vector
: The vectordof::Int
: The degree of freedom
Returns
vector::Vector
: The vector
PeriLab.Solver.synchronise_field
— Methodsynchronise_field(comm, synch_fields::Dict, overlap_map, get_field, synch_field::String, direction::String)
Synchronises field.
Arguments
comm
: The MPI communicatorsynch_fields::Dict
: A dictionary of fieldsoverlap_map
: The overlap mapget_field
: The function to get the fieldsynch_field::String
: The fielddirection::String
: The direction
Returns
nothing
Verlet
PeriLab.Solver.Verlet.calculate_strain
— Methodcalculate_strain(datamanager::Module, nodes::Union{SubArray,Vector{Int64}}, hooke_matrix::Matrix{Float64})
Calculate the von Mises stress.
Arguments
datamanager::Data_manager
: Datamanager.nodes::Union{SubArray,Vector{Int64}}
: The nodes.hooke_matrix::Matrix{Float64}
: The hooke matrix.
Returns
datamanager::Data_manager
: Datamanager.
PeriLab.Solver.Verlet.calculate_stresses
— Methodcalculate_stresses(datamanager::Module, block_nodes::Dict{Int64,Vector{Int64}}, options::Dict{String, Any})
Computes the stresses.
Arguments
datamanager::Data_manager
: Datamanager.block_nodes::Dict{Int64,Vector{Int64}}
: List of block nodes.options::Dict{String, Any}
: List of options.
Returns
datamanager::Data_manager
: Datamanager.
PeriLab.Solver.Verlet.calculate_von_mises_stress
— Methodcalculate_von_mises_stress(datamanager::Module, nodes::Union{SubArray,Vector{Int64}})
Calculate the von Mises stress.
Arguments
datamanager::Data_manager
: Datamanager.nodes::Union{SubArray,Vector{Int64}}
: The nodes.
Returns
datamanager::Data_manager
: Datamanager.
PeriLab.Solver.Verlet.compute_crititical_time_step
— Methodcompute_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
: Iftrue
, mechanical properties are considered in the calculation.thermo::Bool
: Iftrue
, 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 ifthermo
istrue
to calculate thermodynamic critical time steps.compute_mechanical_critical_time_step
: Used ifmechanical
istrue
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.
PeriLab.Solver.Verlet.compute_mechanical_critical_time_step
— Methodcompute_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 [10].
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.
PeriLab.Solver.Verlet.compute_thermodynamic_critical_time_step
— Methodcompute_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 [11].
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.
PeriLab.Solver.Verlet.find_and_set_core_value_avg
— Methodfind_and_set_core_value_avg(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value avg
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.Verlet.find_and_set_core_value_max
— Methodfind_and_set_core_value_max(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value max
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.Verlet.find_and_set_core_value_min
— Methodfind_and_set_core_value_min(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value min
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.Verlet.find_and_set_core_value_sum
— Methodfind_and_set_core_value_sum(comm::MPI.Comm, value::Union{Float64,Int64})
Find and set core value sum
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Union{Float64,Int64}
: The value
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.Verlet.gather_values
— Methodgather_values(comm::MPI.Comm, value::Any)
Gather values
Arguments
comm::MPI.Comm
: The MPI communicatorvalue::Any
: The value
Returns
recv_msg::Any
: The received message
PeriLab.Solver.Verlet.get_cs_denominator
— Methodget_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.
PeriLab.Solver.Verlet.get_forces_from_force_density
— Methodget_forces_from_force_density(datamanager::Module)
Computes the forces from the force densities.
Arguments
datamanager::Data_manager
: Datamanager.
Returns
datamanager::Data_manager
: Datamanager.
PeriLab.Solver.Verlet.get_integration_steps
— Methodget_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.
PeriLab.Solver.Verlet.get_partial_stresses
— Methodget_partial_stresses(datamanager::Module, nodes::Vector{Int64})
Computes the partial stresses.
Arguments
datamanager::Data_manager
: Datamanager.nodes::Vector{Int64}
: List of block nodes.
Returns
datamanager::Data_manager
: Datamanager.
PeriLab.Solver.Verlet.init_solver
— Methodinit_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
: Iftrue
, mechanical properties are considered in the calculation.thermo::Bool
: Iftrue
, 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 ifdt
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
andfind_and_set_core_value_max
: Used to set core values in a distributed computing environment.
PeriLab.Solver.Verlet.run_solver
— Methodrun_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
- Initialize simulation parameters and data fields.
- Perform Verlet integration over a specified number of time steps.
- Update data fields and properties based on the solver options.
- Write simulation results using the
write_results
function.
PeriLab.Solver.Verlet.send_single_value_from_vector
— Methodsend_single_value_from_vector(comm::MPI.Comm, controller::Int64, values::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}, type::Type)
Sends a single value from a vector to a controller
Arguments
comm::MPI.Comm
: The MPI communicatorcontroller::Int64
: The controllervalues::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The valuestype::Type
: The type
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.Verlet.send_value
— Methodsend_value(comm::MPI.Comm, controller, send_msg)
Sends a value to a controller
Arguments
comm::MPI.Comm
: The MPI communicatorcontroller::Int64
: The controllersend_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The send message
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.Verlet.send_vector_from_root_to_core_i
— Methodsend_vector_from_root_to_core_i(comm::MPI.Comm, send_msg, recv_msg, distribution)
Sends a vector from the root to the core i
Arguments
comm::MPI.Comm
: The MPI communicatorsend_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The send messagerecv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The receive messagedistribution::Vector{Int64}
: The distribution
Returns
recv_msg::Union{Int64,Vector{Float64},Vector{Int64},Vector{Bool}}
: The received message
PeriLab.Solver.Verlet.split_vector
— Methodsplit_vector(input, row_nums, dof)
Split a vector into a vector of matrices
Arguments
input::Vector
: The input vectorrow_nums::Vector
: The row numbersdof::Int
: The degree of freedom
Returns
result::Vector
: The result vector
PeriLab.Solver.Verlet.synch_controller_bonds_to_responder
— Methodsynch_controller_bonds_to_responder(comm::MPI.Comm, overlapnodes, array, dof)
Synch the controller bonds to the responder
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesarray::Array
: The arraydof::Int
: The degree of freedom
Returns
array::Array
: The array
PeriLab.Solver.Verlet.synch_controller_bonds_to_responder_flattened
— Methodsynch_controller_bonds_to_responder_flattened(comm::MPI.Comm, overlapnodes, array, dof)
Synch the controller bonds to the responder
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesarray::Array
: The arraydof::Int
: The degree of freedom
Returns
array::Array
: The array
PeriLab.Solver.Verlet.synch_controller_to_responder
— Methodsynch_controller_to_responder(comm::MPI.Comm, overlapnodes, vector, dof)
Synch the controller to the responder
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesvector::Vector
: The vectordof::Int
: The degree of freedom
Returns
vector::Vector
: The vector
PeriLab.Solver.Verlet.synch_responder_to_controller
— Methodsynch_responder_to_controller(comm::MPI.Comm, overlapnodes, vector, dof)
Synch the responder to the controller
Arguments
comm::MPI.Comm
: The MPI communicatoroverlapnodes::Dict
: The overlap nodesvector::Vector
: The vectordof::Int
: The degree of freedom
Returns
vector::Vector
: The vector
PeriLab.Solver.Verlet.test_timestep
— Methodtest_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 withcritical_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 originalcritical_time_step
ort
, whichever is smaller.
Boundary_conditions
PeriLab.Solver.Boundary_conditions.apply_bc_dirichlet
— Methodapply_bc_dirichlet(bcs::Dict, datamanager::Module, time::Float64)
Apply the boundary conditions
Arguments
bcs::Dict{Any,Any}
: The boundary conditionsdatamanager::Module
: Datamanagertime::Float64
: The current time
Returns
datamanager::Module
: Datamanager
PeriLab.Solver.Boundary_conditions.apply_bc_dirichlet_force
— Methodapplybcdirichlet_force(bcs::Dict, datamanager::Module, time::Float64)
Apply the boundary conditions
Arguments
bcs::Dict{Any,Any}
: The boundary conditionsdatamanager::Module
: Datamanagertime::Float64
: The current time
Returns
datamanager::Module
: Datamanager
PeriLab.Solver.Boundary_conditions.apply_bc_neumann
— Methodapply_bc_neumann(bcs::Dict, datamanager::Module, time::Float64)
Apply the boundary conditions
Arguments
bcs::Dict{Any,Any}
: The boundary conditionsdatamanager::Module
: Datamanagertime::Float64
: The current time
Returns
datamanager::Module
: Datamanager
PeriLab.Solver.Boundary_conditions.boundary_condition
— Methodboundary_condition(params::Dict, datamanager)
Initialize the boundary condition
Arguments
params::Dict
: The parametersdatamanager::Module
: Datamanager
Returns
bcs_out::Dict{Any,Any}
: The boundary conditions
PeriLab.Solver.Boundary_conditions.check_valid_bcs
— Methodcheck_valid_bcs(bcs::Dict{String,Any}, datamanager::Module
Check if the boundary conditions are valid
Arguments
bcs::Dict{String,Any}
: The boundary conditionsdatamanager::Module
: The data manager module
Returns
working_bcs::Dict{String,Any}
: The valid boundary conditions
PeriLab.Solver.Boundary_conditions.clean_up
— Methodclean_up(bc::String)
Clean up the boundary condition
Arguments
bc::String
: The boundary condition
Returns
bc::String
: The cleaned up boundary condition
PeriLab.Solver.Boundary_conditions.eval_bc!
— Functioneval_bc!(field_values::Union{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)"
PeriLab.Solver.Boundary_conditions.init_BCs
— Methodinit_BCs(params::Dict, datamanager)
Initialize the boundary conditions
Arguments
params::Dict
: The parametersdatamanager::Module
: Datamanager
Returns
bcs::Dict{Any,Any}
: The boundary conditions