Commit 31f2c259c63cf0ca76232d3a16f53b27036658d1
1 parent
20ce530d4f
Exists in
master
and in
2 other branches
Added test example for fvn_lm
git-svn-id: https://lxsd.femto-st.fr/svn/fvn@79 b657c933-2333-4658-acf2-d3c7c2708721
Showing 3 changed files with 91 additions and 1 deletions Inline Diff
fvn_test/Makefile
1 | 1 | |||
include $(BTREE)/Make.inc | 2 | 2 | include $(BTREE)/Make.inc | |
3 | 3 | |||
programs = test_fac$(exext) test_matinv$(exext) test_specfunc$(exext) \ | 4 | 4 | programs = test_fac$(exext) test_matinv$(exext) test_specfunc$(exext) \ | |
test_det$(exext) test_matcon$(exext) test_matev$(exext) test_inter1d$(exext) \ | 5 | 5 | test_det$(exext) test_matcon$(exext) test_matev$(exext) test_inter1d$(exext) \ | |
test_inter2d$(exext) test_inter3d$(exext) test_akima$(exext) test_lsp$(exext) test_muller$(exext) \ | 6 | 6 | test_inter2d$(exext) test_inter3d$(exext) test_akima$(exext) test_lsp$(exext) test_muller$(exext) \ | |
test_integ$(exext) test_bsyn$(exext) test_bsjn$(exext) test_bskn$(exext) test_bsin$(exext) test_operators$(exext) test_ze1$(exext) \ | 7 | 7 | test_integ$(exext) test_bsyn$(exext) test_bsjn$(exext) test_bskn$(exext) test_bsin$(exext) test_operators$(exext) test_ze1$(exext) \ | |
test_besri$(exext) test_besrj$(exext) test_bestime$(exext) \ | 8 | 8 | test_besri$(exext) test_besrj$(exext) test_bestime$(exext) \ | |
test_sparse_zi$(exext) test_sparse_zl$(exext) test_sparse_di$(exext) test_sparse_dl$(exext) | 9 | 9 | test_sparse_zi$(exext) test_sparse_zl$(exext) test_sparse_di$(exext) test_sparse_dl$(exext) \ | |
10 | test_lm$(exext) | |||
10 | 11 | |||
prog:$(programs) | 11 | 12 | prog:$(programs) | |
12 | 13 | |||
clean: | 13 | 14 | clean: | |
rm -f {*.o,*.oo,*.ipo,*.a,*.mod,*.dat} | 14 | 15 | rm -f {*.o,*.oo,*.ipo,*.a,*.mod,*.dat} | |
rm -f $(programs) | 15 | 16 | rm -f $(programs) | |
16 | 17 | |||
%$(exext): %.o | 17 | 18 | %$(exext): %.o | |
$(LINK) $(LINKFLAGS) $< init_random_seed.o $(LINKFVN) -o $@ | 18 | 19 | $(LINK) $(LINKFLAGS) $< init_random_seed.o $(LINKFVN) -o $@ | |
19 | 20 | |||
%.o: %.f90 | 20 | 21 | %.o: %.f90 | |
$(F95) $(F95FLAGS) -c $< -o $@ | 21 | 22 | $(F95) $(F95FLAGS) -c $< -o $@ | |
22 | 23 |
fvn_test/excursion.f90
File was created | 1 | module excursion | ||
2 | ||||
3 | real(8), dimension(:), pointer :: xm => NULL() | |||
4 | real(8), dimension(:), pointer :: ym => NULL() |
fvn_test/test_lm.f90
File was created | 1 | program gls | ||
2 | use fvn_misc | |||
3 | use excursion | |||
4 | implicit none | |||
5 | ||||
6 | interface | |||
7 | subroutine zef(m, n, x, fvec, iflag) | |||
8 | integer(4), intent(in) :: m,n | |||
9 | real(8), dimension(:), intent(in) :: x | |||
10 | real(8), dimension(:), intent(inout) :: fvec | |||
11 | integer(4), intent(inout) :: iflag | |||
12 | end subroutine | |||
13 | end interface | |||
14 | ||||
15 | ||||
16 | integer,parameter :: npoints=13,deg=3 | |||
17 | integer :: status,i | |||
18 | real(kind=dp_kind) :: xstep,xc,yc | |||
19 | real(kind=dp_kind) :: coeff(deg+1) | |||
20 | real(kind=dp_kind) :: fvec(npoints) | |||
21 | integer(4) :: iwa(deg+1) | |||
22 | integer(4) :: info | |||
23 | real(8) :: tol | |||
24 | ||||
25 | allocate(xm(npoints),ym(npoints)) | |||
26 | ||||
27 | xm = (/ -3.8,-2.7,-2.2,-1.9,-1.1,-0.7,0.5,1.7,2.,2.8,3.2,3.8,4. /) | |||
28 | ym = (/ -3.1,-2.,-0.9,0.8,1.8,0.4,2.1,1.8,3.2,2.8,3.9,5.2,7.5 /) | |||
29 | open(2,file='fvn_lsp_double_mesure.dat') | |||
30 | open(3,file='fvn_lm.dat') | |||
31 | do i=1,npoints | |||
32 | write(2,44) xm(i),ym(i) | |||
33 | end do | |||
34 | close(2) | |||
35 | ||||
36 | coeff=1. | |||
37 | call fvn_lm(zef,npoints,deg+1,coeff,info) | |||
38 | write(*,*) "info : ",info | |||
39 | ||||
40 | xstep=(xm(npoints)-xm(1))/1000. | |||
41 | do i=1,1000 | |||
42 | xc=xm(1)+(i-1)*xstep | |||
43 | yc=poly(xc,coeff) | |||
44 | write(3,44) xc,yc | |||
45 | end do | |||
46 | close(3) | |||
47 | write(*,*) "All done, plot results with gnuplot using command :" | |||
48 | write(*,*) "pl 'fvn_lsp_double_mesure.dat' u 1:2 w p,'fvn_lm.dat' u 1:2 w l" | |||
49 | 44 FORMAT(4(1X,1PE22.14)) | |||
50 | contains | |||
51 | function poly(x,coeff) | |||
52 | use fvn_common | |||
53 | implicit none | |||
54 | real(kind=dp_kind) :: x | |||
55 | real(kind=dp_kind) :: coeff(deg+1) | |||
56 | real(kind=dp_kind) :: poly | |||
57 | integer :: i | |||
58 | poly=0. | |||
59 | do i=1,deg+1 | |||
60 | poly=poly+coeff(i)*x**(i-1) | |||
61 | end do | |||
62 | end function | |||
63 | end program | |||
64 | ||||
65 | subroutine zef(m, n, x, fvec, iflag) | |||
66 | use excursion | |||
67 | implicit none | |||
68 | integer(4), intent(in) :: m,n |