FV3 Bundle
oda_core.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
20 !
21 !<CONTACT EMAIL="matthew.harrison@noaa.gov"> Matthew Harrison
22 !</CONTACT>
23 !
24 !<OVERVIEW>
25 ! core ocean data assimilation subroutines for ocean models. This module
26 ! includes interfaces to :
27 !
28 ! (i) initialize the ODA core with observations and model
29 ! (ii) request ocean state increments using observed or idealized data
30 ! (iii) calculate model guess differences and output to file
31 ! (iv) terminate the DA core.
32 !
33 ! NOTE: Profile files conform to v0.1a profile metadata standard.
34 ! This metadata standard will evolve in time and will maintain
35 ! backward compability.
36 !
37 ! NOTE: This code is still under development. ODA data structures should be
38 ! opaque in future releases.
39 !
40 !</OVERVIEW>
41 !
42 !<NAMELIST NAME="oda_core_nml">
43 ! <DATA NAME="max_misfit" TYPE="real">
44 ! flag measurement if abs(omf) exceeds max_misfit
45 !</DATA>
46 ! <DATA NAME="min_obs_err_t" TYPE="real">
47 ! minumum rms temperature observation error (degC) for
48 ! variable obs error, else nominal temp error
49 !</DATA>
50 ! <DATA NAME="min_obs_err_s" TYPE="real">
51 ! minimum rms salinity observation error (g/kg) for
52 ! variable obs error,else nominal salt error
53 !</DATA>
54 ! <DATA NAME="eta_tide_const" TYPE="real">
55 ! Tidal internal displacement amplitude (meters) for observational error estimation.
56 !</DATA>
57 ! <DATA NAME="min_prof_depth" TYPE="real">
58 ! Minimum profile depth (meters)
59 !</DATA>
60 ! <DATA NAME="max_prof_spacing" TYPE="real">
61 ! Data must be contiguous to this resolution (meters). Otherwise flag is activated.
62 ! </DATA>
63 ! <DATA NAME="data_window" TYPE="integer">
64 ! Half-width of profile time window (days)
65 ! </DATA>
66 ! <DATA NAME="add_tidal_aliasing" TYPE="logical">
67 ! Add tidal aliasing to observation error. Use eta_tide_const * dT/dz to estimate
68 ! observation.
69 ! </DATA>
70 ! <DATA NAME="max_profiles" TYPE="integer">
71 ! Allocation size of profile array
72 ! </DATA>
73 ! <DATA NAME="max_sfc_obs" TYPE="integer">
74 ! Allocation size of surface data array
75 ! </DATA>
76 ! <DATA NAME="temp_obs_rmse" TYPE="real">
77 ! nominal temperature error rms error
78 ! </DATA>
79 ! <DATA NAME="salt_obs_rmse" TYPE="real">
80 ! nominal salinity error rms error
81 ! </DATA>
82 !</NAMELIST>
83  use fms_mod, only : file_exist,read_data
84  use mpp_mod, only : mpp_error, fatal, note, mpp_sum, stdout,&
85  mpp_sync_self, mpp_pe,mpp_npes,mpp_root_pe,&
89  use time_manager_mod, only : time_type, operator( <= ), operator( - ), &
90  operator( > ), operator ( < ), set_time, set_date, &
92  use get_cal_time_mod, only : get_cal_time
93  use axis_utils_mod, only : frac_index
94  use constants_mod, only : radian, pi
95  use oda_types_mod
97 
98  implicit none
99 
100  private
101 
102  real, private :: max_misfit = 5.0 ! reject fg diff gt max_misfit
103 
104 
105  real, parameter, private :: depth_min=0.0, depth_max=10000.
106  real, parameter, private :: temp_min=-3.0, temp_max=40.
107  real, parameter, private :: salt_min=0.0, salt_max=45.
108 
109  integer :: max_profiles = 250000, max_sfc_obs = 1
110  integer, parameter, private :: max_files=100
111  integer, parameter, private :: profile_file = 1,sfc_file= 2,&
113 
114 
115 ! parameters for obs error inflation due to internal tidal displacements
116 
117  real, private :: min_obs_err_t = 0.5, min_obs_err_s=0.1, eta_tide_const = 7.0
118 
119  type(ocean_profile_type), target, save, private, allocatable :: profiles(:)
120  type(ocean_surface_type), target, save, private, allocatable :: sfc_obs(:)
121 
122  integer, private, save :: num_profiles, num_sfc_obs ! total number of observations
123 
124  integer, private, save :: isc, iec, jsc, jec, isd, ied, jsd, jed ! indices for local and global domain
125  integer, private, save :: isg, ieg, jsg, jeg
126  integer, private, save :: nk
127 
128 ! parameters used to restrict profiles
129  real, private, save :: min_prof_depth = 200.0 ! profile data ending
130  ! above this level are
131  ! discarded.
132  real, private, save :: max_prof_spacing = 1.e5 ! reject profile data
133  ! if not contiguous
134  ! at this resolution
135 
136 ! internal arrays for forward and backward observations
137  real, dimension(:,:,:), allocatable, private, save :: sum_wgt, nobs
138 
139  type(time_type) , dimension(0:100), public :: time_window
140 
141 ! the DA module maintains a unique grid,domain association
142 
143  type(grid_type), pointer :: grd
144 
145 ! DA grid (1-d) coordinates. Generalized horizontal coordinates are not supporte.
146 
147  real, allocatable, dimension(:) :: x_grid, y_grid
148 
149  real :: temp_obs_rmse = 0.7071
150  real :: salt_obs_rmse = 0.1
151 
152  logical :: add_tidal_aliasing=.false.
153 
154  logical :: localize_data = .true.
155 
156  logical :: debug=.false.
157 
158  integer :: ndebug=10
159 
160  integer, allocatable, dimension(:) :: nprof_in_comp_domain
161 
162  namelist /oda_core_nml/ max_misfit,add_tidal_aliasing,min_obs_err_t,&
165 
166 
167 ! NOTE: Surface observation type is not yet implemented.
168 
169 ! transform from Grd => Obs
170 
171  interface forward_obs
172  module procedure forward_obs_profile
173  module procedure forward_obs_sfc
174  end interface
175 
176 ! transform from Obs => Grd
177 
178  interface backward_obs
179  module procedure backward_obs_profile
180  module procedure backward_obs_sfc
181  end interface
182 
183 
184 ! one-time association between profile and grid, used for storing
185 ! interpolation weights
186 
188  module procedure assign_forward_model_profile
189  module procedure assign_forward_model_sfc
190  end interface
191 
192 ! duplicate observation array
193 
194  interface copy_obs
195  module procedure copy_obs_prof
196  module procedure copy_obs_sfc
197  end interface
198 
199 ! multiply observation data by inverse error variance
200 
201  interface mult_obs_i_mse
202  module procedure mult_obs_i_mse_profile
203  module procedure mult_obs_i_mse_sfc
204  end interface
205 
206 ! difference between two observations (e.g. model first guess and observations)
207 
208  interface diff_obs
209  module procedure diff_obs_profile
210  module procedure diff_obs_sfc
211  end interface
212 
213 ! inflate observational error based on first guess misfit
214 ! and time window.
215 
217  module procedure adjust_obs_error_profile
218  module procedure adjust_obs_error_sfc
219  end interface
220 
221 
222  interface nullify_obs
223  module procedure nullify_obs_prof
224  end interface
225 
226  public :: forward_obs, backward_obs, &
231 
232  contains
233 
234 
235  subroutine open_profile_dataset(filename, localize)
236 ! <DESCRIPTION>
237 ! open dataset containing profile information in fms station format.
238 ! store internally.
239 ! </DESCRIPTION>
240 
241  use mpp_io_mod, only : mpp_open, mpp_get_atts, mpp_get_info, &
242  mpp_get_fields, mpp_read, mpp_single, mpp_multi, mpp_netcdf,&
243  axistype, atttype, fieldtype, mpp_rdonly, mpp_get_axes, mpp_close,&
244  mpp_get_att_char
245 
246  character(len=*), intent(in) :: filename
247  logical, intent(in), optional :: localize
248  logical :: found_neighbor,continue_looking
249  integer, parameter :: max_levels=1000
250  integer :: unit, ndim, nvar, natt, nstation
251  character(len=32) :: fldname, axisname
252  type(atttype), allocatable, dimension(:), target :: global_atts
253  type(atttype), pointer :: version => null()
254  type(axistype), pointer :: depth_axis => null(), station_axis => null()
255  type(axistype), allocatable, dimension(:), target :: axes
256  type(fieldtype), allocatable, dimension(:), target :: fields
257 
258  type(fieldtype), pointer :: field_lon => null(), field_lat => null(), field_probe => null(),&
259  field_time => null(), field_depth => null()
260  type(fieldtype), pointer :: field_temp => null(), field_salt => null(), &
261  field_error => null(), field_link => null()
262  ! NOTE: fields are restricted to be in separate files
263  real :: lon, lat, time, depth(max_levels), temp(max_levels), salt(max_levels),&
264  error(max_levels), rprobe, profile_error, rlink
265  integer :: nlev, probe, yr, mon, day, hr, minu, sec, kl, outunit
266  integer :: num_levs, num_levs_temp, num_levs_salt,&
267  k, kk, ll, i, i0, j0, k0, nlevs, a, nn, ii, nlinks
268  real :: ri0, rj0, rk0, dx1, dx2, dy1, dy2
269  character(len=128) :: time_units, attname, catt
270  type(time_type) :: profile_time
271  integer :: flag_t(max_levels), flag_s(max_levels), cpe
272  logical :: data_is_local, &
273  continue
274  logical :: found_temp=.false., found_salt=.false.
275  real, dimension(max_links,max_levels) :: temp_bfr, salt_bfr, depth_bfr
276  integer, dimension(max_links,max_levels) :: flag_t_bfr, flag_s_bfr
277  real :: temp_missing=missing_value,salt_missing=missing_value,&
278  depth_missing=missing_value
279  real :: max_prof_depth, zdist
280 
281 
282  ! read timeseries of local observational data from NetCDF files
283  ! and allocate ocean_obs arrays.
284 
285  ! File structure:
286  ! dimensions:
287  ! depth_index;station=UNLIMITED;
288  ! variables:
289  ! depth_index(depth_index);
290  ! station(station);
291  ! longitude(station);
292  ! latitude(station);
293  ! time(station);
294  ! data(station,depth_index);
295  ! depth(station,depth_index);
296  ! probe(station);
297  ! err(station, depth_index);
298 
299  cpe = mpp_pe()
300 
301  dx1 = (x_grid(isc)-x_grid(isc-1))/2.0
302  dx2 = (x_grid(iec+1)-x_grid(iec))/2.0
303  dy1 = (y_grid(jsc)-y_grid(jsc-1))/2.0
304  dy2 = (y_grid(jec+1)-y_grid(jec))/2.0
305 
306  localize_data = .true.
307 
308  if (PRESENT(localize)) localize_data = localize
309 
310  call mpp_open(unit,filename,form=mpp_netcdf,fileset=mpp_single,&
311  threading=mpp_multi,action=mpp_rdonly)
312  call mpp_get_info(unit, ndim, nvar, natt, nstation)
313 
314  outunit = stdout()
315  write(outunit,*) 'Opened profile dataset :',trim(filename)
316 
317  ! get version number of profiles
318 
319  allocate(global_atts(natt))
320  call mpp_get_atts(unit,global_atts)
321 
322  do i=1,natt
323  catt = mpp_get_att_char(global_atts(i))
324  select case (lowercase(trim(catt)))
325  case ('version')
326  version => global_atts(i)
327  end select
328  end do
329 
330  if (.NOT.ASSOCIATED(version)) then
331  call mpp_error(note,'no version number available for profile file, assuming v0.1a ')
332  else
333  write(outunit,*) 'Reading profile dataset version = ',trim(catt)
334  endif
335 
336 
337  ! get axis information
338 
339  allocate (axes(ndim))
340  call mpp_get_axes(unit,axes)
341  do i=1,ndim
342  call mpp_get_atts(axes(i),name=axisname)
343  select case (lowercase(trim(axisname)))
344  case ('depth_index')
345  depth_axis => axes(i)
346  case ('station_index')
347  station_axis => axes(i)
348  end select
349  end do
350 
351  if (.NOT.ASSOCIATED(depth_axis) .or. .NOT.ASSOCIATED(station_axis)) then
352  call mpp_error(fatal,'depth and/or station axes do not exist in input file')
353  endif
354 
355 
356 ! get selected field information.
357 ! NOTE: not checking for all variables here.
358 
359  allocate(fields(nvar))
360  call mpp_get_fields(unit,fields)
361  do i=1,nvar
362  call mpp_get_atts(fields(i),name=fldname)
363  select case (lowercase(trim(fldname)))
364  case ('longitude')
365  field_lon => fields(i)
366  case ('latitude')
367  field_lat => fields(i)
368  case ('probe')
369  field_probe => fields(i)
370  case ('time')
371  field_time => fields(i)
372  case ('temp')
373  field_temp => fields(i)
374  case ('salt')
375  field_salt => fields(i)
376  case ('depth')
377  field_depth => fields(i)
378  case ('link')
379  field_link => fields(i)
380  case ('error')
381  field_error => fields(i)
382  end select
383  enddo
384 
385  call mpp_get_atts(depth_axis,len=nlevs)
386 
387  if (nlevs > max_levels) call mpp_error(fatal,'increase parameter max_levels ')
388 
389  if (nlevs < 1) call mpp_error(fatal)
390 
391  outunit = stdout()
392  write(outunit,*) 'There are ', nstation, ' records in this dataset'
393  write(outunit,*) 'Searching for profiles matching selection criteria ...'
394 
395  if (ASSOCIATED(field_temp)) found_temp=.true.
396  if (ASSOCIATED(field_salt)) found_salt=.true.
397 
398  if (.not. found_temp .and. .not. found_salt) then
399  write(outunit,*) 'temp or salt not found in profile file'
400  call mpp_error(fatal)
401  endif
402 
403  call mpp_get_atts(field_time,units=time_units)
404  if (found_temp) call mpp_get_atts(field_temp,missing=temp_missing)
405  if (found_salt) call mpp_get_atts(field_salt,missing=salt_missing)
406 
407  call mpp_get_atts(field_depth,missing=depth_missing)
408 
409  if (found_salt) then
410  write(outunit,*) 'temperature and salinity where available'
411  else
412  write(outunit,*) 'temperature only records'
413  endif
414 
415 
416  i=1
417  continue=.true.
418 
419  do while (continue)
420 
421  depth(:) = missing_value
422  temp(:) = missing_value
423  salt(:) = missing_value
424 ! required fields
425  call mpp_read(unit,field_lon,lon,tindex=i)
426  call mpp_read(unit,field_lat,lat, tindex=i)
427  call mpp_read(unit,field_time,time,tindex=i)
428  call mpp_read(unit,field_depth,depth(1:nlevs),tindex=i)
429  if (found_temp) call mpp_read(unit,field_temp,temp(1:nlevs),tindex=i)
430  if (found_salt) call mpp_read(unit,field_salt,salt(1:nlevs),tindex=i)
431 ! not required fields
432  if (ASSOCIATED(field_error)) then
433  call mpp_read(unit,field_error,profile_error,tindex=i)
434  endif
435  if (ASSOCIATED(field_probe)) then
436  call mpp_read(unit, field_probe, rprobe,tindex=i)
437  endif
438  if (ASSOCIATED(field_link)) then
439  call mpp_read(unit,field_link,rlink,tindex=i)
440  else
441  rlink = 0.0
442  endif
443  probe=rprobe
444  data_is_local = .false.
445 ! NOTE: assuming grid is modulo 360 here. This needs to be generalized.
446 
447  if (lon .lt. x_grid(isg-1) ) lon = lon + 360.0
448  if (lon .gt. x_grid(ieg+1) ) lon = lon - 360.0
449 
450 ! localized data is within region bounded by halo points
451 ! (halo size = 1) adjacent to boundary points of computational domain
452 
453  if (lon >= x_grid(isc-1) .and. &
454  lon < x_grid(iec+1) .and. &
455  lat >= y_grid(jsc-1) .and. &
456  lat < y_grid(jec+1)) data_is_local = .true.
457 
458 
459  profile_time = get_cal_time(time,time_units,'julian')
460 
461 
462  if ( data_is_local .OR. .NOT.localize_data) then
463 
465  if (num_profiles > max_profiles) then
466  call mpp_error(fatal,'maximum number of profiles exceeded.&
467  &Resize parameter max_profiles in ocean_obs_mod')
468 
469  endif
470 
472 
473  num_levs_temp = 0
474  num_levs_salt = 0
475  do k = 1, nlevs
476 
477 ! flag=0 denotes a valid profile level, anything else
478 ! is invalid. See NODC codes.
479 !================================================================
480 !0 - accepted station
481 !1 - failed annual standard deviation check
482 !2 - two or more density inversions (Levitus, 1982 criteria)
483 !3 - flagged cruise
484 !4 - failed seasonal standard deviation check
485 !5 - failed monthly standard deviation check
486 !6 - flag 1 and flag 4
487 !7 - bullseye from standard level data or failed annual and monthly
488 ! standard deviation check
489 !8 - failed seasonal and monthly standard deviation check
490 !9 - failed annual, seasonal, and monthly standard deviation check
491 !================================================================
492 
493  flag_t(k) = 0;flag_s(k) = 0
494 
495  if (.not.found_salt) then
496  flag_s(k) = 1
497  endif
498 
499  if (depth(k) .eq. depth_missing .or. depth(k) .lt.depth_min&
500  .or. depth(k) .gt. depth_max) then
501  depth(k) = missing_value
502  flag_t(k)=1
503  flag_s(k)=1
504  endif
505 
506  if (found_temp .and. flag_t(k) .eq. 0) then
507  if (temp(k) .eq. temp_missing .or. temp(k) .lt. temp_min&
508  .or. temp(k) .gt. temp_max) then
509  temp(k) = missing_value
510  flag_t(k) = 1
511  flag_s(k) = 1 ! flag salt if no temperature data
512  else
513  num_levs_temp = num_levs_temp+1
514  endif
515  endif
516  if (found_salt .and. flag_s(k) .eq. 0) then
517  if (salt(k) .eq. salt_missing .or. salt(k) .lt. salt_min&
518  .or. salt(k) .gt. salt_max) then
519  salt(k) = missing_value
520  flag_s(k) = 1
521  else
522  num_levs_salt = num_levs_salt+1
523  endif
524  endif
525 
526  enddo
527 
528 ! large profile are stored externally in separate records
529 ! follow the links to get complete profile
530 
531  ii=i+1
532  nlinks = 0
533  do while (rlink > 0.0 .and. nlinks .lt. max_links)
534  nlinks=nlinks+1
535  depth_bfr(nlinks,:) = missing_value
536  temp_bfr(nlinks,:) = missing_value
537  salt_bfr(nlinks,:) = missing_value
538  call mpp_read(unit,field_depth,depth_bfr(nlinks,1:nlevs),tindex=ii)
539  if (found_temp) call mpp_read(unit,field_temp,temp_bfr(nlinks,1:nlevs),tindex=ii)
540  if (found_salt) call mpp_read(unit,field_salt,salt_bfr(nlinks,1:nlevs),tindex=ii)
541  call mpp_read(unit,field_link,rlink,tindex=ii)
542  ii=ii+1
543  enddo
544  i=ii ! set record counter to start of next profile
545 
546  if (nlinks > 0) then
547  do nn = 1, nlinks
548  do k=1, nlevs
549  flag_t_bfr(nn,k) = 0
550  flag_s_bfr(nn,k) = 0
551  if (depth_bfr(nn,k) .eq. depth_missing .or.&
552  depth_bfr(nn,k) .lt. depth_min .or. &
553  depth_bfr(nn,k) .gt. depth_max) then
554  depth_bfr(nn,k) = missing_value
555  flag_t_bfr(nn,k) = 1
556  flag_s_bfr(nn,k) = 1
557  endif
558  if (found_temp .and. flag_t_bfr(nn,k) .eq. 0) then
559  if (temp_bfr(nn,k) .eq. temp_missing .or.&
560  temp_bfr(nn,k) .lt. temp_min .or.&
561  temp_bfr(nn,k) .gt. temp_max) then
562  temp_bfr(nn,k) = missing_value
563  flag_t_bfr(nn,k) = 1
564  flag_s_bfr(nn,k) = 1
565  else
566  num_levs_temp = num_levs_temp+1
567  endif
568  endif
569  if (found_salt .and. flag_s_bfr(nn,k) .eq. 0) then
570  if (salt_bfr(nn,k) .eq. salt_missing .or.&
571  salt_bfr(nn,k) .lt. salt_min .or.&
572  salt_bfr(nn,k) .gt. salt_max) then
573  salt_bfr(nn,k) = missing_value
574  flag_t_bfr(nn,k) = 0
575  flag_s_bfr(nn,k) = 1
576  else
577  num_levs_salt = num_levs_salt+1
578  endif
579  endif
580  enddo
581  enddo
582  endif
583 
584  num_levs = max(num_levs_temp,num_levs_salt)
585 
586  if (num_levs == 0) then
587  if (i .gt. nstation) continue = .false.
588  cycle
589  endif
590 
591  allocate(profiles(num_profiles)%depth(num_levs))
593  if (num_levs_temp .gt. 0) then
594  allocate(profiles(num_profiles)%data_t(num_levs))
596  allocate(profiles(num_profiles)%flag_t(num_levs))
597  profiles(num_profiles)%flag_t= 1
598  profiles(num_profiles)%nvar=1
599  endif
600 
601  if (num_levs_salt .gt. 0) then
602  allocate(profiles(num_profiles)%data_s(num_levs))
604  allocate(profiles(num_profiles)%flag_s(num_levs))
605  profiles(num_profiles)%flag_s= 1
607  endif
608 
609 
610  if (probe < 1 ) probe = 0
611  profiles(num_profiles)%probe = probe
612  profiles(num_profiles)%levels = num_levs
613  profiles(num_profiles)%lat = lat
614  profiles(num_profiles)%lon = lon
615  allocate(profiles(num_profiles)%ms_t(num_levs))
616  profiles(num_profiles)%ms_t(:) = temp_obs_rmse**2.0 ! default error variance for temperature
617 
618  if(num_levs_salt .gt. 0) then
619  allocate(profiles(num_profiles)%ms_s(num_levs))
620  profiles(num_profiles)%ms_s(:) = salt_obs_rmse**2.0 ! default error variance for salinity
621  endif
622 
623  kk= 1
624  do k = 1, nlevs
625  if (flag_t(k) .eq. 0) then
626  if (kk > profiles(num_profiles)%levels) then
627  call mpp_error(fatal)
628  endif
629  profiles(num_profiles)%depth(kk) = depth(k)
630  profiles(num_profiles)%data_t(kk) = temp(k)
631  profiles(num_profiles)%flag_t(kk) = 0
632  if (found_salt .and. flag_s(k) .eq. 0) then
633  profiles(num_profiles)%data_s(kk) = salt(k)
634  profiles(num_profiles)%flag_s(kk) = 0
635  endif
636  kk=kk+1
637  endif
638  enddo
639 
640  do nn = 1, nlinks
641  do k = 1, nlevs
642  if (flag_t_bfr(nn,k) .eq. 0) then
643  if (kk > profiles(num_profiles)%levels) then
644  call mpp_error(fatal)
645  endif
646  profiles(num_profiles)%depth(kk) = depth_bfr(nn,k)
647  profiles(num_profiles)%data_t(kk) = temp_bfr(nn,k)
648  profiles(num_profiles)%flag_t(kk) = 0
649  if (found_salt .and. flag_s_bfr(nn,k) .eq. 0) then
650  profiles(num_profiles)%data_s(kk) = salt_bfr(nn,k)
651  profiles(num_profiles)%flag_s(kk) = 0
652  endif
653  kk=kk+1
654  endif
655  enddo
656  enddo
657 
658  profiles(num_profiles)%time = profile_time
659 
660 ! calculate interpolation coefficients (make sure to account for grid offsets here!)
661 ! NOTE: this only works for lat/lon grids. Lower left indices.
662 !
663 
664  ri0 = frac_index(lon, x_grid(isg-1:ieg+1)) - 1.
665  rj0 = frac_index(lat, y_grid(jsg-1:jeg+1)) - 1.
666  i0 = floor(ri0)
667  j0 = floor(rj0)
668  profiles(num_profiles)%i_index = ri0
669  profiles(num_profiles)%j_index = rj0
670  profiles(num_profiles)%accepted = .true.
671  if (i0 < 0 .or. j0 < 0) then
672  profiles(num_profiles)%accepted = .false.
673  endif
674  if (i0 > ieg .or. j0 > jeg) then
675  call mpp_error(fatal,'grid lat/lon index is out of bounds ')
676  endif
677  if ((i0 < isc-1 .or. i0 > iec) .and. localize_data) then
678  call mpp_error(fatal,'grid lat/lon index is out of bounds ')
679  endif
680  if ((j0 < jsc-1 .or. j0 > jec) .and. localize_data) then
681  call mpp_error(fatal,'grid lat/lon index is out of bounds ')
682  endif
683 !
684 ! flag the profile if it sits on a model land point
685 !
686  if (profiles(num_profiles)%accepted ) then
687  if (grd%mask(i0,j0,1) == 0.0 .or. &
688  grd%mask(i0+1,j0,1) == 0.0 .or. &
689  grd%mask(i0,j0+1,1) == 0.0 .or. &
690  grd%mask(i0+1,j0+1,1) == 0.0) then
691  profiles(num_profiles)%accepted = .false.
692  endif
693  endif
694 
695  if (profiles(num_profiles)%accepted) then
696  allocate(profiles(num_profiles)%k_index(profiles(num_profiles)%levels))
697  max_prof_depth=0.0
698  do k=1, profiles(num_profiles)%levels
699  k0=0
700  if (profiles(num_profiles)%flag_t(k).eq.0) then
701  rk0 = frac_index(profiles(num_profiles)%depth(k), grd%z(:))
702  k0 = floor(rk0)
703  if ( k0 == -1) then
704  if (profiles(num_profiles)%depth(k) .lt. grd%z(1)) then
705  k0 = 1
706  rk0 = 1.0
707  else if (profiles(num_profiles)%depth(k) .gt. grd%z(grd%nk)) then
708  profiles(num_profiles)%flag_t(k) = 1
709  endif
710  endif
711  else
712  cycle
713  endif
714 
715  if (k0 .gt. size(grd%z)-1 ) then
716  write(*,*) 'k0 out of bounds, rk0,k0= ',rk0,k0
717  write(*,*) 'Z_bound= ',grd%z_bound
718  write(*,*) 'Profile%depth= ',profiles(num_profiles)%depth
719 
720  call mpp_error(fatal)
721  endif
722 
723  profiles(num_profiles)%k_index(k) = rk0
724 
725 ! flag depth level if adjacent model grid points are land
726 
727  if (profiles(num_profiles)%flag_t(k) .eq. 0) then
728  if (i0 .lt. 0 .or. j0 .lt. 0 .or. k0 .lt. 0) then
729  write(*,*) 'profile index out of bounds'
730  write(*,*) 'i0,j0,k0=',i0,j0,k0
731  write(*,*) 'lon,lat,depth=',profiles(num_profiles)%lon,&
733  call mpp_error(fatal)
734  endif
735 
736  if (grd%mask(i0,j0,k0) == 0.0 .or. &
737  grd%mask(i0+1,j0,k0) == 0.0 .or. &
738  grd%mask(i0,j0+1,k0) == 0.0 .or. &
739  grd%mask(i0+1,j0+1,k0) == 0.0) then
740  profiles(num_profiles)%flag_t(k) = 1
741  endif
742  if (grd%mask(i0,j0,k0+1) == 0.0 .or. &
743  grd%mask(i0+1,j0,k0+1) == 0.0 .or. &
744  grd%mask(i0,j0+1,k0+1) == 0.0 .or. &
745  grd%mask(i0+1,j0+1,k0+1) == 0.0) then
746  profiles(num_profiles)%flag_t(k) = 1
747  endif
748  if (profiles(num_profiles)%flag_t(k) .eq. 0) then
749  max_prof_depth = profiles(num_profiles)%depth(k)
750  endif
751  endif
752 
753  enddo ! Prof%levels loop
754 
755 ! Flag profile if it is too shallow.
756 
757  if (max_prof_depth .lt. min_prof_depth) then
758  profiles(num_profiles)%accepted = .false.
759  endif
760 
761  found_neighbor=.false.
762 
763  do k=2,profiles(num_profiles)%levels - 1
764  if (profiles(num_profiles)%flag_t(k) .eq. 0) then
765  kk = k-1
766  found_neighbor = .false.
767  continue_looking = .true.
768  do while (continue_looking .and. kk .ge. 1)
769  if (profiles(num_profiles)%flag_t(kk) .eq. 0) then
770  zdist = profiles(num_profiles)%depth(k) - profiles(num_profiles)%depth(kk)
771  if (zdist .gt. max_prof_spacing) then
772  profiles(num_profiles)%accepted = .false.
773  goto 199
774  else
775  continue_looking = .false.
776  found_neighbor = .true.
777  endif
778  else
779  kk = kk - 1
780  endif
781  end do
782  kk = k+1
783  continue_looking = .true.
784  do while (continue_looking .and. kk .le. profiles(num_profiles)%levels)
785  if (profiles(num_profiles)%flag_t(kk).eq. 0) then
786  zdist = profiles(num_profiles)%depth(kk) - profiles(num_profiles)%depth(k)
787  if (zdist .gt. max_prof_spacing) then
788  profiles(num_profiles)%accepted = .false.
789  goto 199
790  else
791  continue_looking = .false.
792  found_neighbor = .true.
793  endif
794  else
795  kk = kk+1
796  endif
797  enddo
798  endif
799  enddo
800 
801  if (.not. found_neighbor) profiles(num_profiles)%accepted = .false.
802 
803 199 continue
804 
805  endif ! if Prof%accept
806  else ! data is not local
807  i = i+1
808  endif ! if data_is_local
809 
810 
811  if (i .gt. nstation) continue = .false.
812 
813  enddo
814 
815 ! a = nprof_in_comp_domain(cpe)
816 
817 
818 
819 ! call mpp_broadcast(nprof_in_comp_domain(cpe),cpe)
820 
821 ! call mpp_sum(a)
822 
823 ! write(stdout(),*) 'A grand total of ',int(a),' profiles satisify acceptance criteria'
824 
825 ! do i=0,mpp_npes()
826 ! write(stdout(),*) 'pe=',i,'profile count=',nprof_in_comp_domain(i)
827 ! enddo
828 
829  call mpp_sync_self()
830  call mpp_close(unit)
831 
832  end subroutine open_profile_dataset
833 
834  subroutine get_obs(model_time, Prof, Sfc, nprof, nsfc)
836 
837  ! get profiles and sfc
838  ! obs relevant to current analysis interval
839 
840  type(time_type), intent(in) :: model_time
841  type(ocean_profile_type), dimension(:) :: prof
842  type(ocean_surface_type), dimension(:) :: sfc
843  integer, intent(inout) :: nprof, nsfc
844 
845  integer :: i,k,yr,mon,day,hr,minu,sec,a,mon_obs,yr_obs, outunit
846  type(time_type) :: tdiff
847  character(len=1) :: cchar
848 
849  nprof=0
850  nsfc=0
851 
852  outunit = stdout()
853  write(outunit,*) 'Gathering profiles for current analysis time'
854  call get_date(model_time,yr,mon,day,hr,minu,sec)
855  write(outunit,'(a,i4,a,i2,a,i2)') 'Current yyyy/mm/dd= ',yr,'/',mon,'/',day
856  write(outunit,*) 'num_profiles=',num_profiles
857 
858 
859  do i=1,num_profiles
860 
861  if (debug .and. i.le.ndebug) then
862  call get_date(profiles(i)%time,yr,mon,day,hr,minu,sec)
863  write(*,*) 'in get_obs prof time: yy/mm/dd= ',yr,mon,day
864  endif
865 
866  if (profiles(i)%time <= model_time) then
867  tdiff = model_time - profiles(i)%time
868  else
869  tdiff = profiles(i)%time - model_time
870  endif
871 
872  if (debug .and. i .le. ndebug) then
873  write(*,*) 'Prof%accepted=',profiles(i)%accepted
874  endif
875 
876  if (tdiff <= time_window(0) .and. &
877  profiles(i)%accepted) then
878  nprof=nprof+1
879  if (nprof > size(prof,1)) &
880  call mpp_error(fatal,'increase size of Prof array before call to get_obs')
881  call copy_obs(profiles(i:i),prof(nprof:nprof))
882  prof(nprof)%tdiff = tdiff
883  if (debug .and. nprof .le. ndebug) then
884  call get_time(tdiff,sec,day)
885  write(*,'(a,i3,a,2f10.5)') 'Accepted profile #',i,' : lon,lat= ',prof(nprof)%lon,prof(nprof)%lat
886  do k=1,prof(nprof)%levels
887  if (prof(nprof)%nvar .eq. 2) then
888  write(*,'(a,i3,a,2f10.5,2i2,2f8.5)') 'Profile #',i,' : temp,salt,flag_t,flag_s,ms_t,ms_s= ',&
889  prof(nprof)%data_t(k),prof(nprof)%data_s(k),prof(nprof)%flag_t(k),prof(nprof)%flag_s(k),&
890  prof(nprof)%ms_t(k),prof(nprof)%ms_s(k)
891  else
892  write(*,'(a,i3,a,2f10.5)') 'Profile #',i,' : temp,flag_t= ',prof(nprof)%data_t(k),prof(nprof)%flag_t(k)
893  endif
894  enddo
895  endif
896 
897  else
898  if (debug .and. i .le. ndebug) then
899  call get_time(tdiff,sec,day)
900  write(*,'(a,i3,a,2f10.5)') 'Rejected profile #',i,' : lon,lat= ',prof(i)%lon,prof(i)%lat
901  do k=1,prof(i)%levels
902  if (prof(i)%nvar .eq. 2) then
903  write(*,'(a,i3,a,2f10.5,2i2)') 'Profile #',i,' : temp,salt,flag_t,flag_s= ',prof(i)%data_t(k),prof(i)%data_s(k),prof(i)%flag_t(k),prof(i)%flag_s(k)
904  else
905  write(*,'(a,i3,a,2f10.5)') 'Profile #',i,' : temp,flag_t= ',prof(i)%data_t(k),prof(i)%flag_t(k)
906  endif
907  enddo
908  endif
909  endif
910 
911  enddo
912 
913  a=nprof
914  call mpp_sum(a)
915  write(outunit,*) 'A total of ',a,' profiles are being used for the current analysis step'
916 
917  return
918 
919  end subroutine get_obs
920 
921  subroutine oda_core_init(Domain, Grid, localize)
923  use fms_mod, only : open_namelist_file, check_nml_error, close_file
924 
925  type(domain2d), intent(in) :: domain
926  type(grid_type), target, intent(in) :: grid
927  logical, intent(in), optional :: localize
928 
929 
930  integer :: ioun, ierr, io_status
931 
932 #ifdef INTERNAL_FILE_NML
933  read (input_nml_file, oda_core_nml, iostat=io_status)
934  ierr = check_nml_error(io_status,'oda_core_nml')
935 #else
936  ioun = open_namelist_file()
937  read(ioun,nml=oda_core_nml,iostat = io_status)
938  ierr = check_nml_error(io_status,'oda_core_nml')
939  call close_file(ioun)
940 #endif
941 
942 ! allocate(nprof_in_comp_domain(0:mpp_npes()-1))
943 ! nprof_in_comp_domain = 0
944 
945  grd => grid
946 
947  call mpp_get_compute_domain(domain,isc,iec,jsc,jec)
948  call mpp_get_data_domain(domain,isd,ied,jsd,jed)
949  call mpp_get_global_domain(domain,isg,ieg,jsg,jeg)
950  nk = size(grid%z)
951 
952  allocate(sum_wgt(isd:ied,jsd:jed,nk))
953  allocate(nobs(isd:ied,jsd:jed,nk))
954 
955  if (PRESENT(localize)) localize_data = localize
956 
958 
959  call write_ocean_data_init()
960 
961  end subroutine oda_core_init
962 
963  subroutine purge_obs()
965  num_profiles=0
966  num_sfc_obs=0
967 
968  end subroutine purge_obs
969 
970  subroutine forward_obs_profile(Model_obs, fg_t, fg_s)
972 ! map first guess to observation locations
973 ! note that forward operator only becomes associated after
974 ! this call
975 
976  type(ocean_profile_type), dimension(:), intent(inout) ::Model_obs
977  real, dimension(isd:ied,jsd:jed,nk) , intent(in) :: fg_t ! first guess for temperature
978  real, dimension(isd:ied,jsd:jed,nk) , intent(in), optional :: fg_s ! first guess for salinity
979 
980  integer :: n, i0, j0, k, num_prof, k0
981  real :: a,b,c
982 
983  character(len=128) :: mesg
984 
985  num_prof = size(model_obs)
986 
987  sum_wgt = 0.0
988 
989  do n = 1, num_prof
990  i0 = floor(model_obs(n)%i_index)
991  j0 = floor(model_obs(n)%j_index)
992  a = model_obs(n)%i_index - i0
993  b = model_obs(n)%j_index - j0
994 
995  if (a >= 1.0 .or. a < 0.0 ) call mpp_error(fatal)
996  if (b >= 1.0 .or. b < 0.0 ) call mpp_error(fatal)
997 
998 
999  if (i0 < isc-1 .or. i0 > iec) then
1000  write(mesg,'(a,i4,2x,i4)') 'out of bounds in forward_obs: i0,j0= ',i0,j0
1001  call mpp_error(fatal,trim(mesg))
1002  endif
1003 
1004  if (j0 < jsc-1 .or. j0 > jec) then
1005  write(mesg,*) 'out of bounds in forward_obs: i0,j0= ',i0,j0
1006  call mpp_error(fatal,trim(mesg))
1007  endif
1008 
1009  if (ASSOCIATED(model_obs(n)%data_t) .and. model_obs(n)%accepted) then
1010 
1011  if (ASSOCIATED(model_obs(n)%Forward_model%wgt))&
1012  model_obs(n)%Forward_model%wgt => null()
1013  allocate(model_obs(n)%Forward_model%wgt(2,2,2,model_obs(n)%levels))
1014 
1015  model_obs(n)%Forward_model%wgt = 0.0
1016 
1017  do k = 1, model_obs(n)%levels
1018  if (model_obs(n)%flag_t(k) .eq. 0) then
1019  k0 = floor(model_obs(n)%k_index(k))
1020  if (k0 < 1 .or. k0 > grd%nk-1) then
1021  write(mesg,*) 'out of bounds in forward_obs: k0= ',k0
1022  call mpp_error(fatal,trim(mesg))
1023  endif
1024  c = model_obs(n)%k_index(k) - k0
1025 
1026  if (c >= 1.0 .or. c < 0.0 ) call mpp_error(fatal)
1027 
1028  model_obs(n)%Forward_model%wgt(1,1,1,k) = (1.0-a)*(1.0-b)*(1.0-c)
1029  model_obs(n)%Forward_model%wgt(2,1,1,k) = a*(1.0-b)*(1.0-c)
1030  model_obs(n)%Forward_model%wgt(1,2,1,k) = (1.0-a)*b*(1.0-c)
1031  model_obs(n)%Forward_model%wgt(2,2,1,k) = a*b*(1.0-c)
1032  model_obs(n)%Forward_model%wgt(1,1,2,k) = (1.0-a)*(1.0-b)*c
1033  model_obs(n)%Forward_model%wgt(2,1,2,k) = a*(1.0-b)*c
1034  model_obs(n)%Forward_model%wgt(1,2,2,k) = (1.0-a)*b*c
1035  model_obs(n)%Forward_model%wgt(2,2,2,k) = a*b*c
1036 
1037  sum_wgt(i0,j0,k0) = sum_wgt(i0,j0,k0)+&
1038  model_obs(n)%Forward_model%wgt(1,1,1,k)
1039  sum_wgt(i0+1,j0,k0) = sum_wgt(i0+1,j0,k0)+&
1040  model_obs(n)%Forward_model%wgt(2,1,1,k)
1041  sum_wgt(i0,j0+1,k0) = sum_wgt(i0,j0+1,k0)+&
1042  model_obs(n)%Forward_model%wgt(1,2,1,k)
1043  sum_wgt(i0+1,j0+1,k0) = sum_wgt(i0+1,j0+1,k0)+&
1044  model_obs(n)%Forward_model%wgt(2,2,1,k)
1045  sum_wgt(i0,j0,k0+1) = sum_wgt(i0,j0,k0+1)+&
1046  model_obs(n)%Forward_model%wgt(1,1,2,k)
1047  sum_wgt(i0+1,j0,k0+1) = sum_wgt(i0+1,j0,k0+1)+&
1048  model_obs(n)%Forward_model%wgt(2,1,2,k)
1049  sum_wgt(i0,j0+1,k0+1) = sum_wgt(i0,j0+1,k0+1)+&
1050  model_obs(n)%Forward_model%wgt(1,2,2,k)
1051  sum_wgt(i0+1,j0+1,k0+1) = sum_wgt(i0+1,j0+1,k0+1)+&
1052  model_obs(n)%Forward_model%wgt(2,2,2,k)
1053 
1054 
1055  model_obs(n)%data_t(k) = &
1056  fg_t(i0,j0,k0)*model_obs(n)%Forward_model%wgt(1,1,1,k) &
1057  + fg_t(i0+1,j0,k0)*model_obs(n)%Forward_model%wgt(2,1,1,k) &
1058  + fg_t(i0,j0+1,k0)*model_obs(n)%Forward_model%wgt(1,2,1,k) &
1059  + fg_t(i0+1,j0+1,k0)*model_obs(n)%Forward_model%wgt(2,2,1,k) &
1060  + fg_t(i0,j0,k0+1)*model_obs(n)%Forward_model%wgt(1,1,2,k) &
1061  + fg_t(i0+1,j0,k0+1)*model_obs(n)%Forward_model%wgt(2,1,2,k) &
1062  + fg_t(i0,j0+1,k0+1)*model_obs(n)%Forward_model%wgt(1,2,2,k) &
1063  + fg_t(i0+1,j0+1,k0+1)*model_obs(n)%Forward_model%wgt(2,2,2,k)
1064 
1065  if (ASSOCIATED(model_obs(n)%data_s) .and. PRESENT(fg_s)) then
1066  model_obs(n)%data_s(k) = &
1067  fg_s(i0,j0,k0)*model_obs(n)%Forward_model%wgt(1,1,1,k) &
1068  + fg_s(i0+1,j0,k0)*model_obs(n)%Forward_model%wgt(2,1,1,k) &
1069  + fg_s(i0,j0+1,k0)*model_obs(n)%Forward_model%wgt(1,2,1,k) &
1070  + fg_s(i0+1,j0+1,k0)*model_obs(n)%Forward_model%wgt(2,2,1,k) &
1071  + fg_s(i0,j0,k0+1)*model_obs(n)%Forward_model%wgt(1,1,2,k) &
1072  + fg_s(i0+1,j0,k0+1)*model_obs(n)%Forward_model%wgt(2,1,2,k) &
1073  + fg_s(i0,j0+1,k0+1)*model_obs(n)%Forward_model%wgt(1,2,2,k) &
1074  + fg_s(i0+1,j0+1,k0+1)*model_obs(n)%Forward_model%wgt(2,2,2,k)
1075  endif
1076  else
1077  if (ASSOCIATED(model_obs(n)%data_t)) then
1078  model_obs(n)%data_t(k) = missing_value
1079  endif
1080  if (ASSOCIATED(model_obs(n)%data_s)) then
1081  model_obs(n)%data_s(k) = missing_value
1082  endif
1083  endif
1084  enddo
1085  endif
1086  enddo
1087 
1088  return
1089 
1090  end subroutine forward_obs_profile
1091 
1092 
1093 
1094 
1095 
1096  subroutine forward_obs_sfc(Sfc, Guess, Diff)
1097  type(ocean_surface_type), intent(in) :: Sfc
1098  type(ocean_dist_type), intent(in) :: Guess
1099  type(ocean_surface_type), intent(inout) :: Diff
1100 
1101 
1102  return
1103 
1104  end subroutine forward_obs_sfc
1105 
1106  subroutine backward_obs_profile(Obs,model_t,model_s)
1108 ! map observations back to model locations
1109 !
1110  type(ocean_profile_type), dimension(:), intent(in) :: Obs
1111  real, dimension(isd:ied,jsd:jed,nk), intent(inout) :: model_t
1112  real, dimension(isd:ied,jsd:jed,nk), intent(inout), optional :: model_s
1113 
1114  real :: tmp
1115  integer :: i0,j0,k0,k,n
1116 
1117  model_t = 0.0
1118  if (PRESENT(model_s)) model_s = 0.0
1119 
1120  do n=1,size(obs)
1121  if (ASSOCIATED(obs(n)%data_t) .and. obs(n)%accepted) then
1122  i0 = floor(obs(n)%i_index)
1123  j0 = floor(obs(n)%j_index)
1124 ! profiles are assumed to lie
1125 ! in domain bounded by first halo points
1126 
1127  if (i0 < isd .or. i0 > ied-1) cycle
1128  if (j0 < jsd .or. j0 > jed-1) cycle
1129 
1130  if (.not.ASSOCIATED(obs(n)%Forward_model%wgt)) call mpp_error(fatal,'forward operator not associated with obs')
1131 
1132  do k=1, obs(n)%levels
1133  if (obs(n)%flag_t(k) .eq. 0) then
1134  k0 = floor(obs(n)%k_index(k))
1135  if (k0 < 1 .or. k0 > grd%nk-1) call mpp_error(fatal,'profile k indx out of bnds')
1136 
1137  tmp = obs(n)%data_t(k)
1138  model_t(i0,j0,k0) = model_t(i0,j0,k0) + tmp*obs(n)%Forward_model%wgt(1,1,1,k)
1139  model_t(i0+1,j0,k0) = model_t(i0+1,j0,k0) + tmp*obs(n)%Forward_model%wgt(2,1,1,k)
1140  model_t(i0,j0+1,k0) = model_t(i0,j0+1,k0) + tmp*obs(n)%Forward_model%wgt(1,2,1,k)
1141  model_t(i0+1,j0+1,k0) = model_t(i0+1,j0+1,k0) + tmp*obs(n)%Forward_model%wgt(2,2,1,k)
1142  model_t(i0,j0,k0+1) = model_t(i0,j0,k0+1) + tmp*obs(n)%Forward_model%wgt(1,1,2,k)
1143  model_t(i0+1,j0,k0+1) = model_t(i0+1,j0,k0+1) + tmp*obs(n)%Forward_model%wgt(2,1,2,k)
1144  model_t(i0,j0+1,k0+1) = model_t(i0,j0+1,k0+1) + tmp*obs(n)%Forward_model%wgt(1,2,2,k)
1145  model_t(i0+1,j0+1,k0+1) = model_t(i0+1,j0+1,k0+1) + tmp*obs(n)%Forward_model%wgt(2,2,2,k)
1146 
1147  if (PRESENT(model_s)) then
1148 
1149  tmp = obs(n)%data_s(k)
1150  model_s(i0,j0,k0) = model_s(i0,j0,k0) + tmp*obs(n)%Forward_model%wgt(1,1,1,k)
1151  model_s(i0+1,j0,k0) = model_s(i0+1,j0,k0) + tmp*obs(n)%Forward_model%wgt(2,1,1,k)
1152  model_s(i0,j0+1,k0) = model_s(i0,j0+1,k0) + tmp*obs(n)%Forward_model%wgt(1,2,1,k)
1153  model_s(i0+1,j0+1,k0) = model_s(i0+1,j0+1,k0) + tmp*obs(n)%Forward_model%wgt(2,2,1,k)
1154  model_s(i0,j0,k0+1) = model_s(i0,j0,k0+1) + tmp*obs(n)%Forward_model%wgt(1,1,2,k)
1155  model_s(i0+1,j0,k0+1) = model_s(i0+1,j0,k0+1) + tmp*obs(n)%Forward_model%wgt(2,1,2,k)
1156  model_s(i0,j0+1,k0+1) = model_s(i0,j0+1,k0+1) + tmp*obs(n)%Forward_model%wgt(1,2,2,k)
1157  model_s(i0+1,j0+1,k0+1) = model_s(i0+1,j0+1,k0+1) + tmp*obs(n)%Forward_model%wgt(2,2,2,k)
1158  endif
1159 
1160  end if
1161 
1162 
1163  end do
1164  end if
1165  end do
1166 
1167 
1168  where(sum_wgt > 0.0)
1169  model_t = model_t /sum_wgt
1170  elsewhere
1171  model_t = 0.0
1172  end where
1173 
1174  if (PRESENT(model_s)) then
1175  where(sum_wgt > 0.0)
1176  model_s = model_s /sum_wgt
1177  elsewhere
1178  model_s = 0.0
1179  end where
1180  endif
1181 
1182  end subroutine backward_obs_profile
1183 
1184 
1185  subroutine backward_obs_sfc(Obs,model)
1187  type(ocean_surface_type), dimension(:), intent(in) :: Obs
1188  real, dimension(isd:ied,jsd:jed,nk), intent(inout) :: model
1189 
1190 
1191  end subroutine backward_obs_sfc
1192 
1193  subroutine assign_forward_model_profile(Obs1,Obs2)
1195  type(ocean_profile_type), dimension(:), target, intent(in) :: Obs1
1196  type(ocean_profile_type), dimension(:), intent(inout) :: Obs2
1197 
1198 
1199  integer :: n
1200 
1201  if (size(obs1) .ne. size(obs2)) call mpp_error(fatal)
1202 
1203  do n=1,size(obs1)
1204 
1205  obs2(n)%Forward_model%wgt => obs1(n)%Forward_model%wgt
1206 
1207  enddo
1208 
1209  end subroutine assign_forward_model_profile
1210 
1211 
1212  subroutine assign_forward_model_sfc(Obs1,Obs2)
1214  type(ocean_surface_type), target, intent(in) :: Obs1
1215  type(ocean_surface_type), intent(inout) :: Obs2
1216 
1217 
1218  return
1219  end subroutine assign_forward_model_sfc
1220 
1221  subroutine diff_obs_profile(prof1, prof2, Diff)
1223  type(ocean_profile_type), dimension(:), intent(in) :: prof1
1224  type(ocean_profile_type), dimension(:), intent(in) :: prof2
1225  type(ocean_profile_type), dimension(:), intent(inout) :: Diff
1226 
1227  integer :: n,k
1228 
1229  if (size(prof1) .ne. size(prof2) ) call mpp_error(fatal)
1230 
1231  if (size(prof1) .ne. size(diff) ) call mpp_error(fatal)
1232 
1233 
1234  do n=1,size(prof1)
1235  do k=1,prof1(n)%levels
1236  if (prof1(n)%flag_t(k) .eq. 0) then
1237  diff(n)%data_t(k) = prof2(n)%data_t(k)-prof1(n)%data_t(k)
1238  else
1239  diff(n)%data_t(k) = missing_value
1240  endif
1241  if (abs(diff(n)%data_t(k)) .gt. max_misfit) then
1242  diff(n)%flag_t(k) = 1
1243  endif
1244  if (ASSOCIATED(prof1(n)%data_s)) then
1245  if (prof1(n)%flag_s(k) .eq. 0) then
1246  diff(n)%data_s(k) = prof2(n)%data_s(k)-prof1(n)%data_s(k)
1247  else
1248  diff(n)%data_s(k) = missing_value
1249  endif
1250  if (abs(diff(n)%data_s(k)) .gt. max_misfit) then
1251  diff(n)%flag_s(k) = 1
1252  endif
1253  endif
1254  enddo
1255  enddo
1256 
1257 
1258  end subroutine diff_obs_profile
1259 
1260  subroutine diff_obs_sfc(prof1,prof2,Diff)
1262  type(ocean_surface_type), dimension(:), intent(in) :: prof1, prof2
1263  type(ocean_surface_type), dimension(:), intent(inout) :: Diff
1264 
1265 
1266  end subroutine diff_obs_sfc
1267 
1268  subroutine copy_obs_prof(obs_in, obs_out)
1270  type(ocean_profile_type), dimension(:), intent(in) :: obs_in
1271  type(ocean_profile_type), dimension(:), intent(inout) :: obs_out
1272 
1273 
1274  integer :: n
1275 
1276  if (size(obs_in) .ne. size(obs_out)) call mpp_error(fatal,&
1277  'size mismatch in call to copy_obs_prof')
1278 
1279 
1280  do n=1,size(obs_in)
1281  call nullify_obs(obs_out(n))
1282  obs_out(n)%nvar = obs_in(n)%nvar
1283  obs_out(n)%project = obs_in(n)%project
1284  obs_out(n)%probe = obs_in(n)%probe
1285  obs_out(n)%ref_inst = obs_in(n)%ref_inst
1286  obs_out(n)%wod_cast_num = obs_in(n)%wod_cast_num
1287  obs_out(n)%fix_depth = obs_in(n)%fix_depth
1288  obs_out(n)%ocn_vehicle = obs_in(n)%ocn_vehicle
1289  obs_out(n)%database_id = obs_in(n)%database_id
1290  obs_out(n)%levels = obs_in(n)%levels
1291  obs_out(n)%profile_flag = obs_in(n)%profile_flag
1292  obs_out(n)%profile_flag_s = obs_in(n)%profile_flag_s
1293  obs_out(n)%lon = obs_in(n)%lon
1294  obs_out(n)%lat = obs_in(n)%lat
1295  obs_out(n)%accepted = obs_in(n)%accepted
1296  ALLOCATE(obs_out(n)%depth(obs_in(n)%levels))
1297  obs_out(n)%depth(:) = obs_in(n)%depth(:)
1298  ALLOCATE(obs_out(n)%data_t(obs_in(n)%levels))
1299  obs_out(n)%data_t(:) = obs_in(n)%data_t(:)
1300  ALLOCATE(obs_out(n)%flag_t(obs_in(n)%levels))
1301  obs_out(n)%flag_t(:) = obs_in(n)%flag_t(:)
1302  if (ASSOCIATED(obs_in(n)%data_s)) then
1303  ALLOCATE(obs_out(n)%data_s(obs_in(n)%levels))
1304  obs_out(n)%data_s(:) = obs_in(n)%data_s(:)
1305  ALLOCATE(obs_out(n)%flag_s(obs_in(n)%levels))
1306  obs_out(n)%flag_s(:) = obs_in(n)%flag_s(:)
1307  endif
1308 
1309  obs_out(n)%time = obs_in(n)%time
1310  obs_out(n)%yyyy = obs_in(n)%yyyy
1311  obs_out(n)%mmdd = obs_in(n)%mmdd
1312  obs_out(n)%i_index = obs_in(n)%i_index
1313  obs_out(n)%j_index = obs_in(n)%j_index
1314  ALLOCATE(obs_out(n)%k_index(obs_in(n)%levels))
1315  obs_out(n)%k_index = obs_in(n)%k_index
1316  ALLOCATE(obs_out(n)%ms_t(obs_in(n)%levels))
1317  obs_out(n)%ms_t = obs_in(n)%ms_t
1318  if (ASSOCIATED(obs_in(n)%ms_s)) then
1319  ALLOCATE(obs_out(n)%ms_s(obs_in(n)%levels))
1320  obs_out(n)%ms_s = obs_in(n)%ms_s
1321  endif
1322 
1323 
1324 
1325  obs_out(n)%tdiff = obs_in(n)%tdiff
1326  obs_out(n)%nbr_index = obs_in(n)%nbr_index
1327  obs_out(n)%nbr_dist = obs_in(n)%nbr_dist
1328  if (ASSOCIATED(obs_in(n)%Model_grid)) &
1329  obs_out(n)%Model_grid => obs_in(n)%Model_Grid
1330  enddo
1331 
1332 end subroutine copy_obs_prof
1333 
1334  subroutine copy_obs_sfc(Obs_in, Obs_out)
1335  type(ocean_surface_type), dimension(:), intent(in) :: Obs_in
1336  type(ocean_surface_type), dimension(:), intent(inout) :: Obs_out
1337 
1338 
1339  return
1340 
1341  end subroutine copy_obs_sfc
1342 
1343  subroutine adjust_obs_error_profile(Prof)
1345  use time_manager_mod, only : get_time
1346 
1347  type(ocean_profile_type), dimension(:), intent(inout) :: Prof
1348  integer :: secs, days, n, k, secs_w, days_w, m
1349  real :: tfac, Ims
1350 
1351  do n=1,size(prof)
1352  call get_time(prof(n)%tdiff,secs, days)
1353  m=prof(n)%probe
1354  call get_time(time_window(m),secs_w,days_w)
1355  tfac = (days + secs/86400.) / days_w
1356  tfac = 1. - min(1.,tfac)
1357  if (tfac > 1.0 ) call mpp_error(fatal)
1358  if (tfac < 0.0 ) call mpp_error(fatal)
1359  do k=1,prof(n)%levels
1360  prof(n)%ms_t(k) = 1.0/ max(1.e-1,tfac) * prof(n)%ms_t(k)
1361  if (ASSOCIATED(prof(n)%data_s)) then
1362  prof(n)%ms_s(k) = 1.0/ max(1.e-1,tfac) * prof(n)%ms_s(k)
1363  endif
1364  end do
1365  end do
1366 
1367  end subroutine adjust_obs_error_profile
1368 
1369  subroutine adjust_obs_error_sfc(Diff)
1371  type(ocean_surface_type), intent(inout) :: Diff
1372 
1373  return
1374 
1375  end subroutine adjust_obs_error_sfc
1376 
1377 
1378  subroutine mult_obs_i_mse_profile(Obs)
1380  type(ocean_profile_type), dimension(:), intent(inout) :: Obs
1381 
1382  integer :: n,k
1383  real :: Ims
1384 
1385  do n=1,size(obs)
1386  do k = 1, obs(n)%levels
1387  ims = 1./obs(n)%ms_t(k)
1388  if (obs(n)%flag_t(k) .eq. 0) obs(n)%data_t(k) = ims*obs(n)%data_t(k)
1389  if (ASSOCIATED(obs(n)%data_s)) then
1390  ims = 1/obs(n)%ms_s(k)
1391  if (obs(n)%flag_s(k) .eq. 0) obs(n)%data_s(k) = ims*obs(n)%data_s(k)
1392  endif
1393  end do
1394  end do
1395 
1396  end subroutine mult_obs_i_mse_profile
1397 
1398  subroutine mult_obs_i_mse_sfc(a, Obs)
1400  real, dimension(:), intent(in) :: a
1401  type(ocean_surface_type), intent(inout) :: Obs
1402 
1403  end subroutine mult_obs_i_mse_sfc
1404 
1405 
1406 
1407 
1408  ! </SUBROUTINE>
1409 ! <FUNCTION NAME="lowercase">
1410 ! <DESCRIPTION>
1411 ! Turn a string from uppercase to lowercase, do nothing if the
1412 ! string is already in lowercase.
1413 ! </DESCRIPTION>
1414  function lowercase (cs)
1415  character(len=*), intent(in) :: cs
1416  character(len=len(cs)) :: lowercase
1417  character :: ca(len(cs))
1418 
1419  integer, parameter :: co=iachar('a')-iachar('A') ! case offset
1420 
1421  ca = transfer(cs,"x",len(cs))
1422  where (ca >= "A" .and. ca <= "Z") ca = achar(iachar(ca)+co)
1423  lowercase = transfer(ca,cs)
1424 
1425  end function lowercase
1426 
1427  subroutine init_observations(localize)
1429  use fms_mod, only : open_namelist_file,close_file,check_nml_error
1430  use mpp_io_mod, only : mpp_open, mpp_ascii, mpp_rdonly, mpp_multi, mpp_single
1431  use mpp_domains_mod, only : mpp_global_field
1432 
1433  logical, intent(in), optional :: localize
1434 
1435  integer :: data_window = 15 ! default data half-window is 15 days
1436 
1437  integer :: i,j
1438 
1439 
1440  type obs_entry_type
1441  character(len=128) :: filename
1442  character(len=16) :: file_type
1443  end type obs_entry_type
1444 
1445  namelist /ocean_obs_nml/ data_window, max_prof_spacing, min_prof_depth
1446 
1447  character(len=128) :: input_files(max_files) = ''
1448  integer :: nfiles, filetype(max_files), ioun, io_status, ierr,&
1449  unit, nrecs, n
1450  character(len=256) :: record
1451  type(obs_entry_type) :: tbl_entry
1452 
1453 #ifdef INTERNAL_FILE_NML
1454  read (input_nml_file, ocean_obs_nml, iostat=io_status)
1455  ierr = check_nml_error(io_status,'ocean_obs_nml')
1456 #else
1457  ioun = open_namelist_file()
1458  read(ioun,nml=ocean_obs_nml,iostat = io_status)
1459  ierr = check_nml_error(io_status,'ocean_obs_nml')
1460  call close_file(ioun)
1461 #endif
1462 
1463  time_window(:) = set_time(0,data_window)
1464 
1465  nfiles=0
1466 
1467  if (file_exist('ocean_obs_table') ) then
1468  call mpp_open(unit,'ocean_obs_table',action=mpp_rdonly)
1469  nfiles = 0;nrecs=0
1470  do while (nfiles <= max_files)
1471  read(unit,'(a)',end=99,err=98) record
1472  nrecs=nrecs+1
1473  if (record(1:1) == '#') cycle
1474  read(record,*,err=98,end=98) tbl_entry
1475  nfiles=nfiles+1
1476  input_files(nfiles) = tbl_entry%filename
1477  select case (trim(tbl_entry%file_type))
1478  case ('profiles')
1479  filetype(nfiles) = profile_file
1480  case ('sfc')
1481  filetype(nfiles) = sfc_file
1482  case ('idealized')
1483  filetype(nfiles) = idealized_profiles
1484  case default
1485  call mpp_error(fatal,'error in obs_table entry format')
1486  end select
1487 98 continue
1488  enddo
1489  call mpp_error(fatal,' number of obs files exceeds max_files parameter')
1490 99 continue
1491 
1492  endif
1493  num_profiles=0
1494  num_sfc_obs=0
1495 
1496 ! get local indices for Model grid
1497 ! Since we currently only support regular grids, the
1498 ! input 2-d grid array is converted to 1-d
1499 ! halo points are added
1500 
1501 ! xhalo=isc-isd
1502 ! yhalo=jsc-jsd
1503 ! if (xhalo.ne.ied-iec) call mpp_error(FATAL)
1504 ! if (yhalo.ne.jed-jec) call mpp_error(FATAL)
1505 
1506  allocate(x_grid(isg-1:ieg+1))
1507  allocate(y_grid(jsg-1:jeg+1))
1508 
1509 
1510  x_grid(isg:ieg) = grd%x(isg:ieg,jsg)
1511  y_grid(jsg:jeg) = grd%y(ieg/4,jsg:jeg)
1512 
1513  allocate(profiles(max_profiles))
1514 
1515  if (grd%cyclic) then
1516  x_grid(isg-1) = x_grid(ieg) - 360. ! assume grid is modulo 360 which is reasonable for data assimilation
1517  x_grid(ieg+1) = x_grid(isg) + 360.
1518  else
1519  x_grid(isg-1) = x_grid(isg) - 1.e-10
1520  x_grid(ieg+1) = x_grid(ieg) + 1.e-10
1521  endif
1522 
1523  y_grid(jsg-1) = y_grid(jsg) - 1.e-10
1524  y_grid(jeg+1) = y_grid(jeg) + 1.e-10
1525 
1526 
1527  do n=1, nfiles
1528  select case (filetype(n))
1529  case (profile_file)
1530  call open_profile_dataset(trim(input_files(n)), localize)
1531  case (idealized_profiles)
1532  call create_ideal_profiles(localize)
1533  case default
1534  call mpp_error(fatal,'filetype not currently supported')
1535  end select
1536  enddo
1537 
1538 
1539  return
1540 
1541  end subroutine init_observations
1542 
1543  subroutine add_tidal_error(Prof)
1544 ! NOT IMPLEMENTED YET !!!
1545  type(ocean_profile_type), intent(inout) :: Prof
1546 
1547  integer :: k
1548  real :: dtdz, err, a1, a2
1549 
1550  if (.not.ASSOCIATED(prof%ms_t)) then
1551  allocate(prof%ms_t(prof%levels))
1552  prof%ms_t(:) = min_obs_err_t
1553  endif
1554 
1555  do k=2,prof%levels - 1
1556  if (prof%flag_t(k-1) .eq. 0 .and. prof%flag_t(k+1) .eq. 0) then
1557  dtdz = (prof%data_t(k+1)-prof%data_t(k-1))/(prof%depth(k+1)-prof%depth(k-1))
1558  a1 = abs(dtdz) * eta_tide_const
1559  err = max(a1,min_obs_err_t)
1560  prof%ms_t(k) = err*err
1561  if (ASSOCIATED(prof%data_s)) then
1562  dtdz = (prof%data_s(k+1)-prof%data_s(k-1))/(prof%depth(k+1)-prof%depth(k-1))
1563  a1 = abs(dtdz) * eta_tide_const
1564  err = max(a1,min_obs_err_s)
1565  prof%ms_s(k) = err*err
1566  endif
1567  endif
1568  enddo
1569 
1570  end subroutine add_tidal_error
1571 
1572  subroutine create_ideal_profiles(localize)
1575 
1576  logical, intent(in), optional :: localize
1577  logical :: localize_data = .true.
1578  integer, parameter :: nlevels = 100 ! number of vertical levels for idealized profiles
1579  real, parameter :: width_trans = 250.0 ! with over which to transition from sfc value to bottom value
1580  real, parameter :: bot_depth = 2000.0 ! bottom depth for idealized profiles
1581  real, allocatable, dimension(:) :: lon,lat, sst, sss, bot_temp, bot_salt, depth
1582  real, allocatable, dimension(:,:) :: temp, salt, temp_error, salt_error
1583  integer, allocatable, dimension(:) :: yr, mon, day, hr, mm, ss
1584  integer :: nstation, unit, n, noobs, i, k
1585  real :: ri0,rj0,rk0, mid_depth, dtdf, temp_cent, depthC_I_trans, dsdf, salt_cent
1586  type(time_type) :: profile_time
1587  logical :: data_is_local
1588  integer :: model, parse_ok, cpe
1589  integer :: i0,j0,k0
1590  real :: dz, a, dx1, dx2, dy1, dy2
1591 
1592  real :: temp_missing=missing_value,salt_missing=missing_value,depth_missing=missing_value
1593  character(len=32) :: fld_type, fld_name
1594  type(method_type), allocatable, dimension(:) :: ocean_obs_methods
1595 
1596  if (PRESENT(localize)) localize_data = localize
1597 
1598 
1599  cpe = mpp_pe()
1600 
1601  dx1 = (x_grid(isc)-x_grid(isc-1))/2.0
1602  dx2 = (x_grid(iec+1)-x_grid(iec))/2.0
1603  dy1 = (y_grid(jsc)-y_grid(jsc-1))/2.0
1604  dy2 = (y_grid(jec+1)-y_grid(jec))/2.0
1605 
1606  model = model_ocean
1607  n = find_field_index(model,'ideal_profiles')
1608  call get_field_info(n,fld_type,fld_name,model,noobs)
1609 
1610  allocate(ocean_obs_methods(noobs))
1611  allocate(lon(noobs),lat(noobs), yr(noobs), mon(noobs), day(noobs), &
1612  hr(noobs), mm(noobs), ss(noobs), &
1613  sst(noobs), sss(noobs), bot_temp(noobs), bot_salt(noobs))
1614  allocate(temp(noobs,nlevels), salt(noobs,nlevels), temp_error(noobs,nlevels), salt_error(noobs,nlevels))
1615  allocate(depth(nlevels))
1616 
1617  call get_field_methods(n,ocean_obs_methods)
1618  do i=1,noobs
1619  parse_ok = parse(ocean_obs_methods(i)%method_control,'lon',lon(i))
1620  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1621  if (lon(i) .lt. x_grid(isg) ) lon(i) = lon(i) + 360.0
1622  if (lon(i) .gt. x_grid(ieg) ) lon(i) = lon(i) - 360.0
1623  parse_ok = parse(ocean_obs_methods(i)%method_control,'lat',lat(i))
1624  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1625  parse_ok = parse(ocean_obs_methods(i)%method_control,'yr',yr(i))
1626  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1627  parse_ok = parse(ocean_obs_methods(i)%method_control,'mon',mon(i))
1628  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1629  parse_ok = parse(ocean_obs_methods(i)%method_control,'day',day(i))
1630  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1631  parse_ok = parse(ocean_obs_methods(i)%method_control,'hr',hr(i))
1632  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1633  parse_ok = parse(ocean_obs_methods(i)%method_control,'mm',mm(i))
1634  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1635  parse_ok = parse(ocean_obs_methods(i)%method_control,'ss',ss(i))
1636  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1637  parse_ok = parse(ocean_obs_methods(i)%method_control,'sst',sst(i))
1638  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1639  parse_ok = parse(ocean_obs_methods(i)%method_control,'sss',sss(i))
1640  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1641  parse_ok = parse(ocean_obs_methods(i)%method_control,'bot_temp',bot_temp(i))
1642  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1643  parse_ok = parse(ocean_obs_methods(i)%method_control,'bot_salt',bot_salt(i))
1644  if (parse_ok == 0) call mpp_error(fatal,'==>Error oda_core_mod: idealized_ocean_profiles table entry error')
1645  enddo
1646 
1647  if (noobs == 0 ) then
1648  call mpp_error(fatal,'==> NOTE from oda_core_mod: no idealized profiles given in field table')
1649  return
1650  endif
1651 
1652  dz = bot_depth/(nlevels-1)
1653 
1654  do k=1,nlevels
1655  depth(k) = (k-1)*dz
1656  enddo
1657 
1658  mid_depth = bot_depth/2.0
1659  depthc_i_trans = mid_depth / width_trans
1660  do i=1,noobs
1661  dtdf = (bot_temp(i) - sst(i)) / (2.0*atan(1.0) + atan(depthc_i_trans))
1662  temp_cent = sst(i) + dtdf * atan(depthc_i_trans)
1663  temp(i,1) = sst(i)
1664  do k=2,nlevels-1
1665  temp(i,k) = temp_cent + dtdf * atan((depth(k) - mid_depth)/width_trans)
1666  enddo
1667  temp(i,nlevels) = bot_temp(i)
1668 
1669  dsdf = (bot_salt(i) - sss(i)) / (2.0*atan(1.0) + atan(depthc_i_trans))
1670  salt_cent = sss(i) + dsdf * atan(depthc_i_trans)
1671  salt(i,1) = sss(i)
1672  do k=2,nlevels-1
1673  salt(i,k) = salt_cent + dsdf * atan((depth(k) - mid_depth)/width_trans)
1674  enddo
1675  salt(i,nlevels) = bot_salt(i)
1676  enddo
1677 
1678 
1679  num_profiles=0
1680  do i=1,noobs
1681 
1682  data_is_local = .false.
1683 
1684 
1685 ! localized data is within region bounded by halo points
1686 ! (halo size = 1) adjacent to boundary points of computational domain
1687 
1688  if (lon(i) >= x_grid(isc-1) .and. &
1689  lon(i) < x_grid(iec+1) .and. &
1690  lat(i) >= y_grid(jsc-1) .and. &
1691  lat(i) < y_grid(jec+1)) data_is_local = .true.
1692 
1693 
1694  profile_time = set_date(yr(i),mon(i),day(i),hr(i),mm(i),ss(i))
1695 
1696  if ( data_is_local .OR. .NOT.localize_data) then
1697  if (lon(i) >= x_grid(isc)-dx1 .and. &
1698  lon(i) < x_grid(iec)+dx2 .and. &
1699  lat(i) >= y_grid(jsc)-dy1 .and. &
1700  lat(i) < y_grid(jec)+dy2) then
1701 ! nprof_in_comp_domain(cpe) = nprof_in_comp_domain(cpe)+1
1702  endif
1704  if (num_profiles > max_profiles) then
1705  call mpp_error(fatal,'maximum number of profiles exceeded. Resize parameter max_profiles in ocean_obs_mod')
1706  endif
1707 
1708 
1709  profiles(num_profiles)%Model_Grid => grd
1710  profiles(num_profiles)%nvar = 2
1711  profiles(num_profiles)%profile_flag = 0
1712  profiles(num_profiles)%profile_flag_s = 0
1713  profiles(num_profiles)%accepted = .true.
1714  allocate(profiles(num_profiles)%depth(nlevels))
1715  profiles(num_profiles)%depth=depth(1:nlevels)
1716  allocate(profiles(num_profiles)%data_t(nlevels))
1717  profiles(num_profiles)%data_t=temp(i,:)
1718  allocate(profiles(num_profiles)%flag_t(nlevels))
1719  profiles(num_profiles)%flag_t= 0
1720 
1721  allocate(profiles(num_profiles)%data_s(nlevels))
1722  profiles(num_profiles)%data_s=salt(i,:)
1723  allocate(profiles(num_profiles)%flag_s(nlevels))
1724  profiles(num_profiles)%flag_s= 0
1725 
1726  profiles(num_profiles)%probe = 0.
1727  profiles(num_profiles)%levels = nlevels
1728  profiles(num_profiles)%lat = lat(i)
1729  profiles(num_profiles)%lon = lon(i)
1730  allocate(profiles(num_profiles)%ms_t(nlevels))
1731  profiles(num_profiles)%ms_t(:) = min_obs_err_t ! default error variance for temperature
1732  allocate(profiles(num_profiles)%ms_s(nlevels))
1733  profiles(num_profiles)%ms_s(:) = min_obs_err_s ! default error variance for salinity
1734 
1735  profiles(num_profiles)%time = profile_time
1736 
1737 ! calculate interpolation coefficients (make sure to account for grid offsets here!)
1738 ! note that this only works for lat/lon grids
1739 
1740  ri0 = frac_index(lon(i), x_grid(isg-1:ieg+1)) - 1.
1741  rj0 = frac_index(lat(i), y_grid(jsg-1:jeg+1)) - 1.
1742  i0 = floor(ri0)
1743  j0 = floor(rj0)
1744  profiles(num_profiles)%i_index = ri0
1745  profiles(num_profiles)%j_index = rj0
1746  profiles(num_profiles)%accepted = .true.
1747  if (i0 < 0 .or. j0 < 0) then
1748  profiles(num_profiles)%accepted = .false.
1749  endif
1750  if (i0 > ieg+1 .or. j0 > jeg+1) then
1751  call mpp_error(fatal,'grid lat/lon index is out of bounds ')
1752  endif
1753  if (i0 < isc-1 .or. i0 > iec) then
1754  call mpp_error(fatal,'grid lat/lon index is out of bounds ')
1755  endif
1756  if (j0 < jsc-1 .or. j0 > jec) then
1757  call mpp_error(fatal,'grid lat/lon index is out of bounds ')
1758  endif
1759  if (profiles(num_profiles)%accepted ) then
1760  if (grd%mask(i0,j0,1) == 0.0 .or. &
1761  grd%mask(i0+1,j0,1) == 0.0 .or. &
1762  grd%mask(i0,j0+1,1) == 0.0 .or. &
1763  grd%mask(i0+1,j0+1,1) == 0.0) then
1764  profiles(num_profiles)%accepted = .false.
1765  endif
1766  endif
1767 
1768 
1769 
1770  if (profiles(num_profiles)%accepted) then
1771  allocate(profiles(num_profiles)%k_index(profiles(num_profiles)%levels))
1772  do k=1, profiles(num_profiles)%levels
1773  rk0 = frac_index(profiles(num_profiles)%depth(k), grd%z(:))
1774  k0 = floor(rk0)
1775  if ( k0 == -1) then
1776  if (profiles(num_profiles)%depth(k) .ne. missing_value .and. &
1777  profiles(num_profiles)%depth(k) .lt. grd%z(1)) then
1778  k0 = 1
1779  rk0 = 1.0
1780  endif
1781  endif
1782 
1783 
1784  if (k0 .gt. size(grd%z)-1 ) then
1785  write(*,*) 'k0 out of bounds, rk0,k0= ',rk0,k0
1786  write(*,*) 'Z_bound= ',grd%z_bound
1787  write(*,*) 'Profile%depth= ',profiles(num_profiles)%depth
1788  call mpp_error(fatal)
1789  endif
1790 
1791  profiles(num_profiles)%k_index(k) = rk0
1792 
1793  if (profiles(num_profiles)%flag_t(k) .eq. 0) then
1794  if (grd%mask(i0,j0,k0) == 0.0 .or. &
1795  grd%mask(i0+1,j0,k0) == 0.0 .or. &
1796  grd%mask(i0,j0+1,k0) == 0.0 .or. &
1797  grd%mask(i0+1,j0+1,k0) == 0.0) then
1798  profiles(num_profiles)%flag_t(k) = 1
1799  endif
1800  if (grd%mask(i0,j0,k0+1) == 0.0 .or. &
1801  grd%mask(i0+1,j0,k0+1) == 0.0 .or. &
1802  grd%mask(i0,j0+1,k0+1) == 0.0 .or. &
1803  grd%mask(i0+1,j0+1,k0+1) == 0.0) then
1804  profiles(num_profiles)%flag_t(k) = 1
1805  endif
1806  if (profiles(num_profiles)%data_t(k) == missing_value &
1807  .or. profiles(num_profiles)%depth(k) == missing_value) then
1808  profiles(num_profiles)%flag_t(k) = 1
1809  endif
1810  endif
1811 
1812  enddo
1813  endif
1814  endif
1815 
1816  enddo
1817 
1818 ! a = nprof_in_comp_domain(cpe)
1819 
1820 ! call mpp_broadcast(nprof_in_comp_domain(cpe),cpe)
1821 
1822 ! call mpp_sum(a)
1823 
1824 ! write(stdout(),*) 'A grand total of ',a,' profiles satisify acceptance criteria'
1825 
1826 ! do i=0,mpp_npes()-1
1827 ! write(stdout(),*) 'pe=',i,'profile count=',nprof_in_comp_domain(i)
1828 ! enddo
1829 
1830  end subroutine create_ideal_profiles
1831 
1832 
1833  subroutine nullify_obs_prof(profile)
1835  type(ocean_profile_type), intent(inout) :: profile
1836 
1837 
1838  profile%nvar = 0
1839  profile%project=0.
1840  profile%probe=0.
1841  profile%ref_inst=0.
1842  profile%wod_cast_num=0
1843  profile%fix_depth=0.
1844  profile%ocn_vehicle=0.
1845  profile%database_id=0.
1846  profile%levels=0
1847  profile%profile_flag=-1
1848  profile%profile_flag_s=-1
1849  profile%lon=-1.0e10
1850  profile%lat=-1.0e10
1851  profile%accepted=.false.
1852  profile%nlinks=0
1853  if (ASSOCIATED(profile%next)) profile%next=>null()
1854  if (ASSOCIATED(profile%depth)) profile%depth=>null()
1855  if (ASSOCIATED(profile%data_t)) profile%data_t=>null()
1856  if (ASSOCIATED(profile%data_s)) profile%data_s=>null()
1857  if (ASSOCIATED(profile%flag_t)) profile%flag_t=>null()
1858  if (ASSOCIATED(profile%flag_s)) profile%flag_s=>null()
1859  profile%temp_err=0.0
1860  profile%salt_err=0.0
1861  if (ASSOCIATED(profile%ms_t)) profile%ms_t=>null()
1862  if (ASSOCIATED(profile%ms_s)) profile%ms_s=>null()
1863  profile%time = set_time(0,0)
1864  profile%yyyy = 0
1865  profile%mmdd = 0
1866  if (ASSOCIATED(profile%model_time)) profile%model_time=>null()
1867  if (ASSOCIATED(profile%model_grid)) profile%model_grid=>null()
1868  if (ASSOCIATED(profile%k_index)) profile%k_index=>null()
1869  profile%i_index=-1.0
1870  profile%j_index=-1.0
1871  profile%tdiff = set_time(0,0)
1872 
1873  end subroutine nullify_obs_prof
1874 
1875 end module oda_core_mod
Definition: fms.F90:20
integer, parameter, private sfc_file
Definition: oda_core.F90:111
subroutine, public get_date(time, year, month, day, hour, minute, second, tick, err_msg)
subroutine nullify_obs_prof(profile)
Definition: oda_core.F90:1834
subroutine backward_obs_profile(Obs, model_t, model_s)
Definition: oda_core.F90:1107
real, private min_obs_err_t
Definition: oda_core.F90:117
type(ocean_profile_type), dimension(:), allocatable, target, save, private profiles
Definition: oda_core.F90:119
subroutine create_ideal_profiles(localize)
Definition: oda_core.F90:1573
subroutine copy_obs_sfc(Obs_in, Obs_out)
Definition: oda_core.F90:1335
integer, save, private iec
Definition: oda_core.F90:124
subroutine adjust_obs_error_sfc(Diff)
Definition: oda_core.F90:1370
real, dimension(:), allocatable x_grid
Definition: oda_core.F90:147
subroutine mult_obs_i_mse_profile(Obs)
Definition: oda_core.F90:1379
type(time_type), dimension(0:100), public time_window
Definition: oda_core.F90:139
integer, save, private isd
Definition: oda_core.F90:124
subroutine add_tidal_error(Prof)
Definition: oda_core.F90:1544
real, private max_misfit
Definition: oda_core.F90:102
subroutine backward_obs_sfc(Obs, model)
Definition: oda_core.F90:1186
integer, save, private nk
Definition: oda_core.F90:126
subroutine adjust_obs_error_profile(Prof)
Definition: oda_core.F90:1344
integer, save, private jed
Definition: oda_core.F90:124
integer, parameter, public model_ocean
subroutine, public get_obs(model_time, Prof, Sfc, nprof, nsfc)
Definition: oda_core.F90:835
subroutine, public get_field_info(n, fld_type, fld_name, model, num_methods)
integer, save, private num_sfc_obs
Definition: oda_core.F90:122
real, private eta_tide_const
Definition: oda_core.F90:117
subroutine, public oda_core_init(Domain, Grid, localize)
Definition: oda_core.F90:922
integer, parameter, public max_links
Maximum number of records per profile for storage for profiles.
Definition: oda_types.F90:45
type(time_type) function, public get_cal_time(time_increment, units, calendar, permit_calendar_conversion)
real, parameter, private depth_max
Definition: oda_core.F90:105
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
real, parameter, private salt_min
Definition: oda_core.F90:107
real, private min_obs_err_s
Definition: oda_core.F90:117
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
real, parameter, private temp_min
Definition: oda_core.F90:106
subroutine assign_forward_model_profile(Obs1, Obs2)
Definition: oda_core.F90:1194
real, save, private max_prof_spacing
Definition: oda_core.F90:132
integer, save, private jeg
Definition: oda_core.F90:125
type(ocean_surface_type), dimension(:), allocatable, target, save, private sfc_obs
Definition: oda_core.F90:120
real, parameter, public radian
Equal to RAD_TO_DEG for backward compatability. [rad/deg].
Definition: constants.F90:121
integer, save, private num_profiles
Definition: oda_core.F90:122
subroutine diff_obs_sfc(prof1, prof2, Diff)
Definition: oda_core.F90:1261
subroutine forward_obs_sfc(Sfc, Guess, Diff)
Definition: oda_core.F90:1097
subroutine, public write_ocean_data_init()
real, parameter, private temp_max
Definition: oda_core.F90:106
integer max_sfc_obs
Definition: oda_core.F90:109
real, dimension(:), allocatable y_grid
Definition: oda_core.F90:147
integer, save, private jsg
Definition: oda_core.F90:125
real salt_obs_rmse
Definition: oda_core.F90:150
type(grid_type), pointer grd
Definition: oda_core.F90:143
integer, save, private isg
Definition: oda_core.F90:125
subroutine assign_forward_model_sfc(Obs1, Obs2)
Definition: oda_core.F90:1213
integer, save, private ied
Definition: oda_core.F90:124
integer ndebug
Definition: oda_core.F90:158
logical debug
Definition: oda_core.F90:156
logical localize_data
Definition: oda_core.F90:154
real, parameter, public missing_value
Definition: oda_types.F90:69
integer, save, private isc
Definition: oda_core.F90:124
real temp_obs_rmse
Definition: oda_core.F90:149
integer, save, private jsc
Definition: oda_core.F90:124
real function, public frac_index(value, array)
Definition: axis_utils.F90:368
#define max(a, b)
Definition: mosaic_util.h:33
subroutine copy_obs_prof(obs_in, obs_out)
Definition: oda_core.F90:1269
subroutine diff_obs_profile(prof1, prof2, Diff)
Definition: oda_core.F90:1222
real, save, private min_prof_depth
Definition: oda_core.F90:129
real, dimension(:,:,:), allocatable, save, private sum_wgt
Definition: oda_core.F90:137
subroutine init_observations(localize)
Definition: oda_core.F90:1428
integer max_profiles
Definition: oda_core.F90:109
integer, save, private jec
Definition: oda_core.F90:124
integer, parameter, private idealized_profiles
Definition: oda_core.F90:111
real, parameter, private depth_min
Definition: oda_core.F90:105
integer, parameter, private profile_file
Definition: oda_core.F90:111
#define min(a, b)
Definition: mosaic_util.h:32
subroutine, public get_time(Time, seconds, days, ticks, err_msg)
integer, dimension(:), allocatable nprof_in_comp_domain
Definition: oda_core.F90:160
subroutine, public purge_obs()
Definition: oda_core.F90:964
subroutine, public get_field_methods(n, methods)
subroutine, public open_profile_dataset(filename, localize)
Definition: oda_core.F90:236
integer, save, private ieg
Definition: oda_core.F90:125
integer, save, private jsd
Definition: oda_core.F90:124
character(len=len(cs)) function lowercase(cs)
Definition: oda_core.F90:1415
Derived type containing the data.
subroutine forward_obs_profile(Model_obs, fg_t, fg_s)
Definition: oda_core.F90:971
subroutine mult_obs_i_mse_sfc(a, Obs)
Definition: oda_core.F90:1399
logical add_tidal_aliasing
Definition: oda_core.F90:152
real, parameter, private salt_max
Definition: oda_core.F90:107
real(fp), parameter, public pi
integer, parameter, private max_files
Definition: oda_core.F90:110