tem_spatial_vector_for_coord Function

private function tem_spatial_vector_for_coord(me, coord, n, ncomp) result(res)

This function invokes different spatial boundary kinds like constant, lua function or predefined Fortran function for given coord

If a spatial block is not defined and a temporal block is defined in the lua file, the return value = ref_value. If both spatial and temporal block are not defined in the lua file, the return value = 1.0_rk. based spatial_kind(kind).

  1. const - set constant value
  2. lua_fun - lua function
  3. gausspulse - fortran gauss pulse function
  4. 2dcrvp - fortran spinning vortex function
  5. parabol - fotran parabolic function

Arguments

Type IntentOptional Attributes Name
type(tem_spatial_type) :: me

spatial type for given boundary state

real(kind=rk), intent(in) :: coord(n,3)

barycentric Ids of an elements. 1st index goes over number of elements and 2nd index goes over x,y,z coordinates

integer, intent(in) :: n

number of return values

integer, intent(in) :: ncomp

number of components per returned value

Return Value real(kind=rk), (n,ncomp)

return value of a function


Calls

proc~~tem_spatial_vector_for_coord~~CallsGraph proc~tem_spatial_vector_for_coord tem_spatial_vector_for_coord interface~tem_spatial_lua_for tem_spatial_lua_for proc~tem_spatial_vector_for_coord->interface~tem_spatial_lua_for interface~tem_spatial_parabol3d_for tem_spatial_parabol3d_for proc~tem_spatial_vector_for_coord->interface~tem_spatial_parabol3d_for interface~tem_spongelayer_radial_for tem_spongeLayer_radial_for proc~tem_spatial_vector_for_coord->interface~tem_spongelayer_radial_for interface~tem_spatial_parabol2d_for tem_spatial_parabol2d_for proc~tem_spatial_vector_for_coord->interface~tem_spatial_parabol2d_for interface~tem_spongelayer_box2d_for tem_spongeLayer_box2d_for proc~tem_spatial_vector_for_coord->interface~tem_spongelayer_box2d_for interface~tem_spongelayer_box_for tem_spongeLayer_box_for proc~tem_spatial_vector_for_coord->interface~tem_spongelayer_box_for interface~tem_spongelayer_plane_for tem_spongeLayer_plane_for proc~tem_spatial_vector_for_coord->interface~tem_spongelayer_plane_for proc~tem_evaluate_pml tem_evaluate_pml proc~tem_spatial_vector_for_coord->proc~tem_evaluate_pml proc~tem_eval_polygon_material tem_eval_polygon_material proc~tem_spatial_vector_for_coord->proc~tem_eval_polygon_material proc~tem_eval_polygon_material_3d tem_eval_polygon_material_3d proc~tem_spatial_vector_for_coord->proc~tem_eval_polygon_material_3d proc~tem_abort tem_abort proc~tem_spatial_vector_for_coord->proc~tem_abort

Called by

proc~~tem_spatial_vector_for_coord~~CalledByGraph proc~tem_spatial_vector_for_coord tem_spatial_vector_for_coord proc~tem_spatial_vector_for_index tem_spatial_vector_for_index proc~tem_spatial_vector_for_index->proc~tem_spatial_vector_for_coord interface~tem_spatial_for tem_spatial_for interface~tem_spatial_for->proc~tem_spatial_vector_for_coord interface~tem_spatial_for->proc~tem_spatial_vector_for_index proc~tem_spacetime_for_treeids tem_spacetime_for_treeIDs proc~tem_spacetime_for_treeids->interface~tem_spatial_for proc~tem_spacetime_vector_for_index tem_spacetime_vector_for_index proc~tem_spacetime_vector_for_index->interface~tem_spatial_for proc~tem_spacetime_vector_for_coord tem_spacetime_vector_for_coord proc~tem_spacetime_vector_for_index->proc~tem_spacetime_vector_for_coord proc~tem_spatial_scalar_storeval tem_spatial_scalar_storeVal proc~tem_spatial_scalar_storeval->interface~tem_spatial_for proc~tem_spacetime_for_coord tem_spacetime_for_coord proc~tem_spacetime_for_coord->interface~tem_spatial_for proc~tem_spacetime_scalar_for_index tem_spacetime_scalar_for_index proc~tem_spacetime_scalar_for_index->interface~tem_spatial_for proc~tem_spacetime_scalar_for_index->proc~tem_spacetime_for_coord proc~tem_spacetime_vector_for_treeids tem_spacetime_vector_for_treeIDs proc~tem_spacetime_vector_for_treeids->interface~tem_spatial_for proc~tem_spacetime_vector_for_coord->interface~tem_spatial_for proc~tem_spatial_vector_storeval tem_spatial_vector_storeVal proc~tem_spatial_vector_storeval->interface~tem_spatial_for interface~tem_spacetime_for tem_spacetime_for interface~tem_spacetime_for->proc~tem_spacetime_for_treeids interface~tem_spacetime_for->proc~tem_spacetime_vector_for_index interface~tem_spacetime_for->proc~tem_spacetime_for_coord interface~tem_spacetime_for->proc~tem_spacetime_scalar_for_index interface~tem_spacetime_for->proc~tem_spacetime_vector_for_treeids interface~tem_spacetime_for->proc~tem_spacetime_vector_for_coord proc~tem_spacetime_for_stcoord tem_spacetime_for_stcoord proc~tem_spacetime_for_stcoord->proc~tem_spacetime_for_coord interface~tem_spatial_storeval tem_spatial_storeVal interface~tem_spatial_storeval->proc~tem_spatial_scalar_storeval interface~tem_spatial_storeval->proc~tem_spatial_vector_storeval

Contents


Source Code

  function tem_spatial_vector_for_coord( me, coord, n, ncomp ) result( res )
    ! -------------------------------------------------------------------- !
    !> spatial type for given boundary state
    type(tem_spatial_type) :: me
    !> number of return values
    integer, intent(in) :: n
    !> number of components per returned value
    integer, intent(in) :: ncomp
    !> barycentric Ids of an elements.
    !! 1st index goes over number of elements and
    !! 2nd index goes over x,y,z coordinates
    real(kind=rk), intent(in) :: coord(n,3)
    !> return value of a function
    real( kind=rk ) :: res(n,ncomp)
    ! -------------------------------------------------------------------- !
    integer :: i
    ! -------------------------------------------------------------------- !

    select case( trim(adjustl(me%kind)) )
    case( 'none', 'const' )
      do i = 1, nComp
        res(:,i) = me%const(i)
      end do

    case( 'lua_fun' )
      res = tem_spatial_lua_for(me%lua_fun_ref, me%conf, coord, n, ncomp)

    case( 'spongelayer_plane', 'spongelayer_plane_2d', 'spongelayer_plane_1d' )
      res = tem_spongeLayer_plane_for(me%spongePlane, nComp, coord, n)

    case( 'spongeLayer_box' )
      res = tem_spongeLayer_box_for(me%spongeBox, nComp, coord, n)

    case( 'spongeLayer_box_2d' )
      res = tem_spongeLayer_box2d_for(me%spongeBox, nComp, coord, n)

    case( 'spongelayer_radial_2d' )
      res = tem_spongeLayer_radial_for( &
        & me    = me%spongeRadial,      &
        & nComp = nComp,                &
        & coord = coord,                &
        & n     = n,                    &
        & nDim  = 2                     )

    case( 'spongelayer_radial' )
      res = tem_spongeLayer_radial_for( &
        & me    = me%spongeRadial,      &
        & nComp = nComp,                &
        & coord = coord,                &
        & n     = n,                    &
        & nDim  = 3                     )

    case( 'radial_spongelayer_2d' )
      res = tem_spongelayer_radial_for( &
        & me    = me%spongeRadial,      &
        & nComp = nComp,                &
        & coord = coord,                &
        & nDim  = 2,                    &
        & n     = n                     )

    case( 'radial_spongelayer' )
      res = tem_spongelayer_radial_for( &
        & me    = me%spongeRadial,      &
        & nComp = nComp,                &
        & coord = coord,                &
        & nDim  = 3,                    &
        & n     = n                     )

    case( 'pml' )
      res = tem_evaluate_pml(me%pml, nComp, coord, n)

    case( 'polygon_material' )
      res = tem_eval_polygon_material( me    = me%polygon_material, &
        &                              coord = coord,               &
        &                              n     = n                    )

    case( 'polygon_material_3d' )
      res = tem_eval_polygon_material_3d( me    = me%polygon_material, &
        &                                 coord = coord,               &
        &                                 n     = n                    )

    case( 'parabol' )
      select case( trim(adjustl(me%parabol%geometry%canoND(1)%kind)) )
      case( 'line' )
        res(:,1) = tem_spatial_parabol2d_for(      &
          & me    = me%parabol%geometry%canoND(1), &
          & coord = coord,                         &
          & n     = n                              )

      case( 'plane' )
        res(:,1) = tem_spatial_parabol3d_for(      &
          & me    = me%parabol%geometry%canoND(1), &
          & coord = coord,                         &
          & n     = n                              )

      end select

      do i = 2, nComp
        res(:,i) = res(:,1)*me%parabol%amplitude(i)
      end do

      res(:,1) = res(:,1)*me%parabol%amplitude(1)

    case('rectangular', 'gate')
      do i = 1, n
        if( (abs(coord(i,2)) < me%rect_ly)       &
          & .and. (abs(coord(i,3)) < me%rect_lz) ) then
          res(i,1) = 1.0_rk
        else
          res(i,1) = 0.0_rk
        end if
      end do

      res(:, 2:nComp) = 0.0_rk

    case default
      write(logUnit(1),*)'ERROR: No vectorial routine for spatial functions of'
      write(logUnit(1),*)'       kind "' // trim(adjustl(me%kind)) // '"'
      write(logUnit(1),*)'       available! Have a look in the ' &
        & // 'tem_spatial_module'
      write(logUnit(1),*)'       for implemented vectorial functions.'
      call tem_abort()

    end select

   end function tem_spatial_vector_for_coord