/[lmdze]/trunk/Sources/phylmd/Mobidic/regr_pr_av.f
ViewVC logotype

Diff of /trunk/Sources/phylmd/Mobidic/regr_pr_av.f

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

trunk/libf/phylmd/o3_mob_ph_m.f90 revision 3 by guez, Wed Feb 27 13:16:39 2008 UTC trunk/Sources/phylmd/Mobidic/regr_pr_av.f revision 168 by guez, Wed Sep 9 10:41:47 2015 UTC
# Line 1  Line 1 
1  module o3_Mob_ph_m  module regr_pr_av_m
2    
3    implicit none    implicit none
4    
5  contains  contains
6    
7    function o3_Mob_ph(ncid, name)    subroutine regr_pr_av(ncid, name, julien, paprs, v3)
8    
9      ! This function reads a single Mobidic ozone parameter from a file and      ! "regr_pr_av" stands for "regrid pressure averaging".
10      ! packs it on the "physics" grid.  
11        ! This procedure reads a 2D latitude-pressure field from a NetCDF
12        ! file, at a given day, regrids this field in pressure to the LMDZ
13        ! vertical grid and packs it to the LMDZ horizontal "physics"
14        ! grid.
15    
16        ! We assume that, in the input file, the field has 3 dimensions:
17        ! latitude, pressure, julian day.
18    
19        ! We assume that the input field is already on the LMDZ "rlatu"
20        ! latitudes, except that latitudes are in ascending order in the
21        ! input file.
22    
23        ! The target vertical LMDZ grid is the grid of layer boundaries.
24        ! Regridding in pressure is done by averaging a step function of pressure.
25    
26      use dimens_m, only: iim, jjm, llm      use dimens_m, only: iim, jjm, llm
27      use dimphy, only: klon      use dimphy, only: klon
28      use netcdf95, only: nf95_inq_varid, nf90_get_var, handle_err      use grid_change, only: gr_dyn_phy
29      use grid_change, only: dyn_phy      use netcdf95, only: nf95_inq_varid, nf95_get_var
30        use nr_util, only: assert
31        use numer_rec_95, only: regr1_step_av
32        use press_coefoz_m, only: press_in_edg
33    
34      integer, intent(in):: ncid ! NetCDF ID of the file      integer, intent(in):: ncid ! NetCDF ID of the file
35      character(len=*), intent(in):: name ! of the NetCDF variable      character(len=*), intent(in):: name ! of the NetCDF variable
36        integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
37    
38      real o3_Mob_ph(klon, llm, 12)      real, intent(in):: paprs(:, :) ! (klon, llm + 1)
39      ! (ozone parameter from Mobidic on the "physics" grid)      ! (pression pour chaque inter-couche, en Pa)
40      ! (Third dimension is the number of the month in the year.  
41      ! "o3_Mob_ph(i, k, month)" is at longitude "xlon(i)", latitude      real, intent(out):: v3(:, :) ! (klon, llm)
42      ! "xlat(i)", middle of layer "k".)      ! regridded field on the partial "physics" grid
43        ! "v3(i, k)" is at longitude "xlon(i)", latitude "xlat(i)", in
44        ! layer "k".
45    
46      ! Variables local to the procedure:      ! Variables local to the procedure:
     integer varid, ncerr  
     integer k, month  
47    
48      real o3_Mob_dyn(iim + 1, jjm + 1, llm, 12)      integer varid ! for NetCDF
49      ! (ozone parameter from Mobidic on the "dynamics" grid)  
50      ! Fourth dimension is the number of the month in the year.      real v1(jjm + 1, size(press_in_edg) - 1)
51      ! "o3_Mob_dyn(i, j, k, month)" is at longitude "rlonv(i)", latitude      ! input field at day "julien"
52      ! "rlatu(j)", middle of layer "k".)      ! "v1(j, k)" is at latitude "rlatu(j)" and for
53        ! pressure interval "[press_in_edg(k), press_in_edg(k+1)]".
54    
55        real v2(klon, size(press_in_edg) - 1)
56        ! Field on the "physics" horizontal grid. "v2(i, k)" is at
57        ! longitude "xlon(i)", latitude "xlat(i)" and for pressure
58        ! interval "[press_in_edg(k), press_in_edg(k+1)]".)
59    
60        integer i
61    
62      !--------------------------------------------      !--------------------------------------------
63    
64        call assert(shape(v3) == (/klon, llm/), "regr_pr_av klon llm")
65        call assert(shape(paprs) == (/klon, llm+1/), "regr_pr_av paprs")
66    
67      call nf95_inq_varid(ncid, name, varid)      call nf95_inq_varid(ncid, name, varid)
     ncerr = nf90_get_var(ncid, varid, o3_Mob_dyn)  
     call handle_err("o3_Mob_ph nf90_get_var " // name, ncerr, ncid)  
68    
69      ! Latitudes are in increasing order in the input file while      ! Get data at the right day from the input file:
70      ! "rlatu" is in decreasing order, so invert:      call nf95_get_var(ncid, varid, v1, start=(/1, 1, julien/))
71      o3_Mob_dyn = o3_Mob_dyn(:, jjm+1:1:-1, :, :)      ! Latitudes are in ascending order in the input file while
72      forall (k = 1:llm, month = 1:12) &      ! "rlatu" is in descending order so we need to invert order:
73           o3_Mob_ph(:, k, month) = pack(o3_Mob_dyn(:, :, k, month), dyn_phy)      v1 = v1(jjm+1:1:-1, :)
74    
75        v2 = gr_dyn_phy(spread(v1, dim = 1, ncopies = iim + 1))
76    
77        ! Regrid in pressure at each horizontal position:
78        do i = 1, klon
79           v3(i, llm:1:-1) = regr1_step_av(v2(i, :), press_in_edg, &
80                paprs(i, llm+1:1:-1))
81           ! (invert order of indices because "paprs" is in descending order)
82        end do
83    
84    end function o3_Mob_ph    end subroutine regr_pr_av
85    
86  end module o3_Mob_ph_m  end module regr_pr_av_m

Legend:
Removed from v.3  
changed lines
  Added in v.168

  ViewVC Help
Powered by ViewVC 1.1.21