tem_read_stlFiles Subroutine

private subroutine tem_read_stlFiles(stl_data)

Read in all STL files from disk which are specified in the config file

Arguments

Type IntentOptional Attributes Name
type(tem_stlData_type), intent(inout) :: stl_data

stl data of current spatial object


Calls

proc~~tem_read_stlfiles~~CallsGraph proc~tem_read_stlfiles tem_read_stlFiles proc~tem_read_stlb tem_read_stlb proc~tem_read_stlfiles->proc~tem_read_stlb proc~tem_size_stlb tem_size_stlb proc~tem_read_stlfiles->proc~tem_size_stlb stla_read stla_read proc~tem_read_stlfiles->stla_read stla_size stla_size proc~tem_read_stlfiles->stla_size proc~tem_open tem_open proc~tem_read_stlb->proc~tem_open proc~tem_abort tem_abort proc~tem_size_stlb->proc~tem_abort proc~tem_size_stlb->proc~tem_open mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~tem_open->proc~tem_abort proc~newunit newunit proc~tem_open->proc~newunit proc~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower

Called by

proc~~tem_read_stlfiles~~CalledByGraph proc~tem_read_stlfiles tem_read_stlFiles proc~tem_load_stl tem_load_stl proc~tem_load_stl->proc~tem_read_stlfiles proc~tem_load_shape_single tem_load_shape_single proc~tem_load_shape_single->proc~tem_load_stl interface~tem_load_shape tem_load_shape interface~tem_load_shape->proc~tem_load_shape_single proc~tem_load_shapes tem_load_shapes interface~tem_load_shape->proc~tem_load_shapes proc~tem_load_shapes->proc~tem_load_shape_single proc~load_spatial_parabol load_spatial_parabol proc~load_spatial_parabol->interface~tem_load_shape proc~tem_load_convergenceheader tem_load_convergenceHeader proc~tem_load_convergenceheader->interface~tem_load_shape proc~tem_load_spacetime_single tem_load_spacetime_single proc~tem_load_spacetime_single->interface~tem_load_shape proc~tem_load_trackingconfig tem_load_trackingConfig proc~tem_load_trackingconfig->interface~tem_load_shape

Source Code

  subroutine tem_read_stlFiles( stl_data )
    ! ---------------------------------------------------------------------------!
    !> stl data of current spatial object
    type( tem_stlData_type ), intent(inout) :: stl_data
    ! ---------------------------------------------------------------------------!
    !local variable
    integer :: nodeOffset !< Offset in the nodelist for multiple STLs
    integer :: triOffset !< Offset in the trilist for multiple STLs
    integer :: iFile, ierr
    integer :: nNodes
    integer :: nTris !< total number of triangles loaded
    ! ---------------------------------------------------------------------------!

    write(logunit(1),*) " Reading in STL Headers ..."
    ! Read headerfiles to determine number of nodes and tris
    do iFile = 1, size(stl_data%head)
      select case( stl_data%head(iFile)%fileformat )
      case( stl_ascii )
        call stla_size( trim(stl_data%head(iFile)%filename),  &
          &                   stl_data%head(iFile)%nSolids,   &
          &                   stl_data%head(iFile)%nNodes,    &
          &                   stl_data%head(iFile)%nTris,     &
          &                   stl_data%head(iFile)%nTexts)
      case( stl_bin )
        call tem_size_stlb( stl_data%head(iFile)%filename,                           &
          &                 stl_data%head(iFile)%nNodes,                             &
          &                 stl_data%head(iFile)%nTris)

      end select
    end do

    !Sum over all read in STL Files to generate one list of elements in the end
    nTris     = 0
    nNodes   = 0
    do iFile = 1, size(stl_data%head)
      nTris = nTris + stl_data%head(iFile)%nTris
      nNodes = nNodes + stl_data%head(iFile)%nNodes
    end do

    write(logunit(5),"(A,I0)") " Total number of triangles: ", nTris
    write(logunit(5),"(A,I0)") " Total number of nodes:     ", nNodes

    !allocate node and triangle arrays
    stl_data%nNodes = nNodes
    stl_data%nTris = nTris
    allocate(stl_data%nodes(1:3,nNodes))
    allocate(stl_data%tri_node(1:3,nTris))

    ! Read in node values from STL files to stl_data%nodes
    ! Store three nodes position of each triangle to stl_date%tri_node
    write(logunit(1),*) " Reading in STL Files ..."

    nodeOffset = 1
    triOffset = 1
    do iFile = 1, size(stl_data%head)
      select case(stl_data%head(iFile)%fileformat)
      case(stl_ascii)
        call stla_read(input_file_name = trim(stl_data%head(iFile)%filename), &
        &    node_num  = stl_data%head(iFile)%nNodes,  &
        &    face_num  = stl_data%head(iFile)%nTris, &
        &    node_xyz  = stl_data%nodes(1:3,nodeOffset:nodeOffset+stl_data%head(iFile)%nNodes-1), &
        &    face_node = stl_data%tri_node(1:3,triOffset:triOffset+stl_data%head(iFile)%nTris-1), &
        &    ierror    = ierr)

      case(stl_bin)
        call tem_read_stlb(filename = stl_data%head(iFile)%filename, &
        &    nNodesRead = stl_data%head(iFile)%nNodes,  &
        &    nTrisRead  = stl_data%head(iFile)%nTris, &
        &    nodes  = stl_data%nodes(1:3,nodeOffset:nodeOffset+stl_data%head(iFile)%nNodes-1), &
        &    tri_node = stl_data%tri_node(1:3,triOffset:triOffset+stl_data%head(iFile)%nTris-1), &
        &    ierror    = ierr)
      end select

      stl_data%tri_node(1:3,triOffset:triOffset+stl_data%head(iFile)%nTris-1) =  &
      & stl_data%tri_node(1:3,triOffset:triOffset+stl_data%head(iFile)%nTris-1) + (nodeOffset-1)

      nodeOffset = nodeOffset + stl_data%head(iFile)%nNodes
      triOffset = triOffset + stl_data%head(iFile)%nTris

    end do

    write(logunit(1),*) " Done."

  end subroutine tem_read_stlFiles