da_transform_vtoy.inc
References to this file elsewhere.
1 subroutine da_transform_vtoy(cv_size, be, ep, cv, iv, vp, vv, xbx, &
2 y, &
3 grid, config_flags)
4
5 !----------------------------------------------------------------------
6 ! Purpose: Does transform of control variable (V) to Obs-space (Y)
7 !----------------------------------------------------------------------
8
9 implicit none
10
11 integer, intent(in) :: cv_size ! Size of cv array.
12 type(be_type), intent(in) :: be ! background error structure.
13 type(ep_type), intent(in) :: ep ! Ensemble perturbation structure.
14 real, intent(in) :: cv(1:cv_size) ! control variables.
15 type(iv_type), intent(inout) :: iv ! innovation vector (o-b).
16 type(vp_type), intent(inout) :: vp ! Grdipt/level CV.
17 type(vp_type), intent(inout) :: vv ! Grdipt/EOF CV.
18 type(xbx_type), intent(in) :: xbx ! For header & non-grid arrays.
19 type(y_type), intent(inout) :: y ! y = H(x_inc).
20 type(domain), intent(inout) :: grid
21 type(grid_config_rec_type), intent(inout) :: config_flags
22 type (xbx_type) :: xbx_tmp
23
24 integer :: nobwin
25 #ifdef DM_PARALLEL
26 integer :: wrf_done_unit
27 #endif
28
29 character(len=4) :: filnam
30
31 if (trace_use) call da_trace_entry("da_transform_vtoy")
32
33 call da_transform_vtox(grid, cv_size, xbx, be, ep, cv, vv, vp)
34
35 if (var4d) then
36 call da_transfer_xatowrftl(grid, config_flags, 'tl01')
37
38 if (trace_use) call da_trace("da_transform_vtoy","Starting da_run_wrfplus_tl.ksh")
39 #ifdef DM_PARALLEL
40 if (var4d_coupling == var4d_coupling_disk_simul) then
41
42 if (rootproc) then
43 call da_system("da_run_wrfplus_tl.ksh pre")
44 call da_system("rm -rf wrf_done")
45 call da_system("touch wrf_go_ahead")
46 call da_get_unit(wrf_done_unit)
47 do while (.true.)
48 open(wrf_done_unit,file="wrf_done",status="old",err=303)
49 close(wrf_done_unit)
50 exit
51 303 continue
52 call da_system("sleep 1")
53 end do
54 call da_free_unit(wrf_done_unit)
55 call da_system("da_run_wrfplus_tl.ksh post")
56 end if
57 ! Wait until PE 0 agrees that TL finished
58 call wrf_get_dm_communicator ( comm )
59 call mpi_barrier(comm, ierr)
60 end if
61 #else
62 call da_system("da_run_wrfplus_tl.ksh")
63 #endif
64 if (trace_use) call da_trace("da_transform_vtoy","Finished da_run_wrfplus_tl.ksh")
65
66 do nobwin=1, num_fgat_time
67 write(unit=filnam, fmt='(a, i2.2)') 'fg', nobwin
68 call da_med_initialdata_input( grid , config_flags, &
69 filnam)
70 call da_setup_firstguess( xbx_tmp, grid)
71
72 iv%time = nobwin
73 iv%info(:)%n1 = iv%info(:)%plocal(iv%time-1) + 1
74 iv%info(:)%n2 = iv%info(:)%plocal(iv%time)
75 write(filnam,'(a2,i2.2)') 'tl',nobwin
76 call da_transfer_wrftltoxa(grid, config_flags, filnam)
77 call da_transform_xtoy(grid, iv, y)
78 end do
79 call da_med_initialdata_input( grid , config_flags, 'fg01')
80 call da_setup_firstguess( xbx_tmp, grid)
81 else ! not var4d
82 call da_transform_xtoy(grid, iv, y)
83 end if ! var4d
84
85 if (trace_use) call da_trace_exit("da_transform_vtoy")
86
87 end subroutine da_transform_vtoy
88
89