get_fio_magnetic_fields Subroutine

public subroutine get_fio_magnetic_fields(prtcls, F, params)

Uses

  • proc~~get_fio_magnetic_fields~~UsesGraph proc~get_fio_magnetic_fields get_fio_magnetic_fields omp_lib omp_lib proc~get_fio_magnetic_fields->omp_lib

Arguments

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

Calls

proc~~get_fio_magnetic_fields~~CallsGraph proc~get_fio_magnetic_fields get_fio_magnetic_fields omp_get_thread_num omp_get_thread_num proc~get_fio_magnetic_fields->omp_get_thread_num interface~fio_eval_field fio_eval_field proc~get_fio_magnetic_fields->interface~fio_eval_field

Called by

proc~~get_fio_magnetic_fields~~CalledByGraph proc~get_fio_magnetic_fields get_fio_magnetic_fields proc~interp_fields interp_fields proc~interp_fields->proc~get_fio_magnetic_fields 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_magnetic_fields(prtcls, F, params)
  USE omp_lib
  IMPLICIT NONE
  
    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)         :: Btmp
    TYPE(C_PTR), DIMENSION(size(prtcls%hint)) :: hint
    INTEGER             :: thread_num

!    write(output_unit_write,*) 'in m3dc1 B'
    
    if (prtcls%cart) then
       !$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(pp,status,x)
       do pp = 1, SIZE(prtcls%hint)
          if (prtcls%flagCon(pp) .EQ. 1_is) then
             x = prtcls%X(pp,:)*params%cpp%length
             status = fio_eval_field(F%FIO_B, x(1),   &
                  prtcls%B(pp,1),                        &
                  prtcls%hint(pp))

             if (status .eq. FIO_NO_DATA) then
                prtcls%B(pp,:) = 0
                prtcls%flagCon(pp) = 0_is
             else if (status .ne. FIO_SUCCESS) then
                prtcls%flagCon(pp) = 0_is
             end if
          end if
       end do
       !$OMP END PARALLEL DO
    else

!       write(output_unit_write,*) 'in cart false'
       !hint=prtcls%hint
       !write(output_unit_write,*) 'hint: ',hint

       Btmp=0._rp
       
       !$OMP PARALLEL DO DEFAULT(none) &
       !$OMP& SHARED(prtcls,params,F) &
       !$OMP& PRIVATE(pp,status,x,thread_num) &
       !$OMP& FIRSTPRIVATE(Btmp)
       do pp = 1, SIZE(prtcls%hint)

          thread_num = OMP_GET_THREAD_NUM()
          
          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(6,*) 'thread',thread_num,'X',x
             
             !             prtcls%hint(pp)=c_null_ptr

             !write(output_unit_write,*) 'thread',thread_num,'before interpolating B'
             
             status = fio_eval_field(F%FIO_B, x(1),                      &
                  Btmp(1),prtcls%hint(pp))

                     
             if (status .eq. FIO_NO_DATA) then
                prtcls%B(pp,:) = 0
                prtcls%flagCon(pp) = 0_is
             else if (status .ne. FIO_SUCCESS) then
                prtcls%flagCon(pp) = 0_is
             end if
            
             
             if (.not.params%GC_coords) then
                
                prtcls%B(pp,1)=(Btmp(1)*cos(x(2))-Btmp(2)*sin(x(2)))/ &
                     params%cpp%Bo
                prtcls%B(pp,2)=(Btmp(1)*sin(x(2))+Btmp(2)*cos(x(2)))/ &
                     params%cpp%Bo
                prtcls%B(pp,3)=Btmp(3)/params%cpp%Bo
                
             else
                
                prtcls%B(pp,1)=Btmp(1)/params%cpp%Bo
                prtcls%B(pp,2)=Btmp(2)/params%cpp%Bo
                prtcls%B(pp,3)=Btmp(3)/params%cpp%Bo
             end if
             
             
          end if
       end do
       !$OMP END PARALLEL DO

    end if
  end subroutine get_fio_magnetic_fields