Commit 048ac7e6 authored by Zolnierczuk, Piotr's avatar Zolnierczuk, Piotr
Browse files

store/use initial phase offset

data_types.f90: physics%estimated_field_integral_asymmetry to store
                initial phase offset in Tesla-meters.
fit_data.f90:   removed optional phase_offset parameter in fit_data

read_data.f90:  set_initial_phase_offset deg=>Tesla-meter

drspine.f90:    allow to specify phase offset for each run on read,
                can be overriden in fit command
parent 6fd6cc37
Loading
Loading
Loading
Loading
+22 −10
Original line number Diff line number Diff line
@@ -1282,15 +1282,19 @@ CONTAINS
    character(len=MAX_PATH_LENGTH) :: filename
    character(len=MAX_LINE_LENGTH) :: chrole
    !character(len=MAX_LINE_LENGTH) :: chgroup
    real(kind=DBL) :: tfac
    integer        :: i, k
    integer        :: inewtfac
    integer        :: inewtfac, inewphoff
    real(kind=DBL) :: tfac
    real(kind=DBL) :: phase_offset

    call msg_info('drspine', '===> reading')
    chrole  = ' '
    tfac    = 1.0
    phase_offset = 0.0
    if ( parse_role_arg(chrole) < 0 ) return
    tfac    =  get_named_value('tfac', tfac, inewtfac)
    phase_offset = get_named_value('phase_offset', phase_offset, inewphoff)

    !
    call configure_instrument_geometry()

@@ -1341,6 +1345,8 @@ CONTAINS
       !endif
       if (inewtfac>0) &
            call set_transmission_factor(data_scan(k), tfac)
       if (inewphoff>0) &
            call set_initial_phase_offset(data_scan(k), phase_offset)
       if (len_trim(program_param%file_checkflag)>0) &
            call check_scan_data(data_scan(k), trim(program_param%file_checkflag))
       program_param%last_run = data_scan(k)%id
@@ -1634,7 +1640,7 @@ CONTAINS
  !-------------------------------------------------------------
  subroutine cmd_fit()
    implicit none
    integer :: i
    integer :: i, inewphoff
    integer :: fit_flag
    integer :: fit_run
    logical :: fit_resolution, fit_sample, fit_background
@@ -1658,7 +1664,7 @@ CONTAINS
       fit_sample     = .false.
       fit_background = .false.
    end if
    phase_offset = get_named_value('phase_offset',phase_offset, ier)
    phase_offset = get_named_value('phase_offset',phase_offset, inewphoff)

    ! now parse "what"
    cwhat = trim(chrnxt('fit', ier)) ! weird way of getting next parameter
@@ -1730,10 +1736,12 @@ CONTAINS
       call msg_info('fit', '===> fitting sample')
       do i=1, data_manager_size()
          if( has_role(data_scan(i)%role, ROLE_SAMPLE) ) then
            if (inewphoff>0) &
              call set_initial_phase_offset(data_scan(i), phase_offset)
            if( fit_flag == -1 ) then
              call fit_echo_data(data_scan(i), PHASE_USE, phase_offset=phase_offset )
              call fit_echo_data(data_scan(i), PHASE_USE)
            else
              call fit_echo_data(data_scan(i),fit_flag, phase_offset=phase_offset)
              call fit_echo_data(data_scan(i),fit_flag)
            endif
          endif
       end do
@@ -1744,10 +1752,12 @@ CONTAINS
       call msg_info('fit', '===> fitting background')
       do i=1, data_manager_size()
          if( has_role(data_scan(i)%role, ROLE_BACKGROUND) ) then
            if (inewphoff>0) &
              call set_initial_phase_offset(data_scan(i), phase_offset)
            if( fit_flag == -1 ) then
              call fit_echo_data(data_scan(i), PHASE_USE, phase_offset=phase_offset)
              call fit_echo_data(data_scan(i), PHASE_USE)
            else
              call fit_echo_data(data_scan(i), fit_flag, phase_offset=phase_offset)
              call fit_echo_data(data_scan(i), fit_flag)
            endif
          endif
       end do
@@ -1757,10 +1767,12 @@ CONTAINS
       call msg_info('fit', '===> fitting run')
       do i=1, data_manager_size()
          if ( data_scan(i)%id /= fit_run ) cycle
          if (inewphoff>0) &
              call set_initial_phase_offset(data_scan(i), phase_offset)
          if  ( fit_flag == -1 ) then
              call fit_echo_data(data_scan(i), PHASE_FIT_DEFAULT, phase_offset=phase_offset) ! FIXME: is this a right default
              call fit_echo_data(data_scan(i), PHASE_FIT_DEFAULT) ! FIXME: is this a right default
          else
              call fit_echo_data(data_scan(i), fit_flag, phase_offset=phase_offset)
              call fit_echo_data(data_scan(i), fit_flag)
          end if
       end do
    end if
+3 −1
Original line number Diff line number Diff line
@@ -272,7 +272,9 @@ contains
    write(out_unit,*) '### physics.phase_per_ampere_per_lambda.[deg/amp/A]:', &
         trim(cformat_value(div_value(physics%phase_per_ampere_per_lambda, &
                                      DEGREE/ANGSTROEM)))

    write(out_unit,*) '### physics.estimated_field_integral_assymetry [uTm]:', &
         trim(cformat_value(div_value(physics%estimated_field_integral_asymmetry,&
                                      UTESLA)))
    write(out_unit,'(1x,a,2(1x,g13.5))') '### physics.lambda(min,max).:', &
         physics%lambda_min, physics%lambda_max

+11 −6
Original line number Diff line number Diff line
@@ -31,14 +31,12 @@ contains
  !!
  !! @param [in,out]      data         - scan_struct data (a @f$ \tau @f$ scan)
  !! @param [in]          fit_flag     - fitting level flag @see "phase_" flags in drspine_parameters
  !! @param [in,optional] phase_offset - initial phase offset
  !!
  !! @description
  !! This is user callable subroutine - public interface to fit_data module.
  subroutine fit_echo_data(data, fit_flag, phase_offset)
  subroutine fit_echo_data(data, fit_flag)
    type(scan_struct), target, intent(inout) :: data
    integer, intent(in)                      :: fit_flag
    real(kind=DBL), intent(in), optional     :: phase_offset
    !
    type(scan_data_struct), pointer :: sig_point, ref_point
    integer :: iscan
@@ -49,11 +47,18 @@ contains
         ', flag='//trim(cformat_fit_flag(fit_flag)))

    djoffset = 0.0d0
    if ( present(phase_offset) .and. data%number_of_points>0 ) &
        djoffset = get_field_integral(deg2rad(phase_offset), data%scan_point(1)%physics%lambda) !just need a wavelength
    if ( data%number_of_points>0 ) then
        associate(phys => data%scan_point(1)%physics )
        if ( phys%estimated_field_integral_asymmetry%state /= VAL_UNDEFINED ) then
            djoffset = phys%estimated_field_integral_asymmetry%value
            call msg_debug('fit_echo_data', &
                msg_fmt("('initial phase offset: ', g13.6)", rad2deg(get_precession_phase(djoffset, phys%lambda))))
        end if
        end associate
    end if
    ref_point => null()
    do iscan=1, data%number_of_points
    call msg_debug('fit_echo_data', msg_fmt("('scan point[0]: ', i4)", iscan))
    do iscan=1, data%number_of_points
       sig_point => data%scan_point(iscan)

       if ( associated(sig_point%matching_ref) ) ref_point => sig_point%matching_ref
+3 −3
Original line number Diff line number Diff line
@@ -311,8 +311,8 @@ MODULE new_com
  public  :: parse
  public  :: errsig
  ! private :: pushn ! (PAZ) not-used
  public  :: get1
  public  :: getvec
  ! public  :: get1  ! (PAZ) not-used
  ! public  :: getvec! (PAZ) not-used
  public  :: getnva
  public  :: rotavc
  public  :: unrot
@@ -320,7 +320,7 @@ MODULE new_com
  public  :: matmux
  public  :: matvec
  public  :: mattxp
  public  :: matcpy
  ! public  :: matcpy! (PAZ) not-used
  public  :: matone
  public  :: axrot
  public  :: cross
+21 −4
Original line number Diff line number Diff line
@@ -14,7 +14,8 @@ module read_data
  private

  public read_echo_data, read_transmission_data !!---, norm_counter
  public set_transmission_factor, check_scan_data, check_magnetic_fields
  public set_initial_phase_offset, set_transmission_factor
  public check_scan_data, check_magnetic_fields


!!---  integer :: norm_counter = -1 ! use default value
@@ -1749,8 +1750,7 @@ contains
    !
    integer :: iscan
    integer :: it, ix, iy!, ip


    !
    do iscan=1, scan_data%number_of_points
       associate(cur_point => scan_data%scan_point(iscan))
       do ix=1, cur_point%no_xpix
@@ -1762,7 +1762,24 @@ contains
       end do
       end associate
    end do

  end subroutine set_transmission_factor

  ! ==============================================================================
  ! set initial phase offset for the scan
  ! and convert it from degrees to Tesla-meters
  ! ==============================================================================
  subroutine set_initial_phase_offset(scan_data, phase_offset)
    type(scan_struct), intent(inout) :: scan_data
    real(kind=DBL), intent(in) :: phase_offset
    !
    integer :: iscan
    !
    do iscan=1, scan_data%number_of_points
       associate(phys => scan_data%scan_point(iscan)%physics)
        call init_value(phys%estimated_field_integral_asymmetry, &
                get_field_integral(deg2rad(phase_offset), phys%lambda))
       end associate
    end do
  end subroutine set_initial_phase_offset

end module read_data