1program midas_genCoeff
2 !
3 !:Purpose: Main program to compute radiance bias correction coefficients by linear regression.
4 !
5 ! ---
6 !
7 !:Algorithm: O-A are computed from input analysis fields resulting from a stand alone
8 ! 3DVar analysis assimilating only anchoring observations (considered as "trials")
9 ! and bgckalt files.
10 !
11 ! --
12 !
13 !:File I/O: The required input files and produced output files are listed as follows.
14 !
15 !============================================== ==============================================================
16 ! Input and Output Files (NWP application) Description of file
17 !============================================== ==============================================================
18 ! ``flnml`` In - Main namelist file with parameters user may modify
19 ! ``flnml_static`` In - The "static" namelist that should not be modified
20 ! ``trlm_$NN`` (e.g. ``trlm_01``) In - Background state (a.k.a. trial) files for each timestep
21 ! ``analysisgrid`` In - File defining grid for computing the analysis increment
22 ! ``obsfiles_$FAM/obs$FAM_$NNNN_$NNNN`` In - Observation file for each "family" and MPI task
23 ! Remainder are files related to radiance obs:
24 ! ``stats_tovs`` In - Observation error file for radiances
25 ! ``stats_tovs_symmetricObsErr`` In - user-defined symmetric TOVS errors for all sky
26 ! ``Cmat_$PLATFORM_$SENSOR.dat`` In - Inter-channel observation-error correlations
27 ! ``rtcoef_$PLATFORM_$SENSOR.H5`` In - RTTOV coefficient file HDF-5 format
28 ! ``rtcoef_$PLATFORM_$SENSOR.dat`` In - RTTOV coefficient file ASCII format
29 ! ``ozoneclim98`` In - ozone climatology standard file (Fortuin and Kelder)
30 !============================================== ==============================================================
31 !
32 ! --
33 !
34 !:Synopsis: Below is a summary of the ``genCoeff`` program calling sequence:
35 !
36 ! - **Initial setups:**
37 !
38 ! - Setup time grid using initial time from trlm_01
39 !
40 ! - Setup horizontal and vertical grid objects for "analysis
41 ! grid" from ``analysisgrid`` file.
42 !
43 ! - Setup ``obsSpaceData`` object and read observations from
44 ! files: ``inn_setupObs``.
45 !
46 ! - Setup ``columnData`` and ``gridStateVector`` modules (read
47 ! list of analysis variables from namelist) and allocate column
48 ! object for storing trial on analysis levels.
49 !
50 ! - Setup the observation error statistics in ``obsSpaceData``
51 ! object: ``oer_setObsErrors``.
52 !
53 ! - Allocate a stateVector object on the trial grid and then
54 ! read the trials: ``gio_readTrials``.
55 !
56 ! - **Coefficients computation**
57 !
58 ! - Horizontally Interpolate "trial" fields to trial columns
59 !
60 ! - if needed substract bias correction from ObsSpaceData to get raw O-F and obs.
61 !
62 ! - Remove outliers
63 !
64 ! - Compute innovation from "trial" fields
65 !
66 ! - Refresh Bias Correction (probably useless)
67 !
68 ! - Perform linear regression
69 !
70 ! - write coefficients to output file.
71 !
72 ! - if requested compute and output to file raw (.i.e without bias correction) O-F statistics
73 !
74 ! - compute and apply bias coorection to obsSpace Data
75 !
76 ! - if requested compute and output to file bias corrected O-F statistics
77 !
78 ! - **Final step**
79 !
80 ! - deallocate memory (radiance bias correction module and obsSpaceData)
81 !
82 ! --
83 !
84 !:Options: `List of namelist blocks <../namelists_in_each_program.html#genCoeff>`_
85 ! that can affect the ``genCoeff`` program.
86 !
87
88 use version_mod
89 use codePrecision_mod
90 use ramDisk_mod
91 use utilities_mod
92 use midasMpi_mod
93 use mathPhysConstants_mod
94 use horizontalCoord_mod
95 use verticalCoord_mod
96 use timeCoord_mod
97 use obsSpaceData_mod
98 use columnData_mod
99 use gridStateVector_mod
100 use gridStateVectorFileIO_mod
101 use obsFiles_mod
102 use innovation_mod
103 use obsErrors_mod
104 use biasCorrectionSat_mod
105
106 implicit none
107
108 integer, external :: exdb,exfin,fnom, fclos, get_max_rss
109 integer :: ierr,istamp
110
111 type(struct_obs), target :: obsSpaceData
112 type(struct_columnData), target :: columnTrlOnAnlIncLev
113 type(struct_gsv) :: stateVectorTrialHighRes
114 type(struct_hco), pointer :: hco_anl => null()
115 type(struct_hco), pointer :: hco_core => null()
116 type(struct_vco), pointer :: vco_anl => null()
117 type(struct_hco), pointer :: hco_trl => null()
118 type(struct_vco), pointer :: vco_trl => null()
119
120 logical :: allocHeightSfc
121
122 character(len=48), parameter :: obsMpiStrategy = 'LIKESPLITFILES', &
123 varMode = 'analysis'
124
125 istamp = exdb('GENCOEFF','DEBUT','NON')
126
127 call ver_printNameAndVersion('genCoeff','Bias Correction Coefficient Computation')
128
129 ! MPI initialization
130 call mmpi_initialize
131
132 call tmg_init(mmpi_myid, 'TMG_INFO')
133
134 call utl_tmg_start(0,'Main')
135
136
137 ! 1. Top level setup
138 call ram_setup()
139
140 ! Do initial set up
141 call gencoeff_setup('VAR') ! obsColumnMode
142
143 ! Reading trials
144 call inn_getHcoVcoFromTrlmFile( hco_trl, vco_trl )
145 allocHeightSfc = ( vco_trl%Vcode /= 0 )
146
147 call gsv_allocate( stateVectorTrialHighRes, tim_nstepobs, hco_trl, vco_trl, &
148 dateStamp_opt=tim_getDateStamp(), mpi_local_opt=.true., &
149 mpi_distribution_opt='Tiles', dataKind_opt=4, &
150 allocHeightSfc_opt=allocHeightSfc, hInterpolateDegree_opt='LINEAR', &
151 beSilent_opt=.false. )
152 call gsv_zero( stateVectorTrialHighRes )
153 call gio_readTrials( stateVectorTrialHighRes )
154 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
155
156 ! Horizontally interpolate trials to trial columns
157 call inn_setupColumnsOnTrlLev( columnTrlOnAnlIncLev, obsSpaceData, hco_core, &
158 stateVectorTrialHighRes )
159 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
160
161 call utl_tmg_start(110,'--BiasCorrection')
162
163 ! Remove bias correction if requested
164 call bcs_removeBiasCorrection(obsSpaceData,"TO")
165
166 call bcs_removeOutliers(obsSpaceData)
167
168 call utl_tmg_stop(110)
169
170 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
171
172 ! Compute observation innovations
173 call inn_computeInnovation(columnTrlOnAnlIncLev,obsSpaceData)
174
175 call utl_tmg_start(110,'--BiasCorrection')
176
177 ! Refresh bias correction if requested
178 call bcs_refreshBiasCorrection(obsSpaceData,columnTrlOnAnlIncLev)
179
180 call utl_tmg_start(111,'----Regression')
181 call bcs_do_regression(columnTrlOnAnlIncLev,obsSpaceData)
182 call utl_tmg_stop(111)
183
184 ! Write coefficients to file
185 call bcs_writebias()
186
187 ! output O-F statistics befor bias correction
188 call bcs_computeResidualsStatistics(obsSpaceData,"_raw")
189
190 ! fill OBS_BCOR with computed bias correction
191 call bcs_calcBias(obsSpaceData,columnTrlOnAnlIncLev)
192
193 ! output O-F statistics after bias coorection
194 call bcs_computeResidualsStatistics(obsSpaceData,"_corrected")
195
196 ! Deallocate internal bias correction structures
197 call bcs_finalize()
198
199 call utl_tmg_stop(110)
200
201 ! Deallocate copied obsSpaceData
202 call obs_finalize(obsSpaceData)
203
204 ! 3. Job termination
205
206 istamp = exfin('GENCOEFF','FIN','NON')
207
208 call utl_tmg_stop(0)
209
210 call tmg_terminate(mmpi_myid, 'TMG_INFO')
211
212 call rpn_comm_finalize(ierr)
213
214contains
215
216 subroutine gencoeff_setup(obsColumnMode)
217 !
218 ! :Purpose: Control of the preprocessing of bais correction coefficient computation
219 !
220 implicit none
221
222 ! Arguments:
223 character(len=*), intent(in) :: obsColumnMode
224
225 ! Locals:
226 integer :: dateStamp, dateStampFromObs
227
228 write(*,*) ''
229 write(*,*) '----------------------------------------'
230 write(*,*) '-- Starting subroutine gencoeff_setup --'
231 write(*,*) '----------------------------------------'
232
233 !
234 !- Initialize the Temporal grid
235 !
236 call tim_setup
237
238 !
239 !- Initialize observation file names and set datestamp from trial file
240 !
241 call obsf_setup(dateStampFromObs, varMode)
242
243 dateStamp = tim_getDatestampFromFile('./trlm_01')
244
245 if (dateStamp > 0) then
246 call tim_setDatestamp(datestamp)
247 else
248 call utl_abort('var_setup: Problem getting dateStamp from first trial field')
249 end if
250
251 !
252 !- Initialize constants
253 !
254 if ( mmpi_myid == 0 ) then
255 call mpc_printConstants(6)
256 call pre_printPrecisions
257 end if
258
259 !
260 !- Initialize variables of the model states
261 !
262 call gsv_setup
263 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
264 !
265 !- Initialize the Analysis grid
266 !
267 if(mmpi_myid == 0) write(*,*)''
268 if(mmpi_myid == 0) write(*,*)'gencoeff_setup: Set hco parameters for analysis grid'
269 call hco_SetupFromFile(hco_anl, './analysisgrid', 'ANALYSIS', 'Analysis' ) ! IN
270
271 if ( hco_anl % global ) then
272 hco_core => hco_anl
273 else
274 !- Initialize the core (Non-Extended) analysis grid
275 if(mmpi_myid == 0) write(*,*)'gencoeff_setup: Set hco parameters for core grid'
276 call hco_SetupFromFile( hco_core, './analysisgrid', 'COREGRID', 'AnalysisCore' ) ! IN
277 end if
278
279 !
280 !- Initialisation of the analysis grid vertical coordinate from analysisgrid file
281 !
282 call vco_SetupFromFile( vco_anl, & ! OUT
283 './analysisgrid') ! IN
284
285 call col_setVco(columnTrlOnAnlIncLev,vco_anl)
286 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
287
288 !
289 !- Setup and read observations
290 !
291 call inn_setupObs(obsSpaceData, hco_anl, obsColumnMode, obsMpiStrategy, varMode) ! IN
292 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
293
294 !
295 !- Basic setup of columnData module
296 !
297 call col_setup
298 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
299
300 !
301 !- Initialize the observation error covariances
302 !
303 if (.not. bcs_mimicSatbcor) call oer_setObsErrors(obsSpaceData, varMode) ! IN
304 write(*,*) 'Memory Used: ',get_max_rss()/1024,'Mb'
305
306 end subroutine gencoeff_setup
307
308
309
310end program midas_genCoeff