calculate_SC_E1D_FS Subroutine

public subroutine calculate_SC_E1D_FS(params, F, dintJphidPSIP)

Arguments

Type IntentOptional AttributesName
type(KORC_PARAMS), intent(in) :: params
type(FIELDS), intent(inout) :: F
real(kind=rp), intent(in), dimension(F%dim_1D):: dintJphidPSIP

Calls

proc~~calculate_sc_e1d_fs~~CallsGraph proc~calculate_sc_e1d_fs calculate_SC_E1D_FS mpi_allreduce mpi_allreduce proc~calculate_sc_e1d_fs->mpi_allreduce proc~initialize_sc1d_field_interpolant_fs initialize_SC1D_field_interpolant_FS proc~calculate_sc_e1d_fs->proc~initialize_sc1d_field_interpolant_fs ezspline_init ezspline_init proc~initialize_sc1d_field_interpolant_fs->ezspline_init ezspline_allocated ezspline_allocated proc~initialize_sc1d_field_interpolant_fs->ezspline_allocated ezspline_error ezspline_error proc~initialize_sc1d_field_interpolant_fs->ezspline_error ezspline_setup ezspline_setup proc~initialize_sc1d_field_interpolant_fs->ezspline_setup

Called by

proc~~calculate_sc_e1d_fs~~CalledByGraph proc~calculate_sc_e1d_fs calculate_SC_E1D_FS proc~adv_gcinterp_psi_top_fs adv_GCinterp_psi_top_FS proc~adv_gcinterp_psi_top_fs->proc~calculate_sc_e1d_fs program~main main program~main->proc~adv_gcinterp_psi_top_fs

Contents

Source Code


Source Code

  subroutine calculate_SC_E1D_FS(params,F,dintJphidPSIP)

    TYPE(FIELDS), INTENT(INOUT)                 :: F
    TYPE(KORC_PARAMS), INTENT(IN) 		:: params
    real(rp),dimension(F%dim_1D),intent(in) :: dintJphidPSIP
    real(rp),dimension(F%dim_1D) :: Jsamall,Jexp,dJdt
    real(rp),dimension(F%dim_1D) :: a,b,c,u,gam,r,alpha,beta,gamma
    real(rp) :: dPSIP,bet
    integer :: ii
    INTEGER 				:: mpierr

!    if (params%mpi_params%rank .EQ. 0) then
!       write(output_unit_write,*) 'Calculating SC_E1D'
!    end if

    !write(output_unit_write,*) 'dintJphidPSIP',dintJphidPSIP(F%dim_1D)
    
    dPSIP=F%PSIP_1D(2)-F%PSIP_1D(1)
    

    ! Add sampled current densities from all MPI processes Jsamone,
    ! and output of total sampled current density Jsamall to each
    ! MPI process.
    
    call MPI_ALLREDUCE(dintJphidPSIP,Jsamall,F%dim_1D,MPI_REAL8,MPI_SUM, &
         MPI_COMM_WORLD,mpierr)

    !write(output_unit_write,*) 'JSamAll',Jsamall(F%dim_1D)
    
    !write(output_unit_write,*) 'Jsam: ',Jsamall(1:5)
    
    Jexp=Jsamall*F%Ip0
    
    F%J3_SC_1D%PHI=F%J2_SC_1D%PHI
    F%J2_SC_1D%PHI=F%J1_SC_1D%PHI
    F%J1_SC_1D%PHI=Jexp    
    
    ! Calculating time-derivative of E_phi
    
    dJdt=(3*F%J1_SC_1D%PHI-4*F%J2_SC_1D%PHI+F%J3_SC_1D%PHI)/ &
         (2*F%dt_E_SC)
    
!    write(output_unit_write,*) params%mpi_params%rank,'J(1)',F%J_SC_1D%PHI(1)
    
    ! Solving 1D Poisson equation with tridiagonal matrix solve

    alpha=F%ddMagPsiSqdPsiPSq
    beta=F%dMagPsiSqdPsiP
    gamma=C_MU*dJdt
   
    
    a=-alpha*dPSIP/2._rp+beta
    b=-2._rp*beta
    c=alpha*dPSIP/2._rp+beta
    u=0._rp
    gam=0._rp
!    r=-2*dr**2*C_MU*Jexp
    r=dPSIP**2*gamma

    c(2)=c(2)-a(2)*a(1)/c(1)
    b(2)=b(2)-a(2)*b(1)/c(1)
    r(2)=r(2)-a(2)*r(1)/c(1)
    
    bet=b(2)
    u(2)=r(2)/bet
    do ii=3_idef,F%dim_1D-1
       gam(ii)=c(ii-1)/bet
       bet=b(ii)-a(ii)*gam(ii)
       if (bet.eq.0) then
          stop 'tridiag failed'         
       end if
       u(ii)=(r(ii)-a(ii)*u(ii-1))/bet
    end do

    do ii=F%dim_1D-2,2,-1
       u(ii)=u(ii)-gam(ii+1)*u(ii+1)
    end do

    u(1)=2*u(2)-u(3)

    F%E_SC_1D%PHI=u

    if (params%mpi_params%rank.eq.0) then
       write(output_unit_write,*) 'J1(1)',F%J1_SC_1D%PHI(1)
       write(output_unit_write,*) 'J2(1)',F%J2_SC_1D%PHI(1)
       write(output_unit_write,*) 'J3(1)',F%J3_SC_1D%PHI(1)
       
       write(output_unit_write,*) 'E(1)',F%E_SC_1D%PHI(1)
    end if
       
    ! Normalizing inductive E_phi
    
    F%E_SC_1D%PHI=F%E_SC_1D%PHI/params%cpp%Eo

#ifdef PSPLINE
    call initialize_SC1D_field_interpolant_FS(params,F)
#endif
    
  end subroutine calculate_SC_E1D_FS