tem_dump_stlb Subroutine

public subroutine tem_dump_stlb(outprefix, nodes, triangles, proc, header, normals, time)

This routine dumps a set of nodes and triangles to disc.

The nodes and their connectivity are passed to the routine. The normals are passed optional or calculated internally. The outputfile name is composed of the $outprefix,$time,'.stl'.

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: outprefix

output prefix for the filename

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

nodes to be dumped (size: 3*nNodes)

integer, intent(in) :: triangles(:,:)

triangles to be dumped (size: 3, nTrias)

type(tem_comm_env_type), intent(in) :: proc

process description to use

character(len=80), intent(in), optional :: header

optional header to be dumped

real(kind=rk), intent(in), optional :: normals(:,:)

optional array of normals, if not passed normals will be calculated internally

type(tem_time_type), intent(in), optional :: time

optional simulation time to be appended to the filename


Calls

proc~~tem_dump_stlb~~CallsGraph proc~tem_dump_stlb tem_dump_stlb mpi_reduce mpi_reduce proc~tem_dump_stlb->mpi_reduce proc~cross_product3d cross_product3D proc~tem_dump_stlb->proc~cross_product3d interface~append~23 append proc~tem_dump_stlb->interface~append~23 proc~tem_abort tem_abort proc~tem_dump_stlb->proc~tem_abort interface~init~22 init proc~tem_dump_stlb->interface~init~22 proc~tem_time_sim_stamp tem_time_sim_stamp proc~tem_dump_stlb->proc~tem_time_sim_stamp proc~tem_open tem_open proc~tem_dump_stlb->proc~tem_open interface~destroy~22 destroy proc~tem_dump_stlb->interface~destroy~22 proc~append_arrayga2d_real append_arrayga2d_real interface~append~23->proc~append_arrayga2d_real proc~append_singlega2d_real append_singlega2d_real interface~append~23->proc~append_singlega2d_real mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~init_ga2d_real init_ga2d_real interface~init~22->proc~init_ga2d_real proc~tem_open->proc~tem_abort proc~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower proc~newunit newunit proc~tem_open->proc~newunit proc~destroy_ga2d_real destroy_ga2d_real interface~destroy~22->proc~destroy_ga2d_real interface~expand~22 expand proc~append_arrayga2d_real->interface~expand~22 proc~append_singlega2d_real->interface~expand~22 proc~expand_ga2d_real expand_ga2d_real interface~expand~22->proc~expand_ga2d_real

Contents

Source Code


Source Code

  subroutine tem_dump_stlb( outprefix, nodes, triangles, proc, header,   &
    &                       normals, time )
    ! --------------------------------------------------------------------------
    !> output prefix for the filename
    character(len=*), intent(in) :: outprefix
    !> nodes to be dumped (size: 3*nNodes)
    real(kind=rk), intent(in) :: nodes(:)
    !> triangles to be dumped (size: 3, nTrias)
    integer, intent(in) :: triangles(:,:)
    !> process description to use
    type(tem_comm_env_type), intent(in) :: proc
    !> optional header to be dumped
    character(len=80), optional, intent(in) :: header
    !> optional array of normals, if not passed normals will be calculated
    !! internally
    real(kind=rk), optional, intent(in) :: normals(:,:)
    !> optional simulation time to be appended to the filename
    type(tem_time_type), optional, intent(in) :: time
    ! --------------------------------------------------------------------------
    real(kind=single_k) :: loc_normals( 3, size( triangles, 2 ))
    integer :: iTria
    ! temporary vectors to calculate the normals
    real(kind=rk) :: a(3), b(3)
    ! filename
    character(len=PathLen) :: filename
    ! output unit
    integer :: outUnit
    ! local header
    character(len=80) :: loc_header
    ! attribute to be dumped
    character(len=2) :: attribute
    ! error variable
    integer :: iError
    ! timestamp for the filename
    character(len=12) :: timeStamp
    ! temporary min and max position for X,Y,Z coordinates in the linearized
    ! array of nodes
    integer :: minPos1, maxPos1, minPos2, maxPos2, minPos3, maxPos3
    ! temporary array of node coordinates
    ! size = size(nodes)
    real(kind=rk), allocatable :: dump_nodes(:)
    ! number of entries in the nodes array
    integer :: nEntries
    real( kind=rk ) :: huge_real
    logical :: validTria
    ! temporary growing array of valid triangles to be dumped
    type(grw_int2darray_type) :: dump_trias
    integer :: iPoint
    ! --------------------------------------------------------------------------

    huge_real = huge(huge_real)

    nEntries = size( nodes )
    allocate( dump_nodes( nEntries ))

    ! first communicate all point coordinates to rank 0
    call mpi_reduce( nodes, dump_nodes, nEntries, rk_mpi, MPI_MIN, proc%root,  &
                     proc%comm, iError )

    ! now only root continues to calculate and dump information
    if( proc%rank .eq. proc%root )then

      if( present( normals ))then
        if( size( triangles, 2 ) .eq. size( normals, 2 ) )then
          loc_normals = real(normals, kind=single_k)
        else
          write(logUnit(0),*) " The number of triangles have to match the " &
            &              // "number of normals!!! This is not the case "  &
            &              // "(nTrias: ", size( triangles, 2 ),            &
            &                 " vs. nNormals: ", size( normals, 2 ),        &
            &                 " Stopping..."
          call tem_abort()
        end if
      else
        ! calculate the normals
        do iTria = 1, size( triangles, 2 )
          minPos1 = (triangles( 1, iTria )-1)*3+1
          maxPos1 = (triangles( 1, iTria )-1)*3+3
          minPos2 = (triangles( 2, iTria )-1)*3+1
          maxPos2 = (triangles( 2, iTria )-1)*3+3
          minPos3 = (triangles( 3, iTria )-1)*3+1
          maxPos3 = (triangles( 3, iTria )-1)*3+3
          a = dump_nodes( minPos2:maxPos2 ) - dump_nodes( minPos1:maxPos1 )
          b = dump_nodes( minPos3:maxPos3 ) - dump_nodes( minPos1:maxPos1 )
          loc_normals( :, iTria ) = real(cross_product3D( a, b ), kind=single_k)
        end do
      end if

      ! initialize the growing array of valid triangles to be dumped
      call init( me = dump_trias, width = 3 )

      ! now check the individual coordinates and remove all triangles that
      ! own a point coordinate which is not .lt. huge(rk)
      do iTria = 1, size( triangles, 2 )
        ! suppose the triangle is valid
        validTria = .true.
        point_loop: do iPoint = 1, 3
          minPos1 = (triangles( iPoint, iTria )-1)*3+1
          maxPos1 = (triangles( iPoint, iTria )-1)*3+3
          ! check if any coordinate for iPoint is not valid
          if( .not. all( dump_nodes( minPos1:maxPos1 ) .lt. huge_real ))then
            ! if this is the case set the validTria to false and
            ! exit the point loop
            validTria = .false.
            exit point_loop
          end if
        end do point_loop
        if( validTria )then
          ! append the valid triangle
          call append( me  = dump_trias,           &
            &          val = triangles( :, iTria ))
        end if
      end do

      if( present( header ))then
        loc_header = header
      else
        loc_header = ''
      end if

      ! initialize the attribute to be empty
      attribute = ''

      if( present( time ))then
        timestamp = adjustl(tem_time_sim_stamp(time))
        ! assemble the filename
        write(filename,'(a)')trim(outprefix)//'_t'//trim(timestamp)//".stl"
      else
        write(filename,'(a)')trim(outprefix)//".stl"
      end if

      ! open the output file and
      call tem_open( newunit = outUnit,        &
        &            file    = trim(filename), &
        &            access  = 'stream',       &
        &            action  = 'write',        &
        &            status  = 'replace'       )

      ! ... dump the header
      write(outUnit) loc_header
      ! ... dump the number of triangles
      write(outUnit) dump_trias%nVals

      ! .. for every triangle dump
      do iTria = 1, dump_trias%nVals
        ! ... calculate the correct positions in the linearized array
        minPos1 = (dump_trias%val( 1, iTria )-1)*3+1
        maxPos1 = (dump_trias%val( 1, iTria )-1)*3+3
        minPos2 = (dump_trias%val( 2, iTria )-1)*3+1
        maxPos2 = (dump_trias%val( 2, iTria )-1)*3+3
        minPos3 = (dump_trias%val( 3, iTria )-1)*3+1
        maxPos3 = (dump_trias%val( 3, iTria )-1)*3+3
        ! ... the normal vector
        write(outUnit)loc_normals( :, iTria )
        ! ... the vertices
        write(outUnit)real( dump_nodes( minPos1:maxPos1 ) ,kind=single_k )
        write(outUnit)real( dump_nodes( minPos2:maxPos2 ) ,kind=single_k )
        write(outUnit)real( dump_nodes( minPos3:maxPos3 ) ,kind=single_k )
        write(outUnit)attribute
      end do

      close(outunit)
      ! destroy the growing 2D array of valid triangles
      call destroy( me = dump_trias )

    end if

    deallocate( dump_nodes )

  end subroutine tem_dump_stlb