SCALE-RM
Functions/Subroutines
mod_da_driver Module Reference

module Data Assimilation driver More...

Functions/Subroutines

subroutine, public da_driver_setup
 Setup. More...
 
subroutine, public da_driver_finalize
 finalize More...
 
subroutine, public da_driver_update
 Data Assimilation. More...
 

Detailed Description

module Data Assimilation driver

Description
Data Assimilation driver
Author
Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
{PREFIX}_DENS DENS for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m3 DENS
{PREFIX}_MOMX MOMX for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s MOMX
{PREFIX}_MOMY MOMY for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s MOMY
{PREFIX}_MOMZ MOMZ for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s MOMZ
{PREFIX}_PRES PRES for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
Pa PRES
{PREFIX}_QC QC for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg QC
{PREFIX}_QG QG for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg QG
{PREFIX}_QI QI for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg QI
{PREFIX}_QR QR for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg QR
{PREFIX}_QS QS for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg QS
{PREFIX}_QV QV for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg QV
{PREFIX}_RH RH for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
% RH
{PREFIX}_RHOT RHOT for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
K*kg/m3 RHOT
{PREFIX}_T TEMP for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
K TEMP
{PREFIX}_U U for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s U
{PREFIX}_V V for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s V
{PREFIX}_W W for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s W
{PREFIX}_mean_DENS Ensemble mean of DENS for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m3 ENS_mean_DENS
{PREFIX}_mean_MOMX Ensemble mean of MOMX for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s ENS_mean_MOMX
{PREFIX}_mean_MOMY Ensemble mean of MOMY for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s ENS_mean_MOMY
{PREFIX}_mean_MOMZ Ensemble mean of MOMZ for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s ENS_mean_MOMZ
{PREFIX}_mean_PRES Ensemble mean of PRES for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
Pa ENS_mean_PRES
{PREFIX}_mean_QC Ensemble mean of QC for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_mean_QC
{PREFIX}_mean_QG Ensemble mean of QG for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_mean_QG
{PREFIX}_mean_QI Ensemble mean of QI for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_mean_QI
{PREFIX}_mean_QR Ensemble mean of QR for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_mean_QR
{PREFIX}_mean_QS Ensemble mean of QS for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_mean_QS
{PREFIX}_mean_QV Ensemble mean of QV for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_mean_QV
{PREFIX}_mean_RH Ensemble mean of RH for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
% ENS_mean_RH
{PREFIX}_mean_RHOT Ensemble mean of RHOT for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
K*kg/m3 ENS_mean_RHOT
{PREFIX}_mean_T Ensemble mean of TEMP for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
K ENS_mean_TEMP
{PREFIX}_mean_U Ensemble mean of U for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s ENS_mean_U
{PREFIX}_mean_V Ensemble mean of V for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s ENS_mean_V
{PREFIX}_mean_W Ensemble mean of W for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s ENS_mean_W
{PREFIX}_sprd_DENS Ensemble spread of DENS for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m3 ENS_sprd_DENS
{PREFIX}_sprd_MOMX Ensemble spread of MOMX for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s ENS_sprd_MOMX
{PREFIX}_sprd_MOMY Ensemble spread of MOMY for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s ENS_sprd_MOMY
{PREFIX}_sprd_MOMZ Ensemble spread of MOMZ for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/m2/s ENS_sprd_MOMZ
{PREFIX}_sprd_PRES Ensemble spread of PRES for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
Pa ENS_sprd_PRES
{PREFIX}_sprd_QC Ensemble spread of QC for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_sprd_QC
{PREFIX}_sprd_QG Ensemble spread of QG for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_sprd_QG
{PREFIX}_sprd_QI Ensemble spread of QI for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_sprd_QI
{PREFIX}_sprd_QR Ensemble spread of QR for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_sprd_QR
{PREFIX}_sprd_QS Ensemble spread of QS for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_sprd_QS
{PREFIX}_sprd_QV Ensemble spread of QV for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
kg/kg ENS_sprd_QV
{PREFIX}_sprd_RH Ensemble spread of RH for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
% ENS_sprd_RH
{PREFIX}_sprd_RHOT Ensemble spread of RHOT for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
K*kg/m3 ENS_sprd_RHOT
{PREFIX}_sprd_T Ensemble spread of TEMP for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
K ENS_sprd_TEMP
{PREFIX}_sprd_U Ensemble spread of U for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s ENS_sprd_U
{PREFIX}_sprd_V Ensemble spread of V for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s ENS_sprd_V
{PREFIX}_sprd_W Ensemble spread of W for {PREFIX};
{PREFIX} depends on the physics schemes, e.g., QV, QC, QR.
m/s ENS_sprd_W

Function/Subroutine Documentation

◆ da_driver_setup()

subroutine, public mod_da_driver::da_driver_setup

Setup.

Definition at line 100 of file mod_da_driver.F90.

100  use mpi
101  use scale_prc, only: &
103  prc_nprocs, &
104  prc_myrank, &
105  prc_abort
106  use scale_prc_cartesc, only: &
107  prc_num_x, &
108  prc_num_y
109  use scale_atmos_grid_cartesc, only: &
110  dx, &
111  dy
112  use scale_atmos_grid_cartesc_index, only: &
113  kmax, &
114  imax, &
115  jmax, &
116  ka, ks, ke, &
117  ia, is, ie, &
118  ja, js, je, &
119  khalo, &
120  ihalo, &
121  jhalo
122  use scale_atmos_hydrometeor, only: &
123  n_hyd
124  use scale_topography, only: &
126  use scale_comm_ensemble, only: &
127  ensemble_world => comm_ensemble_world, &
128  ensemble_nprocs => comm_ensemble_nprocs, &
129  ensemble_myrank => comm_ensemble_myrank, &
131  use scale_letkf, only: &
133  use mod_da_vars, only: &
134  obs_in_num
135  use mod_da_param_estimation, only: &
137  implicit none
138  !---------------------------------------------------------------------------
139 
140  log_newline
141  log_info("DA_driver_setup",*) 'Setup'
142 
143  if ( rp == sp ) then
144  datatype = mpi_real
145  else if( rp == dp ) then
146  datatype = mpi_double_precision
147  else
148  log_error("DA_driver_setup",*) 'The precision has not been implemented yet:', rp
149  call prc_abort
150  endif
151 
153 
154  log_info("DA_driver_setup",*) 'ENSEMBLE myrank/nprocs:', ensemble_myrank, ensemble_nprocs
155 
156  allocate( qhyd( ka, ia, ja, n_hyd ) )
157  qhyd(:,:,:,:) = 0.0_rp
158 
159  allocate( rh( ka, ia, ja ) )
160  allocate( qdry( ka, ia, ja ) )
161  allocate( rtot( ka, ia, ja ) )
162  allocate( cvtot( ka, ia, ja ) )
163  allocate( cptot( ka, ia, ja ) )
164  rh(:,:,:) = 0.0_rp
165  qdry(:,:,:) = 0.0_rp
166  rtot(:,:,:) = 0.0_rp
167  cvtot(:,:,:) = 0.0_rp
168  cptot(:,:,:) = 0.0_rp
169 
170  allocate( prec( ia, ja ) )
171  prec(:,:) = 0.0_rp
172 
173  allocate( ens_mean_dens(ka,ia,ja) )
174  allocate( ens_mean_momx(ka,ia,ja) )
175  allocate( ens_mean_momy(ka,ia,ja) )
176  allocate( ens_mean_momz(ka,ia,ja) )
177  allocate( ens_mean_rhot(ka,ia,ja) )
178  allocate( ens_mean_u(ka,ia,ja) )
179  allocate( ens_mean_v(ka,ia,ja) )
180  allocate( ens_mean_w(ka,ia,ja) )
181  allocate( ens_mean_temp(ka,ia,ja) )
182  allocate( ens_mean_pres(ka,ia,ja) )
183  allocate( ens_mean_qv(ka,ia,ja) )
184  allocate( ens_mean_qc(ka,ia,ja) )
185  allocate( ens_mean_qr(ka,ia,ja) )
186  allocate( ens_mean_qi(ka,ia,ja) )
187  allocate( ens_mean_qs(ka,ia,ja) )
188  allocate( ens_mean_qg(ka,ia,ja) )
189  allocate( ens_mean_rh(ka,ia,ja) )
190  ens_mean_dens(:,:,:) = 0.0_rp
191  ens_mean_momx(:,:,:) = 0.0_rp
192  ens_mean_momy(:,:,:) = 0.0_rp
193  ens_mean_momz(:,:,:) = 0.0_rp
194  ens_mean_rhot(:,:,:) = 0.0_rp
195  ens_mean_u(:,:,:) = 0.0_rp
196  ens_mean_v(:,:,:) = 0.0_rp
197  ens_mean_w(:,:,:) = 0.0_rp
198  ens_mean_temp(:,:,:) = 0.0_rp
199  ens_mean_pres(:,:,:) = 0.0_rp
200  ens_mean_qv(:,:,:) = 0.0_rp
201  ens_mean_qc(:,:,:) = 0.0_rp
202  ens_mean_qr(:,:,:) = 0.0_rp
203  ens_mean_qi(:,:,:) = 0.0_rp
204  ens_mean_qs(:,:,:) = 0.0_rp
205  ens_mean_qg(:,:,:) = 0.0_rp
206  ens_mean_rh(:,:,:) = 0.0_rp
207 
208  allocate( ens_sprd_dens(ka,ia,ja) )
209  allocate( ens_sprd_momx(ka,ia,ja) )
210  allocate( ens_sprd_momy(ka,ia,ja) )
211  allocate( ens_sprd_momz(ka,ia,ja) )
212  allocate( ens_sprd_rhot(ka,ia,ja) )
213  allocate( ens_sprd_u(ka,ia,ja) )
214  allocate( ens_sprd_v(ka,ia,ja) )
215  allocate( ens_sprd_w(ka,ia,ja) )
216  allocate( ens_sprd_temp(ka,ia,ja) )
217  allocate( ens_sprd_pres(ka,ia,ja) )
218  allocate( ens_sprd_qv(ka,ia,ja) )
219  allocate( ens_sprd_qc(ka,ia,ja) )
220  allocate( ens_sprd_qr(ka,ia,ja) )
221  allocate( ens_sprd_qi(ka,ia,ja) )
222  allocate( ens_sprd_qs(ka,ia,ja) )
223  allocate( ens_sprd_qg(ka,ia,ja) )
224  allocate( ens_sprd_rh(ka,ia,ja) )
225  ens_sprd_dens(:,:,:) = 0.0_rp
226  ens_sprd_momx(:,:,:) = 0.0_rp
227  ens_sprd_momy(:,:,:) = 0.0_rp
228  ens_sprd_momz(:,:,:) = 0.0_rp
229  ens_sprd_rhot(:,:,:) = 0.0_rp
230  ens_sprd_u(:,:,:) = 0.0_rp
231  ens_sprd_v(:,:,:) = 0.0_rp
232  ens_sprd_w(:,:,:) = 0.0_rp
233  ens_sprd_temp(:,:,:) = 0.0_rp
234  ens_sprd_pres(:,:,:) = 0.0_rp
235  ens_sprd_qv(:,:,:) = 0.0_rp
236  ens_sprd_qc(:,:,:) = 0.0_rp
237  ens_sprd_qr(:,:,:) = 0.0_rp
238  ens_sprd_qi(:,:,:) = 0.0_rp
239  ens_sprd_qs(:,:,:) = 0.0_rp
240  ens_sprd_qg(:,:,:) = 0.0_rp
241  ens_sprd_rh(:,:,:) = 0.0_rp
242 
243  call letkf_setup( obs_in_num, &
244  ensemble_world, &
245  ensemble_nprocs, &
246  ensemble_myrank, &
248  prc_nprocs, &
249  prc_myrank, &
250  prc_num_x, &
251  prc_num_y, &
252  ka, ks, ke, &
253  ia, is, ie, &
254  ja, js, je, &
255  kmax, &
256  imax, &
257  jmax, &
258  khalo, &
259  ihalo, &
260  jhalo, &
261  dx, &
262  dy, &
263  topography_zsfc(:,:) )
264 
266 
267  return

References scale_comm_ensemble::comm_ensemble_myrank, scale_comm_ensemble::comm_ensemble_nprocs, scale_comm_ensemble::comm_ensemble_setup(), scale_comm_ensemble::comm_ensemble_world, mod_da_param_estimation::da_param_estimation_setup(), scale_precision::dp, scale_atmos_grid_cartesc::dx, scale_atmos_grid_cartesc::dy, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::ihalo, scale_atmos_grid_cartesc_index::imax, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jhalo, scale_atmos_grid_cartesc_index::jmax, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::khalo, scale_atmos_grid_cartesc_index::kmax, scale_atmos_grid_cartesc_index::ks, scale_letkf::letkf_setup(), scale_atmos_hydrometeor::n_hyd, mod_da_vars::obs_in_num, scale_prc::prc_abort(), scale_prc::prc_local_comm_world, scale_prc::prc_myrank, scale_prc::prc_nprocs, scale_prc_cartesc::prc_num_x, scale_prc_cartesc::prc_num_y, scale_precision::sp, and scale_topography::topography_zsfc.

Referenced by mod_rm_driver::rm_driver().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ da_driver_finalize()

subroutine, public mod_da_driver::da_driver_finalize

finalize

Definition at line 273 of file mod_da_driver.F90.

273  use scale_comm_ensemble, only: &
275  use scale_letkf, only: &
277  use mod_da_param_estimation, only: &
279  implicit none
280  !---------------------------------------------------------------------------
281 
282  log_newline
283  log_info("DA_driver_finalize",*) 'Finalize'
284 
285  deallocate( qhyd )
286 
287  deallocate( rh )
288  deallocate( qdry )
289  deallocate( rtot )
290  deallocate( cvtot )
291  deallocate( cptot )
292 
293  deallocate( prec )
294 
295  deallocate( ens_mean_dens )
296  deallocate( ens_mean_momx )
297  deallocate( ens_mean_momy )
298  deallocate( ens_mean_momz )
299  deallocate( ens_mean_rhot )
300  deallocate( ens_mean_u )
301  deallocate( ens_mean_v )
302  deallocate( ens_mean_w )
303  deallocate( ens_mean_temp )
304  deallocate( ens_mean_pres )
305  deallocate( ens_mean_qv )
306  deallocate( ens_mean_qc )
307  deallocate( ens_mean_qr )
308  deallocate( ens_mean_qi )
309  deallocate( ens_mean_qs )
310  deallocate( ens_mean_qg )
311  deallocate( ens_mean_rh )
312 
313  deallocate( ens_sprd_u )
314  deallocate( ens_sprd_v )
315  deallocate( ens_sprd_w )
316  deallocate( ens_sprd_temp )
317  deallocate( ens_sprd_pres )
318  deallocate( ens_sprd_qv )
319  deallocate( ens_sprd_qc )
320  deallocate( ens_sprd_qr )
321  deallocate( ens_sprd_qi )
322  deallocate( ens_sprd_qs )
323  deallocate( ens_sprd_qg )
324  deallocate( ens_sprd_rh )
325 
327  call letkf_finalize
329 
330  return

References scale_comm_ensemble::comm_ensemble_finalize(), mod_da_param_estimation::da_param_estimation_finalize(), and scale_letkf::letkf_finalize().

Referenced by mod_rm_driver::rm_driver().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ da_driver_update()

subroutine, public mod_da_driver::da_driver_update

Data Assimilation.

Definition at line 336 of file mod_da_driver.F90.

336  use scale_const, &
337  pre00 => const_pre00
338  use scale_atmos_grid_cartesc_index, only: &
339  ka, ks, ke, &
340  ia, is, ie, &
341  ja, js, je
342  use scale_atmos_grid_cartesc_real, only: &
344  use scale_atmos_hydrometeor, only: &
345  i_hc, &
346  i_hr, &
347  i_hi, &
348  i_hs, &
349  i_hg
350  use scale_tracer, only: &
351  qa, &
352  tracer_cv, &
353  tracer_cp, &
354  tracer_r, &
356  use scale_atmos_thermodyn, only: &
357  atmos_thermodyn_specific_heat
358  use scale_comm_ensemble, only: &
359  ensemble_world => comm_ensemble_world, &
360  ensemble_nprocs => comm_ensemble_nprocs, &
361  ensemble_myrank => comm_ensemble_myrank
362  use scale_topography, only: &
363  topo => topography_zsfc
364 #ifdef DA
365  use scale_letkf, only: &
367  letkf_obs_clear, &
371 #endif
372  use mod_atmos_vars, only: &
373  atmos_vars_get_diagnostic, &
374  u, &
375  v, &
376  w, &
377  temp, &
378  pres, &
379  pott, &
380  dens, &
381  momz, &
382  momx, &
383  momy, &
384  rhot, &
385  qtrc, &
386  qv, &
387  qc, &
388  qr, &
389  qi, &
390  qs, &
391  qg
392  use mod_atmos_phy_sf_vars, only: &
393  ps => atmos_phy_sf_sfc_pres, &
394  u10m => atmos_phy_sf_u10, &
395  v10m => atmos_phy_sf_v10, &
396  t2m => atmos_phy_sf_t2, &
397  q2m => atmos_phy_sf_q2
398  use mod_atmos_phy_mp_driver, only: &
400  use mod_da_vars, only: &
401  obs_in_num, &
402  obs_in_format, &
403  obs_in_basename, &
404  obs_in_maskfile, &
407  use mod_da_param_estimation, only: &
409  implicit none
410 
411  character(len=H_LONG) :: filename
412 
413  integer :: i, j, k
414  !---------------------------------------------------------------------------
415 
416 #ifdef DA
417  call prof_rapstart('DA_Update', 1)
418 
419  !---------------------------------------------------------------------------
420  ! Read observation data
421  !---------------------------------------------------------------------------
423 
424  !---------------------------------------------------------------------------
425  ! Get observation operator
426  !---------------------------------------------------------------------------
427  !!
428  !! Compute observation operator, return the results in obsda
429  !! with additional space for externally processed observations
430  !!
431  call atmos_vars_get_diagnostic( 'RH', rh(:,:,:) )
432  call atmos_vars_get_diagnostic( 'PREC', prec(:,:) )
433 
434  rh(:,:,:) = rh(:,:,:) * 0.01_rp ! [%] -> [no unit]
435 
437  u, v, w, temp, pres, qv, qc, qr, qi, qs, qg, rh, hgt, &
438  topo, ps, prec, u10m, v10m, t2m, q2m )
439 
440  !---------------------------------------------------------------------------
441  ! Process observation data
442  !---------------------------------------------------------------------------
444 
445  !---------------------------------------------------------------------------
446  ! Output first guess
447  !---------------------------------------------------------------------------
448  call make_output( 'GUES' )
449 
450  !---------------------------------------------------------------------------
451  ! Data Assimilation (LETKF)
452  !---------------------------------------------------------------------------
453  call letkf_system( obs_in_num, & ! [IN]
454  obs_in_format, & ! [IN]
455  u(ks:ke,is:ie,js:je), & ! [INOUT]
456  v(ks:ke,is:ie,js:je), & ! [INOUT]
457  w(ks:ke,is:ie,js:je), & ! [INOUT]
458  temp(ks:ke,is:ie,js:je), & ! [INOUT]
459  pres(ks:ke,is:ie,js:je), & ! [INOUT]
460  qv(ks:ke,is:ie,js:je), & ! [INOUT]
461  qc(ks:ke,is:ie,js:je), & ! [INOUT]
462  qr(ks:ke,is:ie,js:je), & ! [INOUT]
463  qi(ks:ke,is:ie,js:je), & ! [INOUT]
464  qs(ks:ke,is:ie,js:je), & ! [INOUT]
465  qg(ks:ke,is:ie,js:je) ) ! [INOUT]
466 
467  !---------------------------------------------------------------------------
468  ! Parameter Estimation (ETKF)
469  !---------------------------------------------------------------------------
471 
472  !---------------------------------------------------------------------------
473  ! Clearing observation data structure (each step)
474  !---------------------------------------------------------------------------
476 
477  !---------------------------------------------------------------------------
478  ! QV/QHYD negative filter
479  !---------------------------------------------------------------------------
480  if( positive_definite_q ) then
481  qv(:,:,:) = max( qv(:,:,:), 0.0_rp )
482  end if
483  if( positive_definite_qhyd ) then
484  qc(:,:,:) = max( qc(:,:,:), 0.0_rp )
485  qr(:,:,:) = max( qr(:,:,:), 0.0_rp )
486  qi(:,:,:) = max( qi(:,:,:), 0.0_rp )
487  qs(:,:,:) = max( qs(:,:,:), 0.0_rp )
488  qg(:,:,:) = max( qg(:,:,:), 0.0_rp )
489  end if
490 
491  !---------------------------------------------------------------------------
492  ! overwrite the prognostic variables
493  !---------------------------------------------------------------------------
494  do j = js, je
495  do i = is, ie
496  do k = ks, ke
497  qhyd(k,i,j,i_hc) = qc(k,i,j)
498  qhyd(k,i,j,i_hr) = qr(k,i,j)
499  qhyd(k,i,j,i_hi) = qi(k,i,j)
500  qhyd(k,i,j,i_hs) = qs(k,i,j)
501  qhyd(k,i,j,i_hg) = qg(k,i,j)
502  enddo
503  enddo
504  enddo
505 
506  call atmos_phy_mp_driver_qhyd2qtrc( ka, ks, ke, & ! [in]
507  ia, is, ie, & ! [in]
508  ja, js, je, & ! [in]
509  qv, & ! [in]
510  qhyd, & ! [in]
511  qtrc ) ! [out]
512 
513  call atmos_thermodyn_specific_heat( &
514  ka, ks, ke, ia, is, ie, ja, js, je, qa, & ! (in)
515  qtrc(:,:,:,:), & ! (in)
516  tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), & ! (in)
517  qdry(:,:,:), rtot(:,:,:), cvtot(:,:,:), cptot(:,:,:) ) ! (out)
518 
519  do j = js, je
520  do i = is, ie
521  do k = ks, ke
522  dens(k,i,j) = pres(k,i,j) / ( temp(k,i,j) * rtot(k,i,j) )
523  pott(k,i,j) = temp(k,i,j) * ( pre00 / pres(k,i,j) )**(rtot(k,i,j)/cptot(k,i,j))
524  rhot(k,i,j) = dens(k,i,j) * pott(k,i,j)
525  enddo
526  enddo
527  enddo
528 
529  do j = js, je
530  do i = is, ie
531  do k = ks, ke-1
532  momz(k,i,j) = 0.5_rp * ( w(k,i,j)*dens(k,i,j) + w(k+1,i,j)*dens(k+1,i,j) )
533  enddo
534  enddo
535  enddo
536  do j = js, je
537  do i = is, ie
538  momz(ke,i,j) = w(ke,i,j) * dens(ke,i,j)
539  enddo
540  enddo
541 
542  do j = js, je
543  do i = is, ie-1
544  do k = ks, ke
545  momx(k,i,j) = 0.5_rp * ( u(k,i,j)*dens(k,i,j) + u(k,i+1,j)*dens(k,i+1,j) )
546  enddo
547  enddo
548  enddo
549  do j = js, je
550  do k = ks, ke
551  momx(k,ie,j) = u(k,ie,j) * dens(k,ie,j)
552  enddo
553  enddo
554 
555  do j = js, je-1
556  do i = is, ie
557  do k = ks, ke
558  momy(k,i,j) = 0.5_rp * ( v(k,i,j)*dens(k,i,j) + v(k,i,j+1)*dens(k,i,j+1) )
559  end do
560  end do
561  end do
562  do i = is, ie
563  do k = ks, ke
564  momy(k,i,je) = v(k,i,je) * dens(k,i,je)
565  enddo
566  enddo
567 
568  !---------------------------------------------------------------------------
569  ! Output analysis
570  !---------------------------------------------------------------------------
571  call make_output( 'ANLS' )
572 
573  call prof_rapend ('DA_Update', 1)
574 #endif
575 
576  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc(), mod_atmos_phy_sf_vars::atmos_phy_sf_q2, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_pres, mod_atmos_phy_sf_vars::atmos_phy_sf_t2, mod_atmos_phy_sf_vars::atmos_phy_sf_u10, mod_atmos_phy_sf_vars::atmos_phy_sf_v10, scale_comm_ensemble::comm_ensemble_myrank, scale_comm_ensemble::comm_ensemble_nprocs, scale_comm_ensemble::comm_ensemble_world, scale_const::const_pre00, scale_const::const_undef, mod_da_vars::da_compute_ens_history, mod_da_param_estimation::da_param_estimation_update(), mod_atmos_vars::dens, scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hg, scale_atmos_hydrometeor::i_hi, scale_atmos_hydrometeor::i_hr, scale_atmos_hydrometeor::i_hs, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::imax, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jmax, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::kmax, scale_atmos_grid_cartesc_index::ks, scale_letkf::letkf_obs_clear(), scale_letkf::letkf_obs_initialize(), scale_letkf::letkf_obs_operator(), scale_letkf::letkf_obs_readfile(), scale_letkf::letkf_system(), mod_atmos_vars::momx, mod_atmos_vars::momy, mod_atmos_vars::momz, mod_da_vars::obs_in_basename, mod_da_vars::obs_in_format, mod_da_vars::obs_in_maskfile, mod_da_vars::obs_in_num, mod_da_vars::positive_definite_q, mod_da_vars::positive_definite_qhyd, mod_atmos_vars::pott, scale_prc::prc_abort(), mod_atmos_vars::pres, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, mod_atmos_vars::qc, mod_atmos_vars::qg, mod_atmos_vars::qi, mod_atmos_vars::qr, mod_atmos_vars::qs, mod_atmos_vars::qtrc, mod_atmos_vars::qv, mod_atmos_vars::rhot, mod_atmos_vars::temp, scale_topography::topography_zsfc, scale_tracer::tracer_cp, scale_tracer::tracer_cv, scale_tracer::tracer_mass, scale_tracer::tracer_r, mod_atmos_vars::u, mod_atmos_vars::v, and mod_atmos_vars::w.

Referenced by mod_rm_driver::rm_driver().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_comm_ensemble::comm_ensemble_world
integer, public comm_ensemble_world
Definition: scale_comm_ensemble.F90:35
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_hydrometeor::i_hr
integer, parameter, public i_hr
liquid water rain
Definition: scale_atmos_hydrometeor.F90:98
mod_da_vars::obs_in_basename
character(len=h_long), dimension(obs_in_basename_maxsize), public obs_in_basename
Definition: mod_da_vars.F90:46
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
mod_atmos_vars::qr
real(rp), dimension(:,:,:), pointer, public qr
Definition: mod_atmos_vars.F90:99
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:39
mod_atmos_phy_sf_vars::atmos_phy_sf_v10
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_v10
Definition: mod_atmos_phy_sf_vars.F90:96
mod_atmos_vars::pott
real(rp), dimension(:,:,:), allocatable, target, public pott
Definition: mod_atmos_vars.F90:133
scale_atmos_grid_cartesc_index::ihalo
integer, public ihalo
Definition: scale_atmos_grid_cartesC_index.F90:44
scale_atmos_hydrometeor::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:100
mod_da_param_estimation
module Parameter Estimation
Definition: mod_da_param_estimation.F90:11
scale_tracer::tracer_mass
real(rp), dimension(qa_max), public tracer_mass
Definition: scale_tracer.F90:47
scale_comm_ensemble
module Communication for Ensemble system
Definition: scale_comm_ensemble.F90:12
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_prc::prc_local_comm_world
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:89
scale_comm_ensemble::comm_ensemble_setup
subroutine, public comm_ensemble_setup
Setup.
Definition: scale_comm_ensemble.F90:52
scale_letkf::letkf_obs_operator
subroutine, public letkf_obs_operator(OBS_IN_NUM, OBS_IN_FORMAT, U, V, W, TEMP, PRES, QV, QC, QR, QI, QS, QG, RH, HGT, TOPO, PS, RAIN, U10M, V10M, T2M, Q2M, nobs_extern)
Definition: scale_letkf.F90:1061
scale_letkf::letkf_obs_readfile
subroutine, public letkf_obs_readfile(OBS_IN_NUM, OBS_IN_FORMAT, OBS_IN_BASENAME, OBS_IN_MASKFILE)
Definition: scale_letkf.F90:906
scale_letkf
module LETKF for Data-Assimilation
Definition: scale_letkf.F90:12
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
mod_atmos_phy_sf_vars
module ATMOSPHERIC Surface Variables
Definition: mod_atmos_phy_sf_vars.F90:12
mod_da_vars::obs_in_format
character(len=h_long), dimension(obs_in_basename_maxsize), public obs_in_format
Definition: mod_da_vars.F90:47
mod_da_param_estimation::da_param_estimation_finalize
subroutine, public da_param_estimation_finalize
finalize
Definition: mod_da_param_estimation.F90:147
scale_atmos_grid_cartesc_index::khalo
integer, parameter, public khalo
Definition: scale_atmos_grid_cartesC_index.F90:43
scale_topography
module TOPOGRAPHY
Definition: scale_topography.F90:11
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:80
scale_letkf::letkf_obs_initialize
subroutine, public letkf_obs_initialize(OBS_IN_NUM, nobs_extern)
Definition: scale_letkf.F90:1453
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:81
scale_atmos_grid_cartesc::dx
real(rp), public dx
Definition: scale_atmos_grid_cartesC.F90:39
scale_atmos_grid_cartesc_index::imax
integer, public imax
Definition: scale_atmos_grid_cartesC_index.F90:37
mod_da_vars
module Data Assimilation Variables
Definition: mod_da_vars.F90:12
mod_atmos_phy_mp_driver
module atmosphere / physics / cloud microphysics
Definition: mod_atmos_phy_mp_driver.F90:12
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_u10
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_u10
Definition: mod_atmos_phy_sf_vars.F90:95
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:99
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
mod_da_param_estimation::da_param_estimation_update
subroutine, public da_param_estimation_update
Data Assimilation.
Definition: mod_da_param_estimation.F90:158
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:76
scale_comm_ensemble::comm_ensemble_nprocs
integer, public comm_ensemble_nprocs
Definition: scale_comm_ensemble.F90:37
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
mod_da_vars::positive_definite_q
logical, public positive_definite_q
Definition: mod_da_vars.F90:50
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:77
scale_letkf::letkf_finalize
subroutine, public letkf_finalize()
Definition: scale_letkf.F90:881
scale_atmos_grid_cartesc_index::kmax
integer, public kmax
Definition: scale_atmos_grid_cartesC_index.F90:36
scale_comm_ensemble::comm_ensemble_finalize
subroutine, public comm_ensemble_finalize
finalize
Definition: scale_comm_ensemble.F90:116
scale_tracer::tracer_cv
real(rp), dimension(qa_max), public tracer_cv
Definition: scale_tracer.F90:42
mod_da_vars::obs_in_maskfile
character(len=h_long), public obs_in_maskfile
Definition: mod_da_vars.F90:48
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
mod_atmos_vars::v
real(rp), dimension(:,:,:), allocatable, target, public v
Definition: mod_atmos_vars.F90:131
mod_atmos_vars::w
real(rp), dimension(:,:,:), allocatable, target, public w
Definition: mod_atmos_vars.F90:129
mod_atmos_vars::momx
real(rp), dimension(:,:,:), allocatable, target, public momx
Definition: mod_atmos_vars.F90:78
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::jhalo
integer, public jhalo
Definition: scale_atmos_grid_cartesC_index.F90:45
mod_atmos_vars::temp
real(rp), dimension(:,:,:), allocatable, target, public temp
Definition: mod_atmos_vars.F90:134
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_atmos_vars::momy
real(rp), dimension(:,:,:), allocatable, target, public momy
Definition: mod_atmos_vars.F90:79
mod_atmos_vars::qv
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
Definition: mod_atmos_vars.F90:97
scale_atmos_hydrometeor::i_hc
integer, parameter, public i_hc
liquid water cloud
Definition: scale_atmos_hydrometeor.F90:97
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
mod_atmos_vars::pres
real(rp), dimension(:,:,:), allocatable, target, public pres
Definition: mod_atmos_vars.F90:135
scale_atmos_grid_cartesc::dy
real(rp), public dy
Definition: scale_atmos_grid_cartesC.F90:39
mod_atmos_vars::u
real(rp), dimension(:,:,:), allocatable, target, public u
Definition: mod_atmos_vars.F90:130
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_topography::topography_zsfc
real(rp), dimension(:,:), allocatable, public topography_zsfc
absolute ground height [m]
Definition: scale_topography.F90:39
scale_prc::prc_nprocs
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:90
mod_atmos_phy_sf_vars::atmos_phy_sf_t2
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_t2
Definition: mod_atmos_phy_sf_vars.F90:97
scale_tracer::tracer_cp
real(rp), dimension(qa_max), public tracer_cp
Definition: scale_tracer.F90:43
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_pres
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_pres
Definition: mod_atmos_phy_sf_vars.F90:73
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
scale_atmos_grid_cartesc_index::jmax
integer, public jmax
Definition: scale_atmos_grid_cartesC_index.F90:38
mod_atmos_vars::qi
real(rp), dimension(:,:,:), pointer, public qi
Definition: mod_atmos_vars.F90:100
scale_tracer::tracer_r
real(rp), dimension(qa_max), public tracer_r
Definition: scale_tracer.F90:44
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_prc_cartesc::prc_num_y
integer, public prc_num_y
y length of 2D processor topology
Definition: scale_prc_cartesC.F90:43
scale_atmos_thermodyn
module atmosphere / thermodyn
Definition: scale_atmos_thermodyn.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_q2
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_q2
Definition: mod_atmos_phy_sf_vars.F90:98
mod_da_vars::obs_in_num
integer, public obs_in_num
Definition: mod_da_vars.F90:53
mod_atmos_vars::qc
real(rp), dimension(:,:,:), pointer, public qc
Definition: mod_atmos_vars.F90:98
scale_prc_cartesc::prc_num_x
integer, public prc_num_x
x length of 2D processor topology
Definition: scale_prc_cartesC.F90:42
scale_comm_ensemble::comm_ensemble_myrank
integer, public comm_ensemble_myrank
Definition: scale_comm_ensemble.F90:36
scale_letkf::letkf_setup
subroutine, public letkf_setup(OBS_IN_NUM, ensemble_comm, ensemble_nprocs, ensemble_myrank, local_comm, local_nprocs, local_myrank, PRC_NUM_X, PRC_NUM_Y, KA, KS, KE, IA, IS, IE, JA, JS, JE, KMAX, IMAX, JMAX, KHALO, IHALO, JHALO, delta_x, delta_y, Zsfc)
Setup.
Definition: scale_letkf.F90:615
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_letkf::letkf_system
subroutine, public letkf_system(OBS_IN_NUM, OBS_IN_FORMAT, U, V, W, TEMP, PRES, QV, QC, QR, QI, QS, QG)
Definition: scale_letkf.F90:2189
scale_letkf::letkf_obs_clear
subroutine, public letkf_obs_clear(OBS_IN_NUM)
Definition: scale_letkf.F90:1016
mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc
subroutine, public atmos_phy_mp_driver_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
Definition: mod_atmos_phy_mp_driver.F90:1553
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95
mod_da_vars::positive_definite_qhyd
logical, public positive_definite_qhyd
Definition: mod_da_vars.F90:51
mod_da_param_estimation::da_param_estimation_setup
subroutine, public da_param_estimation_setup
Setup.
Definition: mod_da_param_estimation.F90:70
mod_atmos_vars::qg
real(rp), dimension(:,:,:), pointer, public qg
Definition: mod_atmos_vars.F90:102
mod_atmos_vars::qs
real(rp), dimension(:,:,:), pointer, public qs
Definition: mod_atmos_vars.F90:101
scale_atmos_hydrometeor::i_hg
integer, parameter, public i_hg
graupel
Definition: scale_atmos_hydrometeor.F90:101