initialize_m3d_c1 Subroutine

public subroutine initialize_m3d_c1(params, F, P, spp, init)

Arguments

Type IntentOptional AttributesName
type(KORC_PARAMS), intent(inout) :: params
type(FIELDS), intent(inout) :: F
type(PROFILES), intent(inout) :: P
type(SPECIES), intent(inout), DIMENSION(:):: spp
logical, intent(in) :: init

Calls

proc~~initialize_m3d_c1~~CallsGraph proc~initialize_m3d_c1 initialize_m3d_c1 interface~fio_close_field fio_close_field proc~initialize_m3d_c1->interface~fio_close_field interface~fio_set_int_option fio_set_int_option proc~initialize_m3d_c1->interface~fio_set_int_option proc~load_prev_iter load_prev_iter proc~initialize_m3d_c1->proc~load_prev_iter interface~fio_eval_field fio_eval_field proc~initialize_m3d_c1->interface~fio_eval_field interface~fio_get_field fio_get_field proc~initialize_m3d_c1->interface~fio_get_field interface~fio_get_options fio_get_options proc~initialize_m3d_c1->interface~fio_get_options interface~fio_open_source fio_open_source proc~initialize_m3d_c1->interface~fio_open_source interface~fio_get_real_field_parameter fio_get_real_field_parameter proc~initialize_m3d_c1->interface~fio_get_real_field_parameter interface~fio_allocate_search_hint fio_allocate_search_hint proc~initialize_m3d_c1->interface~fio_allocate_search_hint mpi_bcast mpi_bcast proc~load_prev_iter->mpi_bcast h5fclose_f h5fclose_f proc~load_prev_iter->h5fclose_f proc~korc_abort korc_abort proc~load_prev_iter->proc~korc_abort h5fopen_f h5fopen_f proc~load_prev_iter->h5fopen_f interface~load_from_hdf5 load_from_hdf5 proc~load_prev_iter->interface~load_from_hdf5 mpi_abort mpi_abort proc~korc_abort->mpi_abort proc~rload_from_hdf5 rload_from_hdf5 interface~load_from_hdf5->proc~rload_from_hdf5 proc~iload_from_hdf5 iload_from_hdf5 interface~load_from_hdf5->proc~iload_from_hdf5 proc~rload_from_hdf5->proc~korc_abort h5dclose_f h5dclose_f proc~rload_from_hdf5->h5dclose_f h5dread_f h5dread_f proc~rload_from_hdf5->h5dread_f proc~iload_from_hdf5->proc~korc_abort proc~iload_from_hdf5->h5dclose_f proc~iload_from_hdf5->h5dread_f

Called by

proc~~initialize_m3d_c1~~CalledByGraph proc~initialize_m3d_c1 initialize_m3d_c1 program~main main program~main->proc~initialize_m3d_c1

Contents

Source Code


Source Code

  SUBROUTINE initialize_m3d_c1(params, F, P, spp,init)

    TYPE(KORC_PARAMS), INTENT(INOUT)           :: params
    TYPE(FIELDS), INTENT(INOUT)                :: F
    TYPE(PROFILES), INTENT(INOUT)              :: P
    TYPE(SPECIES), DIMENSION(:), INTENT(INOUT) :: spp
    LOGICAL, INTENT(IN)  :: init
    
    INTEGER                                    :: ii
    INTEGER                                    :: pp
    INTEGER                                    :: status
    INTEGER                                    :: isrc
    real(c_double)  ::  time0,time1
    INTEGER (C_INT)                         :: FIO_tmp
    TYPE(C_PTR) :: hint_tmp
    real(rp), DIMENSION(3) :: x
    REAL(rp), DIMENSION(3)         :: Btmp
    
    if (init) then

       status = fio_open_source(FIO_M3DC1_SOURCE,           &
            TRIM(params%magnetic_field_filename)            &
            // C_NULL_CHAR, F%isrc)
       
       F%Efield = Efield
       F%PSIp_lim=PSIp_lim
       F%PSIp_0=PSIp_0
       F%ReInterp_2x1t=ReInterp_2x1t

       if (params%restart.OR.params%proceed) then

          isrc=F%isrc
          status = fio_get_options(isrc)
          status = fio_set_int_option(FIO_TIMESLICE, ind0_2x1t)
          status = fio_get_field(isrc, FIO_MAGNETIC_FIELD, F%FIO_B)
          
          hint_tmp=c_null_ptr
          x(1)=spp(1)%Ro
          x(2)=spp(1)%PHIo
          x(3)=spp(1)%Zo

          status = fio_eval_field(F%FIO_B, x(1),                      &
               Btmp(1),hint_tmp)

          F%Bo = Btmp(2)
          F%Eo = 1.0
          F%Ro = 1.0
          F%Zo = 1.0

          status = fio_close_field(F%FIO_B)
          
          call load_prev_iter(params)
          F%ind0_2x1t=params%prev_iter_2x1t+1
          
       else
          F%ind0_2x1t=ind0_2x1t
       end if

       F%ind_2x1t=F%ind0_2x1t



    else
       status = fio_close_field(F%FIO_B)
       status = fio_close_field(F%FIO_B+1)
       status = fio_close_field(F%FIO_A)

       status = fio_close_field(P%FIO_ne)
       status = fio_close_field(P%FIO_te)
       status = fio_close_field(P%FIO_ni)
       
    end if

    isrc=F%isrc
    
    status = fio_get_options(isrc)
       

    if (.not.F%ReInterp_2x1t) then
       status = fio_set_int_option(FIO_TIMESLICE, params%time_slice)
    else
       status = fio_set_int_option(FIO_TIMESLICE, F%ind_2x1t)
    end if

    
    status = fio_get_field(isrc, FIO_MAGNETIC_FIELD, F%FIO_B)
    status = fio_get_field(isrc, FIO_ELECTRIC_FIELD, F%FIO_E)
    status = fio_get_field(isrc, FIO_VECTOR_POTENTIAL, F%FIO_A)

    status = fio_get_real_field_parameter(F%FIO_B, FIO_TIME, time0)

    write(output_unit_write,*) 'FIO present time index',F%ind_2x1t
    write(output_unit_write,*) 'FIO present time',time0

    if (F%ReInterp_2x1t) then
       status = fio_set_int_option(FIO_TIMESLICE, F%ind_2x1t+1)

       status = fio_get_field(isrc, FIO_MAGNETIC_FIELD, FIO_tmp)
       
       status = fio_get_real_field_parameter(FIO_tmp, FIO_TIME, time1)
       write(output_unit_write,*) 'FIO next time index',F%ind_2x1t+1
       write(output_unit_write,*) 'FIO next time',time1

       status = fio_set_int_option(FIO_TIMESLICE, F%ind_2x1t)

       params%snapshot_frequency=time1-time0

       !write(6,*) 'snapshot_frequency',params%snapshot_frequency
       !write(6,*) 'dt',params%dt
       
       if (.not.init) then
          
          params%t_skip = FLOOR(params%snapshot_frequency/params%cpp%time/ &
               params%dt,ip)
          
       end if
              
    end if    
    
    if (.not.F%Efield) F%FIO_E=-1

    status = fio_set_int_option(FIO_SPECIES, FIO_ELECTRON);
    status = fio_get_field(isrc, FIO_DENSITY, P%FIO_ne);
    status = fio_get_field(isrc, FIO_TEMPERATURE, P%FIO_te);

    status = fio_set_int_option(FIO_SPECIES, FIO_MAIN_ION);
    status = fio_get_field(isrc, FIO_DENSITY, P%FIO_ni);


    if (init) then
       if (.NOT.(params%restart.OR.params%proceed)) then
          hint_tmp=c_null_ptr
          x(1)=spp(1)%Ro
          x(2)=spp(1)%PHIo
          x(3)=spp(1)%Zo

          status = fio_eval_field(F%FIO_B, x(1),                      &
               Btmp(1),hint_tmp)

          F%Bo = Btmp(2)
          F%Eo = 1.0
          F%Ro = 1.0
          F%Zo = 1.0
       end if

       do ii = 1, params%num_species

          do pp = 1, spp(ii)%ppp
             status = fio_allocate_search_hint(isrc, spp(ii)%vars%hint(pp))
             !spp(ii)%vars%hint(pp)=c_null_ptr
          end do

          spp(ii)%vars%cart = .false.
       end do
    end if

    if (params%mpi_params%rank .EQ. 0) then
       write(output_unit_write,*) 'Calculate B',F%FIO_B
       write(output_unit_write,'("Magnetic field: ",E17.10)') F%Bo
       write(output_unit_write,*) 'Calculate E',F%FIO_E
       write(output_unit_write,'("Electric field: ",E17.10)') F%Eo
       write(output_unit_write,*) 'Calculate A',F%FIO_A
       write(output_unit_write,*) 'Calculate ne',P%FIO_ne
       write(output_unit_write,*) 'Calculate Te',P%FIO_te
       write(output_unit_write,*) 'Calculate ni',P%FIO_ni
    end if
       
  END SUBROUTINE initialize_m3d_c1