load_particles_ic Subroutine

public subroutine load_particles_ic(params, spp, F)

Arguments

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

Core KORC simulation parameters.

type(SPECIES), intent(inout), DIMENSION(:), ALLOCATABLE:: spp

An instance of KORC's derived type SPECIES containing all the information of different electron species. See korc_types.f90.

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

Calls

proc~~load_particles_ic~~CallsGraph proc~load_particles_ic load_particles_ic mpi_barrier mpi_barrier proc~load_particles_ic->mpi_barrier h5fclose_f h5fclose_f proc~load_particles_ic->h5fclose_f proc~korc_abort korc_abort proc~load_particles_ic->proc~korc_abort interface~load_array_from_hdf5 load_array_from_hdf5 proc~load_particles_ic->interface~load_array_from_hdf5 mpi_scatter mpi_scatter proc~load_particles_ic->mpi_scatter h5fopen_f h5fopen_f proc~load_particles_ic->h5fopen_f mpi_abort mpi_abort proc~korc_abort->mpi_abort proc~rload_2d_array_from_hdf5 rload_2d_array_from_hdf5 interface~load_array_from_hdf5->proc~rload_2d_array_from_hdf5 proc~rload_3d_array_from_hdf5 rload_3d_array_from_hdf5 interface~load_array_from_hdf5->proc~rload_3d_array_from_hdf5 proc~rload_1d_array_from_hdf5 rload_1d_array_from_hdf5 interface~load_array_from_hdf5->proc~rload_1d_array_from_hdf5 proc~rload_2d_array_from_hdf5->proc~korc_abort h5dclose_f h5dclose_f proc~rload_2d_array_from_hdf5->h5dclose_f h5dread_f h5dread_f proc~rload_2d_array_from_hdf5->h5dread_f proc~rload_3d_array_from_hdf5->proc~korc_abort proc~rload_3d_array_from_hdf5->h5dclose_f proc~rload_3d_array_from_hdf5->h5dread_f proc~rload_1d_array_from_hdf5->proc~korc_abort proc~rload_1d_array_from_hdf5->h5dclose_f proc~rload_1d_array_from_hdf5->h5dread_f

Called by

proc~~load_particles_ic~~CalledByGraph proc~load_particles_ic load_particles_ic proc~set_up_particles_ic set_up_particles_ic proc~set_up_particles_ic->proc~load_particles_ic program~main main program~main->proc~set_up_particles_ic

Contents

Source Code


Source Code

  subroutine load_particles_ic(params,spp,F)
    !! @note Subroutine that loads all the electrons' data from
    !! "restart_file.h5" to restart a simulation.
    TYPE(KORC_PARAMS), INTENT(INOUT) 			:: params
    !! Core KORC simulation parameters.
    TYPE(SPECIES), DIMENSION(:), ALLOCATABLE, INTENT(INOUT) :: spp
    !! An instance of KORC's derived type SPECIES containing all the
    !! information of different electron species. See korc_types.f90.
    TYPE(FIELDS), INTENT(INOUT) 			:: F
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: X_send_buffer,X_send_buffer_tmp,X_send_buffer_tmp1
    !! Temporary buffer used by MPI for scattering the electrons' position
    !! to different MPI processes.
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: X_receive_buffer
    !! Temporary buffer used by MPI for scattering the electrons' position
    !! among MPI processes.
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: V_send_buffer,V_send_buffer_tmp,V_send_buffer_tmp1
    !! Temporary buffer used by MPI for scattering the electrons' velocity
    !! among MPI processes.
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: V1_send_buffer_tmp,V2_send_buffer_tmp,V3_send_buffer_tmp
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: V1_send_buffer_tmp1,V2_send_buffer_tmp1,V3_send_buffer_tmp1
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: V_receive_buffer
    !! Temporary buffer used by MPI for scattering the electrons' velocity
    !! among MPI processes.
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: AUX_send_buffer
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: AUX_send_buffer_tmp
    !!  Temporary buffer used by MPI to scatter various electrons' variables
    !! among MPI processes.
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: AUX_receive_buffer
    !! Temporary buffer used by MPI to scatter various electrons' variables
    !! among MPI processes.
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: JSC0_buffer
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: JSC1_buffer
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: JSC2_buffer
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: JSC3_buffer
    REAL(rp), DIMENSION(:), ALLOCATABLE 		:: ESC_buffer


    CHARACTER(MAX_STRING_LENGTH) 			:: filename
    !! String containing the name of the HDF5 file.
    CHARACTER(MAX_STRING_LENGTH) 			:: dset
    !! Name of data set to be saved to file.
    INTEGER(HID_T) 					:: h5file_id
    !! HDF5 file identifier.
    CHARACTER(19) 					:: tmp_str
    !! Temporary string used to manipulate various strings.
    INTEGER 						:: h5error
    !! HDF5 error status.
    INTEGER 						:: mpierr
    !! Electron species iterator.
    INTEGER 						:: ss,ii,jj
    !! MPI error status.
    INTEGER  :: recieve_num,send_num,nmpi_ratio
    
    do ss=1_idef,params%num_species

       !write(6,*) params%mpi_params%nmpi,params%mpi_params%nmpi_prev
       
       if (mod(spp(ss)%ppp,params%mpi_params%nmpi/ &
            params%mpi_params%nmpi_prev).ne.0) then
          write(output_unit_write,'("ppp must be divisible by factor increase in MPI nodes")')
          call KORC_ABORT(14)
       endif

       nmpi_ratio=params%mpi_params%nmpi/params%mpi_params%nmpi_prev
       
       send_num=spp(ss)%ppp*params%mpi_params%nmpi_prev
       recieve_num=spp(ss)%ppp/nmpi_ratio
       
       ALLOCATE(X_send_buffer(3*send_num))
       ALLOCATE(X_send_buffer_tmp(3*send_num))
       ALLOCATE(X_send_buffer_tmp1(3*send_num))
       ALLOCATE(X_receive_buffer(3*recieve_num))

       ALLOCATE(V_send_buffer(3*send_num))
       ALLOCATE(V_send_buffer_tmp(3*send_num))
       ALLOCATE(V_send_buffer_tmp1(3*send_num))
       ALLOCATE(V_receive_buffer(3*recieve_num))

       ALLOCATE(V1_send_buffer_tmp(send_num))
       ALLOCATE(V2_send_buffer_tmp(send_num))
       ALLOCATE(V3_send_buffer_tmp(send_num))
       ALLOCATE(V1_send_buffer_tmp1(send_num))
       ALLOCATE(V2_send_buffer_tmp1(send_num))
       ALLOCATE(V3_send_buffer_tmp1(send_num))
       
       ALLOCATE(AUX_send_buffer(send_num))
       ALLOCATE(AUX_send_buffer_tmp(send_num))
       ALLOCATE(AUX_receive_buffer(recieve_num))
       
       if (params%mpi_params%rank.EQ.0_idef) then
          
          filename = TRIM(params%path_to_outputs) // "restart_file.h5"
          call h5fopen_f(filename, H5F_ACC_RDONLY_F, h5file_id, h5error)
          if (h5error .EQ. -1) then
             write(output_unit_write,'("KORC ERROR: Something went wrong in: &
                  &load_particles_ic --> h5fopen_f")')
             call KORC_ABORT(14)
          end if

          write(tmp_str,'(I18)') ss

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/X"
          call load_array_from_hdf5(h5file_id,dset,X_send_buffer)

          call h5fclose_f(h5file_id, h5error)

          !write(6,*) 'shape of X_send_buffer',shape(X_send_buffer)
          !write(6,*) 'X_send_buffer',X_send_buffer
          !write(6,*) 'recieve_num',recieve_num

#if DBG_CHECK
          write(6,*) 'X_send_buffer',X_send_buffer*params%cpp%length
#endif
          
          if (params%load_balance) then
             do ii=0,params%mpi_params%nmpi_prev-1

                V1_send_buffer_tmp(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)= &
                     X_send_buffer(3*ii*spp(ss)%ppp+1:(3*ii+1)*spp(ss)%ppp)
                V2_send_buffer_tmp(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)= &
                     X_send_buffer((3*ii+1)*spp(ss)%ppp+1:(3*ii+2)*spp(ss)%ppp)
                V3_send_buffer_tmp(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)= &
                     X_send_buffer((3*ii+2)*spp(ss)%ppp+1:(3*ii+3)*spp(ss)%ppp)
                
             end do
             do ii=1,params%mpi_params%nmpi_prev*spp(ss)%ppp

                V1_send_buffer_tmp1((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     V1_send_buffer_tmp(ii)
                V2_send_buffer_tmp1((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     V2_send_buffer_tmp(ii)
                V3_send_buffer_tmp1((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     V3_send_buffer_tmp(ii)
                                                
             end do
             do ii=0,params%mpi_params%nmpi_prev-1

                X_send_buffer_tmp1(3*ii*spp(ss)%ppp+1:(3*ii+1)*spp(ss)%ppp)= &
                     V1_send_buffer_tmp1(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)
                X_send_buffer_tmp1((3*ii+1)*spp(ss)%ppp+1:(3*ii+2)*spp(ss)%ppp)= &
                     V2_send_buffer_tmp1(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)
                X_send_buffer_tmp1((3*ii+2)*spp(ss)%ppp+1:(3*ii+3)*spp(ss)%ppp)= &
                     V3_send_buffer_tmp1(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)


             end do
          else
             X_send_buffer_tmp1=X_send_buffer
          end if
          
          do jj=0_idef,params%mpi_params%nmpi_prev-1_idef
             do ii=0_idef,nmpi_ratio-1_idef

                X_send_buffer_tmp(jj*spp(ss)%ppp*3+recieve_num*(3*ii)+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(3*ii+1))= &
                     X_send_buffer_tmp1(jj*spp(ss)%ppp*3+recieve_num*ii+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(ii+1))

                X_send_buffer_tmp(jj*spp(ss)%ppp*3+recieve_num*(3*ii+1)+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(3*ii+2))=&
                     X_send_buffer_tmp1(spp(ss)%ppp*(3*jj+1)+recieve_num*ii+1: &
                     spp(ss)%ppp*(3*jj+1)+recieve_num*(ii+1))

                X_send_buffer_tmp(jj*spp(ss)%ppp*3+recieve_num*(3*ii+2)+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(3*ii+3))=&
                     X_send_buffer_tmp1(spp(ss)%ppp*(3*jj+2)+recieve_num*ii+1: &
                     spp(ss)%ppp*(3*jj+2)+recieve_num*(ii+1))
             end do
          end do

#if DBG_CHECK
          write(6,*) 'X_send_buffer_tmp1',X_send_buffer_tmp1*params%cpp%length
          write(6,*) 'X_send_buffer_tmp',X_send_buffer_tmp*params%cpp%length
#endif

          
       end if

       X_receive_buffer = 0.0_rp
       CALL MPI_SCATTER(X_send_buffer_tmp,3*recieve_num,MPI_REAL8, &
            X_receive_buffer,3*recieve_num,MPI_REAL8,0,MPI_COMM_WORLD,mpierr)
       if (params%orbit_model(1:2).EQ.'FO') then  
          spp(ss)%vars%X(1:recieve_num,:) = &
               RESHAPE(X_receive_buffer,(/recieve_num,3/))
          !set dummy values for new particles, so EZspline doesn't throw errors

          do ii=1_idef,spp(ss)%ppp-recieve_num
             spp(ss)%vars%X(recieve_num+ii,:)=spp(ss)%vars%X(1,:)
          end do
          
       else if (params%orbit_model(1:2).EQ.'GC') then
          spp(ss)%vars%Y(1:recieve_num,:) = &
               RESHAPE(X_receive_buffer,(/recieve_num,3/))

          do ii=1_idef,spp(ss)%ppp-recieve_num
             spp(ss)%vars%Y(recieve_num+ii,:)=spp(ss)%vars%Y(1,:)
          end do
          
          spp(ss)%vars%Y(:,2)=modulo(spp(ss)%vars%Y(:,2),2*C_PI)
       end if

       if (params%mpi_params%rank.EQ.0_idef) then
          filename = TRIM(params%path_to_outputs) // "restart_file.h5"
          call h5fopen_f(filename, H5F_ACC_RDONLY_F, h5file_id, h5error)
          if (h5error .EQ. -1) then
             write(output_unit_write,'("KORC ERROR: Something went wrong in: &
                  &load_particles_ic --> h5fopen_f")')
             call KORC_ABORT(14)
          end if

          write(tmp_str,'(I18)') ss

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/V"
          call load_array_from_hdf5(h5file_id,dset,V_send_buffer)

          call h5fclose_f(h5file_id, h5error)

#if DBG_CHECK
          write(6,*) 'V_send_buffer',V_send_buffer
#endif
          
          if (params%load_balance) then
             do ii=0,params%mpi_params%nmpi_prev-1

                V1_send_buffer_tmp(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)= &
                     V_send_buffer(3*ii*spp(ss)%ppp+1:(3*ii+1)*spp(ss)%ppp)
                V2_send_buffer_tmp(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)= &
                     V_send_buffer((3*ii+1)*spp(ss)%ppp+1:(3*ii+2)*spp(ss)%ppp)
                V3_send_buffer_tmp(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)= &
                     V_send_buffer((3*ii+2)*spp(ss)%ppp+1:(3*ii+3)*spp(ss)%ppp)
                
             end do
             do ii=1,params%mpi_params%nmpi_prev*spp(ss)%ppp

                V1_send_buffer_tmp1((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     V1_send_buffer_tmp(ii)
                V2_send_buffer_tmp1((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     V2_send_buffer_tmp(ii)
                V3_send_buffer_tmp1((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     V3_send_buffer_tmp(ii)
                                                
             end do
             do ii=0,params%mpi_params%nmpi_prev-1

                V_send_buffer_tmp1(3*ii*spp(ss)%ppp+1:(3*ii+1)*spp(ss)%ppp)= &
                     V1_send_buffer_tmp1(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)
                V_send_buffer_tmp1((3*ii+1)*spp(ss)%ppp+1:(3*ii+2)*spp(ss)%ppp)= &
                     V2_send_buffer_tmp1(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)
                V_send_buffer_tmp1((3*ii+2)*spp(ss)%ppp+1:(3*ii+3)*spp(ss)%ppp)= &
                     V3_send_buffer_tmp1(ii*spp(ss)%ppp+1:(ii+1)*spp(ss)%ppp)


             end do
          else
             V_send_buffer_tmp1=V_send_buffer
          end if
          
          
          do jj=0_idef,params%mpi_params%nmpi_prev-1_idef
             do ii=0_idef,nmpi_ratio-1_idef

                V_send_buffer_tmp(jj*spp(ss)%ppp*3+recieve_num*(3*ii)+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(3*ii+1))= &
                     V_send_buffer_tmp1(jj*spp(ss)%ppp*3+recieve_num*ii+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(ii+1))

                V_send_buffer_tmp(jj*spp(ss)%ppp*3+recieve_num*(3*ii+1)+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(3*ii+2))=&
                     V_send_buffer_tmp1(spp(ss)%ppp*(3*jj+1)+recieve_num*ii+1: &
                     spp(ss)%ppp*(3*jj+1)+recieve_num*(ii+1))

                V_send_buffer_tmp(jj*spp(ss)%ppp*3+recieve_num*(3*ii+2)+1: &
                     jj*spp(ss)%ppp*3+recieve_num*(3*ii+3))=&
                     V_send_buffer_tmp1(spp(ss)%ppp*(3*jj+2)+recieve_num*ii+1: &
                     spp(ss)%ppp*(3*jj+2)+recieve_num*(ii+1))
             end do
          end do


#if DBG_CHECK
          write(6,*) 'V_send_buffer_tmp1',V_send_buffer_tmp1
          write(6,*) 'V_send_buffer_tmp',V_send_buffer_tmp
#endif
          
       end if

       

       V_receive_buffer = 0.0_rp
       CALL MPI_SCATTER(V_send_buffer_tmp,3*recieve_num,MPI_REAL8, &
            V_receive_buffer,3*recieve_num,MPI_REAL8,0,MPI_COMM_WORLD,mpierr)
       spp(ss)%vars%V(1:recieve_num,:) = &
            RESHAPE(V_receive_buffer,(/recieve_num,3/))

       do ii=1_idef,spp(ss)%ppp-recieve_num
          spp(ss)%vars%V(recieve_num+ii,:)=spp(ss)%vars%V(1,:)
       end do
       
       if (params%mpi_params%rank.EQ.0_idef) then
          filename = TRIM(params%path_to_outputs) // "restart_file.h5"
          call h5fopen_f(filename, H5F_ACC_RDONLY_F, h5file_id, h5error)
          if (h5error .EQ. -1) then
             write(output_unit_write,'("KORC ERROR: Something went wrong in: &
                  &load_particles_ic --> h5fopen_f")')
             call KORC_ABORT(14)
          end if

          write(tmp_str,'(I18)') ss

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/flagCon"
          call load_array_from_hdf5(h5file_id,dset,AUX_send_buffer)

          call h5fclose_f(h5file_id, h5error)

          if (params%load_balance) then
             do ii=1,params%mpi_params%nmpi_prev*spp(ss)%ppp

                !write(6,*) 'send index',ii
                !write(6,*) 'tmp index',(mod(ii-1,params%mpi_params%nmpi))* &
                !     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1

                AUX_send_buffer_tmp((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     AUX_send_buffer(ii)
             end do
          else
             AUX_send_buffer_tmp=AUX_send_buffer
          end if
          
       end if
       
       AUX_receive_buffer = 0.0_rp
       CALL MPI_SCATTER(AUX_send_buffer_tmp,recieve_num,MPI_REAL8, &
            AUX_receive_buffer,recieve_num,MPI_REAL8,0,MPI_COMM_WORLD,mpierr)
       spp(ss)%vars%flagCon(1:recieve_num) = INT(AUX_receive_buffer,is)

       if (params%mpi_params%rank.EQ.0_idef) then
          filename = TRIM(params%path_to_outputs) // "restart_file.h5"
          call h5fopen_f(filename, H5F_ACC_RDONLY_F, h5file_id, h5error)
          if (h5error .EQ. -1) then
             write(output_unit_write,'("KORC ERROR: Something went wrong in: &
                  &load_particles_ic --> h5fopen_f")')
             call KORC_ABORT(14)
          end if

          write(tmp_str,'(I18)') ss

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/flagCol"
          call load_array_from_hdf5(h5file_id,dset,AUX_send_buffer)

          call h5fclose_f(h5file_id, h5error)

          if (params%load_balance) then
             do ii=1,params%mpi_params%nmpi_prev*spp(ss)%ppp

                !write(6,*) 'send index',ii
                !write(6,*) 'tmp index',(mod(ii-1,params%mpi_params%nmpi))* &
                !     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1

                AUX_send_buffer_tmp((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     AUX_send_buffer(ii)
             end do
          else
             AUX_send_buffer_tmp=AUX_send_buffer
          end if
          
       end if
       
       AUX_receive_buffer = 0.0_rp
       CALL MPI_SCATTER(AUX_send_buffer_tmp,recieve_num,MPI_REAL8, &
            AUX_receive_buffer,recieve_num,MPI_REAL8,0,MPI_COMM_WORLD,mpierr)
       spp(ss)%vars%flagCol(1:recieve_num) = INT(AUX_receive_buffer,is)

       if (params%mpi_params%rank.EQ.0_idef) then
          filename = TRIM(params%path_to_outputs) // "restart_file.h5"
          call h5fopen_f(filename, H5F_ACC_RDONLY_F, h5file_id, h5error)
          if (h5error .EQ. -1) then
             write(output_unit_write,'("KORC ERROR: Something went wrong in: &
                  &load_particles_ic --> h5fopen_f")')
             call KORC_ABORT(14)
          end if

          write(tmp_str,'(I18)') ss

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/flagRE"
          call load_array_from_hdf5(h5file_id,dset,AUX_send_buffer)

          call h5fclose_f(h5file_id, h5error)

          !write(6,*) 'ppp',spp(ss)%ppp
          !write(6,*) 'nmpi_prev',params%mpi_params%nmpi_prev
          !write(6,*) 'nmpi',params%mpi_params%nmpi
          !write(6,*) 'nmpi_ratio',nmpi_ratio
          
          if (params%load_balance) then
             do ii=1,params%mpi_params%nmpi_prev*spp(ss)%ppp

                !write(6,*) 'send index',ii
                !write(6,*) 'tmp index',(mod(ii-1,params%mpi_params%nmpi))* &
                !     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1

                AUX_send_buffer_tmp((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     AUX_send_buffer(ii)
             end do
          else
             AUX_send_buffer_tmp=AUX_send_buffer
          end if

#if DBG_CHECK
          write(6,*) 'flagRE_send_buffer',AUX_send_buffer
          write(6,*) 'flagRE_send_buffer_tmp',AUX_send_buffer_tmp
#endif
          
       end if             
       
       AUX_receive_buffer = 0.0_rp
       CALL MPI_SCATTER(AUX_send_buffer_tmp,recieve_num,MPI_REAL8, &
            AUX_receive_buffer,recieve_num,MPI_REAL8,0,MPI_COMM_WORLD,mpierr)
       spp(ss)%vars%flagRE(1:recieve_num) = INT(AUX_receive_buffer,is)
       
       if (params%mpi_params%rank.EQ.0_idef) then
          filename = TRIM(params%path_to_outputs) // "restart_file.h5"
          call h5fopen_f(filename, H5F_ACC_RDONLY_F, h5file_id, h5error)
          if (h5error .EQ. -1) then
             write(output_unit_write,'("KORC ERROR: Something went wrong in: &
                  &load_particles_ic --> h5fopen_f")')
             call KORC_ABORT(14)
          end if

          write(tmp_str,'(I18)') ss

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/g"
          call load_array_from_hdf5(h5file_id,dset,AUX_send_buffer)

          call h5fclose_f(h5file_id, h5error)

          if (params%load_balance) then
             do ii=1,params%mpi_params%nmpi_prev*spp(ss)%ppp

                !write(6,*) 'send index',ii
                !write(6,*) 'tmp index',(mod(ii-1,params%mpi_params%nmpi))* &
                !     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1

                AUX_send_buffer_tmp((mod(ii-1,params%mpi_params%nmpi))* &
                     spp(ss)%ppp/nmpi_ratio+(ii-1)/params%mpi_params%nmpi+1) = &
                     AUX_send_buffer(ii)
             end do
          else
             AUX_send_buffer_tmp=AUX_send_buffer
          end if
          
       end if

       AUX_receive_buffer = 0.0_rp
       CALL MPI_SCATTER(AUX_send_buffer_tmp,recieve_num,MPI_REAL8, &
            AUX_receive_buffer,recieve_num,MPI_REAL8,0,MPI_COMM_WORLD,mpierr)
       spp(ss)%vars%g(1:recieve_num) = AUX_receive_buffer

       if (params%SC_E) then
          
          ALLOCATE(JSC0_buffer(F%dim_1D))
          ALLOCATE(JSC1_buffer(F%dim_1D))
          ALLOCATE(JSC2_buffer(F%dim_1D))
          ALLOCATE(JSC3_buffer(F%dim_1D))
          ALLOCATE(ESC_buffer(F%dim_1D))

          filename = TRIM(params%path_to_outputs) // "restart_file.h5"
          call h5fopen_f(filename, H5F_ACC_RDONLY_F, h5file_id, h5error)
          if (h5error .EQ. -1) then
             write(output_unit_write,'("KORC ERROR: Something went wrong in: &
                  &load_particles_ic --> h5fopen_f")')
             call KORC_ABORT(14)
          end if

          write(tmp_str,'(I18)') ss

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/J0_SC"
          call load_array_from_hdf5(h5file_id,dset,JSC0_buffer)
          
          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/J1_SC"
          call load_array_from_hdf5(h5file_id,dset,JSC1_buffer)

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/J2_SC"
          call load_array_from_hdf5(h5file_id,dset,JSC2_buffer)

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/J3_SC"
          call load_array_from_hdf5(h5file_id,dset,JSC3_buffer)

          dset = "/spp_" // TRIM(ADJUSTL(tmp_str)) // "/E_SC"
          call load_array_from_hdf5(h5file_id,dset,ESC_buffer)
          
          call h5fclose_f(h5file_id, h5error)

          F%J0_SC_1D%PHI=JSC0_buffer
          F%J1_SC_1D%PHI=JSC1_buffer
          F%J2_SC_1D%PHI=JSC2_buffer
          F%J3_SC_1D%PHI=JSC3_buffer
          F%E_SC_1D%PHI=ESC_buffer/params%cpp%Eo
          
          DEALLOCATE(JSC0_buffer)
          DEALLOCATE(JSC1_buffer)
          DEALLOCATE(JSC2_buffer)
          DEALLOCATE(JSC3_buffer)
          DEALLOCATE(ESC_buffer)
          
       end if
       
       DEALLOCATE(X_send_buffer)
       DEALLOCATE(X_receive_buffer)

       DEALLOCATE(V_send_buffer)
       DEALLOCATE(V_receive_buffer)

       DEALLOCATE(AUX_send_buffer)
       DEALLOCATE(AUX_receive_buffer)
    end do

    if (params%orbit_model(1:2).EQ.'GC') then
       params%GC_coords=.TRUE.
    end if

#if DBG_CHECK    
    call MPI_BARRIER(MPI_COMM_WORLD,mpierr)
    
    if (params%mpi_params%rank.eq.0_idef) then
       write(6,*) 'mpi',params%mpi_params%rank
       write(6,*) 'Y_R',spp(1)%vars%Y(:,1)*params%cpp%length
       write(6,*) 'Y_PHI',spp(1)%vars%Y(:,2)*params%cpp%length
       write(6,*) 'Y_Z',spp(1)%vars%Y(:,3)*params%cpp%length
       write(6,*) 'V_PLL',spp(1)%vars%V(:,1)*params%cpp%velocity*params%cpp%mass
       write(6,*) 'V_MU',spp(1)%vars%V(:,2)*params%cpp%energy/params%cpp%Bo
       write(6,*) 'flagCon',spp(1)%vars%flagCon
       write(6,*) 'flagCol',spp(1)%vars%flagCol
       write(6,*) 'flagRE',spp(1)%vars%flagRE
    end if

    call MPI_BARRIER(MPI_COMM_WORLD,mpierr)
    
    if (params%mpi_params%rank.eq.1_idef) then
       write(6,*) 'mpi',params%mpi_params%rank
       write(6,*) 'Y_R',spp(1)%vars%Y(:,1)*params%cpp%length
       write(6,*) 'Y_PHI',spp(1)%vars%Y(:,2)*params%cpp%length
       write(6,*) 'Y_Z',spp(1)%vars%Y(:,3)*params%cpp%length
       write(6,*) 'V_PLL',spp(1)%vars%V(:,1)*params%cpp%velocity*params%cpp%mass
       write(6,*) 'V_MU',spp(1)%vars%V(:,2)*params%cpp%energy/params%cpp%Bo
       write(6,*) 'flagCon',spp(1)%vars%flagCon
       write(6,*) 'flagCol',spp(1)%vars%flagCol
       write(6,*) 'flagRE',spp(1)%vars%flagRE
    end if

    call MPI_BARRIER(MPI_COMM_WORLD,mpierr)
    
    if (params%mpi_params%rank.eq.2_idef) then
       write(6,*) 'mpi',params%mpi_params%rank
       write(6,*) 'Y_R',spp(1)%vars%Y(:,1)*params%cpp%length
       write(6,*) 'Y_PHI',spp(1)%vars%Y(:,2)*params%cpp%length
       write(6,*) 'Y_Z',spp(1)%vars%Y(:,3)*params%cpp%length
       write(6,*) 'V_PLL',spp(1)%vars%V(:,1)*params%cpp%velocity*params%cpp%mass
       write(6,*) 'V_MU',spp(1)%vars%V(:,2)*params%cpp%energy/params%cpp%Bo
       write(6,*) 'flagCon',spp(1)%vars%flagCon
       write(6,*) 'flagCol',spp(1)%vars%flagCol
       write(6,*) 'flagRE',spp(1)%vars%flagRE
    end if

    call MPI_BARRIER(MPI_COMM_WORLD,mpierr)
    
    if (params%mpi_params%rank.eq.3_idef) then
       write(6,*) 'mpi',params%mpi_params%rank
       write(6,*) 'Y_R',spp(1)%vars%Y(:,1)*params%cpp%length
       write(6,*) 'Y_PHI',spp(1)%vars%Y(:,2)*params%cpp%length
       write(6,*) 'Y_Z',spp(1)%vars%Y(:,3)*params%cpp%length
       write(6,*) 'V_PLL',spp(1)%vars%V(:,1)*params%cpp%velocity*params%cpp%mass
       write(6,*) 'V_MU',spp(1)%vars%V(:,2)*params%cpp%energy/params%cpp%Bo
       write(6,*) 'flagCon',spp(1)%vars%flagCon
       write(6,*) 'flagCol',spp(1)%vars%flagCol
       write(6,*) 'flagRE',spp(1)%vars%flagRE
    end if

    call MPI_BARRIER(MPI_COMM_WORLD,mpierr)
#endif

    
  end subroutine load_particles_ic