Loading Makefile.version +1 −1 Original line number Diff line number Diff line Loading @@ -7,7 +7,7 @@ export PROJARCH=$(PROJECT)-$(VERSION_MAJOR).$(VERSION_MINOR) export VERSION_MAJOR=1 export VERSION_MINOR=0 export VERSION_RELEASE=4 export VERSION_RELEASE=5 git_rev=$(shell git rev-parse --short HEAD 2> /dev/null) ifeq "$(git_rev)" "" Loading sources/drspine.f90 +21 −8 Original line number Diff line number Diff line Loading @@ -1033,6 +1033,7 @@ CONTAINS end subroutine cmd_bins !> helper subroutine for cmd_histo and cmd_collect !! FIXME: (paz) unused, obsolete function subroutine auto_histogramming(xbins, role) type(double_bin_struct), intent(inout) :: xbins integer, intent(in) :: role Loading @@ -1044,6 +1045,7 @@ CONTAINS nsize = xbins%nbins nsize = 0 n = 0 do i=1, data_manager_size() if (is_valid_scan(data_scan(i), role=role)) then nsize = nsize + get_scan_struct_size(data_scan(i)) Loading Loading @@ -1256,9 +1258,10 @@ CONTAINS character(len=MAX_LINE_LENGTH) :: chrole !character(len=MAX_LINE_LENGTH) :: chgroup integer :: i, k integer :: inewtfac, inewphoff real(kind=DBL) :: tfac real(kind=DBL) :: phase_offset ! flags and "override parameters integer :: itfac, iphase_offset, isample_temp, iscattering_angle real(kind=DBL) :: tfac, phase_offset, sample_temp, scattering_angle ! logical :: fexists ! flag to check whether file exists call msg_info('drspine', '===> reading') Loading @@ -1266,8 +1269,10 @@ CONTAINS 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) tfac = get_named_value('tfac', tfac, itfac) phase_offset = get_named_value('phase_offset', phase_offset, iphase_offset) sample_temp = get_named_value('sample_temp', sample_temp, isample_temp) scattering_angle = get_named_value('scattering_angle', scattering_angle, iscattering_angle) ! make sure we have all the geometry right call configure_instrument_geometry() Loading Loading @@ -1323,10 +1328,18 @@ CONTAINS !else ! continue !endif if (inewtfac>0) & ! FIXME: perhaps we need to come up with a less laborous scheme ! override parameters, if set if (itfac>0) & call set_transmission_factor(data_scan(k), tfac) if (inewphoff>0) & if (iphase_offset>0) & call set_initial_phase_offset(data_scan(k), phase_offset) if (isample_temp>0) & call set_sample_temperature(data_scan(k), sample_temp) if (iscattering_angle>0) & call set_scattering_angle(data_scan(k), DEG2RAD(scattering_angle), instrument_parameters) ! check data, if required 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 Loading sources/read_data.f90 +35 −0 Original line number Diff line number Diff line Loading @@ -15,6 +15,7 @@ module read_data public read_echo_data, read_transmission_data !!---, norm_counter public set_initial_phase_offset, set_transmission_factor public set_sample_temperature, set_scattering_angle public check_scan_data, check_magnetic_fields Loading Loading @@ -1778,4 +1779,38 @@ contains end do end subroutine set_initial_phase_offset ! ============================================================================== ! set sample temperature ! ============================================================================== subroutine set_sample_temperature(scan_data, sample_temp) type(scan_struct), intent(inout) :: scan_data real(kind=DBL), intent(in) :: sample_temp ! integer :: iscan ! do iscan=1, scan_data%number_of_points call init_value(scan_data%scan_point(iscan)%technical%sample_temperature,& sample_temp) end do call init_value(scan_data%sample%temperature, sample_temp) end subroutine set_sample_temperature ! ============================================================================== ! set sample temperature ! ============================================================================== subroutine set_scattering_angle(scan_data, scattering_angle, instrument_config) type(scan_struct), intent(inout) :: scan_data real(kind=DBL), intent(in) :: scattering_angle type(config_struct),intent(in) :: instrument_config ! integer :: iscan ! do iscan=1, scan_data%number_of_points call init_value(scan_data%scan_point(iscan)%physics%scattering_angle,& scattering_angle) end do ! recalculate pixel angles and q-vectors call fill_geometry_data(scan_data, instrument_config) end subroutine set_scattering_angle end module read_data Loading
Makefile.version +1 −1 Original line number Diff line number Diff line Loading @@ -7,7 +7,7 @@ export PROJARCH=$(PROJECT)-$(VERSION_MAJOR).$(VERSION_MINOR) export VERSION_MAJOR=1 export VERSION_MINOR=0 export VERSION_RELEASE=4 export VERSION_RELEASE=5 git_rev=$(shell git rev-parse --short HEAD 2> /dev/null) ifeq "$(git_rev)" "" Loading
sources/drspine.f90 +21 −8 Original line number Diff line number Diff line Loading @@ -1033,6 +1033,7 @@ CONTAINS end subroutine cmd_bins !> helper subroutine for cmd_histo and cmd_collect !! FIXME: (paz) unused, obsolete function subroutine auto_histogramming(xbins, role) type(double_bin_struct), intent(inout) :: xbins integer, intent(in) :: role Loading @@ -1044,6 +1045,7 @@ CONTAINS nsize = xbins%nbins nsize = 0 n = 0 do i=1, data_manager_size() if (is_valid_scan(data_scan(i), role=role)) then nsize = nsize + get_scan_struct_size(data_scan(i)) Loading Loading @@ -1256,9 +1258,10 @@ CONTAINS character(len=MAX_LINE_LENGTH) :: chrole !character(len=MAX_LINE_LENGTH) :: chgroup integer :: i, k integer :: inewtfac, inewphoff real(kind=DBL) :: tfac real(kind=DBL) :: phase_offset ! flags and "override parameters integer :: itfac, iphase_offset, isample_temp, iscattering_angle real(kind=DBL) :: tfac, phase_offset, sample_temp, scattering_angle ! logical :: fexists ! flag to check whether file exists call msg_info('drspine', '===> reading') Loading @@ -1266,8 +1269,10 @@ CONTAINS 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) tfac = get_named_value('tfac', tfac, itfac) phase_offset = get_named_value('phase_offset', phase_offset, iphase_offset) sample_temp = get_named_value('sample_temp', sample_temp, isample_temp) scattering_angle = get_named_value('scattering_angle', scattering_angle, iscattering_angle) ! make sure we have all the geometry right call configure_instrument_geometry() Loading Loading @@ -1323,10 +1328,18 @@ CONTAINS !else ! continue !endif if (inewtfac>0) & ! FIXME: perhaps we need to come up with a less laborous scheme ! override parameters, if set if (itfac>0) & call set_transmission_factor(data_scan(k), tfac) if (inewphoff>0) & if (iphase_offset>0) & call set_initial_phase_offset(data_scan(k), phase_offset) if (isample_temp>0) & call set_sample_temperature(data_scan(k), sample_temp) if (iscattering_angle>0) & call set_scattering_angle(data_scan(k), DEG2RAD(scattering_angle), instrument_parameters) ! check data, if required 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 Loading
sources/read_data.f90 +35 −0 Original line number Diff line number Diff line Loading @@ -15,6 +15,7 @@ module read_data public read_echo_data, read_transmission_data !!---, norm_counter public set_initial_phase_offset, set_transmission_factor public set_sample_temperature, set_scattering_angle public check_scan_data, check_magnetic_fields Loading Loading @@ -1778,4 +1779,38 @@ contains end do end subroutine set_initial_phase_offset ! ============================================================================== ! set sample temperature ! ============================================================================== subroutine set_sample_temperature(scan_data, sample_temp) type(scan_struct), intent(inout) :: scan_data real(kind=DBL), intent(in) :: sample_temp ! integer :: iscan ! do iscan=1, scan_data%number_of_points call init_value(scan_data%scan_point(iscan)%technical%sample_temperature,& sample_temp) end do call init_value(scan_data%sample%temperature, sample_temp) end subroutine set_sample_temperature ! ============================================================================== ! set sample temperature ! ============================================================================== subroutine set_scattering_angle(scan_data, scattering_angle, instrument_config) type(scan_struct), intent(inout) :: scan_data real(kind=DBL), intent(in) :: scattering_angle type(config_struct),intent(in) :: instrument_config ! integer :: iscan ! do iscan=1, scan_data%number_of_points call init_value(scan_data%scan_point(iscan)%physics%scattering_angle,& scattering_angle) end do ! recalculate pixel angles and q-vectors call fill_geometry_data(scan_data, instrument_config) end subroutine set_scattering_angle end module read_data