da_transform_vtoy_adj.inc

References to this file elsewhere.
1 subroutine da_transform_vtoy_adj(iter,cv_size, be, ep, cv, iv, vp, vv, xa, xb, &
2                                   xbx, xp, y, &
3                                   grid, config_flags, cg_jcdf)
4 
5    !-------------------------------------------------------------------------
6    ! Purpose:  Does Adjoint of control variable (V) transform to Obs-space(Y)
7    !-------------------------------------------------------------------------
8 
9    implicit none
10 
11    integer, intent(in)              :: iter
12    integer, intent(in)              :: cv_size ! Size of cv array.
13    type (be_type), intent(in)       :: be     ! background error structure.
14    type (ep_type), intent(in)       :: ep     ! ensemble perturbation structure.
15    real, intent(out)                :: cv(1:cv_size) ! control variables.
16    type (ob_type), intent(inout)    :: iv     ! innovation vector (o-b).
17    type (vp_type), intent(inout)    :: vp     ! Grdipt/level CV.
18    type (vp_type), intent(inout)    :: vv     ! Grdipt/EOF CV.
19    type (x_type) , intent(inout)    :: xa     ! gridded analy. incs. (local)
20    type (xb_type), intent(in)       :: xb     ! first guess (local).
21    type (xbx_type),intent(in)       :: xbx    ! For header & non-grid arrays.
22    type (xpose_type), intent(inout) :: xp     ! Domain decomposition vars.
23    type (y_type) , intent(inout)    :: y      ! y = H(x_inc).
24 
25    type(domain), intent(inout) :: grid
26    type(grid_config_rec_type), intent(inout) :: config_flags
27    integer, intent(in) :: cg_jcdf
28 
29    integer :: nobwin,ndynopt
30 #ifdef DM_PARALLEL
31    integer :: wrf_done_unit
32 #endif
33 
34    character(len=4) :: filnam
35 
36    call da_trace_entry("da_transform_vtoy_adj")
37 
38    if (var4d) then
39       do nobwin=num_fgat_time,1,-1
40          iv%current_ob_time = nobwin
41          call da_zero_x(xa)
42          call da_transform_xtoy_adj(xb, iv, xp, y, xa)
43          write(unit=filnam,fmt='(a2,i2.2)') 'af',nobwin
44 
45          call da_transfer_wrftltoxa_adj(grid, config_flags, filnam)
46       end do
47 
48       ndynopt      = grid%dyn_opt
49       grid%dyn_opt = DYN_EM_TL
50       call nl_set_dyn_opt (1 , DYN_EM_TL)
51 
52       if (jcdfi_use .AND. iter > 0 .AND. cg_jcdf == 1) then
53          call da_med_initialdata_input(grid , config_flags, 'tldf')
54       else
55          grid%em_g_u_2 = 0.0
56          grid%em_g_v_2 = 0.0
57          grid%em_g_w_2 = 0.0
58          grid%em_g_t_2 = 0.0
59          grid%em_g_ph_2 = 0.0
60          grid%em_g_mu_2 = 0.0
61          grid%g_moist = 0.0
62       end if
63 
64       call med_hist_out(grid , 3 , config_flags)
65 
66       grid%dyn_opt = ndynopt
67       call nl_set_dyn_opt (1 , DYN_EM)
68 
69       call da_trace("da_transform_vtoy_adj","Starting da_run_wrfplus_ad.ksh")
70 
71 #ifdef DM_PARALLEL
72       if (var4d_coupling == var4d_coupling_disk_simul) then
73          call da_system("da_run_wrfplus_ad.ksh pre")
74 
75          if (rootproc) then
76             call da_system("rm -rf wrf_done")
77             call da_system("touch wrf_go_ahead")
78             call da_get_unit(wrf_done_unit)
79             do while (.true.)
80                open(wrf_done_unit,file="wrf_done",status="old",err=303)
81                close(wrf_done_unit)
82                exit
83 303            continue
84                call da_system("sleep 1")
85             end do
86             call da_free_unit(wrf_done_unit)
87          end if
88          ! Wait until PE 0 agrees that AD finished
89          CALL wrf_get_dm_communicator ( comm )
90          call mpi_barrier(comm, ierr)
91          call da_system("da_run_wrfplus_ad.ksh post")
92       end if
93 #else
94       call da_system("da_run_wrfplus_ad.ksh")
95 #endif
96       call da_trace("da_transform_vtoy_adj","Finished da_run_wrfplus_ad.ksh")
97 
98       call da_transfer_xatowrftl_adj(grid, config_flags, 'gr01')
99    else  ! var4d
100       call da_zero_x(xa)
101       call da_transform_xtoy_adj(xb, iv, xp, y, xa)
102    end if ! var4d
103 
104    cv = 0.0
105    call da_transform_vtox_adj(cv_size, xb, xbx, be, ep, xa, xp, vp, vv, cv)
106 
107    call da_trace_exit("da_transform_vtoy_adj")
108 
109 end subroutine da_transform_vtoy_adj
110 
111