IO - Functions

Index

IO

PeriLab.IO.merge_exodus_filesFunction
merge_exodus_files(result_files::Vector{Any}, output_dir::String)

Merges exodus output files

Arguments

  • result_files::Vector{Any}: The result files
  • output_dir::String: The file directory
source
PeriLab.IO.close_result_filesFunction
close_result_files(result_files::Vector{Dict})

Closes the result files

Arguments

  • result_files::Vector{Dict}: The result files

Returns

  • true: File is closed
  • false: File was already closed
source
close_result_files(result_files::Vector{Dict}, outputs::Dict{Int64,Dict{}})

Closes the result files if the flush_file flag is not set

Arguments

  • result_files::Vector{Dict}: The result files
  • outputs::Dict{Int64,Dict{}}: The output settings
source
PeriLab.IO.delete_filesFunction
delete_files(result_files::Vector{Dict})

Deletes the result files

Arguments

  • result_files: The result files
source
PeriLab.IO.get_file_sizeFunction
get_file_size(result_files::Vector{Dict})

Gets the file size of the result files

Arguments

  • result_files: The result files

Returns

  • total_file_size: The total file size
source
PeriLab.IO.clearNP1Function
clearNP1(name::String)

Clears the NP1 from the name

Arguments

  • name::String: The name

Returns

  • name::String: The cleared name
source
PeriLab.IO.get_results_mappingFunction
get_results_mapping(params::Dict, path::String, datamanager::Module)

Gets the results mapping

Arguments

  • params::Dict: The parameters
  • path::String: The path
  • datamanager::Module: The datamanager

Returns

  • output_mapping::Dict{Int64,Dict{}}: The results mapping
source
PeriLab.IO.initialize_dataFunction
initialize_data(filename::String, filedirectory::String, datamanager::Module, comm::MPI.Comm, to::TimerOutputs.TimerOutput)

Initialize data.

Arguments

  • filename::String: The name of the input file.
  • filedirectory::String: The directory of the input file.
  • datamanager::Module: The datamanager
  • comm::MPI.Comm: The MPI communicator
  • to::TimerOutputs.TimerOutput: The TimerOutput

Returns

  • data::Dict: The data
source
PeriLab.IO.init_write_resultsFunction
init_write_results(params::Dict, output_dir::String, path::String, datamanager::Module, nsteps::Int64, PERILAB_VERSION::String)

Initialize write results.

Arguments

  • params::Dict: The parameters
  • output_dir::String: The output directory.
  • path::String: The path
  • datamanager::Module: The datamanager
  • nsteps::Int64: The number of steps

Returns

  • result_files::Array: The result files
  • outputs::Dict: The outputs
source
PeriLab.IO.write_resultsFunction
write_results(result_files::Vector{Any}, time::Float64, max_damage::Float64, outputs::Dict, datamanager::Module)

Write results.

Arguments

  • result_files::Vector{Any}: The result files
  • time::Float64: The time
  • max_damage::Float64: The maximum damage
  • outputs::Dict: The outputs
  • datamanager::Module: The datamanager

Returns

  • result_files::Vector{Any}: The result files
source
PeriLab.IO.get_global_valuesFunction
get_global_values(output::Dict, datamanager::Module)

Get global values.

Arguments

  • output::Dict: The output
  • datamanager::Module: The datamanager

Returns

  • global_values::Vector: The global values
source
PeriLab.IO.find_global_core_value!Function
find_global_core_value!(global_value::Union{Int64,Float64}, calculation_type::String, nnodes::Int64, datamanager::Module)

Find global core value.

Arguments

  • global_value::Union{Int64,Float64}: The global value
  • calculation_type::String: The calculation type
  • nnodes::Int64: The number of nodes
  • datamanager::Module: The datamanager

Returns

  • global_value::Union{Int64,Float64}: The global value
source
PeriLab.IO.show_block_summaryFunction
show_block_summary(solver_options::Dict, params::Dict, log_file::String, silent::Bool, comm::MPI.Comm, datamanager::Module)

Show block summary.

Arguments

  • solver_options::Dict: The solver options
  • params::Dict: The params
  • log_file::String: The log file
  • silent::Bool: The silent flag
  • comm::MPI.Comm: The comm
  • datamanager::Module: The datamanager
source
PeriLab.IO.calculate_nodelistFunction
calculate_nodelist(datamanager::Module, field_key::String, dof::Union{Int64,Vector{Int64}}, calculation_type::String, local_nodes::Vector{Int64})

Calculate the global value of a field for a given set of nodes.

Arguments

  • datamanager::Data_manager: Datamanager.
  • field_key::String: Field key.
  • dof::Union{Int64,Vector{Int64}}: Degree of freedom
  • calculation_type::String: Calculation type.
  • local_nodes::Vector{Int64}: Node set.

Returns

  • value::Vector: Global value.
  • nnodes::Int64: Number of nodes.
source
PeriLab.IO.calculate_blockFunction
calculate_block(datamanager::Module, field_key::String, dof::Int64, calculation_type::String, block::Int64)

Calculate the global value of a field for a given block.

Arguments

  • datamanager::Data_manager: Datamanager.
  • field_key::String: Field key.
  • dof::Union{Int64,Vector{Int64}}: Degree of freedom
  • calculation_type::String: Calculation type.
  • block::Int64: Block number.

Returns

  • value::Float64: Global value.
  • nnodes::Int64: Number of nodes.
source
PeriLab.IO.global_value_sumFunction
global_value_sum(field::SubArray, dof::Union{Int64,Vector{Int64}}, nodes::Union{SubArray,Vector{Int64}})

Calculate the global sum of a field for given nodes.

Arguments

  • field::SubArray: Field.
  • dof::Union{Int64,Vector{Int64}: Degree of freedom
  • nodes::Union{SubArray,Vector{Int64}}: Nodes.

Returns

  • returnValue::Vector: Global value.
source
PeriLab.IO.global_value_maxFunction
global_value_max(field::SubArray, dof::Union{Int64,Vector{Int64}}, nodes::Union{SubArray,Vector{Int64}})

Calculate the global maximum of a field for given nodes.

Arguments

  • field::SubArray: Field.
  • dof::Union{Int64,Vector{Int64}}: Degree of freedom
  • nodes::Union{SubArray,Vector{Int64}}: Nodes.

Returns

  • returnValue::Vector: Global value.
source
PeriLab.IO.global_value_minFunction
global_value_min(field::SubArray, dof::Union{Int64,Vector{Int64}}, nodes::Union{SubArray,Vector{Int64}})

Calculate the global minimum of a field for given nodes.

Arguments

  • field::SubArray: Field.
  • dof::Union{Int64,Vector{Int64}}: Degree of freedom
  • nodes::Union{SubArray,Vector{Int64}}: Nodes.

Returns

  • returnValue::Vector: Global value.
source
PeriLab.IO.global_value_avgFunction
global_value_avg(field::SubArray, dof::Union{Int64,Vector{Int64}}, nodes::Union{SubArray,Vector{Int64}})

Calculate the global average of a field for given nodes.

Arguments

  • field::SubArray: Field.
  • dof::Union{Int64,Vector{Int64}}: Degree of freedom
  • nodes::Union{SubArray,Vector{Int64}}: Nodes.

Returns

  • returnValue::Vector: Global value.
source

Read_Mesh

PeriLab.IO.init_dataFunction
init_data(params::Dict, path::String, datamanager::Module, comm::MPI.Comm, to::TimerOutput)

Initializes the data for the mesh.

Arguments

  • params::Dict: The parameters for the simulation.
  • path::String: The path to the mesh file.
  • datamanager::Data_manager: The data manager.
  • comm::MPI.Comm: The MPI communicator.
  • to::TimerOutput: The timer output.

Returns

  • datamanager::Data_manager: The data manager.
  • params::Dict: The parameters for the simulation.
source
PeriLab.IO.create_and_distribute_bond_normFunction
create_and_distribute_bond_norm(comm::MPI.Comm, datamanager::Module, nlist_filtered_ids::Vector{Vector{Int64}}, distribution::Vector{Int64}, bond_norm::Vector{Float64}, dof::Int64)

Create and distribute the bond norm

Arguments

  • comm::MPI.Comm: MPI communicator
  • datamanager::Module: Data manager
  • nlist_filtered_ids::Vector{Vector{Int64}}: The filtered neighborhood list
  • distribution::Vector{Int64}: The distribution
  • bond_norm::Vector{Float64}: The bond norm
  • dof::Int64: The degree of freedom
source
PeriLab.IO.get_local_element_topologyFunction
get_local_element_topology(datamanager::Module, topology::Vector{Vector{Int64}}, distribution::Vector{Int64})

Get the local element topology

Arguments

  • datamanager::Module: The datamanager
  • topology::Vector{Vector{Int64}}: The topology
  • distribution::Vector{Int64}: The distribution

Returns

  • datamanager::Module: The datamanager
source
PeriLab.IO.get_local_overlap_mapFunction
get_local_overlap_map()

Changes entries in the overlap map from the global numbering to the local computer core one.

Arguments

  • overlap_map::Dict{Int64, Dict{Int64, String}}: overlap map with global nodes.
  • distribution::Vector{Vector{Int64}}: global nodes distribution at cores, needed for the gobal to local mapping
  • ranks Array{Int64} : number of used computer cores

Returns

  • overlap_map::Dict{Int64, Dict{Int64, String}}: returns overlap map with local nodes.

Example:

get_local_overlap_map(overlap_map, distribution, ranks)  # returns local nodes
source
PeriLab.IO.local_nodes_from_dictFunction
local_nodes_from_dict(glob_to_loc::Dict{Int,Int}, global_nodes::Vector{Int64})

Changes entries in the overlap map from the global numbering to the local computer core one.

Arguments

  • glob_to_loc::Dict{Int,Int}: global to local mapping
  • global_nodes::Vector{Int64}: global nodes

Returns

  • overlap_map::Dict{Int64, Dict{Int64, String}}: returns overlap map with local nodes.
source
PeriLab.IO.distribute_neighborhoodlist_to_coresFunction
distribute_neighborhoodlist_to_cores(comm::MPI.Comm, datamanager::Module, nlist, distribution)

Distributes the neighborhood list to the cores.

Arguments

  • comm::MPI.Comm: MPI communicator
  • datamanager::Module: Data manager
  • nlist: neighborhood list
  • distribution Array{Int64}: global nodes distribution at cores

Returns

  • datamanager::Module: data manager
source
PeriLab.IO.get_local_neighborsFunction
get_local_neighbors(mapping, nlist_core)

Gets the local neighborhood list from the global neighborhood list

Arguments

  • mapping: mapping function
  • nlist_core: global neighborhood list

Returns

  • nlist_core: local neighborhood list
source
PeriLab.IO.get_bond_geometryFunction
get_bond_geometry(datamanager::Module)

Gets the bond geometry

Arguments

  • datamanager::Module: Data manager

Returns

  • datamanager::Module: data manager
source
PeriLab.IO.define_nsetsFunction
define_nsets(nsets::Dict{String,Vector{Any,Int64}}, datamanager::Module)

Defines the node sets

Arguments

  • nsets::Dict{String,Vector{Any,Int64}}: Node sets read from files
  • datamanager::Module: Data manager
source
PeriLab.IO.distribution_to_coresFunction
distribution_to_cores(comm::MPI.Comm, datamanager::Module, mesh, distribution, dof::Int64)

Distributes the mesh data to the cores

Arguments

  • comm::MPI.Comm: MPI communicator
  • datamanager::Module: Data manager
  • mesh: Mesh
  • distribution Array{Int64}: global nodes distribution at cores
  • dof::Int64: Degree of freedom

Returns

  • datamanager::Module: data manager
source
PeriLab.IO.check_mesh_elementsFunction
check_mesh_elements(mesh, dof)

Process and analyze mesh data to create an dictionary containing information about mesh elements for further processing.

Arguments

  • mesh::DataFrame: The input mesh data represented as a DataFrame.
  • dof::Int64: The degrees of freedom (DOF) for the mesh elements.

Returns

A dictionary containing information about mesh elements, which can be used for further processing or uploading.

Example

```julia meshdata = DataFrame(x1 = [1.0, 2.0, 3.0], x2 = [4.0, 5.0, 6.0], volume = [10.0, 20.0, 30.0]) dof = 3 result = checkmeshelements(meshdata, dof)

source
PeriLab.IO.read_meshFunction
read_mesh(filename::String, params::Dict)

Read mesh data from a file and return it as a DataFrame.

Arguments

  • filename::String: The path to the mesh file.
  • params::Dict: The input parameters.

Returns

  • mesh::DataFrame: The mesh data as a DataFrame.
source
PeriLab.IO.area_of_polygonFunction
area_of_polygon(vertices)

Calculate the area of a polygon.

Arguments

  • vertices: The vertices of the polygon.

Returns

  • area: The area of the polygon.
source
PeriLab.IO.set_dofFunction
set_dof(mesh::DataFrame)

Set the degrees of freedom (DOF) for the mesh elements.

Arguments

  • mesh::DataFrame: The input mesh data represented as a DataFrame.

Returns

  • dof::Int64: The degrees of freedom (DOF) for the mesh elements.
source
PeriLab.IO.load_and_evaluate_meshFunction
load_and_evaluate_mesh(params::Dict, path::String, ranksize::Int64, to::TimerOutput)

Load and evaluate the mesh data.

Arguments

  • params::Dict: The input parameters.
  • path::String: The path to the mesh file.
  • ranksize::Int64: The number of ranks.
  • to::TimerOutput: The timer output

Returns

  • distribution::Array{Int64,1}: The distribution of the mesh elements.
  • mesh::DataFrame: The mesh data as a DataFrame.
  • ntype::Dict: The type of the mesh elements.
  • overlap_map::Array{Array{Int64,1},1}: The overlap map of the mesh elements.
  • nlist::Array{Array{Int64,1},1}: The neighborhood list of the mesh elements.
  • dof::Int64: The degrees of freedom (DOF) for the mesh elements.
  • nsets::Dict: The node sets
  • topology::Int64::Array{Int64,nelement:nodes}`: The topology of elements.
  • el_distribution::Array{Int64,1}: The distribution of the finite elements.
source
PeriLab.IO.create_neighborhoodlistFunction
create_neighborhoodlist(mesh::DataFrame, params::Dict, dof::Int64)

Create the neighborhood list of the mesh elements.

Arguments

  • mesh::DataFrame: The input mesh data represented as a DataFrame.
  • params::Dict: The input parameters.
  • dof::Int64: The degrees of freedom (DOF) for the mesh elements.

Returns

  • nlist::Array{Array{Int64,1},1}: The neighborhood list of the mesh elements.
source
PeriLab.IO.get_number_of_neighbornodesFunction
get_number_of_neighbornodes(nlist::Vector{Vector{Int64}})

Get the number of neighbors for each node.

Arguments

  • nlist::Vector{Vector{Int64}}: The neighborhood list of the mesh elements.

Returns

  • length_nlist::Vector{Int64}: The number of neighbors for each node.
source
PeriLab.IO.node_distributionFunction
node_distribution(nlist::Vector{Vector{Int64}}, size::Int64)

Create the distribution of the nodes.

Arguments

  • nlist::Vector{Vector{Int64}}: The neighborhood list of the mesh elements.
  • size::Int64: The number of ranks.
  • distribution_type::String: The distribution type.

Returns

  • distribution::Vector{Vector{Int64}}: The distribution of the nodes.
  • ptc::Vector{Int64}: Defines at which core / rank each node lies.
  • ntype::Dict: The type of the nodes.
source
PeriLab.IO._init_overlap_map_Function
_init_overlap_map_(size)

Initialize the overlap map.

Arguments

  • size::Int64: The number of ranks.

Returns

  • overlap_map::Dict{Int64,Dict{Int64,Dict{String,Vector{Int64}}}}: The overlap map.
source
PeriLab.IO.create_overlap_mapFunction
create_overlap_map(distribution, ptc, size)

Create the overlap map.

Arguments

  • distribution::Array{Int64,1}: The distribution of the nodes.
  • ptc::Array{Int64,1}: The number of nodes in each rank.
  • size::Int64: The number of ranks.

Returns

  • overlap_map::Dict{Int64,Dict{Int64,Dict{String,Vector{Int64}}}}: The overlap map.
source
PeriLab.IO.create_distributionFunction
create_distribution(nnodes::Int64, size::Int64)

Calculate the initial size of each chunk for a nearly equal number of nodes vs. cores this algorithm might lead to the problem, that the last core is not equally loaded

Arguments

  • nnodes::Int64: The number of nodes.
  • size::Int64: The number of cores.

Returns

  • distribution::Array{Int64,1}: The distribution of the nodes.
  • point_to_core::Array{Int64,1}: The number of nodes in each rank.
source
PeriLab.IO.create_distribution_node_basedFunction
create_distribution_node_based(nnodes::Int64,nlist::Vector{Vector{Int64}}, size::Int64)

Calculate the initial size of each chunk for a nearly equal number of nodes vs. cores this algorithm might lead to the problem, that the last core is not equally loaded

Arguments

  • nnodes::Int64: The number of nodes.
  • nlist::Vector{Vector{Int64}}: The neighborhood list.
  • size::Int64: The number of cores.

Returns

  • distribution::Array{Int64,1}: The distribution of the nodes.
  • point_to_core::Array{Int64,1}: The number of nodes in each rank.
source
PeriLab.IO.create_distribution_neighbor_basedFunction
create_distribution_neighbor_based(nnodes::Int64,nlist::Vector{Vector{Int64}}, size::Int64)

Calculate the initial size of each chunk for a nearly equal number of nodes vs. cores this algorithm might lead to the problem, that the last core is not equally loaded

Arguments

  • nnodes::Int64: The number of nodes.
  • nlist::Vector{Vector{Int64}}: The neighborhood list.
  • size::Int64: The number of cores.

Returns

  • distribution::Array{Int64,1}: The distribution of the nodes.
  • point_to_core::Array{Int64,1}: The number of nodes in each rank.
source
PeriLab.IO.neighborsFunction
neighbors(mesh, params::Dict, coor)

Compute the neighbor list for each node in a mesh based on their proximity using a BallTree data structure.

Arguments

  • mesh: A mesh data structure containing the coordinates and other information.
  • params: paramss needed for computing the neighbor list.
  • coor: A vector of coordinate names along which to compute the neighbor list.

Returns

An array of neighbor lists, where each element represents the neighbors of a node in the mesh.

source
PeriLab.IO.bond_intersects_discFunction
bond_intersects_disc(p0::Vector{Float64}, p1::Vector{Float64}, center::Vector{Float64}, normal::Vector{Float64}, radius::Float64)

Check if a line segment intersects a disk.

Arguments

  • p0::Vector{Float64}: The start point of the line segment.
  • p1::Vector{Float64}: The end point of the line segment.
  • center::Vector{Float64}: The center of the disk.
  • normal::Vector{Float64}: The normal of the plane.
  • radius::Float64: The radius of the disk.

Returns

  • Bool: True if the line segment intersects the disk, False otherwise.
source
PeriLab.IO.bond_intersect_infinite_planeFunction
bond_intersect_infinite_plane(p0::Vector{Float64}, p1::Vector{Float64}, lower_left_corner::Vector{Float64}, normal::Vector{Float64})

Check if a line segment intersects an infinite plane.

Arguments

  • p0::Vector{Float64}: The start point of the line segment.
  • p1::Vector{Float64}: The end point of the line segment.
  • lower_left_corner::Vector{Float64}: The lower left corner of the plane.
  • normal::Vector{Float64}: The normal of the plane.

Returns

  • Bool: True if the line segment intersects the plane, False otherwise.
source
PeriLab.IO.bond_intersect_rectangle_planeFunction
bond_intersect_rectangle_plane(x::Vector{Float64}, lower_left_corner::Vector{Float64}, bottom_unit_vector::Vector{Float64}, normal::Vector{Float64}, side_length::Float64, bottom_length::Float64)

Check if a bond intersects a rectangle plane.

Arguments

  • x::Vector{Float64}: The point.
  • lower_left_corner::Vector{Float64}: The lower left corner of the rectangle.
  • bottom_unit_vector::Vector{Float64}: The unit vector along the bottom of the rectangle.
  • normal::Vector{Float64}: The normal of the plane.
  • side_length::Float64: The side length of the rectangle.
  • bottom_length::Float64: The bottom length of the rectangle.

Returns

  • Bool: True if the point is inside the rectangle, False otherwise.
source
PeriLab.IO.apply_bond_filtersFunction
apply_bond_filters(nlist::Vector{Vector{Int64}}, mesh::DataFrame, params::Dict, dof::Int64)

Apply the bond filters to the neighborhood list.

Arguments

  • nlist::Vector{Vector{Int64}}: The neighborhood list.
  • mesh::DataFrame: The mesh.
  • params::Dict: The parameters.
  • dof::Int64: The degrees of freedom.

Returns

  • nlist::Vector{Vector{Int64}}: The filtered neighborhood list.
  • nlist_filtered_ids::Vector{Vector{Int64}}: The filtered neighborhood list.
source
PeriLab.IO.disk_filterFunction
disk_filter(nnodes::Int64, data::Matrix{Float64}, filter::Dict, nlist::Vector{Vector{Int64}}, dof::Int64)

Apply the disk filter to the neighborhood list.

Arguments

  • nnodes::Int64: The number of nodes.
  • data::Matrix{Float64}: The data.
  • filter::Dict: The filter.
  • nlist::Vector{Vector{Int64}}: The neighborhood list.
  • dof::Int64: The degrees of freedom.

Returns

  • filter_flag::Vector{Bool}: The filter flag.
  • normal::Vector{Float64}: The normal vector of the disk.
source
PeriLab.IO.rectangular_plane_filterFunction
rectangular_plane_filter(nnodes::Int64, data::Matrix{Float64}, filter::Dict, nlist::Vector{Vector{Int64}}, dof::Int64)

Apply the rectangular plane filter to the neighborhood list.

Arguments

  • nnodes::Int64: The number of nodes.
  • data::Matrix{Float64}: The data.
  • filter::Dict: The filter.
  • nlist::Vector{Vector{Int64}}: The neighborhood list.
  • dof::Int64: The degrees of freedom.

Returns

  • filter_flag::Vector{Bool}: The filter flag.
  • normal::Vector{Float64}: The normal vector of the disk.
source
PeriLab.IO.glob_to_locFunction
glob_to_loc(distribution)

Get the global to local mapping

Arguments

  • distribution: The distribution

Returns

  • glob_to_loc: The global to local mapping
source
PeriLab.IO.check_for_duplicate_in_dataframeFunction
check_for_duplicate_in_dataframe(mesh::DataFrame)

check duplicated entries and throws an error if one is there. If not everything is ok.

Arguments

  • mesh::DataFrame: The input mesh data represented as a DataFrame.
source
PeriLab.IO.check_types_in_dataframeFunction
check_types_in_dataframe(mesh::DataFrame)

check if block_id in mesh contains only int.

Arguments

  • mesh::DataFrame: The input mesh data represented as a DataFrame.
source

Geometry

PeriLab.IO.Geometry.bond_geometryFunction
 bond_geometry(nodes::Union{SubArray,Vector{Int64}}, nlist, coor, undeformed_bond, undeformed_bond_length)

Calculate bond geometries between nodes based on their coordinates.

Arguments

  • nodes::Union{SubArray,Vector{Int64}}: A vector of integers representing node IDs.
  • nlist: A data structure (e.g., a list or array) representing neighboring node IDs for each node.
  • coor: A matrix representing the coordinates of each node.
  • undeformed_bond: A preallocated array or data structure to store bond geometries.
  • undeformed_bond_length: A preallocated array or data structure to store bond distances.

Output

  • undeformed_bond: An updated undeformed_bond array with calculated bond geometries.

Description

This function calculates bond geometries between nodes. For each node in nodes, it computes the bond vector between the node and its neighboring nodes based on their coordinates. It also calculates the distance (magnitude) of each bond vector.

If the distance of any bond vector is found to be zero, indicating identical point coordinates, an error is raised.

Example

```julia nodes = [1, 2, 3] dof = 2 nlist = [[2, 3], [1, 3], [1, 2]] coor = [0.0 0.0; 1.0 0.0; 0.0 1.0] undeformed_bond = zeros(Float64, length(nodes), length(nlist[1]), dof + 1)

undeformedbond(nodes, dof, nlist, coor, undeformedbond)

source
PeriLab.IO.Geometry.shape_tensorFunction
shape_tensor(nodes::Union{SubArray, Vector{Int64}}, dof::Int64, nlist, volume, omega, bond_damage, undeformed_bond, shape_tensor, inverse_shape_tensor)

Calculate the shape tensor and its inverse for a set of nodes in a computational mechanics context.

Arguments

  • nodes::Union{SubArray, Vector{Int64}}: A vector of integers representing node IDs.
  • dof::Int64: An integer representing the degrees of freedom.
  • nlist: A data structure (e.g., a list or array) representing neighboring node IDs for each node.
  • volume: A vector or array containing volume information for each node.
  • omega: A vector or array containing omega information for each node.
  • bond_damage: A data structure representing bond damage for each node.
  • undeformed_bond: A data structure representing bond geometries for each node.
  • shape_tensor: A preallocated 3D array to store the shape tensors for each node.
  • inverse_shape_tensor: A preallocated 3D array to store the inverse shape tensors for each node.

Output

  • shape_tensor: An updated shape_tensor array with calculated shape tensors.
  • inverse_shape_tensor: An updated inverse_shape_tensor array with calculated inverse shape tensors.

Description

This function calculates the shape tensor and its inverse for a set of nodes in a computational mechanics context. The shape tensor is a key quantity used in continuum mechanics to describe material deformation. It is calculated based on bond damage, bond geometries, volume, and omega information for each node.

For each node in nodes, the function iterates through degrees of freedom (dof) and computes elements of the shape tensor. The inverse of the shape tensor is also calculated and stored in inverse_shape_tensor.

Example

```julia nodes = [1, 2, 3] dof = 3 nlist = [[2, 3], [1, 3], [1, 2]] volume = [0.1, 0.2, 0.3] omega = [0.5, 0.4, 0.6] bonddamage = zeros(Float64, length(nodes), length(nlist[1])) undeformedbond = rand(Float64, length(nodes), length(nlist[1]), dof) shapetensor = zeros(Float64, length(nodes), dof, dof) inverseshape_tensor = zeros(Float64, length(nodes), dof, dof)

shapetensor(nodes, dof, nlist, volume, omega, bonddamage, undeformedbond, shapetensor, inverseshapetensor)

source
PeriLab.IO.Geometry.compute_deformation_gradientFunction
compute_deformation_gradient(nodes::Union{SubArray, Vector{Int64}}, dof::Int64, nlist, volume, omega, bond_damage, undeformed_bond, deformed_bond, inverse_shape_tensor, deformation_gradient)

Calculate the deformation gradient tensor for a set of nodes in a computational mechanics context.

Arguments

  • nodes::Union{SubArray, Vector{Int64}}: A vector of integers representing node IDs.
  • dof::Int64: An integer representing the degrees of freedom.
  • nlist: A data structure (e.g., a list or array) representing neighboring node IDs for each node.
  • volume: A vector or array containing volume information for each node.
  • omega: A vector or array containing omega information for each node.
  • bond_damage: A data structure representing bond damage for each node.
  • undeformed_bond: A data structure representing bond geometries for each node.
  • deformed_bond: A data structure representing deformed bond properties for each node.
  • inverse_shape_tensor: A data structure representing the inverse shape tensors for each node.
  • deformation_gradient: A preallocated 3D array to store the deformation gradient tensors for each node.

Output

  • deformation_gradient: An updated deformation_gradient array with calculated deformation gradient tensors.

Description

This function calculates the deformation gradient tensor for a set of nodes in a computational mechanics context. The deformation gradient tensor characterizes the deformation of a material.

For each node in nodes, the function iterates through degrees of freedom (dof) and computes elements of the deformation gradient tensor based on bond damage, deformed bond properties, bond geometries, volume, and omega information. The calculated deformation gradient tensor is stored in deformation_gradient.

Example

```julia nodes = [1, 2, 3] dof = 3 nlist = [[2, 3], [1, 3], [1, 2]] volume = [0.1, 0.2, 0.3] omega = [0.5, 0.4, 0.6] bonddamage = zeros(Float64, length(nodes), length(nlist[1])) undeformedbond = rand(Float64, length(nodes), length(nlist[1]), dof) deformedbond = rand(Float64, length(nodes), length(nlist[1]), dof) inverseshapetensor = rand(Float64, length(nodes), dof, dof) deformationgradient = zeros(Float64, length(nodes), dof, dof)

computedeformationgradient(nodes, dof, nlist, volume, omega, bonddamage, undeformedbond, deformedbond, inverseshapetensor, deformationgradient)

source
PeriLab.IO.Geometry.compute_strainFunction
function compute_strain(nodes::Union{Base.OneTo{Int64},Vector{Int64}, SubArray}, deformation_gradient, strain)

Calculate strains for specified nodes based on deformation gradients.

Arguments

  • nodes::Union{SubArray, Vector{Int64}}: List of nodes
  • deformation_gradient: Deformation gradient at current time step (2D or 3D array).

Returns

  • Updated strain array containing strains.

This function iterates over the specified nodes and computes strain at each node using the given deformation gradients.

source
PeriLab.IO.Geometry.rotation_tensorFunction
function rotation_tensor(angles::Vector{Float64})

Creates the rotation tensor for 2D or 3D applications. Uses Rotations.jl package.

Arguments

  • angles::Vector{Float64}: Vector of angles definede in degrees of length one or three

Returns

  • Rotation tensor
source

ReadInputDeck

PeriLab.IO.read_inputFunction
read_input(filename::String)

Reads the input deck from a yaml file

Arguments

  • filename::String: The name of the yaml file

Returns

  • params::Dict{String,Any}: The parameters read from the yaml file
source
PeriLab.IO.read_input_fileFunction
read_input_file(filename::String)

Reads the input deck from a yaml file

Arguments

  • filename::String: The name of the yaml file

Returns

  • Dict{String,Any}: The validated parameters read from the yaml file.
source

WriteExodusResults

PeriLab.IO.create_result_fileFunction
create_result_file(filename::Union{AbstractString,String}, num_nodes::Int64, num_dim::Int64, num_elem_blks::Int64, num_node_sets::Int64)

Creates a exodus file for the results

Arguments

  • filename::Union{AbstractString,String}: The name of the file to create
  • num_nodes::Int64: The number of nodes
  • num_dim::Int64: The number of dimensions
  • num_elem_blks::Int64: The number of element blocks
  • num_node_sets::Int64: The number of node sets

Returns

  • result_file::Dict{String,Any}: A dictionary containing the filename and the exodus file
source
create_result_file(filename::String, outputs::Dict)

Creates a csv file for the results

Arguments

  • filename::String: The name of the file to create
  • outputs::Dict: The outputs dictionary

Returns

  • Dict: The result file
source
PeriLab.IO.paraview_specificsFunction
paraview_specifics(dof::Int64)

Returns the paraview specific dof

Arguments

  • dof::Int64: The degrees of freedom

Returns

  • paraview_specifics::String: The paraview specific dof
source
PeriLab.IO.get_paraview_coordinatesFunction
get_paraview_coordinates(dof::Int64, refDof::Int64)

Returns the paraview specific dof

Arguments

  • dof::Int64: The degrees of freedom
  • refDof::Int64: The reference degrees of freedom

Returns

  • paraview_specifics::String: The paraview specific dof
source
PeriLab.IO.get_block_nodesFunction
get_block_nodes(block_Id::Union{SubArray,Vector{Int64}}, block::Int64)

Returns the nodes of a block

Arguments

  • block_Id::Union{SubArray,Vector{Int64}}: The block Id
  • block::Int64: The block

Returns

  • nodes::Vector{Int64}: The nodes of the block
source
PeriLab.IO.init_results_in_exodusFunction
init_results_in_exodus(exo::ExodusDatabase, output::Dict{}, coords::Union{Matrix{Int64},Matrix{Float64}}, block_Id::Vector{Int64}, uniqueBlocks::Vector{Int64}, nsets::Dict{String,Vector{Int64}}, global_ids::Vector{Int64}, PERILAB_VERSION::String)

Initializes the results in exodus

Arguments

  • exo::ExodusDatabase: The exodus database
  • output::Dict{String,Any}: The output
  • coords::Union{Matrix{Int64},Matrix{Float64}}: The coordinates
  • block_Id::Vector{Int64}: The block Id
  • uniqueBlocks::Vector{Int64}: The unique blocks
  • nsets::Dict{String,Vector{Int64}}: The node sets
  • global_ids::Vector{Int64}: The global ids

Returns

  • result_file::Dict{String,Any}: The result file
source
PeriLab.IO.write_step_and_timeFunction
write_step_and_time(exo::ExodusDatabase, step::Int64, time::Float64)

Writes the step and time in the exodus file

Arguments

  • exo::ExodusDatabase: The exodus file
  • step::Int64: The step
  • time::Float64: The time

Returns

  • exo::ExodusDatabase: The exodus file
source
PeriLab.IO.write_nodal_results_in_exodusFunction
write_nodal_results_in_exodus(exo::ExodusDatabase, step::Int64, output::Dict, datamanager::Module)

Writes the nodal results in the exodus file

Arguments

  • exo::ExodusDatabase: The exodus file
  • step::Int64: The step
  • output::Dict: The output
  • datamanager::Module: The datamanager

Returns

  • exo::ExodusDatabase: The exodus file
source
PeriLab.IO.write_global_results_in_exodusFunction
write_global_results_in_exodus(exo::ExodusDatabase, step::Int64, global_values)

Writes the global results in the exodus file

Arguments

  • exo::ExodusDatabase: The exodus file
  • step::Int64: The step
  • global_values: The global values

Returns

  • exo::ExodusDatabase: The exodus file
source
PeriLab.IO.merge_exodus_fileFunction
merge_exodus_file(file_name::Union{AbstractString,String})

Merges the exodus file

Arguments

  • file_name::Union{AbstractString,String}: The name of the file to merge

Returns

  • exo::ExodusDatabase: The exodus file
source

WriteCSVResults

PeriLab.IO.write_global_results_in_csvFunction
write_global_results_in_csv(csv_file::IOStream, time::Float64, global_values)

Writes the global results to the csv file

Arguments

  • csv_file::IOStream: The csv file
  • global_values: The global values
source

Helpers

PeriLab.Solver.Helpers.find_indicesFunction
find_indices(vector, what)

Returns the indices of vector that are equal to what.

Arguments

  • vector::Vector: The vector to search in.
  • what: The value to search for.

Returns

  • indices::Vector: The indices of vector that are equal to what.
source
PeriLab.Solver.Helpers.find_activeFunction
find_active(active::Vector{Bool})

Returns the indices of active that are true.

Arguments

  • active::Vector{Bool}: The vector to search in.

Returns

  • indices::Vector: The indices of active that are true.
source
PeriLab.Solver.Helpers.get_active_update_nodesFunction
get_active_update_nodes(active::SubArray, update_list::SubArray, block_nodes::Dict{Int64,Vector{Int64}}, block::Int64)

Returns the active nodes and the update nodes.

Arguments

  • active::SubArray: The active vector.
  • update_list::SubArray: The update vector.
  • block_nodes::Dict{Int64,Vector{Int64}}: The vector to search in.
  • block::Int64: The block_id.

Returns

  • active_nodes::Vector{Int64}: The nodes of active that are true.
  • update_nodes::Vector{Int64}: The nodes of update that are true.
source
PeriLab.Solver.Helpers.find_files_with_endingFunction
find_files_with_ending(folder_path::AbstractString, file_ending::AbstractString)

Returns a list of files in folder_path that end with file_ending.

Arguments

  • folder_path::AbstractString: The path to the folder.
  • file_ending::AbstractString: The ending of the files.

Returns

  • file_list::Vector{String}: The list of files that end with file_ending.
source
PeriLab.Solver.Helpers.check_inf_or_nanFunction
check_inf_or_nan(array, msg)

Checks if the sum of the array is finite. If not, an error is raised.

Arguments

  • array: The array to check.
  • msg: The error message to raise.

Returns

  • Bool: true if the sum of the array is finite, false otherwise.
source
PeriLab.Solver.Helpers.matrix_styleFunction
matrix_style(A)

Include a scalar or an array and reshape it to style needed for LinearAlgebra package

Arguments

  • A: The array or scalar to reshape

Returns

  • Array: The reshaped array
source
PeriLab.Solver.Helpers.progress_barFunction
progress_bar(rank::Int64, nsteps::Int64, silent::Bool)

Create a progress bar if the rank is 0. The progress bar ranges from 1 to nsteps + 1.

Arguments

  • rank::Int64: An integer to determine if the progress bar should be created.
  • nsteps::Int64: The total number of steps in the progress bar.
  • silent::Bool: de/activates the progress bar

Returns

  • ProgressBar or UnitRange: If rank is 0, a ProgressBar object is returned. Otherwise, a range from 1 to nsteps + 1 is returned.
source
PeriLab.Solver.Helpers.get_fourth_orderFunction
get_fourth_order(CVoigt, dof)

Constructs a symmetric fourth-order tensor from a Voigt notation vector. It uses Tensors.jl package.

This function takes a Voigt notation vector CVoigt and the degree of freedom dof to create a symmetric fourth-order tensor. The CVoigt vector contains components that represent the tensor in Voigt notation, and dof specifies the dimension of the tensor.

Arguments

  • CVoigt::Matrix{Float64}: A vector containing components of the tensor in Voigt notation.
  • dof::Int64: The dimension of the resulting symmetric fourth-order tensor.

Returns

  • SymmetricFourthOrderTensor{dof}: A symmetric fourth-order tensor of dimension dof.

Example

```julia CVoigt = [1.0, 2.0, 3.0, 4.0, 5.0, 6.0] dof = 3 result = getfourthorder(CVoigt, dof)

source
PeriLab.Solver.Helpers.find_inverse_bond_idFunction
find_inverse_bond_id(nlist::SubArray)

Finds the inverse of the bond id in the nlist.

Arguments

  • nlist::SubArray: The nlist to find the inverse of.

Returns

  • inverse_nlist::Vector{Dict{Int64,Int64}}: The inverse nlist.
source
PeriLab.Solver.Helpers.rotateFunction
rotate(nodes::Union{SubArray,Vector{Int64}}, dof::Int64, matrix::Union{SubArray,Array{Float64,3}}, angles::SubArray, back::Bool)

Rotates the matrix.

Arguments

  • nodes::Union{SubArray,Vector{Int64}}: List of block nodes.
  • matrix::Union{SubArray,Array{Float64,3}}: Matrix.
  • rot::SubArray: Rotation tensor.
  • back::Bool: Back.

Returns

  • matrix::SubArray: Matrix.
source
PeriLab.Solver.Helpers.rotate_second_order_tensorFunction
rotate_second_order_tensor(angles::Union{Vector{Float64},Vector{Int64}}, tensor::Matrix{Float64}, dof::Int64, back::Bool)

Rotates the second order tensor.

Arguments

  • angles::Union{Vector{Float64},Vector{Int64}}: Angles.
  • tensor::Matrix{Float64}: Second order tensor.
  • dof::Int64: Degree of freedom.
  • back::Bool: Back.

Returns

  • tensor::Matrix{Float64}: Second order tensor.
source

Pre_Calculation

PeriLab.Solver.Model_Factory.Pre_Calculation.init_modelFunction
init_model(datamanager::Module, nodes::Union{SubArray,Vector{Int64}, block::Int64)

Initializes the model.

Arguments

  • datamanager::Data_manager: Datamanager
  • nodes::Union{SubArray,Vector{Int64}}: The nodes.
  • block::Int64: Block.

Returns

  • datamanager::Data_manager: Datamanager.
source
PeriLab.Solver.Model_Factory.Pre_Calculation.compute_modelFunction
compute_model(datamanager::Module, nodes::Union{SubArray,Vector{Int64}}, model_param::Dict, block::Int64, time::Float64, dt::Float64,to::TimerOutput,)

Computes the pre calculation models

Arguments

  • datamanager::Module: The datamanager
  • nodes::Union{SubArray,Vector{Int64}}: The nodes
  • model_param::Dict: The model parameters
  • block::Int64: The block
  • time::Float64: The current time
  • dt::Float64: The time step

Returns

  • datamanager::Module: The datamanager
source
PeriLab.Solver.Model_Factory.Pre_Calculation.check_dependenciesFunction
check_dependencies(datamanager::Module, block_nodes::Dict{Int64,Vector{Int64}}

Check if materials are used which needs a form of pre calculation. If so, the option will be set.

Arguments

  • datamanager::Module: Datamanager.
  • block_nodes::Dict{Int64,Vector{Int64}}: block nodes.

Returns

  • datamanager::Data_manager: Datamanager.
source

parameter_handling

PeriLab.IO.Parameter_Handling.validate_structure_recursiveFunction
validate_structure_recursive(expected::Dict, actual::Dict, validate::Bool, checked_keys::Array, path::String="")

Validates the parameters against the expected structure

Arguments

  • expected::Dict: The expected structure
  • actual::Dict: The actual structure
  • validate::Bool: The validation results
  • checked_keys::Array: The keys that have been checked
  • path::String: The current path

Returns

  • validate::Bool: The validation result
  • checked_keys::Array: The keys that have been checked
source
PeriLab.IO.Parameter_Handling.get_densityFunction
get_density(params::Dict, block_id::Int64)

Get the density of a block.

Arguments

  • params::Dict: The parameters
  • block_id::Int64: The ID of the block

Returns

  • density::Float64: The density of the block
source
PeriLab.IO.Parameter_Handling.get_heat_capacityFunction
get_heat_capacity(params::Dict, block_id::Int64)

Get the heat capacity of a block.

Arguments

  • params::Dict: The parameters
  • block_id::Int64: The ID of the block

Returns

  • heat_capacity::Float64: The heat capacity of the block
source
PeriLab.IO.Parameter_Handling.get_horizonFunction
get_horizon(params::Dict, block_id::Int64)

Get the horizon of a block.

Arguments

  • params::Dict: The parameters
  • block_id::Int64: The ID of the block

Returns

  • horizon::Float64: The horizon of the block
source
PeriLab.IO.Parameter_Handling.get_valuesFunction
get_values(params::Dict, block_id::Int64, valueName::String, defaultValue::Union{Float64,Bool,Nothing})

Get the value of a block.

Arguments

  • params::Dict: The parameters
  • block_id::Int64: The ID of the block
  • valueName::String: The name of the value
  • defaultValue::Union{Float64,Bool,Nothing: The default value

Returns

  • value::Float64: The value of the block
source
PeriLab.IO.Parameter_Handling.get_block_modelsFunction
get_block_models(params::Dict, block_id::Int64)

Get the models of a block.

Arguments

  • params::Dict: The parameters
  • block_id::Int64: The ID of the block

Returns

  • modelDict::Dict{String,String}: The models of the block
source
PeriLab.IO.Parameter_Handling.get_computesFunction
get_computes(params::Dict, variables::Vector{String})

Get the computes.

Arguments

  • params::Dict: The parameters dictionary.
  • variables::Vector{String}: The variables.

Returns

  • computes::Dict{String,Dict{Any,Any}}: The computes.
source
PeriLab.IO.Parameter_Handling.get_bond_filtersFunction
get_bond_filters(params::Dict)

Returns the bond filters from the parameters

Arguments

  • params::Dict: The parameters

Returns

  • check::Bool: Whether the bond filters are defined
  • bfList::Dict{String,Dict{String,Any}}: The bond filters
source
PeriLab.IO.Parameter_Handling.get_node_setsFunction
get_node_sets(params::Dict, path::String)

Returns the node sets from the parameters

Arguments

  • params::Dict: The parameters
  • path::String: The path to the mesh file

Returns

  • nsets::Dict{String,Any}: The node sets
source
PeriLab.IO.Parameter_Handling.get_output_fieldnamesFunction
get_output_fieldnames(outputs::Dict, variables::Vector{String}, computes::Vector{String}, output_type::String)

Gets the output fieldnames.

Arguments

  • outputs::Dict: The outputs
  • variables::Vector{String}: The variables
  • computes::Vector{String}: The computes
  • output_type::String: The output type

Returns

  • output_fieldnames::Vector{String}: The output fieldnames
source
PeriLab.IO.Parameter_Handling.get_outputsFunction
get_outputs(params::Dict, variables::Vector{String}, compute_names::Vector{String})

Gets the outputs.

Arguments

  • params::Dict: The parameters
  • variables::Vector{String}: The variables
  • compute_names::Vector{String}: The compute names

Returns

  • outputs::Dict: The outputs
source
PeriLab.IO.Parameter_Handling.get_model_parameterFunction
get_model_parameter(params, model, id)

Retrieve a model parameter from a dictionary of parameters.

This function retrieves a specific model parameter from a dictionary of parameters based on the provided model and identifier (id).

Arguments

  • params::Dict: A dictionary containing various parameters.

  • model::String: The model type for which the parameter is sought.

  • id::String: The identifier (name) of the specific model parameter.

Returns

  • parameter::Any: The retrieved model parameter, or nothing if the parameter is not found.

Errors

  • If the specified model is defined in blocks but no model definition block exists, an error message is logged, and the function returns nothing.

  • If the model with the given identifier is defined in blocks but missing in the model's definition, an error message is logged, and the function returns nothing.

Example

```julia params = Dict( "Models" => Dict( "Models" => Dict( "ModelA" => 42, "ModelB" => 24 ) ) )

model = "Models" id = "ModelA"

result = getmodelparameter(params, model, id) if result !== nothing println("Parameter id: result") else println("Parameter not found.") end

source
PeriLab.IO.Parameter_Handling.get_models_optionFunction
get_models_option(params, options)

Process models-related options based on the provided parameters.

This function processes models-related options based on the parameters dictionary and updates the options dictionary accordingly.

Arguments

  • params::Dict: A dictionary containing various parameters, including models-related information.

  • options::Dict: A dictionary containing options to be updated based on the models parameters.

Returns

  • updated_options::Dict: A dictionary containing updated options based on the models parameters.

Errors

  • If the 'Pre Calculation' section exists in the 'Models' block but does not contain required options, an error message is logged.

  • If a material model is missing the 'Material Model' specification, an error message is logged.

Example

```julia params = Dict( "Models" => Dict( "Pre Calculation" => Dict( "Option1" => true, "Option2" => false ), "Material Models" => Dict( 1 => Dict( "Material Model" => "Correspondence" ), 2 => Dict( "Material Model" => "Bond Associated" ) ) ) )

options = Dict( "Option1" => false, "Option2" => true )

updatedoptions = getmodelsoption(params, options) println("Updated Options: ", updatedoptions)

source
PeriLab.IO.Parameter_Handling.get_headerFunction
get_header(filename::Union{String,AbstractString})

Returns the header line and the header.

Arguments

  • filename::Union{String,AbstractString}: The filename of the file.

Returns

  • header_line::Int: The header line.
  • header::Vector{String}: The header.
source

Set_modules

PeriLab.Solver.Model_Factory.Material.Set_modules.find_jl_filesFunction
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

```julia jlfiles = findjlfiles("/path/to/modules") for jlfile in jlfiles println("Found Julia file: ", jlfile) end

source
PeriLab.Solver.Model_Factory.Material.Set_modules.find_module_filesFunction
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

```julia result = findmodulefiles("/path/to/modules", "myfunction") for moduleinfo in result println("Function found in module: ", moduleinfo["Module Name"]) println("Module file path: ", moduleinfo["File"]) end

source
PeriLab.Solver.Model_Factory.Material.Set_modules.include_filesFunction
include_files(module_list::Vector{Any})

Include files specified in a list of modules.

This function iterates over a list of modules and includes the files specified in each module's "File" key.

Arguments

  • module_list::Vector{Any}: A list of modules where each module is expected to be a dictionary-like object with a "File" key specifying the file path.

Examples

```julia include_files([Dict("File" => "module1.jl"), Dict("File" => "module2.jl")])

source
PeriLab.Solver.Model_Factory.Material.Set_modules.create_module_specificsFunction
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

```julia modulelist = Dict("Module1" => "Module1Name", "Module2" => "Module2Name") specifics = Dict("Module1Name" => "module1function", "Module2Name" => "module2function") values = (arg1, arg2) createmodulespecifics("Module1Name", modulelist, specifics, values)

source
create_module_specifics(name::String, module_list::Dict{String,AbstractString}(),specifics::Dict{String,String}())
# Returns: the function itself
source

Logging_module

PeriLab.Logging_module.progress_filterFunction
progress_filter(log_args)

Filter progress messages.

Arguments

  • log_args: The log arguments.

Returns

  • true: If the message is not a progress message.
  • false: If the message is a progress message.
source
PeriLab.Logging_module.init_loggingFunction
init_logging(filename::String, debug::Bool, silent::Bool, rank::Int64, size::Int64)

Initialize the logging.

Arguments

  • filename::String: The filename.
  • debug::Bool: If debug is true.
  • silent::Bool: If silent is true.
  • rank::Int64: The rank.
  • size::Int64: The size.
source