initialize_the_solver Subroutine

public subroutine initialize_the_solver(me, config_file_name, x)

Initialize the mission.

This is the "forward-backward" formulation.

Type Bound

my_solver_type

Arguments

Type IntentOptional Attributes Name
class(my_solver_type), intent(inout) :: me
character(len=*), intent(in) :: config_file_name

the config file to read

real(kind=wp), intent(out), optional, dimension(:), allocatable :: x

initial guess


Calls

proc~~initialize_the_solver~~CallsGraph proc~initialize_the_solver my_solver_type%initialize_the_solver initialize initialize proc~initialize_the_solver->initialize proc~define_problem_size mission_type%define_problem_size proc~initialize_the_solver->proc~define_problem_size proc~get_sparsity_pattern mission_type%get_sparsity_pattern proc~initialize_the_solver->proc~get_sparsity_pattern proc~initialize_the_mission mission_type%initialize_the_mission proc~initialize_the_solver->proc~initialize_the_mission proc~read_config_file my_solver_type%read_config_file proc~initialize_the_solver->proc~read_config_file status status proc~initialize_the_solver->status proc~get_sparsity_pattern->proc~define_problem_size proc~initialize_the_mission->initialize proc~initialize_the_mission->proc~define_problem_size proc~initialize_the_mission->proc~get_sparsity_pattern initialize_splinded_ephemeris initialize_splinded_ephemeris proc~initialize_the_mission->initialize_splinded_ephemeris magnitude magnitude proc~initialize_the_mission->magnitude proc~get_problem_arrays mission_type%get_problem_arrays proc~initialize_the_mission->proc~get_problem_arrays proc~get_scales_from_segs mission_type%get_scales_from_segs proc~initialize_the_mission->proc~get_scales_from_segs proc~set_segment_inputs segment%set_segment_inputs proc~initialize_the_mission->proc~set_segment_inputs set_dpert set_dpert proc~initialize_the_mission->set_dpert set_sparsity_pattern set_sparsity_pattern proc~initialize_the_mission->set_sparsity_pattern proc~read_config_file->initialize evaluate evaluate proc~read_config_file->evaluate none~get config_file%get proc~read_config_file->none~get proc~generate_patch_points mission_type%generate_patch_points proc~read_config_file->proc~generate_patch_points proc~read_epoch read_epoch proc~read_config_file->proc~read_epoch proc~update_epoch mission_type%update_epoch proc~read_config_file->proc~update_epoch proc~get_char config_file%get_char none~get->proc~get_char proc~get_int config_file%get_int none~get->proc~get_int proc~get_logical config_file%get_logical none~get->proc~get_logical proc~get_real config_file%get_real none~get->proc~get_real proc~get_real_vec config_file%get_real_vec none~get->proc~get_real_vec proc~generate_patch_points->initialize compute_crtpb_parameter compute_crtpb_parameter proc~generate_patch_points->compute_crtpb_parameter first_call first_call proc~generate_patch_points->first_call integrate integrate proc~generate_patch_points->integrate fill_vector fill_vector proc~get_problem_arrays->fill_vector proc~get_scales_from_segs->proc~define_problem_size proc~get_scales_from_segs->fill_vector proc~read_epoch->none~get jd_to_et jd_to_et proc~read_epoch->jd_to_et julian_date julian_date proc~read_epoch->julian_date icrf_frame icrf_frame proc~set_segment_inputs->icrf_frame transform transform proc~set_segment_inputs->transform two_body_rotating_frame two_body_rotating_frame proc~set_segment_inputs->two_body_rotating_frame et_to_jd et_to_jd proc~update_epoch->et_to_jd proc~update_epoch->jd_to_et proc~update_epoch->julian_date julian_date_to_calendar_date julian_date_to_calendar_date proc~update_epoch->julian_date_to_calendar_date get get proc~get_char->get proc~get_int->get proc~get_logical->get proc~get_real->get proc~get_real_vec->get

Called by

proc~~initialize_the_solver~~CalledByGraph proc~initialize_the_solver my_solver_type%initialize_the_solver proc~halo_solver_main halo_solver_main proc~halo_solver_main->proc~initialize_the_solver

Source Code

    subroutine initialize_the_solver(me,config_file_name,x)

    implicit none

    class(my_solver_type),intent(inout) :: me
    character(len=*),intent(in) :: config_file_name !! the config file to read
    real(wp),dimension(:),allocatable,intent(out),optional :: x !! initial guess

    integer :: n          !! number of opt vars for the solver
    integer :: m          !! number of equality constraints for the solver
    logical :: status_ok  !! status flag for solver initialization
    integer :: istat      !! status code from solver
    integer,dimension(:),allocatable :: irow,icol !! sparsity pattern

    ! note: we don't know the problem size
    ! until we read the config file.

    ! first we have to read the config file:
    ! this will populate some of the mission variables
    call me%read_config_file(config_file_name)

    ! initialize the mission
    call me%mission%init(x)
    call me%mission%define_problem_size(n,m)

    ! initialize the solver:
    select case (me%mission%solver_mode)
    case(1)
        ! dense - uses lapack to solve the linear system
        call me%initialize(     n                = n,            &
                                m                = m,            &
                                max_iter         = 100,          & ! maximum number of iteration
                                func             = halo_func,    &
                                grad             = halo_grad,    &
                                tol              = me%mission%nlesolver_tol,    & ! tolerance
                                step_mode        = 4,            & ! 3-point "line search" (2 intervals)
                                n_intervals      = 2,            & ! number of intervals for step_mode=4
                                alpha_min = 0.2_wp, &
                                alpha_max = 0.8_wp, &
                                use_broyden      = .false.,      & ! broyden update
                                !use_broyden=.true.,broyden_update_n=10, & ! ... test ...
                                export_iteration = halo_export   )

    case(5)
        ! this is using the qr_mumps solver as a user-defined solver to nlesolver-fortran.
        ! the solver is defined in qrm_solver. you must use the WITH_QRMUMPS preprocessor
        ! directive to use this method and link the code with the appropriate libraries.

        call me%mission%get_sparsity_pattern(irow,icol) ! it's already been computed, but for now, just compute it again for this call
        call me%initialize(     n                = n,            &
                                m                = m,            &
                                max_iter         = 100,          & ! maximum number of iteration
                                func             = halo_func,    &
                                grad_sparse      = halo_grad_sparse,    &
                                tol              = me%mission%nlesolver_tol,    & ! tolerance
                                step_mode        = 4,            & ! 3-point "line search" (2 intervals)
                                n_intervals      = 2,            & ! number of intervals for step_mode=4
                                alpha_min        = 0.2_wp, &
                                alpha_max        = 0.8_wp, &
                                use_broyden      = .false.,      & ! broyden update
                                sparsity_mode    = me%mission%solver_mode, &  ! use a sparse solver
                                custom_solver_sparse = qrm_solver, &  ! the qr_mumps solver wrapper
                                irow             = irow, &  ! sparsity pattern
                                icol             = icol, &
                                export_iteration = halo_export   )
    case (2:4)
        ! varions sparse options available in nlesolver-fortran
        call me%mission%get_sparsity_pattern(irow,icol) ! it's already been computed, but for now, just compute it again for this call
        call me%initialize(     n                = n,            &
                                m                = m,            &
                                max_iter         = 100,          & ! maximum number of iteration
                                func             = halo_func,    &
                                grad_sparse      = halo_grad_sparse,    &
                                tol              = me%mission%nlesolver_tol,    & ! tolerance
                                ! step_mode = 1,& ! TEST TEST TEST
                                ! alpha = 1.0_wp,&
                                step_mode        = 4,            & ! 3-point "line search" (2 intervals)
                                n_intervals      = 2,            & ! number of intervals for step_mode=4
                                use_broyden      = .false.,      & ! broyden update
                                !use_broyden=.true.,broyden_update_n=10, & ! ... test ...
                                export_iteration = halo_export,  &
                                sparsity_mode = me%mission%solver_mode, &  ! use a sparse solver
                                atol          = 1.0e-12_wp,&  ! relative error in definition of `A`
                                btol          = 1.0e-12_wp,&  ! relative error in definition of `b`
!                                damp          = 0.00001_wp, & !  TEST: LSQR damp factor !
                                damp          = 0.0_wp, & !  TEST: LSQR damp factor !
   !                             damp          = 0.1_wp, & !  TEST: LSQR damp factor !
                                itnlim        = 1000000, &  ! max iterations
                                irow          = irow, &  ! sparsity pattern
                                icol          = icol, &
                                lusol_method = 0     ) ! test

    case default
        error stop 'invalid solver_mode'
    end select

    call me%status(istat=istat)
    status_ok = istat == 0
    if (.not. status_ok) then
        write(*,*) 'istat = ', istat
        error stop 'error in initialize_the_solver'
    end if

    end subroutine initialize_the_solver