source: trunk/LMDZ.GENERIC/utilities/concatnc.F90 @ 3537

Last change on this file since 3537 was 1905, checked in by mturbet, 7 years ago

import some generic utilities from Mars GCM

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