get_fio_vector_potential Subroutine

public subroutine get_fio_vector_potential(prtcls, F, params)

Arguments

Type IntentOptional AttributesName
type(PARTICLES), intent(inout) :: prtcls
type(FIELDS), intent(in) :: F
type(KORC_PARAMS), intent(in) :: params

Calls

proc~~get_fio_vector_potential~~CallsGraph proc~get_fio_vector_potential get_fio_vector_potential interface~fio_eval_field fio_eval_field proc~get_fio_vector_potential->interface~fio_eval_field

Called by

proc~~get_fio_vector_potential~~CalledByGraph proc~get_fio_vector_potential get_fio_vector_potential proc~interp_fields interp_fields proc~interp_fields->proc~get_fio_vector_potential proc~get_fields get_fields proc~get_fields->proc~interp_fields proc~mh_psi MH_psi proc~mh_psi->proc~get_fields proc~sample_hollmann_distribution_3d_psi sample_Hollmann_distribution_3D_psi proc~sample_hollmann_distribution_3d_psi->proc~get_fields proc~unitvectors unitVectors proc~unitvectors->proc~get_fields proc~get_hollmann_distribution_3d_psi get_Hollmann_distribution_3D_psi proc~get_hollmann_distribution_3d_psi->proc~sample_hollmann_distribution_3d_psi proc~gyro_distribution gyro_distribution proc~gyro_distribution->proc~unitvectors

Contents


Source Code

  subroutine get_fio_vector_potential(prtcls, F, params)
    TYPE(PARTICLES), INTENT(INOUT) :: prtcls
    TYPE(FIELDS), INTENT(IN)       :: F
    TYPE(KORC_PARAMS), INTENT(IN)  :: params
    INTEGER (C_INT)                :: status
    INTEGER                        :: pp
    REAL(rp), DIMENSION(3)         :: x
    REAL(rp), DIMENSION(3)         :: Atmp
    integer(ip)  ::  ss

    if (size(prtcls%Y,1).eq.1) then
       ss = size(prtcls%Y,1)
    else
       if (prtcls%Y(2,1).eq.0) then
          ss=1_idef
       else
          ss = size(prtcls%Y,1)
       end if
    endif

    Atmp=0._rp
    
    !$OMP PARALLEL DO DEFAULT(SHARED) &
    !$OMP& PRIVATE(pp,status,x) &
    !$OMP& FIRSTPRIVATE(Atmp)
    do pp = 1,ss
       if (prtcls%flagCon(pp) .EQ. 1_is) then
          x(1) = prtcls%Y(pp,1)*params%cpp%length
          x(2) = prtcls%Y(pp,2)
          x(3) = prtcls%Y(pp,3)*params%cpp%length

          !prtcls%hint(pp)=c_null_ptr

          !write(output_unit_write,*) F%FIO_A,x,Atmp

          status = fio_eval_field(F%FIO_A, x(1),                      &
               Atmp(1),prtcls%hint(pp))

          if (status .eq. FIO_SUCCESS) then
             prtcls%PSI_P(pp)=-Atmp(2)*x(1)
          else if (status .eq. FIO_NO_DATA) then
             prtcls%PSI_P(pp) = 100._rp
             prtcls%flagCon(pp) = 0_is
          else if (status .ne. FIO_SUCCESS) then
             prtcls%flagCon(pp) = 0_is
             prtcls%PSI_P(pp) = 100._rp
          end if


       end if
    end do
    !$OMP END PARALLEL DO

  end subroutine get_fio_vector_potential