read_job Subroutine

public subroutine read_job(job_file)

Note

Function that reads a script file. Keywords are identified and corresponding actions are triggered.

Arguments

Type IntentOptional Attributes Name
character(len=512), intent(in) :: job_file

job file with macros to execute


Calls

proc~~read_job~~CallsGraph proc~read_job read_job get_unit get_unit proc~read_job->get_unit proc~acf_theo acf_theo proc~read_job->proc~acf_theo proc~add_nois add_nois proc~read_job->proc~add_nois proc~alloc_tabs alloc_tabs proc~read_job->proc~alloc_tabs proc~apod_acf apod_acf proc~read_job->proc~apod_acf proc~apod_sur apod_sur proc~read_job->proc~apod_sur proc~calc_acf calc_acf proc~read_job->proc~calc_acf proc~calc_ffh calc_ffh proc~read_job->proc~calc_ffh proc~calc_ord calc_ord proc~read_job->proc~calc_ord proc~calc_z_f calc_z_f proc~read_job->proc~calc_z_f proc~calc_z_i calc_z_i proc~read_job->proc~calc_z_i proc~def_size def_size proc~read_job->proc~def_size proc~digi_fil digi_fil proc~read_job->proc~digi_fil proc~end_loop end_loop proc~read_job->proc~end_loop proc~end_scri end_scri proc~read_job->proc~end_scri proc~make_msk make_msk proc~read_job->proc~make_msk proc~make_scratches make_scratches proc~read_job->proc~make_scratches proc~make_tex make_tex proc~read_job->proc~make_tex proc~nb_procs nb_procs proc~read_job->proc~nb_procs proc~plt__acf plt__acf proc~read_job->proc~plt__acf proc~read_img read_img proc~read_job->proc~read_img proc~repr_img repr_img proc~read_job->proc~repr_img proc~save_img save_img proc~read_job->proc~save_img proc~smooth__ smooth__ proc~read_job->proc~smooth__ proc~spct_sur spct_sur proc~read_job->proc~spct_sur proc~sta_loop sta_loop proc~read_job->proc~sta_loop proc~sta_scri sta_scri proc~read_job->proc~sta_scri proc~sta_theo sta_theo proc~read_job->proc~sta_theo proc~stat_sur stat_sur proc~read_job->proc~stat_sur proc~sub_surf sub_surf proc~read_job->proc~sub_surf proc~surface_analysis surface_analysis proc~read_job->proc~surface_analysis selectcase selectcase proc~read_job->selectcase str_remove_chars str_remove_chars proc~read_job->str_remove_chars ellipse_acf ellipse_acf proc~acf_theo->ellipse_acf proc~calc_imp_acf calc_imp_acf proc~acf_theo->proc~calc_imp_acf proc~build_heights build_heights proc~add_nois->proc~build_heights scramble scramble proc~add_nois->scramble std_array std_array proc~add_nois->std_array proc~apod_acf->ellipse_acf proc~apod2 apod2 proc~apod_acf->proc~apod2 apod apod proc~apod_sur->apod proc~apod_sur->std_array proc~acf_wiener acf_wiener proc~calc_acf->proc~acf_wiener calc_fftw3_real_bwd calc_fftw3_real_bwd proc~calc_ffh->calc_fftw3_real_bwd calc_fftw3_real_fwd calc_fftw3_real_fwd proc~calc_ffh->calc_fftw3_real_fwd calc_moments calc_moments proc~calc_ffh->calc_moments init_order init_order proc~calc_ord->init_order sort_array2 sort_array2 proc~calc_ord->sort_array2 proc~calc_z_f->proc~build_heights proc~calc_z_f->std_array proc~calc_z_i->proc~build_heights proc~calc_z_i->scramble proc~calc_z_i->std_array proc~digi_fil->calc_fftw3_real_bwd proc~digi_fil->calc_fftw3_real_fwd proc~digi_fil->std_array end_fftw3 end_fftw3 proc~end_scri->end_fftw3 calc_fftw3 calc_fftw3 proc~make_tex->calc_fftw3 proc~make_tex->calc_moments omp_get_num_procs omp_get_num_procs proc~nb_procs->omp_get_num_procs proc~plt__acf->get_unit dir_separator dir_separator proc~plt__acf->dir_separator proc~plt__acf->ellipse_acf mkdir mkdir proc~plt__acf->mkdir read_surf read_surf proc~read_img->read_surf proc~repr_img->proc~alloc_tabs proc~repr_img->proc~surface_analysis proc~repr_img->apod proc~repr_img->calc_moments fft_filter fft_filter proc~repr_img->fft_filter init_scal init_scal proc~repr_img->init_scal proc~repr_img->proc~acf_wiener proc~repr_img->scramble proc~repr_img->sort_array2 proc~repr_img->std_array write_surf write_surf proc~repr_img->write_surf proc~save_img->init_scal proc~save_img->write_surf proc~smooth__->fft_filter proc~smooth__->std_array proc~spct_sur->get_unit proc~spct_sur->proc~apod_sur proc~spct_sur->calc_fftw3 trans_corner2center trans_corner2center proc~spct_sur->trans_corner2center random_init random_init proc~sta_scri->random_init proc~stat_sur->proc~apod_sur proc~stat_sur->proc~surface_analysis proc~stat_sur->proc~acf_wiener proc~stat_sur->sort_array2 proc~stat_sur->std_array proc~sub_surf->calc_moments proc~sub_surf->end_fftw3 fftw_plan_with_nthreads fftw_plan_with_nthreads proc~sub_surf->fftw_plan_with_nthreads proc~sub_surf->init_order omp_get_max_threads omp_get_max_threads proc~sub_surf->omp_get_max_threads proc~sub_surf->proc~acf_wiener proc~sub_surf->proc~build_heights proc~sub_surf->proc~calc_imp_acf proc~calc_res_acf calc_res_acf proc~sub_surf->proc~calc_res_acf progress_bar_terminal progress_bar_terminal proc~sub_surf->progress_bar_terminal proc~sub_surf->sort_array2 tab_end_fftw3_real tab_end_fftw3_real proc~sub_surf->tab_end_fftw3_real tab_init_fftw3_real tab_init_fftw3_real proc~sub_surf->tab_init_fftw3_real proc~surface_analysis->get_unit abbott_param abbott_param proc~surface_analysis->abbott_param calc_median calc_median proc~surface_analysis->calc_median proc~surface_analysis->calc_moments calcul_asfc_hermite calcul_asfc_hermite proc~surface_analysis->calcul_asfc_hermite calcul_normales calcul_normales proc~surface_analysis->calcul_normales proc~surface_analysis->ellipse_acf proc~surface_analysis->fft_filter indice_fractal indice_fractal proc~surface_analysis->indice_fractal multiple_anisotropy multiple_anisotropy proc~surface_analysis->multiple_anisotropy peaks_and_pits_curvatures peaks_and_pits_curvatures proc~surface_analysis->peaks_and_pits_curvatures proc~surface_analysis->proc~acf_wiener surf_area surf_area proc~surface_analysis->surf_area topology topology proc~surface_analysis->topology proc~acf_wiener->calc_fftw3_real_bwd proc~acf_wiener->calc_fftw3_real_fwd proc~acf_wiener->trans_corner2center tab_calc_fftw3_real_bwd tab_calc_fftw3_real_bwd proc~acf_wiener->tab_calc_fftw3_real_bwd tab_calc_fftw3_real_fwd tab_calc_fftw3_real_fwd proc~acf_wiener->tab_calc_fftw3_real_fwd proc~build_heights->sort_array2 proc~build_heights->std_array proc~pikaia_skku_solver pikaia_skku_solver proc~build_heights->proc~pikaia_skku_solver proc~profil_theo_trie_1d profil_theo_trie_1D proc~build_heights->proc~profil_theo_trie_1d proc~calc_imp_acf->proc~apod2 proc~autocov_impo autocov_impo proc~calc_imp_acf->proc~autocov_impo init pikaia_class%init proc~pikaia_skku_solver->init solve pikaia_class%solve proc~pikaia_skku_solver->solve proc~profil_theo_trie_1d->std_array

Called by

proc~~read_job~~CalledByGraph proc~read_job read_job proc~prg_surf prg_surf proc~prg_surf->proc~read_job program~main main program~main->proc~prg_surf

Source Code

   subroutine read_job(job_file)
   !================================================================================================
   !<@note Function that reads a script file. Keywords are identified and corresponding actions are
   !< triggered.
   !<@endnote
   !------------------------------------------------------------------------------------------------
   implicit none
   character(len=512), intent(in) :: job_file   !! *job file with macros to execute*

      integer  (kind=I4) :: vide
      character(len=512) :: keyword
      logical  (kind=I4) :: first

      ! script file
      call get_unit( JOB )
      open(unit = JOB, file = trim(job_file), status = 'old')

      call check_empty_lines()

      ! witness file
      call get_unit( SPY )
      open(unit = SPY, file = "out/spy.txt", status = 'unknown')

      LINE_READ = 0
      first     = .true. ! determine if calc_z_f has been already used

      ! default function used: the tangent function
      PARAM%func_gen = FCT_TANG

      ! default is 2 bounds when generating profiles
      PARAM%nparam = 2

      PARAM%calc_mstt = .true.

      PARAM%calc_zf = .true.

      ! default cutting plane for correlation lengths
      PARAM%orig_surf%cut = 0.5
      PARAM%curr_surf%cut = 0.5

o:    do

         keyword = repeat( ' ', len(keyword) )

         read(JOB, *, iostat = vide) keyword ; LINE_READ = LINE_READ + 1

         ! remove unwanted characters from keyword
         keyword = str_remove_chars(string = trim(keyword), chars = '- # *')

         write(SPY, '(a6,I4.4,a)') "line: ", LINE_READ, ' ', trim(keyword)

         selectcase( keyword(1:8) )

            case('ACF_THEO')
               ! desired acf
               call acf_theo()

            case('ADD_NOIS')
               ! add noise to current surface
               call add_nois()

            case('ALLOCATE')
               ! allocate tabs in PARAM type
               call alloc_tabs()

            case('ANALYSIS')
               call surface_analysis()

            case('APOD_ACF')
               ! apodize acf
               call apod_acf()

            case('APOD_SUR')
               ! apodize surface
               call apod_sur()

            case('CALC_ACF')
               ! determine the surface acf
               call calc_acf()

            case('CALC_FFH')
               ! deigital filter
               call calc_ffh()

            case('CALC_ORD')
               ! determine height order
               call calc_ord()

            case('CALC_Z_F')
               ! final heights
               call calc_z_f()

            case('CALC_Z_I')
               ! starting heights
               call calc_z_i()

            case('DEF_SIZE')
               ! image size
               call def_size()

            case('DIGI_FIL')
               ! apply digital filter
               call digi_fil()

            case('END_LOOP')
               ! loop end
               call end_loop()

            case('END_SCRI')
               ! close the script reading
               call end_scri()

               exit o

            case('MAKE_MSK')
               ! turn image into mask
               call make_msk()

            case('MAKE_SCR')
               ! make secratches
               call make_scratches()

            case('MAKE_TEX')
               ! build a texture
               call make_tex()

            case('NB_PROCS')
               ! number of threads
               call nb_procs()

            case('PLT__ACF')
               ! print the correlation graphs and/or determine
               ! if the stop criterion is reached.
               call plt__acf()

            case('READ_PRF')
               ! read image
               call read_img()

            case('REPR_IMG')
               ! reproduce image
               call repr_img()

            case('SAVE_ACF')
               ! save the acf surface
               if ( sum( PARAM%imp_acf ) == 0 ) then
                  call save_img( tab = PARAM%acf_surf )     ! IN: real acf
               else
                  call save_img( tab = PARAM%imp_acf )      ! IN: wanted acf
               endif

            case('SAVE_PRF')
               ! save image
               call save_img( tab = PARAM%surf )            ! IN: surface to save

            case('SMOOTH__')
               ! low-pass filter
               call smooth__()

            case('SPCT_SUR')
               ! surface spectrum frequencies
               call spct_sur()

            case('STA_LOOP')
               ! loop start
               call sta_loop()

            case('STA_SCRI')
               ! start script
               call sta_scri()

            case('STA_THEO')
               ! desired stat moments
               call sta_theo()

            case('STAT_SUR')
               ! surface moments as reference
               call stat_sur()

            case('SUB_SURF')
               ! extract the best surface
               call sub_surf()

         endselect

      enddo o

      close(JOB)

   contains

      subroutine check_empty_lines()
         !! Check for empty lines in the script
      implicit none

         do

            read(JOB, '(A)', iostat = vide) keyword

            if (vide < 0) then

               rewind(JOB)
               return

            endif

            if (len_trim(keyword) == 0) stop 'Empty line in script file'

         enddo

      return
      endsubroutine check_empty_lines

   endsubroutine read_job