1 | % ================================================================ |
---|
2 | % Chapter observation operator (OBS) |
---|
3 | % ================================================================ |
---|
4 | \chapter{Observation and model comparison (OBS)} |
---|
5 | \label{OBS} |
---|
6 | |
---|
7 | Authors: D. Lea, K. Mogensen, A. Weaver, M. Martin ... % do we keep that ? |
---|
8 | |
---|
9 | \minitoc |
---|
10 | |
---|
11 | |
---|
12 | \newpage |
---|
13 | $\ $\newline % force a new line |
---|
14 | |
---|
15 | The OBS branch is a diagnostic branch which reads in observation files (profile temperature and |
---|
16 | salinity, sea surface temperature, sea level anomaly and sea ice concentration) and calculates |
---|
17 | an interpolated model equivalent value at the observation location and nearest model timestep. |
---|
18 | This is code with was originally developed for use with NEMOVAR. |
---|
19 | |
---|
20 | In the case of temperature data from moored buoys (TAO, TRITON, PIRATA) which in the |
---|
21 | ENACT/ENSEMBLES data-base are available as daily averaged quantities the observation operator |
---|
22 | averages the model temperature fields over one day before interpolating them to the observation |
---|
23 | location. For SST and altimeter observations, a 2D horizontal interpolator is needed. For {\em in |
---|
24 | situ} profiles, a 1D vertical interpolator is needed in addition to the 2D interpolator. |
---|
25 | |
---|
26 | The resulting data is saved in a ``feedback'' file or files which can be used for model validation |
---|
27 | and verification and also to provide information for data assimilation. This code is controlled by |
---|
28 | the namelist \textit{nam\_obs}. To build with the OBS code active \key{diaobs} must be set. |
---|
29 | |
---|
30 | There is a brief description of all the namelist options provided. |
---|
31 | |
---|
32 | Missing information: description of \key{sp}, \key{datetime\_out} |
---|
33 | |
---|
34 | |
---|
35 | % ================================================================ |
---|
36 | % Example |
---|
37 | % ================================================================ |
---|
38 | \section{Running the observation operator code example} |
---|
39 | \label{OBS_example} |
---|
40 | |
---|
41 | This section describes an example of running the observation operator code using |
---|
42 | profile data which can be freely downloaded. This shows how to adapt an existing |
---|
43 | run and build of NEMO to run the observation operator. |
---|
44 | |
---|
45 | First compile NEMO with \key{diaobs} set. |
---|
46 | |
---|
47 | Next download some ENSEMBLES EN3 data from the website http://www.hadobs.org. |
---|
48 | You should choose observations which are valid for the period of your test run |
---|
49 | because the observation operator compares the model and observations for a |
---|
50 | matching date and time. |
---|
51 | |
---|
52 | You will need to add the following to the namelist to run the observation |
---|
53 | operator on this data - set the \np{enactfiles} namelist parameter to the observation |
---|
54 | file name you wish to use (or link in): |
---|
55 | |
---|
56 | %------------------------------------------namobs_example----------------------------------------------------- |
---|
57 | \namdisplay{namobs_example} |
---|
58 | %------------------------------------------------------------------------------------------------------------- |
---|
59 | |
---|
60 | The option \np{ln\_t3d} and \np{ln\_s3d} switch on the temperature and salinity profile |
---|
61 | observation operator code. \np{ln\_ena} switch turns on the reading on ENACT/ENSEMBLES type |
---|
62 | profile data. The filename of the ENACT/ENSEMBLES data is set with enactfiles. Not this can be |
---|
63 | array of multiple files. The observation operator code needs to convert from observation |
---|
64 | latitude and longitude to model grid point. This is done using the grid searching part of the |
---|
65 | code. Setting \np{ln\_grid\_global} means that the code distributes the observations evenly |
---|
66 | between processors in a round-robin. Alternatively each processor will work with observations |
---|
67 | located within the model subdomain. The grid searching can be expensive, particularly for large |
---|
68 | numbers of observations, setting \np{ln\_grid\_search\_lookup} allows the use of a lookup |
---|
69 | table. This will need to generated the first time if it does not exist in the run directory |
---|
70 | however once produced it will significantly speed up future grid searches. |
---|
71 | |
---|
72 | The next step is viewing and working with the data. The NEMOVAR system contains utilities to |
---|
73 | plot the feedback files, convert and recombine the files which are available on request from |
---|
74 | the NEMOVAR team. |
---|
75 | |
---|
76 | \section{Technical details} |
---|
77 | \label{OBS_details} |
---|
78 | |
---|
79 | Here we show a more complete example namelist and also describe the observation files that may |
---|
80 | be used with the observation operator |
---|
81 | |
---|
82 | %------------------------------------------namobs-------------------------------------------------------- |
---|
83 | \namdisplay{namobs} |
---|
84 | %------------------------------------------------------------------------------------------------------------- |
---|
85 | |
---|
86 | This name list uses the "feedback" type observation file input format for profile, sea level |
---|
87 | anomaly and sea surface temperature data |
---|
88 | |
---|
89 | \subsection{Profile feedback type observation file header} |
---|
90 | |
---|
91 | \begin{alltt} |
---|
92 | \tiny |
---|
93 | \begin{verbatim} |
---|
94 | netcdf profiles_01 { |
---|
95 | dimensions: |
---|
96 | N_OBS = 603 ; |
---|
97 | N_LEVELS = 150 ; |
---|
98 | N_VARS = 2 ; |
---|
99 | N_QCF = 2 ; |
---|
100 | N_ENTRIES = 1 ; |
---|
101 | N_EXTRA = 1 ; |
---|
102 | STRINGNAM = 8 ; |
---|
103 | STRINGGRID = 1 ; |
---|
104 | STRINGWMO = 8 ; |
---|
105 | STRINGTYP = 4 ; |
---|
106 | STRINGJULD = 14 ; |
---|
107 | variables: |
---|
108 | char VARIABLES(N_VARS, STRINGNAM) ; |
---|
109 | VARIABLES:long_name = "List of variables in feedback files" ; |
---|
110 | char ENTRIES(N_ENTRIES, STRINGNAM) ; |
---|
111 | ENTRIES:long_name = "List of additional entries for each variable in feedback files" ; |
---|
112 | char EXTRA(N_EXTRA, STRINGNAM) ; |
---|
113 | EXTRA:long_name = "List of extra variables" ; |
---|
114 | char STATION_IDENTIFIER(N_OBS, STRINGWMO) ; |
---|
115 | STATION_IDENTIFIER:long_name = "Station identifier" ; |
---|
116 | char STATION_TYPE(N_OBS, STRINGTYP) ; |
---|
117 | STATION_TYPE:long_name = "Code instrument type" ; |
---|
118 | double LONGITUDE(N_OBS) ; |
---|
119 | LONGITUDE:long_name = "Longitude" ; |
---|
120 | LONGITUDE:units = "degrees_east" ; |
---|
121 | LONGITUDE:_Fillvalue = 99999.f ; |
---|
122 | double LATITUDE(N_OBS) ; |
---|
123 | LATITUDE:long_name = "Latitude" ; |
---|
124 | LATITUDE:units = "degrees_north" ; |
---|
125 | LATITUDE:_Fillvalue = 99999.f ; |
---|
126 | double DEPTH(N_OBS, N_LEVELS) ; |
---|
127 | DEPTH:long_name = "Depth" ; |
---|
128 | DEPTH:units = "metre" ; |
---|
129 | DEPTH:_Fillvalue = 99999.f ; |
---|
130 | int DEPTH_QC(N_OBS, N_LEVELS) ; |
---|
131 | DEPTH_QC:long_name = "Quality on depth" ; |
---|
132 | DEPTH_QC:Conventions = "q where q =[0,9]" ; |
---|
133 | DEPTH_QC:_Fillvalue = 0 ; |
---|
134 | int DEPTH_QC_FLAGS(N_OBS, N_LEVELS, N_QCF) ; |
---|
135 | DEPTH_QC_FLAGS:long_name = "Quality flags on depth" ; |
---|
136 | DEPTH_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
137 | double JULD(N_OBS) ; |
---|
138 | JULD:long_name = "Julian day" ; |
---|
139 | JULD:units = "days since JULD_REFERENCE" ; |
---|
140 | JULD:Conventions = "relative julian days with decimal part (as parts of day)" ; |
---|
141 | JULD:_Fillvalue = 99999.f ; |
---|
142 | char JULD_REFERENCE(STRINGJULD) ; |
---|
143 | JULD_REFERENCE:long_name = "Date of reference for julian days" ; |
---|
144 | JULD_REFERENCE:Conventions = "YYYYMMDDHHMMSS" ; |
---|
145 | int OBSERVATION_QC(N_OBS) ; |
---|
146 | OBSERVATION_QC:long_name = "Quality on observation" ; |
---|
147 | OBSERVATION_QC:Conventions = "q where q =[0,9]" ; |
---|
148 | OBSERVATION_QC:_Fillvalue = 0 ; |
---|
149 | int OBSERVATION_QC_FLAGS(N_OBS, N_QCF) ; |
---|
150 | OBSERVATION_QC_FLAGS:long_name = "Quality flags on observation" ; |
---|
151 | OBSERVATION_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
152 | OBSERVATION_QC_FLAGS:_Fillvalue = 0 ; |
---|
153 | int POSITION_QC(N_OBS) ; |
---|
154 | POSITION_QC:long_name = "Quality on position (latitude and longitude)" ; |
---|
155 | POSITION_QC:Conventions = "q where q =[0,9]" ; |
---|
156 | POSITION_QC:_Fillvalue = 0 ; |
---|
157 | int POSITION_QC_FLAGS(N_OBS, N_QCF) ; |
---|
158 | POSITION_QC_FLAGS:long_name = "Quality flags on position" ; |
---|
159 | POSITION_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
160 | POSITION_QC_FLAGS:_Fillvalue = 0 ; |
---|
161 | int JULD_QC(N_OBS) ; |
---|
162 | JULD_QC:long_name = "Quality on date and time" ; |
---|
163 | JULD_QC:Conventions = "q where q =[0,9]" ; |
---|
164 | JULD_QC:_Fillvalue = 0 ; |
---|
165 | int JULD_QC_FLAGS(N_OBS, N_QCF) ; |
---|
166 | JULD_QC_FLAGS:long_name = "Quality flags on date and time" ; |
---|
167 | JULD_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
168 | JULD_QC_FLAGS:_Fillvalue = 0 ; |
---|
169 | int ORIGINAL_FILE_INDEX(N_OBS) ; |
---|
170 | ORIGINAL_FILE_INDEX:long_name = "Index in original data file" ; |
---|
171 | ORIGINAL_FILE_INDEX:_Fillvalue = -99999 ; |
---|
172 | float POTM_OBS(N_OBS, N_LEVELS) ; |
---|
173 | POTM_OBS:long_name = "Potential temperature" ; |
---|
174 | POTM_OBS:units = "Degrees Celsius" ; |
---|
175 | POTM_OBS:_Fillvalue = 99999.f ; |
---|
176 | float POTM_Hx(N_OBS, N_LEVELS) ; |
---|
177 | POTM_Hx:long_name = "Model interpolated potential temperature" ; |
---|
178 | POTM_Hx:units = "Degrees Celsius" ; |
---|
179 | POTM_Hx:_Fillvalue = 99999.f ; |
---|
180 | int POTM_QC(N_OBS) ; |
---|
181 | POTM_QC:long_name = "Quality on potential temperature" ; |
---|
182 | POTM_QC:Conventions = "q where q =[0,9]" ; |
---|
183 | POTM_QC:_Fillvalue = 0 ; |
---|
184 | int POTM_QC_FLAGS(N_OBS, N_QCF) ; |
---|
185 | POTM_QC_FLAGS:long_name = "Quality flags on potential temperature" ; |
---|
186 | POTM_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
187 | POTM_QC_FLAGS:_Fillvalue = 0 ; |
---|
188 | int POTM_LEVEL_QC(N_OBS, N_LEVELS) ; |
---|
189 | POTM_LEVEL_QC:long_name = "Quality for each level on potential temperature" ; |
---|
190 | POTM_LEVEL_QC:Conventions = "q where q =[0,9]" ; |
---|
191 | POTM_LEVEL_QC:_Fillvalue = 0 ; |
---|
192 | int POTM_LEVEL_QC_FLAGS(N_OBS, N_LEVELS, N_QCF) ; |
---|
193 | POTM_LEVEL_QC_FLAGS:long_name = "Quality flags for each level on potential temperature" ; |
---|
194 | POTM_LEVEL_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
195 | POTM_LEVEL_QC_FLAGS:_Fillvalue = 0 ; |
---|
196 | int POTM_IOBSI(N_OBS) ; |
---|
197 | POTM_IOBSI:long_name = "ORCA grid search I coordinate" ; |
---|
198 | int POTM_IOBSJ(N_OBS) ; |
---|
199 | POTM_IOBSJ:long_name = "ORCA grid search J coordinate" ; |
---|
200 | int POTM_IOBSK(N_OBS, N_LEVELS) ; |
---|
201 | POTM_IOBSK:long_name = "ORCA grid search K coordinate" ; |
---|
202 | char POTM_GRID(STRINGGRID) ; |
---|
203 | POTM_GRID:long_name = "ORCA grid search grid (T,U,V)" ; |
---|
204 | float PSAL_OBS(N_OBS, N_LEVELS) ; |
---|
205 | PSAL_OBS:long_name = "Practical salinity" ; |
---|
206 | PSAL_OBS:units = "PSU" ; |
---|
207 | PSAL_OBS:_Fillvalue = 99999.f ; |
---|
208 | float PSAL_Hx(N_OBS, N_LEVELS) ; |
---|
209 | PSAL_Hx:long_name = "Model interpolated practical salinity" ; |
---|
210 | PSAL_Hx:units = "PSU" ; |
---|
211 | PSAL_Hx:_Fillvalue = 99999.f ; |
---|
212 | int PSAL_QC(N_OBS) ; |
---|
213 | PSAL_QC:long_name = "Quality on practical salinity" ; |
---|
214 | PSAL_QC:Conventions = "q where q =[0,9]" ; |
---|
215 | PSAL_QC:_Fillvalue = 0 ; |
---|
216 | int PSAL_QC_FLAGS(N_OBS, N_QCF) ; |
---|
217 | PSAL_QC_FLAGS:long_name = "Quality flags on practical salinity" ; |
---|
218 | PSAL_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
219 | PSAL_QC_FLAGS:_Fillvalue = 0 ; |
---|
220 | int PSAL_LEVEL_QC(N_OBS, N_LEVELS) ; |
---|
221 | PSAL_LEVEL_QC:long_name = "Quality for each level on practical salinity" ; |
---|
222 | PSAL_LEVEL_QC:Conventions = "q where q =[0,9]" ; |
---|
223 | PSAL_LEVEL_QC:_Fillvalue = 0 ; |
---|
224 | int PSAL_LEVEL_QC_FLAGS(N_OBS, N_LEVELS, N_QCF) ; |
---|
225 | PSAL_LEVEL_QC_FLAGS:long_name = "Quality flags for each level on practical salinity" ; |
---|
226 | PSAL_LEVEL_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
227 | PSAL_LEVEL_QC_FLAGS:_Fillvalue = 0 ; |
---|
228 | int PSAL_IOBSI(N_OBS) ; |
---|
229 | PSAL_IOBSI:long_name = "ORCA grid search I coordinate" ; |
---|
230 | int PSAL_IOBSJ(N_OBS) ; |
---|
231 | PSAL_IOBSJ:long_name = "ORCA grid search J coordinate" ; |
---|
232 | int PSAL_IOBSK(N_OBS, N_LEVELS) ; |
---|
233 | PSAL_IOBSK:long_name = "ORCA grid search K coordinate" ; |
---|
234 | char PSAL_GRID(STRINGGRID) ; |
---|
235 | PSAL_GRID:long_name = "ORCA grid search grid (T,U,V)" ; |
---|
236 | float TEMP(N_OBS, N_LEVELS) ; |
---|
237 | TEMP:long_name = "Insitu temperature" ; |
---|
238 | TEMP:units = "Degrees Celsius" ; |
---|
239 | TEMP:_Fillvalue = 99999.f ; |
---|
240 | |
---|
241 | // global attributes: |
---|
242 | :title = "NEMO observation operator output" ; |
---|
243 | :Convention = "NEMO unified observation operator output" ; |
---|
244 | } |
---|
245 | \end{verbatim} |
---|
246 | \end{alltt} |
---|
247 | |
---|
248 | \subsection{Sea level anomaly feedback type observation file header} |
---|
249 | |
---|
250 | \begin{alltt} |
---|
251 | \tiny |
---|
252 | \begin{verbatim} |
---|
253 | netcdf sla_01 { |
---|
254 | dimensions: |
---|
255 | N_OBS = 41301 ; |
---|
256 | N_LEVELS = 1 ; |
---|
257 | N_VARS = 1 ; |
---|
258 | N_QCF = 2 ; |
---|
259 | N_ENTRIES = 1 ; |
---|
260 | N_EXTRA = 1 ; |
---|
261 | STRINGNAM = 8 ; |
---|
262 | STRINGGRID = 1 ; |
---|
263 | STRINGWMO = 8 ; |
---|
264 | STRINGTYP = 4 ; |
---|
265 | STRINGJULD = 14 ; |
---|
266 | variables: |
---|
267 | char VARIABLES(N_VARS, STRINGNAM) ; |
---|
268 | VARIABLES:long_name = "List of variables in feedback files" ; |
---|
269 | char ENTRIES(N_ENTRIES, STRINGNAM) ; |
---|
270 | ENTRIES:long_name = "List of additional entries for each variable in feedback files" ; |
---|
271 | char EXTRA(N_EXTRA, STRINGNAM) ; |
---|
272 | EXTRA:long_name = "List of extra variables" ; |
---|
273 | char STATION_IDENTIFIER(N_OBS, STRINGWMO) ; |
---|
274 | STATION_IDENTIFIER:long_name = "Station identifier" ; |
---|
275 | char STATION_TYPE(N_OBS, STRINGTYP) ; |
---|
276 | STATION_TYPE:long_name = "Code instrument type" ; |
---|
277 | double LONGITUDE(N_OBS) ; |
---|
278 | LONGITUDE:long_name = "Longitude" ; |
---|
279 | LONGITUDE:units = "degrees_east" ; |
---|
280 | LONGITUDE:_Fillvalue = 99999.f ; |
---|
281 | double LATITUDE(N_OBS) ; |
---|
282 | LATITUDE:long_name = "Latitude" ; |
---|
283 | LATITUDE:units = "degrees_north" ; |
---|
284 | LATITUDE:_Fillvalue = 99999.f ; |
---|
285 | double DEPTH(N_OBS, N_LEVELS) ; |
---|
286 | DEPTH:long_name = "Depth" ; |
---|
287 | DEPTH:units = "metre" ; |
---|
288 | DEPTH:_Fillvalue = 99999.f ; |
---|
289 | int DEPTH_QC(N_OBS, N_LEVELS) ; |
---|
290 | DEPTH_QC:long_name = "Quality on depth" ; |
---|
291 | DEPTH_QC:Conventions = "q where q =[0,9]" ; |
---|
292 | DEPTH_QC:_Fillvalue = 0 ; |
---|
293 | int DEPTH_QC_FLAGS(N_OBS, N_LEVELS, N_QCF) ; |
---|
294 | DEPTH_QC_FLAGS:long_name = "Quality flags on depth" ; |
---|
295 | DEPTH_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
296 | double JULD(N_OBS) ; |
---|
297 | JULD:long_name = "Julian day" ; |
---|
298 | JULD:units = "days since JULD_REFERENCE" ; |
---|
299 | JULD:Conventions = "relative julian days with decimal part (as parts of day)" ; |
---|
300 | JULD:_Fillvalue = 99999.f ; |
---|
301 | char JULD_REFERENCE(STRINGJULD) ; |
---|
302 | JULD_REFERENCE:long_name = "Date of reference for julian days" ; |
---|
303 | JULD_REFERENCE:Conventions = "YYYYMMDDHHMMSS" ; |
---|
304 | int OBSERVATION_QC(N_OBS) ; |
---|
305 | OBSERVATION_QC:long_name = "Quality on observation" ; |
---|
306 | OBSERVATION_QC:Conventions = "q where q =[0,9]" ; |
---|
307 | OBSERVATION_QC:_Fillvalue = 0 ; |
---|
308 | int OBSERVATION_QC_FLAGS(N_OBS, N_QCF) ; |
---|
309 | OBSERVATION_QC_FLAGS:long_name = "Quality flags on observation" ; |
---|
310 | OBSERVATION_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
311 | OBSERVATION_QC_FLAGS:_Fillvalue = 0 ; |
---|
312 | int POSITION_QC(N_OBS) ; |
---|
313 | POSITION_QC:long_name = "Quality on position (latitude and longitude)" ; |
---|
314 | POSITION_QC:Conventions = "q where q =[0,9]" ; |
---|
315 | POSITION_QC:_Fillvalue = 0 ; |
---|
316 | int POSITION_QC_FLAGS(N_OBS, N_QCF) ; |
---|
317 | POSITION_QC_FLAGS:long_name = "Quality flags on position" ; |
---|
318 | POSITION_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
319 | POSITION_QC_FLAGS:_Fillvalue = 0 ; |
---|
320 | int JULD_QC(N_OBS) ; |
---|
321 | JULD_QC:long_name = "Quality on date and time" ; |
---|
322 | JULD_QC:Conventions = "q where q =[0,9]" ; |
---|
323 | JULD_QC:_Fillvalue = 0 ; |
---|
324 | int JULD_QC_FLAGS(N_OBS, N_QCF) ; |
---|
325 | JULD_QC_FLAGS:long_name = "Quality flags on date and time" ; |
---|
326 | JULD_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
327 | JULD_QC_FLAGS:_Fillvalue = 0 ; |
---|
328 | int ORIGINAL_FILE_INDEX(N_OBS) ; |
---|
329 | ORIGINAL_FILE_INDEX:long_name = "Index in original data file" ; |
---|
330 | ORIGINAL_FILE_INDEX:_Fillvalue = -99999 ; |
---|
331 | float SLA_OBS(N_OBS, N_LEVELS) ; |
---|
332 | SLA_OBS:long_name = "Sea level anomaly" ; |
---|
333 | SLA_OBS:units = "metre" ; |
---|
334 | SLA_OBS:_Fillvalue = 99999.f ; |
---|
335 | float SLA_Hx(N_OBS, N_LEVELS) ; |
---|
336 | SLA_Hx:long_name = "Model interpolated sea level anomaly" ; |
---|
337 | SLA_Hx:units = "metre" ; |
---|
338 | SLA_Hx:_Fillvalue = 99999.f ; |
---|
339 | int SLA_QC(N_OBS) ; |
---|
340 | SLA_QC:long_name = "Quality on sea level anomaly" ; |
---|
341 | SLA_QC:Conventions = "q where q =[0,9]" ; |
---|
342 | SLA_QC:_Fillvalue = 0 ; |
---|
343 | int SLA_QC_FLAGS(N_OBS, N_QCF) ; |
---|
344 | SLA_QC_FLAGS:long_name = "Quality flags on sea level anomaly" ; |
---|
345 | SLA_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
346 | SLA_QC_FLAGS:_Fillvalue = 0 ; |
---|
347 | int SLA_LEVEL_QC(N_OBS, N_LEVELS) ; |
---|
348 | SLA_LEVEL_QC:long_name = "Quality for each level on sea level anomaly" ; |
---|
349 | SLA_LEVEL_QC:Conventions = "q where q =[0,9]" ; |
---|
350 | SLA_LEVEL_QC:_Fillvalue = 0 ; |
---|
351 | int SLA_LEVEL_QC_FLAGS(N_OBS, N_LEVELS, N_QCF) ; |
---|
352 | SLA_LEVEL_QC_FLAGS:long_name = "Quality flags for each level on sea level anomaly" ; |
---|
353 | SLA_LEVEL_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
354 | SLA_LEVEL_QC_FLAGS:_Fillvalue = 0 ; |
---|
355 | int SLA_IOBSI(N_OBS) ; |
---|
356 | SLA_IOBSI:long_name = "ORCA grid search I coordinate" ; |
---|
357 | int SLA_IOBSJ(N_OBS) ; |
---|
358 | SLA_IOBSJ:long_name = "ORCA grid search J coordinate" ; |
---|
359 | int SLA_IOBSK(N_OBS, N_LEVELS) ; |
---|
360 | SLA_IOBSK:long_name = "ORCA grid search K coordinate" ; |
---|
361 | char SLA_GRID(STRINGGRID) ; |
---|
362 | SLA_GRID:long_name = "ORCA grid search grid (T,U,V)" ; |
---|
363 | float MDT(N_OBS, N_LEVELS) ; |
---|
364 | MDT:long_name = "Mean Dynamic Topography" ; |
---|
365 | MDT:units = "metre" ; |
---|
366 | MDT:_Fillvalue = 99999.f ; |
---|
367 | |
---|
368 | // global attributes: |
---|
369 | :title = "NEMO observation operator output" ; |
---|
370 | :Convention = "NEMO unified observation operator output" ; |
---|
371 | } |
---|
372 | \end{verbatim} |
---|
373 | \end{alltt} |
---|
374 | |
---|
375 | \subsection{Sea surface temperature feedback type observation file header} |
---|
376 | |
---|
377 | \begin{alltt} |
---|
378 | \tiny |
---|
379 | \begin{verbatim} |
---|
380 | netcdf sst_01 { |
---|
381 | dimensions: |
---|
382 | N_OBS = 33099 ; |
---|
383 | N_LEVELS = 1 ; |
---|
384 | N_VARS = 1 ; |
---|
385 | N_QCF = 2 ; |
---|
386 | N_ENTRIES = 1 ; |
---|
387 | STRINGNAM = 8 ; |
---|
388 | STRINGGRID = 1 ; |
---|
389 | STRINGWMO = 8 ; |
---|
390 | STRINGTYP = 4 ; |
---|
391 | STRINGJULD = 14 ; |
---|
392 | variables: |
---|
393 | char VARIABLES(N_VARS, STRINGNAM) ; |
---|
394 | VARIABLES:long_name = "List of variables in feedback files" ; |
---|
395 | char ENTRIES(N_ENTRIES, STRINGNAM) ; |
---|
396 | ENTRIES:long_name = "List of additional entries for each variable in feedback files" ; |
---|
397 | char STATION_IDENTIFIER(N_OBS, STRINGWMO) ; |
---|
398 | STATION_IDENTIFIER:long_name = "Station identifier" ; |
---|
399 | char STATION_TYPE(N_OBS, STRINGTYP) ; |
---|
400 | STATION_TYPE:long_name = "Code instrument type" ; |
---|
401 | double LONGITUDE(N_OBS) ; |
---|
402 | LONGITUDE:long_name = "Longitude" ; |
---|
403 | LONGITUDE:units = "degrees_east" ; |
---|
404 | LONGITUDE:_Fillvalue = 99999.f ; |
---|
405 | double LATITUDE(N_OBS) ; |
---|
406 | LATITUDE:long_name = "Latitude" ; |
---|
407 | LATITUDE:units = "degrees_north" ; |
---|
408 | LATITUDE:_Fillvalue = 99999.f ; |
---|
409 | double DEPTH(N_OBS, N_LEVELS) ; |
---|
410 | DEPTH:long_name = "Depth" ; |
---|
411 | DEPTH:units = "metre" ; |
---|
412 | DEPTH:_Fillvalue = 99999.f ; |
---|
413 | int DEPTH_QC(N_OBS, N_LEVELS) ; |
---|
414 | DEPTH_QC:long_name = "Quality on depth" ; |
---|
415 | DEPTH_QC:Conventions = "q where q =[0,9]" ; |
---|
416 | DEPTH_QC:_Fillvalue = 0 ; |
---|
417 | int DEPTH_QC_FLAGS(N_OBS, N_LEVELS, N_QCF) ; |
---|
418 | DEPTH_QC_FLAGS:long_name = "Quality flags on depth" ; |
---|
419 | DEPTH_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
420 | double JULD(N_OBS) ; |
---|
421 | JULD:long_name = "Julian day" ; |
---|
422 | JULD:units = "days since JULD_REFERENCE" ; |
---|
423 | JULD:Conventions = "relative julian days with decimal part (as parts of day)" ; |
---|
424 | JULD:_Fillvalue = 99999.f ; |
---|
425 | char JULD_REFERENCE(STRINGJULD) ; |
---|
426 | JULD_REFERENCE:long_name = "Date of reference for julian days" ; |
---|
427 | JULD_REFERENCE:Conventions = "YYYYMMDDHHMMSS" ; |
---|
428 | int OBSERVATION_QC(N_OBS) ; |
---|
429 | OBSERVATION_QC:long_name = "Quality on observation" ; |
---|
430 | OBSERVATION_QC:Conventions = "q where q =[0,9]" ; |
---|
431 | OBSERVATION_QC:_Fillvalue = 0 ; |
---|
432 | int OBSERVATION_QC_FLAGS(N_OBS, N_QCF) ; |
---|
433 | OBSERVATION_QC_FLAGS:long_name = "Quality flags on observation" ; |
---|
434 | OBSERVATION_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
435 | OBSERVATION_QC_FLAGS:_Fillvalue = 0 ; |
---|
436 | int POSITION_QC(N_OBS) ; |
---|
437 | POSITION_QC:long_name = "Quality on position (latitude and longitude)" ; |
---|
438 | POSITION_QC:Conventions = "q where q =[0,9]" ; |
---|
439 | POSITION_QC:_Fillvalue = 0 ; |
---|
440 | int POSITION_QC_FLAGS(N_OBS, N_QCF) ; |
---|
441 | POSITION_QC_FLAGS:long_name = "Quality flags on position" ; |
---|
442 | POSITION_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
443 | POSITION_QC_FLAGS:_Fillvalue = 0 ; |
---|
444 | int JULD_QC(N_OBS) ; |
---|
445 | JULD_QC:long_name = "Quality on date and time" ; |
---|
446 | JULD_QC:Conventions = "q where q =[0,9]" ; |
---|
447 | JULD_QC:_Fillvalue = 0 ; |
---|
448 | int JULD_QC_FLAGS(N_OBS, N_QCF) ; |
---|
449 | JULD_QC_FLAGS:long_name = "Quality flags on date and time" ; |
---|
450 | JULD_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
451 | JULD_QC_FLAGS:_Fillvalue = 0 ; |
---|
452 | int ORIGINAL_FILE_INDEX(N_OBS) ; |
---|
453 | ORIGINAL_FILE_INDEX:long_name = "Index in original data file" ; |
---|
454 | ORIGINAL_FILE_INDEX:_Fillvalue = -99999 ; |
---|
455 | float SST_OBS(N_OBS, N_LEVELS) ; |
---|
456 | SST_OBS:long_name = "Sea surface temperature" ; |
---|
457 | SST_OBS:units = "Degree centigrade" ; |
---|
458 | SST_OBS:_Fillvalue = 99999.f ; |
---|
459 | float SST_Hx(N_OBS, N_LEVELS) ; |
---|
460 | SST_Hx:long_name = "Model interpolated sea surface temperature" ; |
---|
461 | SST_Hx:units = "Degree centigrade" ; |
---|
462 | SST_Hx:_Fillvalue = 99999.f ; |
---|
463 | int SST_QC(N_OBS) ; |
---|
464 | SST_QC:long_name = "Quality on sea surface temperature" ; |
---|
465 | SST_QC:Conventions = "q where q =[0,9]" ; |
---|
466 | SST_QC:_Fillvalue = 0 ; |
---|
467 | int SST_QC_FLAGS(N_OBS, N_QCF) ; |
---|
468 | SST_QC_FLAGS:long_name = "Quality flags on sea surface temperature" ; |
---|
469 | SST_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
470 | SST_QC_FLAGS:_Fillvalue = 0 ; |
---|
471 | int SST_LEVEL_QC(N_OBS, N_LEVELS) ; |
---|
472 | SST_LEVEL_QC:long_name = "Quality for each level on sea surface temperature" ; |
---|
473 | SST_LEVEL_QC:Conventions = "q where q =[0,9]" ; |
---|
474 | SST_LEVEL_QC:_Fillvalue = 0 ; |
---|
475 | int SST_LEVEL_QC_FLAGS(N_OBS, N_LEVELS, N_QCF) ; |
---|
476 | SST_LEVEL_QC_FLAGS:long_name = "Quality flags for each level on sea surface temperature" ; |
---|
477 | SST_LEVEL_QC_FLAGS:Conventions = "NEMOVAR flag conventions" ; |
---|
478 | SST_LEVEL_QC_FLAGS:_Fillvalue = 0 ; |
---|
479 | int SST_IOBSI(N_OBS) ; |
---|
480 | SST_IOBSI:long_name = "ORCA grid search I coordinate" ; |
---|
481 | int SST_IOBSJ(N_OBS) ; |
---|
482 | SST_IOBSJ:long_name = "ORCA grid search J coordinate" ; |
---|
483 | int SST_IOBSK(N_OBS, N_LEVELS) ; |
---|
484 | SST_IOBSK:long_name = "ORCA grid search K coordinate" ; |
---|
485 | char SST_GRID(STRINGGRID) ; |
---|
486 | SST_GRID:long_name = "ORCA grid search grid (T,U,V)" ; |
---|
487 | |
---|
488 | // global attributes: |
---|
489 | :title = "NEMO observation operator output" ; |
---|
490 | :Convention = "NEMO unified observation operator output" ; |
---|
491 | } |
---|
492 | \end{verbatim} |
---|
493 | \end{alltt} |
---|
494 | |
---|
495 | \section{Theoretical details} |
---|
496 | |
---|
497 | \subsection{Horizontal interpolation methods} |
---|
498 | |
---|
499 | Consider an observation point ${\rm P}$ with |
---|
500 | with longitude and latitude $({\lambda_{}}_{\rm P}, \phi_{\rm P})$ and the |
---|
501 | four nearest neighbouring model grid points ${\rm A}$, ${\rm B}$, ${\rm C}$ |
---|
502 | and ${\rm D}$ with longitude and latitude ($\lambda_{\rm A}$, $\phi_{\rm A}$), |
---|
503 | ($\lambda_{\rm B}$, $\phi_{\rm B}$) etc. |
---|
504 | All horizontal interpolation methods implemented in NEMO |
---|
505 | estimate the value of a model variable $x$ at point $P$ as |
---|
506 | a weighted linear combination of the values of the model |
---|
507 | variables at the grid points ${\rm A}$, ${\rm B}$ etc.: |
---|
508 | \begin{eqnarray} |
---|
509 | {x_{}}_{\rm P} & \hspace{-2mm} = \hspace{-2mm} & |
---|
510 | \frac{1}{w} \left( {w_{}}_{\rm A} {x_{}}_{\rm A} + |
---|
511 | {w_{}}_{\rm B} {x_{}}_{\rm B} + |
---|
512 | {w_{}}_{\rm C} {x_{}}_{\rm C} + |
---|
513 | {w_{}}_{\rm D} {x_{}}_{\rm D} \right) |
---|
514 | \end{eqnarray} |
---|
515 | where ${w_{}}_{\rm A}$, ${w_{}}_{\rm B}$ etc. are the respective weights for the |
---|
516 | model field at points ${\rm A}$, ${\rm B}$ etc., and |
---|
517 | $w = {w_{}}_{\rm A} + {w_{}}_{\rm B} + {w_{}}_{\rm C} + {w_{}}_{\rm D}$. |
---|
518 | |
---|
519 | Four different possibilities are available for computing the weights. |
---|
520 | |
---|
521 | \begin{enumerate} |
---|
522 | |
---|
523 | \item[1.] {\bf Great-Circle distance-weighted interpolation.} The weights |
---|
524 | are computed as a function of the great-circle distance $s(P, \cdot)$ |
---|
525 | between $P$ and the model grid points $A$, $B$ etc. For example, |
---|
526 | the weight given to the field ${x_{}}_{\rm A}$ is specified as the |
---|
527 | product of the distances from ${\rm P}$ to the other points: |
---|
528 | \begin{eqnarray} |
---|
529 | {w_{}}_{\rm A} = s({\rm P}, {\rm B}) \, s({\rm P}, {\rm C}) \, s({\rm P}, {\rm D}) |
---|
530 | \nonumber |
---|
531 | \end{eqnarray} |
---|
532 | where |
---|
533 | \begin{eqnarray} |
---|
534 | s\left ({\rm P}, {\rm M} \right ) |
---|
535 | & \hspace{-2mm} = \hspace{-2mm} & |
---|
536 | \cos^{-1} \! \left\{ |
---|
537 | \sin {\phi_{}}_{\rm P} \sin {\phi_{}}_{\rm M} |
---|
538 | + \cos {\phi_{}}_{\rm P} \cos {\phi_{}}_{\rm M} |
---|
539 | \cos ({\lambda_{}}_{\rm M} - {\lambda_{}}_{\rm P}) |
---|
540 | \right\} |
---|
541 | \end{eqnarray} |
---|
542 | and $M$ corresponds to $B$, $C$ or $D$. |
---|
543 | A more stable form of the great-circle distance formula for |
---|
544 | small distances ($x$ near 1) involves the arcsine function |
---|
545 | ($e.g.$ see p.~101 of \citet{Daley_Barker_Bk01}: |
---|
546 | \begin{eqnarray} |
---|
547 | s\left( {\rm P}, {\rm M} \right) |
---|
548 | & \hspace{-2mm} = \hspace{-2mm} & |
---|
549 | \sin^{-1} \! \left\{ \sqrt{ 1 - x^2 } \right\} |
---|
550 | \nonumber |
---|
551 | \end{eqnarray} |
---|
552 | where |
---|
553 | \begin{eqnarray} |
---|
554 | x & \hspace{-2mm} = \hspace{-2mm} & |
---|
555 | {a_{}}_{\rm M} {a_{}}_{\rm P} + {b_{}}_{\rm M} {b_{}}_{\rm P} + {c_{}}_{\rm M} {c_{}}_{\rm P} |
---|
556 | \nonumber |
---|
557 | \end{eqnarray} |
---|
558 | and |
---|
559 | \begin{eqnarray} |
---|
560 | {a_{}}_{\rm M} & \hspace{-2mm} = \hspace{-2mm} & \sin {\phi_{}}_{\rm M}, |
---|
561 | \nonumber \\ |
---|
562 | {a_{}}_{\rm P} & \hspace{-2mm} = \hspace{-2mm} & \sin {\phi_{}}_{\rm P}, |
---|
563 | \nonumber \\ |
---|
564 | {b_{}}_{\rm M} & \hspace{-2mm} = \hspace{-2mm} & \cos {\phi_{}}_{\rm M} \cos {\phi_{}}_{\rm M}, |
---|
565 | \nonumber \\ |
---|
566 | {b_{}}_{\rm P} & \hspace{-2mm} = \hspace{-2mm} & \cos {\phi_{}}_{\rm P} \cos {\phi_{}}_{\rm P}, |
---|
567 | \nonumber \\ |
---|
568 | {c_{}}_{\rm M} & \hspace{-2mm} = \hspace{-2mm} & \cos {\phi_{}}_{\rm M} \sin {\phi_{}}_{\rm M}, |
---|
569 | \nonumber \\ |
---|
570 | {c_{}}_{\rm P} & \hspace{-2mm} = \hspace{-2mm} & \cos {\phi_{}}_{\rm P} \sin {\phi_{}}_{\rm P}. |
---|
571 | \nonumber |
---|
572 | \nonumber |
---|
573 | \end{eqnarray} |
---|
574 | |
---|
575 | \item[2.] {\bf Great-Circle distance-weighted interpolation with small angle |
---|
576 | approximation.} Similar to the previous interpolation but with the |
---|
577 | distance $s$ computed as |
---|
578 | \begin{eqnarray} |
---|
579 | s\left( {\rm P}, {\rm M} \right) |
---|
580 | & \hspace{-2mm} = \hspace{-2mm} & |
---|
581 | \sqrt{ \left( {\phi_{}}_{\rm M} - {\phi_{}}_{\rm P} \right)^{2} |
---|
582 | + \left( {\lambda_{}}_{\rm M} - {\lambda_{}}_{\rm P} \right)^{2} |
---|
583 | \cos^{2} {\phi_{}}_{\rm M} } |
---|
584 | \end{eqnarray} |
---|
585 | where $M$ corresponds to $A$, $B$, $C$ or $D$. |
---|
586 | |
---|
587 | \item[3.] {\bf Bilinear interpolation for a regular spaced grid.} The |
---|
588 | interpolation is split into two 1D interpolations in the longitude |
---|
589 | and latitude directions, respectively. |
---|
590 | |
---|
591 | \item[4.] {\bf Bilinear remapping interpolation for a general grid.} An |
---|
592 | iterative scheme that involves first mapping a quadrilateral cell |
---|
593 | into a cell with coordinates (0,0), (1,0), (0,1) and (1,1). This |
---|
594 | method is based on the SCRIP interpolation package \citep{Jones_MWR99,Jones_Bk01}. |
---|
595 | |
---|
596 | \end{enumerate} |
---|
597 | |
---|
598 | \subsection{Grid search} |
---|
599 | |
---|
600 | For many grids used by the NEMO model, such as the ORCA family, |
---|
601 | the horizontal grid coordinates $i$ and $j$ are not simple functions |
---|
602 | of latitude and longitude. Therefore, it is not always straightforward |
---|
603 | to determine the grid points surrounding any given observational position. |
---|
604 | Before the interpolation can be performed, a search |
---|
605 | algorithm is then required to determine the corner points of |
---|
606 | the quadrilateral cell in which the observation is located. |
---|
607 | This is the most difficult and time consuming part of the |
---|
608 | 2D interpolation procedure. |
---|
609 | A robust test for determining if an observation falls |
---|
610 | within a given quadrilateral cell is as follows. Let |
---|
611 | ${\rm P}({\lambda_{}}_{\rm P} ,{\phi_{}}_{\rm P} )$ denote the observation point, |
---|
612 | and let ${\rm A}({\lambda_{}}_{\rm A} ,{\phi_{}}_{\rm A} )$, |
---|
613 | ${\rm B}({\lambda_{}}_{\rm B} ,{\phi_{}}_{\rm B} )$, |
---|
614 | ${\rm C}({\lambda_{}}_{\rm C} ,{\phi_{}}_{\rm C} )$ |
---|
615 | and |
---|
616 | ${\rm D}({\lambda_{}}_{\rm D} ,{\phi_{}}_{\rm D} )$ denote |
---|
617 | the bottom left, bottom right, top left and top right |
---|
618 | corner points of the cell, respectively. |
---|
619 | To determine if P is inside |
---|
620 | the cell, we verify that the cross-products |
---|
621 | \begin{eqnarray} |
---|
622 | \begin{array}{lllll} |
---|
623 | {{\bf r}_{}}_{\rm PA} \times {{\bf r}_{}}_{\rm PC} |
---|
624 | & = & [({\lambda_{}}_{\rm A}\; -\; {\lambda_{}}_{\rm P} ) |
---|
625 | ({\phi_{}}_{\rm C} \; -\; {\phi_{}}_{\rm P} ) |
---|
626 | - ({\lambda_{}}_{\rm C}\; -\; {\lambda_{}}_{\rm P} ) |
---|
627 | ({\phi_{}}_{\rm A} \; -\; {\phi_{}}_{\rm P} )] \; \widehat{\bf k} \\ |
---|
628 | {{\bf r}_{}}_{\rm PB} \times {{\bf r}_{}}_{\rm PA} |
---|
629 | & = & [({\lambda_{}}_{\rm B}\; -\; {\lambda_{}}_{\rm P} ) |
---|
630 | ({\phi_{}}_{\rm A} \; -\; {\phi_{}}_{\rm P} ) |
---|
631 | - ({\lambda_{}}_{\rm A}\; -\; {\lambda_{}}_{\rm P} ) |
---|
632 | ({\phi_{}}_{\rm B} \; -\; {\phi_{}}_{\rm P} )] \; \widehat{\bf k} \\ |
---|
633 | {{\bf r}_{}}_{\rm PC} \times {{\bf r}_{}}_{\rm PD} |
---|
634 | & = & [({\lambda_{}}_{\rm C}\; -\; {\lambda_{}}_{\rm P} ) |
---|
635 | ({\phi_{}}_{\rm D} \; -\; {\phi_{}}_{\rm P} ) |
---|
636 | - ({\lambda_{}}_{\rm D}\; -\; {\lambda_{}}_{\rm P} ) |
---|
637 | ({\phi_{}}_{\rm C} \; -\; {\phi_{}}_{\rm P} )] \; \widehat{\bf k} \\ |
---|
638 | {{\bf r}_{}}_{\rm PD} \times {{\bf r}_{}}_{\rm PB} |
---|
639 | & = & [({\lambda_{}}_{\rm D}\; -\; {\lambda_{}}_{\rm P} ) |
---|
640 | ({\phi_{}}_{\rm B} \; -\; {\phi_{}}_{\rm P} ) |
---|
641 | - ({\lambda_{}}_{\rm B}\; -\; {\lambda_{}}_{\rm P} ) |
---|
642 | ({\phi_{}}_{\rm D} \; - \; {\phi_{}}_{\rm P} )] \; \widehat{\bf k} \\ |
---|
643 | \end{array} |
---|
644 | \label{eq:cross} |
---|
645 | \end{eqnarray} |
---|
646 | point in the opposite direction to the unit normal |
---|
647 | $\widehat{\bf k}$ (i.e., that the coefficients of |
---|
648 | $\widehat{\bf k}$ are negative), |
---|
649 | where ${{\bf r}_{}}_{\rm PA}$, ${{\bf r}_{}}_{\rm PB}$, |
---|
650 | etc. correspond to the vectors between points P and A, |
---|
651 | P and B, etc.. The method used is |
---|
652 | similar to the method used in |
---|
653 | the Scripp interpolation package \citep{Jones_Bk01}. |
---|
654 | |
---|
655 | In order to speed up the grid search, there is the possibility to construct |
---|
656 | a lookup table for a user specified resolution. This lookup |
---|
657 | table contains the lower and upper bounds on the $i$ and $j$ indices |
---|
658 | to be searched for on a regular grid. For each observation position, |
---|
659 | the closest point on the regular grid of this position is computed and |
---|
660 | the $i$ and $j$ ranges of this point searched to determine the precise |
---|
661 | four points surrounding the observation. |
---|
662 | |
---|
663 | \subsection{Parallel aspects of horizontal interpolation} |
---|
664 | |
---|
665 | For horizontal interpolation, there is the basic problem that the |
---|
666 | observations are unevenly distributed on the globe. In numerical |
---|
667 | models, it is common to divide the model grid into subgrids (or |
---|
668 | domains) where each subgrid is executed on a single processing element |
---|
669 | with explicit message passing for exchange of information along the |
---|
670 | domain boundaries when running on a massively parallel processor (MPP) |
---|
671 | system. This approach is used by the \NEMO. |
---|
672 | |
---|
673 | For observations there is no natural distribution since the |
---|
674 | observations are not equally distributed on the globe. |
---|
675 | Two options have been made available: 1) geographical distribution; |
---|
676 | and 2) round-robin. |
---|
677 | |
---|
678 | \subsubsection{Geographical distribution of observations among processors} |
---|
679 | |
---|
680 | \begin{figure} |
---|
681 | \begin{center} |
---|
682 | \includegraphics[width=10cm,height=12cm,angle=-90.]{./TexFiles/Figures/Fig_ASM_obsdist_local} |
---|
683 | \end{center} |
---|
684 | \caption{Example of the distribution of observations with the geographical |
---|
685 | distribution of observational data.} |
---|
686 | \label{fig:obslocal} |
---|
687 | \end{figure} |
---|
688 | |
---|
689 | This is the simplest option in which the observations are distributed according |
---|
690 | to the domain of the grid-point parallelization. Figure~\ref{fig:obslocal} |
---|
691 | shows an example of the distribution of the {\em in situ} data on processors |
---|
692 | with a different colour for each observation |
---|
693 | on a given processor for a 4 $\times$ 2 decomposition with ORCA2. |
---|
694 | The grid-point domain decomposition is clearly visible on the plot. |
---|
695 | |
---|
696 | The advantage of this approach is that all |
---|
697 | information needed for horizontal interpolation is available without |
---|
698 | any MPP communication. Of course, this is under the assumption that |
---|
699 | we are only using a $2 \times 2$ grid-point stencil for the interpolation |
---|
700 | (e.g., bilinear interpolation). For higher order interpolation schemes this |
---|
701 | is no longer valid. A disadvantage with the above scheme is that the number of |
---|
702 | observations on each processor can be very different. If the cost of |
---|
703 | the actual interpolation is expensive relative to the communication of |
---|
704 | data needed for interpolation, this could lead to load imbalance. |
---|
705 | |
---|
706 | \subsubsection{Round-robin distribution of observations among processors} |
---|
707 | |
---|
708 | \begin{figure} |
---|
709 | \begin{center} |
---|
710 | \includegraphics[width=10cm,height=12cm,angle=-90.]{./TexFiles/Figures/Fig_ASM_obsdist_global} |
---|
711 | \end{center} |
---|
712 | \caption{Example of the distribution of observations with the round-robin |
---|
713 | distribution of observational data.} |
---|
714 | \label{fig:obsglobal} |
---|
715 | \end{figure} |
---|
716 | |
---|
717 | An alternative approach is to distribute the observations equally |
---|
718 | among processors and use message passing in order to retrieve |
---|
719 | the stencil for interpolation. The simplest distribution of the observations |
---|
720 | is to distribute them using a round-robin scheme. Figure~\ref{fig:obsglobal} |
---|
721 | shows the distribution of the {\em in situ} data on processors for the |
---|
722 | round-robin distribution of observations with a different colour for |
---|
723 | each observation on a given processor for a 4 $\times$ 2 decomposition |
---|
724 | with ORCA2 for the same input data as in Fig.~\ref{fig:obslocal}. |
---|
725 | The observations are now clearly randomly distributed on the globe. |
---|
726 | In order to be able to perform horizontal interpolation in this case, |
---|
727 | a subroutine has been developed that retrieves any grid points in the |
---|
728 | global space. |
---|
729 | |
---|
730 | \subsection{Vertical interpolation operator} |
---|
731 | |
---|
732 | The vertical interpolation is achieved using either a cubic spline or |
---|
733 | linear interpolation. For the cubic spline, the top and |
---|
734 | bottom boundary conditions for the second derivative of the |
---|
735 | interpolating polynomial in the spline are set to zero. |
---|
736 | At the bottom boundary, this is done using the land-ocean mask. |
---|