source: trunk/LMDZ.MARS/util/concatnc.F90 @ 2236

Last change on this file since 2236 was 2147, checked in by emillour, 5 years ago

Mars GCM utilities:
Add the information in concat when time axis is set to "Solar longitude", and try to check in localtime that the input time axis in indeed in sols.
EM

  • Property svn:executable set to *
File size: 47.0 KB
Line 
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,axis)
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,axis)
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
664character (len=4),intent(in) :: axis
665! axis: "ls" or "sol"
666
667!==============================================================================
668! Local variables:
669!==============================================================================
670!integer :: latdim,londim,altdim,timedim
671integer :: nvarid,ierr
672integer :: ctldimout
673! nvarid: [netcdf] ID of a variable
674! ierr: [netcdf]  return error code (from called subroutines)
675
676!==============================================================================
677! 1. Create (and open) output file
678!==============================================================================
679write(*,*) "creating "//trim(adjustl(filename))//'...'
680ierr = NF_CREATE(filename,IOR(NF_CLOBBER,NF_64BIT_OFFSET),nout)
681! NB: setting NF_CLOBBER mode means that it's OK to overwrite an existing file
682if (ierr.NE.NF_NOERR) then
683   WRITE(*,*)'ERROR: Impossible to create the file.'
684   stop ""
685endif
686
687!==============================================================================
688! 2. Define/write "dimensions" and get their IDs
689!==============================================================================
690
691ierr = NF_DEF_DIM(nout, "latitude", size(lat), latdimout)
692ierr = NF_DEF_DIM(nout, "longitude", size(lon), londimout)
693ierr = NF_DEF_DIM(nout, "altitude", size(alt), altdimout)
694if (size(ctl).ne.0) ierr = NF_DEF_DIM(nout, "index", size(ctl), ctldimout)
695ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout)
696ierr = NF_DEF_DIM(nout, "GCM_layers", GCM_layers, layerdimout)
697ierr = NF_DEF_DIM(nout, "GCM_interlayers",GCM_layers+1,interlayerdimout)
698
699! End netcdf define mode
700ierr = NF_ENDDEF(nout)
701
702!==============================================================================
703! 3. Write "Time" and its attributes
704!==============================================================================
705if (axis=="sol") then
706  call def_var(nout,"Time","Time","days since 0000-00-0 00:00:00",1,&
707             (/timedimout/),timevarout,ierr)
708else ! Ls
709  call def_var(nout,"Time","Solar longitude","days since 0000-00-0 00:00:00",1,&
710             (/timedimout/),timevarout,ierr)
711endif
712!==============================================================================
713! 4. Write "latitude" (data and attributes)
714!==============================================================================
715
716call def_var(nout,"latitude","latitude","degrees_north",1,&
717             (/latdimout/),nvarid,ierr)
718
719#ifdef NC_DOUBLE
720ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,lat)
721#else
722ierr = NF_PUT_VAR_REAL (nout,nvarid,lat)
723#endif
724
725!==============================================================================
726! 4. Write "longitude" (data and attributes)
727!==============================================================================
728
729call def_var(nout,"longitude","East longitude","degrees_east",1,&
730             (/londimout/),nvarid,ierr)
731
732#ifdef NC_DOUBLE
733ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,lon)
734#else
735ierr = NF_PUT_VAR_REAL (nout,nvarid,lon)
736#endif
737
738!==============================================================================
739! 5. Write "altitude" (data and attributes)
740!==============================================================================
741
742! Switch to netcdf define mode
743ierr = NF_REDEF (nout)
744
745#ifdef NC_DOUBLE
746ierr = NF_DEF_VAR (nout,"altitude",NF_DOUBLE,1,altdimout,nvarid)
747#else
748ierr = NF_DEF_VAR (nout,"altitude",NF_FLOAT,1,altdimout,nvarid)
749#endif
750
751ierr = NF_PUT_ATT_TEXT (nout,nvarid,"long_name",8,"altitude")
752ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',2,"km")
753ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',2,"up")
754
755! End netcdf define mode
756ierr = NF_ENDDEF(nout)
757
758#ifdef NC_DOUBLE
759ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,alt)
760#else
761ierr = NF_PUT_VAR_REAL (nout,nvarid,alt)
762#endif
763
764!==============================================================================
765! 6. Write "controle" (data and attributes)
766!==============================================================================
767
768if (size(ctl).ne.0) then
769   ! Switch to netcdf define mode
770   ierr = NF_REDEF (nout)
771
772#ifdef NC_DOUBLE
773   ierr = NF_DEF_VAR (nout,"controle",NF_DOUBLE,1,ctldimout,nvarid)
774#else
775   ierr = NF_DEF_VAR (nout,"controle",NF_FLOAT,1,ctldimout,nvarid)
776#endif
777
778   ierr = NF_PUT_ATT_TEXT (nout,nvarid,"title",18,"Control parameters")
779
780   ! End netcdf define mode
781   ierr = NF_ENDDEF(nout)
782
783#ifdef NC_DOUBLE
784   ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,ctl)
785#else
786   ierr = NF_PUT_VAR_REAL (nout,nvarid,ctl)
787#endif
788endif
789
790end Subroutine initiate
791!******************************************************************************
792subroutine init2(infid,lonlen,latlen,altlen,GCM_layers, &
793                 outfid,londimout,latdimout,altdimout, &
794                 layerdimout,interlayerdimout)
795!==============================================================================
796! Purpose:
797! Copy ap() , bp(), aps(), bps(), aire() and phisinit()
798! from input file to outpout file
799!==============================================================================
800! Remarks:
801! The NetCDF files must be open
802!==============================================================================
803
804implicit none
805
806include "netcdf.inc" ! NetCDF definitions
807
808!==============================================================================
809! Arguments:
810!==============================================================================
811integer, intent(in) :: infid  ! NetCDF output file ID
812integer, intent(in) :: lonlen ! # of grid points along longitude
813integer, intent(in) :: latlen ! # of grid points along latitude
814integer, intent(in) :: altlen ! # of grid points along latitude
815integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers
816integer, intent(in) :: outfid ! NetCDF output file ID
817integer, intent(in) :: londimout ! longitude dimension ID
818integer, intent(in) :: latdimout ! latitude dimension ID
819integer, intent(in) :: altdimout ! altitude dimension ID
820integer, intent(in) :: layerdimout ! GCM_layers dimension ID
821integer, intent(in) :: interlayerdimout ! GCM_layers+1 dimension ID
822!==============================================================================
823! Local variables:
824!==============================================================================
825real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates
826real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates
827real,dimension(:),allocatable :: sigma ! sigma levels
828real,dimension(:,:),allocatable :: aire ! mesh areas
829real,dimension(:,:),allocatable :: phisinit ! Ground geopotential
830integer :: apsid,bpsid,sigmaid,phisinitid
831integer :: apid,bpid
832integer :: ierr
833integer :: tmpvarid ! temporary variable ID
834logical :: area ! is "aire" available ?
835logical :: phis ! is "phisinit" available ?
836logical :: hybrid ! are "aps" and "bps" available ?
837logical :: apbp ! are "ap" and "bp" available ?
838
839!==============================================================================
840! 1. Read data from input file
841!==============================================================================
842
843! hybrid coordinate aps
844!write(*,*) "aps: altlen=",altlen," GCM_layers=",GCM_layers
845allocate(aps(GCM_layers),stat=ierr)
846if (ierr.ne.0) then
847  write(*,*) "init2: failed to allocate aps!"
848  stop
849endif
850ierr=NF_INQ_VARID(infid,"aps",tmpvarid)
851if (ierr.ne.NF_NOERR) then
852  write(*,*) "Ooops. Failed to get aps ID. OK, will look for sigma coord."
853  hybrid=.false.
854else
855  ierr=NF_GET_VAR_REAL(infid,tmpvarid,aps)
856  hybrid=.true.
857  if (ierr.ne.NF_NOERR) then
858    stop "init2 Error: Failed reading aps"
859  endif
860
861  ! hybrid coordinate bps
862!  write(*,*) "bps: altlen=",altlen," GCM_layers=",GCM_layers
863  allocate(bps(GCM_layers),stat=ierr)
864  if (ierr.ne.0) then
865    write(*,*) "init2: failed to allocate bps!"
866    stop
867  endif
868  ierr=NF_INQ_VARID(infid,"bps",tmpvarid)
869  if (ierr.ne.NF_NOERR) then
870    stop "init2 Error: Failed to get bps ID."
871  endif
872  ierr=NF_GET_VAR_REAL(infid,tmpvarid,bps)
873  if (ierr.ne.NF_NOERR) then
874    stop "init2 Error: Failed reading bps"
875  endif
876endif
877
878! hybrid coordinate ap
879allocate(ap(GCM_layers+1),stat=ierr)
880if (ierr.ne.0) then
881  write(*,*) "init2: failed to allocate ap!"
882  stop
883else
884  ierr=NF_INQ_VARID(infid,"ap",tmpvarid)
885  if (ierr.ne.NF_NOERR) then
886    write(*,*) "Ooops. Failed to get ap ID. OK."
887    apbp=.false.
888  else
889    ierr=NF_GET_VAR_REAL(infid,tmpvarid,ap)
890    apbp=.true.
891    if (ierr.ne.NF_NOERR) then
892      stop "Error: Failed reading ap"
893    endif
894  endif
895endif
896
897! hybrid coordinate bp
898allocate(bp(GCM_layers+1),stat=ierr)
899if (ierr.ne.0) then
900  write(*,*) "init2: failed to allocate bp!"
901  stop
902else
903  ierr=NF_INQ_VARID(infid,"bp",tmpvarid)
904  if (ierr.ne.NF_NOERR) then
905    write(*,*) "Ooops. Failed to get bp ID. OK."
906    apbp=.false.
907  else
908    ierr=NF_GET_VAR_REAL(infid,tmpvarid,bp)
909    apbp=.true.
910    if (ierr.ne.NF_NOERR) then
911      stop "Error: Failed reading bp"
912    endif
913  endif
914endif
915
916! sigma levels (if any)
917if (.not.hybrid) then
918  allocate(sigma(GCM_layers),stat=ierr)
919  if (ierr.ne.0) then
920    write(*,*) "init2: failed to allocate sigma"
921    stop
922  endif
923  ierr=NF_INQ_VARID(infid,"sigma",tmpvarid)
924  ierr=NF_GET_VAR_REAL(infid,tmpvarid,sigma)
925  if (ierr.ne.NF_NOERR) then
926    stop "init2 Error: Failed reading sigma"
927  endif
928endif ! of if (.not.hybrid)
929
930! mesh area
931allocate(aire(lonlen,latlen),stat=ierr)
932if (ierr.ne.0) then
933  write(*,*) "init2: failed to allocate aire!"
934  stop
935endif
936ierr=NF_INQ_VARID(infid,"aire",tmpvarid)
937if (ierr.ne.NF_NOERR) then
938  write(*,*)"init2 warning: Failed to get aire ID."
939  area = .false.
940else
941  ierr=NF_GET_VAR_REAL(infid,tmpvarid,aire)
942  if (ierr.ne.NF_NOERR) then
943    stop "init2 Error: Failed reading aire"
944  endif
945  area = .true.
946endif
947
948! ground geopotential phisinit
949allocate(phisinit(lonlen,latlen),stat=ierr)
950if (ierr.ne.0) then
951  write(*,*) "init2: failed to allocate phisinit!"
952  stop
953endif
954ierr=NF_INQ_VARID(infid,"phisinit",tmpvarid)
955if (ierr.ne.NF_NOERR) then
956  write(*,*)"init2 warning: Failed to get phisinit ID."
957  phis = .false.
958else
959  ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
960  if (ierr.ne.NF_NOERR) then
961    stop "init2 Error: Failed reading phisinit"
962  endif
963  phis = .true.
964endif
965
966!==============================================================================
967! 2. Write
968!==============================================================================
969
970!==============================================================================
971! 2.2. Hybrid coordinates ap() , bp(), aps() and bps()
972!==============================================================================
973if(hybrid) then
974! define aps
975  call def_var(nout,"aps","hybrid pressure at midlayers"," ",1,&
976             (/layerdimout/),apsid,ierr)
977  if (ierr.ne.NF_NOERR) then
978    stop "init2 Error: Failed to def_var aps"
979  endif
980
981! write aps
982#ifdef NC_DOUBLE
983  ierr=NF_PUT_VAR_DOUBLE(outfid,apsid,aps)
984#else
985  ierr=NF_PUT_VAR_REAL(outfid,apsid,aps)
986#endif
987  if (ierr.ne.NF_NOERR) then
988    stop "init2 Error: Failed to write aps"
989  endif
990
991! define bps
992  call def_var(nout,"bps","hybrid sigma at midlayers"," ",1,&
993             (/layerdimout/),bpsid,ierr)
994  if (ierr.ne.NF_NOERR) then
995    stop "init2 Error: Failed to def_var bps"
996  endif
997
998! write bps
999#ifdef NC_DOUBLE
1000  ierr=NF_PUT_VAR_DOUBLE(outfid,bpsid,bps)
1001#else
1002  ierr=NF_PUT_VAR_REAL(outfid,bpsid,bps)
1003#endif
1004  if (ierr.ne.NF_NOERR) then
1005    stop "init2 Error: Failed to write bps"
1006  endif
1007
1008  if (apbp) then
1009!   define ap
1010    call def_var(nout,"ap","hybrid sigma at interlayers"," ",1,&
1011             (/interlayerdimout/),apid,ierr)
1012    if (ierr.ne.NF_NOERR) then
1013      stop "Error: Failed to def_var ap"
1014    endif
1015
1016! write ap
1017#ifdef NC_DOUBLE
1018    ierr=NF_PUT_VAR_DOUBLE(outfid,apid,ap)
1019#else
1020    ierr=NF_PUT_VAR_REAL(outfid,apid,ap)
1021#endif
1022    if (ierr.ne.NF_NOERR) then
1023      stop "Error: Failed to write ap"
1024    endif
1025
1026! define bp
1027    call def_var(nout,"bp","hybrid sigma at interlayers"," ",1,&
1028             (/interlayerdimout/),bpid,ierr)
1029    if (ierr.ne.NF_NOERR) then
1030      stop "Error: Failed to def_var bp"
1031    endif
1032
1033!   write bp
1034#ifdef NC_DOUBLE
1035    ierr=NF_PUT_VAR_DOUBLE(outfid,bpid,bp)
1036#else
1037    ierr=NF_PUT_VAR_REAL(outfid,bpid,bp)
1038#endif
1039    if (ierr.ne.NF_NOERR) then
1040      stop "Error: Failed to write bp"
1041    endif
1042  endif ! of if (apbp)
1043
1044else
1045! define sigma
1046  call def_var(nout,"sigma","sigma at midlayers"," ",1,&
1047             (/layerdimout/),sigmaid,ierr)
1048  if (ierr.ne.NF_NOERR) then
1049    stop "init2 Error: Failed to def_var sigma"
1050  endif
1051! write sigma
1052#ifdef NC_DOUBLE
1053  ierr=NF_PUT_VAR_DOUBLE(outfid,sigmaid,sigma)
1054#else
1055  ierr=NF_PUT_VAR_REAL(outfid,sigmaid,sigma)
1056#endif
1057  if (ierr.ne.NF_NOERR) then
1058    stop "init2 Error: Failed to write sigma"
1059  endif
1060endif ! of if (hybrid)
1061
1062!==============================================================================
1063! 2.2. aire() and phisinit()
1064!==============================================================================
1065
1066if (area) then
1067  ! define aire
1068  call def_var(nout,"aire","Mesh area","m2",2,&
1069           (/londimout,latdimout/),tmpvarid,ierr)
1070  if (ierr.ne.NF_NOERR) then
1071     stop "init2 Error: Failed to def_var aire"
1072  endif
1073 
1074  ! write aire
1075#ifdef NC_DOUBLE
1076  ierr=NF_PUT_VAR_DOUBLE(outfid,tmpvarid,aire)
1077#else
1078  ierr=NF_PUT_VAR_REAL(outfid,tmpvarid,aire)
1079#endif
1080  if (ierr.ne.NF_NOERR) then
1081    stop "init2 Error: Failed to write aire"
1082  endif
1083endif ! of if (area)
1084
1085IF (phis) THEN
1086
1087  !define phisinit
1088   call def_var(nout,"phisinit","Ground level geopotential"," ",2,&
1089            (/londimout,latdimout/),phisinitid,ierr)
1090   if (ierr.ne.NF_NOERR) then
1091     stop "init2 Error: Failed to def_var phisinit"
1092   endif
1093
1094  ! write phisinit
1095#ifdef NC_DOUBLE
1096  ierr=NF_PUT_VAR_DOUBLE(outfid,phisinitid,phisinit)
1097#else
1098  ierr=NF_PUT_VAR_REAL(outfid,phisinitid,phisinit)
1099#endif
1100  if (ierr.ne.NF_NOERR) then
1101    stop "init2 Error: Failed to write phisinit"
1102  endif
1103
1104ENDIF ! of IF (phis)
1105
1106
1107! Cleanup
1108if (allocated(aps)) deallocate(aps)
1109if (allocated(bps)) deallocate(bps)
1110if (allocated(ap)) deallocate(ap)
1111if (allocated(bp)) deallocate(bp)
1112if (allocated(sigma)) deallocate(sigma)
1113if (allocated(phisinit)) deallocate(phisinit)
1114if (allocated(aire)) deallocate(aire)
1115
1116end subroutine init2
1117!******************************************************************************
1118subroutine def_var(nid,name,title,units,nbdim,dim,nvarid,ierr)
1119!==============================================================================
1120! Purpose: Write a variable (i.e: add a variable to a dataset)
1121! called "name"; along with its attributes "title", "units"...
1122! to a file (following the NetCDF format)
1123!==============================================================================
1124! Remarks:
1125! The NetCDF file must be open
1126!==============================================================================
1127
1128implicit none
1129
1130include "netcdf.inc" ! NetCDF definitions
1131
1132!==============================================================================
1133! Arguments:
1134!==============================================================================
1135integer, intent(in) :: nid
1136! nid: [netcdf] file ID #
1137character (len=*), intent(in) :: name
1138! name(): [netcdf] variable's name
1139character (len=*), intent(in) :: title
1140! title(): [netcdf] variable's "title" attribute
1141character (len=*), intent(in) :: units
1142! unit(): [netcdf] variable's "units" attribute
1143integer, intent(in) :: nbdim
1144! nbdim: number of dimensions of the variable
1145integer, dimension(nbdim), intent(in) :: dim
1146! dim(nbdim): [netcdf] dimension(s) ID(s)
1147integer, intent(out) :: nvarid
1148! nvarid: [netcdf] ID # of the variable
1149integer, intent(out) :: ierr
1150! ierr: [netcdf] subroutines returned error code
1151
1152! Switch to netcdf define mode
1153ierr=NF_REDEF(nid)
1154
1155! Insert the definition of the variable
1156#ifdef NC_DOUBLE
1157ierr=NF_DEF_VAR(nid,adjustl(name),NF_DOUBLE,nbdim,dim,nvarid)
1158#else
1159ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid)
1160#endif
1161
1162! Write the attributes
1163ierr=NF_PUT_ATT_TEXT(nid,nvarid,"title",len_trim(adjustl(title)),adjustl(title))
1164ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units))
1165
1166! End netcdf define mode
1167ierr=NF_ENDDEF(nid)
1168
1169end subroutine def_var
1170!******************************************************************************
1171subroutine change_time_axis(nid,ierr)
1172!==============================================================================
1173! Purpose:
1174! Read "time" variable from a dataset, convert it from "sol" (martian
1175! days) to "Ls" (solar longitude) and write the result back in the file
1176!==============================================================================
1177! Remarks:
1178! The NetCDF file must be opened before this subroutine is called
1179!==============================================================================
1180
1181implicit none
1182
1183include "netcdf.inc" ! NetCDF definitions
1184
1185!==============================================================================
1186! Arguments:
1187!==============================================================================
1188integer, intent(in) :: nid
1189! nid: [netcdf] file ID
1190integer, intent(out) :: ierr
1191! ierr: [netcdf]  return error code
1192
1193!==============================================================================
1194! Local variables:
1195!==============================================================================
1196integer :: nvarid
1197! nvarid: ID of the "Time" variable
1198integer :: timelen
1199! timelen: size of the arrays
1200integer :: timedim
1201! timedim: ID of the "Time" dimension
1202integer i
1203
1204real, dimension(:), allocatable :: time,ls
1205! time(): time, given in sols
1206! ls(): time, given in Ls (solar longitude)
1207
1208!==============================================================================
1209! 1. Read
1210!==============================================================================
1211
1212ierr=NF_INQ_DIMID(nid,"Time",timedim)
1213ierr=NF_INQ_VARID(nid,"Time",nvarid)
1214if (ierr.NE.NF_NOERR) then
1215   write(*,*) 'ERROR in change_time_axis: Field <Time> not found'
1216   print*, NF_STRERROR(ierr)
1217   stop ""
1218endif
1219
1220ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
1221allocate(time(timelen),ls(timelen))
1222#ifdef NC_DOUBLE
1223ierr = NF_GET_VAR_DOUBLE(nid,nvarid,time)
1224#else
1225ierr = NF_GET_VAR_REAL(nid,nvarid,time)
1226#endif
1227if (ierr.ne.NF_NOERR) then
1228  write(*,*) "ERROR in change_time_axis: Failed to load Time"
1229  stop
1230endif
1231
1232!==============================================================================
1233! 2. Convert sols to Ls
1234!==============================================================================
1235
1236do i=1,timelen
1237   call sol2ls(time(i),ls(i))
1238enddo
1239
1240!==============================================================================
1241! 2.1 Check if there are not jumps in Ls (rounding problems in sol2ls)
1242!==============================================================================
1243
1244do i=1,timelen-1
1245   if ((ls(i+1)-ls(i)) > 350) then
1246       write(*,*) "+ 360 deg Ls jump solved:", ls(i), ls(i+1), "at timestep", i
1247      ls(i+1) = ls(i+1) - 360
1248       write(*,*) " corrected to now be:", ls(i), ls(i+1)
1249   else if ((ls(i)-ls(i+1)) > 350) then
1250       write(*,*) "- 360 deg Ls jump solved:", ls(i), ls(i+1), "at timestep", i
1251      ls(i+1) = ls(i+1) + 360
1252       write(*,*) " corrected to now be:", ls(i), ls(i+1)
1253   endif
1254enddo
1255
1256!==============================================================================
1257! 3. Write
1258!==============================================================================
1259
1260#ifdef NC_DOUBLE
1261ierr = NF_PUT_VAR_DOUBLE(nid,nvarid,ls)
1262#else
1263ierr = NF_PUT_VAR_REAL(nid,nvarid,ls)
1264#endif
1265if (ierr.ne.NF_NOERR) then
1266  write(*,*) "ERROR in change_time_axis: Failed to write Ls"
1267  stop
1268endif
1269
1270end subroutine change_time_axis
1271!******************************************************************************
1272subroutine sol2ls(sol,Ls)
1273!==============================================================================
1274! Purpose:
1275! Convert a date/time, given in sol (martian day),
1276! into solar longitude date/time, in Ls (in degrees),
1277! where sol=0 is (by definition) the northern hemisphere
1278!  spring equinox (where Ls=0).
1279!==============================================================================
1280! Notes:
1281! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year,
1282! "Ls" will be increased by N*360
1283! Won't work as expected if sol is negative (then again,
1284! why would that ever happen?)
1285!==============================================================================
1286
1287implicit none
1288
1289!==============================================================================
1290! Arguments:
1291!==============================================================================
1292real,intent(in) :: sol
1293real,intent(out) :: Ls
1294
1295!==============================================================================
1296! Local variables:
1297!==============================================================================
1298real year_day,peri_day,timeperi,e_elips,twopi,degrad
1299data year_day /669./            ! # of sols in a martian year
1300data peri_day /485.0/           
1301data timeperi /1.9082314/
1302data e_elips  /0.093358/
1303data twopi       /6.2831853/    ! 2.*pi
1304data degrad   /57.2957795/      ! pi/180
1305
1306real zanom,xref,zx0,zdx,zteta,zz
1307
1308integer count_years
1309integer iter
1310
1311!==============================================================================
1312! 1. Compute Ls
1313!==============================================================================
1314
1315zz=(sol-peri_day)/year_day
1316zanom=twopi*(zz-nint(zz))
1317xref=abs(zanom)
1318
1319!  The equation zx0 - e * sin (zx0) = xref, solved by Newton
1320zx0=xref+e_elips*sin(xref)
1321do iter=1,20 ! typically, 2 or 3 iterations are enough
1322   zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0))
1323   zx0=zx0+zdx
1324   if(abs(zdx).le.(1.e-7)) then
1325!      write(*,*)'iter:',iter,'     |zdx|:',abs(zdx)
1326      exit
1327   endif
1328enddo
1329
1330if(zanom.lt.0.) zx0=-zx0
1331
1332zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
1333Ls=zteta-timeperi
1334
1335if(Ls.lt.0.) then
1336   Ls=Ls+twopi
1337else
1338   if(Ls.gt.twopi) then
1339      Ls=Ls-twopi
1340   endif
1341endif
1342
1343Ls=degrad*Ls
1344! Ls is now in degrees
1345
1346!==============================================================================
1347! 1. Account for (eventual) years included in input date/time sol
1348!==============================================================================
1349
1350count_years=0 ! initialize
1351zz=sol  ! use "zz" to store (and work on) the value of sol
1352do while (zz.ge.year_day)
1353   count_years=count_years+1
1354   zz=zz-year_day
1355enddo
1356
1357! Add 360 degrees to Ls for every year
1358if (count_years.ne.0) then
1359   Ls=Ls+360.*count_years
1360endif
1361
1362
1363end subroutine sol2ls
1364!******************************************************************************
1365subroutine  missing_value(nout,nvarid,valid_range,missing)
1366!==============================================================================
1367! Purpose:
1368! Write "valid_range" and "missing_value" attributes (of a given
1369! variable) to a netcdf file
1370!==============================================================================
1371! Remarks:
1372! NetCDF file must be open
1373! Variable (nvarid) ID must be set
1374!==============================================================================
1375
1376implicit none
1377
1378include "netcdf.inc"  ! NetCDF definitions
1379
1380!==============================================================================
1381! Arguments:
1382!==============================================================================
1383integer, intent(in) :: nout
1384! nout: [netcdf] file ID #
1385integer, intent(in) :: nvarid
1386! varid: [netcdf] variable ID #
1387real, dimension(2), intent(in) :: valid_range
1388! valid_range(2): [netcdf] "valid_range" attribute (min and max)
1389real, intent(in) :: missing
1390! missing: [netcdf] "missing_value" attribute
1391
1392!==============================================================================
1393! Local variables:
1394!==============================================================================
1395integer :: ierr
1396! ierr: [netcdf] subroutine returned error code
1397!      INTEGER nout,nvarid,ierr
1398!      REAL missing
1399!      REAL valid_range(2)
1400
1401! Switch to netcdf dataset define mode
1402ierr = NF_REDEF (nout)
1403if (ierr.ne.NF_NOERR) then
1404   print*,'missing_value: '
1405   print*, NF_STRERROR(ierr)
1406endif
1407
1408! Write valid_range() attribute
1409#ifdef NC_DOUBLE
1410ierr=NF_PUT_ATT_DOUBLE(nout,nvarid,'valid_range',NF_DOUBLE,2,valid_range)
1411#else
1412ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range)
1413#endif
1414
1415if (ierr.ne.NF_NOERR) then
1416   print*,'missing_value: valid_range attribution failed'
1417   print*, NF_STRERROR(ierr)
1418   !write(*,*) 'NF_NOERR', NF_NOERR
1419   !CALL abort
1420   stop ""
1421endif
1422
1423! Write "missing_value" attribute
1424#ifdef NC_DOUBLE
1425ierr= NF_PUT_ATT_DOUBLE(nout,nvarid,'missing_value',NF_DOUBLE,1,missing)
1426#else
1427ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,missing)
1428#endif
1429
1430if (ierr.NE.NF_NOERR) then
1431   print*, 'missing_value: missing value attribution failed'
1432   print*, NF_STRERROR(ierr)
1433!    WRITE(*,*) 'NF_NOERR', NF_NOERR
1434!    CALL abort
1435   stop ""
1436endif
1437
1438! End netcdf dataset define mode
1439ierr = NF_ENDDEF(nout)
1440
1441end subroutine  missing_value
1442!******************************************************************************
1443
1444end program concatnc
Note: See TracBrowser for help on using the repository browser.