tem_read_stlb Subroutine

public subroutine tem_read_stlb(filename, nNodesRead, nTrisRead, nodes, tri_node, iError)

This routine actually reads the data (points, triangles, normals) from the binary file and stores them.

Arguments

Type IntentOptional Attributes Name
character(len=PathLen), intent(in) :: filename

name of the binary stl file

integer, intent(in) :: nNodesRead

Number of nodes read from the file header, to compare against the actual number of nodes read

integer, intent(in) :: nTrisRead

Number of triangles read from the file header, to compare against the actual number of triangles read

real(kind=rk), intent(out) :: nodes(:,:)

point coordinates read from the stl-file size: 3, nPoints_total

integer, intent(out) :: tri_node(:,:)

connectivity array for the triangles size: 3, nTriangles_total

integer, intent(out) :: iError

error while openeing the file, or if the number of nodes/trias do not match to the ones read from the header (if error -> iError > 0)


Calls

proc~~tem_read_stlb~~CallsGraph proc~tem_read_stlb tem_read_stlb proc~tem_open tem_open proc~tem_read_stlb->proc~tem_open proc~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower proc~tem_abort tem_abort proc~tem_open->proc~tem_abort proc~newunit newunit proc~tem_open->proc~newunit mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~tem_read_stlb~~CalledByGraph proc~tem_read_stlb tem_read_stlb proc~tem_read_stlfiles tem_read_stlFiles proc~tem_read_stlfiles->proc~tem_read_stlb proc~tem_readandunify_surfdata tem_readAndUnify_surfData proc~tem_readandunify_surfdata->proc~tem_read_stlb 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 proc~tem_load_shapes tem_load_shapes proc~tem_load_shapes->proc~tem_load_shape_single interface~tem_load_shape tem_load_shape interface~tem_load_shape->proc~tem_load_shape_single

Contents

Source Code


Source Code

  subroutine tem_read_stlb( filename, nNodesRead, nTrisRead, nodes, tri_node,  &
    &                       ierror )
    ! --------------------------------------------------------------------------
    !> name of the binary stl file
    character(len=PathLen),intent(in) :: filename
    !> point coordinates read from the stl-file
    !! size: 3, nPoints_total
    real(kind=rk),intent(out) :: nodes(:,:)
    !> connectivity array for the triangles
    !! size: 3, nTriangles_total
    integer,intent(out) :: tri_node(:,:)
    !> Number of nodes read from the file header, to compare against the actual
    !! number of nodes read
    integer,intent(in)  :: nNodesRead
    !> Number of triangles read from the file header, to compare against the
    !! actual number of triangles read
    integer,intent(in)  :: nTrisRead
    !> error while openeing the file, or if the number of nodes/trias do not
    !! match to the ones read from the header (if error -> iError > 0)
    integer,intent(out) :: iError
    ! --------------------------------------------------------------------------
    ! buffer for reading the data from file
    real(kind=single_k)       :: temp(3) ! has to be single_k
    character(len=80) :: header ! has to be of length 80
    character(len=2) :: attribute
    integer :: nTriangles ! has to be 4 byte integer
    integer :: stlUnit
    integer :: i
    integer :: nNodes
    integer :: nTris
    ! --------------------------------------------------------------------------

    call tem_open( newunit = stlUnit,        &
      &            file    = trim(filename), &
      &            access  = 'stream',       &
      &            action  = 'read',         &
      &            status  = 'old'           )

    read(stlUnit) header
    read(stlUnit) nTriangles

    nNodes = 0
    nTris = 0
    do i=1,nTriangles
      ! assign node coordinates
      ! and nodes to tris
      nTris = nTris+1
      read(stlUnit) temp(:) ! normal vector

      read(stlUnit) temp(:) ! node 1 coords
      nNodes = nNodes+1
      nodes(1:3,nNodes) = real(temp(1:3),kind=rk)
      tri_node(1,nTris) = nNodes
      read(stlUnit) temp(:) ! node 2 coords
      nNodes = nNodes+1
      nodes(1:3,nNodes) = real(temp(1:3),kind=rk)
      tri_node(2,nTris) = nNodes
      read(stlUnit) temp(:) ! node 3 coords
      nNodes = nNodes+1
      nodes(1:3,nNodes) = real(temp(1:3),kind=rk)
      tri_node(3,nTris) = nNodes
      read(stlUnit) attribute
    end do

    close(stlUnit)
    if ( nTris .ne. nTrisRead) then
      write(logUnit(0),*) "Inconsistency found in number of triangles!"
      iError = 1
    endif
    if ( nNodes .ne. nNodesRead) then
      write(logUnit(0),*) "Inconsistency found in number of nodes!"
      iError = 1
    endif

  end subroutine tem_read_stlb