/[lmdze]/trunk/libf/phylmd/Mobidic/regr_pr_comb_coefoz.f90
ViewVC logotype

Diff of /trunk/libf/phylmd/Mobidic/regr_pr_comb_coefoz.f90

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 11 by guez, Thu Jun 5 12:43:08 2008 UTC revision 18 by guez, Thu Aug 7 12:29:13 2008 UTC
# Line 48  contains Line 48  contains
48      ! -- packs the coefficients to the "physics" horizontal grid ;      ! -- packs the coefficients to the "physics" horizontal grid ;
49      ! -- combines the eight coefficients to define the five module variables.      ! -- combines the eight coefficients to define the five module variables.
50    
51      ! We assume that, in "coefoz_LMDZ.nc", the pressure levels are in hPa      use netcdf95, only: nf95_open, nf95_close
     ! and strictly increasing.  
   
     use netcdf95, only: nf95_open, nf95_close, nf95_get_coord  
52      use netcdf, only: nf90_nowrite      use netcdf, only: nf90_nowrite
53      use regr_pr_coefoz, only: regr_pr_av_coefoz, regr_pr_int_coefoz      use regr_pr_coefoz, only: regr_pr_av_coefoz, regr_pr_int_coefoz
54      use phyetat0_m, only: rlat      use phyetat0_m, only: rlat
# Line 59  contains Line 56  contains
56      integer, intent(in):: julien ! jour julien, 1 <= julien <= 360      integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
57    
58      ! Variables local to the procedure:      ! Variables local to the procedure:
     integer ncid ! for NetCDF  
59    
60      real, pointer:: plev(:)      integer ncid ! for NetCDF
     ! (pressure level of input data, converted to Pa, in strictly  
     ! increasing order)  
   
     integer n_plev ! number of pressure levels in the input data  
   
     real, allocatable:: press_in_edg(:)  
     ! (edges of pressure intervals for input data, in Pa, in strictly  
     ! increasing order)  
61    
62      real coefoz(klon, llm)      real coefoz(klon, llm)
63      ! (temporary storage for an ozone coefficient)      ! (temporary storage for an ozone coefficient)
# Line 97  contains Line 85  contains
85    
86      call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)      call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
87    
88      call nf95_get_coord(ncid, "plev", plev)      call regr_pr_av_coefoz(ncid, "a2", julien, a2)
     ! Convert from hPa to Pa because "regr_pr_av" and "regr_pr_int" require so:  
     plev = plev * 100.  
     n_plev = size(plev)  
   
     ! Compute edges of pressure intervals:  
     allocate(press_in_edg(n_plev + 1))  
     press_in_edg(1) = 0.  
     ! We choose edges halfway in logarithm:  
     forall (k = 2:n_plev) press_in_edg(k) = sqrt(plev(k - 1) * plev(k))  
     press_in_edg(n_plev + 1) = huge(0.)  
     ! (infinity, but any value guaranteed to be greater than the  
     ! surface pressure would do)  
   
     call regr_pr_av_coefoz(ncid, "a2", julien, press_in_edg, a2)  
89    
90      call regr_pr_av_coefoz(ncid, "a4", julien, press_in_edg, a4_mass)      call regr_pr_av_coefoz(ncid, "a4", julien, a4_mass)
91      a4_mass = a4_mass * 48. / 29.      a4_mass = a4_mass * 48. / 29.
92    
93      call regr_pr_av_coefoz(ncid, "a6", julien, press_in_edg, a6)      call regr_pr_av_coefoz(ncid, "a6", julien, a6)
94    
95      ! Compute "a6_mass" avoiding underflow, do not divide by 1e4      ! Compute "a6_mass" avoiding underflow, do not divide by 1e4
96      ! before dividing by molecular mass:      ! before dividing by molecular mass:
# Line 127  contains Line 101  contains
101      ! (We use as few local variables as possible, in order to spare      ! (We use as few local variables as possible, in order to spare
102      ! main memory.)      ! main memory.)
103    
104      call regr_pr_av_coefoz(ncid, "P_net_Mob", julien, press_in_edg, c_Mob)      call regr_pr_av_coefoz(ncid, "P_net_Mob", julien, c_Mob)
105    
106      call regr_pr_av_coefoz(ncid, "r_Mob", julien, press_in_edg, coefoz)      call regr_pr_av_coefoz(ncid, "r_Mob", julien, coefoz)
107      c_mob = c_mob - a2 * coefoz      c_mob = c_mob - a2 * coefoz
108    
109      call regr_pr_int_coefoz(ncid, "Sigma_Mob", julien, plev, top_value=0., &      call regr_pr_int_coefoz(ncid, "Sigma_Mob", julien, top_value=0., v3=coefoz)
          v3=coefoz)  
110      c_mob = (c_mob - a6 * coefoz) * 48. / 29.      c_mob = (c_mob - a6 * coefoz) * 48. / 29.
111    
112      call regr_pr_av_coefoz(ncid, "temp_Mob", julien, press_in_edg, coefoz)      call regr_pr_av_coefoz(ncid, "temp_Mob", julien, coefoz)
113      c_mob = c_mob - a4_mass * coefoz      c_mob = c_mob - a4_mass * coefoz
114    
115      call regr_pr_av_coefoz(ncid, "R_Het", julien, press_in_edg, r_het_interm)      call regr_pr_av_coefoz(ncid, "R_Het", julien, r_het_interm)
116      ! Heterogeneous chemistry is only at high latitudes:      ! Heterogeneous chemistry is only at high latitudes:
117      forall (k = 1: llm)      forall (k = 1: llm)
118         where (abs(rlat) <= 45.) r_het_interm(:, k) = 0.         where (abs(rlat) <= 45.) r_het_interm(:, k) = 0.
119      end forall      end forall
120      r_het_interm = r_het_interm * (Clx / 3.8e-9)**2      r_het_interm = r_het_interm * (Clx / 3.8e-9)**2
121    
     deallocate(plev) ! pointer  
122      call nf95_close(ncid)      call nf95_close(ncid)
123    
124    end subroutine regr_pr_comb_coefoz    end subroutine regr_pr_comb_coefoz

Legend:
Removed from v.11  
changed lines
  Added in v.18

  ViewVC Help
Powered by ViewVC 1.1.21