tem_load_stl Subroutine

public subroutine tem_load_stl(stl_data, transform, conf, thandle)

This routine loads STL files from config and reads the triangles from the files into the dynamic array of triangles.

Arguments

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

Array array of triangles in stlData

type(tem_transformation_type), intent(in) :: transform

transformation for spatial object

type(flu_State) :: conf

Lua state

integer, intent(in) :: thandle

Calls

proc~~tem_load_stl~~CallsGraph proc~tem_load_stl tem_load_stl proc~tem_load_stlhead tem_load_stlHead proc~tem_load_stl->proc~tem_load_stlhead proc~tem_read_stlfiles tem_read_stlFiles proc~tem_load_stl->proc~tem_read_stlfiles proc~aot_table_open aot_table_open proc~tem_load_stlhead->proc~aot_table_open proc~tem_load_stlhead_single tem_load_stlHead_single proc~tem_load_stlhead->proc~tem_load_stlhead_single proc~aot_table_close aot_table_close proc~tem_load_stlhead->proc~aot_table_close proc~aot_table_length aot_table_length proc~tem_load_stlhead->proc~aot_table_length stla_size stla_size proc~tem_read_stlfiles->stla_size proc~tem_read_stlb tem_read_stlb proc~tem_read_stlfiles->proc~tem_read_stlb stla_read stla_read proc~tem_read_stlfiles->stla_read proc~tem_size_stlb tem_size_stlb proc~tem_read_stlfiles->proc~tem_size_stlb proc~aot_get_val~2 aot_get_val proc~tem_load_stlhead_single->proc~aot_get_val~2 proc~tem_abort tem_abort proc~tem_load_stlhead_single->proc~tem_abort proc~tem_open tem_open proc~tem_read_stlb->proc~tem_open 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~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower proc~newunit newunit proc~tem_open->proc~newunit

Called by

proc~~tem_load_stl~~CalledByGraph proc~tem_load_stl tem_load_stl 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 interface~tem_load_shape->proc~tem_load_shapes 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_spacetime_single->proc~tem_load_spacetime_single proc~tem_load_trackingconfig tem_load_trackingConfig proc~tem_load_trackingconfig->interface~tem_load_shape proc~load_spatial_predefined load_spatial_predefined proc~load_spatial_predefined->proc~load_spatial_parabol proc~tem_load_spacetime_table tem_load_spacetime_table proc~tem_load_spacetime_table->proc~tem_load_spacetime_single interface~tem_load_spacetime tem_load_spacetime interface~tem_load_spacetime->proc~tem_load_spacetime_single proc~tem_convergence_load tem_convergence_load proc~tem_convergence_load->proc~tem_load_convergenceheader proc~tem_load_tracking tem_load_tracking proc~tem_load_tracking->proc~tem_load_trackingconfig

Contents

Source Code


Source Code

  subroutine tem_load_stl(stl_data, transform, conf, thandle)
    ! --------------------------------------------------------------------------!
    !> Array array of triangles in stlData
    type(tem_stlData_type), intent(out) :: stl_data
    !> transformation for spatial object
    type(tem_transformation_type), intent(in) :: transform
    !> Lua state
    type(flu_state) :: conf
    integer, intent(in) :: thandle !< handle for canonical objects
    ! --------------------------------------------------------------------------!
    integer :: iTri, iVer
    ! --------------------------------------------------------------------------!

    ! Load stl files from config file.
    call tem_load_stlhead(me = stl_data%head, conf = conf, thandle = thandle)

    ! Load triangles and nodes from stl files.
    call tem_read_stlFiles(stl_data = stl_data )

    ! if transformation is active apply transformation to all triangles
    ! first apply deformation and then translation
    if(transform%active) then
      if(transform%deform%active) then
        do iTri=1, stl_data%nTris
          do iVer=1,3
            stl_data%nodes(:, stl_data%tri_node(iVer, iTri)) = &
              &            matmul( transform%deform%matrix, &
              &            stl_data%nodes(:, stl_data%tri_node(iVer, iTri)) )
          enddo
        enddo
      endif
      if(transform%translate%active) then
        do iTri=1, stl_data%nTris
          do iVer=1,3
            stl_data%nodes(:, stl_data%tri_node(iVer, iTri)) = &
              &            transform%translate%vec + &
              &            stl_data%nodes(:, stl_data%tri_node(iVer, iTri))
          enddo
        enddo
      endif
    endif

  end subroutine tem_load_stl