source: LMDZ6/branches/Ocean_skin/libf/phylmd/ecrad/radiation_pdf_sampler.F90 @ 5407

Last change on this file since 5407 was 3908, checked in by idelkadi, 4 years ago

Online implementation of the radiative transfer code ECRAD in the LMDZ model.

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
    • interface routine : radiation_scheme.F90
  • Adaptation of compilation scripts :
    • compilation under CPP key CPP_ECRAD
    • compilation with option "-rad ecard" or "-ecard true"
    • The "-rad old/rtm/ecran" build option will need to replace the "-rrtm true" and "-ecrad true" options in the future.
  • Runing LMDZ simulations with ecrad, you need :
    • logical key iflag_rrtm = 2 in physiq.def
    • namelist_ecrad (DefLists?)
    • the directory "data" containing the configuration files is temporarily placed in ../libfphylmd/ecrad/
  • Compilation and execution are tested in the 1D case. The repository under svn would allow to continue the implementation work: tests, verification of the results, ...
File size: 6.6 KB
Line 
1! radiation_pdf_sampler.F90 - Get samples from a lognormal distribution for McICA
2!
3! (C) Copyright 2015- ECMWF.
4!
5! This software is licensed under the terms of the Apache Licence Version 2.0
6! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7!
8! In applying this licence, ECMWF does not waive the privileges and immunities
9! granted to it by virtue of its status as an intergovernmental organisation
10! nor does it submit to any jurisdiction.
11!
12! Author:  Robin Hogan
13! Email:   r.j.hogan@ecmwf.int
14!
15
16module radiation_pdf_sampler
17
18  use parkind1, only : jprb
19
20  implicit none
21  public
22
23  !---------------------------------------------------------------------
24  ! Derived type for sampling from a lognormal distribution, used to
25  ! generate water content or optical depth scalings for use in the
26  ! Monte Carlo Independent Column Approximation (McICA)
27  type pdf_sampler_type
28    ! Number of points in look-up table for cumulative distribution
29    ! function (CDF) and fractional standard deviation (FSD)
30    ! dimensions
31    integer :: ncdf, nfsd
32
33    ! First value of FSD and the reciprocal of the interval between
34    ! FSD values (which are assumed to be uniformly distributed)
35    real(jprb) :: fsd1, inv_fsd_interval
36
37    ! Value of the distribution for each CDF and FSD bin
38    real(jprb), allocatable, dimension(:,:) :: val
39
40  contains
41
42    procedure :: setup => setup_pdf_sampler
43    procedure :: sample => sample_from_pdf
44    procedure :: masked_sample => sample_from_pdf_masked
45    procedure :: deallocate => deallocate_pdf_sampler
46
47  end type pdf_sampler_type
48
49contains
50
51  !---------------------------------------------------------------------
52  ! Load look-up table from a file
53  subroutine setup_pdf_sampler(this, file_name, iverbose)
54   
55    use yomhook,     only : lhook, dr_hook
56    use easy_netcdf, only : netcdf_file
57
58    class(pdf_sampler_type), intent(inout) :: this
59    character(len=*),        intent(in)    :: file_name
60    integer, optional,       intent(in)    :: iverbose
61
62    type(netcdf_file)  :: file
63    integer            :: iverb
64    real(jprb), allocatable :: fsd(:)
65
66    real(jprb)         :: hook_handle
67
68    if (lhook) call dr_hook('radiation_pdf_sampler:setup',0,hook_handle)
69
70    if (present(iverbose)) then
71      iverb = iverbose
72    else
73      iverb = 2
74    end if
75
76    if (allocated(this%val)) then
77      deallocate(this%val)
78    end if
79
80    call file%open(trim(file_name), iverbose=iverb)
81
82    call file%get('fsd',fsd)
83    call file%get('x',  this%val)
84
85    call file%close()
86
87    this%ncdf = size(this%val,1)
88    this%nfsd = size(this%val,2)
89    this%fsd1 = fsd(1)
90    this%inv_fsd_interval = 1.0_jprb / (fsd(2)-fsd(1))
91
92    deallocate(fsd)
93
94    if (lhook) call dr_hook('radiation_pdf_sampler:setup',1,hook_handle)
95
96  end subroutine setup_pdf_sampler
97
98  !---------------------------------------------------------------------
99  ! Deallocate data in pdf_sampler_type derived type
100  subroutine deallocate_pdf_sampler(this)
101
102    use yomhook,     only : lhook, dr_hook
103
104    class(pdf_sampler_type), intent(inout) :: this
105    real(jprb)         :: hook_handle
106
107    if (lhook) call dr_hook('radiation_pdf_sampler:deallocate',0,hook_handle)
108
109    if (allocated(this%val)) then
110      deallocate(this%val)
111    end if
112
113    if (lhook) call dr_hook('radiation_pdf_sampler:deallocate',1,hook_handle)
114   
115  end subroutine deallocate_pdf_sampler
116
117
118  !---------------------------------------------------------------------
119  ! Extract the value of a lognormal distribution with fractional
120  ! standard deviation "fsd" corresponding to the cumulative
121  ! distribution function value "cdf", and return it in val. Since this
122  ! is an elemental subroutine, fsd, cdf and val may be arrays.
123  elemental subroutine sample_from_pdf(this, fsd, cdf, val)
124   
125    class(pdf_sampler_type), intent(in)  :: this
126
127    ! Fractional standard deviation (0 to 4) and cumulative
128    ! distribution function (0 to 1)
129    real(jprb),              intent(in)  :: fsd, cdf
130
131    ! Sample from distribution
132    real(jprb),              intent(out) :: val
133
134    ! Index to look-up table
135    integer    :: ifsd, icdf
136
137    ! Weights in bilinear interpolation
138    real(jprb) :: wfsd, wcdf
139
140    ! Bilinear interpolation with bounds
141    wcdf = cdf * (this%ncdf-1) + 1.0_jprb
142    icdf = max(1, min(int(wcdf), this%ncdf-1))
143    wcdf = max(0.0_jprb, min(wcdf - icdf, 1.0_jprb))
144
145    wfsd = (fsd-this%fsd1) * this%inv_fsd_interval + 1.0_jprb
146    ifsd = max(1, min(int(wfsd), this%nfsd-1))
147    wfsd = max(0.0_jprb, min(wfsd - ifsd, 1.0_jprb))
148
149    val =    (1.0_jprb-wcdf)*(1.0_jprb-wfsd) * this%val(icdf  ,ifsd)   &
150         & + (1.0_jprb-wcdf)*          wfsd  * this%val(icdf  ,ifsd+1) &
151         & +           wcdf *(1.0_jprb-wfsd) * this%val(icdf+1,ifsd)   &
152         & +           wcdf *          wfsd  * this%val(icdf+1,ifsd+1)
153
154  end subroutine sample_from_pdf
155
156
157  !---------------------------------------------------------------------
158  ! For true elements of mask, extract the values of a lognormal
159  ! distribution with fractional standard deviation "fsd"
160  ! corresponding to the cumulative distribution function values
161  ! "cdf", and return in val. For false elements of mask, return zero
162  ! in val.
163  subroutine sample_from_pdf_masked(this, nsamp, fsd, cdf, val, mask)
164   
165    class(pdf_sampler_type), intent(in)  :: this
166
167    ! Number of samples
168    integer,    intent(in) :: nsamp
169
170    ! Fractional standard deviation (0 to 4) and cumulative
171    ! distribution function (0 to 1)
172    real(jprb), intent(in)  :: fsd(nsamp), cdf(nsamp)
173
174    ! Sample from distribution
175    real(jprb), intent(out) :: val(:)
176
177    ! Mask
178    logical,    intent(in) :: mask(nsamp)
179
180    ! Loop index
181    integer    :: jsamp
182
183    ! Index to look-up table
184    integer    :: ifsd, icdf
185
186    ! Weights in bilinear interpolation
187    real(jprb) :: wfsd, wcdf
188
189    do jsamp = 1,nsamp
190      if (mask(jsamp)) then
191        ! Bilinear interpolation with bounds
192        wcdf = cdf(jsamp) * (this%ncdf-1) + 1.0_jprb
193        icdf = max(1, min(int(wcdf), this%ncdf-1))
194        wcdf = max(0.0_jprb, min(wcdf - icdf, 1.0_jprb))
195       
196        wfsd = (fsd(jsamp)-this%fsd1) * this%inv_fsd_interval + 1.0_jprb
197        ifsd = max(1, min(int(wfsd), this%nfsd-1))
198        wfsd = max(0.0_jprb, min(wfsd - ifsd, 1.0_jprb))
199       
200        val(jsamp)=(1.0_jprb-wcdf)*(1.0_jprb-wfsd) * this%val(icdf  ,ifsd)   &
201             &    +(1.0_jprb-wcdf)*          wfsd  * this%val(icdf  ,ifsd+1) &
202             &    +          wcdf *(1.0_jprb-wfsd) * this%val(icdf+1,ifsd)   &
203             &    +          wcdf *          wfsd  * this%val(icdf+1,ifsd+1)
204      else
205        val(jsamp) = 0.0_jprb
206      end if
207    end do
208  end subroutine sample_from_pdf_masked
209
210end module radiation_pdf_sampler
Note: See TracBrowser for help on using the repository browser.