1 | subroutine MPP_GLOBAL_FIELD_2D_( domain, local, global, flags ) |
---|
2 | type(domain2D), intent(in) :: domain |
---|
3 | MPP_TYPE_, intent(in) :: local(:,:) |
---|
4 | MPP_TYPE_, intent(out) :: global(:,:) |
---|
5 | integer, intent(in), optional :: flags |
---|
6 | MPP_TYPE_ :: local3D (size( local,1),size( local,2),1) |
---|
7 | MPP_TYPE_ :: global3D(size(global,1),size(global,2),1) |
---|
8 | #ifdef use_CRI_pointers |
---|
9 | pointer( lptr, local3D ) |
---|
10 | pointer( gptr, global3D ) |
---|
11 | lptr = LOC( local) |
---|
12 | gptr = LOC(global) |
---|
13 | call mpp_global_field( domain, local3D, global3D, flags ) |
---|
14 | #else |
---|
15 | local3D = RESHAPE( local, SHAPE(local3D) ) |
---|
16 | call mpp_global_field( domain, local3D, global3D, flags ) |
---|
17 | global = RESHAPE( global3D, SHAPE(global) ) |
---|
18 | #endif |
---|
19 | end subroutine MPP_GLOBAL_FIELD_2D_ |
---|
20 | |
---|
21 | subroutine MPP_GLOBAL_FIELD_3D_( domain, local, global, flags ) |
---|
22 | !get a global field from a local field |
---|
23 | !local field may be on compute OR data domain |
---|
24 | type(domain2D), intent(in) :: domain |
---|
25 | MPP_TYPE_, intent(in) :: local(:,:,:) |
---|
26 | MPP_TYPE_, intent(out) :: global(domain%x%global%begin:,domain%y%global%begin:,:) |
---|
27 | integer, intent(in), optional :: flags |
---|
28 | integer :: i, j, k, m, n, nd, nwords, lpos, rpos, ioff, joff |
---|
29 | logical :: xonly, yonly |
---|
30 | MPP_TYPE_ :: clocal (domain%x%compute%size* domain%y%compute%size* size(local,3)) |
---|
31 | MPP_TYPE_ :: cremote(domain%x%compute%max_size*domain%y%compute%max_size*size(local,3)) |
---|
32 | integer :: words_per_long, stackuse |
---|
33 | character(len=8) :: text |
---|
34 | #ifdef use_CRI_pointers |
---|
35 | pointer( ptr_local, clocal ) |
---|
36 | pointer( ptr_remote, cremote ) |
---|
37 | |
---|
38 | ptr_local = LOC(mpp_domains_stack) |
---|
39 | ptr_remote = LOC(mpp_domains_stack(size(clocal)+1)) |
---|
40 | #endif |
---|
41 | |
---|
42 | if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: must first call mpp_domains_init.' ) |
---|
43 | xonly = .FALSE. |
---|
44 | yonly = .FALSE. |
---|
45 | if( PRESENT(flags) )then |
---|
46 | xonly = flags.EQ.XUPDATE |
---|
47 | yonly = flags.EQ.YUPDATE |
---|
48 | if( .NOT.xonly .AND. .NOT.yonly )call mpp_error( WARNING, 'MPP_GLOBAL_FIELD: you must have flags=XUPDATE or YUPDATE.' ) |
---|
49 | end if |
---|
50 | |
---|
51 | #ifdef use_CRI_pointers |
---|
52 | words_per_long = size(transfer(local(1,1,1),mpp_domains_stack)) |
---|
53 | stackuse = (size(clocal)+size(cremote))*words_per_long |
---|
54 | if( stackuse.GT.mpp_domains_stack_size )then |
---|
55 | write( text, '(i8)' )stackuse |
---|
56 | call mpp_error( FATAL, & |
---|
57 | 'MPP_UPDATE_DOMAINS user stack overflow: call mpp_domains0_set_stack_size('//trim(text)//') from all PEs.' ) |
---|
58 | end if |
---|
59 | ptr_local = LOC(mpp_domains_stack) |
---|
60 | ptr_remote = LOC( mpp_domains_stack((size(clocal)+1)*words_per_long) ) |
---|
61 | mpp_domains_stack_hwm = max( mpp_domains_stack_hwm, stackuse ) |
---|
62 | #endif |
---|
63 | if( size(global,1).NE.domain%x%global%size .OR. size(global,2).NE.domain%y%global%size .OR. & |
---|
64 | size(local,3).NE.size(global,3) ) & |
---|
65 | call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: incoming arrays do not match domain.' ) |
---|
66 | if( size(local,1).EQ.domain%x%compute%size .AND. size(local,2).EQ.domain%y%compute%size )then |
---|
67 | !local is on compute domain |
---|
68 | ioff = -domain%x%compute%begin + 1 |
---|
69 | joff = -domain%y%compute%begin + 1 |
---|
70 | else if( size(local,1).EQ.domain%x%data%size .AND. size(local,2).EQ.domain%y%data%size )then |
---|
71 | !local is on data domain |
---|
72 | ioff = -domain%x%data%begin + 1 |
---|
73 | joff = -domain%y%data%begin + 1 |
---|
74 | else |
---|
75 | call mpp_error( FATAL, 'MPP_GLOBAL_FIELD_: incoming field array must match either compute domain or data domain.' ) |
---|
76 | end if |
---|
77 | |
---|
78 | ! make contiguous array from compute domain |
---|
79 | m = 0 |
---|
80 | do k = 1,size(local,3) |
---|
81 | do j = domain%y%compute%begin,domain%y%compute%end |
---|
82 | do i = domain%x%compute%begin,domain%x%compute%end |
---|
83 | m = m + 1 |
---|
84 | clocal(m) = local(i+ioff,j+joff,k) |
---|
85 | global(i,j,k) = clocal(m) !always fill local domain directly |
---|
86 | end do |
---|
87 | end do |
---|
88 | end do |
---|
89 | |
---|
90 | !fill off-domains (note loops begin at an offset of 1) |
---|
91 | if( xonly )then |
---|
92 | nd = size(domain%x%list) |
---|
93 | do n = 1,nd-1 |
---|
94 | lpos = mod(domain%x%pos+nd-n,nd) |
---|
95 | rpos = mod(domain%x%pos +n,nd) |
---|
96 | nwords = domain%x%list(rpos)%compute%size * domain%x%list(rpos)%compute%size * size(local,3) |
---|
97 | call mpp_transmit( clocal, size(clocal), domain%x%list(lpos)%pe, cremote, nwords, domain%x%list(rpos)%pe ) |
---|
98 | m = 0 |
---|
99 | do k = 1,size(global,3) |
---|
100 | do j = domain%y%compute%begin,domain%y%compute%end |
---|
101 | do i = domain%x%list(rpos)%compute%begin,domain%x%list(rpos)%compute%end |
---|
102 | m = m + 1 |
---|
103 | global(i,j,k) = cremote(m) |
---|
104 | end do |
---|
105 | end do |
---|
106 | end do |
---|
107 | end do |
---|
108 | call mpp_sync_self(domain%x%list(:)%pe) |
---|
109 | else if( yonly )then |
---|
110 | nd = size(domain%y%list) |
---|
111 | do n = 1,nd-1 |
---|
112 | lpos = mod(domain%y%pos+nd-n,nd) |
---|
113 | rpos = mod(domain%y%pos +n,nd) |
---|
114 | nwords = domain%y%list(rpos)%compute%size * domain%y%list(rpos)%compute%size * size(local,3) |
---|
115 | call mpp_transmit( clocal, size(clocal), domain%y%list(lpos)%pe, cremote, nwords, domain%y%list(rpos)%pe ) |
---|
116 | m = 0 |
---|
117 | do k = 1,size(global,3) |
---|
118 | do j = domain%y%list(rpos)%compute%begin,domain%y%list(rpos)%compute%end |
---|
119 | do i = domain%x%compute%begin,domain%x%compute%end |
---|
120 | m = m + 1 |
---|
121 | global(i,j,k) = cremote(m) |
---|
122 | end do |
---|
123 | end do |
---|
124 | end do |
---|
125 | end do |
---|
126 | call mpp_sync_self(domain%y%list(:)%pe) |
---|
127 | else |
---|
128 | nd = size(domain%list) |
---|
129 | do n = 1,nd-1 |
---|
130 | lpos = mod(domain%pos+nd-n,nd) |
---|
131 | rpos = mod(domain%pos +n,nd) |
---|
132 | nwords = domain%list(rpos)%x%compute%size * domain%list(rpos)%y%compute%size * size(local,3) |
---|
133 | call mpp_transmit( clocal, size(clocal), domain%list(lpos)%pe, cremote, nwords, domain%list(rpos)%pe ) |
---|
134 | m = 0 |
---|
135 | do k = 1,size(global,3) |
---|
136 | do j = domain%list(rpos)%y%compute%begin,domain%list(rpos)%y%compute%end |
---|
137 | do i = domain%list(rpos)%x%compute%begin,domain%list(rpos)%x%compute%end |
---|
138 | m = m + 1 |
---|
139 | global(i,j,k) = cremote(m) |
---|
140 | end do |
---|
141 | end do |
---|
142 | end do |
---|
143 | end do |
---|
144 | ! call mpp_sync_self(domain%list(:)%pe) |
---|
145 | call mpp_sync_self() |
---|
146 | end if |
---|
147 | |
---|
148 | return |
---|
149 | end subroutine MPP_GLOBAL_FIELD_3D_ |
---|
150 | |
---|
151 | subroutine MPP_GLOBAL_FIELD_4D_( domain, local, global, flags ) |
---|
152 | type(domain2D), intent(in) :: domain |
---|
153 | MPP_TYPE_, intent(in) :: local(:,:,:,:) |
---|
154 | MPP_TYPE_, intent(out) :: global(:,:,:,:) |
---|
155 | integer, intent(in), optional :: flags |
---|
156 | MPP_TYPE_ :: local3D (size( local,1),size( local,2),size( local,3)*size(local,4)) |
---|
157 | MPP_TYPE_ :: global3D(size(global,1),size(global,2),size(global,3)*size(local,4)) |
---|
158 | #ifdef use_CRI_pointers |
---|
159 | pointer( lptr, local3D ) |
---|
160 | pointer( gptr, global3D ) |
---|
161 | lptr = LOC(local) |
---|
162 | gptr = LOC(global) |
---|
163 | call mpp_global_field( domain, local3D, global3D, flags ) |
---|
164 | #else |
---|
165 | local3D = RESHAPE( local, SHAPE(local3D) ) |
---|
166 | call mpp_global_field( domain, local3D, global3D, flags ) |
---|
167 | global = RESHAPE( global3D, SHAPE(global) ) |
---|
168 | #endif |
---|
169 | end subroutine MPP_GLOBAL_FIELD_4D_ |
---|
170 | |
---|
171 | subroutine MPP_GLOBAL_FIELD_5D_( domain, local, global, flags ) |
---|
172 | type(domain2D), intent(in) :: domain |
---|
173 | MPP_TYPE_, intent(in) :: local(:,:,:,:,:) |
---|
174 | MPP_TYPE_, intent(out) :: global(:,:,:,:,:) |
---|
175 | integer, intent(in), optional :: flags |
---|
176 | MPP_TYPE_ :: local3D (size( local,1),size( local,2),size( local,3)*size( local,4)*size(local,5)) |
---|
177 | MPP_TYPE_ :: global3D(size(global,1),size(global,2),size(global,3)*size(global,4)*size(local,5)) |
---|
178 | #ifdef use_CRI_pointers |
---|
179 | pointer( lptr, local3D ) |
---|
180 | pointer( gptr, global3D ) |
---|
181 | lptr = LOC(local) |
---|
182 | gptr = LOC(global) |
---|
183 | call mpp_global_field( domain, local3D, global3D, flags ) |
---|
184 | #else |
---|
185 | local3D = RESHAPE( local, SHAPE(local3D) ) |
---|
186 | call mpp_global_field( domain, local3D, global3D, flags ) |
---|
187 | global = RESHAPE( global3D, SHAPE(global) ) |
---|
188 | #endif |
---|
189 | end subroutine MPP_GLOBAL_FIELD_5D_ |
---|
190 | |
---|
191 | subroutine MPP_GLOBAL1D_FIELD_2D_( domain, local, global ) |
---|
192 | type(domain1D), intent(in) :: domain |
---|
193 | MPP_TYPE_, intent(in) :: local(domain%data%begin: ,:) |
---|
194 | MPP_TYPE_, intent(out) :: global(domain%global%begin:,:) |
---|
195 | integer :: i, n, nwords, lpos, rpos |
---|
196 | MPP_TYPE_, allocatable, dimension(:) :: clocal, cremote |
---|
197 | |
---|
198 | if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: must first call mpp_domains_init.' ) |
---|
199 | |
---|
200 | if( size(local,1).NE.domain%data%size .OR. size(global,1).NE.domain%global%size .OR. & |
---|
201 | size(local,2).NE.size(global,2) ) & |
---|
202 | call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: incoming arrays do not match domain.' ) |
---|
203 | |
---|
204 | allocate( clocal (domain%compute%size *size(local,2)) ) |
---|
205 | allocate( cremote(domain%compute%max_size*size(local,2)) ) |
---|
206 | clocal = TRANSFER( local(domain%compute%begin:domain%compute%end,:), clocal ) |
---|
207 | !always fill local directly |
---|
208 | global(domain%compute%begin:domain%compute%end,:) = local(domain%compute%begin:domain%compute%end,:) |
---|
209 | !fill off-domains (note loops begin at an offset of 1) |
---|
210 | n = size(domain%list) |
---|
211 | do i = 1,n-1 |
---|
212 | lpos = mod(domain%pos+n-i,n) |
---|
213 | rpos = mod(domain%pos +i,n) |
---|
214 | nwords = domain%list(rpos)%compute%size * size(local,2) |
---|
215 | call mpp_transmit( clocal, size(clocal), domain%list(lpos)%pe, cremote, nwords, domain%list(rpos)%pe ) |
---|
216 | global(domain%list(rpos)%compute%begin:domain%list(rpos)%compute%end,:) = & |
---|
217 | RESHAPE( cremote(1:nwords), (/domain%list(rpos)%compute%size,size(local,2)/) ) |
---|
218 | end do |
---|
219 | call mpp_sync_self(domain%list(:)%pe) |
---|
220 | !PV786667: the deallocate stmts can be removed when fixed (7.3.1.3m) |
---|
221 | deallocate( clocal, cremote ) |
---|
222 | |
---|
223 | return |
---|
224 | end subroutine MPP_GLOBAL1D_FIELD_2D_ |
---|