New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Chap_OBS.tex in branches/dev_1784_doc/DOC/TexFiles/Chapters – NEMO

source: branches/dev_1784_doc/DOC/TexFiles/Chapters/Chap_OBS.tex @ 2284

Last change on this file since 2284 was 2284, checked in by djlea, 14 years ago

Add draft documentation for OBS and ASM code

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