IO - Functions
Index
PeriLab.IO._init_overlap_map_PeriLab.IO.apply_bond_filtersPeriLab.IO.area_of_polygonPeriLab.IO.bond_intersect_infinite_planePeriLab.IO.bond_intersect_rectangle_planePeriLab.IO.bond_intersects_discPeriLab.IO.calculate_blockPeriLab.IO.calculate_nodelistPeriLab.IO.calculate_volumePeriLab.IO.check_for_duplicate_in_dataframePeriLab.IO.check_mesh_elementsPeriLab.IO.check_types_in_dataframePeriLab.IO.clearNP1PeriLab.IO.close_result_filePeriLab.IO.close_result_filesPeriLab.IO.close_result_filesPeriLab.IO.create_and_distribute_bond_normPeriLab.IO.create_distributionPeriLab.IO.create_distribution_neighbor_basedPeriLab.IO.create_distribution_node_basedPeriLab.IO.create_global_to_local_mappingPeriLab.IO.create_neighborhoodlistPeriLab.IO.create_overlap_mapPeriLab.IO.create_result_filePeriLab.IO.create_result_filePeriLab.IO.csv_readerPeriLab.IO.define_nsetsPeriLab.IO.delete_filesPeriLab.IO.disk_filterPeriLab.IO.distribute_neighborhoodlist_to_coresPeriLab.IO.distribution_to_coresPeriLab.IO.element_distributionPeriLab.IO.extrudePeriLab.IO.extrude_surface_meshPeriLab.IO.find_global_core_value!PeriLab.IO.get_block_nodesPeriLab.IO.get_bond_geometryPeriLab.IO.get_file_sizePeriLab.IO.get_global_valuesPeriLab.IO.get_local_element_topologyPeriLab.IO.get_local_neighborsPeriLab.IO.get_local_overlap_mapPeriLab.IO.get_mpi_rank_stringPeriLab.IO.get_number_of_neighbornodesPeriLab.IO.get_paraview_coordinatesPeriLab.IO.get_results_mappingPeriLab.IO.global_value_avgPeriLab.IO.global_value_maxPeriLab.IO.global_value_minPeriLab.IO.global_value_sumPeriLab.IO.hex8_volumePeriLab.IO.hex8_volumePeriLab.IO.init_dataPeriLab.IO.init_orientationsPeriLab.IO.init_results_in_exodusPeriLab.IO.init_write_resultsPeriLab.IO.initialize_dataPeriLab.IO.load_and_evaluate_meshPeriLab.IO.local_nodes_from_dictPeriLab.IO.merge_exodus_filePeriLab.IO.merge_exodus_filesPeriLab.IO.movementPeriLab.IO.neighborsPeriLab.IO.node_distributionPeriLab.IO.open_result_filePeriLab.IO.paraview_specificsPeriLab.IO.parseLinePeriLab.IO.read_external_topologyPeriLab.IO.read_inputPeriLab.IO.read_input_filePeriLab.IO.read_meshPeriLab.IO.rectangular_plane_filterPeriLab.IO.set_dofPeriLab.IO.set_output_frequencyPeriLab.IO.show_block_summaryPeriLab.IO.show_mpi_summaryPeriLab.IO.stripCommentsPeriLab.IO.tetrahedron_volumePeriLab.IO.tetrahedron_volumePeriLab.IO.wedge6_volumePeriLab.IO.write_global_results_in_csvPeriLab.IO.write_global_results_in_exodusPeriLab.IO.write_nodal_results_in_exodusPeriLab.IO.write_resultsPeriLab.IO.write_step_and_timePeriLab.Parameter_Handling._get_valuesPeriLab.Parameter_Handling.check_for_duplicatesPeriLab.Parameter_Handling.get_all_keysPeriLab.Parameter_Handling.get_anglesPeriLab.Parameter_Handling.get_bc_definitionsPeriLab.Parameter_Handling.get_block_names_and_idsPeriLab.Parameter_Handling.get_bond_filtersPeriLab.Parameter_Handling.get_calculation_optionsPeriLab.Parameter_Handling.get_computesPeriLab.Parameter_Handling.get_computes_namesPeriLab.Parameter_Handling.get_densityPeriLab.Parameter_Handling.get_end_timePeriLab.Parameter_Handling.get_external_topology_namePeriLab.Parameter_Handling.get_fem_blockPeriLab.Parameter_Handling.get_final_timePeriLab.Parameter_Handling.get_fixed_dtPeriLab.Parameter_Handling.get_flush_filePeriLab.Parameter_Handling.get_headerPeriLab.Parameter_Handling.get_heat_capacityPeriLab.Parameter_Handling.get_horizonPeriLab.Parameter_Handling.get_initial_timePeriLab.Parameter_Handling.get_max_damagePeriLab.Parameter_Handling.get_mesh_namePeriLab.Parameter_Handling.get_model_optionsPeriLab.Parameter_Handling.get_model_parameterPeriLab.Parameter_Handling.get_node_setsPeriLab.Parameter_Handling.get_nstepsPeriLab.Parameter_Handling.get_numerical_dampingPeriLab.Parameter_Handling.get_output_fieldnamesPeriLab.Parameter_Handling.get_output_filenamesPeriLab.Parameter_Handling.get_output_frequenciesPeriLab.Parameter_Handling.get_output_typePeriLab.Parameter_Handling.get_output_variablesPeriLab.Parameter_Handling.get_outputsPeriLab.Parameter_Handling.get_safety_factorPeriLab.Parameter_Handling.get_solver_namePeriLab.Parameter_Handling.get_solver_paramsPeriLab.Parameter_Handling.get_solver_stepsPeriLab.Parameter_Handling.get_start_timePeriLab.Parameter_Handling.get_write_after_damagePeriLab.Parameter_Handling.validate_structure_recursivePeriLab.Parameter_Handling.validate_yaml
IO
PeriLab.IO._init_overlap_map_ — Method
_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.
PeriLab.IO.apply_bond_filters — Method
apply_bond_filters(nlist::BondScalarState{Int64}, mesh::DataFrame, params::Dict, dof::Int64)Apply the bond filters to the neighborhood list.
Arguments
nlist::BondScalarState{Int64}: The neighborhood list.mesh::DataFrame: The mesh.params::Dict: The parameters.dof::Int64: The degrees of freedom.
Returns
nlist::BondScalarState{Int64}: The filtered neighborhood list.nlist_filtered_ids::BondScalarState{Int64}: The filtered neighborhood list.
PeriLab.IO.area_of_polygon — Method
area_of_polygon(vertices)Calculate the area of a polygon.
Arguments
vertices: The vertices of the polygon.
Returns
area: The area of the polygon.
PeriLab.IO.bond_intersect_infinite_plane — Method
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.
PeriLab.IO.bond_intersect_rectangle_plane — Method
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.
PeriLab.IO.bond_intersects_disc — Method
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.
PeriLab.IO.calculate_block — Method
calculate_block( field_key::String, dof::Int64, calculation_type::String, block::Int64)Calculate the global value of a field for a given block.
Arguments
field_key::String: Field key.dof::Union{Int64,Vector{Int64}}: Degree of freedomcalculation_type::String: Calculation type.block::Int64: Block number.
Returns
value::Float64: Global value.nnodes::Int64: Number of nodes.
PeriLab.IO.calculate_nodelist — Method
calculate_nodelist(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
field_key::String: Field key.dof::Union{Int64,Vector{Int64}}: Degree of freedomcalculation_type::String: Calculation type.local_nodes::Vector{Int64}: Node set.
Returns
value::Vector: Global value.nnodes::Int64: Number of nodes.
PeriLab.IO.calculate_volume — Method
calculate_volume(element_type::String, vertices::Vector{Vector{Float64}})Calculate the volume of a element.
Arguments
element_type: The element type of the element.vertices: The vertices of the element.
Returns
volume: The volume of the element.
PeriLab.IO.check_for_duplicate_in_dataframe — Method
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.
PeriLab.IO.check_mesh_elements — Method
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)
PeriLab.IO.check_types_in_dataframe — Method
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.
PeriLab.IO.clearNP1 — Method
clearNP1(name::String)Clears the NP1 from the name
Arguments
name::String: The name
Returns
name::String: The cleared name
PeriLab.IO.close_result_file — Method
close_result_file(result_file::Dict)Closes the result file
Arguments
result_file::Dict: The result file
PeriLab.IO.close_result_files — Method
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 filesoutputs::Dict{Int64,Dict{}}: The output settings
PeriLab.IO.close_result_files — Method
close_result_files(result_files::Vector{Dict})Closes the result files
Arguments
result_files::Vector{Dict}: The result files
Returns
true: File is closedfalse: File was already closed
PeriLab.IO.create_and_distribute_bond_norm — Method
create_and_distribute_bond_norm(comm::MPI.Comm, nlist_filtered_ids::BondScalarState{Int64}, distribution::Vector{Int64}, bond_norm::Vector{Float64}, dof::Int64)Create and distribute the bond norm
Arguments
comm::MPI.Comm: MPI communicatornlist_filtered_ids::BondScalarState{Int64}: The filtered neighborhood listdistribution::Vector{Int64}: The distributionbond_norm::Vector{Float64}: The bond normdof::Int64: The degree of freedom
PeriLab.IO.create_distribution — Method
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.
PeriLab.IO.create_distribution_neighbor_based — Method
create_distribution_neighbor_based(nnodes::Int64,nlist::BondScalarState{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::BondScalarState{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.
PeriLab.IO.create_distribution_node_based — Method
create_distribution_node_based(nnodes::Int64,nlist::BondScalarState{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::BondScalarState{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.
PeriLab.IO.create_global_to_local_mapping — Method
glob_to_loc(distribution)Get the global to local mapping
Arguments
distribution: The distribution
Returns
create_global_to_local_mapping: The global to local mapping
PeriLab.IO.create_neighborhoodlist — Method
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.
PeriLab.IO.create_overlap_map — Method
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.
PeriLab.IO.create_result_file — Function
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 createnum_nodes::Int64: The number of nodesnum_dim::Int64: The number of dimensionsnum_elem_blks::Int64: The number of element blocksnum_node_sets::Int64: The number of node sets
Returns
result_file::Dict{String,Any}: A dictionary containing the filename and the exodus file
PeriLab.IO.create_result_file — Method
create_result_file(filename::String, outputs::Dict)Creates a csv file for the results
Arguments
filename::String: The name of the file to createoutputs::Dict: The outputs dictionary
Returns
Dict: The result file
PeriLab.IO.csv_reader — Method
csv_reader(filename::String)Read csv and return it as a DataFrame.
Arguments
filename::String: The path to the mesh file.
Returns
csvData::DataFrame: The csv data a DataFrame.
PeriLab.IO.define_nsets — Method
define_nsets(nsets::Dict{String,Vector{Int64}})Defines the node sets
Arguments
nsets::Dict{String,Vector{Int64}}: Node sets read from files
PeriLab.IO.delete_files — Method
delete_files(result_files::Vector{Dict})Deletes the result files
Arguments
result_files: The result files
PeriLab.IO.disk_filter — Method
disk_filter(nnodes::Int64, data::Matrix{Float64}, filter::Dict, nlist::BondScalarState{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::BondScalarState{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.
PeriLab.IO.distribute_neighborhoodlist_to_cores — Method
distribute_neighborhoodlist_to_cores(comm::MPI.Comm, nlist, distribution)Distributes the neighborhood list to the cores.
Arguments
comm::MPI.Comm: MPI communicatornlist: neighborhood listdistribution Array{Int64}: global nodes distribution at cores
PeriLab.IO.distribution_to_cores — Method
distribution_to_cores(comm::MPI.Comm, mesh, distribution, dof::Int64)Distributes the mesh data to the cores
Arguments
comm::MPI.Comm: MPI communicatormesh: Meshdistribution Array{Int64}: global nodes distribution at coresdof::Int64: Degree of freedom
PeriLab.IO.element_distribution — Method
element_distribution(topology::Vector{Vector{Int64}}, ptc::Vector{Int64}, size::Int64)Create the distribution of the finite elements. Is needed to avoid multiple element calls. Each element should run only one time at the cores.
Arguments
topology::Vector{Vector{Int64}}: The topology list of the mesh elements.nlist::BondScalarState{Int64}: The neighborhood list of the mesh elements.size::Int64: The number of ranks.
Returns
distribution::Vector{Vector{Int64}}: The distribution of the nodes.etc::Vector{Int64}: The number of nodes in each rank.
PeriLab.IO.extrude — Method
extrude(cmds, dataobject)Example extrusion callback for G1 which calculates total length of filament extruded.
The extruded filament length is obtained by watching the E axis movement in the g-code file.
PeriLab.IO.extrude_surface_mesh — Method
extrude_surface_mesh(mesh::DataFrame)extrude the mesh at the surface of the block
Arguments
mesh::DataFrame: The input mesh data represented as a DataFrame.params::Dict: The input parameters.
PeriLab.IO.find_global_core_value! — Method
find_global_core_value!(global_value::Union{Int64,Float64}, calculation_type::String, nnodes::Int64)Find global core value.
Arguments
global_value::Union{Int64,Float64}: The global valuecalculation_type::String: The calculation typennodes::Int64: The number of nodes
Returns
global_value::Union{Int64,Float64}: The global value
PeriLab.IO.get_block_nodes — Method
get_block_nodes(block_Id::AbstractVector{Int64}, block::Int64)Returns the nodes of a block
Arguments
block_Id::AbstractVector{Int64}: The block Idblock::Int64: The block
Returns
nodes::Vector{Int64}: The nodes of the block
PeriLab.IO.get_bond_geometry — Method
get_bond_geometry()Gets the bond geometry
PeriLab.IO.get_file_size — Method
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
PeriLab.IO.get_global_values — Method
get_global_values(output::Dict,)Get global values.
Arguments
output::Dict: The output
Returns
global_values::Vector: The global values
PeriLab.IO.get_local_element_topology — Method
get_local_element_topology(topology::Vector{Vector{Int64}}, distribution::Vector{Int64})Get the local element topology
Arguments
topology::Vector{Vector{Int64}}: The topologydistribution::Vector{Int64}: The distribution
PeriLab.IO.get_local_neighbors — Method
get_local_neighbors(mapping, nlist_core)Gets the local neighborhood list from the global neighborhood list
Arguments
mapping: mapping functionnlist_core: global neighborhood list
Returns
nlist_core: local neighborhood list
PeriLab.IO.get_local_overlap_map — Method
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 mappingranks 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 nodesPeriLab.IO.get_mpi_rank_string — Method
get_mpi_rank_string(rank::Int64, max_rank::Int64)Get MPI rank string.
Arguments
value::Int64: The rankmax_rank::Int64: The max rank
Returns
result::String: The result
PeriLab.IO.get_number_of_neighbornodes — Method
get_number_of_neighbornodes(nlist::BondScalarState{Int64})Get the number of neighbors for each node.
Arguments
nlist::BondScalarState{Int64}: The neighborhood list of the mesh elements.
Returns
length_nlist::Vector{Int64}: The number of neighbors for each node.
PeriLab.IO.get_paraview_coordinates — Method
get_paraview_coordinates(dof::Int64, refDof::Int64)Returns the paraview specific dof
Arguments
dof::Int64: The degrees of freedomrefDof::Int64: The reference degrees of freedom
Returns
paraview_specifics::String: The paraview specific dof
PeriLab.IO.get_results_mapping — Method
get_results_mapping(params::Dict, path::String)Gets the results mapping
Arguments
params::Dict: The parameterspath::String: The path
Returns
output_mapping::Dict{Int64,Dict{}}: The results mapping
PeriLab.IO.global_value_avg — Method
global_value_avg(field::Union{NodeScalarField{Float64},NodeVectorField{Float64}}, dof::Union{Int64,Vector{Int64}}, nodes::AbstractVector{Int64})Calculate the global average of a field for given nodes.
Arguments
field::Union{NodeScalarField{Float64},NodeVectorField{Float64}}: Field.dof::Union{Int64,Vector{Int64}}: Degree of freedomnodes::AbstractVector{Int64}: Nodes.
Returns
returnValue::Vector: Global value.
PeriLab.IO.global_value_max — Method
global_value_max(field::Union{NodeScalarField{Float64},NodeVectorField{Float64}}, dof::Union{Int64,Vector{Int64}}, nodes::AbstractVector{Int64})Calculate the global maximum of a field for given nodes.
Arguments
field::Union{NodeScalarField{Float64},NodeVectorField{Float64}}: Field.dof::Union{Int64,Vector{Int64}}: Degree of freedomnodes::AbstractVector{Int64}: Nodes.
Returns
returnValue::Vector: Global value.
PeriLab.IO.global_value_min — Method
global_value_min(field::Union{NodeScalarField{Float64},NodeVectorField{Float64}}, dof::Union{Int64,Vector{Int64}}, nodes::AbstractVector{Int64})Calculate the global minimum of a field for given nodes.
Arguments
field::Union{NodeScalarField{Float64},NodeVectorField{Float64}}: Field.dof::Union{Int64,Vector{Int64}}: Degree of freedomnodes::AbstractVector{Int64}: Nodes.
Returns
returnValue::Vector: Global value.
PeriLab.IO.global_value_sum — Method
global_value_sum(field::Union{NodeScalarField{Float64},NodeVectorField{Float64},NodeTensorField{Float64,3}}, dof::Union{Int64,Vector{Int64}}, nodes::AbstractVector{Int64})Calculate the global sum of a field for given nodes.
Arguments
field::Union{NodeScalarField{Float64},NodeVectorField{Float64},NodeTensorField{Float64,3}}: Field.dof::Union{Int64,Vector{Int64}: Degree of freedomnodes::AbstractVector{Int64}: Nodes.
Returns
returnValue::Vector: Global value.
PeriLab.IO.hex8_volume — Method
hex8volume(hexvertices)
Calculate the volume of a hex.
Arguments
hex_vertices: The vertices of the wedge.
Returns
volume: The volume of the wedge.
PeriLab.IO.hex8_volume — Method
hex8volume(hexvertices)
Calculate the volume of a hex.
Arguments
hex_vertices: The vertices of the wedge.
Returns
volume: The volume of the wedge.
PeriLab.IO.init_data — Method
init_data(params::Dict, path::String, comm::MPI.Comm)Initializes the data for the mesh.
Arguments
params::Dict: The parameters for the simulation.path::String: The path to the mesh file.comm::MPI.Comm: The MPI communicator.
Returns
params::Dict: The parameters for the simulation.
PeriLab.IO.init_orientations — Method
init_orientations()Initialize orientations.
PeriLab.IO.init_results_in_exodus — Function
init_results_in_exodus(exo::ExodusDatabase, output::Dict{}, coords::Union{Matrix{Int64},Matrix{Float64}}, block_Id::Vector{Int64}, block_list::Vector{String}, nsets::Dict{String,Vector{Int64}}, global_ids::Vector{Int64}, PERILAB_VERSION::String)Initializes the results in exodus
Arguments
exo::ExodusDatabase: The exodus databaseoutput::Dict{String,Any}: The outputcoords::Union{Matrix{Int64},Matrix{Float64}}: The coordinatesblock_Id::Vector{Int64}: The block Idblock_list::Vector{String}: The unique blocksnsets::Dict{String,Vector{Int64}}: The node setsglobal_ids::Vector{Int64}: The global ids
Returns
result_file::Dict{String,Any}: The result file
PeriLab.IO.init_write_results — Method
init_write_results(params::Dict, output_dir::String, path::String, nsteps::Int64, PERILAB_VERSION::String)Initialize write results.
Arguments
params::Dict: The parametersoutput_dir::String: The output directory.path::String: The pathnsteps::Int64: The number of steps
Returns
result_files::Array: The result filesoutputs::Dict: The outputs
PeriLab.IO.initialize_data — Method
initialize_data(filename::String, filedirectory::String, comm::MPI.Comm)Initialize data.
Arguments
filename::String: The name of the input file.filedirectory::String: The directory of the input file.comm::MPI.Comm: The MPI communicator
Returns
data::Dict: The data
PeriLab.IO.load_and_evaluate_mesh — Method
load_and_evaluate_mesh(params::Dict, path::String, ranksize::Int64)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.
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 setstopology::Int64::Array{Int64,nelement:nodes}`: The topology of elements.el_distribution::Array{Int64,1}: The distribution of the finite elements.
PeriLab.IO.local_nodes_from_dict — Method
local_nodes_from_dict(create_global_to_local_mapping::Dict{Int,Int}, global_nodes::Vector{Int64})Changes entries in the overlap map from the global numbering to the local computer core one.
Arguments
create_global_to_local_mapping::Dict{Int,Int}: global to local mappingglobal_nodes::Vector{Int64}: global nodes
Returns
overlap_map::Dict{Int64, Dict{Int64, String}}: returns overlap map with local nodes.
PeriLab.IO.merge_exodus_file — Method
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
PeriLab.IO.merge_exodus_files — Method
merge_exodus_files(result_files::Vector{Any}, output_dir::String)Merges exodus output files
Arguments
result_files::Vector{Any}: The result filesoutput_dir::String: The file directory
PeriLab.IO.movement — Method
movement(cmds, dataobject)Example movement callback for G0 and G1 which calculates the total distance moved in all axes.
It is calculated by watching the X, Y and Z axes movement.
PeriLab.IO.neighbors — Method
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.
PeriLab.IO.node_distribution — Function
node_distribution(nlist::BondScalarState{Int64}, size::Int64)Create the distribution of the nodes.
Arguments
nlist::BondScalarState{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.
PeriLab.IO.open_result_file — Method
open_result_file(result_file::Dict)Opens the result file
Arguments
result_file::Dict: The result file
PeriLab.IO.paraview_specifics — Method
paraview_specifics(dof::Int64)Returns the paraview specific dof
Arguments
dof::Int64: The degrees of freedom
Returns
paraview_specifics::String: The paraview specific dof
PeriLab.IO.parseLine — Function
parseLine(line::String, returnPair::Bool = true)::Array{Union{String,Pair{String,String}},1}Parse a single line of g-code and return an array of Pair{String,String} or an array of String containing the parsed commands.
The first command usually defines what to do (ie. G01 - linear interpolation) and following commands are the arguments (ie. X 14.312);
Examples
julia> parseLine("G10 X5.Y3. E6.")
4-element Array{Union{Pair{String,String}, String},1}:
"G" => "10"
"X" => "5."
"Y" => "3."
"E" => "6."Return array of strings
julia> parseLine("G10 X5.Y3. E6.", false)
4-element Array{Union{Pair{String,String}, String},1}:
"G10"
"X5."
"Y3."
"E6."PeriLab.IO.read_external_topology — Method
read_external_topology(filename::String)Read external topoloy data from a file and return it as a DataFrame.
Arguments
filename::String: The path to the mesh file.
Returns
external_topology::DataFrame: The external topology data as a DataFrame.
PeriLab.IO.read_input — Method
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
PeriLab.IO.read_input_file — Method
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.
PeriLab.IO.read_mesh — Method
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.
PeriLab.IO.rectangular_plane_filter — Method
rectangular_plane_filter(nnodes::Int64, data::Matrix{Float64}, filter::Dict, nlist::BondScalarState{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::BondScalarState{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.
PeriLab.IO.set_dof — Method
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.
PeriLab.IO.set_output_frequency — Function
set_output_frequency(params::Dict, nsteps::Int64, step_id::Int64)Sets the output frequency.
Arguments
params::Dict: The parametersnsteps::Int64: The number of stepsstep_id::Int64: The step id
PeriLab.IO.show_block_summary — Method
show_block_summary(solver_options::Dict, params::Dict, log_file::String, silent::Bool, comm::MPI.Comm)Show block summary.
Arguments
solver_options::Dict: The solver optionsparams::Dict: The paramslog_file::String: The log filesilent::Bool: The silent flagcomm::MPI.Comm: The Comm_rank
PeriLab.IO.show_mpi_summary — Method
show_mpi_summary(log_file::String, silent::Bool, comm::MPI.Comm)Show MPI summary.
Arguments
log_file::String: The log filesilent::Bool: The silent flagcomm::MPI.Comm: The comm
PeriLab.IO.stripComments — Method
stripComments(line::String)::StringReturn a copy of string line with stripped comments inside parentheses and all characters after a semicolon.
This function also removes whitespace as it it not needed for further parsing.
Examples
julia> stripComments("G92 (G10(aaa)))) ((comment)G) Z0.2 ; this is a comment")
"G92Z0.2"PeriLab.IO.tetrahedron_volume — Method
tetrahedron_volume(tet_vertices)Calculate the volume of a tetrahedron.
Arguments
tet_vertices: The vertices of the tetrahedron.
Returns
volume: The volume of the tetrahedron.
PeriLab.IO.tetrahedron_volume — Method
tetrahedron_volume(tet_vertices)Calculate the volume of a tetrahedron.
Arguments
tet_vertices: The vertices of the tetrahedron.
Returns
volume: The volume of the tetrahedron.
PeriLab.IO.wedge6_volume — Method
wedge6_volume(wedge_vertices)Calculate the volume of a wedge.
Arguments
wedge_vertices: The vertices of the wedge.
Returns
volume: The volume of the wedge.
PeriLab.IO.write_global_results_in_csv — Method
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 fileglobal_values: The global values
PeriLab.IO.write_global_results_in_exodus — Method
write_global_results_in_exodus(exo::ExodusDatabase, step::Int64, global_values)Writes the global results in the exodus file
Arguments
exo::ExodusDatabase: The exodus filestep::Int64: The stepglobal_values: The global values
Returns
exo::ExodusDatabase: The exodus file
PeriLab.IO.write_nodal_results_in_exodus — Method
write_nodal_results_in_exodus(exo::ExodusDatabase, step::Int64, output::Dict)Writes the nodal results in the exodus file
Arguments
exo::ExodusDatabase: The exodus filestep::Int64: The stepoutput::Dict: The output
Returns
exo::ExodusDatabase: The exodus file
PeriLab.IO.write_results — Method
write_results(result_files::Vector{Any}, time::Float64, max_damage::Float64, outputs::Dict)Write results.
Arguments
result_files::Vector{Any}: The result filestime::Float64: The timemax_damage::Float64: The maximum damageoutputs::Dict: The outputs
Returns
result_files::Vector{Any}: The result files
PeriLab.IO.write_step_and_time — Method
write_step_and_time(exo::ExodusDatabase, step::Int64, time::Float64)Writes the step and time in the exodus file
Arguments
exo::ExodusDatabase: The exodus filestep::Int64: The steptime::Float64: The time
Returns
exo::ExodusDatabase: The exodus file
Parameter_Handling
PeriLab.Parameter_Handling._get_values — Function
_get_values(params::Dict, block_id::Int64, valueName::String, defaultValue::Union{Float64,Bool,Nothing})Get the value of a block.
Arguments
params::Dict: The parametersblock_id::Int64: The ID of the blockvalueName::String: The name of the valuedefaultValue::Union{Float64,Bool,Nothing: The default value
Returns
value::Float64: The value of the block
PeriLab.Parameter_Handling.check_for_duplicates — Method
check_for_duplicates(filenames)Check for duplicate filenames.
Arguments
filenames::Vector{String}: The filenames
PeriLab.Parameter_Handling.get_all_keys — Method
get_all_keys(params::Dict)Get all the keys in the parameters
Arguments
params::Dict: The parameters dictionary.
Returns
keys_list::Array: The keys list
PeriLab.Parameter_Handling.get_angles — Method
get_angles(params::Dict, block_id::Int64, dof::Int64)Get the horizon of a block.
Arguments
params::Dict: The parametersblock_id::Int64: The ID of the blockdof::Int64: The dof
Returns
angles::Float64: The angles of the block
PeriLab.Parameter_Handling.get_bc_definitions — Method
get_bc_definitions(params::Dict)Get the boundary condition definitions
Arguments
params::Dict: The parameters
Returns
bcs::Dict{String,Any}: The boundary conditions
PeriLab.Parameter_Handling.get_block_names_and_ids — Method
get_block_names_and_ids(params::Dict, block_ids::Vector{Int64})Get the names of the blocks.
Arguments
params::Dict: The parameters dictionary.block_ids::Vector{Int64}: The IDs of the blocks
Returns
block_names::Vector{String}: The names of the blocks.
PeriLab.Parameter_Handling.get_bond_filters — Method
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 definedbfList::Dict{String,Dict{String,Any}}: The bond filters
PeriLab.Parameter_Handling.get_calculation_options — Method
get_calculation_options(params::Dict)Get the calculation options
Arguments
params::Dict: The parameters dictionary.
Returns
solver_options::Dict: The solver options
PeriLab.Parameter_Handling.get_computes — Method
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.
PeriLab.Parameter_Handling.get_computes_names — Method
get_computes_names(params::Dict)Get the names of the computes.
Arguments
params::Dict: The parameters dictionary.
Returns
computes_names::Vector{String}: The names of the computes.
PeriLab.Parameter_Handling.get_density — Method
get_density(params::Dict, block_id::Int64)Get the density of a block.
Arguments
params::Dict: The parametersblock_id::Int64: The ID of the block
Returns
density::Float64: The density of the block
PeriLab.Parameter_Handling.get_end_time — Method
get_end_time(outputs::Dict, output::String)Get the end_time.
Arguments
outputs::Dict: The outputsoutput::String: The output
Returns
end_time::Float64: The value
PeriLab.Parameter_Handling.get_external_topology_name — Method
get_external_topology_name(params::Dict, path)Returns the name of the mesh file from the parameters
Arguments
params::Dict: The parameterspath::String: Path of the working folder
Returns
String: The name of the finite element topology file
PeriLab.Parameter_Handling.get_fem_block — Method
get_fem_block(params::Dict, block_id::Int64)Get the fem_block of a block.
Arguments
params::Dict: The parametersblock_id::Int64: The ID of the block
Returns
fem_block::Float64: The fem_block of the block
PeriLab.Parameter_Handling.get_final_time — Method
get_final_time(params::Dict)Get the final time
Arguments
params::Dict: The parameters dictionary.
Returns
final_time::Float64: The final time
PeriLab.Parameter_Handling.get_fixed_dt — Method
get_fixed_dt(params::Dict)Get the fixed time step
Arguments
params::Dict: The parameters dictionary.
Returns
fixed_dt::Float64: The fixed time step
PeriLab.Parameter_Handling.get_flush_file — Method
get_flush_file(outputs::Dict, output::String)Gets the flush file.
Arguments
outputs::Dict: The outputsoutput::String: The output
Returns
flush_file::Bool: The flush file
PeriLab.Parameter_Handling.get_header — Method
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.
PeriLab.Parameter_Handling.get_heat_capacity — Method
get_heat_capacity(params::Dict, block_id::Int64)Get the heat capacity of a block.
Arguments
params::Dict: The parametersblock_id::Int64: The ID of the block
Returns
heat_capacity::Float64: The heat capacity of the block
PeriLab.Parameter_Handling.get_horizon — Method
get_horizon(params::Dict, block_id::Int64)Get the horizon of a block.
Arguments
params::Dict: The parametersblock_id::Int64: The ID of the block
Returns
horizon::Float64: The horizon of the block
PeriLab.Parameter_Handling.get_initial_time — Method
get_initial_time(params::Dict)Get the initial time
Arguments
params::Dict: The parameters dictionary.
Returns
initial_time::Float64: The initial time
PeriLab.Parameter_Handling.get_max_damage — Method
getmaxdamage(params::Dict)
Get the maximum damage.
Arguments
params::Dict: The parameters dictionary.
Returns
max_damage::Float64: The value
PeriLab.Parameter_Handling.get_mesh_name — Method
get_mesh_name(params::Dict)Returns the name of the mesh file from the parameters
Arguments
params::Dict: The parameters
Returns
String: The name of the mesh file
PeriLab.Parameter_Handling.get_model_options — Method
get_model_options(params::Dict)Get the solver options
Arguments
params::Dict: The parameters dictionary.
Returns
solver_options::Dict: The solver options
PeriLab.Parameter_Handling.get_model_parameter — Function
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, ornothingif 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
params = Dict(
"Models" => Dict(
"Models" => Dict(
"ModelA" => 42,
"ModelB" => 24
)
)
)
model = "Models"
id = "ModelA"
result = get_model_parameter(params, model, id)
if result !== nothing
println("Parameter id: result")
else
println("Parameter not found.")
endPeriLab.Parameter_Handling.get_node_sets — Method
get_node_sets(params::Dict, path::String)Returns the node sets from the parameters
Arguments
params::Dict: The parameterspath::String: The path to the mesh file
Returns
nsets::Dict{String,Any}: The node sets
PeriLab.Parameter_Handling.get_nsteps — Method
get_nsteps(params::Dict)Get the fixed time step
Arguments
params::Dict: The parameters dictionary.
Returns
nsteps::Int64: The fixed time step
PeriLab.Parameter_Handling.get_numerical_damping — Method
get_numerical_damping(params::Dict)Get the numerical damping
Arguments
params::Dict: The parameters dictionary.
Returns
numerical_damping::Float64: The numerical damping
PeriLab.Parameter_Handling.get_output_fieldnames — Method
get_output_fieldnames(outputs::Dict, variables::Vector{String}, computes::Vector{String}, output_type::String)Gets the output fieldnames.
Arguments
outputs::Dict: The outputsvariables::Vector{String}: The variablescomputes::Vector{String}: The computesoutput_type::String: The output type
Returns
output_fieldnames::Vector{String}: The output fieldnames
PeriLab.Parameter_Handling.get_output_filenames — Method
get_output_filenames(params::Dict, output_dir::String)Gets the output filenames.
Arguments
params::Dict: The parametersoutput_dir::String: The file directory
Returns
filenames::Vector{String}: The filenames
PeriLab.Parameter_Handling.get_output_frequencies — Method
get_output_frequencies(params::Dict, nsteps::Int64)Gets the output frequencies.
Arguments
params::Dict: The parametersnsteps::Int64: The number of steps
Returns
freq::Vector{Int64}: The output frequencies
PeriLab.Parameter_Handling.get_output_type — Method
get_output_type(outputs::Dict, output::String)Gets the output type.
Arguments
outputs::Dict: The outputsoutput::String: The output
Returns
output_type::String: The output type
PeriLab.Parameter_Handling.get_output_variables — Method
get_output_variables(output::String, variables::Vector{String})Get the output variable.
Arguments
output::String: The output variable.variables::Vector{String}: The variables.
Returns
output::String: The output variable.
PeriLab.Parameter_Handling.get_outputs — Method
get_outputs(params::Dict, variables::Vector{String}, compute_names::Vector{String})Gets the outputs.
Arguments
params::Dict: The parametersvariables::Vector{String}: The variablescompute_names::Vector{String}: The compute names
Returns
outputs::Dict: The outputs
PeriLab.Parameter_Handling.get_safety_factor — Method
get_safety_factor(params::Dict)Get the safety factor
Arguments
params::Dict: The parameters dictionary.
Returns
safety_factor::Float64: The safety factor
PeriLab.Parameter_Handling.get_solver_name — Method
get_solver_name(params::Dict)Get the name of the solver
Arguments
params::Dict: The parameters dictionary.
Returns
solver_name::String: The name of the solver
PeriLab.Parameter_Handling.get_solver_params — Method
get_solver_params(params::Dict)Get the solver parameters
Arguments
params::Dict: The parameters dictionary.
Returns
solver_params::Dict: The solver parameters
PeriLab.Parameter_Handling.get_solver_steps — Method
get_solver_steps(params::Dict)Get the solver steps
Arguments
params::Dict: The parameters dictionary.
Returns
solver_steps::List: The solver steps
PeriLab.Parameter_Handling.get_start_time — Method
get_start_time(outputs::Dict, output::String)Get the start_time.
Arguments
outputs::Dict: The outputsoutput::String: The output
Returns
start_time::Float64: The value
PeriLab.Parameter_Handling.get_write_after_damage — Method
get_write_after_damage(outputs::Dict, output::String)Get the write after damage.
Arguments
outputs::Dict: The outputsoutput::String: The output
Returns
write_after_damage::Bool: The value
PeriLab.Parameter_Handling.validate_structure_recursive — Function
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 structureactual::Dict: The actual structurevalidate::Bool: The validation resultschecked_keys::Array: The keys that have been checkedpath::String: The current path
Returns
validate::Bool: The validation resultchecked_keys::Array: The keys that have been checked
PeriLab.Parameter_Handling.validate_yaml — Method
validate_yaml(params::Dict)Validates the parameters against the expected structure
Arguments
params::Dict: The parameters dictionary.
Returns
params::Dict: The parameters dictionary.