read_inputs Subroutine

public subroutine read_inputs()

Reads input namelist file and stores the results into inputs, nbi, and beam_grid 20 for suffixes and seperators e.g. /, _npa.h5, ...

Arguments

None

Calls

proc~~read_inputs~~CallsGraph proc~read_inputs read_inputs proc~rng_seed rng_seed proc~read_inputs->proc~rng_seed proc~my_rank my_rank proc~read_inputs->proc~my_rank proc~identify_transition identify_transition proc~read_inputs->proc~identify_transition

Called by

proc~~read_inputs~~CalledByGraph proc~read_inputs read_inputs program~fidasim fidasim program~fidasim->proc~read_inputs

Contents

Source Code


Source Code

subroutine read_inputs
    !+ Reads input namelist file and stores the results into [[libfida:inputs]],
    !+ [[libfida:nbi]], and [[libfida:beam_grid]]
    character(charlim) :: runid,result_dir, tables_file
    character(charlim) :: distribution_file, equilibrium_file
    character(charlim) :: geometry_file, neutrals_file
    integer            :: pathlen, calc_neutron, seed, calc_cfpd
    integer            :: calc_brems, calc_dcx, calc_halo, calc_cold, calc_bes
    integer            :: calc_fida, calc_pfida, calc_npa, calc_pnpa
    integer            :: calc_birth,calc_fida_wght,calc_npa_wght, calc_res
    integer            :: load_neutrals,verbose,flr,split,stark_components
    integer            :: output_neutral_reservoir
    integer(Int64)     :: n_fida,n_pfida,n_npa,n_pnpa,n_nbi,n_halo,n_dcx,n_birth
    integer(Int32)     :: shot,nlambda,ne_wght,np_wght,nphi_wght,nlambda_wght
    integer(Int32)     :: adaptive, max_cell_splits
    real(Float64)      :: time,lambdamin,lambdamax,emax_wght
    real(Float64)      :: lambdamin_wght,lambdamax_wght
    real(Float64)      :: ab,pinj,einj,current_fractions(3)
    integer(Int32)     :: nx,ny,nz
    real(Float64)      :: xmin,xmax,ymin,ymax,zmin,zmax
    real(Float64)      :: alpha,beta,gamma,origin(3)
    real(Float64)      :: split_tol
    logical            :: exis, error

    NAMELIST /fidasim_inputs/ result_dir, tables_file, distribution_file, &
        geometry_file, equilibrium_file, neutrals_file, shot, time, runid, &
        calc_brems, calc_dcx,calc_halo, calc_cold, calc_fida, calc_bes,&
        calc_pfida, calc_npa, calc_pnpa,calc_birth, calc_res, seed, flr, split, &
        calc_fida_wght, calc_npa_wght, load_neutrals, verbose, stark_components, &
        calc_neutron, calc_cfpd, n_fida, n_pfida, n_npa, n_pnpa, n_nbi, n_halo, n_dcx, n_birth, &
        ab, pinj, einj, current_fractions, output_neutral_reservoir, &
        nx, ny, nz, xmin, xmax, ymin, ymax, zmin, zmax, &
        origin, alpha, beta, gamma, &
        ne_wght, np_wght, nphi_wght, &
        nlambda, lambdamin,lambdamax,emax_wght, &
        nlambda_wght,lambdamin_wght,lambdamax_wght, &
        adaptive, max_cell_splits, split_tol

    inquire(file=namelist_file,exist=exis)
    if(.not.exis) then
        write(*,'(a,a)') 'READ_INPUTS: Input file does not exist: ', trim(namelist_file)
        stop
    endif

    ! variables that are not changed not be auto-initalized
    ! provide reasonable defaults here
    result_dir="."
    tables_file="."
    distribution_file="."
    geometry_file ="."
    equilibrium_file="."
    neutrals_file="."
    shot=0
    time=0
    runid="0"
    seed = -1
    calc_brems=0
    calc_bes=0
    calc_dcx=0
    calc_halo=0
    calc_cold=0
    calc_fida=0
    calc_pfida=0
    calc_npa=0
    calc_pnpa=0
    calc_birth=0
    calc_res=0
    flr=2
    stark_components=0
    split=1
    calc_fida_wght=0
    calc_npa_wght=0
    load_neutrals=0
    output_neutral_reservoir=1
    verbose=0
    calc_neutron=0
    calc_cfpd=0
    n_fida=0
    n_pfida=0
    n_npa=0
    n_pnpa=0
    n_nbi=0
    n_halo=0
    n_dcx=0
    n_birth=0
    ab=0
    pinj=0
    einj=0
    current_fractions=0
    nx=0
    ny=0
    nz=0
    xmin=0
    xmax=0
    ymin=0
    ymax=0
    zmin=0
    zmax=0
    origin=0
    alpha=0
    beta=0
    gamma=0
    ne_wght=0
    np_wght=0
    nphi_wght=0
    nlambda=0
    lambdamin=0
    lambdamax=0
    emax_wght=0
    nlambda_wght=0
    lambdamin_wght=0
    lambdamax_wght=0
    adaptive=0
    max_cell_splits=1
    split_tol=0

    open(13,file=namelist_file)
    read(13,NML=fidasim_inputs)
    close(13)

    !!General Information
    inputs%shot_number=shot
    inputs%time=time
    inputs%runid=runid
    inputs%result_dir=result_dir

    !!Input Files
    inputs%tables_file=tables_file
    inputs%geometry_file=geometry_file
    inputs%equilibrium_file=equilibrium_file
    inputs%distribution_file=distribution_file
    inputs%neutrals_file=neutrals_file

    !! RNG seed
    inputs%seed = seed
    if(inputs%seed.lt.0) inputs%seed = rng_seed()

    !!Simulation Switches
    if((calc_brems+calc_bes+calc_dcx+calc_halo+&
        calc_cold+calc_fida+calc_pfida).gt.0) then
        inputs%calc_spec=1
        inputs%tot_spectra=calc_brems+calc_bes+calc_dcx+calc_halo+&
                           calc_cold+calc_fida+calc_pfida
    else
        inputs%calc_spec=0
        inputs%tot_spectra=0
    endif

    inputs%calc_beam = 0
    if((calc_bes+calc_birth+calc_dcx+&
        calc_halo+calc_fida+calc_npa+&
        calc_fida_wght+calc_npa_wght).gt.0) then
        inputs%calc_nbi_dens=1
        inputs%calc_beam=1
    else
        inputs%calc_nbi_dens=0
    endif

    if((calc_dcx+calc_halo+calc_fida+calc_npa+&
        calc_fida_wght+calc_npa_wght).gt.0) then
        inputs%calc_dcx_dens=1
        inputs%calc_beam=1
    else
        inputs%calc_dcx_dens=0
    endif

    if((calc_halo+calc_fida+calc_npa+&
        calc_fida_wght+calc_npa_wght).gt.0) then
        inputs%calc_halo_dens=1
        inputs%calc_beam=1
    else
        inputs%calc_halo_dens=0
    endif

    inputs%calc_brems=calc_brems
    inputs%calc_bes=calc_bes
    inputs%calc_dcx=calc_dcx
    inputs%calc_halo=calc_halo
    inputs%calc_cold=calc_cold
    inputs%calc_fida=calc_fida
    inputs%calc_pfida=calc_pfida
    inputs%calc_npa=calc_npa
    inputs%calc_pnpa=calc_pnpa
    inputs%calc_birth=calc_birth
    inputs%calc_fida_wght=calc_fida_wght
    inputs%calc_npa_wght=calc_npa_wght
    inputs%calc_neutron=calc_neutron
    inputs%calc_cfpd=calc_cfpd
    inputs%calc_res = calc_res

    !! Misc. Settings
    inputs%load_neutrals=load_neutrals
    inputs%output_neutral_reservoir=output_neutral_reservoir

    inputs%verbose=verbose
    inputs%flr = flr
    inputs%stark_components = stark_components
    inputs%split = split

    !!Monte Carlo Settings
    inputs%n_fida=max(10,n_fida)
    inputs%n_pfida=max(10,n_pfida)
    inputs%n_npa=max(10,n_npa)
    inputs%n_pnpa=max(10,n_pnpa)
    inputs%n_nbi=max(10,n_nbi)
    inputs%n_halo=max(10,n_halo)
    inputs%n_dcx=max(10,n_dcx)
    inputs%n_birth= max(1,nint(n_birth/real(n_nbi)))

    !!Neutral Beam Settings
    beam_mass=ab
    nbi%current_fractions=current_fractions
    nbi%einj=einj
    nbi%pinj=pinj

    !!Weight Function Settings
    inputs%ne_wght=ne_wght
    inputs%np_wght=np_wght
    inputs%nphi_wght=nphi_wght
    inputs%emax_wght=emax_wght
    inputs%nlambda_wght = nlambda_wght
    inputs%lambdamin_wght=lambdamin_wght
    inputs%lambdamax_wght=lambdamax_wght

    !!Wavelength Grid Settings
    inputs%nlambda=nlambda
    inputs%lambdamin=lambdamin
    inputs%lambdamax=lambdamax
    inputs%dlambda=(inputs%lambdamax-inputs%lambdamin)/inputs%nlambda

    !!Adaptive Time Step Settings
    inputs%adaptive=adaptive
    if(inputs%adaptive.eq.0) then
        inputs%max_cell_splits=1
        inputs%split_tol=0.0
    else
        inputs%max_cell_splits=max_cell_splits
        inputs%split_tol=split_tol
    endif

    !!Beam Grid Settings
    beam_grid%nx=nx
    beam_grid%ny=ny
    beam_grid%nz=nz
    beam_grid%xmin=xmin
    beam_grid%xmax=xmax
    beam_grid%ymin=ymin
    beam_grid%ymax=ymax
    beam_grid%zmin=zmin
    beam_grid%zmax=zmax
    beam_grid%alpha=alpha
    beam_grid%beta=beta
    beam_grid%gamma=gamma
    beam_grid%origin=origin

#ifdef _MPI
    if(my_rank().ne.0) inputs%verbose=0
#endif

    if(inputs%verbose.ge.1) then
        write(*,'(a)') "---- Shot settings ----"
        write(*,'(T2,"Shot: ",i8)') inputs%shot_number
        write(*,'(T2,"Time: ",i4," [ms]")') int(inputs%time*1.d3)
        write(*,'(T2,"Runid: ",a)') trim(adjustl(inputs%runid))
        write(*,*) ''
        write(*,'(a)') "---- Input files ----"
    endif

    error = .False.

    inquire(file=inputs%tables_file,exist=exis)
    if(exis) then
        if(inputs%verbose.ge.1) then
            write(*,'(T2,"Tables file: ",a)') trim(inputs%tables_file)
        endif
    else
        if(inputs%verbose.ge.0) then
            write(*,'(a,a)') 'READ_INPUTS: Tables file does not exist: ', &
                             trim(inputs%tables_file)
        endif
        error = .True.
    endif

    inquire(file=inputs%geometry_file,exist=exis)
    if(exis) then
        if(inputs%verbose.ge.1) then
            write(*,'(T2,"Geometry file: ",a)') trim(inputs%geometry_file)
        endif
    else
        if(inputs%verbose.ge.0) then
            write(*,'(a,a)') 'READ_INPUTS: Geometry file does not exist: ', &
                             trim(inputs%geometry_file)
        endif
        error = .True.
    endif

    inquire(file=inputs%equilibrium_file,exist=exis)
    if(exis) then
        if(inputs%verbose.ge.1) then
            write(*,'(T2,"Equilibrium file: ",a)') trim(inputs%equilibrium_file)
        endif
    else
        if(inputs%verbose.ge.0) then
            write(*,'(a,a)') 'READ_INPUTS: Equilibrium file does not exist: ', &
                              trim(inputs%equilibrium_file)
        endif
        error = .True.
    endif

    inquire(file=inputs%distribution_file,exist=exis)
    if(exis) then
        if(inputs%verbose.ge.1) then
            write(*,'(T2,"Distribution file: ",a)') trim(inputs%distribution_file)
        endif
    else
        if(inputs%verbose.ge.0) then
            write(*,'(a,a)') 'READ_INPUTS: Distribution file does not exist: ', &
                             trim(inputs%distribution_file)
        endif
        error = .True.
    endif

    pathlen = len_trim(inputs%result_dir)+len_trim(inputs%runid) + 20
    !+20 for suffixes and seperators e.g. /, _npa.h5, ...
    if(pathlen.gt.charlim) then
        if(inputs%verbose.ge.0) then
            write(*,'(a,i3,a,i3)') 'READ_INPUTS: Result directory path + runID use too many characters: ', &
                                   pathlen-20,'>', charlim-20
        endif
        error = .True.
    endif

    if(inputs%verbose.ge.1) then
        write(*,*) ''
    endif

    if(inputs%adaptive.ge.0.and.inputs%adaptive.le.7) then
        write(*,'(a)') "---- Adaptive time step settings ----"
        write(*,'(T2,"Adaptive: ",i2)') inputs%adaptive
        if(inputs%adaptive.gt.0) then
            if(inputs%max_cell_splits.gt.1) then
                write(*,'(T2,"Max cell splits: ",i4)') inputs%max_cell_splits
            else
                write(*,'(a,i4)') 'READ_INPUTS: max_cell_splits must be greater than 1 if adaptive time stepping is on: ', &
                                  inputs%max_cell_splits
                error=.True.
            endif
            if(inputs%split_tol.gt.0.0.and.inputs%split_tol.lt.1.0) then
                write(*, '(T2,"Split tolerance: ",1f9.6)') inputs%split_tol
            else
                write(*,'(a,f2.6)') 'READ_INPUTS: split_tol must be positive if adaptive time stepping is on: ', &
                                inputs%split_tol
                error=.True.
            endif
        else
            write(*,'(a)') ' Adaptive time stepping is off'
        endif
        write(*,*) ''
    else
        write(*,'(a,i2)') 'READ_INPUTS: Invalid adaptive switch setting, must be within [0,7]: ', &
                          inputs%adaptive
        error=.True.
    endif

    if(error) then
        stop
    endif

    ! Identify the transition from lambdamax and lambdamin
    call identify_transition(n_stark, stark_pi, stark_sigma, &
                                stark_intens, stark_wavel, line_lambda0)

end subroutine read_inputs