cmhamel / Exodus.jl

A julia interface for accessing the ExodusII data format
MIT License
9 stars 0 forks source link

Merging (EPU) returns "There were at least 2 elements with duplicated ids detected" #156

Closed JTHesse closed 11 months ago

JTHesse commented 12 months ago

I'm trying to create multiple Exodus files and merge them afterwards. However, epu won't run correctly as the element ids are equal in both files.

I think the best solution is to add an id offset.

Here you can find the example i'm referring to:

using Exodus

# data to write
coords1 = [
    1.0 0.5 0.5 1.0 0.0 0.0 0.5 1.0 0.0
    1.0 1.0 0.5 0.5 1.0 0.5 0.0 0.0 0.0
]
coords2 = [
    2.0 1.5 1.5 2.0 1.0 1.0 1.5 2.0 1.0
    1.0 1.0 0.5 0.5 1.0 0.5 0.0 0.0 0.0
]

conn = [
    1 2 4 3
    2 5 3 6
    3 6 7 9
    4 3 8 7
]

# make some hack variables to write
v_nodal_1 = rand(9)
v_nodal_2 = rand(9)

v_elem_1 = rand(4)
v_elem_2 = rand(4)

# set the types
maps_int_type = Int32
ids_int_type = Int32
bulk_int_type = Int32
float_type = Float64

# initialization parameters
num_dim, num_nodes = size(coords1)
num_elems = size(conn, 2)
num_elem_blks = 1
num_side_sets = 0
num_node_sets = 0

# make init
init = Initialization{bulk_int_type}(
    num_dim, num_nodes, num_elems,
    num_elem_blks, num_side_sets, num_node_sets
)

# finally make empty exo database
exo1 = ExodusDatabase(
    "test_write.e.2.0", "w", init,
    maps_int_type, ids_int_type, bulk_int_type, float_type
)
exo2 = ExodusDatabase(
    "test_write.e.2.1", "w", init,
    maps_int_type, ids_int_type, bulk_int_type, float_type
)

# how to write coordinates
write_coordinates(exo1, coords1)
write_coordinates(exo2, coords2)
# how to write a block
write_block(exo1, 1, "QUAD4", conn)
write_block(exo2, 1, "QUAD4", conn)
# need at least one timestep to output variables
write_time(exo1, 1, 0.0)
write_time(exo2, 1, 0.0)
# write number of variables and their names
write_names(exo1, NodalVariable, ["v_nodal_1"])
write_names(exo2, NodalVariable, ["v_nodal_1"])
write_names(exo1, ElementVariable, ["v_elem_1"])
write_names(exo2, ElementVariable, ["v_elem_1"])
# write variable values the 1 is for the time step
write_values(exo1, NodalVariable, 1, "v_nodal_1", v_nodal_1)
write_values(exo2, NodalVariable, 1, "v_nodal_1", v_nodal_1)
# the first 1 is for the time step 
# and the second 1 is for the block number
write_values(exo1, ElementVariable, 1, 1, "v_elem_1", v_elem_1)
write_values(exo2, ElementVariable, 1, 1, "v_elem_1", v_elem_1)
# don't forget to close the exodusdatabase, it can get corrupted otherwise if you're writing
close(exo1)
close(exo2)

epu("test_write.e.2.0")
cmhamel commented 12 months ago

@JTHesse

I think the more appropriate tool to use here is ejoin. Unfortunately, I have not wrapped that executable in Exodus_jll yet. I can add it to the to-do list though.

JTHesse commented 11 months ago

I forgot to mention that we are utilizing MPI and therefore I thought epu is the right function for that.

cmhamel commented 11 months ago

Are you trying to have different exodus databases for each MPI rank?

If so, consider first making the full mesh and splitting it with the decomp method. Then write data on each exodus database and epu at the end. Decomp should handle the element offset for you automatically (I think?).

JTHesse commented 11 months ago

Yes exactly, unfortunately we don't always have access to a binary .g mesh file. I just had a look at our old Project which is utilizing the exodus format in a similar way. Based on that, the two functions _ex_put_node_nummap and _ex_put_elem_nummap should be exactly what we are looking for, maybe it is possible to add those functionalities.

cmhamel commented 11 months ago

Got it.

I'm almost positive I have that functionality supported in some capacity but it's likely not well tested or documented.

I'll take a look this evening and get back to you.

cmhamel commented 11 months ago

@JTHesse update to Exodus v0.10.0 and try the following code snippet which I added to the testing coverage

# data to write
  coords1 = [
    1.0 0.5 0.5 1.0 0.0 0.0 0.5 1.0 0.0
    1.0 1.0 0.5 0.5 1.0 0.5 0.0 0.0 0.0
  ]
  coords2 = [
    2.0 1.5 1.5 2.0 1.0 1.0 1.5 2.0 1.0
    1.0 1.0 0.5 0.5 1.0 0.5 0.0 0.0 0.0
  ]

  conn = [
    1 2 4 3
    2 5 3 6
    3 6 7 9
    4 3 8 7
  ]

  node_map_1 = Int32[1, 2, 3, 4, 5, 6, 7, 8, 9]
  node_map_2 = node_map_1 .+ Int32(1000)
  elem_map_1 = Int32[1, 2, 3, 4]
  elem_map_2 = Int32[1000, 1001, 1002, 1003]

  # make some hack variables to write
  v_nodal_1 = rand(9)
  v_nodal_2 = rand(9)

  v_elem_1 = rand(4)
  v_elem_2 = rand(4)

  # set the types
  maps_int_type = Int32
  ids_int_type = Int32
  bulk_int_type = Int32
  float_type = Float64

  # initialization parameters
  num_dim, num_nodes = size(coords1)
  num_elems = size(conn, 2)
  num_elem_blks = 1
  num_side_sets = 0
  num_node_sets = 0

  # make init
  init = Initialization{bulk_int_type}(
    num_dim, num_nodes, num_elems,
    num_elem_blks, num_side_sets, num_node_sets
  )

  # finally make empty exo database
  exo1 = ExodusDatabase(
    "test_write.e.2.0", "w", init,
    maps_int_type, ids_int_type, bulk_int_type, float_type
  )
  exo2 = ExodusDatabase(
    "test_write.e.2.1", "w", init,
    maps_int_type, ids_int_type, bulk_int_type, float_type
  )

  # how to write coordinates
  write_coordinates(exo1, coords1)
  write_coordinates(exo2, coords2)
  # how to write a block
  write_block(exo1, 1, "QUAD4", conn)
  write_block(exo2, 1, "QUAD4", conn)
  # write element id map
  write_id_map(exo1, NodeMap, node_map_1)
  write_id_map(exo2, NodeMap, node_map_2)
  write_id_map(exo1, ElementMap, elem_map_1)
  write_id_map(exo2, ElementMap, elem_map_2)
  # need at least one timestep to output variables
  write_time(exo1, 1, 0.0)
  write_time(exo2, 1, 0.0)
  # write number of variables and their names
  write_names(exo1, NodalVariable, ["v_nodal_1"])
  write_names(exo2, NodalVariable, ["v_nodal_1"])
  write_names(exo1, ElementVariable, ["v_elem_1"])
  write_names(exo2, ElementVariable, ["v_elem_1"])
  # write variable values the 1 is for the time step
  write_values(exo1, NodalVariable, 1, "v_nodal_1", v_nodal_1)
  write_values(exo2, NodalVariable, 1, "v_nodal_1", v_nodal_1)
  # the first 1 is for the time step 
  # and the second 1 is for the block number
  write_values(exo1, ElementVariable, 1, 1, "v_elem_1", v_elem_1)
  write_values(exo2, ElementVariable, 1, 1, "v_elem_1", v_elem_1)
  # don't forget to close the exodusdatabase, it can get corrupted otherwise if you're writing
  close(exo1)
  close(exo2)

  epu("test_write.e.2.0")

If this works for you please close the issue.

JTHesse commented 11 months ago

Works perfect, thank you very much!