CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
rrtmgp_lw_cloud_optics.F90
1
3
10 use mo_rte_kind, only: wl, wp
11 use mo_cloud_optics_rrtmgp, only: ty_cloud_optics => ty_cloud_optics_rrtmgp
12 use rrtmgp_lw_gas_optics, only: lw_gas_props
13 use radiation_tools, only: check_error_msg
14 use netcdf
15 use mpi_f08
16
17 implicit none
18
19 type(ty_cloud_optics) :: lw_cloud_props
20 integer :: &
21 nrghice_fromfilelw, nbandlw, nsize_liqlw, nsize_icelw, nsizereglw, &
22 ncoeff_extlw, ncoeff_ssa_glw, nboundlw, npairslw
23 real(wp), dimension(:,:), allocatable :: &
24 lut_extliqlw, & !< LUT shortwave liquid extinction coefficient
25 lut_ssaliqlw, & !< LUT shortwave liquid single scattering albedo
26 lut_asyliqlw, & !< LUT shortwave liquid asymmetry parameter
27 band_limscldlw
28 real(wp), dimension(:,:,:), allocatable :: &
29 lut_exticelw, & !< LUT shortwave ice extinction coefficient
30 lut_ssaicelw, & !< LUT shortwave ice single scattering albedo
31 lut_asyicelw
32
33 ! Parameters used for rain and snow(+groupel) RRTMGP cloud-optics
34 real(wp), parameter :: &
35 absrain = 0.33e-3, &
36 abssnow0 = 1.5, &
37 abssnow1 = 2.34e-3
38 real(wp) :: &
39 radliq_lwrlw, & !< Liquid particle size lower bound for LUT interpolation
40 radliq_uprlw, & !< Liquid particle size upper bound for LUT interpolation
41 radice_lwrlw, & !< Ice particle size upper bound for LUT interpolation
42 radice_uprlw
43
44contains
45
46 ! ######################################################################################
47 ! SUBROUTINE rrtmgp_lw_cloud_optics_init()
48 ! ######################################################################################
50 subroutine rrtmgp_lw_cloud_optics_init(rrtmgp_root_dir, rrtmgp_lw_file_clouds, &
51 nrghice, mpicomm, mpirank, mpiroot, errmsg, errflg)
52
53 ! Inputs
54 character(len=128),intent(in) :: &
55 rrtmgp_root_dir, & !< RTE-RRTMGP root directory
56 rrtmgp_lw_file_clouds
57
58 integer, intent(inout) :: &
59 nrghice
60 type(mpi_comm), intent(in) :: &
61 mpicomm
62 integer, intent(in) :: &
63 mpirank, & !< Current MPI rank
64 mpiroot
65
66 ! Outputs
67 character(len=*), intent(out) :: &
68 errmsg
69 integer, intent(out) :: &
70 errflg
71
72 ! Local variables
73 integer :: dimID,varID,status,ncid,mpierr
74 character(len=264) :: lw_cloud_props_file
75
76 ! Initialize
77 errmsg = ''
78 errflg = 0
79
80 ! Filenames are set in the physics_nml
81 lw_cloud_props_file = trim(rrtmgp_root_dir)//trim(rrtmgp_lw_file_clouds)
82
83 ! #######################################################################################
84 !
85 ! Read dimensions for longwave cloud-optics fields...
86 ! (ONLY master processor(0), if MPI enabled)
87 !
88 ! #######################################################################################
89 if (mpirank .eq. mpiroot) then
90 write (*,*) 'Reading RRTMGP longwave cloud-optics metadata ... '
91
92 ! Open file
93 status = nf90_open(trim(lw_cloud_props_file), nf90_nowrite, ncid)
94
95 ! Read dimensions
96 status = nf90_inq_dimid(ncid, 'nband', dimid)
97 status = nf90_inquire_dimension(ncid, dimid, len=nbandlw)
98 status = nf90_inq_dimid(ncid, 'nrghice', dimid)
99 status = nf90_inquire_dimension(ncid, dimid, len=nrghice_fromfilelw)
100 status = nf90_inq_dimid(ncid, 'nsize_liq', dimid)
101 status = nf90_inquire_dimension(ncid, dimid, len=nsize_liqlw)
102 status = nf90_inq_dimid(ncid, 'nsize_ice', dimid)
103 status = nf90_inquire_dimension(ncid, dimid, len=nsize_icelw)
104 status = nf90_inq_dimid(ncid, 'nsizereg', dimid)
105 status = nf90_inquire_dimension(ncid, dimid, len=nsizereglw)
106 status = nf90_inq_dimid(ncid, 'ncoeff_ext', dimid)
107 status = nf90_inquire_dimension(ncid, dimid, len=ncoeff_extlw)
108 status = nf90_inq_dimid(ncid, 'ncoeff_ssa_g', dimid)
109 status = nf90_inquire_dimension(ncid, dimid, len=ncoeff_ssa_glw)
110 status = nf90_inq_dimid(ncid, 'nbound', dimid)
111 status = nf90_inquire_dimension(ncid, dimid, len=nboundlw)
112 status = nf90_inq_dimid(ncid, 'pair', dimid)
113 status = nf90_inquire_dimension(ncid, dimid, len=npairslw)
114
115 endif ! On master processor
116
117 ! Other processors waiting...
118 call mpi_barrier(mpicomm, mpierr)
119
120 ! #######################################################################################
121 !
122 ! Broadcast dimensions...
123 ! (ALL processors)
124 !
125 ! #######################################################################################
126 call mpi_bcast(nbandlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
127 call mpi_bcast(nsize_liqlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
128 call mpi_bcast(nsize_icelw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
129 call mpi_bcast(nsizereglw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
130 call mpi_bcast(ncoeff_extlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
131 call mpi_bcast(ncoeff_ssa_glw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
132 call mpi_bcast(nboundlw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
133 call mpi_bcast(npairslw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
134
135 ! Has the number of ice-roughnesses to use been provided from the namelist?
136 ! If so, override nrghice from cloud-optics file
137 if (nrghice .ne. 0) nrghice_fromfilelw = nrghice
138 call mpi_bcast(nrghice_fromfilelw, 1, mpi_integer, mpiroot, mpicomm, mpierr)
139
140 ! #######################################################################################
141 !
142 ! Allocate space for arrays...
143 ! (ALL processors)
144 !
145 ! #######################################################################################
146 if (.not. allocated(lut_extliqlw)) allocate(lut_extliqlw(nsize_liqlw, nbandlw))
147 if (.not. allocated(lut_ssaliqlw)) allocate(lut_ssaliqlw(nsize_liqlw, nbandlw))
148 if (.not. allocated(lut_asyliqlw)) allocate(lut_asyliqlw(nsize_liqlw, nbandlw))
149 if (.not. allocated(lut_exticelw)) allocate(lut_exticelw(nsize_icelw, nbandlw, nrghice_fromfilelw))
150 if (.not. allocated(lut_ssaicelw)) allocate(lut_ssaicelw(nsize_icelw, nbandlw, nrghice_fromfilelw))
151 if (.not. allocated(lut_asyicelw)) allocate(lut_asyicelw(nsize_icelw, nbandlw, nrghice_fromfilelw))
152 if (.not. allocated(band_limscldlw)) allocate(band_limscldlw(2,nbandlw))
153
154 ! #######################################################################################
155 !
156 ! Read in data ...
157 ! (ONLY master processor(0), if MPI enabled)
158 !
159 ! #######################################################################################
160 if (mpirank .eq. mpiroot) then
161 ! Read in fields from file
162 write (*,*) 'Reading RRTMGP longwave cloud data (LUT) ... '
163 status = nf90_inq_varid(ncid,'radliq_lwr',varid)
164 status = nf90_get_var(ncid,varid,radliq_lwrlw)
165 status = nf90_inq_varid(ncid,'radliq_upr',varid)
166 status = nf90_get_var(ncid,varid,radliq_uprlw)
167 status = nf90_inq_varid(ncid,'radice_lwr',varid)
168 status = nf90_get_var(ncid,varid,radice_lwrlw)
169 status = nf90_inq_varid(ncid,'radice_upr',varid)
170 status = nf90_get_var(ncid,varid,radice_uprlw)
171 status = nf90_inq_varid(ncid,'lut_extliq',varid)
172 status = nf90_get_var(ncid,varid,lut_extliqlw)
173 status = nf90_inq_varid(ncid,'lut_ssaliq',varid)
174 status = nf90_get_var(ncid,varid,lut_ssaliqlw)
175 status = nf90_inq_varid(ncid,'lut_asyliq',varid)
176 status = nf90_get_var(ncid,varid,lut_asyliqlw)
177 status = nf90_inq_varid(ncid,'lut_extice',varid)
178 status = nf90_get_var(ncid,varid,lut_exticelw)
179 status = nf90_inq_varid(ncid,'lut_ssaice',varid)
180 status = nf90_get_var(ncid,varid,lut_ssaicelw)
181 status = nf90_inq_varid(ncid,'lut_asyice',varid)
182 status = nf90_get_var(ncid,varid,lut_asyicelw)
183 status = nf90_inq_varid(ncid,'bnd_limits_wavenumber',varid)
184 status = nf90_get_var(ncid,varid,band_limscldlw)
185
186 ! Close file
187 status = nf90_close(ncid)
188 endif ! Master process
189
190 ! Other processors waiting...
191 call mpi_barrier(mpicomm, mpierr)
192
193 ! #######################################################################################
194 !
195 ! Broadcast data...
196 ! (ALL processors)
197 !
198 ! #######################################################################################
199
200 ! Real scalars
201#ifdef RTE_USE_SP
202 call mpi_bcast(radliq_lwrlw, 1, mpi_real, mpiroot, mpicomm, mpierr)
203 call mpi_bcast(radliq_uprlw, 1, mpi_real, mpiroot, mpicomm, mpierr)
204 call mpi_bcast(radice_lwrlw, 1, mpi_real, mpiroot, mpicomm, mpierr)
205 call mpi_bcast(radice_uprlw, 1, mpi_real, mpiroot, mpicomm, mpierr)
206#else
207 call mpi_bcast(radliq_lwrlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
208 call mpi_bcast(radliq_uprlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
209 call mpi_bcast(radice_lwrlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
210 call mpi_bcast(radice_uprlw, 1, mpi_double_precision, mpiroot, mpicomm, mpierr)
211#endif
212
213 ! Real arrays
214#ifdef RTE_USE_SP
215 call mpi_bcast(band_limscldlw, size(band_limscldlw), mpi_real, mpiroot, mpicomm, mpierr)
216 call mpi_bcast(lut_extliqlw, size(lut_extliqlw), mpi_real, mpiroot, mpicomm, mpierr)
217 call mpi_bcast(lut_ssaliqlw, size(lut_ssaliqlw), mpi_real, mpiroot, mpicomm, mpierr)
218 call mpi_bcast(lut_asyliqlw, size(lut_asyliqlw), mpi_real, mpiroot, mpicomm, mpierr)
219 call mpi_bcast(lut_exticelw, size(lut_exticelw), mpi_real, mpiroot, mpicomm, mpierr)
220 call mpi_bcast(lut_ssaicelw, size(lut_ssaicelw), mpi_real, mpiroot, mpicomm, mpierr)
221 call mpi_bcast(lut_asyicelw, size(lut_asyicelw), mpi_real, mpiroot, mpicomm, mpierr)
222#else
223 call mpi_bcast(band_limscldlw, size(band_limscldlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
224 call mpi_bcast(lut_extliqlw, size(lut_extliqlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
225 call mpi_bcast(lut_ssaliqlw, size(lut_ssaliqlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
226 call mpi_bcast(lut_asyliqlw, size(lut_asyliqlw), mpi_double_precision, mpiroot, mpicomm, mpierr)
227 call mpi_bcast(lut_exticelw, size(lut_exticelw), mpi_double_precision, mpiroot, mpicomm, mpierr)
228 call mpi_bcast(lut_ssaicelw, size(lut_ssaicelw), mpi_double_precision, mpiroot, mpicomm, mpierr)
229 call mpi_bcast(lut_asyicelw, size(lut_asyicelw), mpi_double_precision, mpiroot, mpicomm, mpierr)
230#endif
231
232 ! #######################################################################################
233 !
234 ! Initialize RRTMGP DDT's...
235 !
236 ! #######################################################################################
237 call check_error_msg('lw_cloud_optics_init',lw_cloud_props%load(band_limscldlw, &
238 radliq_lwrlw, radliq_uprlw, radice_lwrlw, radice_uprlw, &
239 lut_extliqlw, lut_ssaliqlw, lut_asyliqlw, &
240 lut_exticelw, lut_ssaicelw, lut_asyicelw))
241
242 call check_error_msg('lw_cloud_optics_init',lw_cloud_props%set_ice_roughness(nrghice))
243
244 end subroutine rrtmgp_lw_cloud_optics_init
245end module rrtmgp_lw_cloud_optics
This module contains tools for radiation.
This module contains two routines: The first initializes data and functions needed to compute the lon...
This module contains two routines: One to initialize the k-distribution data and functions needed to ...