source: LMDZ5/trunk/libf/phylmd/cosp/cosp_stats.F90 @ 2824

Last change on this file since 2824 was 2428, checked in by idelkadi, 9 years ago

Mise a jour du simulateur COSP (passage de la version v3.2 a la version v1.4) :

  • mise a jour des sources pour ISCCP, CALIPSO et PARASOL
  • prise en compte des changements de phases pour les nuages (Calipso)
  • rajout de plusieurs diagnostiques (fraction nuageuse en fonction de la temperature, ...)

http://lmdz.lmd.jussieu.fr/Members/aidelkadi/cosp

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 13.4 KB
Line 
1! (c) British Crown Copyright 2008, the Met Office.
2! All rights reserved.
3! $Revision: 88 $, $Date: 2013-11-13 15:08:38 +0100 (mer. 13 nov. 2013) $
4! $URL: http://cfmip-obs-sim.googlecode.com/svn/stable/v1.4.0/cosp_stats.F90 $
5!
6! Redistribution and use in source and binary forms, with or without modification, are permitted
7! provided that the following conditions are met:
8!
9!     * Redistributions of source code must retain the above copyright notice, this list
10!       of conditions and the following disclaimer.
11!     * Redistributions in binary form must reproduce the above copyright notice, this list
12!       of conditions and the following disclaimer in the documentation and/or other materials
13!       provided with the distribution.
14!     * Neither the name of the Met Office nor the names of its contributors may be used
15!       to endorse or promote products derived from this software without specific prior written
16!       permission.
17!
18! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR
19! IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND
20! FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR
21! CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
22! DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
23! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER
24! IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
25! OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26
27!
28! History:
29! Jul 2007 - A. Bodas-Salcedo - Initial version
30! Jul 2008 - A. Bodas-Salcedo - Added capability of producing outputs in standard grid
31! Oct 2008 - J.-L. Dufresne   - Bug fixed. Assignment of Npoints,Nlevels,Nhydro,Ncolumns in COSP_STATS
32! Oct 2008 - H. Chepfer       - Added PARASOL reflectance arguments
33! Jun 2010 - T. Yokohata, T. Nishimura and K. Ogochi - Added NEC SXs optimisations
34! Jan 2013 - G. Cesana        - Added betaperp and temperature arguments
35!                             - Added phase 3D/3Dtemperature/Map output variables in diag_lidar
36!
37!
38#include "cosp_defs.h"
39MODULE MOD_COSP_STATS
40  USE MOD_COSP_CONSTANTS
41  USE MOD_COSP_TYPES
42  USE MOD_LLNL_STATS
43  USE MOD_LMD_IPSL_STATS
44  IMPLICIT NONE
45
46CONTAINS
47
48!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
49!------------------- SUBROUTINE COSP_STATS ------------------------
50!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
51SUBROUTINE COSP_STATS(gbx,sgx,cfg,sgradar,sglidar,vgrid,stradar,stlidar)
52
53   ! Input arguments
54   type(cosp_gridbox),intent(in) :: gbx
55   type(cosp_subgrid),intent(in) :: sgx
56   type(cosp_config),intent(in)  :: cfg
57   type(cosp_sgradar),intent(in) :: sgradar
58   type(cosp_sglidar),intent(in) :: sglidar
59   type(cosp_vgrid),intent(in)   :: vgrid
60   ! Output arguments
61   type(cosp_radarstats),intent(inout) :: stradar ! Summary statistics for radar
62   type(cosp_lidarstats),intent(inout) :: stlidar ! Summary statistics for lidar
63
64   ! Local variables
65   integer :: Npoints  !# of grid points
66   integer :: Nlevels  !# of levels
67   integer :: Nhydro   !# of hydrometeors
68   integer :: Ncolumns !# of columns
69   integer :: Nlr
70   logical :: ok_lidar_cfad = .false.
71   real,dimension(:,:,:),allocatable :: Ze_out,betatot_out,betamol_in,betamol_out,ph_in,ph_out
72   real,dimension(:,:),allocatable :: ph_c,betamol_c
73   real,dimension(:,:,:),allocatable ::  betaperptot_out, temp_in, temp_out
74   real,dimension(:,:),allocatable :: temp_c
75
76   Npoints  = gbx%Npoints
77   Nlevels  = gbx%Nlevels
78   Nhydro   = gbx%Nhydro
79   Ncolumns = gbx%Ncolumns
80   Nlr      = vgrid%Nlvgrid
81
82   if (cfg%LcfadLidarsr532) ok_lidar_cfad=.true.
83
84   if (vgrid%use_vgrid) then ! Statistics in a different vertical grid
85        allocate(Ze_out(Npoints,Ncolumns,Nlr),betatot_out(Npoints,Ncolumns,Nlr), &
86                 betamol_in(Npoints,1,Nlevels),betamol_out(Npoints,1,Nlr),betamol_c(Npoints,Nlr), &
87                 ph_in(Npoints,1,Nlevels),ph_out(Npoints,1,Nlr),ph_c(Npoints,Nlr))
88        Ze_out = 0.0
89        betatot_out  = 0.0
90        betamol_out= 0.0
91        betamol_c  = 0.0
92        ph_in(:,1,:)  = gbx%ph(:,:)
93        ph_out  = 0.0
94        ph_c    = 0.0
95        allocate(betaperptot_out(Npoints,Ncolumns,Nlr),temp_in(Npoints,1,Nlevels),temp_out(Npoints,1,Nlr), &
96                 temp_c(Npoints,Nlr))
97        betaperptot_out = 0.0
98        temp_in = 0.0
99        temp_out = 0.0
100        temp_c = 0.0
101
102        !++++++++++++ Radar CFAD ++++++++++++++++
103        if (cfg%Lradar_sim) then
104            call cosp_change_vertical_grid(Npoints,Ncolumns,Nlevels,gbx%zlev,gbx%zlev_half,sgradar%Ze_tot, &
105                                           Nlr,vgrid%zl,vgrid%zu,Ze_out,log_units=.true.)
106            stradar%cfad_ze = cosp_cfad(Npoints,Ncolumns,Nlr,DBZE_BINS,Ze_out, &
107                                        DBZE_MIN,DBZE_MAX,CFAD_ZE_MIN,CFAD_ZE_WIDTH)
108        endif
109        !++++++++++++ Lidar CFAD ++++++++++++++++
110        if (cfg%Llidar_sim) then
111            betamol_in(:,1,:) = sglidar%beta_mol(:,:)
112            call cosp_change_vertical_grid(Npoints,1,Nlevels,gbx%zlev,gbx%zlev_half,betamol_in, &
113                                           Nlr,vgrid%zl,vgrid%zu,betamol_out)
114            call cosp_change_vertical_grid(Npoints,Ncolumns,Nlevels,gbx%zlev,gbx%zlev_half,sglidar%beta_tot, &
115                                           Nlr,vgrid%zl,vgrid%zu,betatot_out)
116
117            temp_in(:,1,:) = gbx%T(:,:)
118            call cosp_change_vertical_grid(Npoints,Ncolumns,Nlevels,gbx%zlev,gbx%zlev_half,sglidar%betaperp_tot, &
119                                           Nlr,vgrid%zl,vgrid%zu,betaperptot_out)
120            call cosp_change_vertical_grid(Npoints,1,Nlevels,gbx%zlev,gbx%zlev_half,temp_in, &
121                                           Nlr,vgrid%zl,vgrid%zu,temp_out)
122            temp_c(:,:) = temp_out(:,1,:)
123
124            call cosp_change_vertical_grid(Npoints,1,Nlevels,gbx%zlev,gbx%zlev_half,ph_in, &
125                                           Nlr,vgrid%zl,vgrid%zu,ph_out)
126            ph_c(:,:) = ph_out(:,1,:)
127            betamol_c(:,:) = betamol_out(:,1,:)
128            ! Stats from lidar_stat_summary
129            call diag_lidar(Npoints,Ncolumns,Nlr,SR_BINS,PARASOL_NREFL &
130                            ,temp_c,betatot_out,betaperptot_out,betamol_c,sglidar%refl,gbx%land,ph_c &
131                            ,LIDAR_UNDEF,ok_lidar_cfad &
132                            ,stlidar%cfad_sr,stlidar%srbval &
133                            ,LIDAR_NCAT,stlidar%lidarcld,stlidar%lidarcldphase &
134                            ,stlidar%cldlayer,stlidar%cldlayerphase,stlidar%lidarcldtmp &
135                            ,stlidar%parasolrefl)
136        endif
137
138        !++++++++++++ Lidar-only cloud amount and lidar&radar total cloud mount ++++++++++++++++
139        if (cfg%Lradar_sim.and.cfg%Llidar_sim) call cosp_lidar_only_cloud(Npoints,Ncolumns,Nlr, &
140                                    temp_c,betatot_out,betaperptot_out,betamol_c,Ze_out, &
141                                    stradar%lidar_only_freq_cloud,stradar%radar_lidar_tcc)
142        deallocate(temp_out,temp_c,betaperptot_out)
143
144        ! Deallocate arrays at coarse resolution
145        deallocate(Ze_out,betatot_out,betamol_in,betamol_out,betamol_c,ph_in,ph_out,ph_c)
146   else ! Statistics in model levels
147        !++++++++++++ Radar CFAD ++++++++++++++++
148        if (cfg%Lradar_sim) stradar%cfad_ze = cosp_cfad(Npoints,Ncolumns,Nlr,DBZE_BINS,sgradar%Ze_tot, &
149                                        DBZE_MIN,DBZE_MAX,CFAD_ZE_MIN,CFAD_ZE_WIDTH)
150        !++++++++++++ Lidar CFAD ++++++++++++++++
151        ! Stats from lidar_stat_summary
152        if (cfg%Llidar_sim) call diag_lidar(Npoints,Ncolumns,Nlr,SR_BINS,PARASOL_NREFL &
153                        ,sglidar%temp_tot,sglidar%beta_tot,sglidar%betaperp_tot,sglidar%beta_mol,sglidar%refl,gbx%land,gbx%ph &
154                        ,LIDAR_UNDEF,ok_lidar_cfad &
155                        ,stlidar%cfad_sr,stlidar%srbval &
156                        ,LIDAR_NCAT,stlidar%lidarcld,stlidar%lidarcldphase,stlidar%cldlayer,stlidar%cldlayerphase &
157                        ,stlidar%lidarcldtmp,stlidar%parasolrefl)
158        !++++++++++++ Lidar-only cloud amount and lidar&radar total cloud mount ++++++++++++++++
159        if (cfg%Lradar_sim.and.cfg%Llidar_sim) call cosp_lidar_only_cloud(Npoints,Ncolumns,Nlr, &
160                                    sglidar%temp_tot,sglidar%beta_tot,sglidar%betaperp_tot,sglidar%beta_mol,sgradar%Ze_tot, &
161                                    stradar%lidar_only_freq_cloud,stradar%radar_lidar_tcc)
162   endif
163   ! Replace undef
164   where (stlidar%cfad_sr   == LIDAR_UNDEF) stlidar%cfad_sr   = R_UNDEF
165   where (stlidar%lidarcld  == LIDAR_UNDEF) stlidar%lidarcld  = R_UNDEF
166   where (stlidar%cldlayer  == LIDAR_UNDEF) stlidar%cldlayer  = R_UNDEF
167   where (stlidar%parasolrefl == LIDAR_UNDEF) stlidar%parasolrefl = R_UNDEF
168   where (stlidar%cldlayerphase  == LIDAR_UNDEF) stlidar%cldlayerphase  = R_UNDEF
169   where (stlidar%lidarcldphase  == LIDAR_UNDEF) stlidar%lidarcldphase  = R_UNDEF
170   where (stlidar%lidarcldtmp  == LIDAR_UNDEF) stlidar%lidarcldtmp  = R_UNDEF
171
172END SUBROUTINE COSP_STATS
173
174!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
175!---------- SUBROUTINE COSP_CHANGE_VERTICAL_GRID ----------------
176!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
177SUBROUTINE COSP_CHANGE_VERTICAL_GRID(Npoints,Ncolumns,Nlevels,zfull,zhalf,y,Nglevels,newgrid_bot,newgrid_top,r,log_units)
178   implicit none
179   ! Input arguments
180   integer,intent(in) :: Npoints  !# of grid points
181   integer,intent(in) :: Nlevels  !# of levels
182   integer,intent(in) :: Ncolumns !# of columns
183   real,dimension(Npoints,Nlevels),intent(in) :: zfull ! Height at model levels [m] (Bottom of model layer)
184   real,dimension(Npoints,Nlevels),intent(in) :: zhalf ! Height at half model levels [m] (Bottom of model layer)
185   real,dimension(Npoints,Ncolumns,Nlevels),intent(in) :: y     ! Variable to be changed to a different grid
186   integer,intent(in) :: Nglevels  !# levels in the new grid
187   real,dimension(Nglevels),intent(in) :: newgrid_bot ! Lower boundary of new levels  [m]
188   real,dimension(Nglevels),intent(in) :: newgrid_top ! Upper boundary of new levels  [m]
189   logical,optional,intent(in) :: log_units ! log units, need to convert to linear units
190   ! Output
191   real,dimension(Npoints,Ncolumns,Nglevels),intent(out) :: r ! Variable on new grid
192
193   ! Local variables
194   integer :: i,j,k
195   logical :: lunits
196   integer :: l
197   real :: w ! Weight
198   real :: dbb, dtb, dbt, dtt ! Distances between edges of both grids
199   integer :: Nw  ! Number of weights
200   real :: wt  ! Sum of weights
201   real,dimension(Nlevels) :: oldgrid_bot,oldgrid_top ! Lower and upper boundaries of model grid
202   real :: yp ! Local copy of y at a particular point.
203              ! This allows for change of units.
204
205   lunits=.false.
206   if (present(log_units)) lunits=log_units
207
208   r = 0.0
209
210   do i=1,Npoints
211     ! Calculate tops and bottoms of new and old grids
212     oldgrid_bot = zhalf(i,:)
213     oldgrid_top(1:Nlevels-1) = oldgrid_bot(2:Nlevels)
214     oldgrid_top(Nlevels) = zfull(i,Nlevels) +  zfull(i,Nlevels) - zhalf(i,Nlevels) ! Top level symmetric
215     l = 0 ! Index of level in the old grid
216     ! Loop over levels in the new grid
217     do k = 1,Nglevels
218       Nw = 0 ! Number of weigths
219       wt = 0.0 ! Sum of weights
220       ! Loop over levels in the old grid and accumulate total for weighted average
221       do
222         l = l + 1
223         w = 0.0 ! Initialise weight to 0
224         ! Distances between edges of both grids
225         dbb = oldgrid_bot(l) - newgrid_bot(k)
226         dtb = oldgrid_top(l) - newgrid_bot(k)
227         dbt = oldgrid_bot(l) - newgrid_top(k)
228         dtt = oldgrid_top(l) - newgrid_top(k)
229         if (dbt >= 0.0) exit ! Do next level in the new grid
230         if (dtb > 0.0) then
231           if (dbb <= 0.0) then
232             if (dtt <= 0) then
233               w = dtb
234             else
235               w = newgrid_top(k) - newgrid_bot(k)
236             endif
237           else
238             if (dtt <= 0) then
239               w = oldgrid_top(l) - oldgrid_bot(l)
240             else
241               w = -dbt
242             endif
243           endif
244           ! If layers overlap (w/=0), then accumulate
245           if (w /= 0.0) then
246             Nw = Nw + 1
247             wt = wt + w
248             do j=1,Ncolumns
249               if (lunits) then
250                 if (y(i,j,l) /= R_UNDEF) then
251                   yp = 10.0**(y(i,j,l)/10.0)
252                 else
253                   yp = 0.0
254                 endif
255               else
256                 yp = y(i,j,l)
257               endif
258               r(i,j,k) = r(i,j,k) + w*yp
259             enddo
260           endif
261         endif
262       enddo
263       l = l - 2
264       if (l < 1) l = 0
265       ! Calculate average in new grid
266       if (Nw > 0) then
267         do j=1,Ncolumns
268           r(i,j,k) = r(i,j,k)/wt
269         enddo
270       endif
271     enddo
272   enddo
273
274   ! Set points under surface to R_UNDEF, and change to dBZ if necessary
275   do k=1,Nglevels
276     do j=1,Ncolumns
277       do i=1,Npoints
278         if (newgrid_top(k) > zhalf(i,1)) then ! Level above model bottom level
279           if (lunits) then
280             if (r(i,j,k) <= 0.0) then
281               r(i,j,k) = R_UNDEF
282             else
283               r(i,j,k) = 10.0*log10(r(i,j,k))
284             endif
285           endif
286         else ! Level below surface
287           r(i,j,k) = R_GROUND
288         endif
289       enddo
290     enddo
291   enddo
292
293END SUBROUTINE COSP_CHANGE_VERTICAL_GRID
294
295END MODULE MOD_COSP_STATS
Note: See TracBrowser for help on using the repository browser.