Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_setup_structures / da_write_kma_increments.inc
blob0e65de4de01c1fb6e1e0d5ae4e63002d56adf325
1 subroutine da_write_kma_increments(grid)
3    !---------------------------------------------------------------------------
4    ! Purpose: Gathers KMA analysis increments and writes 
5    !           on "anl_inc_unit" unit 
6    !---------------------------------------------------------------------------
8    implicit none
10    type (domain), intent(in) :: grid
12    ! Arrays for write out increments:
13    integer                                     :: ix, jy, kz
15 #ifdef DM_PARALLEL
16    real, dimension(1:grid%xb%mix,1:grid%xb%mjy)          :: gbuf_2d
17    real, dimension(1:grid%xb%mix,1:grid%xb%mjy,1:grid%xb%mkz) :: gbuf
18    real, dimension(:,:)  , allocatable :: psfc_g
19    real, dimension(:,:,:), allocatable :: u_g, v_g, t_g, q_g, p_g
20 #endif
22    integer                                     :: i, j, k,anl_inc_unit
24    if (trace_use) call da_trace_entry("da_write_kma_increments")
26    ! Dimension of the domain:
27    ix = grid%xb%mix
28    jy = grid%xb%mjy
29    kz = grid%xb%mkz
31 #ifdef DM_PARALLEL
33    ! 3-d and 2-d increments:
35    allocate (psfc_g (1:ix,1:jy))
36    allocate (   u_g (1:ix,1:jy,1:kz))
37    allocate (   v_g (1:ix,1:jy,1:kz))
38    allocate (   t_g (1:ix,1:jy,1:kz))
39    allocate (   q_g (1:ix,1:jy,1:kz))
40    allocate (   p_g (1:ix,1:jy,1:kz))
42    call da_patch_to_global(grid, grid%xa%psfc, gbuf_2d) 
43    if (rootproc) then 
44       psfc_g(1:ix,1:jy) = gbuf_2d(1:ix,1:jy) 
45    end if 
47    call da_patch_to_global(grid, grid%xa%u, gbuf) 
48    if (rootproc) then 
49       u_g(1:ix,1:jy,1:kz) = gbuf(1:ix,1:jy,1:kz) 
50    end if 
52    call da_patch_to_global(grid, grid%xa%v, gbuf) 
53    if (rootproc) then 
54       v_g(1:ix,1:jy,1:kz) = gbuf(1:ix,1:jy,1:kz) 
55    end if 
57    call da_patch_to_global(grid, grid%xa%t, gbuf) 
58    if (rootproc) then 
59       t_g(1:ix,1:jy,1:kz) = gbuf(1:ix,1:jy,1:kz) 
60    end if 
62    call da_patch_to_global(grid, grid%xa%q, gbuf) 
63    if (rootproc) then 
64       q_g(1:ix,1:jy,1:kz) = gbuf(1:ix,1:jy,1:kz) 
65    end if 
67    call da_patch_to_global(grid, grid%xa%p, gbuf) 
68    if (rootproc) then 
69       p_g(1:ix,1:jy,1:kz) = gbuf(1:ix,1:jy,1:kz) 
70    end if 
71 #endif
73    if (rootproc) then
74       ! 3d- and 2d-increments:
76       call da_get_unit(anl_inc_unit)
77       open(unit=anl_inc_unit,file="analysis_increments_kma",status="replace", &
78          form="unformatted")
79 #ifdef DM_PARALLEL
80       write(anl_inc_unit) ((psfc_g(i,j),i=ids,ide),j=jds,jde)
81       write(anl_inc_unit) (((u_g(i,j,k),i=ids,ide),j=ids,jde),k=kds,kde)
82       write(anl_inc_unit) (((v_g(i,j,k),i=ids,ide),j=ids,jde),k=kds,kde)
83       write(anl_inc_unit) (((t_g(i,j,k),i=ids,ide),j=ids,jde),k=kds,kde)
84       write(anl_inc_unit) (((q_g(i,j,k),i=ids,ide),j=ids,jde),k=kds,kde)
85       write(anl_inc_unit) (((p_g(i,j,k),i=ids,ide),j=ids,jde),k=kds,kde)
86 #else
87       write(anl_inc_unit) ((grid%xa%psfc(i,j),i=ids,ide),j=jds,jde)
88       write(anl_inc_unit) (((grid%xa%u(i,j,k),i=ids,ide),j=jds,jde),k=kds,kde)
89       write(anl_inc_unit) (((grid%xa%v(i,j,k),i=ids,ide),j=jds,jde),k=kds,kde)
90       write(anl_inc_unit) (((grid%xa%t(i,j,k),i=ids,ide),j=jds,jde),k=kds,kde)
91       write(anl_inc_unit) (((grid%xa%q(i,j,k),i=ids,ide),j=jds,jde),k=kds,kde)
92       write(anl_inc_unit) (((grid%xa%p(i,j,k),i=ids,ide),j=jds,jde),k=kds,kde)
93 #endif
94       close(anl_inc_unit)
95       call da_free_unit(anl_inc_unit)
96    end if
98    if (trace_use) call da_trace_exit("da_write_kma_increments")
100 end subroutine da_write_kma_increments