ic_2dcrvpY_for Function

public function ic_2dcrvpY_for(me, coord, n) result(res)

This function defines the y-velocity component of the spinning (= co-rotating) vortex pair Source: complex velocity potential of both vortices complex coordinates: z = x+i*y Gamma ... circulation

 b = r0*exp(i*omega*t)
 w(z,t) = Gamma/(2Pi*i)*ln(z^2-b^2)
 dw/dz = Gamma/(2Pi*i)*z/(z^2-b^2)
 u =  Re( dw/dz( z, t=0 )
 v = -Im( dw/dz( z, t=0 )

see ic_2dcrvpX_for as a documentation reference. This routine is exactly the same, except for that in the end, instead of evaluating the Re of the potential function, we have to evalute -Im

Arguments

Type IntentOptional Attributes Name
type(ic_2dcrvp_type), intent(in) :: me

global gauss pulse data

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

coordinate of an element

integer, intent(in) :: n

number of return values

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

return value which is sent to state variable


Calls

proc~~ic_2dcrvpy_for~~CallsGraph proc~ic_2dcrvpy_for ic_2dcrvpY_for proc~cutoff_factor cutoff_factor proc~ic_2dcrvpy_for->proc~cutoff_factor

Called by

proc~~ic_2dcrvpy_for~~CalledByGraph proc~ic_2dcrvpy_for ic_2dcrvpY_for proc~ic_2dcrvppressure_for ic_2dcrvpPressure_for proc~ic_2dcrvppressure_for->proc~ic_2dcrvpy_for proc~tem_spatial_for_coord tem_spatial_for_coord proc~tem_spatial_for_coord->proc~ic_2dcrvpy_for proc~tem_spatial_for_coord->proc~ic_2dcrvppressure_for interface~tem_spatial_for tem_spatial_for interface~tem_spatial_for->proc~tem_spatial_for_coord proc~tem_spatial_scalar_for_index tem_spatial_scalar_for_index interface~tem_spatial_for->proc~tem_spatial_scalar_for_index proc~tem_spatial_scalar_for_index->proc~tem_spatial_for_coord proc~tem_spacetime_for_coord tem_spacetime_for_coord proc~tem_spacetime_for_coord->interface~tem_spatial_for proc~tem_spacetime_for_treeids tem_spacetime_for_treeIDs proc~tem_spacetime_for_treeids->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_coord tem_spacetime_vector_for_coord proc~tem_spacetime_vector_for_coord->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_index->proc~tem_spacetime_vector_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_spatial_scalar_storeval tem_spatial_scalar_storeVal proc~tem_spatial_scalar_storeval->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_coord interface~tem_spacetime_for->proc~tem_spacetime_for_treeids interface~tem_spacetime_for->proc~tem_spacetime_scalar_for_index interface~tem_spacetime_for->proc~tem_spacetime_vector_for_coord interface~tem_spacetime_for->proc~tem_spacetime_vector_for_index interface~tem_spacetime_for->proc~tem_spacetime_vector_for_treeids 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 proc~tem_spacetime_for_stcoord tem_spacetime_for_stcoord proc~tem_spacetime_for_stcoord->proc~tem_spacetime_for_coord

Source Code

  function ic_2dcrvpY_for(me, coord, n) result(res)
    ! ---------------------------------------------------------------------------
    !> number of return values
    integer, intent(in) :: n
    !> global gauss pulse data
    type(ic_2dcrvp_type), intent(in) :: me
    !> coordinate of an element
    real(kind=rk), intent(in) :: coord(n, 3)
    !> return value which is sent to state variable
    real(kind=rk) :: res(n)
    ! ---------------------------------------------------------------------------
    real(kind=rk) :: r      ! radius from center of rotation
    real(kind=rk) :: rC, x_max ! core radius
    real(kind=rk) :: omega     ! angular speed
    complex(kind=rk) :: z_coord, b_coord, z_compl, vortex1, vortex2, umax
    complex(kind=rk),parameter :: i_complex = (0._rk,1._rk)
    integer :: iElem
    ! ---------------------------------------------------------------------------
    do iElem = 1, n
      z_coord    = ( coord(iElem,1) - me%center(1))                            &
        &        + i_complex*(coord(iElem,2) - me%center(2))
      z_compl    = ( coord(iElem,1) - me%center(1))                            &
        &        - i_complex*(coord(iElem,2) - me%center(2))

      r = abs(z_coord)
      rC = me%radius_C

      omega = me%circulation/(Pi*me%radius_rot*me%radius_rot*4._rk)
      b_coord    = me%radius_rot*( cos( omega*me%t )                           &
        &        + i_complex*sin( omega*me%t))

      ! left vortex. compare potential_single_w_z
      vortex1 =  me%circulation/(2._rk*PI*i_complex*(z_coord+b_coord))
      ! right vortex
      vortex2 =  me%circulation/(2._rk*PI*i_complex*(z_coord-b_coord))
      ! Get the velocity on the core radius for the core model
      x_max   = -real(b_coord + rC, kind=rk )
      umax    =  me%circulation/(2._rk*PI*i_complex*(x_max + b_coord))
      if( me%rankineModel ) then
      ! Rankine vortex model, if inside the core radius
      ! right vortex
        if( real(((z_coord) - b_coord)*(z_compl-b_coord), kind=rk)             &
          &                                                   .lt. rC*rC ) then
          vortex2 = -(umax*z_compl/rC - umax*b_coord/rC)
        elseif( real((z_coord + b_coord)*(z_compl + b_coord), kind=rk)         &
          &                                                   .lt. rC**2 ) then
          vortex1 = -(umax*z_compl/rC + umax*b_coord/rC)
        endif
      endif
      res( iElem ) = cutoff_factor(me%cutoff, r)                               &
        &          * (-aimag( vortex1 + vortex2))
    end do

  end function ic_2dcrvpY_for