da_transform_vtoy.inc
References to this file elsewhere.
1 subroutine da_transform_vtoy(cv_size, be, ep, cv, iv, vp, vv, xa, xb, xbx, &
2 xp, 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 (ob_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 (x_type) , intent(inout) :: xa ! gridded analy. incs. (local)
19 type (xb_type), intent(in) :: xb ! first guess (local).
20 type (xbx_type),intent(in) :: xbx ! For header & non-grid arrays.
21 type (xpose_type), intent(inout) :: xp ! Domain decomposition vars.
22 type (y_type) , intent(inout) :: y ! y = H(x_inc).
23
24 type(domain), intent(inout) :: grid
25 type(grid_config_rec_type), intent(inout) :: config_flags
26
27 integer :: nobwin
28 #ifdef DM_PARALLEL
29 integer :: wrf_done_unit
30 #endif
31
32 character(len=4) :: filnam
33
34 call da_trace_entry("da_transform_vtoy")
35
36 call da_transform_vtox(cv_size, xb, xbx, be, ep, cv, vv, vp, xp, xa)
37
38 if (var4d) then
39 call da_transfer_xatowrftl(grid, config_flags, 'tl01')
40
41 call da_trace("da_transform_vtoy","Starting da_run_wrfplus_tl.ksh")
42 #ifdef DM_PARALLEL
43 if (var4d_coupling == var4d_coupling_disk_simul) then
44 call da_system("da_run_wrfplus_tl.ksh pre")
45
46 if (rootproc) then
47 call da_system("rm -rf wrf_done")
48 call da_system("touch wrf_go_ahead")
49 call da_get_unit(wrf_done_unit)
50 do while (.true.)
51 open(wrf_done_unit,file="wrf_done",status="old",err=303)
52 close(wrf_done_unit)
53 exit
54 303 continue
55 call da_system("sleep 1")
56 end do
57 call da_free_unit(wrf_done_unit)
58 end if
59 ! Wait until PE 0 agrees that TL finished
60 CALL wrf_get_dm_communicator ( comm )
61 call mpi_barrier(comm, ierr)
62 call da_system("da_run_wrfplus_tl.ksh post")
63 end if
64 #else
65 call da_system("da_run_wrfplus_tl.ksh")
66 #endif
67 call da_trace("da_transform_vtoy","Finished da_run_wrfplus_tl.ksh")
68
69 do nobwin=1, num_fgat_time
70 iv%current_ob_time = nobwin
71 write(filnam,'(a2,i2.2)') 'tl',nobwin
72 call da_transfer_wrftltoxa(grid, config_flags, filnam)
73 call da_transform_xtoy(xb, iv, xa, xp, y)
74 end do
75 else ! not var4d
76 call da_transform_xtoy(xb, iv, xa, xp, y)
77 end if ! var4d
78
79 call da_trace_exit("da_transform_vtoy")
80
81 end subroutine da_transform_vtoy
82
83