read_parameters_solver_settings_nonlinear_convergence Subroutine

subroutine read_parameters_solver_settings_nonlinear_convergence(convergences, json, key_base)

Arguments

Type IntentOptional Attributes Name
type(type_convergence_criteria) :: convergences
type(json_file), intent(inout) :: json
character(len=*), intent(in) :: key_base

Calls

proc~~read_parameters_solver_settings_nonlinear_convergence~~CallsGraph proc~read_parameters_solver_settings_nonlinear_convergence read_parameters_solver_settings_nonlinear_convergence destroy destroy proc~read_parameters_solver_settings_nonlinear_convergence->destroy get get proc~read_parameters_solver_settings_nonlinear_convergence->get log_warning log_warning proc~read_parameters_solver_settings_nonlinear_convergence->log_warning print_error_message print_error_message proc~read_parameters_solver_settings_nonlinear_convergence->print_error_message proc~error_message error_message proc~read_parameters_solver_settings_nonlinear_convergence->proc~error_message proc~join join proc~read_parameters_solver_settings_nonlinear_convergence->proc~join log_error log_error proc~error_message->log_error

Called by

proc~~read_parameters_solver_settings_nonlinear_convergence~~CalledByGraph proc~read_parameters_solver_settings_nonlinear_convergence read_parameters_solver_settings_nonlinear_convergence proc~read_parameters_solver_settings_nonlinear read_parameters_solver_settings_nonlinear proc~read_parameters_solver_settings_nonlinear->proc~read_parameters_solver_settings_nonlinear_convergence proc~read_parameters_solver_settings read_parameters_solver_settings proc~read_parameters_solver_settings->proc~read_parameters_solver_settings_nonlinear proc~inout_read_basic_parameters inout_read_basic_parameters proc~inout_read_basic_parameters->proc~read_parameters_solver_settings interface~inout_read_basic_parameters type_input%inout_read_basic_parameters interface~inout_read_basic_parameters->proc~inout_read_basic_parameters proc~initialize_type_input type_input%initialize_type_input proc~initialize_type_input->interface~inout_read_basic_parameters

Source Code

    subroutine read_parameters_solver_settings_nonlinear_convergence(convergences, json, key_base)
        implicit none
        type(type_convergence_criteria) :: convergences
        type(json_file), intent(inout) :: json
        character(*), intent(in) :: key_base

        logical :: found
        character(:), allocatable :: key

        key = join([key_base, criteria])
        call json%get(key, convergences%criteria, found)
        call json%print_error_message(output_unit)
        if (.not. found) then
            call json%destroy()
            call error_message(904, c_opt=key)
        else if (.not. any(valid_local_criteria_types(:) == convergences%criteria)) then
            call json%destroy()
            call error_message(905, c_opt=key)
        end if

        if (convergences%criteria == valid_local_criteria_types(3)) then
            key = join([key_base, logic])
            call json%get(key, convergences%logic, found)
            call json%print_error_message(output_unit)
            if (.not. found) then
                call global_logger%log_warning(message="Using default convergences logic 'and'.")
                convergences%logic = "and"
            else if (.not. any(valid_logic_types(:) == convergences%logic)) then
                call json%destroy()
                call error_message(905, c_opt=key)
            end if
        end if

        select case (convergences%criteria)
        case (valid_local_criteria_types(1)) ! absolute

            key = join([key_base, absolute_tolerance])
            call json%get(key, convergences%absolute_tolerance, found)
            call json%print_error_message(output_unit)
            if (.not. found) then
                call global_logger%log_warning(message="Using default absolute convergences of 1.0d-6.")
                convergences%absolute_tolerance = 1.0d-6
            else if (convergences%absolute_tolerance < 0.0d0) then
                call json%destroy()
                call error_message(905, c_opt=key)
            end if
        case (valid_local_criteria_types(2)) ! relative

            key = join([key_base, relative_tolerance])
            call json%get(key, convergences%relative_tolerance, found)
            call json%print_error_message(output_unit)
            if (.not. found) then
                call global_logger%log_warning(message="Using default relative convergences of 1.0d-6.")
                convergences%relative_tolerance = 1.0d-6
            else if (convergences%relative_tolerance < 0.0d0) then
                call json%destroy()
                call error_message(905, c_opt=key)
            end if
        case (valid_local_criteria_types(3)) ! absolute and relative
            key = join([key_base, absolute_tolerance])
            call json%get(key, convergences%absolute_tolerance, found)
            call json%print_error_message(output_unit)
            if (.not. found) then
                call global_logger%log_warning(message="Using default absolute convergences of 1.0d-6.")
                convergences%absolute_tolerance = 1.0d-6
            else if (convergences%absolute_tolerance < 0.0d0) then
                call json%destroy()
                call error_message(905, c_opt=key)
            end if
            key = join([key_base, relative_tolerance])
            call json%get(key, convergences%relative_tolerance, found)
            call json%print_error_message(output_unit)
            if (.not. found) then
                call global_logger%log_warning(message="Using default relative convergences of 1.0d-6.")
                convergences%relative_tolerance = 1.0d-6
            else if (convergences%relative_tolerance < 0.0d0) then
                call json%destroy()
                call error_message(905, c_opt=key)
            end if
        end select
    end subroutine read_parameters_solver_settings_nonlinear_convergence