tem_polygon_material_load Subroutine

public subroutine tem_polygon_material_load(me, conf, thandle)

read list of vertices

Arguments

Type IntentOptional Attributes Name
type(tem_polygon_material_type), intent(out) :: me

Polygon data structure to fill with information provided by the user in config.

type(flu_State) :: conf

Handle to the Lua script containing the polygon definition

integer, intent(in), optional :: thandle

Handle for the table containing the polygon definition.


Calls

proc~~tem_polygon_material_load~~CallsGraph proc~tem_polygon_material_load tem_polygon_material_load aot_get_val aot_get_val proc~tem_polygon_material_load->aot_get_val aot_table_close aot_table_close proc~tem_polygon_material_load->aot_table_close aot_table_length aot_table_length proc~tem_polygon_material_load->aot_table_length aot_table_open aot_table_open proc~tem_polygon_material_load->aot_table_open proc~tem_abort tem_abort proc~tem_polygon_material_load->proc~tem_abort mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~tem_polygon_material_load~~CalledByGraph proc~tem_polygon_material_load tem_polygon_material_load proc~load_spatial_predefined load_spatial_predefined proc~load_spatial_predefined->proc~tem_polygon_material_load proc~tem_load_spatial tem_load_spatial proc~tem_load_spatial->proc~load_spatial_predefined proc~load_spacetime_predefined load_spacetime_predefined proc~load_spacetime_predefined->proc~tem_load_spatial proc~tem_load_ic tem_load_ic proc~tem_load_ic->proc~tem_load_spatial proc~tem_load_spacetime_single tem_load_spacetime_single proc~tem_load_spacetime_single->proc~load_spacetime_predefined

Source Code

  subroutine tem_polygon_material_load(me, conf, thandle)
    ! ----------------------------------------------------------------------
    !> Polygon data structure to fill with information provided
    !! by the user in config.
    type(tem_polygon_material_type), intent(out) :: me

    !> Handle to the Lua script containing the polygon definition
    type(flu_state) :: conf

    !> Handle for the table containing the polygon definition.
    integer, intent(in), optional :: thandle
    ! ----------------------------------------------------------------------
    real(kind=rk), allocatable :: defout(:)
    integer :: vertex_table, valtable, vertices_table
    integer :: iVertex
    integer :: iError
    integer, allocatable :: vError(:)
    integer :: iError_v(2)
    integer :: iPoly
    ! ----------------------------------------------------------------------

    write(logUnit(1),*) 'Loading predefined function polygonal material:'
    valtable = 0
    call aot_get_val( L       = conf,    &
      &               thandle = thandle, &
      &               key     = 'zmin',  &
      &               val     = me%zmin, &
      &               default = 0.0_rk,  &
      &               ErrCode = iError   )

    if (btest(iError, aoterr_fatal)) then
      write(logunit(1),*) 'ERROR in tem_polygon_material_load: Not able to'
      write(logunit(1),*) '      get a value for zmin!'
      write(logunit(1),*) '      This also applies if no inval is provided.'
      call tem_abort()
    end if
    call aot_get_val( L       = conf,    &
      &               thandle = thandle, &
      &               key     = 'zmax',  &
      &               val     = me%zmax, &
      &               default = 0.0_rk,  &
      &               ErrCode = iError   )

    if (btest(iError, aoterr_fatal)) then
      write(logunit(1),*) 'ERROR in tem_polygon_material_load: Not able to'
      write(logunit(1),*) '      get a value for zmax!'
      write(logunit(1),*) '      This also applies if no inval is provided.'
      call tem_abort()
    end if
    call aot_table_open( L       = conf,    &
      &                  parent  = thandle, &
      &                  key     = 'inval', &
      &                  thandle = valtable )
    if (valtable  == 0) then
      ! inval not provided as a table, try to read it as a scalar.
      allocate(me%inval(1))
      call aot_get_val( L       = conf,        &
        &               thandle = thandle,     &
        &               key     = 'inval',     &
        &               val     = me%inval(1), &
        &               default = 1.0_rk,      &
        &               ErrCode = iError       )

      if (btest(iError, aoterr_Fatal)) then
        write(logunit(1),*) 'ERROR in tem_polygon_material_load: Not able to'
        write(logunit(1),*) '      get value for inval!'
        call tem_abort()
      end if

      me%nComponents = 1
      allocate(me%outval(me%nComponents))
      allocate(vError(me%nComponents))

      ! Outval needs to be consistent with the inval definition, if inval was
      ! defined as a scalar, outval also has to be a scalar!
      ! We do not check for tables with single entries in this case.
      call aot_get_val( L       = conf,         &
        &               thandle = thandle,      &
        &               key     = 'outval',     &
        &               val     = me%outval(1), &
        &               default = 0.0_rk,       &
        &               ErrCode = iError        )

      if (btest(iError, aoterr_fatal)) then
        write(logunit(1),*) 'ERROR in tem_polygon_material_load: Not able to'
        write(logunit(1),*) '      get a value for outval!'
        write(logunit(1),*) '      Note, that outval needs be a scalar, as'
        write(logunit(1),*) '      inval is provided as a scalar.'
        write(logunit(1),*) '      This also applies if no inval is provided.'
        call tem_abort()
      end if

    else

      ! Intable is a table, close it an read it into an array.
      call aot_table_close(L = conf, thandle = valtable)

      ! Value to use inside the polygon
      call aot_get_val( L         = conf,     &
        &               thandle   = thandle,  &
        &               key       = 'inval',  &
        &               val       = me%inval, &
        &               maxlength = 20,       &
        &               default   = [1.0_rk], &
        &               ErrCode   = vError    )

      if (any(btest(vError, aoterr_Fatal))) then
        write(logunit(1),*) 'ERROR in tem_polygon_material_load: Not able to'
        write(logunit(1),*) '      get values for inval!'
        call tem_abort()
      end if
      me%nComponents = size(me%inval)
      deallocate(vError)

      ! Definition of outval needs to be consistent with inval, it has to have
      ! the same number of components, and also needs to be a vector.
      ! However, we define a default of all zeroes, so if outval is 0 for all
      ! components, this definition can be omitted in the user definition.
      allocate(me%outval(me%nComponents))
      allocate(vError(me%nComponents))
      allocate(defout(me%nComponents))

      defout = 0.0_rk

      ! Value to use outside the polygon
      call aot_get_val( L       = conf,      &
        &               thandle = thandle,   &
        &               key     = 'outval',  &
        &               val     = me%outval, &
        &               default = defout,    &
        &               ErrCode = vError     )

      if (any(btest(vError,aoterr_Fatal))) then
        write(logunit(1),*) 'ERROR in tem_polygon_material_load: Not able to'
        write(logunit(1),*) '      get a value for outval!'
        write(logunit(1),*) '      Note, that outval needs to have the same'
        write(logunit(1),*) '      length as inval.'
        call tem_abort()
      end if

      deallocate(vError)
      deallocate(defout)

    end if

    !> read list of vertices
    call aot_table_open( L       = conf,          &
      &                  parent  = thandle,       &
      &                  key     = 'vertices',    &
      &                  thandle = vertices_table )

    if (vertices_table == 0) then
      write(logunit(1),*) 'ERROR in tem_polygon_material_load: No vertices'
      write(logunit(1),*) '      defined, unable to set up a polynomial!'
      write(logunit(1),*) '      Please define vertices tables with lists of'
      write(logunit(1),*) '      2D points.'
      call tem_abort()
    end if
    me%nPoly = aot_table_length( L       = conf,          &
      &                          thandle = vertices_table )
    
    allocate(me%poly_list(me%npoly))
   
   ! me%npoly = 1
    do iPoly = 1, me%npoly
      call aot_table_open( L       = conf,           &
        &                  parent  = vertices_table, &
        &                  pos     = ipoly,          &
        &                  thandle = vertex_table    )
      if (vertex_table == 0) then
        write(logunit(1),*) 'ERROR in tem_polygon_material_load:'
        write(logunit(1),*) '      Please define vertex tables with lists of'
        write(logunit(1),*) '      2D points.'
        call tem_abort()
      end if
      me%poly_list(ipoly)%nVertices =              &
        & aot_table_length( L       = conf,        &
        &                   thandle = vertex_table )

      allocate(me%poly_list(ipoly)%vertex(me%poly_list(ipoly)%nVertices, 2))
      do iVertex=1,me%poly_list(ipoly)%nVertices
        call aot_get_val( val     = me%poly_list(ipoly)%vertex(iVertex, :), &
          &               ErrCode = iError_v,                                  &
          &               L       = conf,                                      &
          &               thandle = vertex_table,                              &
          &               pos     = iVertex                                    )
        if (any(iError_v /= 0)) then
          write(logunit(1),*) 'ERROR in tem_polygon_material_load: Not able to'
          write(logunit(1),*) '      obtain vertex ', iVertex, '!'
          write(logunit(1),*) '      Vertices have to be vectors of length 2,'
          write(logunit(1),*) '      with real numbers as entries.'
          call tem_abort()
        end if
      end do
      call aot_table_close( L       = conf,        &
        &                   thandle = vertex_table )

    end do
    call aot_table_close( L       = conf,          &
      &                   thandle = vertices_table )

  end subroutine tem_polygon_material_load