tem_spatial_for_coord Function

private function tem_spatial_for_coord(me, coord, n) 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 is either 1.0 or default value provided in tem_load_spatial. 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), intent(in) :: 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

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

return value of a function


Calls

proc~~tem_spatial_for_coord~~CallsGraph proc~tem_spatial_for_coord tem_spatial_for_coord proc~tem_spatial_random_for tem_spatial_random_for proc~tem_spatial_for_coord->proc~tem_spatial_random_for proc~ic_tgv_uy_for ic_tgv_uy_for proc~tem_spatial_for_coord->proc~ic_tgv_uy_for proc~ic_tgv_syz_for ic_tgv_Syz_for proc~tem_spatial_for_coord->proc~ic_tgv_syz_for proc~tem_eval_polygon_material_scal tem_eval_polygon_material_scal proc~tem_spatial_for_coord->proc~tem_eval_polygon_material_scal proc~tem_eval_heaviside_gibbs tem_eval_heaviside_gibbs proc~tem_spatial_for_coord->proc~tem_eval_heaviside_gibbs interface~tem_spatial_lua_for tem_spatial_lua_for proc~tem_spatial_for_coord->interface~tem_spatial_lua_for interface~tem_spatial_parabol3d_for tem_spatial_parabol3d_for proc~tem_spatial_for_coord->interface~tem_spatial_parabol3d_for proc~tem_eval_miescatter_magny tem_eval_miescatter_magny proc~tem_spatial_for_coord->proc~tem_eval_miescatter_magny proc~ic_tgv_sxz_for ic_tgv_Sxz_for proc~tem_spatial_for_coord->proc~ic_tgv_sxz_for interface~tem_viscspongelayer_plane_for tem_viscSpongeLayer_plane_for proc~tem_spatial_for_coord->interface~tem_viscspongelayer_plane_for proc~tem_abort tem_abort proc~tem_spatial_for_coord->proc~tem_abort interface~tem_spatial_parabol2d_for tem_spatial_parabol2d_for proc~tem_spatial_for_coord->interface~tem_spatial_parabol2d_for interface~tem_viscspongelayer_radial_for tem_viscSpongeLayer_radial_for proc~tem_spatial_for_coord->interface~tem_viscspongelayer_radial_for proc~tem_eval_cylindricalwave tem_eval_cylindricalWave proc~tem_spatial_for_coord->proc~tem_eval_cylindricalwave interface~tem_spongelayer_box2d_for tem_spongeLayer_box2d_for proc~tem_spatial_for_coord->interface~tem_spongelayer_box2d_for interface~tem_viscspongelayer_box2d_for tem_viscSpongeLayer_box2d_for proc~tem_spatial_for_coord->interface~tem_viscspongelayer_box2d_for proc~ic_2dcrvpx_for ic_2dcrvpX_for proc~tem_spatial_for_coord->proc~ic_2dcrvpx_for proc~tem_eval_miescatter_displz tem_eval_miescatter_displz proc~tem_spatial_for_coord->proc~tem_eval_miescatter_displz proc~ic_tgv_pressure_for ic_tgv_pressure_for proc~tem_spatial_for_coord->proc~ic_tgv_pressure_for proc~ic_tgv_ux_for ic_tgv_ux_for proc~tem_spatial_for_coord->proc~ic_tgv_ux_for proc~tem_eval_polygon_material_scal_3d tem_eval_polygon_material_scal_3d proc~tem_spatial_for_coord->proc~tem_eval_polygon_material_scal_3d interface~tem_spongelayer_plane_for tem_spongeLayer_plane_for proc~tem_spatial_for_coord->interface~tem_spongelayer_plane_for proc~ic_gausspulse_for ic_gausspulse_for proc~tem_spatial_for_coord->proc~ic_gausspulse_for proc~ic_tgv_sxx_for ic_tgv_Sxx_for proc~tem_spatial_for_coord->proc~ic_tgv_sxx_for interface~tem_spongelayer_radial_for tem_spongeLayer_radial_for proc~tem_spatial_for_coord->interface~tem_spongelayer_radial_for interface~tem_viscspongelayer_box_for tem_viscSpongeLayer_box_for proc~tem_spatial_for_coord->interface~tem_viscspongelayer_box_for proc~ic_2dcrvpy_for ic_2dcrvpY_for proc~tem_spatial_for_coord->proc~ic_2dcrvpy_for proc~ic_2dcrvppressure_for ic_2dcrvpPressure_for proc~tem_spatial_for_coord->proc~ic_2dcrvppressure_for proc~tem_eval_miescatter_magnx tem_eval_miescatter_magnx proc~tem_spatial_for_coord->proc~tem_eval_miescatter_magnx proc~ic_tgv_syy_for ic_tgv_Syy_for proc~tem_spatial_for_coord->proc~ic_tgv_syy_for interface~tem_spongelayer_box_for tem_spongeLayer_box_for proc~tem_spatial_for_coord->interface~tem_spongelayer_box_for

Called by

proc~~tem_spatial_for_coord~~CalledByGraph proc~tem_spatial_for_coord tem_spatial_for_coord proc~tem_spatial_scalar_for_index tem_spatial_scalar_for_index proc~tem_spatial_scalar_for_index->proc~tem_spatial_for_coord interface~tem_spatial_for tem_spatial_for interface~tem_spatial_for->proc~tem_spatial_for_coord interface~tem_spatial_for->proc~tem_spatial_scalar_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


Source Code

  function tem_spatial_for_coord( me, coord, n ) result( res )
    ! -------------------------------------------------------------------- !
    !> spatial type for given boundary state
    type(tem_spatial_type), intent(in) :: me
    !> number of return values
    integer, intent(in) :: n
    !> 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)
    ! -------------------------------------------------------------------- !

    select case( trim(adjustl(me%kind)) )
    case( 'none', 'const' )
      res = me%const(1)

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

    case ('random')
      res = tem_spatial_random_for(me%random, n)

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

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

      end select

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

    case( 'viscous_spongelayer_plane' )
      res = tem_viscSpongelayer_plane_for( &
        & me    = me%spongePlane,          &
        & coord = coord,                   &
        & n     = n                        )

    case( 'viscous_spongelayer_box' )
      res = tem_viscSpongelayer_box_for( &
        & me    = me%spongeBox,          &
        & coord = coord,                 &
        & n     = n                      )

    case( 'viscous_spongelayer_box_2d' )
      res = tem_viscSpongelayer_box2d_for( &
        & me    = me%spongeBox,            &
        & coord = coord,                   &
        & n     = n                        )

    case( 'viscous_spongelayer_radial_2d' )
      res = tem_viscSpongelayer_radial_for( &
        & me    = me%spongeRadial,          &
        & coord = coord,                    &
        & nDim  = 2,                        &
        & n     = n                         )

    case( 'viscous_spongelayer_radial' )
      res = tem_viscSpongelayer_radial_for( &
        & me    = me%spongeRadial,          &
        & coord = coord,                    &
        & nDim  = 3,                        &
        & n     = n                         )

    case( 'gausspulse' )
      res = ic_gausspulse_for(me%gausspulse, coord, n)

    case( 'crvpX' )
      res = ic_2dcrvpX_for(me%crvp, coord, n)

    case( 'crvpY' )
      res = ic_2dcrvpY_for(me%crvp, coord, n)

    case( 'crvpPressure' )
      res = ic_2dcrvpPressure_for(me%crvp, coord, n)

    case('miescatter_displacementfieldz')
      res = tem_eval_miescatter_displz( me    = me%mie_fun, &
        &                               coord = coord,      &
        &                               time  = 0.0_rk,     &
        &                               n     = n           )

    case('miescatter_magneticfieldx')
      res = tem_eval_miescatter_magnx( me    = me%mie_fun, &
        &                              coord = coord,      &
        &                              time  = 0.0_rk,     &
        &                              n     = n           )

    case('miescatter_magneticfieldy')
      res = tem_eval_miescatter_magny( me    = me%mie_fun, &
        &                              coord = coord,      &
        &                              time  = 0.0_rk,     &
        &                              n     = n           )

    case('heaviside_gibbs')
      res = tem_eval_heaviside_gibbs( me    = me%heaviside_gibbs_fun, &
        &                             coord = coord,                  &
        &                             n     = n                       )

    case('cylindricalwave')
      res = tem_eval_cylindricalWave( me    = me%cylindricalWave, &
        &                             coord = coord,              &
        &                             time  = 0.0_rk,             &
        &                             n     = n                   )

    case('tgv_p')
      res = ic_tgv_pressure_for( me%tgv, coord, n )

    case('tgv_ux')
      res = ic_tgv_ux_for( me%tgv, coord, n )

    case('tgv_uy')
      res = ic_tgv_uy_for( me%tgv, coord, n )

    case('tgv_sxx')
      res = ic_tgv_sxx_for( me%tgv, coord, n )

    case('tgv_syy')
      res = ic_tgv_syy_for( me%tgv, coord, n )

    case('tgv_sxz')
      res = ic_tgv_sxz_for( me%tgv, coord, n )

    case('tgv_syz')
      res = ic_tgv_syz_for( me%tgv, coord, n )

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

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

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

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

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

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

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

    case default
      call tem_abort(                                                 &
        & 'ERROR: Unknown spatial function in tem_spatial_for_coord.' )

    end select

   end function tem_spatial_for_coord