modified_quasi_fd_newton_method_T0 Subroutine

private impure subroutine modified_quasi_fd_newton_method_T0(this, F, x0, x_sol)

Arguments

Type IntentOptional Attributes Name
class(nlsolver), intent(inout) :: this
procedure(Fun8) :: F
real(kind=rk), intent(in) :: x0
real(kind=rk), intent(out) :: x_sol

Calls

proc~~modified_quasi_fd_newton_method_t0~~CallsGraph proc~modified_quasi_fd_newton_method_t0 forsolver::modified_quasi_fd_newton_method_T0 derivative derivative proc~modified_quasi_fd_newton_method_t0->derivative

Called by

proc~~modified_quasi_fd_newton_method_t0~~CalledByGraph proc~modified_quasi_fd_newton_method_t0 forsolver::modified_quasi_fd_newton_method_T0 proc~newton_rel_t0 forsolver::nlsolver%newton_rel_T0 proc~newton_rel_t0->proc~modified_quasi_fd_newton_method_t0 none~solve forsolver::nlsolver%solve none~solve->proc~newton_rel_t0 program~test_solver10 test_solver10 program~test_solver10->none~solve program~test_solver11 test_solver11 program~test_solver11->none~solve program~test_solver12 test_solver12 program~test_solver12->none~solve program~test_solver13 test_solver13 program~test_solver13->none~solve program~test_solver14 test_solver14 program~test_solver14->none~solve program~test_solver3 test_solver3 program~test_solver3->none~solve program~test_solver4 test_solver4 program~test_solver4->none~solve program~test_solver5 test_solver5 program~test_solver5->none~solve program~test_solver6 test_solver6 program~test_solver6->none~solve program~test_solver7 test_solver7 program~test_solver7->none~solve program~test_solver8 test_solver8 program~test_solver8->none~solve program~test_solver9 test_solver9 program~test_solver9->none~solve

Source Code

   impure subroutine modified_quasi_fd_newton_method_T0(this, F, x0,  x_sol)
      interface
         impure function Fun8(x) result(res)
            import rk
            real(rk), intent(in) :: x
            real(rk)             :: res
         end function Fun8
      end interface

      procedure(Fun8) :: F

      class(nlsolver), intent(inout) :: this
      real(rk),        intent(in)    :: x0
      real(rk),        intent(out)   :: x_sol
      real(rk)                       :: xk
      real(rk)                       :: F_val
      real(rk)                       :: dFdx_val
      real(rk)                       :: criteriaFun
      integer                        :: k
      logical                        :: convergenz
      real(rk)                       :: pk
      real(rk)                       :: qk
      real(rk)                       :: dk
      real(rk)                       :: alphak

      k          = 0
      xk         = x0
      convergenz = .false.

      do while (.not. convergenz .and. k < this%maxit)
         F_val = F(xk)
         if (mod(k, this%nmp) == 0) dFdx_val = derivative(f=F, x=xk, h=this%fdm_tol, method=this%fdm_method)

         criteriaFun = abs(F_val)

         if (this%verbosity == 1) then
            print '(g0, f12.4, 4x, e12.4, 4x, e12.4)', k, xk, F_val, dFdx_val
         end if

         if (criteriaFun <= this%TolFun) then
            convergenz = .true.
            x_sol      = xk
            return
         else
            dk = - F_val / dFdx_val
            alphak = 1.0_rk
            xk = xk + alphak*dk
            k = k + 1
         end if
      end do
   end subroutine modified_quasi_fd_newton_method_T0