Oasis3 4.0.2
|
00001 subroutine MPP_GLOBAL_FIELD_2D_( domain, local, global, flags ) 00002 type(domain2D), intent(in) :: domain 00003 MPP_TYPE_, intent(in) :: local(:,:) 00004 MPP_TYPE_, intent(out) :: global(:,:) 00005 integer, intent(in), optional :: flags 00006 MPP_TYPE_ :: local3D (size( local,1),size( local,2),1) 00007 MPP_TYPE_ :: global3D(size(global,1),size(global,2),1) 00008 #ifdef use_CRI_pointers 00009 pointer( lptr, local3D ) 00010 pointer( gptr, global3D ) 00011 lptr = LOC( local) 00012 gptr = LOC(global) 00013 call mpp_global_field( domain, local3D, global3D, flags ) 00014 #else 00015 local3D = RESHAPE( local, SHAPE(local3D) ) 00016 call mpp_global_field( domain, local3D, global3D, flags ) 00017 global = RESHAPE( global3D, SHAPE(global) ) 00018 #endif 00019 end subroutine MPP_GLOBAL_FIELD_2D_ 00020 00021 subroutine MPP_GLOBAL_FIELD_3D_( domain, local, global, flags ) 00022 !get a global field from a local field 00023 !local field may be on compute OR data domain 00024 type(domain2D), intent(in) :: domain 00025 MPP_TYPE_, intent(in) :: local(:,:,:) 00026 MPP_TYPE_, intent(out) :: global(domain%x%global%begin:,domain%y%global%begin:,:) 00027 integer, intent(in), optional :: flags 00028 integer :: i, j, k, m, n, nd, nwords, lpos, rpos, ioff, joff 00029 logical :: xonly, yonly 00030 MPP_TYPE_ :: clocal (domain%x%compute%size* domain%y%compute%size* size(local,3)) 00031 MPP_TYPE_ :: cremote(domain%x%compute%max_size*domain%y%compute%max_size*size(local,3)) 00032 integer :: words_per_long, stackuse 00033 character(len=8) :: text 00034 #ifdef use_CRI_pointers 00035 pointer( ptr_local, clocal ) 00036 pointer( ptr_remote, cremote ) 00037 00038 ptr_local = LOC(mpp_domains_stack) 00039 ptr_remote = LOC(mpp_domains_stack(size(clocal)+1)) 00040 #endif 00041 00042 if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: must first call mpp_domains_init.' ) 00043 xonly = .FALSE. 00044 yonly = .FALSE. 00045 if( PRESENT(flags) )then 00046 xonly = flags.EQ.XUPDATE 00047 yonly = flags.EQ.YUPDATE 00048 if( .NOT.xonly .AND. .NOT.yonly )call mpp_error( WARNING, 'MPP_GLOBAL_FIELD: you must have flags=XUPDATE or YUPDATE.' ) 00049 end if 00050 00051 #ifdef use_CRI_pointers 00052 words_per_long = size(transfer(local(1,1,1),mpp_domains_stack)) 00053 stackuse = (size(clocal)+size(cremote))*words_per_long 00054 if( stackuse.GT.mpp_domains_stack_size )then 00055 write( text, '(i8)' )stackuse 00056 call mpp_error( FATAL, & 00057 'MPP_UPDATE_DOMAINS user stack overflow: call mpp_domains0_set_stack_size('//trim(text)//') from all PEs.' ) 00058 end if 00059 ptr_local = LOC(mpp_domains_stack) 00060 ptr_remote = LOC( mpp_domains_stack((size(clocal)+1)*words_per_long) ) 00061 mpp_domains_stack_hwm = max( mpp_domains_stack_hwm, stackuse ) 00062 #endif 00063 if( size(global,1).NE.domain%x%global%size .OR. size(global,2).NE.domain%y%global%size .OR. & 00064 size(local,3).NE.size(global,3) ) & 00065 call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: incoming arrays do not match domain.' ) 00066 if( size(local,1).EQ.domain%x%compute%size .AND. size(local,2).EQ.domain%y%compute%size )then 00067 !local is on compute domain 00068 ioff = -domain%x%compute%begin + 1 00069 joff = -domain%y%compute%begin + 1 00070 else if( size(local,1).EQ.domain%x%data%size .AND. size(local,2).EQ.domain%y%data%size )then 00071 !local is on data domain 00072 ioff = -domain%x%data%begin + 1 00073 joff = -domain%y%data%begin + 1 00074 else 00075 call mpp_error( FATAL, 'MPP_GLOBAL_FIELD_: incoming field array must match either compute domain or data domain.' ) 00076 end if 00077 00078 ! make contiguous array from compute domain 00079 m = 0 00080 do k = 1,size(local,3) 00081 do j = domain%y%compute%begin,domain%y%compute%end 00082 do i = domain%x%compute%begin,domain%x%compute%end 00083 m = m + 1 00084 clocal(m) = local(i+ioff,j+joff,k) 00085 global(i,j,k) = clocal(m) !always fill local domain directly 00086 end do 00087 end do 00088 end do 00089 00090 !fill off-domains (note loops begin at an offset of 1) 00091 if( xonly )then 00092 nd = size(domain%x%list) 00093 do n = 1,nd-1 00094 lpos = mod(domain%x%pos+nd-n,nd) 00095 rpos = mod(domain%x%pos +n,nd) 00096 nwords = domain%x%list(rpos)%compute%size * domain%x%list(rpos)%compute%size * size(local,3) 00097 call mpp_transmit( clocal, size(clocal), domain%x%list(lpos)%pe, cremote, nwords, domain%x%list(rpos)%pe ) 00098 m = 0 00099 do k = 1,size(global,3) 00100 do j = domain%y%compute%begin,domain%y%compute%end 00101 do i = domain%x%list(rpos)%compute%begin,domain%x%list(rpos)%compute%end 00102 m = m + 1 00103 global(i,j,k) = cremote(m) 00104 end do 00105 end do 00106 end do 00107 end do 00108 call mpp_sync_self(domain%x%list(:)%pe) 00109 else if( yonly )then 00110 nd = size(domain%y%list) 00111 do n = 1,nd-1 00112 lpos = mod(domain%y%pos+nd-n,nd) 00113 rpos = mod(domain%y%pos +n,nd) 00114 nwords = domain%y%list(rpos)%compute%size * domain%y%list(rpos)%compute%size * size(local,3) 00115 call mpp_transmit( clocal, size(clocal), domain%y%list(lpos)%pe, cremote, nwords, domain%y%list(rpos)%pe ) 00116 m = 0 00117 do k = 1,size(global,3) 00118 do j = domain%y%list(rpos)%compute%begin,domain%y%list(rpos)%compute%end 00119 do i = domain%x%compute%begin,domain%x%compute%end 00120 m = m + 1 00121 global(i,j,k) = cremote(m) 00122 end do 00123 end do 00124 end do 00125 end do 00126 call mpp_sync_self(domain%y%list(:)%pe) 00127 else 00128 nd = size(domain%list) 00129 do n = 1,nd-1 00130 lpos = mod(domain%pos+nd-n,nd) 00131 rpos = mod(domain%pos +n,nd) 00132 nwords = domain%list(rpos)%x%compute%size * domain%list(rpos)%y%compute%size * size(local,3) 00133 call mpp_transmit( clocal, size(clocal), domain%list(lpos)%pe, cremote, nwords, domain%list(rpos)%pe ) 00134 m = 0 00135 do k = 1,size(global,3) 00136 do j = domain%list(rpos)%y%compute%begin,domain%list(rpos)%y%compute%end 00137 do i = domain%list(rpos)%x%compute%begin,domain%list(rpos)%x%compute%end 00138 m = m + 1 00139 global(i,j,k) = cremote(m) 00140 end do 00141 end do 00142 end do 00143 end do 00144 ! call mpp_sync_self(domain%list(:)%pe) 00145 call mpp_sync_self() 00146 end if 00147 00148 return 00149 end subroutine MPP_GLOBAL_FIELD_3D_ 00150 00151 subroutine MPP_GLOBAL_FIELD_4D_( domain, local, global, flags ) 00152 type(domain2D), intent(in) :: domain 00153 MPP_TYPE_, intent(in) :: local(:,:,:,:) 00154 MPP_TYPE_, intent(out) :: global(:,:,:,:) 00155 integer, intent(in), optional :: flags 00156 MPP_TYPE_ :: local3D (size( local,1),size( local,2),size( local,3)*size(local,4)) 00157 MPP_TYPE_ :: global3D(size(global,1),size(global,2),size(global,3)*size(local,4)) 00158 #ifdef use_CRI_pointers 00159 pointer( lptr, local3D ) 00160 pointer( gptr, global3D ) 00161 lptr = LOC(local) 00162 gptr = LOC(global) 00163 call mpp_global_field( domain, local3D, global3D, flags ) 00164 #else 00165 local3D = RESHAPE( local, SHAPE(local3D) ) 00166 call mpp_global_field( domain, local3D, global3D, flags ) 00167 global = RESHAPE( global3D, SHAPE(global) ) 00168 #endif 00169 end subroutine MPP_GLOBAL_FIELD_4D_ 00170 00171 subroutine MPP_GLOBAL_FIELD_5D_( domain, local, global, flags ) 00172 type(domain2D), intent(in) :: domain 00173 MPP_TYPE_, intent(in) :: local(:,:,:,:,:) 00174 MPP_TYPE_, intent(out) :: global(:,:,:,:,:) 00175 integer, intent(in), optional :: flags 00176 MPP_TYPE_ :: local3D (size( local,1),size( local,2),size( local,3)*size( local,4)*size(local,5)) 00177 MPP_TYPE_ :: global3D(size(global,1),size(global,2),size(global,3)*size(global,4)*size(local,5)) 00178 #ifdef use_CRI_pointers 00179 pointer( lptr, local3D ) 00180 pointer( gptr, global3D ) 00181 lptr = LOC(local) 00182 gptr = LOC(global) 00183 call mpp_global_field( domain, local3D, global3D, flags ) 00184 #else 00185 local3D = RESHAPE( local, SHAPE(local3D) ) 00186 call mpp_global_field( domain, local3D, global3D, flags ) 00187 global = RESHAPE( global3D, SHAPE(global) ) 00188 #endif 00189 end subroutine MPP_GLOBAL_FIELD_5D_ 00190 00191 subroutine MPP_GLOBAL1D_FIELD_2D_( domain, local, global ) 00192 type(domain1D), intent(in) :: domain 00193 MPP_TYPE_, intent(in) :: local(domain%data%begin: ,:) 00194 MPP_TYPE_, intent(out) :: global(domain%global%begin:,:) 00195 integer :: i, n, nwords, lpos, rpos 00196 MPP_TYPE_, allocatable, dimension(:) :: clocal, cremote 00197 00198 if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: must first call mpp_domains_init.' ) 00199 00200 if( size(local,1).NE.domain%data%size .OR. size(global,1).NE.domain%global%size .OR. & 00201 size(local,2).NE.size(global,2) ) & 00202 call mpp_error( FATAL, 'MPP_GLOBAL_FIELD: incoming arrays do not match domain.' ) 00203 00204 allocate( clocal (domain%compute%size *size(local,2)) ) 00205 allocate( cremote(domain%compute%max_size*size(local,2)) ) 00206 clocal = TRANSFER( local(domain%compute%begin:domain%compute%end,:), clocal ) 00207 !always fill local directly 00208 global(domain%compute%begin:domain%compute%end,:) = local(domain%compute%begin:domain%compute%end,:) 00209 !fill off-domains (note loops begin at an offset of 1) 00210 n = size(domain%list) 00211 do i = 1,n-1 00212 lpos = mod(domain%pos+n-i,n) 00213 rpos = mod(domain%pos +i,n) 00214 nwords = domain%list(rpos)%compute%size * size(local,2) 00215 call mpp_transmit( clocal, size(clocal), domain%list(lpos)%pe, cremote, nwords, domain%list(rpos)%pe ) 00216 global(domain%list(rpos)%compute%begin:domain%list(rpos)%compute%end,:) = & 00217 RESHAPE( cremote(1:nwords), (/domain%list(rpos)%compute%size,size(local,2)/) ) 00218 end do 00219 call mpp_sync_self(domain%list(:)%pe) 00220 !PV786667: the deallocate stmts can be removed when fixed (7.3.1.3m) 00221 deallocate( clocal, cremote ) 00222 00223 return 00224 end subroutine MPP_GLOBAL1D_FIELD_2D_