Solver - Functions
Index
PeriLab.Solver_Manager.Boundary_Conditions.apply_bc_dirichletPeriLab.Solver_Manager.Boundary_Conditions.apply_bc_neumannPeriLab.Solver_Manager.Boundary_Conditions.boundary_conditionPeriLab.Solver_Manager.Boundary_Conditions.check_valid_bcsPeriLab.Solver_Manager.Boundary_Conditions.clean_upPeriLab.Solver_Manager.Boundary_Conditions.eval_bc!PeriLab.Solver_Manager.Boundary_Conditions.find_bc_free_dofPeriLab.Solver_Manager.Boundary_Conditions.init_BCsPeriLab.Solver_Manager.Static_Solver.init_solverPeriLab.Solver_Manager.Verlet_Solver.calculate_strainPeriLab.Solver_Manager.Verlet_Solver.calculate_stressesPeriLab.Solver_Manager.Verlet_Solver.calculate_von_mises_stressPeriLab.Solver_Manager.Verlet_Solver.compute_crititical_time_stepPeriLab.Solver_Manager.Verlet_Solver.compute_mechanical_critical_time_stepPeriLab.Solver_Manager.Verlet_Solver.compute_thermodynamic_critical_time_stepPeriLab.Solver_Manager.Verlet_Solver.get_cs_denominatorPeriLab.Solver_Manager.Verlet_Solver.get_forces_from_force_densityPeriLab.Solver_Manager.Verlet_Solver.get_integration_stepsPeriLab.Solver_Manager.Verlet_Solver.get_partial_stressesPeriLab.Solver_Manager.Verlet_Solver.init_solverPeriLab.Solver_Manager.Verlet_Solver.run_solverPeriLab.Solver_Manager.Verlet_Solver.test_timestepPeriLab.Solver_Manager.create_module_specificsPeriLab.Solver_Manager.create_module_specificsPeriLab.Solver_Manager.find_jl_filesPeriLab.Solver_Manager.find_module_filesPeriLab.Solver_Manager.initPeriLab.Solver_Manager.remove_modelsPeriLab.Solver_Manager.set_anglesPeriLab.Solver_Manager.set_densityPeriLab.Solver_Manager.set_fem_blockPeriLab.Solver_Manager.set_horizonPeriLab.Solver_Manager.solverPeriLab.Solver_Manager.synchronise_field
Solver
PeriLab.Solver_Manager.create_module_specifics — Method
create_module_specifics(name::String, module_list::Dict{String,AbstractString}(),specifics::Dict{String,String}(), values::Tuple)Searches for a specific function within a list of modules and calls that function if found.
This function iterates over a list of modules specified in module_list and looks for a module-specific function specified in the specifics dictionary. If the module and function are found, it calls that function with the provided values tuple.
Arguments
name::String: The name to match against the module names.module_list::Dict{String, AbstractString}: A dictionary of module names mapped to abstract strings.specifics::Dict{String, String}: A dictionary specifying the module-specific function to call for each module.values::Tuple: A tuple of values to be passed as arguments to the module-specific function.
Example
module_list = Dict("Module1" => "Module1Name", "Module2" => "Module2Name")
specifics = Dict("Module1Name" => "module1_function", "Module2Name" => "module2_function")
values = (arg1, arg2)
create_module_specifics("Module1Name", module_list, specifics, values)PeriLab.Solver_Manager.create_module_specifics — Method
create_module_specifics(name::String, module_list::Dict{String,AbstractString}(),specifics::Dict{String,String}())
# Returns: the function itselfPeriLab.Solver_Manager.find_jl_files — Method
find_jl_files(directory::AbstractString)Recursively find Julia files (.jl) in a directory.
This function recursively searches for Julia source files with the ".jl" extension in the specified directory and its subdirectories. It returns a vector of file paths for all the found .jl files.
Arguments
directory::AbstractString: The directory in which to search for .jl files.
Returns
A vector of strings, where each string is a file path to a .jl file found in the specified directory and its subdirectories.
Example
jl_files = find_jl_files("/path/to/modules")
for jl_file in jl_files
println("Found Julia file: ", jl_file)
endPeriLab.Solver_Manager.find_module_files — Method
find_module_files(directory::AbstractString, specific::String)Search for Julia modules containing a specific function in a given directory.
This function searches for Julia modules (files with .jl extension) in the specified directory and checks if they contain a specific function. It returns a list of dictionaries where each dictionary contains the file path and the name of the module where the specific function is found.
Arguments
directory::AbstractString: The directory to search for Julia modules.specific::String: The name of the specific function to search for.
Returns
An array of dictionaries, where each dictionary has the following keys:
"File": The file path to the module where the specific function is found."Module Name": The name of the module where the specific function is found.
Example
result = find_module_files("/path/to/modules", "my_function")
for module_info in result
println("Function found in module: ", module_info["Module Name"])
println("Module file path: ", module_info["File"])
endPeriLab.Solver_Manager.init — Function
init(params::Dict)Initialize the solver
Arguments
params::Dict: The parameters
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.solver_options::Dict{String,Any}: A dictionary containing solver options.
PeriLab.Solver_Manager.remove_models — Method
remove_models(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
solver_options::Vector{String}: A dictionary of fields
PeriLab.Solver_Manager.set_angles — Method
set_angles(params::Dict, block_nodes::Dict)Sets the density of the nodes in the dictionary.
Arguments
params::Dict: The parametersblock_nodes::Dict: A dictionary mapping block IDs to collections of nodes
PeriLab.Solver_Manager.set_density — Method
set_density(params::Dict, block_nodes::Dict, density::NodeScalarField{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::NodeScalarField{Float64}: The density
Returns
density::NodeScalarField{Float64}: The density
PeriLab.Solver_Manager.set_fem_block — Method
set_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_Manager.set_horizon — Method
set_horizon(params::Dict, block_nodes::Dict, horizon::NodeScalarField{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::NodeScalarField{Float64}: The horizon
Returns
horizon::NodeScalarField{Float64}: The horizon
PeriLab.Solver_Manager.solver — Method
solver(solver_options::Dict{String,Any}, block_nodes::Dict{Int64,Vector{Int64}}, bcs::Dict{Any,Any}, outputs::Dict{Int64,Dict{}}, result_files::Vector{Any}, write_results, 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 conditionsoutputs::Dict{Int64,Dict{}}: A dictionary for output settingsresult_files::Vector{Any}: A vector of result fileswrite_results: A function to write simulation resultssilent::Bool: A boolean flag to suppress progress bars
Returns
result_files: A vector of updated result files
PeriLab.Solver_Manager.synchronise_field — Method
synchronise_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_Solver
PeriLab.Solver_Manager.Verlet_Solver.calculate_strain — Method
calculate_strain(nodes::AbstractVector{Int64}, hooke_matrix::Matrix{Float64})Calculate the von Mises stress.
Arguments
nodes::AbstractVector{Int64}: The nodes.hooke_matrix::Matrix{Float64}: The hooke matrix.
PeriLab.Solver_Manager.Verlet_Solver.calculate_stresses — Method
calculate_stresses(block_nodes::Dict{Int64,Vector{Int64}}, options::Dict{String, Any})Computes the stresses.
Arguments
block_nodes::Dict{Int64,Vector{Int64}}: List of block nodes.options::Dict{String, Any}: List of options.
PeriLab.Solver_Manager.Verlet_Solver.calculate_von_mises_stress — Method
calculate_von_mises_stress(nodes::AbstractVector{Int64})Calculate the von Mises stress.
Arguments
nodes::AbstractVector{Int64}: The nodes.
PeriLab.Solver_Manager.Verlet_Solver.compute_crititical_time_step — Method
compute_crititical_time_step(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
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 ifthermoistrueto calculate thermodynamic critical time steps.compute_mechanical_critical_time_step: Used ifmechanicalistrueto 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_Manager.Verlet_Solver.compute_mechanical_critical_time_step — Method
compute_mechanical_critical_time_step(nodes::AbstractVector{Int64}, bulk_modulus::Float64)Calculate the critical time step for a mechanical simulation using a bond-based approximation [38].
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::AbstractVector{Int64}: The collection of nodes to calculate the critical time step for.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 Data_Manager 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_Manager.Verlet_Solver.compute_thermodynamic_critical_time_step — Method
compute_thermodynamic_critical_time_step(nodes::AbstractVector{Int64}, lambda::Float64, Cv::Float64)Calculate the critical time step for a thermodynamic simulation based on [24].
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::AbstractVector{Int64}: The collection of nodes to calculate the critical time step for.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 Data_Manager 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_Manager.Verlet_Solver.get_cs_denominator — Method
get_cs_denominator(volume::AbstractVector{Float64}, undeformed_bond::AbstractVector{Float64})Calculate the denominator for the critical time step calculation.
Arguments
volume::AbstractVector{Float64}: 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_Manager.Verlet_Solver.get_forces_from_force_density — Method
get_forces_from_force_density()Computes the forces from the force densities.
PeriLab.Solver_Manager.Verlet_Solver.get_integration_steps — Method
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
dtis less than or equal to zero.
PeriLab.Solver_Manager.Verlet_Solver.get_partial_stresses — Method
get_partial_stresses(nodes::Vector{Int64})Computes the partial stresses.
Arguments
nodes::Vector{Int64}: List of block nodes.
PeriLab.Solver_Manager.Verlet_Solver.init_solver — Method
init_solver(params::Dict, bcs::Dict{Any,Any}, 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.bcs::Dict{Any,Any}: Boundary conditionsblock_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 ifdtis not fixed.get_integration_steps: Used to determine the number of integration steps and adjust the time step.find_and_set_core_value_minandfind_and_set_core_value_max: Used to set core values in a distributed computing environment.
PeriLab.Solver_Manager.Verlet_Solver.run_solver — Method
run_solver(
solver_options::Dict{Any,Any},
block_nodes::Dict{Int64,Vector{Int64}},
bcs::Dict{Any,Any},
outputs::Dict{Int64,Dict{}},
result_files::Vector{Any},
synchronise_field,
write_results,
silent::Bool
)Run the Verlet solver for a simulation based on the strategy provided in [1] and [3].
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.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.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 Data_Manager 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_resultsfunction.
PeriLab.Solver_Manager.Verlet_Solver.test_timestep — Method
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 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_steport, whichever is smaller.
Static_Solver
PeriLab.Solver_Manager.Static_Solver.init_solver — Method
init_solver(params::Dict, bcs::Dict{Any,Any}, block_nodes::Dict{Int64,Vector{Int64}}, mechanical::Bool, thermo::Bool)Initialize the Static solver for a simulation.
This function sets up the Static solver for a simulation by initializing various parameters.
Arguments
params::Dict: A dictionary containing simulation parameters.bcs::Dict{Any,Any}: Boundary conditionsblock_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.xtol::Float64: solution tolerance; mininum step size value between two iterationsftol::Float64: residual tolerance; mininum residual value of the functioniterations::Int64: maximum number of iterations per time stepshow_trace::Bool: show the iteration steps of the solver
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.get_integration_steps: Used to determine the number of integration steps and adjust the time step.
Boundary_Conditions
PeriLab.Solver_Manager.Boundary_Conditions.apply_bc_dirichlet — Method
apply_bc_dirichlet(bcs::Dict, time::Float64)Apply the boundary conditions
Arguments
bcs::Dict{Any,Any}: The boundary conditionstime::Float64: The current time
PeriLab.Solver_Manager.Boundary_Conditions.apply_bc_neumann — Method
apply_bc_neumann(bcs::Dict, time::Float64)Apply the boundary conditions
Arguments
bcs::Dict{Any,Any}: The boundary conditionstime::Float64: The current time
PeriLab.Solver_Manager.Boundary_Conditions.boundary_condition — Method
boundary_condition(params::Dict)Initialize the boundary condition
Arguments
params::Dict: The parameters
Returns
bcs_out::Dict{Any,Any}: The boundary conditions
PeriLab.Solver_Manager.Boundary_Conditions.check_valid_bcs — Method
check_valid_bcs(bcs::Dict{String,Any})Check if the boundary conditions are valid
Arguments
bcs::Dict{String,Any}: The boundary conditions
Returns
working_bcs::Dict{String,Any}: The valid boundary conditions
PeriLab.Solver_Manager.Boundary_Conditions.clean_up — Method
clean_up(bc::String)Clean up the boundary condition
Arguments
bc::String: The boundary condition
Returns
bc::String: The cleaned up boundary condition
PeriLab.Solver_Manager.Boundary_Conditions.eval_bc! — Function
eval_bc!(field_values::Union{NodeScalarField{Float64},NodeScalarField{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_Manager.Boundary_Conditions.find_bc_free_dof — Method
find_bc_free_dof(bcs::Dict{String,Any})Finds all dof without a displacement boundary condition. This tuple vector is stored in the Data_Manager.
Arguments
bcs::Dict{String,Any}: The boundary conditions
Returns
PeriLab.Solver_Manager.Boundary_Conditions.init_BCs — Method
init_BCs(params::Dict)Initialize the boundary conditions
Arguments
params::Dict: The parameters
Returns
bcs::Dict{Any,Any}: The boundary conditions