1 | program concatnc |
---|
2 | |
---|
3 | |
---|
4 | ! ******************************************************** |
---|
5 | ! Program to concatenate data from Netcdf "diagfi" files, |
---|
6 | ! outputs from the martian GCM |
---|
7 | ! input : diagfi.nc / concat.nc / stats.nc kind of files |
---|
8 | ! author: Y. Wanherdrick |
---|
9 | ! + aps(), bps() and phisinit() are now also written |
---|
10 | ! to output file (E. Millour, september 2006) |
---|
11 | ! if available (F. Forget, october 2006) |
---|
12 | ! + handle 1D data (EM, January 2007) |
---|
13 | ! + ap(), bp() (F. Forget, February 2008) |
---|
14 | ! + handle the possibility that number of GCM layers (aps, bps |
---|
15 | ! or sigma) may be different from number of vertical levels |
---|
16 | ! of data (which is the case for outputs from 'zrecast') |
---|
17 | ! (EM, April 2010) |
---|
18 | ! + handle absence of ap() and bp() if aps and bps are available |
---|
19 | ! (case of stats file) FF, November 2011 |
---|
20 | ! + read and write controle field, if available. TN, October 2013 |
---|
21 | ! + possibility to concatenate concat files with a coherent |
---|
22 | ! time axis in the output. AB, April 2020 |
---|
23 | ! ******************************************************** |
---|
24 | |
---|
25 | implicit none |
---|
26 | |
---|
27 | include "netcdf.inc" ! NetCDF definitions |
---|
28 | |
---|
29 | character (len=80), dimension(1000) :: file |
---|
30 | ! file(): input file(s) names(s) |
---|
31 | character (len=30), dimension(16) :: notconcat |
---|
32 | ! notconcat(): names of the (16) variables that won't be concatenated |
---|
33 | character (len=50), dimension(:), allocatable :: var |
---|
34 | ! var(): name(s) of variable(s) that will be concatenated |
---|
35 | character (len=50) :: tmpvar,tmpfile,title,units |
---|
36 | ! tmpvar(): used to temporarily store a variable name |
---|
37 | ! tmpfile(): used to temporarily store a file name |
---|
38 | ! title(): [netcdf] title attribute |
---|
39 | ! units(): [netcdf] units attribute |
---|
40 | character (len=100) :: filename,vartmp |
---|
41 | ! filename(): output file name |
---|
42 | ! vartmp(): temporary variable name (read from netcdf input file) |
---|
43 | !character (len=1) :: ccopy |
---|
44 | ! ccpy: 'y' or 'n' answer |
---|
45 | character (len=4) :: axis |
---|
46 | ! axis: "ls" or "sol" |
---|
47 | integer :: nid,ierr,miss |
---|
48 | ! nid: [netcdf] file ID # |
---|
49 | ! ierr: [netcdf] subroutine returned error code |
---|
50 | ! miss: [netcdf] subroutine returned error code |
---|
51 | integer :: i,j,k,inter |
---|
52 | ! for various loops |
---|
53 | integer :: varid |
---|
54 | ! varid: [netcdf] variable ID # |
---|
55 | integer :: memolatlen=0,memolonlen=0,memoaltlen=0,memoctllen=0 |
---|
56 | ! memolatlen: # of elements of lat(), read from the first input file |
---|
57 | ! memolonlen: # of elements of lon(), read from the first input file |
---|
58 | ! memoaltlen: # of elements of alt(), read from the first input file |
---|
59 | ! memoctllen: # of elements of ctl(), read from the first input file |
---|
60 | real, dimension(:), allocatable:: lat,lon,alt,ctl,time |
---|
61 | ! lat(): array, stores latitude coordinates |
---|
62 | ! lon(): array, stores longitude coordinates |
---|
63 | ! alt(): array, stores altitude coordinates |
---|
64 | ! ctl(): array, stores controle coordinates |
---|
65 | ! time(): array, stores time coordinates |
---|
66 | integer :: nbvar,nbfile,nbvarfile,ndim |
---|
67 | ! nbvar: # of variables to concatenate |
---|
68 | ! nbfile: # number of input file(s) |
---|
69 | ! nbvarfile: total # of variables in an input file |
---|
70 | ! ndim: [netcdf] # (3 or 4) of dimensions (for variables) |
---|
71 | integer :: latdim,londim,altdim,ctldim,timedim |
---|
72 | ! latdim: [netcdf] "latitude" dim ID |
---|
73 | ! londim: [netcdf] "longitude" dim ID |
---|
74 | ! altdim: [netcdf] "altdim" dim ID |
---|
75 | ! ctldim: [netcdf] "ctldim" dim ID |
---|
76 | ! timedim: [netcdf] "timedim" dim ID |
---|
77 | integer :: gcmlayerdim ! NetCDF dimension ID for # of layers in GCM |
---|
78 | integer :: latvar,lonvar,altvar,ctlvar,timevar |
---|
79 | ! latvar: [netcdf] ID of "latitude" variable |
---|
80 | ! lonvar: [netcdf] ID of "longitude" variable |
---|
81 | ! altvar: [netcdf] ID of "altitude" variable |
---|
82 | ! ctlvar: [netcdf] ID of "controle" variable |
---|
83 | ! timevar: [netcdf] ID of "Time" variable |
---|
84 | integer :: latlen,lonlen,altlen,ctllen,timelen |
---|
85 | ! latlen: # of elements of lat() array |
---|
86 | ! lonlen: # of elements of lon() array |
---|
87 | ! altlen: # of elements of alt() array |
---|
88 | ! ctllen: # of elements of ctl() array |
---|
89 | ! timelen: # of elemnets of time() array |
---|
90 | integer :: GCM_layers ! number of GCM atmospheric layers (may not be |
---|
91 | ! same as altlen if file is an output of zrecast) |
---|
92 | integer :: nout,latdimout,londimout,altdimout,timedimout,timevarout |
---|
93 | ! nout: [netcdf] output file ID |
---|
94 | ! latdimout: [netcdf] output latitude (dimension) ID |
---|
95 | ! londimout: [netcdf] output longitude (dimension) ID |
---|
96 | ! altdimout: [netcdf] output altitude (dimension) ID |
---|
97 | ! timedimout: [netcdf] output time (dimension) ID |
---|
98 | ! timevarout: [netcdf] ID of output "Time" variable |
---|
99 | integer :: layerdimout ! NetCDF dimension ID for # of layers in GCM |
---|
100 | integer :: interlayerdimout ! dimension ID for # of interlayers in GCM |
---|
101 | integer :: reptime,rep,varidout |
---|
102 | ! reptime: total length of concatenated time() arrays |
---|
103 | ! rep: # number of elements of a time() array to write to the output file |
---|
104 | ! varidout: [netcdf] variable ID # (of a variable to write to the output file) |
---|
105 | integer :: Nnotconcat,var_ok |
---|
106 | ! Nnotconcat: # of (leading)variables that won't be concatenated |
---|
107 | ! var_ok: flag (0 or 1) |
---|
108 | integer, dimension(4) :: corner,edges,dim |
---|
109 | ! corner: [netcdf] |
---|
110 | ! edges: [netcdf] |
---|
111 | ! dim: [netcdf] |
---|
112 | real, dimension(:,:,:,:), allocatable :: var3d |
---|
113 | ! var3D(,,,): 4D array to store a field |
---|
114 | real :: ctlsol |
---|
115 | ! ctlsol: starting sol value of the file, stored in the variable ctl() |
---|
116 | real :: memotime |
---|
117 | ! memotime: (cumulative) time value, in martian days (sols) |
---|
118 | real :: missing |
---|
119 | !PARAMETER(missing=1E+20) |
---|
120 | ! missing: [netcdf] to handle "missing" values when reading/writing files |
---|
121 | real, dimension(2) :: valid_range |
---|
122 | ! valid_range(2): [netcdf] interval in which a value is considered valid |
---|
123 | |
---|
124 | !============================================================================== |
---|
125 | ! 1.1. Get input file name(s) |
---|
126 | !============================================================================== |
---|
127 | write(*,*) |
---|
128 | write(*,*) "-- Program written specifically for NetCDF 'diagfi'/'concat' files from the martian GCM --" |
---|
129 | write(*,*) |
---|
130 | write(*,*) "which files do you want to use?" |
---|
131 | write(*,*) "<Enter> when list ok" |
---|
132 | |
---|
133 | nbfile=0 |
---|
134 | read(*,'(a50)') tmpfile |
---|
135 | do While (len_trim(tmpfile).ne.0) |
---|
136 | nbfile=nbfile+1 |
---|
137 | file(nbfile)=tmpfile |
---|
138 | read(*,'(a50)') tmpfile |
---|
139 | enddo |
---|
140 | |
---|
141 | if(nbfile==0) then |
---|
142 | write(*,*) "no file... game over" |
---|
143 | stop "" |
---|
144 | endif |
---|
145 | |
---|
146 | !============================================================================== |
---|
147 | ! 1.2. Ask for starting day value (memotime) |
---|
148 | !============================================================================== |
---|
149 | |
---|
150 | write(*,*) |
---|
151 | !write(*,*) "Beginning day of the first specified file?" |
---|
152 | write(*,*) "Starting day of the run stored in the first input file?" |
---|
153 | write(*,*) " (Obsolete if the controle field is present, answer any number)" |
---|
154 | write(*,*) "(e.g.: 100 if that run started at time=100 sols)" |
---|
155 | read(*,*) memotime |
---|
156 | |
---|
157 | !============================================================================== |
---|
158 | ! 1.3. Open the first input file |
---|
159 | !============================================================================== |
---|
160 | |
---|
161 | write(*,*) |
---|
162 | write(*,*) "opening "//trim(file(1))//"..." |
---|
163 | ierr = NF_OPEN(file(1),NF_NOWRITE,nid) |
---|
164 | if (ierr.NE.NF_NOERR) then |
---|
165 | write(*,*) 'ERROR: Pb opening file '//file(1) |
---|
166 | stop "" |
---|
167 | endif |
---|
168 | |
---|
169 | ierr=NF_INQ_NVARS(nid,nbvarfile) |
---|
170 | ! nbvarfile now set to be the (total) number of variables in file |
---|
171 | |
---|
172 | !============================================================================== |
---|
173 | ! 1.4. Ask for (output) "Time" axis type |
---|
174 | !============================================================================== |
---|
175 | |
---|
176 | write(*,*) "Warning: to read the result with grads, choose sol" |
---|
177 | write(*,*) "Warning: use ferret to read the non linear scale ls" |
---|
178 | write(*,*) "Which time axis should be given in the output? (sol/ls)" |
---|
179 | read(*,*) axis |
---|
180 | ! loop as long as axis is neither "sol" nor "ls" |
---|
181 | do while ((axis/="sol").AND.(axis/="ls")) |
---|
182 | read(*,*) axis |
---|
183 | enddo |
---|
184 | |
---|
185 | ! The following variables don't need to be concatenated |
---|
186 | notconcat(1)='Time' |
---|
187 | notconcat(2)='controle' |
---|
188 | notconcat(3)='rlonu' |
---|
189 | notconcat(4)='latitude' |
---|
190 | notconcat(5)='longitude' |
---|
191 | notconcat(6)='altitude' |
---|
192 | notconcat(7)='rlatv' |
---|
193 | notconcat(8)='aps' |
---|
194 | notconcat(9)='bps' |
---|
195 | notconcat(10)='ap' |
---|
196 | notconcat(11)='bp' |
---|
197 | notconcat(12)='cu' |
---|
198 | notconcat(13)='cv' |
---|
199 | notconcat(14)='aire' |
---|
200 | notconcat(15)='phisinit' |
---|
201 | notconcat(16)='soildepth' |
---|
202 | |
---|
203 | |
---|
204 | !============================================================================== |
---|
205 | ! 1.5. Get (and check) list of variables to concatenate |
---|
206 | !============================================================================== |
---|
207 | write(*,*) |
---|
208 | Nnotconcat=0 |
---|
209 | do i=1,nbvarfile |
---|
210 | ierr=NF_INQ_VARNAME(nid,i,vartmp) |
---|
211 | ! vartmp now contains the "name" of variable of ID # i |
---|
212 | var_ok=0 |
---|
213 | do inter=1,size(notconcat) |
---|
214 | if (vartmp.eq.notconcat(inter)) then |
---|
215 | var_ok=1 |
---|
216 | Nnotconcat=Nnotconcat+1 |
---|
217 | endif |
---|
218 | enddo |
---|
219 | if (var_ok.eq.0) write(*,*) trim(vartmp) |
---|
220 | enddo |
---|
221 | |
---|
222 | ! Nnotconcat: # of variables that won't be concatenated |
---|
223 | ! nbvarfile: total # of variables in file |
---|
224 | allocate(var(nbvarfile-Nnotconcat)) |
---|
225 | |
---|
226 | |
---|
227 | write(*,*) |
---|
228 | write(*,*) "which variables do you want to concatenate?" |
---|
229 | write(*,*) "all / list of <variables> (separated by <Enter>s)" |
---|
230 | write(*,*) "(an empty line , i.e: just <Enter>, implies end of list)" |
---|
231 | nbvar=0 |
---|
232 | read(*,'(a50)') tmpvar |
---|
233 | do while ((tmpvar/=' ').AND.(trim(tmpvar)/='all')) |
---|
234 | nbvar=nbvar+1 |
---|
235 | var(nbvar)=tmpvar |
---|
236 | read(*,'(a50)') tmpvar |
---|
237 | enddo |
---|
238 | |
---|
239 | if (tmpvar=="all") then |
---|
240 | if (axis=="ls") then |
---|
241 | ! write(*,*) "Do you want to keep the original file? (y/n)" |
---|
242 | ! read(*,*) ccopy |
---|
243 | ! if ((ccopy=="n").or.(ccopy=="N")) then |
---|
244 | ! do i=1,nbfile |
---|
245 | ! ierr=NF_CLOSE(nid) |
---|
246 | ! ierr = NF_OPEN(file(1),NF_WRITE,nid) |
---|
247 | ! call change_time_axis(nid,ierr) |
---|
248 | ! ierr=NF_CLOSE(nid) |
---|
249 | ! STOP "" |
---|
250 | ! enddo |
---|
251 | ! else |
---|
252 | nbvar=nbvarfile-Nnotconcat |
---|
253 | do j=Nnotconcat+1,nbvarfile |
---|
254 | ierr=nf_inq_varname(nid,j,var(j-Nnotconcat)) |
---|
255 | enddo |
---|
256 | ! endif ! of if ((ccopy=="n").or.(ccopy=="N")) |
---|
257 | endif ! of if (axis=="ls") |
---|
258 | ! Variables names from the file are catched |
---|
259 | nbvar=nbvarfile-Nnotconcat |
---|
260 | j=1 |
---|
261 | do i=1,nbvarfile |
---|
262 | ierr=nf_inq_varname(nid,i,vartmp) |
---|
263 | var_ok=0 |
---|
264 | do inter=1,size(notconcat) |
---|
265 | if (vartmp.eq.notconcat(inter)) then |
---|
266 | var_ok=1 |
---|
267 | endif |
---|
268 | enddo |
---|
269 | if (var_ok.eq.0) then |
---|
270 | if (j .gt. nbvar) then |
---|
271 | write(*,*) "PROBLEM HERE !", var |
---|
272 | stop |
---|
273 | endif |
---|
274 | var(j) = vartmp |
---|
275 | write(*,'(a9,1x,i2,1x,a1,1x,a50)') "variable ",j,":",var(j) |
---|
276 | j=j+1 |
---|
277 | endif |
---|
278 | enddo |
---|
279 | else if(nbvar==0) then |
---|
280 | write(*,*) "no variable... game over" |
---|
281 | stop "" |
---|
282 | endif ! of if (tmpvar=="all") |
---|
283 | |
---|
284 | ! Name of the new file |
---|
285 | !========================================================== |
---|
286 | !filename=var(1) |
---|
287 | !do i=2, nbvar |
---|
288 | ! filename=trim(adjustl(filename))//"_"//var(i) |
---|
289 | !enddo |
---|
290 | !filename=trim(adjustl(filename))//".nc" |
---|
291 | |
---|
292 | !============================================================================== |
---|
293 | ! 1.6. Get output file name |
---|
294 | !============================================================================== |
---|
295 | filename="concat.nc" |
---|
296 | |
---|
297 | |
---|
298 | !============================================================================== |
---|
299 | ! 2. Concatenate input file(s) into output file |
---|
300 | !============================================================================== |
---|
301 | |
---|
302 | reptime=0 |
---|
303 | |
---|
304 | do i=1,nbfile |
---|
305 | |
---|
306 | !============================================================================== |
---|
307 | ! 2.1. Open input file |
---|
308 | !============================================================================== |
---|
309 | |
---|
310 | if (i/=1) then |
---|
311 | write(*,*) |
---|
312 | write(*,*) "opening "//trim(file(i))//"..." |
---|
313 | ierr = NF_OPEN(file(i),NF_NOWRITE,nid) |
---|
314 | if (ierr.NE.NF_NOERR) then |
---|
315 | write(*,*) 'ERROR: Pb opening file '//file(i) |
---|
316 | write(*,*) NF_STRERROR(ierr) |
---|
317 | stop "" |
---|
318 | endif |
---|
319 | endif |
---|
320 | |
---|
321 | !============================================================================== |
---|
322 | ! 2.2. Read (and check) dimensions of variables from input file |
---|
323 | !============================================================================== |
---|
324 | |
---|
325 | ierr=NF_INQ_DIMID(nid,"latitude",latdim) |
---|
326 | ierr=NF_INQ_VARID(nid,"latitude",latvar) |
---|
327 | if (ierr.NE.NF_NOERR) then |
---|
328 | write(*,*) 'ERROR: Field <latitude> is missing in file '//file(i) |
---|
329 | stop "" |
---|
330 | endif |
---|
331 | ierr=NF_INQ_DIMLEN(nid,latdim,latlen) |
---|
332 | ! write(*,*) "latlen: ",latlen |
---|
333 | |
---|
334 | ierr=NF_INQ_DIMID(nid,"longitude",londim) |
---|
335 | ierr=NF_INQ_VARID(nid,"longitude",lonvar) |
---|
336 | if (ierr.NE.NF_NOERR) then |
---|
337 | write(*,*) 'ERROR: Field <longitude> is missing in file '//file(i) |
---|
338 | stop "" |
---|
339 | endif |
---|
340 | ierr=NF_INQ_DIMLEN(nid,londim,lonlen) |
---|
341 | ! write(*,*) "lonlen: ",lonlen |
---|
342 | |
---|
343 | ierr=NF_INQ_DIMID(nid,"altitude",altdim) |
---|
344 | ierr=NF_INQ_VARID(nid,"altitude",altvar) |
---|
345 | if (ierr.NE.NF_NOERR) then |
---|
346 | write(*,*) 'ERROR: Field <altitude> is missing in file '//file(i) |
---|
347 | stop "" |
---|
348 | endif |
---|
349 | ierr=NF_INQ_DIMLEN(nid,altdim,altlen) |
---|
350 | ! write(*,*) "altlen: ",altlen |
---|
351 | |
---|
352 | ierr=NF_INQ_DIMID(nid,"index",ctldim) |
---|
353 | ierr=NF_INQ_VARID(nid,"controle",ctlvar) |
---|
354 | if (ierr.NE.NF_NOERR) then |
---|
355 | write(*,*) 'Field <controle> is missing in file '//file(i) |
---|
356 | ctllen=0 |
---|
357 | !stop "" |
---|
358 | else |
---|
359 | ierr=NF_INQ_DIMLEN(nid,ctldim,ctllen) |
---|
360 | endif |
---|
361 | ! write(*,*) "controle: ",controle |
---|
362 | |
---|
363 | ! load size of aps() or sigma() (in case it is not altlen) |
---|
364 | ! default is that GCM_layers=altlen |
---|
365 | ! but for outputs of zrecast, it may be a different value |
---|
366 | ierr=NF_INQ_DIMID(nid,"GCM_layers",gcmlayerdim) |
---|
367 | if (ierr.ne.NF_NOERR) then |
---|
368 | ! didn't find a GCM_layers dimension; therefore we have: |
---|
369 | GCM_layers=altlen |
---|
370 | else |
---|
371 | ! load value of GCM_layers |
---|
372 | ierr=NF_INQ_DIMLEN(nid,gcmlayerdim,GCM_layers) |
---|
373 | endif |
---|
374 | ! write(*,*) "GCM_layers=",GCM_layers |
---|
375 | |
---|
376 | !============================================================================== |
---|
377 | ! 2.3. Read (and check compatibility of) dimensions of |
---|
378 | ! variables from input file |
---|
379 | !============================================================================== |
---|
380 | |
---|
381 | if (ctllen .ne. 0) then ! variable controle |
---|
382 | allocate(ctl(ctllen)) |
---|
383 | #ifdef NC_DOUBLE |
---|
384 | ierr = NF_GET_VAR_DOUBLE(nid,ctlvar,ctl) |
---|
385 | #else |
---|
386 | ierr = NF_GET_VAR_REAL(nid,ctlvar,ctl) |
---|
387 | #endif |
---|
388 | ctlsol = ctl(4) |
---|
389 | endif |
---|
390 | |
---|
391 | if (i==1) then ! First call; initialize/allocate |
---|
392 | memolatlen=latlen |
---|
393 | memolonlen=lonlen |
---|
394 | memoaltlen=altlen |
---|
395 | memoctllen=ctllen |
---|
396 | allocate(lat(latlen)) |
---|
397 | allocate(lon(lonlen)) |
---|
398 | allocate(alt(altlen)) |
---|
399 | #ifdef NC_DOUBLE |
---|
400 | ierr = NF_GET_VAR_DOUBLE(nid,latvar,lat) |
---|
401 | ierr = NF_GET_VAR_DOUBLE(nid,lonvar,lon) |
---|
402 | ierr = NF_GET_VAR_DOUBLE(nid,altvar,alt) |
---|
403 | #else |
---|
404 | ierr = NF_GET_VAR_REAL(nid,latvar,lat) |
---|
405 | ierr = NF_GET_VAR_REAL(nid,lonvar,lon) |
---|
406 | ierr = NF_GET_VAR_REAL(nid,altvar,alt) |
---|
407 | #endif |
---|
408 | if (ctllen .ne. 0) then |
---|
409 | if (modulo(int(memotime),669)/=modulo(int(ctlsol),669)) then |
---|
410 | write(*,*) "WARNING: Starting day of the run is not ",& |
---|
411 | modulo(int(memotime),669)," but ",modulo(int(ctlsol),669),"!!" |
---|
412 | write(*,*) "Starting day of the run has been corrected." |
---|
413 | ! memotime=float(modulo(int(ctl(4)),669)) + ctl(27) |
---|
414 | memotime=float(int(ctlsol)) + ctl(27) |
---|
415 | ctl(4) = 0. ! values written in the output |
---|
416 | ctl(27) = 0. ! file controle variable (resp. day_ini and hour_ini) |
---|
417 | endif |
---|
418 | endif |
---|
419 | ! Initialize output file's lat,lon,alt and time dimensions |
---|
420 | call initiate (filename,lat,lon,alt,ctl,GCM_layers,nout,& |
---|
421 | latdimout,londimout,altdimout,timedimout,& |
---|
422 | layerdimout,interlayerdimout,timevarout,axis) |
---|
423 | ! Initialize output file's aps,bps,ap,bp and phisinit variables |
---|
424 | call init2(nid,lonlen,latlen,altlen,GCM_layers,& |
---|
425 | nout,londimout,latdimout,altdimout,& |
---|
426 | layerdimout,interlayerdimout) |
---|
427 | |
---|
428 | else ! Not a first call, |
---|
429 | ! Check that latitude,longitude and altitude of current input file |
---|
430 | ! are identical to those of the output file |
---|
431 | if (memolatlen/=latlen) then |
---|
432 | write(*,*) "ERROR: Not the same latitude axis" |
---|
433 | stop "" |
---|
434 | else if (memolonlen/=lonlen) then |
---|
435 | write(*,*) "ERROR: Not the same longitude axis" |
---|
436 | stop "" |
---|
437 | else if (memoaltlen/=altlen) then |
---|
438 | write(*,*) "ERROR: Not the same altitude axis" |
---|
439 | stop "" |
---|
440 | else if (memoctllen/=ctllen) then |
---|
441 | write(*,*) "ERROR: Not the same controle axis" |
---|
442 | stop "" |
---|
443 | endif |
---|
444 | endif ! of if (i==1) |
---|
445 | |
---|
446 | !============================================================================== |
---|
447 | ! 2.4. Handle "Time" dimension from input file |
---|
448 | !============================================================================== |
---|
449 | |
---|
450 | !============================================================================== |
---|
451 | ! 2.4.1 Read "Time" dimension from input file |
---|
452 | !============================================================================== |
---|
453 | ierr=NF_INQ_DIMID(nid,"Time",timedim) |
---|
454 | if (ierr.NE.NF_NOERR) then |
---|
455 | write(*,*) 'ERROR: Dimension <Time> is missing in file '//file(i) |
---|
456 | stop "" |
---|
457 | endif |
---|
458 | ierr=NF_INQ_VARID(nid,"Time",timevar) |
---|
459 | if (ierr.NE.NF_NOERR) then |
---|
460 | write(*,*) 'ERROR: Field <Time> is missing in file '//file(i) |
---|
461 | stop "" |
---|
462 | endif |
---|
463 | ierr=NF_INQ_DIMLEN(nid,timedim,timelen) |
---|
464 | ! write(*,*) "timelen: ",timelen |
---|
465 | |
---|
466 | ! allocate time() array and fill it with values from input file |
---|
467 | allocate(time(timelen)) |
---|
468 | |
---|
469 | #ifdef NC_DOUBLE |
---|
470 | ierr = NF_GET_VAR_DOUBLE(nid,timevar,time) |
---|
471 | #else |
---|
472 | ierr = NF_GET_VAR_REAL(nid,timevar,time) |
---|
473 | #endif |
---|
474 | |
---|
475 | !============================================================================== |
---|
476 | ! 2.4.2 Write/extend "Time" dimension/values in output file |
---|
477 | !============================================================================== |
---|
478 | |
---|
479 | rep=0 |
---|
480 | write(*,*) |
---|
481 | write(*,'(a3,1x,f6.1)') "Sol",memotime |
---|
482 | ! if (ctllen.ne.0) write(*,'(a30,1x,f6.1)') "In the current file's ctl: Sol",ctlsol |
---|
483 | |
---|
484 | ! Add (memotime)/(ctlsol) offset and write "concatenated" time values to output file |
---|
485 | if (ctllen.ne.0) then |
---|
486 | do k=reptime+1,reptime+timelen |
---|
487 | rep=rep+1 |
---|
488 | if ((time(rep)+ctlsol).le.memotime)then |
---|
489 | write(*,*) "ERROR : the time intervals of the files are not dissociated" |
---|
490 | stop "" |
---|
491 | endif |
---|
492 | #ifdef NC_DOUBLE |
---|
493 | ierr= NF_PUT_VARA_DOUBLE(nout,timevarout,k,1,ctlsol+time(rep)) |
---|
494 | #else |
---|
495 | ierr= NF_PUT_VARA_REAL(nout,timevarout,k,1,ctlsol+time(rep)) |
---|
496 | #endif |
---|
497 | enddo |
---|
498 | |
---|
499 | else |
---|
500 | do k=reptime+1,reptime+timelen |
---|
501 | rep=rep+1 |
---|
502 | #ifdef NC_DOUBLE |
---|
503 | ierr= NF_PUT_VARA_DOUBLE(nout,timevarout,k,1,memotime+time(rep)) |
---|
504 | #else |
---|
505 | ierr= NF_PUT_VARA_REAL(nout,timevarout,k,1,memotime+time(rep)) |
---|
506 | #endif |
---|
507 | enddo |
---|
508 | endif ! if ctllen.ne.0 |
---|
509 | |
---|
510 | ! Compute new time offset (for further concatenations) |
---|
511 | memotime=memotime+time(timelen) |
---|
512 | |
---|
513 | !============================================================================== |
---|
514 | ! 2.5 Read/write variables |
---|
515 | !============================================================================== |
---|
516 | |
---|
517 | do j=1,nbvar ! loop on variables to read/write |
---|
518 | |
---|
519 | !============================================================================== |
---|
520 | ! 2.5.1 Check that variable to be read existe in input file |
---|
521 | !============================================================================== |
---|
522 | |
---|
523 | write(*,*) "variable ",var(j) |
---|
524 | ierr=nf_inq_varid(nid,var(j),varid) |
---|
525 | if (ierr.NE.NF_NOERR) then |
---|
526 | write(*,*) 'ERROR: Field <',var(j),'> not found in file'//file(i) |
---|
527 | stop "" |
---|
528 | endif |
---|
529 | ierr=nf_inq_varndims(nid,varid,ndim) |
---|
530 | |
---|
531 | !============================================================================== |
---|
532 | ! 2.5.2 Prepare things in order to read/write the variable |
---|
533 | !============================================================================== |
---|
534 | |
---|
535 | ! build dim(),corner() and edges() arrays |
---|
536 | ! and allocate var3d() array |
---|
537 | if (ndim==1) then |
---|
538 | allocate(var3d(timelen,1,1,1)) |
---|
539 | dim(1)=timedimout |
---|
540 | |
---|
541 | ! start indexes (where data values will be written) |
---|
542 | corner(1)=reptime+1 |
---|
543 | corner(2)=1 |
---|
544 | corner(3)=1 |
---|
545 | corner(4)=1 |
---|
546 | |
---|
547 | ! length (along dimensions) of block of data to be written |
---|
548 | edges(1)=timelen |
---|
549 | edges(2)=1 |
---|
550 | edges(3)=1 |
---|
551 | edges(4)=1 |
---|
552 | |
---|
553 | else if (ndim==3) then |
---|
554 | allocate(var3d(lonlen,latlen,timelen,1)) |
---|
555 | dim(1)=londimout |
---|
556 | dim(2)=latdimout |
---|
557 | dim(3)=timedimout |
---|
558 | |
---|
559 | ! start indexes (where data values will be written) |
---|
560 | corner(1)=1 |
---|
561 | corner(2)=1 |
---|
562 | corner(3)=reptime+1 |
---|
563 | corner(4)=1 |
---|
564 | |
---|
565 | ! length (along dimensions) of block of data to be written |
---|
566 | edges(1)=lonlen |
---|
567 | edges(2)=latlen |
---|
568 | edges(3)=timelen |
---|
569 | edges(4)=1 |
---|
570 | |
---|
571 | else if (ndim==4) then |
---|
572 | allocate(var3d(lonlen,latlen,altlen,timelen)) |
---|
573 | dim(1)=londimout |
---|
574 | dim(2)=latdimout |
---|
575 | dim(3)=altdimout |
---|
576 | dim(4)=timedimout |
---|
577 | |
---|
578 | ! start indexes (where data values will be written) |
---|
579 | corner(1)=1 |
---|
580 | corner(2)=1 |
---|
581 | corner(3)=1 |
---|
582 | corner(4)=reptime+1 |
---|
583 | |
---|
584 | ! length (along dimensions) of block of data to be written |
---|
585 | edges(1)=lonlen |
---|
586 | edges(2)=latlen |
---|
587 | edges(3)=altlen |
---|
588 | edges(4)=timelen |
---|
589 | endif |
---|
590 | |
---|
591 | if (i==1) then ! First call: write some definitions to output file |
---|
592 | units=" " |
---|
593 | title=" " |
---|
594 | ierr=nf_get_att_text(nid,varid,"title",title) |
---|
595 | ierr=nf_get_att_text(nid,varid,"units",units) |
---|
596 | call def_var(nout,var(j),title,units,ndim,dim,varidout,ierr) |
---|
597 | else |
---|
598 | ierr=NF_INQ_VARID(nout,var(j),varidout) |
---|
599 | endif |
---|
600 | |
---|
601 | !============================================================================== |
---|
602 | ! 2.5.3 Read from input file and write (append) to the output file |
---|
603 | !============================================================================== |
---|
604 | |
---|
605 | #ifdef NC_DOUBLE |
---|
606 | ierr = NF_GET_VAR_DOUBLE(nid,varid,var3d) |
---|
607 | ierr= NF_PUT_VARA_DOUBLE(nout,varidout,corner,edges,var3d) |
---|
608 | miss=NF_GET_ATT_DOUBLE(nid,varid,"missing_value",missing) |
---|
609 | miss=NF_GET_ATT_DOUBLE(nid,varid,"valid_range",valid_range) |
---|
610 | #else |
---|
611 | ierr = NF_GET_VAR_REAL(nid,varid,var3d) |
---|
612 | ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3d) |
---|
613 | miss=NF_GET_ATT_REAL(nid,varid,"missing_value",missing) |
---|
614 | miss=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range) |
---|
615 | #endif |
---|
616 | |
---|
617 | if (ierr.ne.NF_NOERR) then |
---|
618 | write(*,*) 'PUT_VAR ERROR: ',NF_STRERROR(ierr) |
---|
619 | stop "" |
---|
620 | endif |
---|
621 | |
---|
622 | ! In case there is a "valid_range" attribute |
---|
623 | ! Write "valid_range" and "missing_value" attributes in output file |
---|
624 | if (miss.eq.NF_NOERR) then |
---|
625 | call missing_value(nout,varidout,valid_range,missing) |
---|
626 | endif |
---|
627 | |
---|
628 | ! free var3d() array |
---|
629 | deallocate(var3d) |
---|
630 | |
---|
631 | enddo ! of do j=1,nbvar |
---|
632 | |
---|
633 | ! Free time() and compute/store array length (for further concatenations) |
---|
634 | deallocate(time) |
---|
635 | reptime=reptime+timelen |
---|
636 | |
---|
637 | ! Deallocate controle if needed |
---|
638 | if (ctllen.ne.0) deallocate(ctl) |
---|
639 | |
---|
640 | ! Close input file |
---|
641 | ierr=nf_close(nid) |
---|
642 | |
---|
643 | enddo ! of i=1,nbfile |
---|
644 | |
---|
645 | !============================================================================== |
---|
646 | ! 3. If required, change time axis (from sols to Ls) |
---|
647 | !============================================================================== |
---|
648 | |
---|
649 | if (axis=="ls") then |
---|
650 | call change_time_axis(nout,ierr) |
---|
651 | endif |
---|
652 | |
---|
653 | ! Close output file |
---|
654 | ierr=NF_CLOSE(nout) |
---|
655 | |
---|
656 | contains |
---|
657 | |
---|
658 | !****************************************************************************** |
---|
659 | subroutine initiate (filename,lat,lon,alt,ctl,GCM_layers,nout,& |
---|
660 | latdimout,londimout,altdimout,timedimout,& |
---|
661 | layerdimout,interlayerdimout,timevarout,axis) |
---|
662 | !============================================================================== |
---|
663 | ! Purpose: |
---|
664 | ! Create and initialize a data file (NetCDF format) |
---|
665 | !============================================================================== |
---|
666 | ! Remarks: |
---|
667 | ! The NetCDF file (created in this subroutine) remains open |
---|
668 | !============================================================================== |
---|
669 | |
---|
670 | implicit none |
---|
671 | |
---|
672 | include "netcdf.inc" ! NetCDF definitions |
---|
673 | |
---|
674 | !============================================================================== |
---|
675 | ! Arguments: |
---|
676 | !============================================================================== |
---|
677 | character (len=*), intent(in):: filename |
---|
678 | ! filename(): the file's name |
---|
679 | real, dimension(:), intent(in):: lat |
---|
680 | ! lat(): latitude |
---|
681 | real, dimension(:), intent(in):: lon |
---|
682 | ! lon(): longitude |
---|
683 | real, dimension(:), intent(in):: alt |
---|
684 | ! alt(): altitude |
---|
685 | real, dimension(:), intent(in):: ctl |
---|
686 | ! ctl(): controle |
---|
687 | integer,intent(in) :: GCM_layers ! number of GCM layers |
---|
688 | integer, intent(out):: nout |
---|
689 | ! nout: [netcdf] file ID |
---|
690 | integer, intent(out):: latdimout |
---|
691 | ! latdimout: [netcdf] lat() (i.e.: latitude) ID |
---|
692 | integer, intent(out):: londimout |
---|
693 | ! londimout: [netcdf] lon() ID |
---|
694 | integer, intent(out):: altdimout |
---|
695 | ! altdimout: [netcdf] alt() ID |
---|
696 | integer, intent(out):: timedimout |
---|
697 | ! timedimout: [netcdf] "Time" ID |
---|
698 | integer,intent(out) :: layerdimout |
---|
699 | ! layerdimout: [netcdf] "GCM_layers" ID |
---|
700 | integer,intent(out) :: interlayerdimout |
---|
701 | ! layerdimout: [netcdf] "GCM_layers+1" ID |
---|
702 | integer, intent(out):: timevarout |
---|
703 | ! timevarout: [netcdf] "Time" (considered as a variable) ID |
---|
704 | character (len=4),intent(in) :: axis |
---|
705 | ! axis: "ls" or "sol" |
---|
706 | |
---|
707 | !============================================================================== |
---|
708 | ! Local variables: |
---|
709 | !============================================================================== |
---|
710 | !integer :: latdim,londim,altdim,timedim |
---|
711 | integer :: nvarid,ierr |
---|
712 | integer :: ctldimout |
---|
713 | ! nvarid: [netcdf] ID of a variable |
---|
714 | ! ierr: [netcdf] return error code (from called subroutines) |
---|
715 | |
---|
716 | !============================================================================== |
---|
717 | ! 1. Create (and open) output file |
---|
718 | !============================================================================== |
---|
719 | write(*,*) "creating "//trim(adjustl(filename))//'...' |
---|
720 | ierr = NF_CREATE(filename,IOR(NF_CLOBBER,NF_64BIT_OFFSET),nout) |
---|
721 | ! NB: setting NF_CLOBBER mode means that it's OK to overwrite an existing file |
---|
722 | if (ierr.NE.NF_NOERR) then |
---|
723 | WRITE(*,*)'ERROR: Impossible to create the file.' |
---|
724 | stop "" |
---|
725 | endif |
---|
726 | |
---|
727 | !============================================================================== |
---|
728 | ! 2. Define/write "dimensions" and get their IDs |
---|
729 | !============================================================================== |
---|
730 | |
---|
731 | ierr = NF_DEF_DIM(nout, "latitude", size(lat), latdimout) |
---|
732 | ierr = NF_DEF_DIM(nout, "longitude", size(lon), londimout) |
---|
733 | ierr = NF_DEF_DIM(nout, "altitude", size(alt), altdimout) |
---|
734 | if (size(ctl).ne.0) ierr = NF_DEF_DIM(nout, "index", size(ctl), ctldimout) |
---|
735 | ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout) |
---|
736 | ierr = NF_DEF_DIM(nout, "GCM_layers", GCM_layers, layerdimout) |
---|
737 | ierr = NF_DEF_DIM(nout, "GCM_interlayers",GCM_layers+1,interlayerdimout) |
---|
738 | |
---|
739 | ! End netcdf define mode |
---|
740 | ierr = NF_ENDDEF(nout) |
---|
741 | |
---|
742 | !============================================================================== |
---|
743 | ! 3. Write "Time" and its attributes |
---|
744 | !============================================================================== |
---|
745 | if (axis=="sol") then |
---|
746 | call def_var(nout,"Time","Time","days since 0000-00-0 00:00:00",1,& |
---|
747 | (/timedimout/),timevarout,ierr) |
---|
748 | else ! Ls |
---|
749 | call def_var(nout,"Time","Solar longitude","days since 0000-00-0 00:00:00",1,& |
---|
750 | (/timedimout/),timevarout,ierr) |
---|
751 | endif |
---|
752 | !============================================================================== |
---|
753 | ! 4. Write "latitude" (data and attributes) |
---|
754 | !============================================================================== |
---|
755 | |
---|
756 | call def_var(nout,"latitude","latitude","degrees_north",1,& |
---|
757 | (/latdimout/),nvarid,ierr) |
---|
758 | |
---|
759 | #ifdef NC_DOUBLE |
---|
760 | ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,lat) |
---|
761 | #else |
---|
762 | ierr = NF_PUT_VAR_REAL (nout,nvarid,lat) |
---|
763 | #endif |
---|
764 | |
---|
765 | !============================================================================== |
---|
766 | ! 4. Write "longitude" (data and attributes) |
---|
767 | !============================================================================== |
---|
768 | |
---|
769 | call def_var(nout,"longitude","East longitude","degrees_east",1,& |
---|
770 | (/londimout/),nvarid,ierr) |
---|
771 | |
---|
772 | #ifdef NC_DOUBLE |
---|
773 | ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,lon) |
---|
774 | #else |
---|
775 | ierr = NF_PUT_VAR_REAL (nout,nvarid,lon) |
---|
776 | #endif |
---|
777 | |
---|
778 | !============================================================================== |
---|
779 | ! 5. Write "altitude" (data and attributes) |
---|
780 | !============================================================================== |
---|
781 | |
---|
782 | ! Switch to netcdf define mode |
---|
783 | ierr = NF_REDEF (nout) |
---|
784 | |
---|
785 | #ifdef NC_DOUBLE |
---|
786 | ierr = NF_DEF_VAR (nout,"altitude",NF_DOUBLE,1,altdimout,nvarid) |
---|
787 | #else |
---|
788 | ierr = NF_DEF_VAR (nout,"altitude",NF_FLOAT,1,altdimout,nvarid) |
---|
789 | #endif |
---|
790 | |
---|
791 | ierr = NF_PUT_ATT_TEXT (nout,nvarid,"long_name",8,"altitude") |
---|
792 | ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',2,"km") |
---|
793 | ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',2,"up") |
---|
794 | |
---|
795 | ! End netcdf define mode |
---|
796 | ierr = NF_ENDDEF(nout) |
---|
797 | |
---|
798 | #ifdef NC_DOUBLE |
---|
799 | ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,alt) |
---|
800 | #else |
---|
801 | ierr = NF_PUT_VAR_REAL (nout,nvarid,alt) |
---|
802 | #endif |
---|
803 | |
---|
804 | !============================================================================== |
---|
805 | ! 6. Write "controle" (data and attributes) |
---|
806 | !============================================================================== |
---|
807 | |
---|
808 | if (size(ctl).ne.0) then |
---|
809 | ! Switch to netcdf define mode |
---|
810 | ierr = NF_REDEF (nout) |
---|
811 | |
---|
812 | #ifdef NC_DOUBLE |
---|
813 | ierr = NF_DEF_VAR (nout,"controle",NF_DOUBLE,1,ctldimout,nvarid) |
---|
814 | #else |
---|
815 | ierr = NF_DEF_VAR (nout,"controle",NF_FLOAT,1,ctldimout,nvarid) |
---|
816 | #endif |
---|
817 | |
---|
818 | ierr = NF_PUT_ATT_TEXT (nout,nvarid,"title",18,"Control parameters") |
---|
819 | |
---|
820 | ! End netcdf define mode |
---|
821 | ierr = NF_ENDDEF(nout) |
---|
822 | |
---|
823 | #ifdef NC_DOUBLE |
---|
824 | ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,ctl) |
---|
825 | #else |
---|
826 | ierr = NF_PUT_VAR_REAL (nout,nvarid,ctl) |
---|
827 | #endif |
---|
828 | endif |
---|
829 | |
---|
830 | end Subroutine initiate |
---|
831 | !****************************************************************************** |
---|
832 | subroutine init2(infid,lonlen,latlen,altlen,GCM_layers, & |
---|
833 | outfid,londimout,latdimout,altdimout, & |
---|
834 | layerdimout,interlayerdimout) |
---|
835 | !============================================================================== |
---|
836 | ! Purpose: |
---|
837 | ! Copy ap() , bp(), aps(), bps(), aire() and phisinit() |
---|
838 | ! from input file to outpout file |
---|
839 | !============================================================================== |
---|
840 | ! Remarks: |
---|
841 | ! The NetCDF files must be open |
---|
842 | !============================================================================== |
---|
843 | |
---|
844 | implicit none |
---|
845 | |
---|
846 | include "netcdf.inc" ! NetCDF definitions |
---|
847 | |
---|
848 | !============================================================================== |
---|
849 | ! Arguments: |
---|
850 | !============================================================================== |
---|
851 | integer, intent(in) :: infid ! NetCDF output file ID |
---|
852 | integer, intent(in) :: lonlen ! # of grid points along longitude |
---|
853 | integer, intent(in) :: latlen ! # of grid points along latitude |
---|
854 | integer, intent(in) :: altlen ! # of grid points along latitude |
---|
855 | integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers |
---|
856 | integer, intent(in) :: outfid ! NetCDF output file ID |
---|
857 | integer, intent(in) :: londimout ! longitude dimension ID |
---|
858 | integer, intent(in) :: latdimout ! latitude dimension ID |
---|
859 | integer, intent(in) :: altdimout ! altitude dimension ID |
---|
860 | integer, intent(in) :: layerdimout ! GCM_layers dimension ID |
---|
861 | integer, intent(in) :: interlayerdimout ! GCM_layers+1 dimension ID |
---|
862 | !============================================================================== |
---|
863 | ! Local variables: |
---|
864 | !============================================================================== |
---|
865 | real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates |
---|
866 | real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates |
---|
867 | real,dimension(:),allocatable :: sigma ! sigma levels |
---|
868 | real,dimension(:,:),allocatable :: aire ! mesh areas |
---|
869 | real,dimension(:,:),allocatable :: phisinit ! Ground geopotential |
---|
870 | integer :: apsid,bpsid,sigmaid,phisinitid |
---|
871 | integer :: apid,bpid |
---|
872 | integer :: ierr |
---|
873 | integer :: tmpvarid ! temporary variable ID |
---|
874 | logical :: area ! is "aire" available ? |
---|
875 | logical :: phis ! is "phisinit" available ? |
---|
876 | logical :: hybrid ! are "aps" and "bps" available ? |
---|
877 | logical :: apbp ! are "ap" and "bp" available ? |
---|
878 | |
---|
879 | !============================================================================== |
---|
880 | ! 1. Read data from input file |
---|
881 | !============================================================================== |
---|
882 | |
---|
883 | ! hybrid coordinate aps |
---|
884 | !write(*,*) "aps: altlen=",altlen," GCM_layers=",GCM_layers |
---|
885 | allocate(aps(GCM_layers),stat=ierr) |
---|
886 | if (ierr.ne.0) then |
---|
887 | write(*,*) "init2: failed to allocate aps!" |
---|
888 | stop |
---|
889 | endif |
---|
890 | ierr=NF_INQ_VARID(infid,"aps",tmpvarid) |
---|
891 | if (ierr.ne.NF_NOERR) then |
---|
892 | write(*,*) "Ooops. Failed to get aps ID. OK, will look for sigma coord." |
---|
893 | hybrid=.false. |
---|
894 | else |
---|
895 | ierr=NF_GET_VAR_REAL(infid,tmpvarid,aps) |
---|
896 | hybrid=.true. |
---|
897 | if (ierr.ne.NF_NOERR) then |
---|
898 | stop "init2 Error: Failed reading aps" |
---|
899 | endif |
---|
900 | |
---|
901 | ! hybrid coordinate bps |
---|
902 | ! write(*,*) "bps: altlen=",altlen," GCM_layers=",GCM_layers |
---|
903 | allocate(bps(GCM_layers),stat=ierr) |
---|
904 | if (ierr.ne.0) then |
---|
905 | write(*,*) "init2: failed to allocate bps!" |
---|
906 | stop |
---|
907 | endif |
---|
908 | ierr=NF_INQ_VARID(infid,"bps",tmpvarid) |
---|
909 | if (ierr.ne.NF_NOERR) then |
---|
910 | stop "init2 Error: Failed to get bps ID." |
---|
911 | endif |
---|
912 | ierr=NF_GET_VAR_REAL(infid,tmpvarid,bps) |
---|
913 | if (ierr.ne.NF_NOERR) then |
---|
914 | stop "init2 Error: Failed reading bps" |
---|
915 | endif |
---|
916 | endif |
---|
917 | |
---|
918 | ! hybrid coordinate ap |
---|
919 | allocate(ap(GCM_layers+1),stat=ierr) |
---|
920 | if (ierr.ne.0) then |
---|
921 | write(*,*) "init2: failed to allocate ap!" |
---|
922 | stop |
---|
923 | else |
---|
924 | ierr=NF_INQ_VARID(infid,"ap",tmpvarid) |
---|
925 | if (ierr.ne.NF_NOERR) then |
---|
926 | write(*,*) "Ooops. Failed to get ap ID. OK." |
---|
927 | apbp=.false. |
---|
928 | else |
---|
929 | ierr=NF_GET_VAR_REAL(infid,tmpvarid,ap) |
---|
930 | apbp=.true. |
---|
931 | if (ierr.ne.NF_NOERR) then |
---|
932 | stop "Error: Failed reading ap" |
---|
933 | endif |
---|
934 | endif |
---|
935 | endif |
---|
936 | |
---|
937 | ! hybrid coordinate bp |
---|
938 | allocate(bp(GCM_layers+1),stat=ierr) |
---|
939 | if (ierr.ne.0) then |
---|
940 | write(*,*) "init2: failed to allocate bp!" |
---|
941 | stop |
---|
942 | else |
---|
943 | ierr=NF_INQ_VARID(infid,"bp",tmpvarid) |
---|
944 | if (ierr.ne.NF_NOERR) then |
---|
945 | write(*,*) "Ooops. Failed to get bp ID. OK." |
---|
946 | apbp=.false. |
---|
947 | else |
---|
948 | ierr=NF_GET_VAR_REAL(infid,tmpvarid,bp) |
---|
949 | apbp=.true. |
---|
950 | if (ierr.ne.NF_NOERR) then |
---|
951 | stop "Error: Failed reading bp" |
---|
952 | endif |
---|
953 | endif |
---|
954 | endif |
---|
955 | |
---|
956 | ! sigma levels (if any) |
---|
957 | if (.not.hybrid) then |
---|
958 | allocate(sigma(GCM_layers),stat=ierr) |
---|
959 | if (ierr.ne.0) then |
---|
960 | write(*,*) "init2: failed to allocate sigma" |
---|
961 | stop |
---|
962 | endif |
---|
963 | ierr=NF_INQ_VARID(infid,"sigma",tmpvarid) |
---|
964 | ierr=NF_GET_VAR_REAL(infid,tmpvarid,sigma) |
---|
965 | if (ierr.ne.NF_NOERR) then |
---|
966 | stop "init2 Error: Failed reading sigma" |
---|
967 | endif |
---|
968 | endif ! of if (.not.hybrid) |
---|
969 | |
---|
970 | ! mesh area |
---|
971 | allocate(aire(lonlen,latlen),stat=ierr) |
---|
972 | if (ierr.ne.0) then |
---|
973 | write(*,*) "init2: failed to allocate aire!" |
---|
974 | stop |
---|
975 | endif |
---|
976 | ierr=NF_INQ_VARID(infid,"aire",tmpvarid) |
---|
977 | if (ierr.ne.NF_NOERR) then |
---|
978 | write(*,*)"init2 warning: Failed to get aire ID." |
---|
979 | area = .false. |
---|
980 | else |
---|
981 | ierr=NF_GET_VAR_REAL(infid,tmpvarid,aire) |
---|
982 | if (ierr.ne.NF_NOERR) then |
---|
983 | stop "init2 Error: Failed reading aire" |
---|
984 | endif |
---|
985 | area = .true. |
---|
986 | endif |
---|
987 | |
---|
988 | ! ground geopotential phisinit |
---|
989 | allocate(phisinit(lonlen,latlen),stat=ierr) |
---|
990 | if (ierr.ne.0) then |
---|
991 | write(*,*) "init2: failed to allocate phisinit!" |
---|
992 | stop |
---|
993 | endif |
---|
994 | ierr=NF_INQ_VARID(infid,"phisinit",tmpvarid) |
---|
995 | if (ierr.ne.NF_NOERR) then |
---|
996 | write(*,*)"init2 warning: Failed to get phisinit ID." |
---|
997 | phis = .false. |
---|
998 | else |
---|
999 | ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit) |
---|
1000 | if (ierr.ne.NF_NOERR) then |
---|
1001 | stop "init2 Error: Failed reading phisinit" |
---|
1002 | endif |
---|
1003 | phis = .true. |
---|
1004 | endif |
---|
1005 | |
---|
1006 | !============================================================================== |
---|
1007 | ! 2. Write |
---|
1008 | !============================================================================== |
---|
1009 | |
---|
1010 | !============================================================================== |
---|
1011 | ! 2.2. Hybrid coordinates ap() , bp(), aps() and bps() |
---|
1012 | !============================================================================== |
---|
1013 | if(hybrid) then |
---|
1014 | ! define aps |
---|
1015 | call def_var(nout,"aps","hybrid pressure at midlayers"," ",1,& |
---|
1016 | (/layerdimout/),apsid,ierr) |
---|
1017 | if (ierr.ne.NF_NOERR) then |
---|
1018 | stop "init2 Error: Failed to def_var aps" |
---|
1019 | endif |
---|
1020 | |
---|
1021 | ! write aps |
---|
1022 | #ifdef NC_DOUBLE |
---|
1023 | ierr=NF_PUT_VAR_DOUBLE(outfid,apsid,aps) |
---|
1024 | #else |
---|
1025 | ierr=NF_PUT_VAR_REAL(outfid,apsid,aps) |
---|
1026 | #endif |
---|
1027 | if (ierr.ne.NF_NOERR) then |
---|
1028 | stop "init2 Error: Failed to write aps" |
---|
1029 | endif |
---|
1030 | |
---|
1031 | ! define bps |
---|
1032 | call def_var(nout,"bps","hybrid sigma at midlayers"," ",1,& |
---|
1033 | (/layerdimout/),bpsid,ierr) |
---|
1034 | if (ierr.ne.NF_NOERR) then |
---|
1035 | stop "init2 Error: Failed to def_var bps" |
---|
1036 | endif |
---|
1037 | |
---|
1038 | ! write bps |
---|
1039 | #ifdef NC_DOUBLE |
---|
1040 | ierr=NF_PUT_VAR_DOUBLE(outfid,bpsid,bps) |
---|
1041 | #else |
---|
1042 | ierr=NF_PUT_VAR_REAL(outfid,bpsid,bps) |
---|
1043 | #endif |
---|
1044 | if (ierr.ne.NF_NOERR) then |
---|
1045 | stop "init2 Error: Failed to write bps" |
---|
1046 | endif |
---|
1047 | |
---|
1048 | if (apbp) then |
---|
1049 | ! define ap |
---|
1050 | call def_var(nout,"ap","hybrid sigma at interlayers"," ",1,& |
---|
1051 | (/interlayerdimout/),apid,ierr) |
---|
1052 | if (ierr.ne.NF_NOERR) then |
---|
1053 | stop "Error: Failed to def_var ap" |
---|
1054 | endif |
---|
1055 | |
---|
1056 | ! write ap |
---|
1057 | #ifdef NC_DOUBLE |
---|
1058 | ierr=NF_PUT_VAR_DOUBLE(outfid,apid,ap) |
---|
1059 | #else |
---|
1060 | ierr=NF_PUT_VAR_REAL(outfid,apid,ap) |
---|
1061 | #endif |
---|
1062 | if (ierr.ne.NF_NOERR) then |
---|
1063 | stop "Error: Failed to write ap" |
---|
1064 | endif |
---|
1065 | |
---|
1066 | ! define bp |
---|
1067 | call def_var(nout,"bp","hybrid sigma at interlayers"," ",1,& |
---|
1068 | (/interlayerdimout/),bpid,ierr) |
---|
1069 | if (ierr.ne.NF_NOERR) then |
---|
1070 | stop "Error: Failed to def_var bp" |
---|
1071 | endif |
---|
1072 | |
---|
1073 | ! write bp |
---|
1074 | #ifdef NC_DOUBLE |
---|
1075 | ierr=NF_PUT_VAR_DOUBLE(outfid,bpid,bp) |
---|
1076 | #else |
---|
1077 | ierr=NF_PUT_VAR_REAL(outfid,bpid,bp) |
---|
1078 | #endif |
---|
1079 | if (ierr.ne.NF_NOERR) then |
---|
1080 | stop "Error: Failed to write bp" |
---|
1081 | endif |
---|
1082 | endif ! of if (apbp) |
---|
1083 | |
---|
1084 | else |
---|
1085 | ! define sigma |
---|
1086 | call def_var(nout,"sigma","sigma at midlayers"," ",1,& |
---|
1087 | (/layerdimout/),sigmaid,ierr) |
---|
1088 | if (ierr.ne.NF_NOERR) then |
---|
1089 | stop "init2 Error: Failed to def_var sigma" |
---|
1090 | endif |
---|
1091 | ! write sigma |
---|
1092 | #ifdef NC_DOUBLE |
---|
1093 | ierr=NF_PUT_VAR_DOUBLE(outfid,sigmaid,sigma) |
---|
1094 | #else |
---|
1095 | ierr=NF_PUT_VAR_REAL(outfid,sigmaid,sigma) |
---|
1096 | #endif |
---|
1097 | if (ierr.ne.NF_NOERR) then |
---|
1098 | stop "init2 Error: Failed to write sigma" |
---|
1099 | endif |
---|
1100 | endif ! of if (hybrid) |
---|
1101 | |
---|
1102 | !============================================================================== |
---|
1103 | ! 2.2. aire() and phisinit() |
---|
1104 | !============================================================================== |
---|
1105 | |
---|
1106 | if (area) then |
---|
1107 | ! define aire |
---|
1108 | call def_var(nout,"aire","Mesh area","m2",2,& |
---|
1109 | (/londimout,latdimout/),tmpvarid,ierr) |
---|
1110 | if (ierr.ne.NF_NOERR) then |
---|
1111 | stop "init2 Error: Failed to def_var aire" |
---|
1112 | endif |
---|
1113 | |
---|
1114 | ! write aire |
---|
1115 | #ifdef NC_DOUBLE |
---|
1116 | ierr=NF_PUT_VAR_DOUBLE(outfid,tmpvarid,aire) |
---|
1117 | #else |
---|
1118 | ierr=NF_PUT_VAR_REAL(outfid,tmpvarid,aire) |
---|
1119 | #endif |
---|
1120 | if (ierr.ne.NF_NOERR) then |
---|
1121 | stop "init2 Error: Failed to write aire" |
---|
1122 | endif |
---|
1123 | endif ! of if (area) |
---|
1124 | |
---|
1125 | IF (phis) THEN |
---|
1126 | |
---|
1127 | !define phisinit |
---|
1128 | call def_var(nout,"phisinit","Ground level geopotential"," ",2,& |
---|
1129 | (/londimout,latdimout/),phisinitid,ierr) |
---|
1130 | if (ierr.ne.NF_NOERR) then |
---|
1131 | stop "init2 Error: Failed to def_var phisinit" |
---|
1132 | endif |
---|
1133 | |
---|
1134 | ! write phisinit |
---|
1135 | #ifdef NC_DOUBLE |
---|
1136 | ierr=NF_PUT_VAR_DOUBLE(outfid,phisinitid,phisinit) |
---|
1137 | #else |
---|
1138 | ierr=NF_PUT_VAR_REAL(outfid,phisinitid,phisinit) |
---|
1139 | #endif |
---|
1140 | if (ierr.ne.NF_NOERR) then |
---|
1141 | stop "init2 Error: Failed to write phisinit" |
---|
1142 | endif |
---|
1143 | |
---|
1144 | ENDIF ! of IF (phis) |
---|
1145 | |
---|
1146 | |
---|
1147 | ! Cleanup |
---|
1148 | if (allocated(aps)) deallocate(aps) |
---|
1149 | if (allocated(bps)) deallocate(bps) |
---|
1150 | if (allocated(ap)) deallocate(ap) |
---|
1151 | if (allocated(bp)) deallocate(bp) |
---|
1152 | if (allocated(sigma)) deallocate(sigma) |
---|
1153 | if (allocated(phisinit)) deallocate(phisinit) |
---|
1154 | if (allocated(aire)) deallocate(aire) |
---|
1155 | |
---|
1156 | end subroutine init2 |
---|
1157 | !****************************************************************************** |
---|
1158 | subroutine def_var(nid,name,title,units,nbdim,dim,nvarid,ierr) |
---|
1159 | !============================================================================== |
---|
1160 | ! Purpose: Write a variable (i.e: add a variable to a dataset) |
---|
1161 | ! called "name"; along with its attributes "title", "units"... |
---|
1162 | ! to a file (following the NetCDF format) |
---|
1163 | !============================================================================== |
---|
1164 | ! Remarks: |
---|
1165 | ! The NetCDF file must be open |
---|
1166 | !============================================================================== |
---|
1167 | |
---|
1168 | implicit none |
---|
1169 | |
---|
1170 | include "netcdf.inc" ! NetCDF definitions |
---|
1171 | |
---|
1172 | !============================================================================== |
---|
1173 | ! Arguments: |
---|
1174 | !============================================================================== |
---|
1175 | integer, intent(in) :: nid |
---|
1176 | ! nid: [netcdf] file ID # |
---|
1177 | character (len=*), intent(in) :: name |
---|
1178 | ! name(): [netcdf] variable's name |
---|
1179 | character (len=*), intent(in) :: title |
---|
1180 | ! title(): [netcdf] variable's "title" attribute |
---|
1181 | character (len=*), intent(in) :: units |
---|
1182 | ! unit(): [netcdf] variable's "units" attribute |
---|
1183 | integer, intent(in) :: nbdim |
---|
1184 | ! nbdim: number of dimensions of the variable |
---|
1185 | integer, dimension(nbdim), intent(in) :: dim |
---|
1186 | ! dim(nbdim): [netcdf] dimension(s) ID(s) |
---|
1187 | integer, intent(out) :: nvarid |
---|
1188 | ! nvarid: [netcdf] ID # of the variable |
---|
1189 | integer, intent(out) :: ierr |
---|
1190 | ! ierr: [netcdf] subroutines returned error code |
---|
1191 | |
---|
1192 | ! Switch to netcdf define mode |
---|
1193 | ierr=NF_REDEF(nid) |
---|
1194 | |
---|
1195 | ! Insert the definition of the variable |
---|
1196 | #ifdef NC_DOUBLE |
---|
1197 | ierr=NF_DEF_VAR(nid,adjustl(name),NF_DOUBLE,nbdim,dim,nvarid) |
---|
1198 | #else |
---|
1199 | ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid) |
---|
1200 | #endif |
---|
1201 | |
---|
1202 | ! Write the attributes |
---|
1203 | ierr=NF_PUT_ATT_TEXT(nid,nvarid,"title",len_trim(adjustl(title)),adjustl(title)) |
---|
1204 | ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units)) |
---|
1205 | |
---|
1206 | ! End netcdf define mode |
---|
1207 | ierr=NF_ENDDEF(nid) |
---|
1208 | |
---|
1209 | end subroutine def_var |
---|
1210 | !****************************************************************************** |
---|
1211 | subroutine change_time_axis(nid,ierr) |
---|
1212 | !============================================================================== |
---|
1213 | ! Purpose: |
---|
1214 | ! Read "time" variable from a dataset, convert it from "sol" (martian |
---|
1215 | ! days) to "Ls" (solar longitude) and write the result back in the file |
---|
1216 | !============================================================================== |
---|
1217 | ! Remarks: |
---|
1218 | ! The NetCDF file must be opened before this subroutine is called |
---|
1219 | !============================================================================== |
---|
1220 | |
---|
1221 | implicit none |
---|
1222 | |
---|
1223 | include "netcdf.inc" ! NetCDF definitions |
---|
1224 | |
---|
1225 | !============================================================================== |
---|
1226 | ! Arguments: |
---|
1227 | !============================================================================== |
---|
1228 | integer, intent(in) :: nid |
---|
1229 | ! nid: [netcdf] file ID |
---|
1230 | integer, intent(out) :: ierr |
---|
1231 | ! ierr: [netcdf] return error code |
---|
1232 | |
---|
1233 | !============================================================================== |
---|
1234 | ! Local variables: |
---|
1235 | !============================================================================== |
---|
1236 | integer :: nvarid |
---|
1237 | ! nvarid: ID of the "Time" variable |
---|
1238 | integer :: timelen |
---|
1239 | ! timelen: size of the arrays |
---|
1240 | integer :: timedim |
---|
1241 | ! timedim: ID of the "Time" dimension |
---|
1242 | integer i |
---|
1243 | |
---|
1244 | real, dimension(:), allocatable :: time,ls |
---|
1245 | ! time(): time, given in sols |
---|
1246 | ! ls(): time, given in Ls (solar longitude) |
---|
1247 | |
---|
1248 | !============================================================================== |
---|
1249 | ! 1. Read |
---|
1250 | !============================================================================== |
---|
1251 | |
---|
1252 | ierr=NF_INQ_DIMID(nid,"Time",timedim) |
---|
1253 | ierr=NF_INQ_VARID(nid,"Time",nvarid) |
---|
1254 | if (ierr.NE.NF_NOERR) then |
---|
1255 | write(*,*) 'ERROR in change_time_axis: Field <Time> not found' |
---|
1256 | print*, NF_STRERROR(ierr) |
---|
1257 | stop "" |
---|
1258 | endif |
---|
1259 | |
---|
1260 | ierr=NF_INQ_DIMLEN(nid,timedim,timelen) |
---|
1261 | allocate(time(timelen),ls(timelen)) |
---|
1262 | #ifdef NC_DOUBLE |
---|
1263 | ierr = NF_GET_VAR_DOUBLE(nid,nvarid,time) |
---|
1264 | #else |
---|
1265 | ierr = NF_GET_VAR_REAL(nid,nvarid,time) |
---|
1266 | #endif |
---|
1267 | if (ierr.ne.NF_NOERR) then |
---|
1268 | write(*,*) "ERROR in change_time_axis: Failed to load Time" |
---|
1269 | stop |
---|
1270 | endif |
---|
1271 | |
---|
1272 | !============================================================================== |
---|
1273 | ! 2. Convert sols to Ls |
---|
1274 | !============================================================================== |
---|
1275 | |
---|
1276 | do i=1,timelen |
---|
1277 | call sol2ls(time(i),ls(i)) |
---|
1278 | enddo |
---|
1279 | |
---|
1280 | !============================================================================== |
---|
1281 | ! 2.1 Check if there are not jumps in Ls (rounding problems in sol2ls) |
---|
1282 | !============================================================================== |
---|
1283 | |
---|
1284 | do i=1,timelen-1 |
---|
1285 | if ((ls(i+1)-ls(i)) > 350) then |
---|
1286 | write(*,*) "+ 360 deg Ls jump solved:", ls(i), ls(i+1), "at timestep", i |
---|
1287 | ls(i+1) = ls(i+1) - 360 |
---|
1288 | write(*,*) " corrected to now be:", ls(i), ls(i+1) |
---|
1289 | else if ((ls(i)-ls(i+1)) > 350) then |
---|
1290 | write(*,*) "- 360 deg Ls jump solved:", ls(i), ls(i+1), "at timestep", i |
---|
1291 | ls(i+1) = ls(i+1) + 360 |
---|
1292 | write(*,*) " corrected to now be:", ls(i), ls(i+1) |
---|
1293 | endif |
---|
1294 | enddo |
---|
1295 | |
---|
1296 | !============================================================================== |
---|
1297 | ! 3. Write |
---|
1298 | !============================================================================== |
---|
1299 | |
---|
1300 | #ifdef NC_DOUBLE |
---|
1301 | ierr = NF_PUT_VAR_DOUBLE(nid,nvarid,ls) |
---|
1302 | #else |
---|
1303 | ierr = NF_PUT_VAR_REAL(nid,nvarid,ls) |
---|
1304 | #endif |
---|
1305 | if (ierr.ne.NF_NOERR) then |
---|
1306 | write(*,*) "ERROR in change_time_axis: Failed to write Ls" |
---|
1307 | stop |
---|
1308 | endif |
---|
1309 | |
---|
1310 | end subroutine change_time_axis |
---|
1311 | !****************************************************************************** |
---|
1312 | subroutine sol2ls(sol,Ls) |
---|
1313 | !============================================================================== |
---|
1314 | ! Purpose: |
---|
1315 | ! Convert a date/time, given in sol (martian day), |
---|
1316 | ! into solar longitude date/time, in Ls (in degrees), |
---|
1317 | ! where sol=0 is (by definition) the northern hemisphere |
---|
1318 | ! spring equinox (where Ls=0). |
---|
1319 | !============================================================================== |
---|
1320 | ! Notes: |
---|
1321 | ! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year, |
---|
1322 | ! "Ls" will be increased by N*360 |
---|
1323 | ! Won't work as expected if sol is negative (then again, |
---|
1324 | ! why would that ever happen?) |
---|
1325 | !============================================================================== |
---|
1326 | |
---|
1327 | implicit none |
---|
1328 | |
---|
1329 | !============================================================================== |
---|
1330 | ! Arguments: |
---|
1331 | !============================================================================== |
---|
1332 | real,intent(in) :: sol |
---|
1333 | real,intent(out) :: Ls |
---|
1334 | |
---|
1335 | !============================================================================== |
---|
1336 | ! Local variables: |
---|
1337 | !============================================================================== |
---|
1338 | real year_day,peri_day,timeperi,e_elips,twopi,degrad |
---|
1339 | data year_day /669./ ! # of sols in a martian year |
---|
1340 | data peri_day /485.0/ |
---|
1341 | data timeperi /1.9082314/ |
---|
1342 | data e_elips /0.093358/ |
---|
1343 | data twopi /6.2831853/ ! 2.*pi |
---|
1344 | data degrad /57.2957795/ ! pi/180 |
---|
1345 | |
---|
1346 | real zanom,xref,zx0,zdx,zteta,zz |
---|
1347 | |
---|
1348 | integer count_years |
---|
1349 | integer iter |
---|
1350 | |
---|
1351 | !============================================================================== |
---|
1352 | ! 1. Compute Ls |
---|
1353 | !============================================================================== |
---|
1354 | |
---|
1355 | zz=(sol-peri_day)/year_day |
---|
1356 | zanom=twopi*(zz-nint(zz)) |
---|
1357 | xref=abs(zanom) |
---|
1358 | |
---|
1359 | ! The equation zx0 - e * sin (zx0) = xref, solved by Newton |
---|
1360 | zx0=xref+e_elips*sin(xref) |
---|
1361 | do iter=1,20 ! typically, 2 or 3 iterations are enough |
---|
1362 | zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0)) |
---|
1363 | zx0=zx0+zdx |
---|
1364 | if(abs(zdx).le.(1.e-7)) then |
---|
1365 | ! write(*,*)'iter:',iter,' |zdx|:',abs(zdx) |
---|
1366 | exit |
---|
1367 | endif |
---|
1368 | enddo |
---|
1369 | |
---|
1370 | if(zanom.lt.0.) zx0=-zx0 |
---|
1371 | |
---|
1372 | zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.)) |
---|
1373 | Ls=zteta-timeperi |
---|
1374 | |
---|
1375 | if(Ls.lt.0.) then |
---|
1376 | Ls=Ls+twopi |
---|
1377 | else |
---|
1378 | if(Ls.gt.twopi) then |
---|
1379 | Ls=Ls-twopi |
---|
1380 | endif |
---|
1381 | endif |
---|
1382 | |
---|
1383 | Ls=degrad*Ls |
---|
1384 | ! Ls is now in degrees |
---|
1385 | |
---|
1386 | !============================================================================== |
---|
1387 | ! 1. Account for (eventual) years included in input date/time sol |
---|
1388 | !============================================================================== |
---|
1389 | |
---|
1390 | count_years=0 ! initialize |
---|
1391 | zz=sol ! use "zz" to store (and work on) the value of sol |
---|
1392 | do while (zz.ge.year_day) |
---|
1393 | count_years=count_years+1 |
---|
1394 | zz=zz-year_day |
---|
1395 | enddo |
---|
1396 | |
---|
1397 | ! Add 360 degrees to Ls for every year |
---|
1398 | if (count_years.ne.0) then |
---|
1399 | Ls=Ls+360.*count_years |
---|
1400 | endif |
---|
1401 | |
---|
1402 | |
---|
1403 | end subroutine sol2ls |
---|
1404 | !****************************************************************************** |
---|
1405 | subroutine missing_value(nout,nvarid,valid_range,missing) |
---|
1406 | !============================================================================== |
---|
1407 | ! Purpose: |
---|
1408 | ! Write "valid_range" and "missing_value" attributes (of a given |
---|
1409 | ! variable) to a netcdf file |
---|
1410 | !============================================================================== |
---|
1411 | ! Remarks: |
---|
1412 | ! NetCDF file must be open |
---|
1413 | ! Variable (nvarid) ID must be set |
---|
1414 | !============================================================================== |
---|
1415 | |
---|
1416 | implicit none |
---|
1417 | |
---|
1418 | include "netcdf.inc" ! NetCDF definitions |
---|
1419 | |
---|
1420 | !============================================================================== |
---|
1421 | ! Arguments: |
---|
1422 | !============================================================================== |
---|
1423 | integer, intent(in) :: nout |
---|
1424 | ! nout: [netcdf] file ID # |
---|
1425 | integer, intent(in) :: nvarid |
---|
1426 | ! varid: [netcdf] variable ID # |
---|
1427 | real, dimension(2), intent(in) :: valid_range |
---|
1428 | ! valid_range(2): [netcdf] "valid_range" attribute (min and max) |
---|
1429 | real, intent(in) :: missing |
---|
1430 | ! missing: [netcdf] "missing_value" attribute |
---|
1431 | |
---|
1432 | !============================================================================== |
---|
1433 | ! Local variables: |
---|
1434 | !============================================================================== |
---|
1435 | integer :: ierr |
---|
1436 | ! ierr: [netcdf] subroutine returned error code |
---|
1437 | ! INTEGER nout,nvarid,ierr |
---|
1438 | ! REAL missing |
---|
1439 | ! REAL valid_range(2) |
---|
1440 | |
---|
1441 | ! Switch to netcdf dataset define mode |
---|
1442 | ierr = NF_REDEF (nout) |
---|
1443 | if (ierr.ne.NF_NOERR) then |
---|
1444 | print*,'missing_value: ' |
---|
1445 | print*, NF_STRERROR(ierr) |
---|
1446 | endif |
---|
1447 | |
---|
1448 | ! Write valid_range() attribute |
---|
1449 | #ifdef NC_DOUBLE |
---|
1450 | ierr=NF_PUT_ATT_DOUBLE(nout,nvarid,'valid_range',NF_DOUBLE,2,valid_range) |
---|
1451 | #else |
---|
1452 | ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range) |
---|
1453 | #endif |
---|
1454 | |
---|
1455 | if (ierr.ne.NF_NOERR) then |
---|
1456 | print*,'missing_value: valid_range attribution failed' |
---|
1457 | print*, NF_STRERROR(ierr) |
---|
1458 | !write(*,*) 'NF_NOERR', NF_NOERR |
---|
1459 | !CALL abort |
---|
1460 | stop "" |
---|
1461 | endif |
---|
1462 | |
---|
1463 | ! Write "missing_value" attribute |
---|
1464 | #ifdef NC_DOUBLE |
---|
1465 | ierr= NF_PUT_ATT_DOUBLE(nout,nvarid,'missing_value',NF_DOUBLE,1,missing) |
---|
1466 | #else |
---|
1467 | ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,missing) |
---|
1468 | #endif |
---|
1469 | |
---|
1470 | if (ierr.NE.NF_NOERR) then |
---|
1471 | print*, 'missing_value: missing value attribution failed' |
---|
1472 | print*, NF_STRERROR(ierr) |
---|
1473 | ! WRITE(*,*) 'NF_NOERR', NF_NOERR |
---|
1474 | ! CALL abort |
---|
1475 | stop "" |
---|
1476 | endif |
---|
1477 | |
---|
1478 | ! End netcdf dataset define mode |
---|
1479 | ierr = NF_ENDDEF(nout) |
---|
1480 | |
---|
1481 | end subroutine missing_value |
---|
1482 | !****************************************************************************** |
---|
1483 | |
---|
1484 | end program concatnc |
---|