wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_control / da_control.f90
blob9d71e5c62a088a2850f4b50cb2dc3e2e9f943499
1 module da_control
3 !--------------------------------------------------------------------------
4 ! Purpose: Common reference point for WRFVAR control.
5 !--------------------------------------------------------------------------
7 use module_driver_constants, only : max_domains, max_eta, max_moves, max_bogus, &
8 max_outer_iterations, max_instruments
10 implicit none
12 #include "namelist_defines.inc"
14 ! switches set from other namelist options
15 logical :: use_obsgts
16 logical :: use_rad
18 !---------------------------------------------------------------------------
19 ! [1.0] Physical parameter constants (all NIST standard values):
20 !---------------------------------------------------------------------------
22 ! Fundamental constants:
23 real, parameter :: pi = 3.1415926 ! Value used in WRF.
24 real, parameter :: radian = pi / 180.0
25 real, parameter :: gas_constant = 287.0 ! Value used in WRF.
26 real, parameter :: gas_constant_v = 461.6 ! Value used in WRF.
27 real, parameter :: cp = 7.0*gas_constant/2.0 ! Value used in WRF.
28 real, parameter :: t_kelvin = 273.15
29 real, parameter :: t_triple = 273.16 ! triple point of water
30 ! The imported code for ssmi and radiance uses 273.0 in a way that suggests
31 ! it may not be a lazy definition of the melting point of water, so keep the
32 ! value separate for the moment
33 real, parameter :: t_roughem = 273.0
34 real, parameter :: t_landem = 273.0
36 real, parameter :: kappa = gas_constant / cp
37 real, parameter :: rd_over_rv = gas_constant / gas_constant_v
38 real, parameter :: rd_over_rv1 = 1.0 - rd_over_rv
39 real, parameter :: L_over_Rv = 5418.12
41 real, parameter :: gamma = 1.4
43 ! Earth constants:
44 real, parameter :: gravity = 9.81 ! m/s - value used in WRF.
45 ! real, parameter :: earth_radius = 6378.15
46 real, parameter :: earth_radius = 6370.0 ! Be consistant with WRF
47 ! real, parameter :: earth_omega = 2.0*pi/86400.0 ! Omega
48 real, parameter :: earth_omega = 0.000072921 ! Omega 7.2921*10**-5
50 ! Saturation Vapour Pressure Constants(Rogers & Yau, 1989)
51 real, parameter :: es_alpha = 611.2
52 real, parameter :: es_beta = 17.67
53 real, parameter :: es_gamma = 243.5
54 real, parameter :: es_gammabeta = es_gamma * es_beta
55 real, parameter :: es_gammakelvin = es_gamma - t_kelvin
57 ! Explicit moist constants:
58 real, parameter :: SVP1=0.6112, SVP2=17.67, SVP3=29.65
59 real, parameter :: SVPT0=t_kelvin, TO=t_kelvin
60 real, parameter :: N0R=8.0E6, N0S=2.0E7, RHOS=0.1
61 real, parameter :: AVT=841.99667, BVT=0.8, BVT2=2.5+0.5*BVT, BVT3=3.0+BVT
62 real, parameter :: PPI=1.0/(pi*N0R), PPIS=1.0/(pi*N0S*RHOS)
63 real, parameter :: XLV1=2370.0, XLF0=0.3337E6, XLV0=3.15E6
64 real, parameter :: XLS=XLV0-XLV1*t_triple+XLF0
66 ! Planetary boundary physics constants
67 real, parameter :: k_kar = 0.4 ! Von Karman constant
69 ! Zenith Total Delay:
70 ! Hydrostatic delay:
71 real, parameter :: zdk1 = 2.2768e-5
72 real, parameter :: zdk2 = 2.66e-3
73 real, parameter :: zdk3 = 2.8e-7
74 ! Wet delay:
75 real, parameter :: wdk1 = 2.21e-7
76 real, parameter :: wdk2 = 3.73e-3
78 ! GPS Refractivity constant
79 real, parameter :: a_ew = 0.622
80 real, parameter :: b_ew = 0.378
82 ! GPS Refractivity constant
83 real, parameter :: coeff = (wdk2*1.e8) / 77.6
85 #if RWORDSIZE==8
86 real, parameter :: da_zero = 0D0
87 #else
88 real, parameter :: da_zero = 0.0
89 #endif
91 complex, parameter :: da_zero_complex = (da_zero,da_zero)
93 !---------------------------------------------------------------------------
94 ! [2.0] WRF-Var parameter constants:
95 !---------------------------------------------------------------------------
97 ! Missing values and the index number of the quality control
99 integer, parameter :: missing = -888888
100 real , parameter :: missing_r = -888888.0
101 real , parameter :: Max_StHeight_Diff = 100.0
103 integer, parameter :: cv_options_hum_specific_humidity = 1
104 integer, parameter :: cv_options_hum_relative_humidity = 2
106 ! No-one explains what these options means anywhere
107 integer, parameter :: vert_corr_1 = 1
108 integer, parameter :: vert_corr_2 = 2
110 integer, parameter :: vertical_ip_0 = 0
111 integer, parameter :: vertical_ip_sqrt_delta_p = 1
112 integer, parameter :: vertical_ip_delta_p = 2
114 integer, parameter :: vert_evalue_global = 1
115 integer, parameter :: vert_evalue_local = 2
117 integer, parameter :: alphacv_method_vp = 1
118 integer, parameter :: alphacv_method_xa = 2
120 integer, parameter :: sfc_assi_options_1 = 1
121 integer, parameter :: sfc_assi_options_2 = 2
123 integer, parameter :: check_rh_simple = 1
124 integer, parameter :: check_rh_tpw = 2
126 logical :: anal_type_verify=.false.
127 logical :: anal_type_randomcv=.false.
128 logical :: anal_type_qcobs=.false.
130 integer,parameter :: monitor_on = 1
131 integer,parameter :: monitor_off = 0
133 integer,parameter :: qc_good = 1
134 integer,parameter :: qc_bad = -1
135 integer,parameter :: qc_varbc_bad = -1
137 integer, parameter :: bufr_satellite_id = 1
138 integer, parameter :: bufr_ifov = 2
139 integer, parameter :: bufr_year = 3
140 integer, parameter :: bufr_month = 4
141 integer, parameter :: bufr_day = 5
142 integer, parameter :: bufr_hour = 6
143 integer, parameter :: bufr_minute = 7
144 integer, parameter :: bufr_second = 8
145 integer, parameter :: bufr_lat = 9
146 integer, parameter :: bufr_lon = 10
147 integer, parameter :: bufr_satzen = 11
148 integer, parameter :: bufr_solzen = 12
149 integer, parameter :: bufr_station_height = 13
150 integer, parameter :: bufr_landsea_mask = 14
152 integer, parameter :: nchan_amsua = 15
153 integer, parameter :: nchan_amsub = 5
154 integer, parameter :: nchan_mhs = 5
155 integer, parameter :: nchan_msu = 4
156 integer, parameter :: nchan_hirs2 = 19
157 integer, parameter :: nchan_hirs3 = 19
158 integer, parameter :: nchan_hirs4 = 19
159 integer, parameter :: nchan_ssmis = 24
160 integer, parameter :: nchan_airs = 281
162 ! WRFVAR Minimisation:
164 integer :: iter
165 integer :: cv_size
166 integer, parameter :: MP = 6
167 integer, parameter :: LP = 6
168 integer, parameter :: MAXFEV = 10
169 real, parameter :: FTOL = 1.0E-4
170 real, parameter :: GTOL = 0.9
171 real, parameter :: XTOL = 1.0E-17
172 real, parameter :: STPMIN = 1.0E-20
173 real, parameter :: STPMAX = 1.0E+20
175 ! Background errors:
176 real, parameter :: pplow = 1.0e-8 ! Machine lowest number?
177 real, parameter :: pp_umin = 1.0e-2 ! Minimum u back. error (m/s).
178 real, parameter :: pp_vmin = 1.0e-2 ! Minimum v back. error (m/s).
179 real, parameter :: pp_tmin = 1.0e-2 ! Minimum t back. error (K).
180 real, parameter :: pp_qmin = 1.0e-6 ! Minimum q back. error (kg/kg)
181 real, parameter :: pp_pmin= 1.0e+1 ! Minimum pp back. error (Pa).
183 ! FFTs:
184 integer, parameter :: Forward_FFT = -1 ! Grid to spectral
185 integer, parameter :: Inverse_FFT = 1 ! Spectral to grid.
186 integer, parameter :: num_fft_factors = 10 ! Max number of factors.
188 ! Balance:
189 integer, parameter :: balance_geo = 1 ! Geostrophic balance.
190 integer, parameter :: balance_cyc = 2 ! Cyclostrophic balance.
191 integer, parameter :: balance_geocyc = 3 ! Geo/cyclostrophic balance.
193 ! Adjoint tests:
194 real, parameter :: typical_u_rms = 2.0 ! m/s
195 real, parameter :: typical_v_rms = 2.0 ! m/s
196 real, parameter :: typical_speed_rms = 2.0 ! m/s
197 real, parameter :: typical_tb19v_rms = 1.0 ! K
198 real, parameter :: typical_tb19h_rms = 1.0 ! K
199 real, parameter :: typical_tb22v_rms = 1.0 ! K
200 real, parameter :: typical_tb37v_rms = 1.0 ! K
201 real, parameter :: typical_tb37h_rms = 1.0 ! K
202 real, parameter :: typical_tb85v_rms = 1.0 ! K
203 real, parameter :: typical_tb85h_rms = 1.0 ! K
204 real, parameter :: typical_t_rms = 1.0 ! K
205 real, parameter :: typical_p_rms = 100.0 ! Pa
206 real, parameter :: typical_q_rms = 0.00001 ! g/kg
207 real, parameter :: typical_rho_rms = 0.01 ! kg/m^3
208 real, parameter :: typical_tpw_rms = 0.2 ! cm
209 real, parameter :: typical_ref_rms = 5.0 ! N unit
210 real, parameter :: typical_rh_rms = 20.0 ! %
211 real, parameter :: typical_thickness_rms = 50.0 ! m
212 real, parameter :: typical_qrn_rms = 0.00001 ! g/kg
213 real, parameter :: typical_qcw_rms = 0.00001 ! g/kg
214 real, parameter :: typical_w_rms = 0.1 ! m/s
215 real, parameter :: typical_rv_rms = 1.0 ! m/s
216 real, parameter :: typical_rf_rms = 1.0 ! dBZ
218 ! The following typical mean squared values depend on control variable. They
219 ! are calculated in da_setup_background_errors and used in the VvToVp adjoint
220 ! test:
222 real, parameter :: inv_typ_vp1_sumsq = 0.00001 ! 1/sum(psi**2)
223 real, parameter :: inv_typ_vp2_sumsq = 0.00001 ! 1/sum(chi**2)
224 real, parameter :: inv_typ_vp3_sumsq = 0.00001 ! 1/sum(phi_u**2)
225 real, parameter :: inv_typ_vp4_sumsq = 10000.0 ! 1/sum(q**2)
226 real, parameter :: inv_typ_vp5_sumsq = 0.00001 ! 1/sum(?**2)
227 real, parameter :: inv_typ_vpalpha_sumsq = 1.0 ! 1/sum(?**2)
229 #include "version_decl"
231 integer, parameter :: fg_format_wrf_arw_regional = 1
232 integer, parameter :: fg_format_wrf_nmm_regional = 2
233 integer, parameter :: fg_format_wrf_arw_global = 3
234 integer, parameter :: fg_format_kma_global = 4
236 integer, parameter :: ob_format_bufr = 1
237 integer, parameter :: ob_format_ascii = 2
238 integer, parameter :: ob_format_madis = 3
240 integer, parameter :: convert_fd2uv = 1
241 integer, parameter :: convert_uv2fd = -1
243 ! Fortran unit parameters:
245 ! stdout, stderr, trace_unit all controlled from namelist
247 ! Units 9,10 are used for reading and writing namelist.input/output in WRF
249 ! Do not use get_unit/free_unit because tracing is too low level
250 integer, parameter :: trace_csv_unit = 8
252 integer :: y_unit, yp_unit, cost_unit, grad_unit, stats_unit, jo_unit
253 integer :: check_max_iv_unit, check_buddy_unit, rand_unit, omb_unit, &
254 filtered_obs_unit
255 integer :: biasprep_unit, qcstat_conv_unit
257 integer,parameter :: filename_len = 200
259 integer, parameter :: num_alpha_corr_types = 3
261 integer, parameter :: alpha_corr_type_exp = 1
262 integer, parameter :: alpha_corr_type_soar = 2
263 integer, parameter :: alpha_corr_type_gaussian = 3
265 integer :: alpha_corr_unit1(num_alpha_corr_types)
266 integer :: alpha_corr_unit2(num_alpha_corr_types)
268 integer, parameter :: max_num_of_var = 200 ! Maximum # of stored fields.
270 integer, parameter :: unit_start = 20
271 integer, parameter :: unit_end = 500
272 logical :: unit_used(unit_start:unit_end) = .false.
274 ! grid properties
276 character(len=3), parameter :: grid_ordering = "xyz"
277 character(len=3), parameter :: grid_stagger = "xyz"
279 !---------------------------------------------------------------------------
280 ! [3.0] Variables used in MM5 part of code:
281 !---------------------------------------------------------------------------
283 integer :: map_projection ! 1=LamConf/2=PolarSte/3=Mercator
284 real :: ycntr
285 integer :: coarse_ix ! coarse domain dim in i direction.
286 integer :: coarse_jy ! coarse domain dim in y direction.
287 real :: coarse_ds ! Coarse domain gridlength (km)
288 real :: start_x ! i posn. of (1,1) in coarse domain.
289 real :: start_y ! j posn. of (1,1) in coarse domain.
290 real :: start_lat ! Latitude coresponds to start_(x,y)
291 real :: start_lon ! Longitude coresponds to start_(x,y)
292 real :: delt_lat ! Latitude increments for global grids
293 real :: delt_lon ! Longitude increments for global grids
295 real :: phic ! coarse domain central lat(degree)
296 real :: xlonc ! coarse domain central lon(degree)
297 real :: cone_factor ! Cone Factor
298 real :: truelat1_3dv ! True latitude 1 (degrees)
299 real :: truelat2_3dv ! True latitude 2 (degrees)
300 real :: pole ! Pole latitude (degrees)
301 real :: dsm ! Current domain gridlength (km)
302 real :: psi1 ! ?
303 real :: c2 ! earth_radius * COS(psi1)
305 real :: ptop
306 real, parameter :: t0 = 300.0
308 !------------------------------------------------------------------------------
309 ! 4.0 vertical interpolation options
310 !------------------------------------------------------------------------------
312 integer, parameter :: v_interp_not_specified = missing, &
313 v_interp_p = 1, &
314 v_interp_h = 2
316 !------------------------------------------------------------------------------
317 ! WRFVAR scalar constants:
318 !------------------------------------------------------------------------------
320 integer :: Anal_Space ! Space of analysis
321 ! ( 1 = Full model,
322 ! 2 = Transformed grid,
323 ! 3 = Ob space (PSAS) )
325 integer :: mix ! 1st dimension of analysis grid.
326 integer :: mjy ! 2nd dimension of analysis grid.
327 integer :: mkz ! 3rd dimension of analysis grid.
329 ! Recursive filter:
331 real, allocatable :: rf_turnconds(:) ! RF turning conditions.
333 integer, parameter :: max_ob_levels = 1001 ! Maximum levels for single ob
334 integer, parameter :: max_fgat_time = 100 ! Maximum levels for FGAT.
336 integer :: time
338 logical :: gaussian_lats
341 integer :: cv_size_domain_jb ! Total jb cv size.
342 integer :: cv_size_domain_je ! Total je cv size.
343 integer :: cv_size_domain_jp ! Total jp cv size.
344 integer :: cv_size_domain_js ! Total js cv size.
345 integer :: cv_size_domain ! Total cv size.
347 ! Hybrid:
348 real :: sigma_alpha ! Alpha standard deviation.
349 real :: jb_factor ! Weighting for Background Error Cov.
351 ! Namelist variables in future?:
352 real, parameter :: maximum_rh = 100.0
353 real, parameter :: minimum_rh = 10.0
355 ! other
357 integer, parameter :: jperr = 6
359 ! NCEP errors (U in m/s, V in m/s, T in K, H in %, P in Pa)
360 ! rh has been divided by 2
362 real, parameter :: err_k(0:jperr+1) = &
363 (/200000.0, 100100.0,70000.0,50000.0,30000.0,10000.0,5000.0, 1.0/)
364 real, parameter :: err_u(0:jperr+1) = &
365 (/ 1.4, 1.4, 2.4, 2.8, 3.4, 2.5, 2.7, 2.7/)
366 real, parameter :: err_v(0:jperr+1) = &
367 (/ 1.4, 1.4, 2.4, 2.8, 3.4, 2.5, 2.7 , 2.7 /)
368 real, parameter :: err_t(0:jperr+1) = &
369 (/ 1.8, 1.8, 1.3, 1.3, 2.0, 3.1, 4.0 , 4.0 /)
370 real, parameter :: err_rh(0:jperr+1) = &
371 (/ 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0/)
372 real, parameter :: err_p(0:jperr+1) = &
373 (/ 100.0,100.0, 100.0, 100.0, 100.0, 100.0,100.0,100.0 /)
375 ! Buddy check parameters (YRG, 10/3/2008):
377 real, parameter :: max_buddy_t = 8.0, &
378 max_buddy_uv = 8.0, &
379 max_buddy_z = 8.0, &
380 max_buddy_rh = 40.0, &
381 max_buddy_p = 350.0, &
382 buddy_weight = 1.0, &
383 bin_p_width = 5000.0, &
384 bin_z_width = 500.0
386 ! Define various ways for bad data to be flagged.
388 integer, parameter :: &
389 missing_data = -88, & ! Data is missing with the value of
390 ! missing_r
391 outside_of_domain = -77, & ! Data outside horizontal domain
392 ! or time window, data set to missing_r
393 wrong_direction = -15, & ! Wind vector direction <0 or> 360
394 ! => direction set to missing_r
395 negative_spd = -14, & ! Wind vector norm is negative
396 ! => norm set to missing_r
397 zero_spd = -13, & ! Wind vector norm is zero
398 ! => norm set to missing_r
399 wrong_wind_data = -12, & ! Spike in wind profile
400 ! =>direction and norm set to missing_r
401 zero_t_td = -11, & ! t or td = 0 => t or td, rh and qv
402 ! are set to missing_r,
403 t_fail_supa_inver = -10, & ! superadiabatic temperature
405 wrong_t_sign = - 9, & ! Spike in Temperature profile
407 above_model_lid = - 8, & ! heigh above model lid
408 ! => no action
409 far_below_model_surface = - 7, & ! heigh far below model surface
410 ! => no action
411 below_model_surface = - 6, & ! height below model surface
412 ! => no action
413 standard_atmosphere = - 5, & ! Missing h, p or t
414 ! =>Datum interpolated from standard atm
415 from_background = - 4, & ! Missing h, p or t
416 ! =>Datum interpolated from model
417 fails_error_max = - 3, & ! Datum Fails error max check
418 ! => no action
419 fails_buddy_check = - 2, & ! Datum Fails buddy check
420 ! => no action
421 no_buddies = - 1, & ! Datum has no buddies
422 ! => no action
423 good_quality = 0, & ! OBS datum has good quality
425 convective_adjustment = 1, & ! convective adjustement check
426 ! =>apply correction on t, td, rh and qv
427 surface_correction = 2, & ! Surface datum
428 ! => apply correction on datum
429 Hydrostatic_recover = 3, & ! Height from hydrostaic assumption with
430 ! the OBS data calibration
431 Reference_OBS_recover = 4, & ! Height from reference state with OBS
432 ! data calibration
433 Other_check = 88 ! passed other quality check
435 ! Observations:
437 integer :: num_procs ! Number of total processors.
438 integer :: myproc ! My processor ID.
439 integer, parameter :: root = 0 ! Number of root processor
440 logical :: rootproc ! Am I the root processor
442 integer, parameter :: var4d_coupling_disk_linear = 1
443 integer, parameter :: var4d_coupling_disk_simul = 2
445 integer, parameter :: rtm_option_rttov = 1
446 integer, parameter :: rtm_option_crtm = 2
448 ! rtm_init setup parameter
450 integer, parameter :: maxsensor = 30
452 integer, parameter :: num_ob_indexes = 27
453 integer, parameter :: npres_print = 12
456 ! Tracing
458 integer :: trace_start_points=0 ! Number of routines to initiate trace
460 integer, parameter :: sound = 1
461 integer, parameter :: synop = 2
462 integer, parameter :: pilot = 3
463 integer, parameter :: satem = 4
464 integer, parameter :: geoamv = 5
465 integer, parameter :: polaramv = 6
466 integer, parameter :: airep = 7
467 integer, parameter :: gpspw = 8
468 integer, parameter :: gpsref = 9
469 integer, parameter :: metar = 10
470 integer, parameter :: ships = 11
471 integer, parameter :: ssmi_rv = 12
472 integer, parameter :: ssmi_tb = 13
473 integer, parameter :: ssmt1 = 14
474 integer, parameter :: ssmt2 = 15
475 integer, parameter :: qscat = 16
476 integer, parameter :: profiler = 17
477 integer, parameter :: buoy = 18
478 integer, parameter :: bogus = 19
479 integer, parameter :: pseudo = 20
480 integer, parameter :: radar = 21
481 integer, parameter :: radiance = 22
482 integer, parameter :: airsr = 23
483 integer, parameter :: sonde_sfc = 24
484 integer, parameter :: mtgirs = 25
485 integer, parameter :: tamdar = 26
486 integer, parameter :: tamdar_sfc = 27
488 character(len=14), parameter :: obs_names(num_ob_indexes) = (/ &
489 "sound ", &
490 "synop ", &
491 "pilot ", &
492 "satem ", &
493 "geoamv ", &
494 "polaramv ", &
495 "airep ", &
496 "gpspw ", &
497 "gpsrf ", &
498 "metar ", &
499 "ships ", &
500 "ssmi_rv ", &
501 "ssmi_tb ", &
502 "ssmt1 ", &
503 "ssmt2 ", &
504 "qscat ", &
505 "profiler ", &
506 "buoy ", &
507 "bogus ", &
508 "pseudo ", &
509 "radar ", &
510 "radiance ", &
511 "airs retrieval", &
512 "sonde_sfc ", &
513 "mtgirs ", &
514 "tamdar ", &
515 "tamdar_sfc " &
519 integer, parameter :: max_no_fm = 290
521 integer, parameter :: num_ob_vars=9
523 logical, parameter :: in_report(num_ob_vars,2) = reshape((/&
524 .false.,.false.,.false.,.false.,.false.,.false.,.false.,.false.,.false., & ! sound
525 .true.,.true.,.true.,.true.,.true.,.true.,.false.,.false.,.false./), &
526 (/num_ob_vars,2/))
528 integer, parameter :: report_h = 1
529 integer, parameter :: report_u = 2
530 integer, parameter :: report_v = 3
531 integer, parameter :: report_t = 4
532 integer, parameter :: report_q = 5
533 integer, parameter :: report_p = 6
534 integer, parameter :: report_rh = 7
535 integer, parameter :: report_slp = 8
536 integer, parameter :: report_zk = 9
538 logical :: obs_use(num_ob_indexes) = .false.
540 ! Special cases
542 integer, parameter :: fm_satem = 86
543 integer, parameter :: fm_amv = 88
545 integer, parameter :: fm_index(max_no_fm) = (/ &
546 0,0,0,0,0,0,0,0,0,0, & ! 1-10
547 0,Synop,Ships,0,Metar, & ! 11-15
548 Metar,Ships,buoy,buoy,0, & ! 16-20
549 0,0,0,0,0,0,0,0,0,0, & ! 21-30
550 0,pilot,pilot,pilot,sound, & ! 31-35
551 sound,sound,sound,0,0, & ! 36-40
552 0,airep,0,0,0,0,0,0,0,0, & ! 41-50
553 0,0,0,0,0,0,0,0,0,0, & ! 51-60
554 0,0,0,0,0,0,0,0,0,0, & ! 61-70
555 0,0,0,0,0,0,0,0,0,0, & ! 71-80
556 0,0,0,0,0,satem,0,geoamv,0,0, & ! 81-90
557 0,0,0,0,0,airep,airep,0,0,0, & ! 91-100
558 tamdar,0,0,0,0,0,0,0,0,0, & ! 101-110
559 gpspw,0,0,gpspw,0,gpsref,0,0,0,0, & ! 111-120
560 ssmt1,ssmt2,0,0,ssmi_rv,0,0,0,0,0, & ! 121-130
561 0,profiler,airsr,0,bogus,0,0,0,0,0, & ! 131-140
562 0,0,0,0,0,0,0,0,0,0, & ! 141-150
563 0,0,0,0,0,0,0,0,0,0, & ! 151-160
564 mtgirs,0,0,0,0,0,0,0,0,0, & ! 161-170
565 0,0,0,0,0,0,0,0,0,0, & ! 171-180
566 0,0,0,0,0,0,0,0,0,0, & ! 181-190
567 0,0,0,0,0,0,0,0,0,0, & ! 191-200
568 0,0,0,0,0,0,0,0,0,0, & ! 201-210
569 0,0,0,0,0,0,0,0,0,0, & ! 211-220
570 0,0,0,0,0,0,0,0,0,0, & ! 231-230
571 0,0,0,0,0,0,0,0,0,0, & ! 231-240
572 0,0,0,0,0,0,0,0,0,0, & ! 241-250
573 0,0,0,0,0,0,0,0,0,0, & ! 251-260
574 0,0,0,0,0,0,0,0,0,0, & ! 261-270
575 0,0,0,0,0,0,0,0,0,0, & ! 271-280
576 qscat,0,0,0,0,0,0,0,0,0 /) ! 281-290
578 character(len=120) :: fmt_info ='(a12,1x,a19,1x,a40,1x,i6,3(f12.3,11x),6x,a5)'
579 character(len=120) :: fmt_srfc = '(7(:,f12.3,i4,f7.2))'
580 ! character(len=120) :: fmt_srfc = '(f12.3,i4,f7.2,F12.3,I4,F7.3)'
581 character(len=120) :: fmt_each = &
582 '(3(f12.3,i4,f7.2),11x,3(f12.3,i4,f7.2),11x,3(f12.3,i4,f7.2))'
584 ! lat/long information calculated in da_setup_firstguess_wrf
586 real, parameter :: deg_to_rad = pi/180.0
587 real, parameter :: rad_to_deg = 1.0/deg_to_rad
589 real, allocatable :: cos_xls(:)
590 real, allocatable :: sin_xls(:)
591 real, allocatable :: cos_xle(:)
592 real, allocatable :: sin_xle(:)
594 integer :: ierr ! General error code
595 integer :: comm ! MPI communicator
597 integer :: ids,ide,jds,jde,kds,kde
598 integer :: ims,ime,jms,jme,kms,kme
599 integer :: its,ite,jts,jte,kts,kte
600 integer :: ips,ipe,jps,jpe,kps,kpe
601 integer :: itsy,itey,jtsy,jtey,ktsy,ktey
602 integer :: itsx,itex,jtsx,jtex,ktsx,ktex
604 integer :: num_qcstat_conv(2,num_ob_indexes,num_ob_vars,npres_print)
605 character*4, parameter :: ob_vars(num_ob_vars) = (/'U ','V ','T ',&
606 'Q ','Ps ','Spd ',&
607 'Tpw ','GpsR','Thic'/)
608 real, parameter :: pptop(1:npres_print) = (/ 1000.0, 900.0, 800.0, 600.0, 400.0, 300.0, &
609 250.0, 200.0, 150.0, 100.0, 50.0, 0./)
611 real, parameter :: ppbot(npres_print) = (/ 1200.0, 999.9, 899.9, 799.0, 599.9, 399.9, &
612 299.9, 249.9, 199.9, 149.9, 99.9, 2000./)
614 real, allocatable :: time_slots(:)
616 logical :: global
618 end module da_control