get_fio_FOelectric_fields_p Subroutine

public subroutine get_fio_FOelectric_fields_p(params, F, Y_R, Y_PHI, Y_Z, E_X, E_Y, E_Z, flag, hint)

Arguments

Type IntentOptional AttributesName
type(KORC_PARAMS), intent(in) :: params
type(FIELDS), intent(in) :: F
real(kind=rp), intent(in), DIMENSION(params%pchunk):: Y_R
real(kind=rp), intent(in), DIMENSION(params%pchunk):: Y_PHI
real(kind=rp), intent(in), DIMENSION(params%pchunk):: Y_Z
real(kind=rp), intent(out), DIMENSION(params%pchunk):: E_X
real(kind=rp), intent(out), DIMENSION(params%pchunk):: E_Y
real(kind=rp), intent(out), DIMENSION(params%pchunk):: E_Z
integer(kind=is), intent(inout), DIMENSION(params%pchunk):: flag
type(C_PTR), intent(inout), DIMENSION(params%pchunk):: hint

Calls

proc~~get_fio_foelectric_fields_p~~CallsGraph proc~get_fio_foelectric_fields_p get_fio_FOelectric_fields_p interface~fio_eval_field fio_eval_field proc~get_fio_foelectric_fields_p->interface~fio_eval_field

Called by

proc~~get_fio_foelectric_fields_p~~CalledByGraph proc~get_fio_foelectric_fields_p get_fio_FOelectric_fields_p proc~adv_fofio_top adv_FOfio_top proc~adv_fofio_top->proc~get_fio_foelectric_fields_p proc~fo_init FO_init proc~fo_init->proc~get_fio_foelectric_fields_p program~main main program~main->proc~adv_fofio_top program~main->proc~fo_init

Contents


Source Code

  subroutine get_fio_FOelectric_fields_p(params,F,Y_R,Y_PHI,Y_Z, &
       E_X,E_Y,E_Z,flag,hint)
    TYPE(FIELDS), INTENT(IN)       :: F
    TYPE(KORC_PARAMS), INTENT(IN)  :: params
    REAL(rp), DIMENSION(params%pchunk), INTENT(IN)  :: Y_R,Y_PHI,Y_Z
    REAL(rp), DIMENSION(params%pchunk), INTENT(OUT)  :: E_X,E_Y,E_Z
    INTEGER(is), DIMENSION(params%pchunk), INTENT(INOUT)  :: flag
    TYPE(C_PTR), DIMENSION(params%pchunk), INTENT(INOUT)  :: hint
    INTEGER (C_INT)                :: status
    INTEGER                        :: pp,pchunk
    REAL(rp), DIMENSION(3)         :: x
    REAL(rp), DIMENSION(3)         :: Etmp

    pchunk=params%pchunk

    do pp = 1,pchunk
       if (flag(pp) .EQ. 1_is) then
          x(1) = Y_R(pp)*params%cpp%length
          x(2) = Y_PHI(pp)
          x(3) = Y_Z(pp)*params%cpp%length

          !             prtcls%hint(pp)=c_null_ptr


          status = fio_eval_field(F%FIO_E, x(1),                      &
               Etmp(1),hint(pp))

          if (status .eq. FIO_SUCCESS) then
             E_X(pp)=(Etmp(1)*cos(x(2))-Etmp(2)*sin(x(2)))/ &
                  params%cpp%Eo
             E_Y(pp)=(Etmp(1)*sin(x(2))+Etmp(2)*cos(x(2)))/ &
                  params%cpp%Eo
             E_Z(pp)=Etmp(3)/params%cpp%Eo     
          else if (status .eq. FIO_NO_DATA) then
             E_X(pp) = 0
             E_Y(pp) = 0
             E_Z(pp) = 0                
             flag(pp) = 0_is
          else if (status .ne. FIO_SUCCESS) then
             flag(pp) = 0_is
          end if
         
       end if
    end do

  end subroutine get_fio_FOelectric_fields_p