1 | PROGRAM read_palm_netcdf_data |
---|
2 | |
---|
3 | !------------------------------------------------------------------------------! |
---|
4 | ! This is an example program for reading PALM 2d/3d NetCDF datasets |
---|
5 | ! |
---|
6 | ! The user has to add his own code for further storing and analyzing of |
---|
7 | ! these data!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
8 | ! |
---|
9 | ! The NetCDF include file and library has to be given with the respective |
---|
10 | ! compiler options. Please find out the respective paths on your system and |
---|
11 | ! set them appropriately. |
---|
12 | ! |
---|
13 | ! Here are some examples how this routine should be compiled: |
---|
14 | ! |
---|
15 | ! decalpha: |
---|
16 | ! f95 -fast -r8 -I/usr/local/netcdf-3.5.1/include |
---|
17 | ! -L/usr/local/netcdf-3.5.1/lib -lnetcdf |
---|
18 | ! IBM-Regatta: |
---|
19 | ! xlf95 -qrealsize=8 -q64 -qmaxmem=-1 -Q |
---|
20 | ! -I /aws/dataformats/netcdf-3.6.0-p1/64-32/include |
---|
21 | ! -L/aws/dataformats/netcdf-3.6.0-p1/64-32/lib -lnetcdf -O3 |
---|
22 | ! IBM-Regatta KISTI: |
---|
23 | ! xlf95 -qrealsize=8 -q64 -qmaxmem=-1 -Q |
---|
24 | ! -I /applic/netcdf64/src/f90 |
---|
25 | ! -L/applic/lib/NETCDF64 -lnetcdf -O3 |
---|
26 | ! IBM-Regatta Yonsei (gfdl5): |
---|
27 | ! xlf95 -qrealsize=8 -q64 -qmaxmem=-1 -Q |
---|
28 | ! -I /usr1/users/raasch/pub/netcdf-3.6.0-p1/include |
---|
29 | ! -L/usr1/users/raasch/pub/netcdf-3.6.0-p1/lib -lnetcdf -O3 |
---|
30 | ! IMUK: |
---|
31 | ! ifort read_palm...f90 -o read_palm...x |
---|
32 | ! -I /muksoft/packages/netcdf/linux/include -axW -r8 -nbs |
---|
33 | ! -Vaxlib -L /muksoft/packages/netcdf/linux/lib -lnetcdf |
---|
34 | ! NEC-SX6: |
---|
35 | ! sxf90 read_palm...f90 -o read_palm...x |
---|
36 | ! -I /pool/SX-6/netcdf/netcdf-3.6.0-p1/include -C hopt -Wf '-A idbl4' |
---|
37 | ! -L/pool/SX-6/netcdf/netcdf-3.6.0-p1/lib -lnetcdf |
---|
38 | !------------------------------------------------------------------------------! |
---|
39 | |
---|
40 | USE netcdf |
---|
41 | |
---|
42 | IMPLICIT NONE |
---|
43 | |
---|
44 | ! |
---|
45 | !-- Local variables |
---|
46 | CHARACTER (LEN=10) :: dimname(4), var_name |
---|
47 | CHARACTER (LEN=40) :: filename |
---|
48 | |
---|
49 | CHARACTER (LEN=2000) :: title, var_list |
---|
50 | |
---|
51 | INTEGER :: i, j, k, nc_stat, pos, time_step |
---|
52 | |
---|
53 | INTEGER :: current_level, current_var, id_set, id_var_time, num_var |
---|
54 | |
---|
55 | INTEGER, DIMENSION(4) :: id_dims, id_dims_loc, levels |
---|
56 | |
---|
57 | INTEGER, DIMENSION(1000) :: id_var |
---|
58 | |
---|
59 | REAL :: time(1) |
---|
60 | |
---|
61 | REAL, DIMENSION(:,:,:), ALLOCATABLE :: data_array |
---|
62 | |
---|
63 | |
---|
64 | PRINT*, '*** Please type NetCDF filename to be read:' |
---|
65 | READ*, filename |
---|
66 | |
---|
67 | nc_stat = NF90_OPEN( filename, NF90_WRITE, id_set ) |
---|
68 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 1 ) |
---|
69 | |
---|
70 | ! |
---|
71 | !-- Get the run description header and output |
---|
72 | title = ' ' |
---|
73 | nc_stat = NF90_GET_ATT( id_set, NF90_GLOBAL, 'title', title ) |
---|
74 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 2 ) |
---|
75 | WRITE (*,'(/A/A)') '*** file created by:', TRIM( title ) |
---|
76 | |
---|
77 | ! |
---|
78 | !-- Get the list of variables (order of variables corresponds with the |
---|
79 | !-- order of data on the binary file) |
---|
80 | var_list = ' ' ! GET_ATT does not assign trailing blanks |
---|
81 | nc_stat = NF90_GET_ATT( id_set, NF90_GLOBAL, 'VAR_LIST', var_list ) |
---|
82 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 3 ) |
---|
83 | |
---|
84 | ! |
---|
85 | !-- Inquire id of the time coordinate variable |
---|
86 | nc_stat = NF90_INQ_VARID( id_set, 'time', id_var_time ) |
---|
87 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 4 ) |
---|
88 | |
---|
89 | ! |
---|
90 | !-- Count number of variables; there is one more semicolon in the |
---|
91 | !-- string than variable names |
---|
92 | num_var = -1 |
---|
93 | DO i = 1, LEN( var_list ) |
---|
94 | IF ( var_list(i:i) == ';' ) num_var = num_var + 1 |
---|
95 | ENDDO |
---|
96 | WRITE (*,'(/A,I3,A/)') '*** file contains ', num_var, ' variable(s)' |
---|
97 | |
---|
98 | |
---|
99 | pos = INDEX( var_list, ';' ) |
---|
100 | ! |
---|
101 | !-- Loop over all variables |
---|
102 | DO i = 1, num_var |
---|
103 | |
---|
104 | ! |
---|
105 | !-- Extract variable name from list |
---|
106 | var_list = var_list(pos+1:) |
---|
107 | pos = INDEX( var_list, ';' ) |
---|
108 | var_name = var_list(1:pos-1) |
---|
109 | |
---|
110 | ! |
---|
111 | !-- Get variable ID from name |
---|
112 | nc_stat = NF90_INQ_VARID( id_set, TRIM( var_name ), id_var(i) ) |
---|
113 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 5 ) |
---|
114 | |
---|
115 | ! |
---|
116 | !-- Inquire the dimension IDs |
---|
117 | nc_stat = NF90_INQUIRE_VARIABLE( id_set, id_var(i), & |
---|
118 | dimids = id_dims_loc ) |
---|
119 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 6 ) |
---|
120 | id_dims = id_dims_loc |
---|
121 | |
---|
122 | ! |
---|
123 | !-- Get number of x/y/z/time levels(gridpoints) for that variable |
---|
124 | DO j = 1, 4 |
---|
125 | nc_stat = NF90_INQUIRE_DIMENSION( id_set, id_dims(j),& |
---|
126 | dimname(j), levels(j) ) |
---|
127 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 7 ) |
---|
128 | ENDDO |
---|
129 | |
---|
130 | WRITE (*,100) '*** reading variable "', TRIM(var_name), & |
---|
131 | '", dimensioned as', TRIM(var_name), levels(1)-1, & |
---|
132 | levels(2)-1, levels(3)-1 |
---|
133 | 100 FORMAT (A,A,A/4X,A,'(0:',I4,',0:',I4,',0:',I4,') (x/y/z)'/) |
---|
134 | |
---|
135 | ! |
---|
136 | !-- Allocate the data array to be read |
---|
137 | ALLOCATE( data_array(0:levels(1)-1,0:levels(2)-1,0:levels(3)-1) ) |
---|
138 | |
---|
139 | ! |
---|
140 | !-- Read the data from file for each timestep |
---|
141 | DO j = 1, levels(4) |
---|
142 | |
---|
143 | ! |
---|
144 | !-- Get the time of the current timelevel and output |
---|
145 | nc_stat = NF90_GET_VAR( id_set, id_var_time, time, start = (/ j /), & |
---|
146 | count = (/ 1 /) ) |
---|
147 | |
---|
148 | IF ( nc_stat /= NF90_NOERR ) CALL handle_netcdf_error( 7+i ) |
---|
149 | |
---|
150 | WRITE (*,'(A,I3,A,F8.1,A)') ' reading timelevel ', i, & |
---|
151 | ' time = ', time(1), ' s' |
---|
152 | |
---|
153 | nc_stat = NF90_GET_VAR( id_set, id_var(i), & |
---|
154 | data_array, start = (/ 1, 1, 1, j /), & |
---|
155 | count = (/ levels(1), levels(2), levels(3), 1 /) ) |
---|
156 | |
---|
157 | IF ( nc_stat /= NF90_NOERR ) & |
---|
158 | CALL handle_netcdf_error( 7+levels(4)+i ) |
---|
159 | ! |
---|
160 | !-- ADD YOUR OWN CODE FOR FURTHER STORING OF THESE DATA HERE |
---|
161 | !-- -------------------------------------------------------- |
---|
162 | |
---|
163 | |
---|
164 | ENDDO |
---|
165 | |
---|
166 | WRITE (*,'(/)') |
---|
167 | |
---|
168 | DEALLOCATE( data_array ) |
---|
169 | |
---|
170 | ENDDO |
---|
171 | |
---|
172 | |
---|
173 | |
---|
174 | CONTAINS |
---|
175 | |
---|
176 | |
---|
177 | SUBROUTINE handle_netcdf_error( errno ) |
---|
178 | ! |
---|
179 | !-- Prints out a text message corresponding to the current NetCDF status |
---|
180 | |
---|
181 | IMPLICIT NONE |
---|
182 | |
---|
183 | INTEGER, INTENT(IN) :: errno |
---|
184 | |
---|
185 | IF ( nc_stat /= NF90_NOERR ) THEN |
---|
186 | WRITE (*,'(A,1X,I3/4X,A)') & |
---|
187 | '+++ read_palm_netcdf_data error handle:', & |
---|
188 | errno, TRIM( nf90_strerror( nc_stat ) ) |
---|
189 | STOP |
---|
190 | ENDIF |
---|
191 | |
---|
192 | END SUBROUTINE handle_netcdf_error |
---|
193 | |
---|
194 | |
---|
195 | END PROGRAM read_palm_netcdf_data |
---|
196 | |
---|
197 | |
---|
198 | |
---|