wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_par_util / da_cv_to_global.inc
blob82ca3c5e227a083535c47454890df85064dc8b7d
1 subroutine da_cv_to_global(cv_size, cv_size_global, x, grid, mzs, xg)
4    !-----------------------------------------------------------------------
5    ! Purpose: Gathers local cv-array x into domain cv-array xg(:).  
6    ! Global cv-array xg will only be valid on the "monitor" task.  
7    !
8    ! Must be called by all MPI tasks.  
9    !-----------------------------------------------------------------------
11    implicit none
13    integer,          intent(in)    :: cv_size        ! Size of local cv-array
14    integer,          intent(in)    :: cv_size_global ! Size of domain cv-array
15    real,             intent(in)    :: x(1:cv_size)   ! local cv-array
16    type(domain),     intent(in)    :: grid
17    integer,          intent(in)    :: mzs(:)  ! mz for each variable
18                                          ! (to identify empty or 2D arrays)
19    real,             intent(inout) :: xg(1:cv_size_global) ! global cv-array
21 #ifdef DM_PARALLEL
23    ! Local declarations
24    type (vp_type) :: vv_x    ! Grdipt/EOF cv_array (local)
25    type (vp_type) :: vv_xg   ! Grdipt/EOF cv_array (global)
26    type (xpose_type) :: xpg  ! global dimensions
27    integer :: ids, ide, jds, jde, kds, &
28               ims, ime, jms, jme, kms, &
29               ips, ipe, jps, jpe, kps 
31    integer :: n
33    if (trace_use) call da_trace_entry("da_cv_to_global")      
35    !
36    ! Gather to mimic serial summation order.  
37    !
39    ! k?e varies with variable v1 - v5
40    ids = ids; ide = ide; jds = jds; jde = jde; kds = kds
41    ims = ims; ime = ime; jms = jms; jme = jme; kms = grid%xp%kms
42    ips = grid%xp%ips; ipe = grid%xp%ipe; jps = grid%xp%jps; jpe = grid%xp%jpe; kps = grid%xp%kps
44    ! TOdo:  encapsulate this crap!  
45    ! allocate vv_x
46    allocate(vv_x%v1(ims:ime,jms:jme,mzs(1)))
47    allocate(vv_x%v2(ims:ime,jms:jme,mzs(2)))
48    allocate(vv_x%v3(ims:ime,jms:jme,mzs(3)))
49    allocate(vv_x%v4(ims:ime,jms:jme,mzs(4)))
50    allocate(vv_x%v5(ims:ime,jms:jme,mzs(5)))
51    allocate(vv_x%alpha(ims:ime,jms:jme,kms:kme,mzs(7)))
53    call da_cv_to_vv (cv_size, x, mzs, vv_x)
55    if (rootproc) then
56       ! allocate vv_xg
57       allocate(vv_xg%v1(ids:ide,jds:jde,mzs(1)))
58       allocate(vv_xg%v2(ids:ide,jds:jde,mzs(2)))
59       allocate(vv_xg%v3(ids:ide,jds:jde,mzs(3)))
60       allocate(vv_xg%v4(ids:ide,jds:jde,mzs(4)))
61       allocate(vv_xg%v5(ids:ide,jds:jde,mzs(5)))
62       allocate(vv_xg%alpha(ids:ide,jds:jde,kds:kde,mzs(7)))
63    else
64       ! Allocate dummy array for non-monitor process to keep Fortran
65       ! CICO happy...
66       allocate(vv_xg%v1(1,1,1))
67       allocate(vv_xg%v2(1,1,1))
68       allocate(vv_xg%v3(1,1,1))
69       allocate(vv_xg%v4(1,1,1))
70       allocate(vv_xg%v5(1,1,1))
71       allocate(vv_xg%alpha(1,1,1,1))
72    end if
74    ! TOdo:  encapsulate this crap!  
75    ! gather to global data structures
76    ! it is possible to gather straight into a globally-sized cv-array
77    ! TOdo:  add this optimization later
78    call da_patch_to_global(grid, vv_x%v1,    vv_xg%v1,    mzs(1))
79    call da_patch_to_global(grid, vv_x%v2,    vv_xg%v2,    mzs(2))
80    call da_patch_to_global(grid, vv_x%v3,    vv_xg%v3,    mzs(3))
81    call da_patch_to_global(grid, vv_x%v4,    vv_xg%v4,    mzs(4))
82    call da_patch_to_global(grid, vv_x%v5,    vv_xg%v5,    mzs(5))
84    if ( mzs(7) > 0 ) then
85       do n = 1, mzs(7) ! Ensemble size
86          call da_patch_to_global(grid, vv_x%alpha(:,:,:,n), vv_xg%alpha(:,:,:,n), mzs(6))
87       end do
88    end if
90    ! deallocate vv_x
91    deallocate (vv_x%v1, vv_x%v2, vv_x%v3, vv_x%v4, vv_x%v5, vv_x%alpha)
93    if (rootproc) then
94       ! finally, collapse data back into a globally-sized cv-array
95       xpg%ids = ids; xpg%ide = ide
96       xpg%ims = ids; xpg%ime = ide
97       xpg%its = ids; xpg%ite = ide
98       xpg%jds = jds; xpg%jde = jde
99       xpg%jms = jds; xpg%jme = jde
100       xpg%jts = jds; xpg%jte = jde
101       xpg%kds = kds; xpg%kde = kde
102       xpg%kms = kds; xpg%kme = kde
103       xpg%kts = kds; xpg%kte = kde
104       call da_vv_to_cv(vv_xg, xpg, mzs, cv_size_global, xg)
105    end if
107    ! deallocate vv_xg
108    deallocate(vv_xg%v1, vv_xg%v2, vv_xg%v3, vv_xg%v4, vv_xg%v5, vv_xg%alpha)
110    if (trace_use) call da_trace_exit("da_cv_to_global")    
112 #endif
114 end subroutine da_cv_to_global