Update version info for release v4.6.1 (#2122)
[WRF.git] / var / da / da_par_util / da_patch_to_global_3d.inc
blobfe7b05733e260bb41baef769ac55b6f5badcaee2
1 subroutine da_patch_to_global_3d(grid, vlocal, vglobal, mz)
3    !----------------------------------------------------------------------
4    ! Purpose: Gathers local 3D array vlocal into global array vglobal.  
5    ! Assumes that "k" is not decomposed.  End indices in the "k" dimension 
6    ! are inferred from mz, which can be less than kde.  
7    !
8    ! Must be called by all MPI tasks.  
9    !----------------------------------------------------------------------
11    implicit none
13    type(domain),      intent(in)  :: grid
14    real,              intent(in)  :: vlocal(ims:ime,jms:jme,kms:kme)
15    real,              intent(out) :: vglobal(ids:ide,jds:jde,kds:kde)
16    integer, optional, intent(in)  :: mz
18 #ifdef DM_PARALLEL
20    integer :: local_kde
21    integer :: local_kme
22    integer :: local_kpe
23 #endif
25    if (trace_use_frequent) call da_trace_entry("da_patch_to_global_3d")
27 #ifdef DM_PARALLEL
29    if (present(mz)) then
30       local_kde = kds + mz - 1
31       local_kme = local_kde
32       local_kpe = local_kde
33    else
34       local_kde = kde
35       local_kme = kme
36       local_kpe = kpe
37    end if
39    if (local_kde > 0) then
40       call wrf_patch_to_global_real (vlocal, vglobal, grid%xp%domdesc, &
41          trim(grid_stagger), trim(grid_ordering), &
42          ids, ide, jds, jde, kds, local_kde,  &
43          ims, ime, jms, jme, kms, local_kme,  &
44          ips, ipe, jps, jpe, kps, local_kpe)
45    end if
46 #else
47    vglobal(ids:ide,jds:jde,kds:kde) = vlocal(ids:ide,jds:jde,kds:kde)
48 #endif
50    if (trace_use_frequent) call da_trace_exit("da_patch_to_global_3d")
52 end subroutine da_patch_to_global_3d