source: LMDZ5/tags/proto-testing-20131015/tools/Max_diff_nc_with_lib/Max_diff_nc/max_diff_nc.f90 @ 1893

Last change on this file since 1893 was 1795, checked in by Ehouarn Millour, 11 years ago

Version testing basee sur la r1794


Testing release based on r1794

File size: 4.8 KB
Line 
1program max_diff_nc
2
3  ! This is a program in Fortran 95.
4  ! Author: Lionel GUEZ
5  ! See description in wrapper script.
6
7  ! Maximum memory used will normally be:
8
9  ! -- without computation of average order of magnitude: about 5
10  ! times the memory occupied by the largest variable;
11
12  ! -- with computation of average order of magnitude: about 7 times
13  ! the memory occupied by the largest variable.
14
15  ! This program is meant to be used with a wrapper script so input
16  ! statements do not have prompts.
17
18  use netcdf95, only: nf95_close, nf95_gw_var, nf95_inq_varid, nf95_inquire, &
19       nf95_inquire_variable, nf95_open
20  use netcdf, only: nf90_noerr, nf90_nowrite, nf90_max_name, NF90_FLOAT, &
21       NF90_double
22  use jumble, only: compare
23
24  implicit none
25
26  integer, parameter:: wp = kind(0d0) ! working precision
27
28  integer ncid1, ncid2, ncerr, xtype1, ndims
29  integer nvariables ! number of variables in the first file
30  integer nvar_comp ! number of variables which will be compared
31  integer, allocatable:: varid1(:), varid2(:)
32  real(wp), pointer:: v1_1d(:), v2_1d(:)
33  real(wp), pointer:: v1_2d(:, :), v2_2d(:, :)
34  real(wp), pointer:: v1_3d(:, :, :), v2_3d(:, :, :)
35  real(wp), pointer:: v1_4d(:, :, :, :), v2_4d(:, :, :, :)
36  character(len=nf90_max_name) name1
37  logical same_varid ! compare variables with same varid
38  logical report_id ! report identical variables
39  logical comp_mag ! compute avergage order of magnitude
40  logical quiet
41  character(len=30+nf90_max_name), allocatable:: tag(:)
42  integer i
43
44  !----------------------
45
46  read *, same_varid
47  read *, report_id
48  read *, quiet
49  read *, comp_mag
50  read "(a)", name1
51
52  call nf95_open("in1.nc", nf90_nowrite, ncid1)
53  call nf95_open("in2.nc", nf90_nowrite, ncid2)
54
55  ! Define "nvar_comp", "varid1(:nvar_comp)", "varid2(:nvar_comp)" and
56  ! "tag(:nvar_comp)":
57  if (name1 == "") then
58     ! We want to compare all the variables
59     call nf95_inquire(ncid1, nvariables=nvariables)
60     print *, "Found ", nvariables, " variable(s) in the first file."
61     allocate(varid1(nvariables), varid2(nvariables), tag(nvariables))
62     if (same_varid) then
63        nvar_comp = nvariables
64        varid1 = (/(i, i = 1, nvariables)/)
65        varid2 = varid1
66        do i = 1, nvariables
67           call nf95_inquire_variable(ncid1, varid1(i), name1)
68           tag(i) = 'Variable "' // trim(name1) // '" (name in the first file)'
69        end do
70     else
71        nvar_comp = 0
72        do i = 1, nvariables
73           call nf95_inquire_variable(ncid1, i, name1)
74           call nf95_inq_varid(ncid2, trim(name1), varid2(nvar_comp + 1), &
75                ncerr)
76           if (ncerr == nf90_noerr) then
77              varid1(nvar_comp + 1) = i
78              tag(nvar_comp + 1) = 'Variable "' // trim(name1) // '"'
79              nvar_comp = nvar_comp + 1
80           else
81              print *, 'Could not find "' // trim(name1) &
82                   // '" in the second file. Comparison will be skipped.'
83           end if
84        end do
85     end if
86  else
87     nvar_comp = 1
88     allocate(varid1(1), varid2(1), tag(1))
89     call nf95_inq_varid(ncid1, trim(name1), varid1(1))
90     call nf95_inq_varid(ncid2, trim(name1), varid2(1))
91     tag(1) = 'Variable "' // trim(name1) // '"'
92  end if
93
94  do i = 1, nvar_comp
95     call nf95_inquire_variable(ncid1, varid1(i), xtype=xtype1, ndims=ndims)
96     if (xtype1 == nf90_float .or. xtype1 == nf90_double) then
97        select case (ndims)
98        case (1)
99           call nf95_gw_var(ncid1, varid1(i), v1_1d)
100           call nf95_gw_var(ncid2, varid2(i), v2_1d)
101           call compare(v1_1d, v2_1d, trim(tag(i)), comp_mag, report_id, quiet)
102           deallocate(v1_1d, v2_1d)
103        case (2)
104           call nf95_gw_var(ncid1, varid1(i), v1_2d)
105           call nf95_gw_var(ncid2, varid2(i), v2_2d)
106           call compare(v1_2d, v2_2d, trim(tag(i)), comp_mag, report_id, quiet)
107           deallocate(v1_2d, v2_2d)
108        case (3)
109           call nf95_gw_var(ncid1, varid1(i), v1_3d)
110           call nf95_gw_var(ncid2, varid2(i), v2_3d)
111           call compare(v1_3d, v2_3d, trim(tag(i)), comp_mag, report_id, quiet)
112           deallocate(v1_3d, v2_3d)
113        case (4)
114           call nf95_gw_var(ncid1, varid1(i), v1_4d)
115           call nf95_gw_var(ncid2, varid2(i), v2_4d)
116           call compare(v1_4d, v2_4d, trim(tag(i)), comp_mag, report_id, quiet)
117           deallocate(v1_4d, v2_4d)
118        case default
119           print *
120           print *, "******************"
121           print *, trim(tag(i)) // ":"
122           print *, "Rank not supported."
123           print *, "ndims = ", ndims
124        end select
125     else
126        print *
127        print *, "******************"
128        print *, trim(tag(i)) // ":"
129        print *, 'Not of type "nf90_float or "nf90_double".'
130     end if
131  end do
132
133  call nf95_close(ncid1)
134  call nf95_close(ncid2)
135
136end program max_diff_nc
Note: See TracBrowser for help on using the repository browser.