MOM6
MOM_oda_driver.F90
1 !> Interfaces for MOM6 ensembles and data assimilation.
3 
4  ! This file is part of MOM6. see LICENSE.md for the license.
5 
6 use mpp_mod, only : stdout, stdlog, mpp_error, npes=>mpp_npes,pe=>mpp_pe
7 use mpp_mod, only : set_current_pelist => mpp_set_current_pelist
8 use mpp_mod, only : set_root_pe => mpp_set_root_pe
9 use mpp_mod, only : mpp_sync_self, mpp_sum, get_pelist=>mpp_get_current_pelist, mpp_root_pe
10 use mpp_mod, only : set_stack_size=>mpp_set_stack_size, broadcast=>mpp_broadcast
11 use mpp_io_mod, only : io_set_stack_size=>mpp_io_set_stack_size
12 use mpp_io_mod, only : mpp_single,mpp_multi
13 use mpp_domains_mod, only : domain2d, mpp_global_field
14 use mpp_domains_mod, only : mpp_get_compute_domain, mpp_get_data_domain
15 use mpp_domains_mod, only : mpp_redistribute, mpp_broadcast_domain
16 use mpp_domains_mod, only : set_domains_stack_size=>mpp_domains_set_stack_size
17 use diag_manager_mod, only : register_diag_field, diag_axis_init, send_data
18 use ensemble_manager_mod, only : get_ensemble_id, get_ensemble_size
19 use ensemble_manager_mod, only : get_ensemble_pelist, get_ensemble_filter_pelist
20 use time_manager_mod, only : time_type, decrement_time, increment_time
21 use time_manager_mod, only : get_date, operator(>=),operator(/=),operator(==),operator(<)
22 use constants_mod, only : radius, epsln
23 ! ODA Modules
25 use ocean_da_core_mod, only : ocean_da_core_init, get_profiles
26 !use eakf_oda_mod, only : ensemble_filter
27 use write_ocean_obs_mod, only : open_profile_file
28 use write_ocean_obs_mod, only : write_profile,close_profile_file
29 use kdtree, only : kd_root !# JEDI
30 ! MOM Modules
31 use mom_io, only : slasher, mom_read_data
32 use mom_diag_mediator, only : diag_ctrl, set_axes_info
33 use mom_error_handler, only : fatal, warning, mom_error, mom_mesg, is_root_pe
34 use mom_get_input, only : get_mom_input, directories
35 use mom_grid, only : ocean_grid_type, mom_grid_init
36 use mom_grid_initialize, only : set_grid_metrics
37 use mom_hor_index, only : hor_index_type, hor_index_init
38 use mom_dyn_horgrid, only : dyn_horgrid_type, create_dyn_horgrid, destroy_dyn_horgrid
39 use mom_transcribe_grid, only : copy_dyngrid_to_mom_grid, copy_mom_grid_to_dyngrid
40 use mom_fixed_initialization, only : mom_initialize_fixed, mom_initialize_topography
41 use mom_coord_initialization, only : mom_initialize_coord
43 use mom_string_functions, only : lowercase
44 use mom_ale, only : ale_cs, ale_initthicknesstocoord, ale_init, ale_updateverticalgridtype
45 use mom_domains, only : mom_domains_init, mom_domain_type, clone_mom_domain
46 use mom_remapping, only : remapping_cs, initialize_remapping, remapping_core_h
47 use mom_regridding, only : regridding_cs, initialize_regridding
48 use mom_regridding, only : regridding_main, set_regrid_params
49 use mom_unit_scaling, only : unit_scale_type, unit_scaling_init
51 use mom_verticalgrid, only : verticalgrid_type, verticalgridinit
52 
53 implicit none ; private
54 
55 public :: init_oda, oda_end, set_prior_tracer, get_posterior_tracer
56 public :: set_analysis_time, oda, save_obs_diff, apply_oda_tracer_increments
57 
58 #include <MOM_memory.h>
59 
60 !> Control structure that contains a transpose of the ocean state across ensemble members.
61 type, public :: oda_cs ; private
62  type(ocean_control_struct), pointer :: ocean_prior=> null() !< ensemble ocean prior states in DA space
63  type(ocean_control_struct), pointer :: ocean_posterior=> null() !< ensemble ocean posterior states
64  !! or increments to prior in DA space
65  integer :: nk !< number of vertical layers used for DA
66  type(ocean_grid_type), pointer :: grid => null() !< MOM6 grid type and decomposition for the DA
67  type(ptr_mpp_domain), pointer, dimension(:) :: domains => null() !< Pointer to mpp_domain objects
68  !! for ensemble members
69  type(verticalgrid_type), pointer :: gv => null() !< vertical grid for DA
70  type(unit_scale_type), pointer :: &
71  us => null() !< structure containing various unit conversion factors for DA
72 
73  type(domain2d), pointer :: mpp_domain => null() !< Pointer to a mpp domain object for DA
74  type(grid_type), pointer :: oda_grid !< local tracer grid
75  real, pointer, dimension(:,:,:) :: h => null() !<layer thicknesses [H ~> m or kg m-2] for DA
76  type(thermo_var_ptrs), pointer :: tv => null() !< pointer to thermodynamic variables
77  integer :: ni !< global i-direction grid size
78  integer :: nj !< global j-direction grid size
79  logical :: reentrant_x !< grid is reentrant in the x direction
80  logical :: reentrant_y !< grid is reentrant in the y direction
81  logical :: tripolar_n !< grid is folded at its north edge
82  logical :: symmetric !< Values at C-grid locations are symmetric
83  integer :: assim_method !< Method: NO_ASSIM,EAKF_ASSIM or OI_ASSIM
84  integer :: ensemble_size !< Size of the ensemble
85  integer :: ensemble_id = 0 !< id of the current ensemble member
86  integer, pointer, dimension(:,:) :: ensemble_pelist !< PE list for ensemble members
87  integer, pointer, dimension(:) :: filter_pelist !< PE list for ensemble members
88  integer :: assim_frequency !< analysis interval in hours
89  ! Profiles local to the analysis domain
90  type(ocean_profile_type), pointer :: profiles => null() !< pointer to linked list of all available profiles
91  type(ocean_profile_type), pointer :: cprofiles => null()!< pointer to linked list of current profiles
92  type(kd_root), pointer :: kdroot => null() !< A structure for storing nearest neighbors
93  type(ale_cs), pointer :: ale_cs=>null() !< ALE control structure for DA
94  logical :: use_ale_algorithm !< true is using ALE remapping
95  type(regridding_cs) :: regridcs !< ALE control structure for regridding
96  type(remapping_cs) :: remapcs !< ALE control structure for remapping
97  type(time_type) :: time !< Current Analysis time
98  type(diag_ctrl) :: diag_cs !<Diagnostics control structure
99 end type oda_cs
100 
101 !> A structure with a pointer to a domain2d, to allow for the creation of arrays of pointers.
103  type(domain2d), pointer :: mpp_domain => null() !< pointer to an mpp domain2d
104 end type ptr_mpp_domain
105 
106 !>@{ DA parameters
107 integer, parameter :: no_assim = 0, oi_assim=1, eakf_assim=2
108 !>@}
109 
110 contains
111 
112 !> initialize First_guess (prior) and Analysis grid
113 !! information for all ensemble members
114 subroutine init_oda(Time, G, GV, CS)
116  type(time_type), intent(in) :: Time !< The current model time.
117  type(ocean_grid_type), pointer :: G !< domain and grid information for ocean model
118  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
119  type(oda_cs), pointer, intent(inout) :: CS !< The DA control structure
120 
121 ! Local variables
122  type(thermo_var_ptrs) :: tv_dummy
123  type(dyn_horgrid_type), pointer :: dG=> null()
124  type(hor_index_type), pointer :: HI=> null()
125  type(directories) :: dirs
126 
127  type(grid_type), pointer :: T_grid !< global tracer grid
128  real, dimension(:,:), allocatable :: global2D, global2D_old
129  real, dimension(:), allocatable :: lon1D, lat1D, glon1D, glat1D
130  type(param_file_type) :: PF
131  integer :: n, m, k, i, j, nk
132  integer :: is,ie,js,je,isd,ied,jsd,jed
133  integer :: stdout_unit
134  character(len=32) :: assim_method
135  integer :: npes_pm, ens_info(6), ni, nj
136  character(len=128) :: mesg
137  character(len=32) :: fldnam
138  character(len=30) :: coord_mode
139  character(len=200) :: inputdir, basin_file
140  logical :: reentrant_x, reentrant_y, tripolar_N, symmetric
141 
142  if (associated(cs)) call mpp_error(fatal,'Calling oda_init with associated control structure')
143  allocate(cs)
144 ! Use ens1 parameters , this could be changed at a later time
145 ! if it were desirable to have alternate parameters, e.g. for the grid
146 ! for the analysis
147  call get_mom_input(pf,dirs,ensemble_num=0)
148 
149  call unit_scaling_init(pf, cs%US)
150 
151  call get_param(pf, "MOM", "ASSIM_METHOD", assim_method, &
152  "String which determines the data assimilation method "//&
153  "Valid methods are: \'EAKF\',\'OI\', and \'NO_ASSIM\'", default='NO_ASSIM')
154  call get_param(pf, "MOM", "ASSIM_FREQUENCY", cs%assim_frequency, &
155  "data assimilation frequency in hours")
156  call get_param(pf, "MOM", "USE_REGRIDDING", cs%use_ALE_algorithm , &
157  "If True, use the ALE algorithm (regridding/remapping).\n"//&
158  "If False, use the layered isopycnal algorithm.", default=.false. )
159  call get_param(pf, "MOM", "REENTRANT_X", cs%reentrant_x, &
160  "If true, the domain is zonally reentrant.", default=.true.)
161  call get_param(pf, "MOM", "REENTRANT_Y", cs%reentrant_y, &
162  "If true, the domain is meridionally reentrant.", &
163  default=.false.)
164  call get_param(pf,"MOM", "TRIPOLAR_N", cs%tripolar_N, &
165  "Use tripolar connectivity at the northern edge of the "//&
166  "domain. With TRIPOLAR_N, NIGLOBAL must be even.", &
167  default=.false.)
168  call get_param(pf,"MOM", "NIGLOBAL", cs%ni, &
169  "The total number of thickness grid points in the "//&
170  "x-direction in the physical domain.")
171  call get_param(pf,"MOM", "NJGLOBAL", cs%nj, &
172  "The total number of thickness grid points in the "//&
173  "y-direction in the physical domain.")
174  call get_param(pf, 'MOM', "INPUTDIR", inputdir)
175  inputdir = slasher(inputdir)
176 
177  select case(lowercase(trim(assim_method)))
178  case('eakf')
179  cs%assim_method = eakf_assim
180  case('oi')
181  cs%assim_method = oi_assim
182  case('no_assim')
183  cs%assim_method = no_assim
184  case default
185  call mpp_error(fatal,'Invalid assimilation method provided')
186  end select
187 
188  ens_info = get_ensemble_size()
189  cs%ensemble_size = ens_info(1)
190  npes_pm=ens_info(3)
191  cs%ensemble_id = get_ensemble_id()
192  !! Switch to global pelist
193  allocate(cs%ensemble_pelist(cs%ensemble_size,npes_pm))
194  allocate(cs%filter_pelist(cs%ensemble_size*npes_pm))
195  call get_ensemble_pelist(cs%ensemble_pelist,'ocean')
196  call get_ensemble_filter_pelist(cs%filter_pelist,'ocean')
197 
198  call set_current_pelist(cs%filter_pelist)
199 
200  allocate(cs%domains(cs%ensemble_size))
201  cs%domains(cs%ensemble_id)%mpp_domain => g%Domain%mpp_domain
202  do n=1,cs%ensemble_size
203  if (.not. associated(cs%domains(n)%mpp_domain)) allocate(cs%domains(n)%mpp_domain)
204  call set_root_pe(cs%ensemble_pelist(n,1))
205  call mpp_broadcast_domain(cs%domains(n)%mpp_domain)
206  enddo
207  call set_root_pe(cs%filter_pelist(1))
208  allocate(cs%Grid)
209  ! params NIHALO_ODA, NJHALO_ODA set the DA halo size
210  call mom_domains_init(cs%Grid%Domain,pf,param_suffix='_ODA')
211  allocate(hi)
212  call hor_index_init(cs%Grid%Domain, hi, pf, &
213  local_indexing=.false.) ! Use global indexing for DA
214  call verticalgridinit( pf, cs%GV, cs%US )
215  allocate(dg)
216  call create_dyn_horgrid(dg, hi)
217  call clone_mom_domain(cs%Grid%Domain, dg%Domain,symmetric=.false.)
218  call set_grid_metrics(dg,pf)
219  call mom_initialize_topography(dg%bathyT,dg%max_depth,dg,pf)
220  call mom_initialize_coord(cs%GV, cs%US, pf, .false., &
221  dirs%output_directory, tv_dummy, dg%max_depth)
222  call ale_init(pf, cs%GV, cs%US, dg%max_depth, cs%ALE_CS)
223  call mom_grid_init(cs%Grid, pf, global_indexing=.true.)
224  call ale_updateverticalgridtype(cs%ALE_CS, cs%GV)
225  call copy_dyngrid_to_mom_grid(dg, cs%Grid, cs%US)
226  cs%mpp_domain => cs%Grid%Domain%mpp_domain
227  cs%Grid%ke = cs%GV%ke
228  cs%nk = cs%GV%ke
229  ! initialize storage for prior and posterior
230  allocate(cs%Ocean_prior)
231  call init_ocean_ensemble(cs%Ocean_prior,cs%Grid,cs%GV,cs%ensemble_size)
232  allocate(cs%Ocean_posterior)
233  call init_ocean_ensemble(cs%Ocean_posterior,cs%Grid,cs%GV,cs%ensemble_size)
234  allocate(cs%tv)
235 
236  call get_param(pf, 'oda_driver', "REGRIDDING_COORDINATE_MODE", coord_mode, &
237  "Coordinate mode for vertical regridding.", &
238  default="ZSTAR", fail_if_missing=.false.)
239  call initialize_regridding(cs%regridCS, cs%GV, cs%US, dg%max_depth,pf,'oda_driver',coord_mode,'','')
240  call initialize_remapping(cs%remapCS,'PLM')
241  call set_regrid_params(cs%regridCS, min_thickness=0.)
242  call mpp_get_data_domain(g%Domain%mpp_domain,isd,ied,jsd,jed)
243  if (.not. associated(cs%h)) then
244  allocate(cs%h(isd:ied,jsd:jed,cs%GV%ke)); cs%h(:,:,:)=0.0
245  ! assign thicknesses
246  call ale_initthicknesstocoord(cs%ALE_CS,g,cs%GV,cs%h)
247  endif
248  allocate(cs%tv%T(isd:ied,jsd:jed,cs%GV%ke)); cs%tv%T(:,:,:)=0.0
249  allocate(cs%tv%S(isd:ied,jsd:jed,cs%GV%ke)); cs%tv%S(:,:,:)=0.0
250 
251  call set_axes_info(cs%Grid, cs%GV, cs%US, pf, cs%diag_cs, set_vertical=.true.)
252 
253  call mpp_get_data_domain(cs%mpp_domain,isd,ied,jsd,jed)
254  allocate(cs%oda_grid)
255  cs%oda_grid%x => cs%Grid%geolonT
256  cs%oda_grid%y => cs%Grid%geolatT
257 
258  call get_param(pf, 'oda_driver', "BASIN_FILE", basin_file, &
259  "A file in which to find the basin masks, in variable 'basin'.", &
260  default="basin.nc")
261  basin_file = trim(inputdir) // trim(basin_file)
262  allocate(cs%oda_grid%basin_mask(isd:ied,jsd:jed))
263  cs%oda_grid%basin_mask(:,:) = 0.0
264  call mom_read_data(basin_file,'basin',cs%oda_grid%basin_mask,cs%Grid%domain, timelevel=1)
265 
266 ! get global grid information from ocean_model
267  allocate(t_grid)
268  allocate(t_grid%x(cs%ni,cs%nj))
269  allocate(t_grid%y(cs%ni,cs%nj))
270  allocate(t_grid%basin_mask(cs%ni,cs%nj))
271  call mpp_global_field(cs%mpp_domain, cs%Grid%geolonT, t_grid%x)
272  call mpp_global_field(cs%mpp_domain, cs%Grid%geolatT, t_grid%y)
273  call mpp_global_field(cs%mpp_domain, cs%oda_grid%basin_mask, t_grid%basin_mask)
274  t_grid%ni = cs%ni
275  t_grid%nj = cs%nj
276  t_grid%nk = cs%nk
277  allocate(t_grid%mask(cs%ni,cs%nj,cs%nk))
278  allocate(t_grid%z(cs%ni,cs%nj,cs%nk))
279  allocate(global2d(cs%ni,cs%nj))
280  allocate(global2d_old(cs%ni,cs%nj))
281  t_grid%mask(:,:,:) = 0.0
282  t_grid%z(:,:,:) = 0.0
283 
284  do k = 1, cs%nk
285  call mpp_global_field(g%Domain%mpp_domain, cs%h(:,:,k), global2d)
286  do i=1, cs%ni; do j=1, cs%nj
287  if ( global2d(i,j) > 1 ) then
288  t_grid%mask(i,j,k) = 1.0
289  endif
290  enddo ; enddo
291  if (k == 1) then
292  t_grid%z(:,:,k) = global2d/2
293  else
294  t_grid%z(:,:,k) = t_grid%z(:,:,k-1) + (global2d + global2d_old)/2
295  endif
296  global2d_old = global2d
297  enddo
298 
299  call ocean_da_core_init(cs%mpp_domain, t_grid, cs%Profiles, time)
300 
301  cs%Time=time
302  !! switch back to ensemble member pelist
303  call set_current_pelist(cs%ensemble_pelist(cs%ensemble_id,:))
304 end subroutine init_oda
305 
306 !> Copy ensemble member tracers to ensemble vector.
307 subroutine set_prior_tracer(Time, G, GV, h, tv, CS)
308  type(time_type), intent(in) :: Time !< The current model time
309  type(ocean_grid_type), pointer :: G !< domain and grid information for ocean model
310  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
311  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
312  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
313 
314  type(oda_cs), pointer :: CS !< ocean DA control structure
315  real, dimension(:,:,:), allocatable :: T, S
316  type(ocean_grid_type), pointer :: Grid=>null()
317  integer :: i,j, m, n, ss
318  integer :: is, ie, js, je
319  integer :: isc, iec, jsc, jec
320  integer :: isd, ied, jsd, jed
321  integer :: id
322  logical :: used
323 
324  ! return if not time for analysis
325  if (time < cs%Time) return
326 
327  if (.not. associated(cs%Grid)) call mom_error(fatal,'ODA_CS ensemble horizontal grid not associated')
328  if (.not. associated(cs%GV)) call mom_error(fatal,'ODA_CS ensemble vertical grid not associated')
329 
330  !! switch to global pelist
331  call set_current_pelist(cs%filter_pelist)
332  call mom_mesg('Setting prior')
333 
334  isc=cs%Grid%isc;iec=cs%Grid%iec;jsc=cs%Grid%jsc;jec=cs%Grid%jec
335  call mpp_get_compute_domain(cs%domains(cs%ensemble_id)%mpp_domain,is,ie,js,je)
336  call mpp_get_data_domain(cs%domains(cs%ensemble_id)%mpp_domain,isd,ied,jsd,jed)
337  allocate(t(isd:ied,jsd:jed,cs%nk))
338  allocate(s(isd:ied,jsd:jed,cs%nk))
339 
340  do j=js,je; do i=is,ie
341  call remapping_core_h(cs%remapCS, gv%ke, h(i,j,:), tv%T(i,j,:), &
342  cs%nk, cs%h(i,j,:), t(i,j,:))
343  call remapping_core_h(cs%remapCS, gv%ke, h(i,j,:), tv%S(i,j,:), &
344  cs%nk, cs%h(i,j,:), s(i,j,:))
345  enddo ; enddo
346 
347  do m=1,cs%ensemble_size
348  call mpp_redistribute(cs%domains(m)%mpp_domain, t,&
349  cs%mpp_domain, cs%Ocean_prior%T(:,:,:,m), complete=.true.)
350  call mpp_redistribute(cs%domains(m)%mpp_domain, s,&
351  cs%mpp_domain, cs%Ocean_prior%S(:,:,:,m), complete=.true.)
352  enddo
353  deallocate(t,s)
354 
355  !! switch back to ensemble member pelist
356  call set_current_pelist(cs%ensemble_pelist(cs%ensemble_id,:))
357 
358  return
359 
360 end subroutine set_prior_tracer
361 
362 !> Returns posterior adjustments or full state
363 !!Note that only those PEs associated with an ensemble member receive data
364 subroutine get_posterior_tracer(Time, CS, h, tv, increment)
365  type(time_type), intent(in) :: Time !< the current model time
366  type(oda_cs), pointer :: CS !< ocean DA control structure
367  real, dimension(:,:,:), pointer :: h !< Layer thicknesses [H ~> m or kg m-2]
368  type(thermo_var_ptrs), pointer :: tv !< A structure pointing to various thermodynamic variables
369  logical, optional, intent(in) :: increment !< True if returning increment only
370 
371  type(ocean_control_struct), pointer :: Ocean_increment=>null()
372  integer :: i, j, m
373  logical :: used, get_inc
374 
375  ! return if not analysis time (retain pointers for h and tv)
376  if (time < cs%Time) return
377 
378 
379  !! switch to global pelist
380  call set_current_pelist(cs%filter_pelist)
381  call mom_mesg('Getting posterior')
382 
383  get_inc = .true.
384  if (present(increment)) get_inc = increment
385 
386  if (get_inc) then
387  allocate(ocean_increment)
388  call init_ocean_ensemble(ocean_increment,cs%Grid,cs%GV,cs%ensemble_size)
389  ocean_increment%T = cs%Ocean_posterior%T - cs%Ocean_prior%T
390  ocean_increment%S = cs%Ocean_posterior%S - cs%Ocean_prior%S
391  endif
392  do m=1,cs%ensemble_size
393  if (get_inc) then
394  call mpp_redistribute(cs%mpp_domain, ocean_increment%T(:,:,:,m), &
395  cs%domains(m)%mpp_domain, cs%tv%T, complete=.true.)
396  call mpp_redistribute(cs%mpp_domain, ocean_increment%S(:,:,:,m), &
397  cs%domains(m)%mpp_domain, cs%tv%S, complete=.true.)
398  else
399  call mpp_redistribute(cs%mpp_domain, cs%Ocean_posterior%T(:,:,:,m), &
400  cs%domains(m)%mpp_domain, cs%tv%T, complete=.true.)
401  call mpp_redistribute(cs%mpp_domain, cs%Ocean_posterior%S(:,:,:,m), &
402  cs%domains(m)%mpp_domain, cs%tv%S, complete=.true.)
403  endif
404  enddo
405 
406  tv => cs%tv
407  h => cs%h
408  !! switch back to ensemble member pelist
409  call set_current_pelist(cs%ensemble_pelist(cs%ensemble_id,:))
410 
411 end subroutine get_posterior_tracer
412 
413 !> Gather observations and sall ODA routines
414 subroutine oda(Time, CS)
415  type(time_type), intent(in) :: Time !< the current model time
416  type(oda_cs), intent(inout) :: CS !< the ocean DA control structure
417 
418  integer :: i, j
419  integer :: m
420  integer :: yr, mon, day, hr, min, sec
421 
422  if ( time >= cs%Time ) then
423 
424  !! switch to global pelist
425  call set_current_pelist(cs%filter_pelist)
426 
427  call get_profiles(time, cs%Profiles, cs%CProfiles)
428 #ifdef ENABLE_ECDA
429  call ensemble_filter(cs%Ocean_prior, cs%Ocean_posterior, cs%CProfiles, cs%kdroot, cs%mpp_domain, cs%oda_grid)
430 #endif
431 
432  !! switch back to ensemble member pelist
433  call set_current_pelist(cs%ensemble_pelist(cs%ensemble_id,:))
434 
435  endif
436 
437  return
438 end subroutine oda
439 
440 !> Finalize DA module
441 subroutine oda_end(CS)
442  type(oda_cs), intent(inout) :: CS !< the ocean DA control structure
443 
444 end subroutine oda_end
445 
446 !> Initialize DA module
447 subroutine init_ocean_ensemble(CS,Grid,GV,ens_size)
448  type(ocean_control_struct), pointer :: CS !< Pointer to ODA control structure
449  type(ocean_grid_type), pointer :: Grid !< Pointer to ocean analysis grid
450  type(verticalgrid_type), pointer :: GV !< Pointer to DA vertical grid
451  integer, intent(in) :: ens_size !< ensemble size
452 
453  integer :: n,is,ie,js,je,nk
454 
455  nk=gv%ke
456  is=grid%isd;ie=grid%ied
457  js=grid%jsd;je=grid%jed
458  cs%ensemble_size=ens_size
459  allocate(cs%T(is:ie,js:je,nk,ens_size))
460  allocate(cs%S(is:ie,js:je,nk,ens_size))
461  allocate(cs%SSH(is:ie,js:je,ens_size))
462 ! allocate(CS%id_t(ens_size));CS%id_t(:)=-1
463 ! allocate(CS%id_s(ens_size));CS%id_s(:)=-1
464 ! allocate(CS%U(is:ie,js:je,nk,ens_size))
465 ! allocate(CS%V(is:ie,js:je,nk,ens_size))
466 ! allocate(CS%id_u(ens_size));CS%id_u(:)=-1
467 ! allocate(CS%id_v(ens_size));CS%id_v(:)=-1
468 ! allocate(CS%id_ssh(ens_size));CS%id_ssh(:)=-1
469 
470  return
471 end subroutine init_ocean_ensemble
472 
473 !> Set the next analysis time
474 subroutine set_analysis_time(Time,CS)
475  type(time_type), intent(in) :: Time !< the current model time
476  type(oda_cs), pointer, intent(inout) :: CS !< the DA control structure
477 
478  character(len=160) :: mesg ! The text of an error message
479  integer :: yr, mon, day, hr, min, sec
480 
481  if (time >= cs%Time) then
482  cs%Time=increment_time(cs%Time,cs%assim_frequency*3600)
483 
484  call get_date(time, yr, mon, day, hr, min, sec)
485  write(mesg,*) 'Model Time: ', yr, mon, day, hr, min, sec
486  call mom_mesg("set_analysis_time: "//trim(mesg))
487  call get_date(cs%time, yr, mon, day, hr, min, sec)
488  write(mesg,*) 'Assimilation Time: ', yr, mon, day, hr, min, sec
489  call mom_mesg("set_analysis_time: "//trim(mesg))
490  endif
491  if (cs%Time < time) then
492  call mom_error(fatal, " set_analysis_time: " // &
493  "assimilation interval appears to be shorter than " // &
494  "the model timestep")
495  endif
496  return
497 
498 end subroutine set_analysis_time
499 
500 !> Write observation differences to a file
501 subroutine save_obs_diff(filename,CS)
502  character(len=*), intent(in) :: filename !< name of output file
503  type(oda_cs), pointer, intent(in) :: CS !< pointer to DA control structure
504 
505  integer :: fid ! profile file handle
506  type(ocean_profile_type), pointer :: Prof=>null()
507 
508  fid = open_profile_file(trim(filename), nvar=2, thread=mpp_single, fset=mpp_single)
509  prof=>cs%CProfiles
510 
511  !! switch to global pelist
512  !call set_current_pelist(CS%filter_pelist)
513 
514  do while (associated(prof))
515  call write_profile(fid,prof)
516  prof=>prof%cnext
517  enddo
518  call close_profile_file(fid)
519 
520  !! switch back to ensemble member pelist
521  !call set_current_pelist(CS%ensemble_pelist(CS%ensemble_id,:))
522 
523  return
524 end subroutine save_obs_diff
525 
526 
527 !> Apply increments to tracers
528 subroutine apply_oda_tracer_increments(dt,G,tv,h,CS)
529  real, intent(in) :: dt !< The tracer timestep [s]
530  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
531  type(thermo_var_ptrs), intent(inout) :: tv !< A structure pointing to various thermodynamic variables
532  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
533  intent(in) :: h !< layer thickness [H ~> m or kg m-2]
534  type(oda_cs), intent(inout) :: CS !< the data assimilation structure
536 end subroutine apply_oda_tracer_increments
537 
538 !> \namespace MOM_oda_driver_mod
539 !!
540 !! \section section_ODA The Ocean data assimilation (DA) and Ensemble Framework
541 !!
542 !! The DA framework implements ensemble capability in MOM6. Currently, this framework
543 !! is enabled using the cpp directive ENSEMBLE_OCEAN. The ensembles need to be generated
544 !! at the level of the calling routine for oda_init or above. The ensemble instances may
545 !! exist on overlapping or non-overlapping processors. The ensemble information is accessed
546 !! via the FMS ensemble manager. An independent PE layout is used to gather (prior) ensemble
547 !! member information where this information is stored in the ODA control structure. This
548 !! module was developed in collaboration with Feiyu Lu and Tony Rosati in the GFDL prediction
549 !! group for use in their coupled ensemble framework. These interfaces should be suitable for
550 !! interfacing MOM6 to other data assimilation packages as well.
551 
552 end module mom_oda_driver_mod
A structure with a pointer to a domain2d, to allow for the creation of arrays of pointers.
Interfaces for MOM6 ensembles and data assimilation.
A K-d tree tpe.
Definition: kdtree.f90:9
This module contains the main regridding routines.
Definition: MOM_ALE.F90:9
A set of dummy interfaces for compiling the MOM6 DA driver code.
Initializes fixed aspects of the model, such as horizontal grid metrics, topography and Coriolis...
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Reads the only Fortran name list needed to boot-strap the model.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
Example type for ocean ensemble DA state.
Generates vertical grids as part of the ALE algorithm.
Provides column-wise vertical remapping functions.
Initializes horizontal grid.
Describes the horizontal ocean grid with only dynamic memory arrays.
This module contains I/O framework code.
Definition: MOM_io.F90:2
The MOM6 facility to parse input files for runtime parameters.
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Container for horizontal index ranges for data, computational and global domains. ...
Container for remapping parameters.
Control structure that contains a transpose of the ocean state across ensemble members.
Describes various unit conversion factors.
Copy one MOM_domain_type into another.
Definition: MOM_domains.F90:99
Regridding control structure.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Dummy interfaces for writing ODA data.
Routines for error handling and I/O management.
Describes the vertical ocean grid, including unit conversion factors.
The MOM_domain_type contains information about the domain decompositoin.
Initializes fixed aspects of the related to its vertical coordinate.
An overloaded interface to read various types of parameters.
A null version of K-d tree from geoKdTree.
Definition: kdtree.f90:2
Container for paths and parameter file names.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Module with routines for copying information from a shared dynamic horizontal grid to an ocean-specif...
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Handy functions for manipulating strings.
Provides a transparent vertical ocean grid type and supporting routines.
ALE control structure.
Definition: MOM_ALE.F90:63
Provides transparent structures with groups of MOM6 variables and supporting routines.
Read a data field from a file.
Definition: MOM_io.F90:74
Dummy aata structures and methods for ocean data assimilation.
An overloaded interface to read and log the values of various types of parameters.