1 subroutine da_check(grid, config_flags, cv_size, xbx, be, ep, iv, vv, vp, y &
8 use module_domain, only : xchem_type
9 use da_control, only : test_transformsch
12 !-----------------------------------------------------------------------
14 !-----------------------------------------------------------------------
18 type(domain), intent(inout) :: grid
19 type(grid_config_rec_type), intent(inout) :: config_flags
21 integer, intent(in) :: cv_size ! Size of cv array.
22 type (xbx_type), intent(inout) :: xbx ! Header & non-gridded vars.
23 type (be_type), intent(in) :: be ! background error structure.
24 type (ep_type), intent(in) :: ep ! Ensemble perturbation structure.
25 type (iv_type), intent(inout) :: iv ! ob. increment vector.
26 type (vp_type), intent(inout) :: vv ! Grdipt/EOF CV.
27 type (vp_type), intent(inout) :: vp ! Grdipt/level CV.
28 type (y_type), intent(inout) :: y ! y = h (xa)
31 type(xchem_type), optional, intent(out) :: vchem ! CHEM_IC CV
35 real :: cvtest(1:cv_size) ! background error structure.
36 real :: field(its:ite,jts:jte) ! Field for spectral transform test.
38 call da_trace_entry("da_check")
40 !----------------------------------------------------------------------------
41 ! [1] Set up test data:
42 !----------------------------------------------------------------------------
44 ! Initialize cv values with random data:
45 call random_number(cvtest(:))
46 cvtest(:) = cvtest(:) - 0.5
48 ! vv arrays initialized already.
49 ! vp arrays initialized already.
51 !----------------------------------------------------------------------------
52 ! [2] Perform vtox adjoint tests:
53 !----------------------------------------------------------------------------
55 call da_message((/"Performing vtox adjoint tests"/))
57 ! v_to_vv adjoint test:
59 call da_check_cvtovv_adjoint(grid, cv_size, xbx, be, cvtest, vv)
62 if (test_transformsch) then
63 call da_check_cvtovv_adjoint_chem(grid, cv_size, xbx, be, cvtest, vv, vchem=grid%vchem)
67 !-------------------------------------------------------------------------
68 ! vv_to_vp adjoint test:
69 !-------------------------------------------------------------------------
71 call da_check_vvtovp_adjoint(grid, be % ne, grid%xb, be, vv, vp)
73 !-------------------------------------------------------------------------
75 !-------------------------------------------------------------------------
77 call da_check_vptox_adjoint(grid, be % ne, be, ep, vp, cv_size)
79 !-------------------------------------------------------------------------
80 ! vtox adjoint test: <x,x> = <v_adj,v>
81 !-------------------------------------------------------------------------
83 call da_check_vtox_adjoint(grid, cv_size, xbx, be, ep, cvtest, vv, vp)
86 if (test_transformsch) then
87 call da_check_vchemtox_adjoint(grid, grid%vchem, be, grid%vchem2)
88 call da_check_vtox_adjoint_chem(grid, cv_size, xbx, be, ep, cvtest, vv, vp)
92 !----------------------------------------------------------------------------
93 ! [2] Perform xtoy adjoint tests:
94 !----------------------------------------------------------------------------
96 call da_message((/"Performing xtoy adjoint tests"/))
98 call da_allocate_y(iv, y)
99 call da_zero_x(grid%xa)
101 call da_setup_testfield(grid)
104 ! Make cv_array random.
106 ! call random_number(cvtest(1:cv_size))
107 ! cvtest(1:cv_size) = cvtest(1:cv_size) - 0.5
109 ! call da_transform_vtox(grid, cv_size, xbx, be, ep, cvtest, vv, vp)
111 call da_check_xtoy_adjoint(cv_size, cvtest, xbx, be, grid, config_flags, iv, y)
113 call da_deallocate_y(y)
115 !----------------------------------------------------------------------------
116 ! [3] Perform dynamical constraint test:
117 !----------------------------------------------------------------------------
118 call da_message((/"Performing dynamical constraint adjoint tests"/))
119 call da_zero_x(grid%xa)
120 call da_setup_testfield(grid)
121 call da_check_dynamics_adjoint(cv_size, cvtest, xbx, be, grid, config_flags)
123 !----------------------------------------------------------------------------
124 ! [4] Perform spectral test:
125 !----------------------------------------------------------------------------
129 call da_message((/"Performing spectral tests"/))
131 call random_number(field(:,:))
132 field(:,:) = field(:,:) - 0.5
134 sizec = (be % max_wave+1) * (be % max_wave+2)/2
135 call da_test_spectral(be % max_wave, sizec, xbx, field)
139 call da_trace_exit("da_check")
141 end subroutine da_check