Commit 81a0e791 authored by Zolnierczuk, Piotr's avatar Zolnierczuk, Piotr
Browse files

allow to change echo sign

use command eshape, for example:

	eshape triangular negative
	eshape flux_weighted positive
parent 0562bab4
Loading
Loading
Loading
Loading
+47 −14
Original line number Diff line number Diff line
@@ -237,9 +237,15 @@ program drspine
  istat = add_tab_expansion('eshape gaussian')
  istat = add_tab_expansion('eshape flux_weighted')
  istat = add_tab_expansion('eshape triangular')
  istat = add_tab_expansion('eshape negative')
  istat = add_tab_expansion('eshape positive')
  !
  istat = add_tab_expansion('read')
  istat = add_tab_expansion('read as')
  !istat = add_tab_expansion('read tfac')
  !istat = add_tab_expansion('read phase_offset')
  !istat = add_tab_expansion('read sample_temp')
  !istat = add_tab_expansion('read scattering_angle')
  !
  istat = add_tab_expansion('role')
  istat = add_tab_expansion('role as')
@@ -445,12 +451,17 @@ program drspine
     !> COMMAND: eshape
     !-------------------------------------------------------------
     if(command('eshape ', &
          ' eshape [basic|gaussian|flux_weighted|triangular] '//LF//&
          ' eshape [basic|gaussian|flux_weighted|triangular] [negative|positive]'//LF//&
          '   - set echo shape function '//LF//&
          ''//LF//&
          '   options:'//LF//&
          '     basic|b          - "basic" echo shape (single wavelength) '//LF//&
          '     gaussian|g       - gaussian wavelength spectrum'//LF//&
          '     flux_weighted|fw - flux weighted wavelength spectrum (SNS-NSE)'//LF//&
          '     triangular|t     - triangular gaussian smeared spectrum (J-NSE)$')) then
          '     triangular|t     - triangular gaussian smeared spectrum (J-NSE)'//LF//&
          ''//LF//&
          '     positive         - positive echo sign (maximum at the symmetry phase, default)'//LF//&
          '     negative         - negative echo sign (minimum at the symmetry phase)$')) then
        !              ===============
        call cmd_eshape()
        cycle commandloop
@@ -461,16 +472,21 @@ program drspine
     !-------------------------------------------------------------
     !> COMMAND: read
     if(command('read  ', &
          " read <filename(s)> [as (ref|sample|background|auto)] [tfac <val>] [phase_offset <val>]"//LF//&
          "   - read filename(s) and assign a role to them (default auto)"//LF//&
          "     transmission_ratio <val>  - set the relative transmission factor to <val>"//LF//&
          "     tfac <val>                - set the relative transmission factor (short) to <val>"//LF//&
          "                           <val> = transmission(sample)/transmission(bgr)"//LF//&
          "     phase_offset <val>  - set the initial phase offset"//LF//&
          " read <filename(s)> [as (ref|sample|background|auto)] [options]"//LF//&
          "   - read filename(s) and assign a 'role' to each file"//LF//&
          "                      default 'role' is auto, i.e inferred from the data file"//LF//&
          "     options:"//LF//&
          "       transmission_ratio <tfac>  - set the relative transmission factor to <val>"//LF//&
          "                          <tfac>  = transmission(sample)/transmission(bgr)"//LF//&
          "       tfac <val>                 - same as transmission_ratio"//LF//&
          "       phase_offset <val>         - set the initial phase offset to <val>"//LF//&
          "       sample_temp <val>          - override the sample temperature and set it to <val>"//LF//&
          "       scattering_angle <val>     - override the scattering angle and set it to <val>"//LF//&
          ""//LF//&
          " read "//LF//&
          "   - read filenames given at the OS command line"//LF//&
          ""//LF//&
          " NOTE: transmission factor is only valid for sample data"//LF//&
          " NOTE: transmission factor option is only valid for sample data"//LF//&
          "       and is ignored for other data roles$")) then
        call cmd_read()
        cycle commandloop
@@ -1286,22 +1302,37 @@ CONTAINS
  !!    eshape flux_weighted (fw)
  !-------------------------------------------------------------
  subroutine cmd_eshape
    integer :: nselect, new_shape
    integer :: nselect, new_shape, echo_sign

    new_shape = -1
    nselect   = 0
    echo_sign = 0

    nselect = parse_command_flags('basic,b,gaussian,g,flux_weighted,fw,triangular,t')

    if ( found('negative') ) echo_sign = -1
    if ( found('positive') ) echo_sign = +1 

    if (nselect < 0 ) then
       call msg_error('drspine', 'command eshape: echo shape keywords are mutually exclusive', ERROR_OPTION_SYNTAX)
       return
    else if (nselect==0) then
       write(output_unit,'(a, i0, a)') ' current echo shape: '//trim(cformat_eshape())//' [', eshape_get(), ']'
    else if (nselect==0 .and. echo_sign==0) then
       write(output_unit,'(a, i0, a)', advance='no') ' current echo shape: '//trim(cformat_eshape())//' [', eshape_get(), ']'
       select case(eshape_sign)
       case (-1)
          write(output_unit,'(a,i0,a)') ', sign=negative (', eshape_sign, ')'
       case ( 1)
          write(output_unit,'(a,i0,a)') ', sign=positive (', eshape_sign, ')'
       case default
          call msg_error('drspine', 'command eshape: invalid echo shape sign')
       end select
       return
    end if


    select case(nselect)
    case(0)
        continue ! do not change
    case( 1,  2)
        new_shape = E_BASIC_ESHAPE
    case( 4,  8)
@@ -1312,8 +1343,10 @@ CONTAINS
        new_shape = E_TG_ESHAPE
    end select

    if (echo_sign /= 0) eshape_sign = int(sign(1,echo_sign))

    ! set the echo shape
    call eshape_select(new_shape)
    if (new_shape>=0) call eshape_select(new_shape)

  end subroutine cmd_eshape

+4 −1
Original line number Diff line number Diff line
@@ -54,6 +54,9 @@ module echo_shapes
  !> current echo shape type
  integer :: eshape_type = -1

  !> current echo sign (positive or negative)
  integer :: eshape_sign =  1

  !> triangular echo parameters
  real(kind=8) :: selector_width_fwhm = 0.19_DBL    !!### >> make this a parameter that is set using dlambda
                                                    !!### >> and/or by command
@@ -91,7 +94,7 @@ contains
    else
       this%nbin = [1,1]
    end if
    this%e_sign  = 1
    this%e_sign  = eshape_sign
    this%phase_sens_corr = 1
    this%spectrum_struct = spectrum
  end subroutine init_lambda_struct