co2_sw_rtsafe.F90 coverage: 0.00 %func 0.00 %block
1) module co2_sw_rtsafe_module
2)
3) use PFLOTRAN_Constants_module
4)
5) implicit none
6)
7) #include "petsc/finclude/petscsys.h"
8)
9) contains
10)
11) ! ************************************************************************** !
12)
13) FUNCTION rtsafe(funcd,x1,x2,xacc)
14)
15) IMPLICIT NONE
16) PetscReal, INTENT(IN) :: x1,x2,xacc
17) PetscReal :: rtsafe
18)
19) ! INTERFACE
20) ! SUBROUTINE funcd(x,fval,fderiv)
21) ! IMPLICIT NONE
22) ! PetscReal, INTENT(IN) :: x
23) ! PetscReal, INTENT(OUT) :: fval,fderiv
24) ! END SUBROUTINE funcd
25) ! END INTERFACE
26)
27) external funcd
28) INTEGER, PARAMETER :: MAXIT=100
29) INTEGER :: j
30) PetscReal :: df,dx,dxold,f,fh,fl,temp,xh,xl
31)
32) call funcd(x1,fl,df)
33) call funcd(x2,fh,df)
34) if ((fl > 0.0 .and. fh > 0.0) .or. &
35) (fl < 0.0 .and. fh < 0.0)) &
36) print *, 'root must be bracketed in rtsafe'
37) if (fl == 0.0) then
38) rtsafe=x1
39) RETURN
40) else if (fh == 0.0) then
41) rtsafe=x2
42) RETURN
43) else if (fl < 0.0) then
44) xl=x1
45) xh=x2
46) else
47) xh=x1
48) xl=x2
49) end if
50) rtsafe=0.5*(x1+x2)
51) dxold=dabs(x2-x1)
52) dx=dxold
53) call funcd(rtsafe,f,df)
54) do j=1,MAXIT
55) if (((rtsafe-xh)*df-f)*((rtsafe-xl)*df-f) >= 0.0 .or. &
56) dabs(2.0*f) > dabs(dxold*df) ) then
57) dxold=dx
58) dx=0.5*(xh-xl)
59) rtsafe=xl+dx
60) if (xl == rtsafe) RETURN
61) else
62) dxold=dx
63) dx=f/df
64) temp=rtsafe
65) rtsafe=rtsafe-dx
66) if (temp == rtsafe) RETURN
67) end if
68) if (dabs(dx) < xacc) RETURN
69) call funcd(rtsafe,f,df)
70) if (f < 0.0) then
71) xl=rtsafe
72) else
73) xh=rtsafe
74) end if
75) end do
76) print *,'rtsafe: exceeded maximum iterations'
77) END FUNCTION rtsafe
78)
79) ! ************************************************************************** !
80)
81) subroutine bracket(func,x1,x2)
82)
83) implicit none
84)
85) PetscInt :: i,ifind
86) PetscReal :: fac,f1,f2,x1,x2,df
87)
88) external func
89)
90) fac = 1.2d0
91) call func(x1,f1,df)
92) call func(x2,f2,df)
93) ifind = 1
94) do i = 1, 200
95) if (f1*f2 < 0.d0) return
96) if (dabs(f1) < dabs(f2)) then
97) x1 = x1+fac*(x1-x2)
98) call func(x1,f1,df)
99) else
100) x2 = x2+fac*(x2-x1)
101) call func(x2,f2,df)
102) endif
103) enddo
104) ifind = 0
105) print *,'root bracket failed',x1,x2,f1,f2
106) return
107) end subroutine bracket
108) end module co2_sw_rtsafe_module