interp_fields Subroutine

public subroutine interp_fields(params, prtcls, F)

Arguments

Type IntentOptional AttributesName
type(KORC_PARAMS), intent(in) :: params

Core KORC simulation parameters.

type(PARTICLES), intent(inout) :: prtcls

An instance of PARTICLES containing the variables of a given species.

type(FIELDS), intent(in) :: F

An instance of KORC's derived type FIELDS containing all the information about the fields used in the simulation. See korc_types and korc_fields.


Calls

proc~~interp_fields~~CallsGraph proc~interp_fields interp_fields proc~interp_3d_bfields interp_3D_bfields proc~interp_fields->proc~interp_3d_bfields proc~interp_2d_efields interp_2D_efields proc~interp_fields->proc~interp_2d_efields proc~interp_fofields_mars interp_FOfields_mars proc~interp_fields->proc~interp_fofields_mars proc~check_if_in_lcfs check_if_in_LCFS proc~interp_fields->proc~check_if_in_lcfs proc~interp_2d_curlbfields interp_2D_curlbfields proc~interp_fields->proc~interp_2d_curlbfields proc~calculate_magnetic_field calculate_magnetic_field proc~interp_fields->proc~calculate_magnetic_field proc~get_fio_magnetic_fields get_fio_magnetic_fields proc~interp_fields->proc~get_fio_magnetic_fields proc~korc_abort korc_abort proc~interp_fields->proc~korc_abort proc~interp_2d_bfields interp_2D_bfields proc~interp_fields->proc~interp_2d_bfields proc~get_fio_electric_fields get_fio_electric_fields proc~interp_fields->proc~get_fio_electric_fields proc~get_fio_vector_potential get_fio_vector_potential proc~interp_fields->proc~get_fio_vector_potential proc~gradient_2d_bfields gradient_2D_Bfields proc~interp_fields->proc~gradient_2d_bfields proc~check_if_in_fields_domain check_if_in_fields_domain proc~interp_fields->proc~check_if_in_fields_domain proc~interp_3d_efields interp_3D_efields proc~interp_fields->proc~interp_3d_efields proc~interp_2d_gradbfields interp_2D_gradBfields proc~interp_fields->proc~interp_2d_gradbfields proc~interp_fofields_aorsa interp_FOfields_aorsa proc~interp_fields->proc~interp_fofields_aorsa ezspline_error ezspline_error proc~interp_3d_bfields->ezspline_error proc~interp_2d_efields->ezspline_error proc~interp_fofields_mars->ezspline_error ezspline_interp ezspline_interp proc~interp_fofields_mars->ezspline_interp proc~interp_2d_curlbfields->ezspline_error proc~calculate_magnetic_field->ezspline_error interface~fio_eval_field fio_eval_field proc~get_fio_magnetic_fields->interface~fio_eval_field omp_get_thread_num omp_get_thread_num proc~get_fio_magnetic_fields->omp_get_thread_num mpi_abort mpi_abort proc~korc_abort->mpi_abort proc~interp_2d_bfields->ezspline_error proc~get_fio_electric_fields->interface~fio_eval_field proc~get_fio_vector_potential->interface~fio_eval_field proc~gradient_2d_bfields->ezspline_error proc~interp_3d_efields->ezspline_error proc~interp_2d_gradbfields->ezspline_error proc~interp_fofields_aorsa->ezspline_error proc~interp_fofields_aorsa->ezspline_interp

Called by

proc~~interp_fields~~CalledByGraph proc~interp_fields interp_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


Source Code

subroutine interp_fields(params,prtcls,F)
  !! @note Subroutine that works as an interface for calling the 
  !! appropriate subroutines for interpolating or calculating the 
  !! electric and magnetic fields. @endnote
  TYPE(KORC_PARAMS), INTENT(IN)      :: params
    !! Core KORC simulation parameters.
  TYPE(PARTICLES), INTENT(INOUT) :: prtcls
    !! An instance of PARTICLES containing the variables of a given species.
  TYPE(FIELDS), INTENT(IN)       :: F
    !! An instance of KORC's derived type FIELDS containing all the 
    !! information about the fields used in the simulation.
  !! See [[korc_types]] and [[korc_fields]].
  integer :: pp
  
  if (.not.params%GC_coords) call cart_to_cyl(prtcls%X,prtcls%Y)
  
!  write(output_unit_write,'("BR: ",E17.10)') prtcls%BR(:,1)
!  write(output_unit_write,'("Y: ",E17.10)') prtcls%X(2,1)
!  write(output_unit_write,'("Z: ",E17.10)') prtcls%X(3,1)

#ifdef PSPLINE
  call check_if_in_fields_domain(F,prtcls%Y, prtcls%flagCon)

  call check_if_in_LCFS(F,prtcls%Y, prtcls%initLCFS)
#endif
  
  !write(output_unit_write,*) 'checked domain'

#ifdef FIO
  if (TRIM(params%field_model) .eq. 'M3D_C1'.or. &
       TRIM(params%field_model) .eq. 'NIMROD') then

     if (F%FIO_B .ge. 0) then

        !write(6,*) 'interp_fields'

        call get_fio_magnetic_fields(prtcls, F, params)

     end if

     if (F%FIO_E .ge. 0) then
        call get_fio_electric_fields(prtcls, F, params)
     end if

     if (F%FIO_A .ge. 0) then
        call get_fio_vector_potential(prtcls, F, params)
     end if

     do pp=1,sizeof(prtcls%flagCon)
        if (prtcls%flagCon(pp)==0.and. &
             (.not.(params%restart.OR.params%proceed))) then
           write(6,*) 'RE initialized outside of computational domain!!!'
           call KORC_ABORT(15)
        end if
     end do

  end if
#endif

#ifdef PSPLINE
  if (params%field_model(10:13).eq.'MARS') then
     call interp_FOfields_mars(prtcls, F, params)
  end if

  if (params%field_model(10:14).eq.'AORSA') then
     call interp_FOfields_aorsa(prtcls, F, params)
  end if

  if ((ALLOCATED(F%PSIp).and.F%Bflux).or. &
       (F%ReInterp_2x1t.and.(.not.((TRIM(params%field_model).eq.'M3D_C1').or. &
       (params%field_model(10:13).eq.'MARS').or. &
       (TRIM(params%field_model).eq.'NIMROD'))))) then

     !     write(output_unit_write,'("3 size of PSI_P: ",I16)') size(prtcls%PSI_P)

     !     write(output_unit_write,'("B_X: ",E17.10)') prtcls%B(:,1)
     !     write(output_unit_write,'("B_Z: ",E17.10)') prtcls%B(:,3)
     !     write(output_unit_write,'("B_Y: ",E17.10)') prtcls%B(:,2)
     !     write(output_unit_write,'("PSI_P: ",E17.10)') prtcls%PSI_P

     call calculate_magnetic_field(params,prtcls%Y,F,prtcls%B,prtcls%E, &
          prtcls%PSI_P,prtcls%flagCon)

     !write(output_unit_write,*) 'interp PSIp'

     !     write(output_unit_write,'("interp_fields")')
     !     write(output_unit_write,'("B_X: ",E17.10)') prtcls%B(:,1)
     !     write(output_unit_write,'("B_Z: ",E17.10)') prtcls%B(:,3)
     !     write(output_unit_write,'("B_Y: ",E17.10)') prtcls%B(:,2)

  end if

  if (ALLOCATED(F%PSIp3D).and.F%Bflux3D) then
     
!     write(output_unit_write,'("3 size of PSI_P: ",I16)') size(prtcls%PSI_P)

!     write(output_unit_write,'("B_X: ",E17.10)') prtcls%B(:,1)
!     write(output_unit_write,'("B_Z: ",E17.10)') prtcls%B(:,3)
!     write(output_unit_write,'("B_Y: ",E17.10)') prtcls%B(:,2)
!     write(output_unit_write,'("PSI_P: ",E17.10)') prtcls%PSI_P
     
     call calculate_magnetic_field(params,prtcls%Y,F,prtcls%B,prtcls%E, &
          prtcls%PSI_P,prtcls%flagCon)

!     write(output_unit_write,'("interp_fields")')
!     write(output_unit_write,'("B_X: ",E17.10)') prtcls%B(:,1)
!     write(output_unit_write,'("B_Z: ",E17.10)') prtcls%B(:,3)
!     write(output_unit_write,'("B_Y: ",E17.10)') prtcls%B(:,2)

  end if
  
  if (ALLOCATED(F%B_2D%R).and.F%Bfield) then
     call interp_2D_bfields(params,prtcls%Y,prtcls%B,prtcls%flagCon)
  end if

  if (ALLOCATED(F%B_3D%R).and.F%Bfield) then
     call interp_3D_bfields(params,prtcls%Y,prtcls%B,prtcls%flagCon)
  end if

!  if (ALLOCATED(F%E_2D%R).and.F%Efield) then
!     call interp_2D_efields(params,prtcls%Y,prtcls%E,prtcls%flagCon)
!  end if

  if (ALLOCATED(F%E_3D%R).and.F%Efield.and.(.not.F%ReInterp_2x1t)) then
     call interp_3D_efields(params,prtcls%Y,prtcls%E,prtcls%flagCon)
     
  end if

  if (ALLOCATED(F%E_3D%R).and.F%Efield.and.F%ReInterp_2x1t) then
     call interp_2D_efields(params,prtcls%Y,prtcls%E,prtcls%flagCon)

!     write(output_unit_write,*) 'interpolated efield'
     
  end if
  
  if (params%GC_coords.and.ALLOCATED(F%gradB_2D%R).and.F%Bfield) then
     call interp_2D_gradBfields(prtcls%Y,prtcls%gradB,prtcls%flagCon)

  end if

  if (params%GC_coords.and.ALLOCATED(F%gradB_2D%R).and.F%Bfield) then
     call interp_2D_curlbfields(prtcls%Y,prtcls%curlb,prtcls%flagCon)
  end if

  if(params%GC_coords.and.params%orbit_model(3:6)=='grad') then
     call gradient_2D_bfields(prtcls%Y,prtcls%BR,prtcls%BPHI, &
          prtcls%BZ,prtcls%flagCon)
  end if
#endif
  
end subroutine interp_fields