source: trunk/LMDZ.MARS/util/extract.F90 @ 2140

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

Mars GCM utilities:
Very minor syntax fixes for picky recent gfortran compilers.
EM

File size: 23.2 KB
Line 
1program extract
2
3! program to extract (ie: interpolates) pointwise values of an atmospheric
4! variable from a 'zrecast'ed diagfi file (works if altitude is geometrical
5! height or a pressure vertical coordinates)
6! user has to specify:
7! - name of input file
8! - date (in sols) offset wrt the input file (e.g. if the input file "begins"
9!   at Ls=0, then the offset is 0; if the input file begins at Ls=30, the
10!   offset date corresponding to the first 3 months is 61+66+66=193 sols, etc.)
11! - the "extraction mode":
12!      1: extract individual values; user will specify values of
13!         lon lat alt Ls LT (all on a same line)
14!         on as many lines as there are sought values
15!      2: extract a profile: user will specify on a first line the values of
16!         lon lat Ls LT (all on a same line)
17!         and then only specify values of altitudes (m or Pa depending on the
18!         coordinate in the input file), one per line, at which values are
19!         sought
20! - output values are sent to (ASCII) output file 'infile_var_.dat', where
21!   'infile' is the input file name (without trailing '.nc') and
22!   'var' is the sought variable, for extraction mode 1 as
23!   lines of "lon lat alt Ls LT value" and for a profile (extraction mode 2)
24!   as lines of "alt value"
25!
26!  NB: If there is no data to do an appropriate interpolation to extract
27!      the sought value, then a "missing_value" (taken from the variable's
28!      attribute in the input file, most likely -9.99E33) is returned.
29!
30! EM. Sept. 2011
31
32use netcdf
33implicit none
34
35! Input file:
36character(len=256) :: infile
37character(len=256) :: outfile
38
39character (len=256) :: text ! to store some text
40character (len=64) :: varname ! to store the name of the variable to retreive
41
42! NetCDF stuff
43integer :: status ! NetCDF return code
44integer :: inid ! NetCDF file IDs
45integer :: varid ! to store the ID of a variable
46integer :: lat_dimid,lon_dimid,alt_dimid,time_dimid
47integer :: datashape(4)
48
49real,dimension(:),allocatable :: longitude ! longitude
50integer lonlen ! # of grid points along longitude
51real,dimension(:),allocatable :: latitude ! latitude
52integer latlen ! # of grid points along latitude
53real,dimension(:),allocatable :: altitude ! can be geometric heights or pressure
54integer altlen ! # of grid point along altitude
55real,dimension(:),allocatable :: time ! time
56integer timelen ! # of points along time
57character :: alttype ! altitude coord. type:'z' (altitude, m) 'p' (pressure, Pa)
58real,dimension(:,:,:,:),allocatable :: field
59real :: missing_value ! value to denote non-existant data
60real :: starttimeoffset ! offset (in sols) wrt Ls=0 of sol 0 in file
61integer :: extract_mode ! 1: point-by-point extraction 2: extract a profile 
62
63! point at which data is sought:
64real :: lon,lat,alt,Ls,LT,value
65real :: sol ! sol GCM date corresponding to sought Ls and LT
66
67integer :: nb
68
69!===============================================================================
70! 1.1 Input file
71!===============================================================================
72
73write(*,*) ""
74write(*,*) " Program valid for diagfi.nc or concatnc.nc files"
75write(*,*) "  processed by zrecast "
76write(*,*) " Enter input file name:"
77
78read(*,'(a)') infile
79write(*,*) ""
80
81! open input file
82status=nf90_open(infile,NF90_NOWRITE,inid)
83if (status.ne.nf90_noerr) then
84  write(*,*)"Failed to open datafile ",trim(infile)
85  write(*,*)trim(nf90_strerror(status))
86 stop
87endif
88
89write(*,*) " Beginning date of the file?"
90write(*,*) " (i.e. number of sols since Ls=0 that the Time=0.0 in the input"
91write(*,*) "  file corresponds to)"
92read(*,*) starttimeoffset
93
94write(*,*) " Extraction mode?"
95write(*,*) " ( 1: pointwise extraction , 2: profile extraction)"
96read(*,*,iostat=status) extract_mode
97if ((status.ne.0).or.(extract_mode.lt.1) &
98                 .or.(extract_mode.gt.2)) then
99  write(*,*) "Error: invalid extraction mode:",extract_mode
100  stop
101endif
102
103!===============================================================================
104! 1.2 Input variable to extract
105!===============================================================================
106
107write(*,*) "Enter variable to extract:"
108read(*,*) varname
109! check that input file contains that variable
110status=nf90_inq_varid(inid,trim(varname),varid)
111if (status.ne.nf90_noerr) then
112  write(*,*) "Failed to find variable ",trim(varname)," in ",trim(infile)
113  write(*,*)trim(nf90_strerror(status))
114  write(*,*) " Might as well stop here."
115  stop
116endif
117
118!===============================================================================
119! 1.3 Get grids for dimensions lon,lat,alt,time
120!===============================================================================
121
122! latitude
123status=nf90_inq_dimid(inid,"latitude",lat_dimid)
124if (status.ne.nf90_noerr) then
125  write(*,*)"Failed to find latitude dimension"
126  write(*,*)trim(nf90_strerror(status))
127 stop
128endif
129status=nf90_inquire_dimension(inid,lat_dimid,len=latlen)
130if (status.ne.nf90_noerr) then
131  write(*,*)"Failed to find latitude length"
132  write(*,*)trim(nf90_strerror(status))
133endif
134allocate(latitude(latlen))
135status=nf90_inq_varid(inid,"latitude",varid)
136if (status.ne.nf90_noerr) then
137  write(*,*) "Failed to find latitude ID"
138  write(*,*)trim(nf90_strerror(status))
139  stop
140endif
141! read latitude
142status=nf90_get_var(inid,varid,latitude)
143if (status.ne.nf90_noerr) then
144  write(*,*) "Failed to load latitude"
145  write(*,*)trim(nf90_strerror(status))
146  stop
147endif
148
149!longitude
150status=nf90_inq_dimid(inid,"longitude",lon_dimid)
151if (status.ne.nf90_noerr) then
152  write(*,*)"Failed to find longitude dimension"
153  write(*,*)trim(nf90_strerror(status))
154 stop
155endif
156status=nf90_inquire_dimension(inid,lon_dimid,len=lonlen)
157if (status.ne.nf90_noerr) then
158  write(*,*)"Failed to find longitude length"
159  write(*,*)trim(nf90_strerror(status))
160endif
161allocate(longitude(lonlen))
162status=nf90_inq_varid(inid,"longitude",varid)
163if (status.ne.nf90_noerr) then
164  write(*,*) "Failed to find longitude ID"
165  write(*,*)trim(nf90_strerror(status))
166  stop
167endif
168! read longitude
169status=nf90_get_var(inid,varid,longitude)
170if (status.ne.nf90_noerr) then
171  write(*,*) "Failed to load longitude"
172  write(*,*)trim(nf90_strerror(status))
173  stop
174endif
175
176!time
177status=nf90_inq_dimid(inid,"Time",time_dimid)
178if (status.ne.nf90_noerr) then
179  write(*,*)"Failed to find Time dimension"
180  write(*,*)trim(nf90_strerror(status))
181 stop
182endif
183status=nf90_inquire_dimension(inid,time_dimid,len=timelen)
184if (status.ne.nf90_noerr) then
185  write(*,*)"Failed to find Time length"
186  write(*,*)trim(nf90_strerror(status))
187endif
188allocate(time(timelen))
189status=nf90_inq_varid(inid,"Time",varid)
190if (status.ne.nf90_noerr) then
191  write(*,*) "Failed to find Time ID"
192  write(*,*)trim(nf90_strerror(status))
193  stop
194endif
195! read Time
196status=nf90_get_var(inid,varid,time)
197if (status.ne.nf90_noerr) then
198  write(*,*) "Failed to load Time"
199  write(*,*)trim(nf90_strerror(status))
200  stop
201endif
202! add the offset to time(:)
203time(:)=time(:)+starttimeoffset
204
205!altitude
206status=nf90_inq_dimid(inid,"altitude",alt_dimid)
207if (status.ne.nf90_noerr) then
208  write(*,*)"Failed to find altitude dimension"
209  write(*,*)trim(nf90_strerror(status))
210 stop
211endif
212status=nf90_inquire_dimension(inid,alt_dimid,len=altlen)
213if (status.ne.nf90_noerr) then
214  write(*,*)"Failed to find altitude length"
215  write(*,*)trim(nf90_strerror(status))
216endif
217allocate(altitude(altlen))
218status=nf90_inq_varid(inid,"altitude",varid)
219if (status.ne.nf90_noerr) then
220  write(*,*) "Failed to find altitude ID"
221  write(*,*)trim(nf90_strerror(status))
222  stop
223endif
224! read altitude
225status=nf90_get_var(inid,varid,altitude)
226if (status.ne.nf90_noerr) then
227  write(*,*) "Failed to load altitude"
228  write(*,*)trim(nf90_strerror(status))
229  stop
230endif
231! check altitude attribute "units" to find out altitude type
232status=nf90_get_att(inid,varid,"units",text)
233if (status.ne.nf90_noerr) then
234  write(*,*) "Failed to load altitude units attribute"
235  write(*,*)trim(nf90_strerror(status))
236  stop
237else
238  if (trim(text).eq."Pa") then
239    alttype="p" ! pressure coordinate
240  else if (trim(text).eq."m") then
241    alttype="z" ! altitude coordinate
242  else
243    write(*,*)" I do not understand this unit ",trim(text)," for altitude!"
244    stop
245  endif
246endif
247
248!===============================================================================
249! 1.3 Get input dataset
250!===============================================================================
251status=nf90_inq_varid(inid,trim(varname),varid)
252if (status.ne.nf90_noerr) then
253  write(*,*) "Failed to find variable ",trim(varname)," in ",trim(infile)
254  write(*,*)trim(nf90_strerror(status))
255  write(*,*) " Might as well stop here."
256  stop
257endif
258
259! sanity checks on the variable
260status=nf90_inquire_variable(inid,varid,ndims=nb,dimids=datashape)
261if (status.ne.nf90_noerr) then
262  write(*,*) "Failed to obtain information on variable ",trim(varname)
263  write(*,*)trim(nf90_strerror(status))
264  write(*,*) " Might as well stop here."
265  stop
266else
267  ! check that it is a 4D variable
268  if (nb.ne.4) then
269    write(*,*) "Error, expected a 4D (lon-lat-alt-time) variable!"
270    stop
271  endif
272  ! check that its dimensions are indeed lon,lat,alt,time
273  if (datashape(1).ne.lon_dimid) then
274    write(*,*) "Error, expected first dimension to be longitude!"
275    stop
276  endif
277  if (datashape(2).ne.lat_dimid) then
278    write(*,*) "Error, expected second dimension to be latitude!"
279    stop
280  endif
281  if (datashape(3).ne.alt_dimid) then
282    write(*,*) "Error, expected third dimension to be altitude!"
283    stop
284  endif
285  if (datashape(4).ne.time_dimid) then
286    write(*,*) "Error, expected fourth dimension to be time!"
287    stop
288  endif
289endif
290
291allocate(field(lonlen,latlen,altlen,timelen))
292
293! load dataset
294status=nf90_get_var(inid,varid,field)
295if (status.ne.nf90_noerr) then
296  write(*,*) "Failed to load ",trim(varname)
297  write(*,*)trim(nf90_strerror(status))
298  stop
299else
300  write(*,*) "Loaded ",trim(varname)
301endif
302! get dataset's missing_value attribute
303status=nf90_get_att(inid,varid,"missing_value",missing_value)
304if (status.ne.nf90_noerr) then
305  write(*,*) "Failed to load missing_value attribute"
306  write(*,*)trim(nf90_strerror(status))
307  stop
308else
309  write(*,'("  with missing_value attribute : ",1pe12.5)') missing_value
310endif
311
312!===============================================================================
313! 2. Process and interpolate
314!===============================================================================
315
316!===============================================================================
317! 2.1 Create output file
318!===============================================================================
319
320outfile=trim(infile(1:index(infile,".nc",back=.true.)-1))//"_"//&
321        trim(varname)//".dat"
322open(42,file=outfile,form="formatted")
323write(*,*) "Output file is: ",trim(outfile)
324
325!===============================================================================
326! 2.2 Extract values
327!===============================================================================
328
329if (extract_mode==1) then ! pointwise extraction
330 write(*,*) " Enter values of: lon lat alt Ls LT (all on the same line,"
331 write(*,*) "                    using as many lines as sought data values)"
332 write(*,*) "    (enter anything else, e.g. blank line, nonesense,... to quit)"
333else ! extract_mode==2 , profile extraction
334 write(*,*) " Enter values of: lon lat Ls LT (all on the same line)"
335 read(*,*) lon,lat,Ls,LT
336 write(*,*) " Enter values of altitude (one per line)"
337 write(*,*) "    (enter anything else, e.g. blank line, nonesense,... to quit)"
338endif
339do
340  ! 2.1 read coordinates and do some sanity checks
341  text="" !initialize text to empty string
342  read(*,'(a)') text ! store input as text (to spot empty input line)
343  if (len_trim(adjustl(text)).eq.0) exit
344 
345  if (extract_mode==1) then
346    read(text,*,iostat=status) lon,lat,alt,Ls,LT
347    ! Ls in degrees, LT (local true solar time) in hours, i.e.  in [0:24]
348  else ! extract_mode==2 , read alt only
349    read(text,*,iostat=status) alt
350  endif
351  if (status.ne.0) exit
352
353  if ((lon.lt.-360.).or.(lon.gt.360.)) then
354    write(*,*) "Unexpected value for lon: ",lon
355    stop
356  endif
357  ! we want lon in [-180:180]
358  if (lon.lt.-180.) lon=lon+360.
359  if (lon.gt.180.) lon=lon-360.
360 
361  if ((lat.lt.-90.).or.(lat.gt.90.)) then
362    write(*,*) "Unexpected value for lat: ",lat
363    stop
364  endif
365 
366  if ((Ls.lt.0.).or.(Ls.gt.360.)) then
367    write(*,*) "Unexpected value for Ls: ",Ls
368    stop
369  endif
370 
371  if ((LT.lt.0.).or.(LT.gt.24.)) then
372    write(*,*) "Unexpected value for LT: ",LT
373    stop
374  endif
375 
376  ! 2.2 compute GCM sol date corresponding to sought  Ls and LT
377  call ls2sol(Ls,sol)
378  !shift 'sol' decimal part to ensure a compatible LT with the desired one
379  sol=floor(sol)+(LT-lon/15.)/24.
380  ! handle bordeline cases:
381  if (sol.gt.669) sol=sol-669
382  if (sol.lt.0) sol=sol+669
383 
384!  write(*,*) " Ls=",Ls," LT=",LT," => sol=",sol
385 
386  ! 2.3 do the interpolation
387  call extraction(lon,lat,alt,sol,&
388                  lonlen,latlen,altlen,timelen,&
389                  longitude,latitude,altitude,time,&
390                  field,missing_value,alttype,varname,value)
391  ! 2.4 Write value to output
392  if (extract_mode==1) then ! pointwise extraction
393    write(42,'(6(1x,1pe12.5))')lon,lat,alt,Ls,LT,value
394  else ! profile
395    write(42,'(6(1x,1pe12.5))')alt,value
396  endif
397enddo ! of do while
398
399! close output file
400close(42)
401
402end program extract
403
404!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
405subroutine extraction(lon,lat,alt,sol,&
406                  lonlen,latlen,altlen,timelen,&
407                  longitude,latitude,altitude,time,&
408                  field,missing_value,alttype,varname,value)
409
410implicit none
411! Arguments:
412real,intent(in) :: lon  ! sought longitude (deg, in [-180:180])
413real,intent(in) :: lat  ! sought latitude (deg, in [-90:90])
414real,intent(in) :: alt  ! sought altitude (m or Pa)
415real,intent(in) :: sol  ! sought date (sols)
416integer,intent(in) :: lonlen
417integer,intent(in) :: latlen
418integer,intent(in) :: altlen
419integer,intent(in) :: timelen
420real,intent(in) :: longitude(lonlen)
421real,intent(in) :: latitude(latlen)
422real,intent(in) :: altitude(altlen)
423real,intent(in) :: time(timelen)
424real,intent(in) :: field(lonlen,latlen,altlen,timelen)
425real,intent(in) :: missing_value ! default value in GCM file for "no data"
426character,intent(in) :: alttype ! altitude coord. type:'z' (altitude, m)
427                                !                      'p' (pressure, Pa)
428character(len=*),intent(in) :: varname ! variable name (in GCM file)
429real,intent(out) :: value
430
431! Local variables:
432real,save :: prev_lon=-99 ! previous value of 'lon' routine was called with
433real,save :: prev_lat=-99 ! previous value of 'lat' routine was called with
434real,save :: prev_alt=-9.e20 ! ! previous value of 'alt'
435real,save :: prev_sol=-99 ! previous value of 'sol' routine was called with
436
437! encompasing indexes:
438integer,save :: ilon_inf=-1,ilon_sup=-1 ! along longitude
439integer,save :: ilat_inf=-1,ilat_sup=-1 ! along latitude
440integer,save :: ialt_inf=-1,ialt_sup=-1 ! along altitude
441integer,save :: itim_inf=-1,itim_sup=-1 ! along time
442
443! intermediate interpolated values
444real :: t_interp(2,2,2) ! after time interpolation
445real :: zt_interp(2,2) ! after altitude interpolation
446real :: yzt_interp(2) ! after latitude interpolation
447real :: coeff ! interpolation coefficient
448
449integer :: i
450
451! By default, set value to missing_value
452value=missing_value
453
454! 1. Find encompassing indexes
455if (lon.ne.prev_lon) then
456  do i=1,lonlen-1
457    if (longitude(i).le.lon) then
458      ilon_inf=i
459    else
460      exit
461    endif
462  enddo
463  ilon_sup=ilon_inf+1
464endif ! of if (lon.ne.prev_lon)
465!write(*,*) 'ilon_inf=',ilon_inf," longitude(ilon_inf)=",longitude(ilon_inf)
466
467if (lat.ne.prev_lat) then
468  ! recall that latitudes start from north pole
469  do i=1,latlen-1
470    if (latitude(i).ge.lat) then
471      ilat_inf=i
472    else
473      exit
474    endif
475  enddo
476  ilat_sup=ilat_inf+1
477endif ! of if (lat.ne.prev_lat)
478!write(*,*) 'ilat_inf=',ilat_inf," latitude(ilat_inf)=",latitude(ilat_inf)
479
480if (alt.ne.prev_alt) then
481  if (alttype.eq.'p') then ! pressures are ordered from max to min
482    !handle special case for alt not in the altitude(1:altlen) interval
483    if ((alt.gt.altitude(1)).or.(alt.lt.altitude(altlen))) then
484      ialt_inf=-1
485      ialt_sup=-1
486      ! return to main program (with value=missing_value)
487      return
488    else ! general case
489      do i=1,altlen-1
490        if (altitude(i).ge.alt) then
491          ialt_inf=i
492        else
493          exit
494        endif
495      enddo
496      ialt_sup=ialt_inf+1
497    endif ! of if ((alt.gt.altitude(1)).or.(alt.lt.altitude(altlen)))
498  else ! altitudes (m) are ordered from min to max
499    !handle special case for alt not in the altitude(1:altlen) interval
500    if ((alt.lt.altitude(1)).or.(alt.gt.altitude(altlen))) then
501      ialt_inf=-1
502      ialt_sup=-1
503      ! return to main program (with value=missing_value)
504      return
505    else ! general case
506      do i=1,altlen-1
507        if (altitude(i).le.alt) then
508          ialt_inf=i
509        else
510          exit
511        endif
512      enddo
513      ialt_sup=ialt_inf+1
514    endif ! of if ((alt.lt.altitude(1)).or.(alt.gt.altitude(altlen)))
515  endif ! of if (alttype.eq.'p')
516endif ! of if (alt.ne.prev_alt)
517!write(*,*) 'ialt_inf=',ialt_inf," altitude(ialt_inf)=",altitude(ialt_inf)
518
519if (sol.ne.prev_sol) then
520  !handle special case for sol not in the time(1):time(timenlen) interval
521  if ((sol.lt.time(1)).or.(sol.gt.time(timelen))) then
522    itim_inf=-1
523    itim_sup=-1
524    ! return to main program (with value=missing_value)
525    return
526  else ! general case
527    do i=1,timelen-1
528      if (time(i).le.sol) then
529        itim_inf=i
530      else
531        exit
532      endif
533    enddo
534    itim_sup=itim_inf+1
535  endif ! of if ((sol.lt.time(1)).or.(sol.gt.time(timenlen)))
536endif ! of if (sol.ne.prev_sol)
537!write(*,*) 'itim_inf=',itim_inf," time(itim_inf)=",time(itim_inf)
538!write(*,*) 'itim_sup=',itim_sup," time(itim_sup)=",time(itim_sup)
539
540! 2. Interpolate
541! check that there are no "missing_value" in the field() elements we need
542! otherwise return to main program (with value=missing_value)
543if (field(ilon_inf,ilat_inf,ialt_inf,itim_inf).eq.missing_value) return
544if (field(ilon_inf,ilat_inf,ialt_inf,itim_sup).eq.missing_value) return
545if (field(ilon_inf,ilat_inf,ialt_sup,itim_inf).eq.missing_value) return
546if (field(ilon_inf,ilat_inf,ialt_sup,itim_sup).eq.missing_value) return
547if (field(ilon_inf,ilat_sup,ialt_inf,itim_inf).eq.missing_value) return
548if (field(ilon_inf,ilat_sup,ialt_inf,itim_sup).eq.missing_value) return
549if (field(ilon_inf,ilat_sup,ialt_sup,itim_inf).eq.missing_value) return
550if (field(ilon_inf,ilat_sup,ialt_sup,itim_sup).eq.missing_value) return
551if (field(ilon_sup,ilat_inf,ialt_inf,itim_inf).eq.missing_value) return
552if (field(ilon_sup,ilat_inf,ialt_inf,itim_sup).eq.missing_value) return
553if (field(ilon_sup,ilat_inf,ialt_sup,itim_inf).eq.missing_value) return
554if (field(ilon_sup,ilat_inf,ialt_sup,itim_sup).eq.missing_value) return
555if (field(ilon_sup,ilat_sup,ialt_inf,itim_inf).eq.missing_value) return
556if (field(ilon_sup,ilat_sup,ialt_inf,itim_sup).eq.missing_value) return
557if (field(ilon_sup,ilat_sup,ialt_sup,itim_inf).eq.missing_value) return
558if (field(ilon_sup,ilat_sup,ialt_sup,itim_sup).eq.missing_value) return
559
560! 2.1 Linear interpolation in time
561coeff=(sol-time(itim_inf))/(time(itim_sup)-time(itim_inf))
562t_interp(1,1,1)=field(ilon_inf,ilat_inf,ialt_inf,itim_inf)+ &
563                coeff*(field(ilon_inf,ilat_inf,ialt_inf,itim_sup)- &
564                       field(ilon_inf,ilat_inf,ialt_inf,itim_inf))
565t_interp(1,1,2)=field(ilon_inf,ilat_inf,ialt_sup,itim_inf)+ &
566                coeff*(field(ilon_inf,ilat_inf,ialt_sup,itim_sup)- &
567                       field(ilon_inf,ilat_inf,ialt_sup,itim_inf))
568t_interp(1,2,1)=field(ilon_inf,ilat_sup,ialt_inf,itim_inf)+ &
569                coeff*(field(ilon_inf,ilat_sup,ialt_inf,itim_sup)- &
570                       field(ilon_inf,ilat_sup,ialt_inf,itim_inf))
571t_interp(1,2,2)=field(ilon_inf,ilat_sup,ialt_sup,itim_inf)+ &
572                coeff*(field(ilon_inf,ilat_sup,ialt_sup,itim_sup)- &
573                       field(ilon_inf,ilat_sup,ialt_sup,itim_inf))
574t_interp(2,1,1)=field(ilon_sup,ilat_inf,ialt_inf,itim_inf)+ &
575                coeff*(field(ilon_sup,ilat_inf,ialt_inf,itim_sup)- &
576                       field(ilon_sup,ilat_inf,ialt_inf,itim_inf))
577t_interp(2,1,2)=field(ilon_sup,ilat_inf,ialt_sup,itim_inf)+ &
578                coeff*(field(ilon_sup,ilat_inf,ialt_sup,itim_sup)- &
579                       field(ilon_sup,ilat_inf,ialt_sup,itim_inf))
580t_interp(2,2,1)=field(ilon_sup,ilat_sup,ialt_inf,itim_inf)+ &
581                coeff*(field(ilon_sup,ilat_sup,ialt_inf,itim_sup)- &
582                       field(ilon_sup,ilat_sup,ialt_inf,itim_inf))
583t_interp(2,2,2)=field(ilon_sup,ilat_sup,ialt_sup,itim_inf)+ &
584                coeff*(field(ilon_sup,ilat_sup,ialt_sup,itim_sup)- &
585                       field(ilon_sup,ilat_sup,ialt_sup,itim_inf))
586
587! 2.2 Vertical interpolation
588if (((varname=='rho').or.(varname=='pressure')).and.(alttype=='z')) then
589  ! do the interpolation on the log of the quantity
590  coeff=(alt-altitude(ialt_inf))/(altitude(ialt_sup)-altitude(ialt_inf))
591  zt_interp(1,1)=log(t_interp(1,1,1))+coeff* &
592                             (log(t_interp(1,1,2))-log(t_interp(1,1,1)))
593  zt_interp(1,2)=log(t_interp(1,2,1))+coeff* &
594                             (log(t_interp(1,2,2))-log(t_interp(1,2,1)))
595  zt_interp(2,1)=log(t_interp(2,1,1))+coeff* &
596                             (log(t_interp(2,1,2))-log(t_interp(2,1,1)))
597  zt_interp(2,2)=log(t_interp(2,2,1))+coeff* &
598                             (log(t_interp(2,2,2))-log(t_interp(2,2,1)))
599  zt_interp(1:2,1:2)=exp(zt_interp(1:2,1:2))
600else ! general case
601  coeff=(alt-altitude(ialt_inf))/(altitude(ialt_sup)-altitude(ialt_inf))
602  zt_interp(1,1)=t_interp(1,1,1)+coeff*(t_interp(1,1,2)-t_interp(1,1,1))
603  zt_interp(1,2)=t_interp(1,2,1)+coeff*(t_interp(1,2,2)-t_interp(1,2,1))
604  zt_interp(2,1)=t_interp(2,1,1)+coeff*(t_interp(2,1,2)-t_interp(2,1,1))
605  zt_interp(2,2)=t_interp(2,2,1)+coeff*(t_interp(2,2,2)-t_interp(2,2,1))
606endif
607
608! 2.3 Latitudinal interpolation
609coeff=(lat-latitude(ilat_inf))/(latitude(ilat_sup)-latitude(ilat_inf))
610yzt_interp(1)=zt_interp(1,1)+coeff*(zt_interp(1,2)-zt_interp(1,1))
611yzt_interp(2)=zt_interp(2,1)+coeff*(zt_interp(2,2)-zt_interp(2,1))
612
613! 2.4 longitudinal interpolation
614coeff=(lon-longitude(ilon_inf))/(longitude(ilon_sup)-longitude(ilon_inf))
615value=yzt_interp(1)+coeff*(yzt_interp(2)-yzt_interp(1))
616
617end subroutine extraction
618
619!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
620subroutine ls2sol(ls,sol)
621
622implicit none
623!  Arguments:
624real,intent(in) :: ls
625real,intent(out) :: sol
626
627!  Local:
628double precision xref,zx0,zteta,zz
629!xref: mean anomaly, zteta: true anomaly, zx0: eccentric anomaly
630double precision year_day
631double precision peri_day,timeperi,e_elips
632double precision pi,degrad
633parameter (year_day=668.6d0) ! number of sols in a martian year
634parameter (peri_day=485.35d0) ! date (in sols) of perihelion
635!  timeperi: 2*pi*( 1 - Ls(perihelion)/ 360 ); Ls(perihelion)=250.99
636parameter (timeperi=1.90258341759902d0)
637parameter (e_elips=0.0934d0)  ! eccentricity of orbit
638parameter (pi=3.14159265358979d0)
639parameter (degrad=57.2957795130823d0)
640
641      if (abs(ls).lt.1.0e-5) then
642         if (ls.ge.0.0) then
643            sol = 0.0
644         else
645            sol = real(year_day)
646         end if
647         return
648      end if
649
650      zteta = ls/degrad + timeperi
651      zx0 = 2.0*datan(dtan(0.5*zteta)/dsqrt((1.+e_elips)/(1.-e_elips)))
652      xref = zx0-e_elips*dsin(zx0)
653      zz = xref/(2.*pi)
654      sol = real(zz*year_day + peri_day)
655      if (sol.lt.0.0) sol = sol + real(year_day)
656      if (sol.ge.year_day) sol = sol - real(year_day)
657
658
659end subroutine ls2sol
660
Note: See TracBrowser for help on using the repository browser.