da_setup_testfield.inc

References to this file elsewhere.
1 subroutine da_setup_testfield(xb, xa, xp)
2 
3    !----------------------------------------------------------------------------
4    ! Purpose: produce test increment field based on xb field.
5    !
6    ! Method:  pass through x=uv transfom to ensure satisfies boundary conditions
7    !----------------------------------------------------------------------------
8 
9    implicit none
10 
11    type (xb_type),    intent(in)  :: xb  ! first guess (local).
12    type (x_type) ,    intent(inout) :: xa  ! gridded analy. incs. (local)
13    type (xpose_type), intent(in)  :: xp  ! Transpose variables.
14 
15    integer :: i, j
16 
17    !-------------------------------------------------------------------------
18    ! [1.0]: initialise:
19    !-------------------------------------------------------------------------
20 
21    write(unit=stdout, fmt='(/a/)') &
22       'Starting da_setup_testfield ...'
23 
24    !-------------------------------------------------------------------------
25    ! [2.0]: set up test increment field structure:
26    !-------------------------------------------------------------------------
27 
28    ! [2.1] test wind, temperature, pressure, and humidity.
29 
30    call da_set_tst_trnsf_fld(xp, xa%u, xb%u, typical_u_rms)
31    call da_set_tst_trnsf_fld(xp, xa%v, xb%v, typical_v_rms)
32    call da_set_tst_trnsf_fld(xp, xa%w, xb%w, typical_w_rms)
33    call da_set_tst_trnsf_fld(xp, xa%t, xb%t, typical_t_rms)
34    call da_set_tst_trnsf_fld(xp, xa%p, xb%p, typical_p_rms)
35    call da_set_tst_trnsf_fld(xp, xa%q, xb%q, typical_q_rms)
36    call da_set_tst_trnsf_fld(xp, xa%qcw, xb%qcw, typical_qcw_rms)
37    call da_set_tst_trnsf_fld(xp, xa%qrn, xb%qrn, typical_qrn_rms)
38 
39    ! [2.5] get test density increment from linearised ideal gas law:
40 
41    call da_pt_to_rho_lin(xb, xa, xp)
42 
43    xa%psfc(xp%its:xp%ite, xp%jts:xp%jte) = &
44    xa%p   (xp%its:xp%ite, xp%jts:xp%jte, xp%kts)
45 
46    if (print_detail_testing) then
47       write(unit=stdout, fmt='(2a,4x,a,i8)') &
48          'file:', __FILE__, 'line:', __LINE__
49 
50       write(unit=stdout, fmt=*) 'xp%its, xp%ite, xp%jts, xp%jte) =', &
51          xp%its, xp%ite, xp%jts, xp%jte
52    
53       do j=xp%jts, xp%jte
54          do i=xp%its, xp%ite
55             if (i == j) then
56                write(unit=stdout, fmt='(2(a,i4),a,f14.6)') &
57                   'xa%psfc(', i, ',', j, ') =', xa%psfc(i, j)
58             end if
59          end do
60       end do
61    end if
62 
63    ! Exchange XA halo region.
64    call wrf_dm_halo(xp%domdesc,xp%comms,xp%halo_id4)
65 
66    if (sfc_assi_options == 2) then
67       ! Exchange XA (surface variable) halo region.
68       call wrf_dm_halo(xp%domdesc,xp%comms,xp%halo_id6)
69    end if
70 
71    if (use_ssmt1obs .or. use_ssmt2obs .or. Use_GpspwObs .or. &
72         Use_SsmiTbObs .or. Use_SsmiRetrievalObs) then
73 
74       ! Now do something for PW
75       call da_transform_xtotpw(xa, xb)
76 
77       ! GPS Refractivity:
78       if (use_GpsrefObs) &
79          call da_transform_xtogpsref_lin(xa, xb, xp)
80 
81       if (use_ssmt1obs .or. use_ssmt2obs .or. &
82            Use_SsmiTbObs .or. Use_SsmiRetrievalObs) then
83          call da_transform_xtoseasfcwind_lin(xa, xb)
84       end if
85 
86       ! Exchange XA halo region.
87       call wrf_dm_halo(xp%domdesc,xp%comms,xp%halo_id8)
88    end if
89 
90    write(unit=stdout, fmt='(/a/)') &
91       'End of da_setup_testfield.'
92 
93 end subroutine da_setup_testfield
94 
95