nonlin 1.5.2
A library that provides routines to compute the solutions to systems of nonlinear equations.
Loading...
Searching...
No Matches
nonlin_fcnnvar_helper.f90
1! nonlin_fcnnvar_helper.f90
2
3submodule(nonlin_core) nonlin_fcnnvar_helper
4 implicit none
5contains
6! ------------------------------------------------------------------------------
7 module function fnh_fcn(this, x) result(f)
8 class(fcnnvar_helper), intent(in) :: this
9 real(real64), intent(in), dimension(:) :: x
10 real(real64) :: f
11 if (associated(this%m_fcn)) then
12 f = this%m_fcn(x)
13 end if
14 end function
15
16! ------------------------------------------------------------------------------
17 module function fnh_is_fcn_defined(this) result(x)
18 class(fcnnvar_helper), intent(in) :: this
19 logical :: x
20 x = associated(this%m_fcn)
21 end function
22
23! ------------------------------------------------------------------------------
24 module subroutine fnh_set_fcn(this, fcn, nvar)
25 class(fcnnvar_helper), intent(inout) :: this
26 procedure(fcnnvar), intent(in), pointer :: fcn
27 integer(int32), intent(in) :: nvar
28 this%m_fcn => fcn
29 this%m_nvar = nvar
30 end subroutine
31
32! ------------------------------------------------------------------------------
33 module function fnh_get_nvar(this) result(n)
34 class(fcnnvar_helper), intent(in) :: this
35 integer(int32) :: n
36 n = this%m_nvar
37 end function
38
39! ------------------------------------------------------------------------------
40 module subroutine fnh_set_grad(this, fcn)
41 class(fcnnvar_helper), intent(inout) :: this
42 procedure(gradientfcn), pointer, intent(in) :: fcn
43 this%m_grad => fcn
44 end subroutine
45
46! ------------------------------------------------------------------------------
47 module function fnh_is_grad_defined(this) result(x)
48 class(fcnnvar_helper), intent(in) :: this
49 logical :: x
50 x = associated(this%m_grad)
51 end function
52
53! ------------------------------------------------------------------------------
54 module subroutine fnh_grad_fcn(this, x, g, fv, err)
55 ! Arguments
56 class(fcnnvar_helper), intent(in) :: this
57 real(real64), intent(inout), dimension(:) :: x
58 real(real64), intent(out), dimension(:) :: g
59 real(real64), intent(in), optional :: fv
60 integer(int32), intent(out), optional :: err
61
62 ! Parameters
63 real(real64), parameter :: zero = 0.0d0
64
65 ! Local Variables
66 integer(int32) :: j, n, flag
67 real(real64) :: eps, epsmch, h, temp, f, f1
68
69 ! Initialization
70 if (present(err)) err = 0
71 ! n = this%m_nvar
72 n = this%get_variable_count()
73
74 ! Input Checking
75 flag = 0
76 if (size(x) /= n) then
77 flag = 2
78 else if (size(g) /= n) then
79 flag = 3
80 end if
81 if (flag /= 0) then
82 ! ERROR: Incorrectly sized input arrays
83 if (present(err)) err = flag
84 return
85 end if
86
87 ! Process
88 if (.not.this%is_fcn_defined()) return
89 if (this%is_gradient_defined()) then
90 ! Call the user-defined gradient routine
91 call this%m_grad(x, g)
92 else
93 ! Compute the gradient via finite differences
94 if (present(fv)) then
95 f = fv
96 else
97 f = this%fcn(x)
98 end if
99
100 ! Establish step size factors
101 epsmch = epsilon(epsmch)
102 eps = sqrt(epsmch)
103
104 ! Compute the derivatives
105 do j = 1, n
106 temp = x(j)
107 h = eps * abs(temp)
108 if (h == zero) h = eps
109 x(j) = temp + h
110 f1 = this%fcn(x)
111 x(j) = temp
112 g(j) = (f1 - f) / h
113 end do
114 end if
115 end subroutine
116
117
118end submodule
nonlin_core