/[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/Mobidic/regr_pr_coefoz.f90 revision 12 by guez, Mon Jul 21 16:05:07 2008 UTC trunk/Sources/phylmd/Mobidic/regr_pr_av.f revision 162 by guez, Fri Jul 24 16:54:30 2015 UTC
# Line 1  Line 1 
1  module regr_pr_coefoz  module regr_pr_av_m
2    
3    implicit none    implicit none
4    
5  contains  contains
6    
7    subroutine regr_pr_av_coefoz(ncid, name, julien, press_in_edg, v3)    subroutine regr_pr_av(ncid, name, julien, paprs, v3)
8    
9      ! "regr_pr_av_coefoz" stands for "regrid pressure averaging      ! "regr_pr_av" stands for "regrid pressure averaging".
     ! coefficient ozone".  
     ! This procedure reads a single Mobidic ozone coefficient from  
     ! "coefoz_LMDZ.nc", at the current day, regrids this parameter in  
     ! pressure to the LMDZ vertical grid and packs it to the LMDZ  
     ! horizontal "physics" grid.  
     ! Regridding in pressure is done by averaging a step function.  
10    
11      use dimens_m, only: iim, jjm, llm      ! This procedure reads a 2D latitude-pressure field from a NetCDF
12      use dimphy, only: klon      ! file, at a given day, regrids this field in pressure to the LMDZ
13      use netcdf95, only: nf95_inq_varid, handle_err      ! vertical grid and packs it to the LMDZ horizontal "physics"
14      use netcdf, only: nf90_get_var      ! grid.
     use grid_change, only: dyn_phy  
     use regr_pr, only: regr_pr_av  
     use nrutil, only: assert  
   
     integer, intent(in):: ncid ! NetCDF ID of the file  
     character(len=*), intent(in):: name ! of the NetCDF variable  
     integer, intent(in):: julien ! jour julien, 1 <= julien <= 360  
   
     real, intent(in):: press_in_edg(:)  
     ! (edges of pressure intervals for Mobidic data, in Pa, in  
     ! strictly increasing order)  
   
     real, intent(out):: v3(:, :) ! (klon, llm)  
     ! (ozone coefficient from Mobidic on the "physics" grid)  
     ! ("v3(i, k)" is at longitude "xlon(i)", latitude  
     ! "xlat(i)", middle of layer "k".)  
   
     ! Variables local to the procedure:  
     integer varid, ncerr  
     integer k  
   
     real  v1(jjm + 1, size(press_in_edg) - 1)  
     ! (ozone coefficient from "coefoz_LMDZ.nc" at day "julien")  
     ! ("v1(j, k)" is at latitude "rlatu(j)" and for  
     ! pressure interval "[press_in_edg(k), press_in_edg(k+1)]".)  
   
     real v2(iim + 1, jjm + 1, llm)  
     ! (ozone parameter from Mobidic on the "dynamics" grid)  
     ! "v2(i, j, k)" is at longitude "rlonv(i)", latitude  
     ! "rlatu(j)", middle of layer "k".)  
   
     !--------------------------------------------  
15    
16      call assert(shape(v3) == (/klon, llm/), "regr_pr_av_coefoz")      ! We assume that, in the input file, the field has 3 dimensions:
17        ! latitude, pressure, julian day.
18    
19      call nf95_inq_varid(ncid, name, varid)      ! We assume that the input field is already on the LMDZ "rlatu"
20        ! latitudes, except that latitudes are in ascending order in the
21      ! Get data at the right day from the input file:      ! input file.
     ncerr = nf90_get_var(ncid, varid, v1, start=(/1, 1, julien/))  
     call handle_err("regr_pr_av_coefoz nf90_get_var " // name, ncerr, ncid)  
     ! Latitudes are in increasing order in the input file while  
     ! "rlatu" is in decreasing order so we need to invert order:  
     v1 = v1(jjm+1:1:-1, :)  
   
     ! Regrid in pressure at each horizontal position:  
     v2 = regr_pr_av(v1, press_in_edg)  
   
     forall (k = 1:llm) v3(:, k) = pack(v2(:, :, k), dyn_phy)  
   
   end subroutine regr_pr_av_coefoz  
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.
   subroutine regr_pr_int_coefoz(ncid, name, julien, plev, top_value, v3)  
   
     ! This procedure reads a single Mobidic ozone coefficient from  
     !"coefoz_LMDZ.nc", at the current day, regrids this parameter in  
     ! pressure to the LMDZ vertical grid and packs it to the LMDZ  
     ! horizontal "physics" grid.  
     ! Regridding is by linear interpolation.  
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
     use netcdf95, only: nf95_inq_varid, handle_err  
     use netcdf, only: nf90_get_var  
28      use grid_change, only: dyn_phy      use grid_change, only: dyn_phy
29      use regr_pr, only: regr_pr_int      use netcdf95, only: nf95_inq_varid, nf95_get_var
30      use nrutil, only: assert      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      integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
37    
38      real, intent(in):: plev(:)      real, intent(in):: paprs(:, :) ! (klon, llm + 1)
39      ! (pressure levels of Mobidic data, in Pa, in strictly increasing order)      ! (pression pour chaque inter-couche, en Pa)
   
     real, intent(in):: top_value  
     ! (extra value of ozone coefficient at 0 pressure)  
40    
41      real, intent(out):: v3(:, :) ! (klon, llm)      real, intent(out):: v3(:, :) ! (klon, llm)
42      ! (ozone parameter from Mobidic on the "physics" grid)      ! regridded field on the partial "physics" grid
43      ! ("v3(i, k)" is at longitude "xlon(i)", latitude      ! "v3(i, k)" is at longitude "xlon(i)", latitude "xlat(i)", in
44      ! "xlat(i)", middle of layer "k".)      ! layer "k".
45    
46      ! Variables local to the procedure:      ! Variables local to the procedure:
     integer varid, ncerr  
     integer k  
47    
48      real  v1(jjm + 1, 0:size(plev))      integer varid ! for NetCDF
49      ! (ozone coefficient from "coefoz_LMDZ.nc" at day "julien")  
50      ! ("v1(j, k >=1)" is at latitude "rlatu(j)" and pressure "plev(k)".)      real v1(jjm + 1, size(press_in_edg) - 1)
51        ! input field at day "julien"
52      real v2(iim + 1, jjm + 1, llm)      ! "v1(j, k)" is at latitude "rlatu(j)" and for
53      ! (ozone parameter from Mobidic on the "dynamics" grid)      ! pressure interval "[press_in_edg(k), press_in_edg(k+1)]".
54      ! "v2(i, j, k)" is at longitude "rlonv(i)", latitude  
55      ! "rlatu(j)", middle of layer "k".)      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, k
61    
62      !--------------------------------------------      !--------------------------------------------
63    
64      call assert(shape(v3) == (/klon, llm/), "regr_pr_int_coefoz")      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)
68    
69      ! Get data at the right day from the input file:      ! Get data at the right day from the input file:
70      ncerr = nf90_get_var(ncid, varid, v1(:, 1:), start=(/1, 1, julien/))      call nf95_get_var(ncid, varid, v1, start=(/1, 1, julien/))
71      call handle_err("regr_pr_int_coefoz nf90_get_var " // name, ncerr, ncid)      ! Latitudes are in ascending order in the input file while
72      ! Latitudes are in increasing order in the input file while      ! "rlatu" is in descending order so we need to invert order:
73      ! "rlatu" is in decreasing order so we need to invert order:      v1 = v1(jjm+1:1:-1, :)
     v1(:, 1:) = v1(jjm+1:1:-1, 1:)  
74    
75      ! Complete "v1" with the value at 0 pressure:      forall (k = 1:size(press_in_edg) - 1) v2(:, k) = pack(spread(v1(:, k), &
76      v1(:, 0) = top_value           dim = 1, ncopies = iim + 1), dyn_phy)
77    
78      ! Regrid in pressure at each horizontal position:      ! Regrid in pressure at each horizontal position:
79      v2 = regr_pr_int(v1, (/0., plev/))      do i = 1, klon
80           v3(i, llm:1:-1) = regr1_step_av(v2(i, :), press_in_edg, &
81      forall (k = 1:llm) v3(:, k) = pack(v2(:, :, k), dyn_phy)              paprs(i, llm+1:1:-1))
82           ! (invert order of indices because "paprs" is in descending order)
83        end do
84    
85    end subroutine regr_pr_int_coefoz    end subroutine regr_pr_av
86    
87  end module regr_pr_coefoz  end module regr_pr_av_m

Legend:
Removed from v.12  
changed lines
  Added in v.162

  ViewVC Help
Powered by ViewVC 1.1.21