1program midas_analysisErrorOI
2 !
3 !:Purpose: Calculate analysis-error standard deviation given
4 ! new assimilated observations. It only works for sea ice variables and
5 ! uses a simple OI approach.
6 !
7 ! ---
8 !
9 !:Algorithm: The Optimal Interpolation (OI) is a data assimilation approach
10 ! where both the state and its estimated error are computed using
11 ! observations while taking into account the specified
12 ! uncertainties for both the observations and the background
13 ! state (i.e. the R and B covariance matrices, respectively).
14 ! The computations are done independently at each analysis grid
15 ! point with only local observations, those that have a
16 ! significant influence on the analysis at the grid point
17 ! location. Here the code only implements the calculation to
18 ! update the diagonal of the B covariance matrix.
19 !
20 ! --
21 !
22 !============================================== ==============================================================
23 ! Input and Output Files Description of file
24 !============================================== ==============================================================
25 ! ``flnml`` In - Main namelist file with parameters user may modify
26 ! ``trlm_01`` In - Background-error standard deviation
27 ! ``analysisgrid`` In - File defining grid for computing the analysis error
28 ! ``sea_ice_obs-err`` In - Observation error statistics
29 ! ``bgstddev`` In - Static background-error statistics
30 ! ``bgSeaIceConc`` In - Background sea ice concentration
31 ! ``obsfiles_$FAM/obs$FAM_0001_0001`` In - Observation file for each "family" (only 1 MPI task)
32 ! ``anlm_000m`` Out - Analysis-error on the analysis grid
33 ! ``obsfiles_$FAM.updated/obs$FAM_0001_0001`` Out - Updated obs file for each "family" (only 1 MPI task)
34 !============================================== ==============================================================
35 !
36 ! --
37 !
38 !:Synopsis: Below is a summary of the ``analysisErrorOI`` program calling sequence:
39 !
40 ! - **Initial setups:**
41 !
42 ! - Setup horizontal and vertical grid objects for "analysis
43 ! grid" from ``trlm_01`` file.
44 !
45 ! - Setup ``obsSpaceData`` object and read observations from
46 ! files: ``inn_setupObs``.
47 !
48 ! - Setup ``columnData`` and ``gridStateVector`` modules (read
49 ! list of analysis variables from namelist) and allocate column
50 ! object for storing trial on analysis levels.
51 !
52 ! - Setup the observation error statistics in ``obsSpaceData``
53 ! object: ``oer_setObsErrors``.
54 !
55 ! - Filter out observations from satellites
56 ! not specified in the name list: ``filt_iceConcentration``.
57 !
58 ! - Filter scatterometer backscatter anisotropy observations
59 ! where wind speed is too small: ``filt_backScatAnisIce``.
60 !
61 ! - Setup observation-error for scatterometer backscatter
62 ! anisotropy observations: ``oer_setErrBackScatAnisIce``.
63 !
64 ! - **Main calculation:**
65 !
66 ! - Compute the analysis-error: ``aer_analysisError``.
67 !
68 ! - Update the Days Since Last Obs: ``aer_daysSinceLastObs``.
69 !
70 ! - Update the observation files: ``obsf_writeFiles``.
71 !
72 ! --
73 !
74 !:Options: `List of namelist blocks <../namelists_in_each_program.html#analysiserroroi>`_
75 ! that can affect the ``analysisErrorOI`` program.
76 !
77 ! * The relevant namelist blocks used to configure the
78 ! analysis-error calculation are listed in the following table:
79 !
80 !======================== ================== ==============================================================
81 ! Module Namelist Description of what is controlled
82 !======================== ================== ==============================================================
83 ! ``timeCoord_mod`` ``NAMTIME`` assimilation time window length, temporal resolution of
84 ! the background state
85 ! ``columndata_mod`` ``NAMSTATE`` name of the analysis variable (RPN nomvar, 4-character long),
86 ! only sea ice concentration (GL) is allowed for now
87 ! ``gridstatevector_mod`` " "
88 ! ``obsspacedata_mod`` ``NAMDIMO`` specify the maximum number of header and body elements
89 ! ``obsfilter_mod`` ``NAMFILT`` list of varno to use and bit flags (13-bit#) for filtering
90 ! ``obsfilter_mod`` ``namPlatformIce`` list of observation platforms to assimilate
91 ! ``sqliteread_mod`` ``NAMSQLgl`` list of varno to read from SQLite files for sea ice
92 ! ``sqliteread_mod`` ``namSQLUpdate`` list of elements to update in SQLite files
93 ! ``sqliteread_mod`` ``namSQLInsert`` place holder, could be empty
94 ! ``analysisErrorOI_mod`` ``NAMAER`` set the maximum analysis-error std dev. allowed
95 !======================== ================== ==============================================================
96 !
97 use version_mod
98 use ramDisk_mod
99 use utilities_mod
100 use midasMpi_mod
101 use message_mod
102 use mathPhysConstants_mod
103 use horizontalCoord_mod
104 use verticalCoord_mod
105 use timeCoord_mod
106 use obsSpaceData_mod
107 use columnData_mod
108 use gridStateVector_mod
109 use obsFiles_mod
110 use innovation_mod
111 use obsErrors_mod
112 use obsFilter_mod
113 use analysisErrorOI_mod
114
115 implicit none
116
117 integer :: istamp, exdb, exfin
118 integer :: ierr, dateStampFromObs
119 character(len=48) :: obsMpiStrategy, varMode
120 character(len=*), parameter :: myName = 'analysisErrorOI'
121
122 type(struct_obs) , target :: obsSpaceData
123 type(struct_columnData), target :: trlColumnOnAnlLev
124 type(struct_hco) , pointer :: hco_anl => null()
125 type(struct_vco) , pointer :: vco_anl => null()
126
127 istamp = exdb('ANALYSISERROROI','DEBUT','NON')
128
129 call ver_printNameAndVersion('analysisErrorOI', 'Program to calculate the analysis-error standard deviation using OI.')
130
131 ! MPI initialization
132 call mmpi_initialize
133
134 call tmg_init(mmpi_myid, 'TMG_INFO')
135
136 call utl_tmg_start(0,'Main')
137
138 varMode='analysis'
139
140 ! Setup the ram disk
141 call ram_setup
142
143 obsMpiStrategy = 'LIKESPLITFILES'
144
145 !
146 !- Initialize the Temporal grid and set dateStamp from env variable
147 !
148 call tim_setup()
149
150 if (tim_nstepobs > 1 .or. tim_nstepobsinc > 1) then
151 call utl_abort('analysisErrorOI: The program assumes only one time step.')
152 end if
153
154 !
155 !- Initialize observation file names and set datestamp if not already
156 !
157 call obsf_setup(dateStampFromObs, varMode)
158 if (tim_getDateStamp() == 0) then
159 if (dateStampFromObs > 0) then
160 call tim_setDateStamp(dateStampFromObs)
161 else
162 call utl_abort('analysisErrorOI: DateStamp was not set')
163 end if
164 end if
165
166 !
167 !- Initialize constants
168 !
169 if (mmpi_myid == 0) call mpc_printConstants(6)
170
171 !
172 !- Initialize list of analyzed variables.
173 !
174 call gsv_setup
175 call msg_memUsage(myName)
176
177 !
178 !- Initialize the Analysis grid
179 !
180 call msg(myName,'Set hco parameters for analysis grid', mpiAll_opt=.false.)
181 call hco_SetupFromFile(hco_anl, './analysisgrid', '') ! IN
182
183 !
184 !- Initialisation of the analysis grid vertical coordinate from analysisgrid file
185 !
186 call vco_SetupFromFile(vco_anl, & ! OUT
187 './analysisgrid') ! IN
188
189 call col_setVco(trlColumnOnAnlLev,vco_anl)
190 call msg_memUsage(myName)
191
192 !
193 !- Setup and read observations
194 !
195 call inn_setupObs(obsSpaceData, hco_anl, 'VAR', obsMpiStrategy, varMode) ! IN
196 call msg_memUsage(myName)
197
198 !
199 !- Basic setup of columnData module
200 !
201 call col_setup
202 call msg_memUsage(myName)
203
204 !
205 !- Memory allocation for background column data
206 !
207 call col_allocate(trlColumnOnAnlLev, obs_numheader(obsSpaceData), mpiLocal_opt = .true.)
208
209 !
210 !- Initialize the observation error covariances
211 !
212 call oer_setObsErrors(obsSpaceData, varMode) ! IN
213 call msg_memUsage(myName)
214
215 ! Sea ice concentration
216 if (obs_famExist(obsSpaceData, 'GL')) then
217 call filt_iceConcentration(obsSpaceData, beSilent = .false.)
218 call filt_backScatAnisIce(obsSpaceData, beSilent = .false.)
219 call oer_setErrBackScatAnisIce(obsSpaceData, beSilent = .false.)
220 end if
221
222 ! Compute the analysis error
223 call aer_analysisError(obsSpaceData, hco_anl, vco_anl)
224
225 ! Now write out the observation data files
226 if ( .not. obsf_filesSplit() ) then
227 call msg(myName,'reading/writing global observation files')
228 call obs_expandToMpiGlobal(obsSpaceData)
229 if (mmpi_myid == 0) call obsf_writeFiles(obsSpaceData)
230 else
231 ! redistribute obs data to how it was just after reading the files
232 call obs_MpiRedistribute(obsSpaceData,OBS_IPF)
233 call obsf_writeFiles(obsSpaceData)
234 end if
235
236 ! Deallocate copied obsSpaceData
237 call obs_finalize(obsSpaceData)
238
239 !
240 ! 3. Job termination
241 !
242 istamp = exfin('ANALYSISERROROI','FIN','NON')
243
244 call utl_tmg_stop(0)
245
246 call tmg_terminate(mmpi_myid, 'TMG_INFO')
247
248 call rpn_comm_finalize(ierr)
249
250end program midas_analysisErrorOI