39 use mpp_mod,
only : mpp_root_pe, mpp_pe
41 use fms_mod,
only : write_version_number, file_exist, close_file
82 #include<file_version.h> 97 integer :: unit, ierr, io
101 call write_version_number(
"horiz_interp_spherical_mod", version)
102 #ifdef INTERNAL_FILE_NML 106 if (file_exist(
'input.nml'))
then 107 unit = open_namelist_file( )
108 ierr=1;
do while (ierr /= 0)
109 read (unit, nml=horiz_interp_spherical_nml, iostat=io, end=10)
112 10
call close_file (unit)
178 num_nbrs, max_dist, src_modulo)
180 real,
intent(in),
dimension(:,:) :: lon_in, lat_in, lon_out, lat_out
181 integer,
intent(in),
optional :: num_nbrs
182 real,
optional,
intent(in) :: max_dist
183 logical,
intent(in),
optional :: src_modulo
187 integer :: map_dst_xsize, map_dst_ysize, map_src_xsize, map_src_ysize
188 integer :: map_src_size, num_neighbors
189 real :: max_src_dist, tpi, hpi
190 logical :: src_is_modulo
191 real :: min_theta_dst, max_theta_dst, min_phi_dst, max_phi_dst
192 real :: min_theta_src, max_theta_src, min_phi_src, max_phi_src
195 integer :: num_found(
size(lon_out,1),
size(lon_out,2))
197 real,
dimension(size(lon_out,1),size(lon_out,2)) :: theta_dst, phi_dst
198 real,
dimension(size(lon_in,1)*size(lon_in,2)) :: theta_src, phi_src
205 tpi = 2.0*
pi; hpi = 0.5*
pi 208 if(
present(num_nbrs)) num_neighbors = num_nbrs
209 if (num_neighbors <= 0)
call mpp_error(fatal,
'horiz_interp_spherical_mod: num_neighbors must be > 0')
212 if (
PRESENT(max_dist)) max_src_dist = max_dist
213 interp%max_src_dist = max_src_dist
215 src_is_modulo = .true.
216 if (
PRESENT(src_modulo)) src_is_modulo = src_modulo
219 map_dst_xsize=
size(lon_out,1);map_dst_ysize=
size(lon_out,2)
220 map_src_xsize=
size(lon_in,1); map_src_ysize=
size(lon_in,2)
221 map_src_size = map_src_xsize*map_src_ysize
223 if (map_dst_xsize /=
size(lat_out,1) .or. map_dst_ysize /=
size(lat_out,2)) &
224 call mpp_error(fatal,
'horiz_interp_spherical_mod: destination grids not conformable')
225 if (map_src_xsize /=
size(lat_in,1) .or. map_src_ysize /=
size(lat_in,2)) &
226 call mpp_error(fatal,
'horiz_interp_spherical_mod: source grids not conformable')
228 theta_src = reshape(lon_in,(/map_src_size/))
229 phi_src = reshape(lat_in,(/map_src_size/))
230 theta_dst(:,:) = lon_out(:,:)
231 phi_dst(:,:) = lat_out(:,:)
233 min_theta_dst=tpi;max_theta_dst=0.0;min_phi_dst=
pi;max_phi_dst=-
pi 234 min_theta_src=tpi;max_theta_src=0.0;min_phi_src=
pi;max_phi_src=-
pi 236 where(theta_dst<0.0) theta_dst = theta_dst+tpi
237 where(theta_dst>tpi) theta_dst = theta_dst-tpi
238 where(theta_src<0.0) theta_src = theta_src+tpi
239 where(theta_src>tpi) theta_src = theta_src-tpi
241 where(phi_dst < -hpi) phi_dst = -hpi
242 where(phi_dst > hpi) phi_dst = hpi
243 where(phi_src < -hpi) phi_src = -hpi
244 where(phi_src > hpi) phi_src = hpi
248 min_theta_dst =
min(min_theta_dst,theta_dst(i,j))
249 max_theta_dst =
max(max_theta_dst,theta_dst(i,j))
250 min_phi_dst =
min(min_phi_dst,phi_dst(i,j))
251 max_phi_dst =
max(max_phi_dst,phi_dst(i,j))
256 min_theta_src =
min(min_theta_src,theta_src(i))
257 max_theta_src =
max(max_theta_src,theta_src(i))
258 min_phi_src =
min(min_phi_src,phi_src(i))
259 max_phi_src =
max(max_phi_src,phi_src(i))
262 if (min_phi_dst < min_phi_src) print *,
'=> WARNING: latitute of dest grid exceeds src' 263 if (max_phi_dst > max_phi_src) print *,
'=> WARNING: latitute of dest grid exceeds src' 265 if(.not. src_is_modulo)
then 266 if (min_theta_dst < min_theta_src) print *,
'=> WARNING : longitude of dest grid exceeds src' 267 if (max_theta_dst > max_theta_src) print *,
'=> WARNING : longitude of dest grid exceeds src' 271 if(
ASSOCIATED(interp%i_lon))
then 272 if(
size(interp%i_lon,1) .NE. map_dst_xsize .OR. &
273 size(interp%i_lon,2) .NE. map_dst_ysize )
call mpp_error(fatal, &
274 .NE..OR.
'horiz_interp_spherical_mod: size(Interp%i_lon(:),1) map_dst_xsize '// &
275 .NE.
'size(Interp%i_lon(:),2) map_dst_ysize')
277 allocate(interp%i_lon(map_dst_xsize,map_dst_ysize,
max_neighbors), &
279 interp%src_dist(map_dst_xsize,map_dst_ysize,
max_neighbors), &
280 interp%num_found(map_dst_xsize,map_dst_ysize) )
290 case (
"radial_search")
291 call radial_search(theta_src, phi_src, theta_dst, phi_dst, map_src_xsize, map_src_ysize, &
292 map_src_add, map_src_dist, num_found, num_neighbors,max_src_dist,src_is_modulo)
294 call full_search(theta_src, phi_src, theta_dst, phi_dst, map_src_add, map_src_dist, &
295 num_found, num_neighbors,max_src_dist )
297 call mpp_error(fatal,
"horiz_interp_spherical_new: nml search_method = "// &
303 do n=1,num_found(i,j)
304 if(map_src_add(i,j,n) == 0)
then 305 jlat(n) = 0; ilon(n) = 0
307 jlat(n) = map_src_add(i,j,n)/map_src_xsize + 1
308 ilon(n) = map_src_add(i,j,n) - (jlat(n)-1)*map_src_xsize
309 if(ilon(n) == 0)
then 310 jlat(n) = jlat(n) - 1
311 ilon(n) = map_src_xsize
315 interp%i_lon(i,j,:) = ilon(:)
316 interp%j_lat(i,j,:) = jlat(:)
317 interp%num_found(i,j) = num_found(i,j)
318 interp%src_dist(i,j,:) = map_src_dist(i,j,:)
322 interp%nlon_src = map_src_xsize; interp%nlat_src = map_src_ysize
323 interp%nlon_dst = map_dst_xsize; interp%nlat_dst = map_dst_ysize
372 real,
intent(in),
dimension(:,:) :: data_in
373 real,
intent(out),
dimension(:,:) :: data_out
374 integer,
intent(in),
optional :: verbose
375 real,
intent(in),
dimension(:,:),
optional :: mask_in
376 real,
intent(out),
dimension(:,:),
optional :: mask_out
377 real,
intent(in),
optional :: missing_value
380 real,
dimension(Interp%nlon_dst, Interp%nlat_dst,size(Interp%src_dist,3)) :: wt
381 real,
dimension(Interp%nlon_src, Interp%nlat_src) :: mask_src
382 real,
dimension(Interp%nlon_dst, Interp%nlat_dst) :: mask_dst
383 integer :: nlon_in, nlat_in, nlon_out, nlat_out, num_found
384 integer :: m, n, i, j, k, miss_in, miss_out, i1, i2, j1, j2, iverbose
385 real :: min_in, max_in, avg_in, min_out, max_out, avg_out, sum
388 iverbose = 0;
if (
present(verbose)) iverbose = verbose
390 nlon_in = interp%nlon_src; nlat_in = interp%nlat_src
391 nlon_out = interp%nlon_dst; nlat_out = interp%nlat_dst
393 if(
size(data_in,1) .ne. nlon_in .or.
size(data_in,2) .ne. nlat_in ) &
394 call mpp_error(fatal,
'horiz_interp_spherical_mod: size of input array incorrect')
396 if(
size(data_out,1) .ne. nlon_out .or.
size(data_out,2) .ne. nlat_out ) &
397 call mpp_error(fatal,
'horiz_interp_spherical_mod: size of output array incorrect')
399 mask_src = 1.0; mask_dst = 1.0
400 if(
present(mask_in)) mask_src = mask_in
406 num_found = interp%num_found(m,n)
407 if(num_found == 0 )
then 410 i1 = interp%i_lon(m,n,1); j1 = interp%j_lat(m,n,1)
411 if (mask_src(i1,j1) .lt. 0.5)
then 415 if(num_found .gt. 1 )
then 416 i2 = interp%i_lon(m,n,2); j2 = interp%j_lat(m,n,2)
419 if(abs(interp%src_dist(m,n,2)-interp%src_dist(m,n,1)) .lt.
epsln)
then 420 if((mask_src(i1,j1) .lt. 0.5)) mask_dst(m,n) = 0.0
426 if(mask_src(interp%i_lon(m,n,k),interp%j_lat(m,n,k)) .lt. 0.5 )
then 429 if (interp%src_dist(m,n,k) <=
epsln)
then 432 else if(interp%src_dist(m,n,k) <= interp%max_src_dist )
then 433 wt(m,n,k) = 1.0/interp%src_dist(m,n,k)
440 if (sum >
epsln)
then 442 wt(m,n,k) = wt(m,n,k)/sum
454 if(mask_dst(m,n) .gt. 0.5)
then 455 do k=1, interp%num_found(m,n)
456 i = interp%i_lon(m,n,k)
457 j = interp%j_lat(m,n,k)
458 data_out(m,n) = data_out(m,n)+data_in(i,j)*wt(m,n,k)
461 if(
present(missing_value))
then 462 data_out(m,n) = missing_value
470 if(
present(mask_out)) mask_out = mask_dst
476 if (iverbose > 0)
then 480 call stats (data_in, min_in, max_in, avg_in, miss_in, missing_value, mask=mask_src)
483 call stats (data_out, min_out, max_out, avg_out, miss_out, missing_value, mask=mask_dst)
489 write (*,901) min_in ,max_in, avg_in
490 if (
present(mask_in))
write (*,903) miss_in
491 write (*,902) min_out,max_out,avg_out
492 if (
present(mask_out))
write (*,903) miss_out
494 900
format (/,1x,10(
'-'),
' output from horiz_interp ',10(
'-'))
495 901
format (
' input: min=',f16.9,
' max=',f16.9,
' avg=',f22.15)
496 902
format (
' output: min=',f16.9,
' max=',f16.9,
' avg=',f22.15)
497 903
format (
' number of missing points = ',i6)
508 real,
intent(out),
dimension(:,:,:) :: wt
509 integer,
intent(in),
optional :: verbose
510 real,
intent(in),
dimension(:,:),
optional :: mask_in
511 real,
intent(inout),
dimension(:,:),
optional :: mask_out
512 real,
intent(in),
optional :: missing_value
515 real,
dimension(Interp%nlon_src, Interp%nlat_src) :: mask_src
516 real,
dimension(Interp%nlon_dst, Interp%nlat_dst) :: mask_dst
517 integer :: nlon_in, nlat_in, nlon_out, nlat_out, num_found
518 integer :: m, n, i, j, k, miss_in, miss_out, i1, i2, j1, j2, iverbose
519 real :: min_in, max_in, avg_in, min_out, max_out, avg_out, sum
522 iverbose = 0;
if (
present(verbose)) iverbose = verbose
524 nlon_in = interp%nlon_src; nlat_in = interp%nlat_src
525 nlon_out = interp%nlon_dst; nlat_out = interp%nlat_dst
527 mask_src = 1.0; mask_dst = 1.0
528 if(
present(mask_in)) mask_src = mask_in
534 num_found = interp%num_found(m,n)
537 print*,
'pe=',mpp_pe(),
'num_found=',num_found
541 if(num_found == 0 )
then 544 i1 = interp%i_lon(m,n,1); j1 = interp%j_lat(m,n,1)
545 if (mask_src(i1,j1) .lt. 0.5)
then 549 if(num_found .gt. 1 )
then 550 i2 = interp%i_lon(m,n,2); j2 = interp%j_lat(m,n,2)
553 if(abs(interp%src_dist(m,n,2)-interp%src_dist(m,n,1)) .lt.
epsln)
then 554 if((mask_src(i1,j1) .lt. 0.5)) mask_dst(m,n) = 0.0
560 if(mask_src(interp%i_lon(m,n,k),interp%j_lat(m,n,k)) .lt. 0.5 )
then 563 if (interp%src_dist(m,n,k) <=
epsln)
then 566 else if(interp%src_dist(m,n,k) <= interp%max_src_dist )
then 567 wt(m,n,k) = 1.0/interp%src_dist(m,n,k)
574 if (sum >
epsln)
then 576 wt(m,n,k) = wt(m,n,k)/sum
616 if(
associated(interp%src_dist))
deallocate(interp%src_dist)
617 if(
associated(interp%num_found))
deallocate(interp%num_found)
618 if(
associated(interp%i_lon))
deallocate(interp%i_lon)
619 if(
associated(interp%j_lat))
deallocate(interp%j_lat)
627 subroutine radial_search(theta_src,phi_src,theta_dst,phi_dst, map_src_xsize, map_src_ysize, &
628 map_src_add, map_src_dist, num_found, num_neighbors,max_src_dist,src_is_modulo)
629 real,
intent(in),
dimension(:) :: theta_src, phi_src
630 real,
intent(in),
dimension(:,:) :: theta_dst, phi_dst
631 integer,
intent(in) :: map_src_xsize, map_src_ysize
632 integer,
intent(out),
dimension(:,:,:) :: map_src_add
633 real,
intent(out),
dimension(:,:,:) :: map_src_dist
634 integer,
intent(inout),
dimension(:,:) :: num_found
635 integer,
intent(in) :: num_neighbors
636 real,
intent(in) :: max_src_dist
637 logical,
intent(in) :: src_is_modulo
640 integer,
parameter :: max_nbrs = 50
641 integer :: i, j, jj, i0, j0, n, l,i_left, i_right
642 integer :: map_dst_xsize, map_dst_ysize
643 integer :: i_left1, i_left2, i_right1, i_right2
644 integer :: map_src_size, step, step_size, bound, bound_start, bound_end
645 logical :: continue_search, result, continue_radial_search
648 map_dst_xsize=
size(theta_dst,1);map_dst_ysize=
size(theta_dst,2)
649 map_src_size = map_src_xsize*map_src_ysize
653 continue_search=.true.
655 step_size = sqrt(
real(map_src_size) )
656 do while (continue_search .and. step_size > 0)
657 do while (step <= map_src_size .and. continue_search)
660 if (d <= max_src_dist)
then 662 step,d, num_found(i,j), num_neighbors )
665 i0 = mod(step,map_src_xsize)
667 if (i0 == 0) i0 = map_src_xsize
668 res = float(step)/float(map_src_xsize)
670 continue_radial_search = .true.
671 do while (continue_radial_search)
672 continue_radial_search = .false.
674 if(n > max_nbrs)
exit 677 if (i_left <= 0)
then 678 if (src_is_modulo)
then 679 i_left = map_src_xsize + i_left
688 bound = ( 1 - jj )*map_src_xsize - i_left
689 else if ( jj >= map_src_ysize )
then 690 bound = ( 2*map_src_ysize - jj ) * map_src_xsize - i_left
692 bound = jj * map_src_xsize + i_left
696 if(d<=max_src_dist)
then 698 bound,d, num_found(i,j), num_neighbors)
699 if (result) continue_radial_search = .true.
705 if (i_right > map_src_xsize)
then 706 if (src_is_modulo)
then 707 i_right = i_right - map_src_xsize
709 i_right = map_src_xsize
716 bound = ( 1 - jj )*map_src_xsize - i_right
717 else if ( jj >= map_src_ysize )
then 718 bound = ( 2*map_src_ysize - jj) * map_src_xsize - i_right
721 bound = jj * map_src_xsize + i_right
725 if(d<=max_src_dist)
then 727 bound,d, num_found(i,j), num_neighbors)
728 if (result) continue_radial_search = .true.
734 if( i_left > i_right)
then 738 i_right2 = map_src_xsize
746 bound_start = ( 1 - jj)*map_src_xsize - i_right1
747 bound_end = ( 1 - jj)*map_src_xsize - i_left1
749 bound_start = jj * map_src_xsize + i_left1
750 bound_end = jj * map_src_xsize + i_right1
754 do while (bound <= bound_end)
756 if(d<=max_src_dist)
then 758 bound,d, num_found(i,j), num_neighbors)
759 if (result) continue_radial_search = .true.
765 if(i_left2 > 0 )
then 767 bound_start = ( 1 - jj)*map_src_xsize - i_right2
768 bound_end = ( 1 - jj)*map_src_xsize - i_left2
770 bound_start = jj * map_src_xsize + i_left2
771 bound_end = jj * map_src_xsize + i_right2
775 do while (bound <= bound_end)
777 if(d<=max_src_dist)
then 779 bound,d, num_found(i,j), num_neighbors)
780 if (result) continue_radial_search = .true.
788 if( jj >= map_src_ysize)
then 789 bound_start = ( 2*map_src_ysize - jj ) * map_src_xsize - i_right1
790 bound_end = ( 2*map_src_ysize - jj ) * map_src_xsize - i_left1
792 bound_start = jj * map_src_xsize + i_left1
793 bound_end = jj * map_src_xsize + i_right1
797 do while (bound <= bound_end)
799 if(d<=max_src_dist)
then 801 bound,d, num_found(i,j), num_neighbors)
802 if (result) continue_radial_search = .true.
808 if( jj >= map_src_ysize)
then 809 bound_start = ( 2*map_src_ysize - jj ) * map_src_xsize - i_right2
810 bound_end = ( 2*map_src_ysize - jj ) * map_src_xsize - i_left2
812 bound_start = jj * map_src_xsize + i_left2
813 bound_end = jj * map_src_xsize + i_right2
817 do while (bound <= bound_end)
819 if(d<=max_src_dist)
then 821 bound,d, num_found(i,j), num_neighbors)
822 if (result) continue_radial_search = .true.
829 continue_search = .false.
835 step_size = step_size/2
849 integer,
intent(inout),
dimension(:) :: map_src_add
850 real,
intent(inout),
dimension(:) :: map_src_dist
851 integer,
intent(in) :: src_add
852 real,
intent(in) :: d
853 integer,
intent(inout) :: num_found
854 integer,
intent(in) :: min_nbrs
863 nloop :
do while ( n .le. num_found )
865 dist_chk :
if (d .le. map_src_dist(n))
then 867 if (src_add == map_src_add(m))
then 868 already_exist = .true.
873 num_found = num_found + 1
875 call mpp_error(fatal,
'update_dest_neighbors: '// &
876 'number of neighbor points found is greated than maxium neighbor points' )
878 do m=num_found,n+1,-1
879 map_src_add(m) = map_src_add(m-1)
880 map_src_dist(m) = map_src_dist(m-1)
882 map_src_add(n) = src_add
885 if( num_found > min_nbrs )
then 886 if( map_src_dist(num_found) > map_src_dist(num_found-1) )
then 887 num_found = num_found - 1
889 if( map_src_dist(min_nbrs+1) > map_src_dist(min_nbrs) )
then 896 if(already_exist)
return 899 if( num_found < min_nbrs )
then 900 num_found = num_found + 1
902 map_src_add(num_found) = src_add
903 map_src_dist(num_found) = d
959 real,
intent(in) :: theta1, phi1, theta2, phi2
962 if(theta1 == theta2 .and. phi1 == phi2)
then 967 dot = cos(phi1)*cos(phi2)*cos(theta1-theta2) + sin(phi1)*sin(phi2)
979 subroutine full_search(theta_src,phi_src,theta_dst,phi_dst,map_src_add, map_src_dist,num_found, &
980 num_neighbors,max_src_dist)
981 real,
intent(in),
dimension(:) :: theta_src, phi_src
982 real,
intent(in),
dimension(:,:) :: theta_dst, phi_dst
983 integer,
intent(out),
dimension(:,:,:) :: map_src_add
984 real,
intent(out),
dimension(:,:,:) :: map_src_dist
985 integer,
intent(out),
dimension(:,:) :: num_found
986 integer,
intent(in) :: num_neighbors
987 real,
intent(in) :: max_src_dist
989 integer :: i,j,map_src_size, step
990 integer :: map_dst_xsize,map_dst_ysize
994 map_dst_xsize=
size(theta_dst,1);map_dst_ysize=
size(theta_dst,2)
995 map_src_size =
size(theta_src(:))
999 do step = 1, map_src_size
1001 if( d <= max_src_dist)
then 1003 step,d,num_found(i,j), num_neighbors )
real function spherical_distance(theta1, phi1, theta2, phi2)
logical module_is_initialized
integer root_pe
NAMELIST NAME="horiz_interp_spherical_nml">
logical function update_dest_neighbors(map_src_add, map_src_dist, src_add, d, num_found, min_nbrs)
subroutine, public horiz_interp_spherical(Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value)
subroutine, public stats(dat, low, high, avg, miss, missing_value, mask)
subroutine full_search(theta_src, phi_src, theta_dst, phi_dst, map_src_add, map_src_dist, num_found, num_neighbors, max_src_dist)
integer function, public check_nml_error(IOSTAT, NML_NAME)
character(len=32) search_method
subroutine radial_search(theta_src, phi_src, theta_dst, phi_dst, map_src_xsize, map_src_ysize, map_src_add, map_src_dist, num_found, num_neighbors, max_src_dist, src_is_modulo)
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
subroutine, public horiz_interp_spherical_init
integer, parameter max_neighbors
integer, parameter num_nbrs_default
subroutine, public horiz_interp_spherical_new(Interp, lon_in, lat_in, lon_out, lat_out, num_nbrs, max_dist, src_modulo)
real, parameter max_dist_default
subroutine, public horiz_interp_spherical_del(Interp)
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************subroutine READ_RECORD_CORE_(unit, field, nwords, data, start, axsiz) integer, intent(in) ::unit type(fieldtype), intent(in) ::field integer, intent(in) ::nwords MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in) ::start(:), axsiz(:) integer(SHORT_KIND) ::i2vals(nwords)!rab used in conjunction with transfer intrinsic to determine size of a variable integer(KIND=1) ::one_byte(8) integer ::word_sz!#ifdef __sgi integer(INT_KIND) ::ivals(nwords) real(FLOAT_KIND) ::rvals(nwords)!#else! integer ::ivals(nwords)! real ::rvals(nwords)!#endif real(DOUBLE_KIND) ::r8vals(nwords) pointer(ptr1, i2vals) pointer(ptr2, ivals) pointer(ptr3, rvals) pointer(ptr4, r8vals) if(mpp_io_stack_size< nwords) call mpp_io_set_stack_size(nwords) call mpp_error(FATAL, 'MPP_READ currently requires use_netCDF option') end subroutine READ_RECORD_CORE_ subroutine READ_RECORD_(unit, field, nwords, data, time_level, domain, position, tile_count, start_in, axsiz_in)!routine that is finally called by all mpp_read routines to perform the read!a non-netCDF record contains:! field ID! a set of 4 coordinates(is:ie, js:je) giving the data subdomain! a timelevel and a timestamp(=NULLTIME if field is static)! 3D real data(stored as 1D)!if you are using direct access I/O, the RECL argument to OPEN must be large enough for the above!in a global direct access file, record position on PE is given by %record.!Treatment of timestamp:! We assume that static fields have been passed without a timestamp.! Here that is converted into a timestamp of NULLTIME.! For non-netCDF fields, field is treated no differently, but is written! with a timestamp of NULLTIME. There is no check in the code to prevent! the user from repeatedly writing a static field. integer, intent(in) ::unit, nwords type(fieldtype), intent(in) ::field MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in), optional ::time_level type(domain2D), intent(in), optional ::domain integer, intent(in), optional ::position, tile_count integer, intent(in), optional ::start_in(:), axsiz_in(:) integer, dimension(size(field%axes(:))) ::start, axsiz integer ::tlevel !, subdomain(4) integer ::i, error, is, ie, js, je, isg, ieg, jsg, jeg type(domain2d), pointer ::io_domain=> tlevel if(PRESENT(start_in) .AND. PRESENT(axsiz_in)) then if(size(start(! the data domain and compute domain must refer to the subdomain being passed ! In this ! since that attempts to gather all data on PE size(field%axes(:)) axsiz(i)
subroutine, public horiz_interp_spherical_wght(Interp, wt, verbose, mask_in, mask_out, missing_value)
double dot(const double *p1, const double *p2)
real(fp), parameter, public pi