1C> \ingroup nwad_tests 2C> @{ 3C> 4C> \brief test the implementation of the TANH function 5C> 6C> This is an NWAD unit test. The derivatives of TANH function are compared 7C> against analytic derivatives. The function is 1-dimensional as that is 8C> sufficient for this test. The input data set is randomly generated. 9C> 10 program test_tanh 11 use nwad3 12 implicit none 13 integer :: npt, i 14 parameter (npt = 100) 15 type(nwad_dble) :: x(npt), f(npt) 16 double precision fa(npt),dfa(npt),dfa2(npt),dfa3(npt), tol 17 double precision tmp 18 parameter( tol = 1.0d-10) 19 call random_seed 20 do i = 1, npt 21 call random_number(tmp) 22 tmp = tmp*0.5d0*acos(-1.0d0) 23 x(i) = active(tmp) 24 call submaxima_tanh(x(i)%d0,fa(i),dfa(i),dfa2(i),dfa3(i)) 25 call subad_tanh(x(i),f(i)) 26 if (abs((fa(i)-f(i)%d0)/(fa(i)+f(i)%d0)).gt.tol) then 27 write(*,*)"F : fail:",i,x(i)%d0,fa(i) 28 write(*,*)"F : fail:",i,x(i)%d0,f(i)%d0 29 write(*,*) 30 endif 31 if (abs((dfa(i)-f(i)%d1)/(dfa(i)+f(i)%d1)).gt.tol) then 32 write(*,*)"DF : fail:",i,x(i)%d0,dfa(i) 33 write(*,*)"DF : fail:",i,x(i)%d0,f(i)%d1 34 write(*,*) 35 endif 36 if (abs((dfa2(i)-f(i)%d2)/(dfa2(i)+f(i)%d2)).gt.tol) then 37 write(*,*)"DF2: fail:",i,x(i)%d0,dfa2(i) 38 write(*,*)"DF2: fail:",i,x(i)%d0,f(i)%d2 39 write(*,*) 40 endif 41 if (abs((dfa3(i)-f(i)%d3)/(dfa3(i)+f(i)%d3)).gt.tol) then 42 write(*,*)"DF3: fail:",i,x(i)%d0,dfa3(i) 43 write(*,*)"DF3: fail:",i,x(i)%d0,f(i)%d3 44 write(*,*) 45 endif 46 enddo 47 end 48C> 49C> \brief The test routine 50C> 51 subroutine subad_tanh(x,f) 52 use nwad3 53 implicit none 54 type(nwad_dble) :: x, f 55 f = tanh(sin(x)) 56 end 57C> @} 58c $Id$ 59