example3 Program

Uses

  • program~~example3~~UsesGraph program~example3 example3 module~example3_model example3_model program~example3->module~example3_model module~odrpack odrpack program~example3->module~odrpack module~odrpack_kinds odrpack_kinds program~example3->module~odrpack_kinds module~example3_model->module~odrpack_kinds module~odrpack->module~odrpack_kinds iso_fortran_env iso_fortran_env module~odrpack->iso_fortran_env module~odrpack_kinds->iso_fortran_env

Explicit ODR job, with parameter bounds, weights, delta initialized by user, and central difference derivatives.


Calls

program~~example3~~CallsGraph program~example3 example3 proc~odr odr program~example3->proc~odr interface~ddot ddot proc~odr->interface~ddot proc~access_workspace access_workspace proc~odr->proc~access_workspace proc~check_inputs check_inputs proc~odr->proc~check_inputs proc~check_jac check_jac proc~odr->proc~check_jac proc~derstep derstep proc~odr->proc~derstep proc~eta_fcn eta_fcn proc~odr->proc~eta_fcn proc~eval_jac eval_jac proc~odr->proc~eval_jac proc~fctrw fctrw proc~odr->proc~fctrw proc~init_work init_work proc~odr->proc~init_work proc~loc_iwork loc_iwork proc~odr->proc~loc_iwork proc~loc_rwork loc_rwork proc~odr->proc~loc_rwork proc~move_beta move_beta proc~odr->proc~move_beta proc~pack_vec pack_vec proc~odr->proc~pack_vec proc~print_error_inputs print_error_inputs proc~odr->proc~print_error_inputs proc~print_errors print_errors proc~odr->proc~print_errors proc~print_header print_header proc~odr->proc~print_header proc~print_reports print_reports proc~odr->proc~print_reports proc~scale_mat scale_mat proc~odr->proc~scale_mat proc~select_row select_row proc~odr->proc~select_row proc~set_flags set_flags proc~odr->proc~set_flags proc~trust_region trust_region proc~odr->proc~trust_region proc~unpack_vec unpack_vec proc~odr->proc~unpack_vec proc~vcv_beta vcv_beta proc~odr->proc~vcv_beta proc~workspace_dimensions workspace_dimensions proc~odr->proc~workspace_dimensions proc~access_workspace->proc~loc_iwork proc~access_workspace->proc~loc_rwork proc~check_jac_value check_jac_value proc~check_jac->proc~check_jac_value proc~hstep hstep proc~check_jac->proc~hstep proc~derstep->proc~hstep proc~eval_jac->proc~scale_mat proc~eval_jac->proc~unpack_vec proc~jac_cdiff jac_cdiff proc~eval_jac->proc~jac_cdiff proc~jac_fwdiff jac_fwdiff proc~eval_jac->proc~jac_fwdiff proc~set_ifix set_ifix proc~eval_jac->proc~set_ifix proc~fctr fctr proc~fctrw->proc~fctr proc~init_work->proc~set_flags proc~scale_beta scale_beta proc~init_work->proc~scale_beta proc~scale_delta scale_delta proc~init_work->proc~scale_delta proc~move_beta->proc~hstep proc~print_errors->proc~print_error_inputs proc~print_errors->proc~print_header proc~print_error_derivative print_error_derivative proc~print_errors->proc~print_error_derivative proc~print_error_fcn print_error_fcn proc~print_errors->proc~print_error_fcn proc~print_reports->proc~print_header proc~print_reports->proc~set_flags proc~print_report_final print_report_final proc~print_reports->proc~print_report_final proc~print_report_initial print_report_initial proc~print_reports->proc~print_report_initial proc~print_report_iter print_report_iter proc~print_reports->proc~print_report_iter proc~trust_region->proc~scale_mat proc~lcstep lcstep proc~trust_region->proc~lcstep proc~scale_vec scale_vec proc~trust_region->proc~scale_vec dpodi dpodi proc~vcv_beta->dpodi proc~vcv_beta->proc~lcstep proc~check_jac_curv check_jac_curv proc~check_jac_value->proc~check_jac_curv proc~check_jac_zero check_jac_zero proc~check_jac_value->proc~check_jac_zero proc~fpvb fpvb proc~check_jac_value->proc~fpvb proc~fpvd fpvd proc~check_jac_value->proc~fpvd proc~jac_cdiff->proc~derstep proc~jac_cdiff->proc~hstep proc~jac_fwdiff->proc~derstep proc~jac_fwdiff->proc~hstep proc~lcstep->proc~scale_mat proc~lcstep->proc~fctr dchex dchex proc~lcstep->dchex dqrdc dqrdc proc~lcstep->dqrdc dqrsl dqrsl proc~lcstep->dqrsl dtrco dtrco proc~lcstep->dtrco dtrsl dtrsl proc~lcstep->dtrsl interface~drot drot proc~lcstep->interface~drot interface~drotg drotg proc~lcstep->interface~drotg proc~esubi esubi proc~lcstep->proc~esubi proc~solve_trl solve_trl proc~lcstep->proc~solve_trl proc~vevtr vevtr proc~lcstep->proc~vevtr proc~ppf_tstudent ppf_tstudent proc~print_report_final->proc~ppf_tstudent proc~print_report_initial->proc~hstep proc~check_jac_curv->proc~fpvb proc~check_jac_curv->proc~fpvd proc~check_jac_fp check_jac_fp proc~check_jac_curv->proc~check_jac_fp proc~check_jac_zero->proc~fpvb proc~check_jac_zero->proc~fpvd proc~ppf_normal ppf_normal proc~ppf_tstudent->proc~ppf_normal proc~vevtr->proc~solve_trl proc~check_jac_fp->proc~fpvb proc~check_jac_fp->proc~fpvd

Variables

Type Attributes Name Initial
integer :: i
integer :: iprint
integer :: j
integer :: job
integer :: lundata
integer :: lunrpt
integer :: m
integer :: n
integer :: np
integer :: q
integer, allocatable :: ifixx(:,:)
real(kind=wp), allocatable :: beta(:)
real(kind=wp), allocatable :: x(:,:)
real(kind=wp), allocatable :: y(:,:)
real(kind=wp), allocatable :: wd(:,:,:)
real(kind=wp), allocatable :: we(:,:,:)
real(kind=wp), allocatable :: delta(:,:)

Source Code

program example3
!! Explicit ODR job, with parameter bounds, weights, `delta` initialized by user, and central
!! difference derivatives.

   use odrpack_kinds, only: wp
   use odrpack, only: odr
   use example3_model, only: fcn
   implicit none

   ! Variable declarations
   integer :: i, iprint, j, job, lundata, lunrpt, m, n, np, q
   integer, allocatable :: ifixx(:, :)
   real(kind=wp), allocatable :: beta(:), x(:, :), y(:, :), wd(:, :, :), we(:, :, :), &
                                 delta(:, :)

   ! Set up report files
   open (newunit=lunrpt, file='./example/report3.dat')

   ! Read problem dimensions
   open (newunit=lundata, file='./example/data3.dat')
   read (lundata, fmt=*) n, m, np, q

   ! Allocate arrays
   allocate (beta(np), x(n, m), y(n, q), delta(n, m), we(n, q, q), wd(n, m, m), ifixx(n, m))

   ! Read problem data
   read (lundata, fmt=*) (beta(i), i=1, np)
   do i = 1, n
      read (lundata, fmt=*) (x(i, j), j=1, m), (y(i, j), j=1, q)
   end do
   close (lundata)

   ! Specify task as explicit orthogonal distance regression
   !       With central difference derivatives
   !       Covariance matrix constructed with recomputed derivatives
   !       `delta` initialized by user
   !       Not a restart
   ! And indicate long initial report
   !       No iteration reports
   !       Long final report
   job = 01010
   iprint = 2002

   ! Initialize `delta`, and specify first decade of frequencies as fixed
   do i = 1, n
      if (x(i, 1) < 100.0E0_wp) then
         delta(i, 1) = 0.0E0_wp
         ifixx(i, 1) = 0
      elseif (x(i, 1) <= 150.0E0_wp) then
         delta(i, 1) = 0.0E0_wp
         ifixx(i, 1) = 1
      elseif (x(i, 1) <= 1000.0E0_wp) then
         delta(i, 1) = 25.0E0_wp
         ifixx(i, 1) = 1
      elseif (x(i, 1) <= 10000.0E0_wp) then
         delta(i, 1) = 560.0E0_wp
         ifixx(i, 1) = 1
      elseif (x(i, 1) <= 100000.0E0_wp) then
         delta(i, 1) = 9500.0E0_wp
         ifixx(i, 1) = 1
      else
         delta(i, 1) = 144000.0E0_wp
         ifixx(i, 1) = 1
      end if
   end do

   ! Set weights
   do i = 1, n
      if (x(i, 1) == 100.0E0_wp .or. x(i, 1) == 150.0E0_wp) then
         we(i, 1, 1) = 0.0E0_wp
         we(i, 1, 2) = 0.0E0_wp
         we(i, 2, 1) = 0.0E0_wp
         we(i, 2, 2) = 0.0E0_wp
      else
         we(i, 1, 1) = 559.6E0_wp
         we(i, 1, 2) = -1634.0E0_wp
         we(i, 2, 1) = -1634.0E0_wp
         we(i, 2, 2) = 8397.0E0_wp
      end if
      wd(i, 1, 1) = (1.0E-4_wp)/(x(i, 1)**2)
   end do

   ! Compute solution
   call odr(fcn, n, m, q, np, beta, y, x, &
            delta=delta, we=we, wd=wd, ifixx=ifixx, &
            job=job, iprint=iprint, lunerr=lunrpt, lunrpt=lunrpt)

   close (lunrpt)

end program example3