/[lmdze]/trunk/libf/Test_ozonecm/test_ozonecm.f90
ViewVC logotype

Contents of /trunk/libf/Test_ozonecm/test_ozonecm.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 22 - (show annotations)
Fri Jul 31 15:18:47 2009 UTC (14 years, 11 months ago) by guez
File size: 3302 byte(s)
Superficial modifications
1 program test_ozonecm
2
3 ! This is a program in Fortran 95.
4
5 ! This program computes values of ozone abundance from Royer, on a
6 ! latitude-pressure grid, and writes the values to a NetCDF file.
7 ! The pressure grid is "presnivs" from "disvert".
8
9 use dimens_m, only: jjm, llm
10 use comvert, only: pa, disvert, ap, bp, preff, presnivs
11 use ozonecm_m, only: ozonecm
12 use phyetat0_m, only: rlat
13 use numer_rec, only: arth
14 use netcdf95, only: nf95_create, nf95_def_dim, nf95_def_var, nf95_put_att, &
15 nf95_enddef, nf95_put_var, nf95_close
16 use netcdf, only: nf90_clobber, nf90_float, nf90_global
17
18 implicit none
19
20 real p(llm + 1) ! pressure at LMDZ layer interface, in Pa
21 real wo(jjm + 1, llm, 360) ! column density of ozone in a cell, in kDU
22 real tro3(jjm + 1, llm, 360) ! mole fraction of ozone
23 integer julien, k
24 real, parameter:: RG = 9.80665 ! acceleration of gravity, in m s-2
25 real, parameter:: dobson_u = 2.1415e-05 ! Dobson unit, in kg m-2
26
27 ! For NetCDF:
28 integer ncid, dimid_time, dimid_plev, dimid_latitude
29 integer varid_time, varid_plev, varid_latitude, varid_tro3
30
31 !---------------------
32
33 pa = 5e4
34 call disvert
35 p = ap + bp * preff
36 rlat = arth(-90., 180. / jjm, jjm + 1)
37
38 do julien = 1, 360
39 wo(:, :, julien) = ozonecm(REAL(julien), spread(p, dim=1, ncopies=jjm+1))
40 end do
41
42 forall (k=1: llm)
43 tro3(:, k, :) &
44 = wo(:, k, :) * 1e3 * dobson_u / (p(k)-p(k+1)) * rg / 48. * 29.
45 end forall
46
47 call nf95_create("O3_Royer.nc", nf90_clobber, ncid)
48 call nf95_def_dim(ncid, "time", 360, dimid_time)
49 call nf95_def_dim(ncid, "plev", llm, dimid_plev)
50 call nf95_def_dim(ncid, "latitude", jjm + 1, dimid_latitude)
51 call nf95_def_var(ncid, "time", nf90_float, dimid_time, varid_time)
52 call nf95_def_var(ncid, "plev", nf90_float, dimid_plev, varid_plev)
53 call nf95_def_var(ncid, "latitude", nf90_float, dimid_latitude, &
54 varid_latitude)
55 call nf95_def_var(ncid, "tro3", nf90_float, &
56 (/dimid_latitude, dimid_plev, dimid_time/), varid_tro3)
57
58 call nf95_put_att(ncid, varid_time, "units", "days since 1976-1-1")
59 call nf95_put_att(ncid, varid_time, "calendar", "360_day")
60 call nf95_put_att(ncid, varid_time, "standard_name", "time")
61
62 call nf95_put_att(ncid, varid_plev, "units", "mbar")
63 call nf95_put_att(ncid, varid_plev, "long_name", "pressure")
64 call nf95_put_att(ncid, varid_plev, "standard_name", "air_pressure")
65
66 call nf95_put_att(ncid, varid_latitude, "units", "degrees_north")
67 call nf95_put_att(ncid, varid_latitude, "standard_name", "latitude")
68
69 call nf95_put_att(ncid, varid_tro3, "long_name", "O3 mole fraction")
70 call nf95_put_att(ncid, varid_tro3, "standard_name", &
71 "mole_fraction_of_ozone_in_air")
72
73 call nf95_put_att(ncid, nf90_global, "title", &
74 "ozone climatology from J.-F. Royer")
75
76 call nf95_enddef(ncid)
77
78 call nf95_put_var(ncid, varid_time, (/(julien - 0.5, julien = 1, 360)/))
79
80 call nf95_put_var(ncid, varid_plev, presnivs(llm:1:-1) / 100.)
81 ! sort in ascending order and convert to hPa
82
83 call nf95_put_var(ncid, varid_latitude, rlat)
84
85 call nf95_put_var(ncid, varid_tro3, tro3(:, llm:1:-1, :))
86 ! "plev" is in ascending order whereas "p" is descending
87
88 call nf95_close(ncid)
89
90 print *, 'The file "O3_Royer.nc" has been created.'
91
92 end program test_ozonecm

  ViewVC Help
Powered by ViewVC 1.1.21