complex_step_derivative_T0_T1 Function

private impure function complex_step_derivative_T0_T1(f, x, h) result(dfdx)

Calculates the derivative of a scalar-valued function w.r.t. a vector-valued variable using complex step differentiation.

Arguments

Type IntentOptional Attributes Name
private impure function f(z) result(fz)
Arguments
Type IntentOptional Attributes Name
complex(kind=rk), intent(in), dimension(:) :: z

vector complex variable

Return Value complex(kind=rk)

scalar complex function

real(kind=rk), intent(in), dimension(:) :: x

vector variable

real(kind=rk), intent(in) :: h

perturbation for complex step differentiation

Return Value real(kind=rk), dimension(size(x))

derivative of w.r.t.


Called by

proc~~complex_step_derivative_t0_t1~~CalledByGraph proc~complex_step_derivative_t0_t1 fordiff::complex_step_derivative_T0_T1 interface~derivative fordiff::derivative interface~derivative->proc~complex_step_derivative_t0_t1 program~test1 test1 program~test1->interface~derivative program~test2 test2 program~test2->interface~derivative program~test3 test3 program~test3->interface~derivative program~test4 test4 program~test4->interface~derivative program~test5 test5 program~test5->interface~derivative program~test6 test6 program~test6->interface~derivative

Source Code

   impure function complex_step_derivative_T0_T1(f, x, h) result(dfdx)
      real(rk), dimension(:), intent(in) :: x      !! vector variable
      real(rk),               intent(in) :: h      !! perturbation for complex step differentiation
      real(rk), dimension(size(x))       :: dfdx   !! derivative of \(f\) w.r.t. \(\mathbf{x}\)
      real(rk), dimension(size(x))       :: temp_x !! temporary vector variable
      integer                            :: i      !! loop index

      interface
         !! scalar-valued function to differentiate
         impure function f(z) result(fz)
            use kinds
            complex(rk), dimension(:), intent(in) :: z  !! vector complex variable
            complex(rk)                           :: fz !! scalar complex function
         end function f
      end interface

      if (abs(h)<tiny(0.0_rk)) error stop 'Division by zero. Please provide a non-zero value for h.'

      do i = 1, size(x)
         temp_x    = 0.0_rk
         temp_x(i) = x(i)
         dfdx(i) = aimag(f(cmplx(temp_x, h, rk))) / h
      end do

   end function complex_step_derivative_T0_T1