da_setup_firstguess.inc

References to this file elsewhere.
1 subroutine da_setup_firstguess(xbx, grid)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Allocate and read in components of first guess state.
5    !---------------------------------------------------------------------------
6 
7    implicit none
8 
9    type (xbx_type),intent(out)  :: xbx   ! Header & non-gridded vars.
10 
11    type(domain),intent(inout)   :: grid
12 
13    integer :: is, ie, js, je
14    real    :: ddx , ddy    
15 
16    if (trace_use) call da_trace_entry("da_setup_firstguess")
17 
18    is = grid%xp % its
19    ie = grid%xp % ite
20    js = grid%xp % jts
21    je = grid%xp % jte
22 
23    ! Calculate sin and cosine values used in da_get_avpoles
24 
25    if (global) then
26       if (grid%xp%jts == grid%xp%jds) then
27 
28          allocate(cos_xls(grid%xp%its:grid%xp%ite))
29          allocate(sin_xls(grid%xp%its:grid%xp%ite))
30          cos_xls(grid%xp%its:grid%xp%ite) = & 
31             cos(deg_to_rad*grid%xlong(grid%xp%its:grid%xp%ite,grid%xp%jts))
32          sin_xls(grid%xp%its:grid%xp%ite) = &
33             sin(deg_to_rad*grid%xlong(grid%xp%its:grid%xp%ite,grid%xp%jts))
34       end if
35 
36       if (grid%xp%jte == grid%xp%jde) then 
37          allocate(cos_xle(grid%xp%its:grid%xp%ite))
38          allocate(sin_xle(grid%xp%its:grid%xp%ite))
39          cos_xle(grid%xp%its:grid%xp%ite) = &
40             cos(deg_to_rad*grid%xlong(grid%xp%its:grid%xp%ite,grid%xp%jte))
41          sin_xle(grid%xp%its:grid%xp%ite) = &
42             sin(deg_to_rad*grid%xlong(grid%xp%its:grid%xp%ite,grid%xp%jte))
43       end if
44    end if
45 
46    !---------------------------------------------------------------------------      
47    ! [1.0] Setup and read in fields from first guess:
48    !---------------------------------------------------------------------------      
49 
50    if (fg_format == fg_format_wrf) then
51 
52       ! First guess is a WRF C-grid format file:
53 
54       call da_setup_firstguess_wrf(xbx, grid)
55    else if (fg_format == fg_format_kma_global) then
56       ! First guess is an KMA format file:
57       call da_setup_firstguess_kma(xbx, grid)
58    end if
59 
60    !---------------------------------------------------------------------------
61    ! Exchange halo region for XB arrays.
62    !---------------------------------------------------------------------------
63 
64    if (fg_format == fg_format_wrf) then
65       ! Calculate multiplicative constants for PsiChi_TO_UV 
66       grid%xb%coefx(is:ie,js:je) = 0.5 * grid%xb%map_factor(is:ie,js:je)/grid%xb%ds
67       grid%xb%coefy(is:ie,js:je) = grid%xb%coefx(is:ie,js:je)
68       grid%xb%coefz(is:ie,js:je) = 0.5 / (grid%xb%map_factor(is:ie,js:je)*grid%xb%ds)
69    else if (fg_format == fg_format_kma_global) then
70       ! Calculate multiplicative constants for PsiChi_TO_UV 
71       ddx =  earth_radius*1000 * 2.0 * pi / (grid%xb%ide-grid%xb%ids+1)
72       ddy =  earth_radius*1000       * pi / (grid%xb%jde-grid%xb%jds)
73       grid%xb% coefx(is:ie,js:je) = 0.5 / (ddx * cos(grid%xlat(is:ie,js:je)*pi/180.))
74       grid%xb% coefy(is:ie,js:je) = 0.5 /  ddy
75    else
76       write(unit=message(1),fmt='(A,I5)') &
77          "Wrong choice for fg_format = ",fg_format
78       call da_error(__FILE__,__LINE__,message(1:1))
79    end if
80 
81 #ifdef DM_PARALLEL
82 #include "HALO_INIT.inc"
83 #endif
84    periodic_x = grid%periodic_x
85 
86    if (global) then     
87       ! Set East-West boundary for Xb-array 
88       call da_set_boundary_xb(grid)
89    end if
90 
91    !---------------------------------------------------------------------------      
92    ! [2.0] Setup grid-dependent constants used:
93    !---------------------------------------------------------------------------
94 
95    ! [2.1] Set up fast Fourier & Legendre transform constants:
96 
97    call da_setup_runconstants(grid, xbx)
98 
99    if (trace_use) call da_trace_exit("da_setup_firstguess")
100 
101 end subroutine da_setup_firstguess
102 
103