source: trunk/mesoscale/LMD_MM_MARS/SRC/PREP_MARS/readmeteo.F90 @ 16

Last change on this file since 16 was 11, checked in by aslmd, 14 years ago

spiga@svn-planeto:ajoute le modele meso-echelle martien

File size: 34.2 KB
Line 
1program readmeteo
2
3implicit none
4include "netcdf.inc"
5
6!------------------------------------------------------------------------!
7! Readmeteo prepares an initial state for the WPS pre-processor of WRF   !
8!                                                                        !
9! Input is a diagfi.nc NETCDF file from a LMD/GCM simulation             !
10!                                                                        !
11! Outputs are WRFSI intermediate format files ready for metgrid.exe      !
12!       http://www.mmm.ucar.edu/wrf/OnLineTutorial/WPS/IM_si.htm         !
13!                                                                        !
14! **** Please create a WPSFEED folder (or a symbolic link) ****          !       
15!                                                                        !
16! A. Spiga - 16/04/2007                                                  !
17!            22/05/2007 : need to run zrecast as a preliminary           !
18!            07/06/2007 : external parameter file 'readmeteo.def         !
19!            30/07/2007 : manage additional variables (tsoil, emiss,...) !     
20!            19/10/2007 : no need to run zrecast at all (eta levels)     !   
21!               12/2007 : include co2ice and emissivity                  !
22!               02/2008 : include water vapor and ice                    !
23!                                                                        !
24! spiga@lmd.jussieu.fr                                                   !     
25!------------------------------------------------------------------------!
26
27
28
29!***************************************************************************
30!***************************************************************************
31REAL, PARAMETER :: emiss_prescribed=0.95
32REAL, PARAMETER :: co2ice_prescribed=0.
33CHARACTER (LEN=3), PARAMETER :: ident='LMD'
34
35
36
37!***************************************************************************
38!***************************************************************************
39
40REAL :: ptop
41REAL, PARAMETER :: grav=3.72
42LOGICAL, PARAMETER :: blank=.false.
43
44
45
46!! Variables to be written in the output file
47!! *** NB: conformity with metgrid/src/read_met_module.F90
48CHARACTER (LEN=33) :: output
49INTEGER :: IFV
50CHARACTER (LEN=24) :: HDATE
51REAL :: XFCST
52CHARACTER (LEN=32) :: SOURCE
53CHARACTER (LEN=9) :: FIELD
54CHARACTER (LEN=25) :: UNITS
55CHARACTER (LEN=46) :: DESC
56REAL :: XLVL
57INTEGER :: NX,NY,IPROJ
58CHARACTER (LEN=8) :: STARTLOC
59REAL :: STARTLAT,STARTLON,DELTALAT,DELTALON
60REAL, POINTER, DIMENSION(:,:) :: SLAB
61
62!! Variables related to NETCDF format
63integer :: nid,nid2,nid3,nvarid,ierr,ierr2
64integer :: timedim,londim,latdim,altdim,altdim2
65integer :: rlatvdim,rlonudim
66integer :: timelen,lonlen,latlen,altlen,altlen2
67integer :: rlatvlen,rlonulen
68integer :: soillen,soildim
69integer :: nphys
70integer, dimension(4) :: corner,edges
71
72!! Intermediate data arrays
73integer :: k,l
74real, dimension(:), allocatable :: lat,lon,time,alt,aps,bps,levels
75real, dimension(:,:), allocatable :: vide,ones,ghtsfile
76real, dimension(:,:), allocatable :: interm
77real, dimension(:,:,:), allocatable :: gwparam
78real, dimension(:,:,:), allocatable :: psfile,tsfile
79real, dimension(:,:,:), allocatable :: emissfile,co2icefile
80real, dimension(:,:,:), allocatable :: tnsfile,unsfile,vnsfile
81real, dimension(:,:,:,:), allocatable :: tfile,ufile,vfile,rfile,hfile
82real, dimension(:,:,:,:), allocatable :: eta_gcm
83real, dimension(:,:,:,:), allocatable :: tfileorig,ufileorig,vfileorig
84real, dimension(:,:,:,:), allocatable :: tsoilfile
85real, dimension(:,:,:,:), allocatable :: waterfile, watericefile
86
87!! Reading the parameter file
88integer :: tmp,increment,FILES
89integer :: tmp2,tmp3,tmp4
90character*1 :: answer
91character*4 :: str
92character*2 :: str2, str3, str4
93integer, dimension(:), allocatable :: time_out
94character*13, dimension(:), allocatable :: date_out
95character*19, dimension(:), allocatable :: date_out2
96!***************************************************************************
97!***************************************************************************
98
99
100!!---------------------
101!! Open the datafile
102!!---------------------
103ierr=NF_OPEN ("input_diagfi.nc",NF_NOWRITE,nid)
104IF (ierr.NE.NF_NOERR) THEN
105   write(*,*)'**** Please create a symbolic link called input_diagfi.nc'
106   CALL ABORT
107ENDIF
108
109
110!!--------------------------
111!! Ask for data references
112!!--------------------------
113!write(*,*) "Create n files. What is n ?"
114read(*,*) FILES
115allocate(time_out(FILES))
116allocate(date_out(FILES))
117allocate(date_out2(FILES))
118
119!write(*,*) ""
120!write(*,*) "INPUT: Diagfi time"
121!write(*,*) "list of n subscripts, separated by <Return>s"
122increment=1
123do while (increment.ne.FILES+1)
124    read(*,*) tmp
125    time_out(increment)=tmp
126    increment=increment+1
127enddo
128
129!write(*,*) ""
130!write(*,*) "OUTPUT: WRF time"
131!write(*,*) "fill 4 lines indicating"
132!write(*,*) "-year (4 digits)"
133!write(*,*) "-month (2 digits)"
134!write(*,*) "-day (2 digits)"
135!write(*,*) "-hour (2 digits)"
136increment=1
137do while (increment.ne.FILES+1)
138    read(*,*) str
139    read(*,*) str2
140    read(*,*) str3
141    read(*,*) str4
142    date_out(increment)=str//'-'//str2//'-'//str3//'_'//str4
143    date_out2(increment)=str//'-'//str2//'-'//str3//'_'//str4//':00:00'
144    !print *,date_out(increment)
145    !write(*,*) "ok? (y/n)"
146    read(*,*) answer
147    if (answer.eq.'n') then
148        !write(*,*) "please write again"
149    else
150        !write(*,*) "next one, please"   
151        increment=increment+1
152    endif
153enddo
154
155
156!!-------------------
157!! Get array sizes
158!!-------------------
159SELECT CASE(ident)
160CASE('LMD')
161ierr=NF_INQ_DIMID(nid,"Time",timedim)
162CASE('OXF')
163ierr=NF_INQ_DIMID(nid,"time",timedim)
164END SELECT
165ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
166  write(*,*) "timelen: ",timelen
167
168SELECT CASE(ident)
169CASE('LMD')
170ierr=NF_INQ_DIMID(nid,"latitude",latdim)
171CASE('OXF')
172ierr=NF_INQ_DIMID(nid,"lat",latdim)
173END SELECT
174ierr=NF_INQ_DIMLEN(nid,latdim,latlen)
175  write(*,*) "latlen: ",latlen
176
177SELECT CASE(ident)
178CASE('LMD')
179ierr=NF_INQ_DIMID(nid,"longitude",londim)
180CASE('OXF')
181ierr=NF_INQ_DIMID(nid,"lon",londim)
182END SELECT
183ierr=NF_INQ_DIMLEN(nid,londim,lonlen)
184  write(*,*) "lonlen: ",lonlen
185
186SELECT CASE(ident)
187CASE('LMD')
188ierr=NF_INQ_DIMID(nid,"altitude",altdim)
189CASE('OXF')
190ierr=NF_INQ_DIMID(nid,"sigma",altdim)
191END SELECT
192ierr=NF_INQ_DIMLEN(nid,altdim,altlen)
193  write(*,*) "altlen: ",altlen
194
195
196
197!!-------------------------
198!! Allocate local arrays
199!!-------------------------
200allocate(eta_gcm(lonlen,latlen,altlen,timelen))
201allocate(tfile(lonlen,latlen,altlen,timelen))
202allocate(tsoilfile(lonlen,latlen,altlen,timelen))
203allocate(tfileorig(lonlen,latlen,altlen,timelen))
204allocate(ufile(lonlen,latlen,altlen,timelen))
205allocate(ufileorig(lonlen,latlen,altlen,timelen))
206allocate(vfile(lonlen,latlen,altlen,timelen))
207allocate(vfileorig(lonlen,latlen,altlen,timelen))
208allocate(rfile(lonlen,latlen,altlen,timelen))
209allocate(hfile(lonlen,latlen,altlen,timelen)) 
210allocate(waterfile(lonlen,latlen,altlen,timelen))
211allocate(watericefile(lonlen,latlen,altlen,timelen))
212allocate(psfile(lonlen,latlen,timelen))
213allocate(tsfile(lonlen,latlen,timelen))
214allocate(tnsfile(lonlen,latlen,timelen))
215allocate(unsfile(lonlen,latlen,timelen))
216allocate(vnsfile(lonlen,latlen,timelen))
217allocate(emissfile(lonlen,latlen,timelen))
218allocate(co2icefile(lonlen,latlen,timelen))
219allocate(interm(lonlen,latlen))
220allocate(gwparam(lonlen,latlen,5))
221allocate(ghtsfile(lonlen,latlen))       !! no scan axis
222allocate(vide(lonlen,latlen))
223allocate(ones(lonlen,latlen))
224allocate(lat(latlen), lon(lonlen), alt(altlen), time(timelen))
225allocate(aps(altlen),bps(altlen),levels(altlen))
226
227tfile(:,:,:,:)=0
228tsoilfile(:,:,:,:)=0
229tfileorig(:,:,:,:)=0
230ufileorig(:,:,:,:)=0
231vfileorig(:,:,:,:)=0
232ufile(:,:,:,:)=0
233vfile(:,:,:,:)=0
234rfile(:,:,:,:)=0
235hfile(:,:,:,:)=0
236waterfile(:,:,:,:)=0
237watericefile(:,:,:,:)=0
238psfile(:,:,:)=0
239tsfile(:,:,:)=0
240emissfile(:,:,:)=0
241co2icefile(:,:,:)=0
242tnsfile(:,:,:)=0
243unsfile(:,:,:)=0
244vnsfile(:,:,:)=0
245interm(:,:)=0
246gwparam(:,:,:)=0
247ghtsfile(:,:)=0
248vide(:,:)=0
249ones(:,:)=0
250
251
252!!-------------------
253!! Read dimensions
254!!-------------------
255
256print *,'>>> Read dimensions !'
257
258print *,'Time'
259SELECT CASE(ident)
260CASE('LMD')
261   ierr = NF_INQ_VARID (nid, "Time",nvarid)
262CASE('OXF')
263   ierr = NF_INQ_VARID (nid, "time",nvarid)
264END SELECT
265   IF (ierr .NE. NF_NOERR) THEN
266      PRINT *, "Error: Readmeteo <Time> not found"
267      stop
268   ENDIF
269#ifdef NC_DOUBLE
270   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, time)
271#else
272   ierr = NF_GET_VAR_REAL(nid, nvarid, time)
273#endif
274   print *,time(1),' ... to ... ',time(timelen)
275
276print *,'Latitude'
277SELECT CASE(ident)
278CASE('LMD')
279   ierr = NF_INQ_VARID (nid, "latitude",nvarid)
280CASE('OXF')
281   ierr = NF_INQ_VARID (nid, "lat",nvarid)
282END SELECT
283   IF (ierr .NE. NF_NOERR) THEN
284      PRINT *, "Error: Readmeteo <latitude> not found"
285      stop
286   ENDIF
287#ifdef NC_DOUBLE
288   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, lat)
289#else
290   ierr = NF_GET_VAR_REAL(nid, nvarid, lat)
291#endif
292   print *,lat(1),' ... to ... ',lat(latlen),' ... step: ',lat(latlen)-lat(latlen-1)
293
294print *,'Longitude'
295SELECT CASE(ident)
296CASE('LMD')
297   ierr = NF_INQ_VARID (nid, "longitude",nvarid)
298CASE('OXF')
299   ierr = NF_INQ_VARID (nid, "lon",nvarid)
300END SELECT
301   IF (ierr .NE. NF_NOERR) THEN
302      PRINT *, "Error: Readmeteo <longitude> not found"
303      stop
304   ENDIF
305#ifdef NC_DOUBLE
306   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, lon)
307#else
308   ierr = NF_GET_VAR_REAL(nid, nvarid, lon)
309#endif
310   print *,lon(1),' ... to ... ',lon(lonlen),' ... step: ',lon(lonlen)-lon(lonlen-1)
311
312SELECT CASE(ident)
313CASE('LMD')
314print *, 'Hybrid coordinates'
315   ierr = NF_INQ_VARID (nid, "aps", nvarid)
316   IF (ierr .NE. NF_NOERR) THEN
317      PRINT *, "Error: Readmeteo <aps> not found"
318      stop
319   ENDIF
320#ifdef NC_DOUBLE
321   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, aps)
322#else
323   ierr = NF_GET_VAR_REAL(nid, nvarid, aps)
324#endif
325   ierr = NF_INQ_VARID (nid, "bps", nvarid)
326   IF (ierr .NE. NF_NOERR) THEN
327      PRINT *, "Error: Readmeteo <bps> not found"
328      stop
329   ENDIF
330#ifdef NC_DOUBLE
331   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, bps)
332#else
333   ierr = NF_GET_VAR_REAL(nid, nvarid, bps)
334#endif
335   print *,aps(1),' ... to ... ',aps(altlen)
336   print *,bps(1),' ... to ... ',bps(altlen)
337CASE('OXF')
338   ierr = NF_INQ_VARID (nid, "sigma", nvarid)
339   IF (ierr .NE. NF_NOERR) THEN
340      PRINT *, "Error: Readmeteo <sigma> not found"
341      stop
342   ENDIF
343#ifdef NC_DOUBLE
344   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, bps)
345#else
346   ierr = NF_GET_VAR_REAL(nid, nvarid, bps)
347#endif
348   print *,bps(1),' ... to ... ',bps(altlen)
349END SELECT
350
351
352!!-------------------------------------------
353!! Reading 2D and 3D meteorological fields
354!!-------------------------------------------
355
356IF (blank .EQV. .false.) THEN
357
358print *,'>>> Read 2D optional fields !'
359
360print *,'Emissivity'
361   ierr = NF_INQ_VARID (nid,"emis",nvarid)
362IF (ierr .NE. NF_NOERR) THEN
363        PRINT *, '...warning: not found in diagfi !'
364        PRINT *, '...will be filled with a prescribed value', emiss_prescribed
365        emissfile(:,:,:)=emiss_prescribed       
366ELSE
367#ifdef NC_DOUBLE
368        ierr = NF_GET_VAR_DOUBLE(nid, nvarid, emissfile)
369#else
370        ierr = NF_GET_VAR_REAL(nid, nvarid, emissfile)
371#endif
372ENDIF   
373
374print *,'CO2 ice'
375   ierr = NF_INQ_VARID (nid,"co2ice",nvarid)
376IF (ierr .NE. NF_NOERR) THEN
377        PRINT *, '...warning: not found in diagfi !'
378        PRINT *, '...will be filled with a prescribed value', co2ice_prescribed
379        co2icefile(:,:,:)=co2ice_prescribed     
380ELSE
381#ifdef NC_DOUBLE
382        ierr = NF_GET_VAR_DOUBLE(nid, nvarid, co2icefile)
383#else
384        ierr = NF_GET_VAR_REAL(nid, nvarid, co2icefile)
385#endif
386ENDIF
387
388print *,'>>> Read 2D surface fields !'
389
390print *,'Surface Pressure'
391   ierr = NF_INQ_VARID (nid,"ps",nvarid)
392   IF (ierr .NE. NF_NOERR) THEN
393      PRINT *, "Error: Readmeteo <ps> not found"
394      stop
395   ENDIF
396#ifdef NC_DOUBLE
397   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, psfile)
398#else
399   ierr = NF_GET_VAR_REAL(nid, nvarid, psfile)
400#endif
401
402print *,'Ground Temperature'
403   ierr = NF_INQ_VARID (nid,"tsurf",nvarid)
404   IF (ierr .NE. NF_NOERR) THEN
405      PRINT *, "Error: Readmeteo <tsurf> not found"
406      stop
407   ENDIF
408#ifdef NC_DOUBLE
409   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tsfile)
410#else
411   ierr = NF_GET_VAR_REAL(nid, nvarid, tsfile)
412#endif
413
414
415!!"atmospheric" surface temperature is taken
416!!from original diagfi.nc first level
417!!... level is ~3-5 meters
418print *,'Near-Surface Temperature'
419   ierr = NF_INQ_VARID (nid,"temp",nvarid)
420   IF (ierr .NE. NF_NOERR) THEN
421      ierr = NF_INQ_VARID (nid,"t",nvarid)
422        IF (ierr .NE. NF_NOERR) THEN
423           PRINT *, "Error: Readmeteo <temp> not found"
424           stop
425        ENDIF
426   ENDIF
427#ifdef NC_DOUBLE
428   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tfileorig)
429#else
430   ierr = NF_GET_VAR_REAL(nid, nvarid, tfileorig)
431#endif
432   tnsfile=tfileorig(:,:,1,:)
433
434!!"atmospheric" surface u is taken
435!!from original diagfi.nc first level
436!!... level is ~3-5 meters
437print *,'Near-Surface Zonal Wind'
438   ierr = NF_INQ_VARID (nid,"u",nvarid)
439   IF (ierr .NE. NF_NOERR) THEN
440     PRINT *, "Error: Readmeteo <u> not found"
441     stop
442   ENDIF
443#ifdef NC_DOUBLE
444   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, ufileorig)
445#else
446   ierr = NF_GET_VAR_REAL(nid, nvarid, ufileorig)
447#endif
448   unsfile=ufileorig(:,:,1,:)
449
450!!"atmospheric" surface v is taken
451!!from original diagfi.nc first level
452!!... level is ~3-5 meters
453print *,'Near-Surface Meridional Wind'
454   ierr = NF_INQ_VARID (nid,"v",nvarid)
455   IF (ierr .NE. NF_NOERR) THEN
456     PRINT *, "Error: Readmeteo <v> not found"
457     stop
458   ENDIF
459#ifdef NC_DOUBLE
460   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, vfileorig)
461#else
462   ierr = NF_GET_VAR_REAL(nid, nvarid, vfileorig)
463#endif
464   vnsfile=vfileorig(:,:,1,:)
465
466SELECT CASE(ident)
467CASE('LMD')
468   print *,'Geopotential height at the ground'
469   ierr = NF_INQ_VARID (nid,"phisinit",nvarid)
470   IF (ierr .NE. NF_NOERR) THEN
471      PRINT *, "Error: Readmeteo <phisinit> not found"
472      stop
473   ENDIF
474#ifdef NC_DOUBLE
475   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, ghtsfile)
476#else
477   ierr = NF_GET_VAR_REAL(nid, nvarid, ghtsfile)
478#endif
479!**** from geopotential to geopotential height
480ghtsfile=ghtsfile/grav
481!****
482CASE('OXF')
483!
484! geopotential height ~ altimetry
485!
486print *,'Geopotential height at the ground from file mountain_new.nc'
487ierr=NF_OPEN("mountain_new.nc",NF_NOWRITE,nid3)
488    if (ierr.ne.NF_NOERR) then
489      write(*,*) "Error: Could not open that file either"
490      stop "Might as well stop here"
491    endif
492ierr=NF_INQ_VARID(nid3,"orography",nvarid)
493ierr=NF_GET_VAR_REAL(nid3,nvarid,ghtsfile)
494if (ierr.ne.NF_NOERR) then
495    stop "Error: Failed reading phisinit"
496endif
497ierr=NF_CLOSE(nid3)
498END SELECT
499
500
501print *,'>>> Read 3D meteorological fields ! - This may take some time ...'
502
503print *,'Temperature'
504   ierr = NF_INQ_VARID (nid,"temp",nvarid)
505   IF (ierr .NE. NF_NOERR) THEN
506      ierr = NF_INQ_VARID (nid,"t",nvarid)
507        IF (ierr .NE. NF_NOERR) THEN
508          PRINT *, "Error: Readmeteo <t> not found"
509          stop
510        ENDIF
511   ENDIF
512#ifdef NC_DOUBLE
513   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tfile)
514#else
515   ierr = NF_GET_VAR_REAL(nid, nvarid, tfile)
516#endif
517
518print *,'Zonal wind'   
519   ierr = NF_INQ_VARID (nid,"u",nvarid)
520   IF (ierr .NE. NF_NOERR) THEN
521      PRINT *, "Error: Readmeteo <u> not found"
522      stop
523   ENDIF
524#ifdef NC_DOUBLE
525   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, ufile)
526#else
527   ierr = NF_GET_VAR_REAL(nid, nvarid, ufile)
528#endif
529
530print *,'Meridional wind'
531   ierr = NF_INQ_VARID (nid,"v",nvarid)
532   IF (ierr .NE. NF_NOERR) THEN
533      PRINT *, "Error: Readmeteo <v> not found"
534      stop
535   ENDIF
536#ifdef NC_DOUBLE
537   ierr = NF_GET_VAR_DOUBLE(nid, nvarid, vfile)
538#else
539   ierr = NF_GET_VAR_REAL(nid, nvarid, vfile)
540#endif
541
542
543!!------------------------
544!! special water stuff
545!!------------------------
546    print *,'Water vapor'
547    ierr=NF_INQ_VARID(nid,"q02",nvarid)
548    if (ierr.ne.NF_NOERR) then
549      write(*,*) "...No q02 - Water vapor set to 0"
550      waterfile(:,:,:,:)=0. 
551    endif
552    ierr=NF_GET_VAR_REAL(nid,nvarid,waterfile)
553
554    print *,'Water ice'
555    ierr=NF_INQ_VARID(nid,"q01",nvarid)
556    if (ierr.ne.NF_NOERR) then
557      write(*,*) "...No q01 - Water ice set to 0" 
558      watericefile(:,:,:,:)=0.
559    endif
560    ierr=NF_GET_VAR_REAL(nid,nvarid,watericefile)
561!!------------------------
562!! special water stuff
563!!------------------------
564
565
566!SELECT CASE(ident)
567!CASE('LMD')
568
569    print *,'Soil temperatures'
570    ierr=NF_INQ_VARID(nid,"tsoil",nvarid)
571    if (ierr.ne.NF_NOERR) then
572        write(*,*) "...No tsoil - Soil temperatures set isothermal with surface temperature"
573        DO l=1,altlen
574                tsoilfile(:,:,l,:)=tsfile(:,:,:)
575        ENDDO
576    endif
577    ierr=NF_GET_VAR_REAL(nid,nvarid,tsoilfile)
578
579
580!END SELECT
581
582ierr=NF_CLOSE(nid)
583
584!!!!lott stuff
585!print *,'get lott'
586!ierr=NF_OPEN("init_lott.nc",NF_NOWRITE,nid3)
587!ierr=NF_INQ_VARID(nid3,"ZMEA",nvarid)
588!ierr=NF_GET_VAR_REAL(nid3,nvarid,interm)
589!gwparam(:,:,1)=interm(:,:)
590!ierr=NF_INQ_VARID(nid3,"ZSTD",nvarid)
591!ierr=NF_GET_VAR_REAL(nid3,nvarid,interm)
592!gwparam(:,:,2)=interm(:,:)
593!ierr=NF_INQ_VARID(nid3,"ZSIG",nvarid)
594!ierr=NF_GET_VAR_REAL(nid3,nvarid,interm)
595!gwparam(:,:,3)=interm(:,:)
596!ierr=NF_INQ_VARID(nid3,"ZGAM",nvarid)
597!ierr=NF_GET_VAR_REAL(nid3,nvarid,interm)
598!gwparam(:,:,4)=interm(:,:)
599!ierr=NF_INQ_VARID(nid3,"ZTHE",nvarid)
600!ierr=NF_GET_VAR_REAL(nid3,nvarid,interm)
601!gwparam(:,:,5)=interm(:,:)
602!ierr=NF_CLOSE(nid3)
603
604ENDIF
605
606
607!!-----------------------------
608!! Loop on the written files
609!! >>> what follows is pretty repetitive ...
610!!        gotta do something (one more loop?)
611!!-----------------------------
612
613!!! Equivalent eta levels for WRF interpolation in initialize_real
614!print *,'Computing equivalent eta levels'
615!       
616!       !ptop will be passed through RH(surface)
617!       ptop=MINVAL(aps(altlen)+bps(altlen)*psfile(:,:,:))
618!       print *, 'ptop', ptop
619!
620!print *, 'sample: equivalent eta levels at i=1,j=1'
621!DO k = 1,altlen
622!        levels(k)=-k
623!        eta_gcm(:,:,k,:)=(aps(k)+bps(k)*psfile(:,:,:)-ptop)/(psfile(:,:,:)-ptop)
624!        print *, eta_gcm(1,1,k,1)
625!END DO
626
627!!---better to pass pressure values
628!!---the eta treatment is done in module_initialize_real
629DO k = 1,altlen
630        levels(k)=-k
631        !!dummy decreasing levels
632END DO
633
634
635
636
637!! Dummy surface level is XLVL=200100.
638
639
640DO l=1,FILES
641
642
643!!---------------------------------------------
644!! Write file in intermediate format for WPS
645!! 1. Surface data
646!!---------------------------------------------
647
648!
649! These variables remain the same in the "loop"
650!
651HDATE=date_out(l)
652IFV=4
653XFCST=0.
654SOURCE=ident
655NX=lonlen
656NY=latlen
657ALLOCATE(SLAB(NX,NY))
658IPROJ=0
659STARTLOC='SWCORNER'
660STARTLAT=lat(1)
661STARTLON=lon(1)
662DELTALAT=lat(latlen)-lat(latlen-1)
663DELTALON=lon(lonlen)-lon(lonlen-1)
664!
665! Open the data destination file
666!
667output="./WPSFEED/"//ident//":"//date_out2(l)       
668open(UNIT=1,                    &
669        FILE=output,            &
670        STATUS='REPLACE',       &
671        FORM='UNFORMATTED',     &
672        ACCESS='SEQUENTIAL')
673
674!------------------------!   
675! >>> Write a variable   !
676!    ... Copy&Paste part !
677!------------------------!
678FIELD='T'
679UNITS='K'
680DESC='Atmospheric temperature'
681XLVL=200100.
682!SLAB=tsfile(:,:,time_out(l))
683!SLAB=tfileorig(:,:,1,time_out(l))
684SLAB=tnsfile(:,:,time_out(l))
685        ! And now put everything in the destination file
686        ! ... Header
687        write(1) IFV
688        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
689        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
690        ! ... Data
691        write(1) SLAB
692!print *,'The field '//DESC//' was written to '//output
693
694!------------------------!   
695! >>> Write a variable   !
696!    ... Copy&Paste part !
697!------------------------!
698FIELD='U'
699UNITS='m s{-1}'
700DESC='Zonal wind'
701XLVL=200100.
702!SLAB=ufile(:,:,1,time_out(l))
703!SLAB=ufileorig(:,:,1,time_out(l))
704SLAB=unsfile(:,:,time_out(l))
705        ! And now put everything in the destination file
706        ! ... Header
707        write(1) IFV
708        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
709        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
710        ! ... Data
711        write(1) SLAB
712!print *,'The field '//DESC//' was written to '//output
713
714!------------------------!   
715! >>> Write a variable   !
716!    ... Copy&Paste part !
717!------------------------!
718FIELD='V'
719UNITS='m s{-1}'
720DESC='Meridional wind'
721XLVL=200100.
722!SLAB=vfile(:,:,1,time_out(l))
723!SLAB=vfileorig(:,:,1,time_out(l))
724SLAB=vnsfile(:,:,time_out(l))
725        ! And now put everything in the destination file
726        ! ... Header
727        write(1) IFV
728        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
729        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
730        ! ... Data
731        write(1) SLAB
732!print *,'The field '//DESC//' was written to '//output
733
734!------------------------!   
735! >>> Write a variable   !
736!    ... Copy&Paste part !
737!------------------------!
738FIELD='RH'
739UNITS=''
740DESC='Customized 2D static field'
741XLVL=200100.
742!SLAB=co2icefile(:,:,time_out(l))
743SLAB=vide(:,:)
744!SLAB=vide(:,:)+ptop
745        ! And now put everything in the destination file
746        ! ... Header
747        write(1) IFV
748        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
749        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
750        ! ... Data
751        write(1) SLAB
752!print *,'The field '//DESC//' was written to '//output
753
754!------------------------!   
755! >>> Write a variable   !
756!    ... Copy&Paste part !
757!------------------------!
758FIELD='SOILHGT'
759UNITS='m'
760DESC='Terrain field of source analysis'  !Ground geopotential height
761XLVL=200100.
762SLAB=ghtsfile(:,:)
763        ! And now put everything in the destination file
764        ! ... Header
765        write(1) IFV
766        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
767        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
768        ! ... Data
769        write(1) SLAB
770!print *,'The field '//DESC//' was written to '//output
771
772!------------------------!   
773! >>> Write a variable   !
774!    ... Copy&Paste part !
775!------------------------!
776FIELD='PSFC'
777UNITS='Pa'
778DESC='Surface Pressure'
779XLVL=200100.
780SLAB=psfile(:,:,time_out(l))
781        ! And now put everything in the destination file
782        ! ... Header
783        write(1) IFV
784        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
785        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
786        ! ... Data
787        write(1) SLAB
788!print *,'The field '//DESC//' was written to '//output
789
790!------------------------!   
791! >>> Write a variable   !
792!    ... Copy&Paste part !
793!------------------------!
794FIELD='ST000010'
795UNITS=''
796DESC='Emissivity'
797XLVL=200100.
798SLAB=emissfile(:,:,time_out(l))
799        ! And now put everything in the destination file
800        ! ... Header
801        write(1) IFV
802        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
803        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
804        ! ... Data
805        write(1) SLAB
806!print *,'The field '//DESC//' was written to '//output
807
808!------------------------!   
809! >>> Write a variable   !
810!    ... Copy&Paste part !
811!------------------------!
812FIELD='ST010040'
813UNITS=''
814DESC='CO2 ice'
815XLVL=200100.
816SLAB=co2icefile(:,:,time_out(l))
817        ! And now put everything in the destination file
818        ! ... Header
819        write(1) IFV
820        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
821        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
822        ! ... Data
823        write(1) SLAB
824!print *,'The field '//DESC//' was written to '//output
825
826!------------------------!   
827! >>> Write a variable   !
828!    ... Copy&Paste part !
829!------------------------!
830FIELD='ST040100'
831UNITS=''
832DESC='ZMEA'
833XLVL=200100.
834!SLAB=vide(:,:)
835SLAB=gwparam(:,:,1)
836        ! And now put everything in the destination file
837        ! ... Header
838        write(1) IFV
839        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
840        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
841        ! ... Data
842        write(1) SLAB
843!print *,'The field '//DESC//' was written to '//output
844
845!------------------------!   
846! >>> Write a variable   !
847!    ... Copy&Paste part !
848!------------------------!
849FIELD='ST100200'
850UNITS=''
851DESC='ZSTD'
852XLVL=200100.
853!SLAB=vide(:,:)
854SLAB=gwparam(:,:,2)
855        ! And now put everything in the destination file
856        ! ... Header
857        write(1) IFV
858        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
859        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
860        ! ... Data
861        write(1) SLAB
862!print *,'The field '//DESC//' was written to '//output
863
864!------------------------!   
865! >>> Write a variable   !
866!    ... Copy&Paste part !
867!------------------------!
868FIELD='LANDSEA'
869UNITS='0/1 Flag'
870DESC='Land/Sea flag'
871XLVL=200100.
872SLAB=ones(:,:)
873        ! And now put everything in the destination file
874        ! ... Header
875        write(1) IFV
876        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
877        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
878        ! ... Data
879        write(1) SLAB
880!print *,'The field '//DESC//' was written to '//output
881
882!------------------------!   
883! >>> Write a variable   !
884!    ... Copy&Paste part !
885!------------------------!
886FIELD='SKINTEMP'
887UNITS='K'
888DESC='Ground temperature'
889XLVL=200100.
890!SLAB=vide(:,:)
891SLAB=tsfile(:,:,time_out(l))
892        ! And now put everything in the destination file
893        ! ... Header
894        write(1) IFV
895        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
896        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
897        ! ... Data
898        write(1) SLAB
899!print *,'The field '//DESC//' was written to '//output
900
901!------------------------!   
902! >>> Write a variable   !
903!    ... Copy&Paste part !
904!------------------------!
905FIELD='SM000010'
906UNITS=''
907DESC='ZSIG'
908XLVL=200100.
909!SLAB=vide(:,:)
910SLAB=gwparam(:,:,3)
911        ! And now put everything in the destination file
912        ! ... Header
913        write(1) IFV
914        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
915        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
916        ! ... Data
917        write(1) SLAB
918!print *,'The field '//DESC//' was written to '//output
919
920!------------------------!   
921! >>> Write a variable   !
922!    ... Copy&Paste part !
923!------------------------!
924FIELD='SM010040'
925UNITS=''
926DESC='ZGAM'
927XLVL=200100.
928!SLAB=vide(:,:)
929SLAB=gwparam(:,:,4)
930        ! And now put everything in the destination file
931        ! ... Header
932        write(1) IFV
933        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
934        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
935        ! ... Data
936        write(1) SLAB
937!print *,'The field '//DESC//' was written to '//output
938
939!------------------------!   
940! >>> Write a variable   !
941!    ... Copy&Paste part !
942!------------------------!
943FIELD='SM040100'
944UNITS=''
945DESC='ZTHE'
946XLVL=200100.
947!SLAB=vide(:,:)
948SLAB=gwparam(:,:,5)
949        ! And now put everything in the destination file
950        ! ... Header
951        write(1) IFV
952        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
953        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
954        ! ... Data
955        write(1) SLAB
956!print *,'The field '//DESC//' was written to '//output
957
958!------------------------!   
959! >>> Write a variable   !
960!    ... Copy&Paste part !
961!------------------------!
962FIELD='SM100200'
963UNITS='fraction'
964DESC='Relative humidity'
965XLVL=200100.
966SLAB=vide(:,:)
967        ! And now put everything in the destination file
968        ! ... Header
969        write(1) IFV
970        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
971        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
972        ! ... Data
973        write(1) SLAB
974!print *,'The field '//DESC//' was written to '//output
975
976
977!------------------------!
978! >>> Write a variable   !
979!    ... Copy&Paste part !
980!------------------------!
981FIELD='HV'
982UNITS='kg/kg'
983DESC='Water vapor'
984XLVL=200100.
985SLAB=waterfile(:,:,1,time_out(l))
986        ! And now put everything in the destination file
987        ! ... Header
988        write(1) IFV
989        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
990        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
991        ! ... Data
992        write(1) SLAB
993!print *,'The field '//DESC//' was written to '//output
994
995!------------------------!
996! >>> Write a variable   !
997!    ... Copy&Paste part !
998!------------------------!
999FIELD='HI'
1000UNITS='kg/kg'
1001DESC='Water ice'
1002XLVL=200100.
1003SLAB=watericefile(:,:,1,time_out(l))
1004        ! And now put everything in the destination file
1005        ! ... Header
1006        write(1) IFV
1007        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1008        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
1009        ! ... Data
1010        write(1) SLAB
1011!print *,'The field '//DESC//' was written to '//output
1012
1013!------------------------!
1014! >>> Write a variable   !
1015!    ... Copy&Paste part !
1016!------------------------!
1017FIELD='TSOIL'
1018UNITS='K'
1019DESC='Soil temperature'
1020XLVL=200100.
1021SLAB=tsoilfile(:,:,1,time_out(l))
1022        ! And now put everything in the destination file
1023        ! ... Header
1024        write(1) IFV
1025        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1026        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
1027        ! ... Data
1028        write(1) SLAB
1029!print *,'The field '//DESC//' was written to '//output
1030
1031
1032
1033!!----------------------------------------------------
1034!! Write file in intermediate format for WPS
1035!! 2. 3D meteorological data
1036!! >>> same stuff, but handle vertical profile
1037!! >>> !! XLVL must be decreasing (pressure levels)
1038!!----------------------------------------------------
1039
1040!------------------------!   
1041! >>> Write a variable   !
1042!    ... Copy&Paste part !
1043!------------------------!
1044FIELD='T'
1045UNITS='K'
1046DESC='Atmospheric temperature'
1047DO k = 1,altlen
1048        XLVL=levels(k)
1049        SLAB=tfile(:,:,k,time_out(l))
1050                ! And now put everything in the destination file
1051                ! ... Header
1052        write(1) IFV
1053        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1054        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
1055                ! ... Data
1056                write(1) SLAB
1057END DO
1058!print *,'The field '//DESC//' was written to '//output
1059
1060!------------------------!   
1061! >>> Write a variable   !
1062!    ... Copy&Paste part !
1063!------------------------!
1064FIELD='U'
1065UNITS='m s{-1}'
1066DESC='Zonal wind'
1067DO k = 1,altlen
1068        XLVL=levels(k)
1069        SLAB=ufile(:,:,k,time_out(l))
1070                ! And now put everything in the destination file
1071                ! ... Header
1072        write(1) IFV
1073        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1074        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
1075                ! ... Data
1076                write(1) SLAB
1077END DO
1078!print *,'The field '//DESC//' was written to '//output
1079
1080!------------------------!   
1081! >>> Write a variable   !
1082!    ... Copy&Paste part !
1083!------------------------!
1084FIELD='V'
1085UNITS='m s{-1}'
1086DESC='Meridional wind'
1087DO k = 1,altlen
1088        XLVL=levels(k)
1089        SLAB=vfile(:,:,k,time_out(l))
1090                ! And now put everything in the destination file
1091                ! ... Header
1092        write(1) IFV
1093        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1094        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
1095                ! ... Data
1096                write(1) SLAB
1097END DO
1098!print *,'The field '//DESC//' was written to '//output
1099
1100!------------------------!   
1101! >>> Write a variable   !
1102!    ... Copy&Paste part !
1103!------------------------!
1104FIELD='RH'
1105UNITS=''
1106DESC='Customized 2D static field'
1107DESC='Eta-levels from the GCM'
1108DESC='Pressure values from the GCM'
1109DO k = 1,altlen
1110        XLVL=levels(k)
1111SELECT CASE(ident)
1112CASE('LMD')
1113        SLAB=aps(k)+bps(k)*psfile(:,:,time_out(l))
1114CASE('OXF')
1115        SLAB=bps(k)*psfile(:,:,time_out(l))
1116END SELECT
1117                ! And now put everything in the destination file
1118                ! ... Header
1119        write(1) IFV
1120        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1121        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
1122                ! ... Data
1123                write(1) SLAB
1124END DO
1125!print *,'The field '//DESC//' was written to '//output
1126
1127!------------------------!   
1128! >>> Write a variable   !
1129!    ... Copy&Paste part !
1130!------------------------!
1131FIELD='HGT'
1132UNITS='m'
1133DESC='Height'
1134DO k = 1,altlen
1135        XLVL=levels(k)
1136!!*******
1137!! PROVISOIRE: normalement, il faudrait la hauteur geopotentielle
1138!!*******
1139!!however, not used by initialize_real
1140SELECT CASE(ident)
1141CASE('LMD')
1142        SLAB=10000.*alog(610./(aps(k)+bps(k)*psfile(:,:,time_out(l))))
1143CASE('OXF')
1144        SLAB=10000.*alog(610./(bps(k)*psfile(:,:,time_out(l))))
1145END SELECT
1146                ! And now put everything in the destination file
1147                ! ... Header
1148        write(1) IFV
1149        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1150        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON 
1151                ! ... Data
1152                write(1) SLAB
1153END DO
1154!print *,'The field '//DESC//' was written to '//output
1155
1156!------------------------!
1157! >>> Write a variable   !
1158!    ... Copy&Paste part !
1159!------------------------!
1160FIELD='HV'
1161UNITS='kg/kg'
1162DESC='Water vapor'
1163DO k = 1,altlen
1164        XLVL=levels(k)
1165        SLAB=waterfile(:,:,k,time_out(l))
1166                ! And now put everything in the destination file
1167                ! ... Header
1168        write(1) IFV
1169        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1170        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
1171                ! ... Data
1172                write(1) SLAB
1173END DO
1174!print *,'The field '//DESC//' was written to '//output
1175
1176!------------------------!
1177! >>> Write a variable   !
1178!    ... Copy&Paste part !
1179!------------------------!
1180FIELD='HI'
1181UNITS='kg/kg'
1182DESC='Water ice'
1183DO k = 1,altlen
1184        XLVL=levels(k)
1185        SLAB=watericefile(:,:,k,time_out(l))
1186                ! And now put everything in the destination file
1187                ! ... Header
1188        write(1) IFV
1189        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1190        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
1191                ! ... Data
1192                write(1) SLAB
1193END DO
1194!print *,'The field '//DESC//' was written to '//output
1195
1196!------------------------!
1197! >>> Write a variable   !
1198!    ... Copy&Paste part !
1199!------------------------!
1200FIELD='TSOIL'
1201UNITS='K'
1202DESC='Soil temperature'
1203DO k = 1,altlen
1204        XLVL=levels(k)
1205        SLAB=tsoilfile(:,:,k,time_out(l))
1206                ! And now put everything in the destination file
1207                ! ... Header
1208        write(1) IFV
1209        write(1) HDATE,XFCST,SOURCE,FIELD,UNITS,DESC,XLVL,NX,NY,IPROJ
1210        write(1) STARTLOC,STARTLAT,STARTLON,DELTALAT,DELTALON
1211                ! ... Data
1212                write(1) SLAB
1213END DO
1214!print *,'The field '//DESC//' was written to '//output
1215
1216
1217print *,'****done file '//output, int(100.*float(l)/float(FILES)), ' % '
1218close(1)
1219
1220DEALLOCATE(SLAB)
1221
1222END DO !! end of the files loop
1223
1224
1225!!-------------------------
1226!! DeAllocate local arrays
1227!!-------------------------
1228deallocate(eta_gcm)
1229deallocate(tfile)
1230deallocate(tsoilfile)
1231deallocate(tfileorig)
1232deallocate(ufile)
1233deallocate(ufileorig)
1234deallocate(vfile)
1235deallocate(vfileorig)
1236deallocate(rfile)
1237deallocate(hfile) 
1238deallocate(waterfile)
1239deallocate(watericefile)
1240deallocate(psfile)
1241deallocate(tsfile)
1242deallocate(tnsfile)
1243deallocate(unsfile)
1244deallocate(vnsfile)
1245deallocate(emissfile)
1246deallocate(co2icefile)
1247deallocate(ghtsfile)    !! no scan axis
1248deallocate(vide)
1249deallocate(ones)
1250deallocate(lat, lon, alt, time)
1251deallocate(aps,bps,levels)
1252
1253
1254print *, '------------------------'
1255print *, 'End of pre-WPS process !'
1256print *, '------------------------'
1257print *, 'Here is what you should set in namelist.wps:'
1258print *, " start_date = '"//date_out2(1)//"'"
1259print *, " end_date = '"//date_out2(FILES)//"'"
1260print *, '------------------------'
1261print *, 'Any date between those is ok'
1262print *, 'for example, second available date is ... '//date_out2(2)
1263print *, '------------------------'
1264
1265end
Note: See TracBrowser for help on using the repository browser.