get_fio_profile Subroutine

public subroutine get_fio_profile(prtcls, P, params)

Arguments

Type IntentOptional AttributesName
type(PARTICLES), intent(inout) :: prtcls
type(PROFILES), intent(in) :: P
type(KORC_PARAMS), intent(in) :: params

Calls

proc~~get_fio_profile~~CallsGraph proc~get_fio_profile get_fio_profile interface~fio_eval_field fio_eval_field proc~get_fio_profile->interface~fio_eval_field

Called by

proc~~get_fio_profile~~CalledByGraph proc~get_fio_profile get_fio_profile proc~interp_profiles interp_profiles proc~interp_profiles->proc~get_fio_profile proc~get_profiles get_profiles proc~get_profiles->proc~interp_profiles proc~avalanche_4d Avalanche_4D proc~avalanche_4d->proc~get_profiles proc~get_avalanche_4d get_Avalanche_4D proc~get_avalanche_4d->proc~avalanche_4d

Contents

Source Code


Source Code

  subroutine get_fio_profile(prtcls, P, params)
    TYPE(PARTICLES), INTENT(INOUT) :: prtcls
    TYPE(PROFILES), INTENT(IN)     :: P
    TYPE(KORC_PARAMS), INTENT(IN)  :: params
    INTEGER (C_INT)                :: status
    INTEGER                        :: pp
    REAL(rp), DIMENSION(3)         :: x
    integer(ip)  ::  ss
    REAL(rp)        :: netmp=-1._rp,Tetmp=-1._rp

    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

    if (prtcls%cart) then
       !$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(pp,status,x)
       do pp = 1,ss
          if (prtcls%flagCon(pp) .EQ. 1_is) then
             x = prtcls%X(pp,:)*params%cpp%length
             status = fio_eval_field(P%FIO_ne, x(1),                     &
                  netmp,                         &
                  prtcls%hint(pp))

             if (status .eq. FIO_SUCCESS) then

                if(netmp.le.0) netmp=1._rp

                prtcls%ne(pp) = netmp/params%cpp%density
                
             else if (status .eq. FIO_NO_DATA) then
                prtcls%ne(pp) = 0
             else if (status .ne. FIO_SUCCESS) then
                prtcls%flagCon(pp) = 0_is
                CYCLE
             end if

             status = fio_eval_field(P%FIO_te, x(1),                     &
                  Tetmp,                         &
                  prtcls%hint(pp))

             if (status .eq. FIO_SUCCESS) then
                if(Tetmp.le.0) Tetmp=0.1_rp
             
                prtcls%te(pp) = Tetmp/(params%cpp%temperature/C_E)
             else if (status .eq. FIO_NO_DATA) then
                prtcls%te(pp) = 0
             end if

             status = fio_eval_field(P%FIO_zeff, x(1), &
                  prtcls%Zeff(pp),prtcls%hint(pp))

             if (status .eq. FIO_NO_DATA) then
                prtcls%Zeff(pp) = 1
             end if
          end if
       end do
       !$OMP END PARALLEL DO
    else
       !$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(pp,status,x)
       do pp = 1,ss
          if (prtcls%flagCon(pp) .EQ. 1_is) then
             x(1) = prtcls%Y(1,pp)*params%cpp%length
             x(2) = prtcls%Y(2,pp)
             x(3) = prtcls%Y(3,pp)*params%cpp%length
             status = fio_eval_field(P%FIO_ne, x(1), &
                  netmp,prtcls%hint(pp))

             if (status .eq. FIO_SUCCESS) then

                if(netmp.le.0) netmp=1._rp

                prtcls%ne(pp) = netmp/params%cpp%density
                
             else if (status .eq. FIO_NO_DATA) then
                prtcls%ne(pp) = 0
             else if (status .ne. FIO_SUCCESS) then
                prtcls%flagCon(pp) = 0_is
                CYCLE
             end if

             status = fio_eval_field(P%FIO_te, x(1),&
                  tetmp,prtcls%hint(pp))

             if (status .eq. FIO_SUCCESS) then
                if(Tetmp.le.0) Tetmp=0.1_rp
             
                prtcls%te(pp) = Tetmp/(params%cpp%temperature/C_E)
             else if (status .eq. FIO_NO_DATA) then
                prtcls%te(pp) = 0
             end if

             status = fio_eval_field(P%FIO_zeff, x(1), &
                  prtcls%Zeff(pp),prtcls%hint(pp))

             if (status .eq. FIO_NO_DATA) then
                prtcls%Zeff(pp) = 1
             end if
          end if
       end do
       !$OMP END PARALLEL DO
    end if
  end subroutine get_fio_profile