tem_readAndUnify_surfData Subroutine

public subroutine tem_readAndUnify_surfData(me, useInitPos)

This routine reads the surface data from a set of stl files and stores it in the surfaceData_type.

Arguments

Type IntentOptional Attributes Name
type(tem_surfData_type), intent(inout) :: me

datatype to store the surface information

logical, intent(in), optional :: useInitPos

shall the initial points be stored and used for updating the points later on ???


Calls

proc~~tem_readandunify_surfdata~~CallsGraph proc~tem_readandunify_surfdata tem_readAndUnify_surfData proc~tem_size_stlb tem_size_stlb proc~tem_readandunify_surfdata->proc~tem_size_stlb proc~tem_read_stlb tem_read_stlb proc~tem_readandunify_surfdata->proc~tem_read_stlb proc~tem_abort tem_abort proc~tem_readandunify_surfdata->proc~tem_abort proc~tem_unify_surfacedata tem_unify_surfaceData proc~tem_readandunify_surfdata->proc~tem_unify_surfacedata proc~tem_size_stlb->proc~tem_abort proc~tem_open tem_open proc~tem_size_stlb->proc~tem_open proc~tem_read_stlb->proc~tem_open mpi_abort mpi_abort proc~tem_abort->mpi_abort interface~placeat~35 placeat proc~tem_unify_surfacedata->interface~placeat~35 proc~tem_idofcoord tem_IdOfCoord proc~tem_unify_surfacedata->proc~tem_idofcoord interface~init~22 init proc~tem_unify_surfacedata->interface~init~22 interface~append~23 append proc~tem_unify_surfacedata->interface~append~23 proc~tem_coordofreal tem_CoordOfReal proc~tem_unify_surfacedata->proc~tem_coordofreal interface~destroy~22 destroy proc~tem_unify_surfacedata->interface~destroy~22 proc~placeat_singlega2d_real placeat_singlega2d_real interface~placeat~35->proc~placeat_singlega2d_real proc~placeat_arrayga2d_real placeat_arrayga2d_real interface~placeat~35->proc~placeat_arrayga2d_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~init_ga2d_real init_ga2d_real interface~init~22->proc~init_ga2d_real 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 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~placeat_singlega2d_real->interface~expand~22 proc~placeat_arrayga2d_real->interface~expand~22 proc~append_singlega2d_real->interface~expand~22

Contents


Source Code

  subroutine tem_readAndUnify_surfData( me, useInitPos )
    ! ---------------------------------------------------------------------------
    !> datatype to store the surface information
    type( tem_surfData_type ), intent(inout) :: me
    !> shall the initial points be stored and used for updating the points
    !! later on ???
    logical, optional, intent(in) :: useInitPos
    ! ---------------------------------------------------------------------------
    ! number of points to be read from the mesh
    integer :: nPoints( size( me%stlHead ))
    ! number of triangles to be read from the mesh
    integer :: nTrias( size( me%stlHead ))
    ! total number of points to be read from file
    integer :: nPoints_total
    ! tmp point coordinates
    real(kind=rk), allocatable :: tmp_pointCoords(:,:)
    ! Offset in the temporary pointCoords list for multiple STLs
    integer :: nodeOffset
    ! Offset in the triangle list for multiple STLs
    integer :: triOffset
    integer :: iFile
    integer :: ierr
    logical :: tmp_useInitPos
    ! min and max of all X,Y,Z coordinates
    real(kind=rk) :: minX, minY, minZ, maxX, maxY, maxZ
    ! ---------------------------------------------------------------------------
    if( present( useInitPos ))then
      tmp_useInitPos = useInitPos
    else
      tmp_useInitPos = .false.
    end if

    write(logUnit(10),*) " Reading in STL Headers ..."
    ! Read headerfiles to determine number of nodes and tris
    do iFile = 1, size(me%stlHead)
      call tem_size_stlb( filename = me%stlHead(iFile)%filename,               &
        &                 nNodes   = nPoints( iFile ),                         &
        &                 nTris    = nTrias( iFile ))
    end do

    nPoints_total = sum( nPoints( : ))
    me%nTrias     = sum(  nTrias( : ))

    write(logUnit(2),*) " Total number of triangles: ", me%nTrias
    write(logUnit(2),*) " Total number of points:    ", nPoints_total

    ! allocate the temporary array of coordinates
    allocate(tmp_pointCoords( 3, nPoints_total ))
    ! allocate the global number of triangles (this is constant regardless of
    ! unification)
    allocate(me%trias( 3, me%nTrias ))

    write(logUnit(10),*) " Reading in STL Files ..."

    ! Read in node values from STL files
    nodeOffset = 1
    triOffset = 1
    do iFile = 1, size(me%stlHead)
      ! read in the nodes and triangles to the right positions in the temporary
      ! pointCoords array and in the global triangle array
      call tem_read_stlb( filename   = me%stlHead( iFile )%filename,           &
        &                 nNodesRead = nPoints( iFile ),                       &
        &                 nTrisRead  = nTrias( iFile ),                        &
        &                 nodes      = tmp_pointCoords( 1:3,                   &
        &                           nodeOffset:nodeOffset+nPoints( iFile )-1), &
        &                 tri_node   = me%trias(1:3,                           &
        &                              triOffset:triOffset+nTrias( iFile )-1), &
        &                 ierror     = ierr)

      if( iErr /= 0) then
        write(logUnit(0),*) "An error appeared when reading the surface data " &
          &               //"file "//trim( me%stlHead( iFile )%filename )      &
          &               //". Stopping!!!"
        call tem_abort()
      end if

      ! update the positions in the triangle array according to the nodes read
      me%trias(1:3,triOffset:triOffset+nTrias( iFile )-1) =                    &
        &  me%trias(1:3,triOffset:triOffset+nTrias( iFile )-1) + (nodeOffset-1)

      ! update the triangle and node offsets
      nodeOffset = nodeOffset + nPoints( iFile )
      triOffset = triOffset + nTrias( iFile )

    end do

    write(logUnit(10),*) " Done Reading. Starting to unify ..."

    ! unify the coordinates
    call tem_unify_surfaceData( me, tmp_pointCoords )

    write(logUnit(10),*) " Done Unifying. Starting to calulate the initial "//  &
      &                  "surface areas ..."

    ! allocate the surface areas array with the number of unique surface points
    allocate( me%surfArea( me%nUniquePoints_total ))

    deallocate( tmp_pointCoords )

    write(logUnit(10),*) " Done calculating the surface areas."
    write(logUnit(10),*) " Number of unique surface points: ",                 &
      &                  me%nUniquePoints_total
    write(logUnit(10),*) " Number of triangles: ",                             &
      &                  me%nTrias

    if( tmp_useInitPos )then
      allocate( me%backPointCoords( 3*me%nUniquePoints_total))
      me%backPointCoords = me%pointCoords
    end if

    ! get the min and max coordinates
    minX = minval(me%pointCoords(1:((me%nUniquePoints_total-1))*3+1:3))
    minY = minval(me%pointCoords(2:((me%nUniquePoints_total-1)*3+2):3))
    minZ = minval(me%pointCoords(3:me%nUniquePoints_total*3:3))
    maxX = maxval(me%pointCoords(1:((me%nUniquePoints_total-1)*3+1):3))
    maxY = maxval(me%pointCoords(2:((me%nUniquePoints_total-1)*3+2):3))
    maxZ = maxval(me%pointCoords(3:me%nUniquePoints_total*3:3))

    write(logUnit(3),*)'minX: ', minX
    write(logUnit(3),*)'minY: ', minY
    write(logUnit(3),*)'minZ: ', minZ
    write(logUnit(3),*)'maxX: ', maxX
    write(logUnit(3),*)'maxY: ', maxY
    write(logUnit(3),*)'maxZ: ', maxZ

  end subroutine tem_readAndUnify_surfData