WAVEWATCH III  beta 0.0.1
scrip_interface Module Reference

Data Types

type  weight_data
 

Functions/Subroutines

subroutine scrip_clear
 
subroutine scrip (src_num, dst_num, l_master, l_read, l_test)
 
subroutine scrip_driverexit (errorCode, errormsg)
 
subroutine sort_add_v2 (add1, add2, weights)
 

Variables

type(weight_data), dimension(:), allocatable wgtdata
 

Function/Subroutine Documentation

◆ scrip()

subroutine scrip_interface::scrip ( integer (scrip_i4), intent(in)  src_num,
integer (scrip_i4), intent(in)  dst_num,
logical (scrip_logical), intent(in)  l_master,
logical( scrip_logical), intent(in)  l_read,
logical(scrip_logical), intent(in)  l_test 
)

Definition at line 330 of file scrip_interface.F90.

330 
331  !-----------------------------------------------------------------------
332 
333  use scrip_kindsmod ! module defining data types
334  use scrip_constants ! module for common constants
335  use scrip_iounitsmod ! I/O unit manager
336  use scrip_timers ! CPU timers
337  use scrip_grids ! module with grid information
338  use scrip_remap_vars ! common remapping variables
339  use scrip_remap_conservative ! routines for conservative remap
340 
341 #ifdef W3_SCRIPNC
342  use scrip_remap_write ! routines for remap output
343  use scrip_remap_read ! routines for remap input
344 #endif
345 
346  use scrip_errormod
347 
348  implicit none
349 
350  !-----------------------------------------------------------------------
351  !
352  ! input variables formerly part of namelist
353  !
354  !-----------------------------------------------------------------------
355 
356  character (SCRIP_charLength) :: &
357  interp_file1,& ! filename for output remap data (map1)
358  interp_file2,& ! filename for output remap data (map2)
359  map1_name, & ! name for mapping from grid1 to grid2
360  map2_name, & ! name for mapping from grid2 to grid1
361  map_method, & ! choice for mapping method
362  normalize_opt,&! option for normalizing weights
363  output_opt ! option for output conventions
364 
365  integer (SCRIP_i4) :: &
366  nmap ! number of mappings to compute (1 or 2)
367 
368  !-----------------------------------------------------------------------
369  !
370  ! input variables not part of namelist
371  !
372  !-----------------------------------------------------------------------
373 
374  integer (SCRIP_i4), intent(in) :: dst_num ! number of destination grid GDST
375  integer (SCRIP_i4), intent(in) :: src_num ! number of source grid GSRC
376  logical (SCRIP_Logical), intent(in) :: l_master ! Am I the master processor (do I/O)?
377  logical( SCRIP_Logical), intent(in) :: l_read ! Do I read the remap file?
378  logical(SCRIP_Logical), intent(in) :: l_test ! Whether to include test output
379  ! in subroutines
380 
381  !-----------------------------------------------------------------------
382  !
383  ! local variables
384  !
385  !-----------------------------------------------------------------------
386 
387  integer (SCRIP_i4) :: n, & ! dummy counter
388  iunit ! unit number for namelist file
389 
390  integer (SCRIP_i4) :: &
391  errorCode ! error flag
392 
393  character (12), parameter :: &
394  rtnName = 'SCRIP_driver'
395 
396 #ifdef W3_SCRIPNC
397  character (LEN=3) :: cdst ! 3 character number of destination map
398  character (LEN=3) :: csrc ! 3 character number of source map
399 #endif
400 #ifdef W3_T38
401  CHARACTER (LEN=10) :: CDATE_TIME(3)
402  INTEGER :: DATE_TIME(8)
403  INTEGER :: ELAPSED_TIME, BEG_TIME, END_TIME
404 #endif
405 
406  !-----------------------------------------------------------------------
407  !
408  ! initialize timers and errorcode
409  !
410  !-----------------------------------------------------------------------
411 
412 #ifdef W3_T38
413  if(l_master)write(scrip_stdout,*)'subroutine scrip'
414 #endif
415 
416  call timers_init
417  do n=1,max_timers
418  call timer_clear(n)
419  end do
420 
421  errorcode = scrip_success
422 
423  !-----------------------------------------------------------------------
424  !
425  ! set variables that were previously read in as a namelist
426  !
427  !-----------------------------------------------------------------------
428 
429  num_maps = 1
430 #ifdef W3_SCRIPNC
431  ! Note: Only master does I/O, but all processors need to know about
432  ! file existence
433  interp_file1 = "rmp_src_to_dst_conserv_XXX_XXX.nc"
434  interp_file2 = 'not_used.nc'
435  map1_name = 'source to destination Conservative Mapping'
436  map2_name = 'map not used'
437  write(cdst, "(i3.3)") dst_num
438  write(csrc, "(i3.3)") src_num
439  interp_file1(24:26) = csrc
440  interp_file1(28:30) = cdst
441 #endif
442 
443  map_method = 'conservative'
444  normalize_opt = 'fracarea'
445  output_opt = 'scrip'
446  restrict_type = 'latitude'
447  num_srch_bins = 90
448  luse_grid1_area = .false.
449  luse_grid2_area = .false.
450  npseg=11 ! or num_polar_segs
451  north_thresh=1.5_scrip_r8 ! or npole_threshold
452  south_thresh=-1.5_scrip_r8 ! or spole_threshold
453  nthreads=2 ! or num_threads
454 
455  select case(map_method)
456  case ('conservative')
458  luse_grid_centers = .false.
459  case ('bilinear')
461  luse_grid_centers = .true.
462  case ('bicubic')
464  luse_grid_centers = .true.
465  case ('distwgt')
467  luse_grid_centers = .true.
468  case ('particle')
470  luse_grid_centers = .false.
471  case default
472  call scrip_errorset(errorcode, rtnname, 'unknown mapping method')
473  call scrip_driverexit(errorcode, 'unknown mapping method')
474  end select
475 
476  select case(normalize_opt(1:4))
477  case ('none')
479  case ('frac')
481  case ('dest')
483  case default
484  call scrip_errorset(errorcode, rtnname, 'unknown normalization option')
485  call scrip_driverexit(errorcode, 'unknown normalization option')
486  end select
487 
488  !-----------------------------------------------------------------------
489  !
490  ! initialize grid information for both grids
491  !
492  !-----------------------------------------------------------------------
493 
494 #ifdef W3_T38
495  if(l_master)write(scrip_stdout,*)'calling grid_init'
496 #endif
497 
498  call grid_init( errorcode,l_master,l_test)
499 
500 #ifdef W3_T38
501  if(l_master)write(scrip_stdout, *) 'Computing remappings between: ',grid1_name
502  if(l_master)write(scrip_stdout, *) ' and ',grid2_name
503 #endif
504 
505  !-----------------------------------------------------------------------
506  !
507  ! initialize some remapping variables.
508  !
509  !-----------------------------------------------------------------------
510 
511  call init_remap_vars
512 
513  !-----------------------------------------------------------------------
514  !
515  ! call appropriate interpolation setup routine based on type of
516  ! remapping requested. or read in remapping data.
517 #ifdef W3_SCRIPNC
518  !or, read in remapping data.
519 #endif
520  !
521  !-----------------------------------------------------------------------
522 
523 #ifdef W3_T38
524  call date_and_time (cdate_time(1), cdate_time(2), cdate_time(3), date_time)
525  beg_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
526 #endif
527 
528 #ifdef W3_SCRIPNC
529  if (l_read) then
530  if(l_master)write(scrip_stdout, *) 'Reading remapping data from ', interp_file1
531  call read_remap_ww3(map1_name, interp_file1, errorcode)
532 #endif
533 #ifdef W3_T38
534  call date_and_time (cdate_time(1), cdate_time(2), cdate_time(3), date_time)
535  end_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
536  elapsed_time = end_time - beg_time
537  write(0,*) "SCRIP: READING ", elapsed_time, " MSEC"
538 #endif
539 #ifdef W3_SCRIPNC
540  else
541 #endif
542  select case(map_type)
543  case(map_type_conserv)
544 #ifdef W3_T38
545  if(l_master)write(scrip_stdout,*)'calling remap_conserv'
546 #endif
547  call remap_conserv(l_master,l_test)
548 #ifdef W3_T38
549  if(l_master)write(scrip_stdout,*)'back from remap_conserv'
550  call date_and_time (cdate_time(1), cdate_time(2), cdate_time(3), date_time)
551  end_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
552  elapsed_time = end_time - beg_time
553  write(0,*) "SCRIP: CALCULATING ", elapsed_time, " MSEC"
554 #endif
555  case default
556  call scrip_errorset(errorcode, rtnname, 'Invalid Map Type')
557  call scrip_driverexit(errorcode, 'Invalid Map Type')
558  end select
559 
560  !-----------------------------------------------------------------------
561  !
562  ! reduce size of remapping arrays
563  !
564  !-----------------------------------------------------------------------
565 
566 #ifdef W3_T38
567  call date_and_time (cdate_time(1), cdate_time(2), cdate_time(3), date_time)
568  beg_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
569 #endif
570 
571  if (num_links_map1 /= max_links_map1) then
573  endif
574  if ((num_maps > 1) .and. (num_links_map2 /= max_links_map2)) then
576  endif
577 
578  call sort_add_v2(grid2_add_map1, grid1_add_map1, wts_map1)
579 
580  !-----------------------------------------------------------------------
581  !
582 #ifdef W3_SCRIPNC
583  ! write remapping info to a file.
584 #endif
585  !
586  !-----------------------------------------------------------------------
587 
588 #ifdef W3_SCRIPNC
589  if (l_master) then
590  write(scrip_stdout, *) 'Writing remapping data to ', interp_file1
591  endif
592 #endif
593 
594 #ifdef W3_T38
595  call date_and_time (cdate_time(1), cdate_time(2), cdate_time(3), date_time)
596  end_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
597  elapsed_time = end_time - beg_time
598  write(0,*) "SCRIP: RESIZING ", elapsed_time, " MSEC"
599  call date_and_time (cdate_time(1), cdate_time(2), cdate_time(3), date_time)
600  beg_time = ((date_time(5)*60 + date_time(6))*60 +date_time(7))*1000 + date_time(8)
601 #endif
602 
603  ! Use write_remap if you want the extra variables in the .nc files for diagnostics
604  ! Use write_remap_ww3 if you don't want any extra variables in the .nc files
605 
606 #ifdef W3_SCRIPNC
607  if(l_test)then
608  call write_remap(map1_name, map2_name, interp_file1, interp_file2, &
609  output_opt, l_master, errorcode)
610  else
611  call write_remap_ww3(map1_name, interp_file1, output_opt, &
612  l_master, errorcode)
613  endif
614 #endif
615 
616 #ifdef W3_SCRIPNC
617  end if
618 #endif
619 
620  !-----------------------------------------------------------------------
621 

References scrip_remap_vars::grid1_add_map1, scrip_grids::grid1_name, scrip_remap_vars::grid2_add_map1, scrip_grids::grid2_name, scrip_grids::grid_init(), scrip_remap_vars::init_remap_vars(), scrip_grids::luse_grid1_area, scrip_grids::luse_grid2_area, scrip_grids::luse_grid_centers, scrip_remap_vars::map_type, scrip_remap_vars::map_type_bicubic, scrip_remap_vars::map_type_bilinear, scrip_remap_vars::map_type_conserv, scrip_remap_vars::map_type_distwgt, scrip_remap_vars::map_type_particle, scrip_remap_vars::max_links_map1, scrip_remap_vars::max_links_map2, scrip_timers::max_timers, scrip_remap_vars::norm_opt, scrip_remap_vars::norm_opt_dstarea, scrip_remap_vars::norm_opt_frcarea, scrip_remap_vars::norm_opt_none, scrip_grids::north_thresh, scrip_grids::npseg, scrip_remap_conservative::nthreads, scrip_remap_vars::num_links_map1, scrip_remap_vars::num_links_map2, scrip_remap_vars::num_maps, scrip_grids::num_srch_bins, scrip_remap_read::read_remap_ww3(), scrip_remap_conservative::remap_conserv(), scrip_remap_vars::resize_remap_vars(), scrip_grids::restrict_type, scrip_driverexit(), scrip_errormod::scrip_errorset(), scrip_iounitsmod::scrip_stdout, scrip_errormod::scrip_success, sort_add_v2(), scrip_grids::south_thresh, scrip_timers::timer_clear(), scrip_timers::timers_init(), scrip_remap_write::write_remap(), scrip_remap_write::write_remap_ww3(), and scrip_remap_vars::wts_map1.

Referenced by wmscrpmd::scrip_wrapper().

◆ scrip_clear()

subroutine scrip_interface::scrip_clear

Definition at line 82 of file scrip_interface.F90.

82  !#######################################################################
83 
84  ! 1. Original author :
85  !
86  ! Erick Rogers, NRL
87  !
88  ! 2. Last update :
89  !
90  ! See revisions.
91  !
92  ! 3. Revisions :
93  !
94  ! 5-May-2011 : Origination ( version 4.01 )
95  !
96  ! 4. Copyright :
97  !
98  ! 5. Purpose :
99  !
100  ! "Clear" all variables declared at module level of SCRIP routines
101  ! (clear "common block" equivalent)
102  !
103  ! 6. Method :
104  !
105  ! rules:
106  ! - if not an array (scalar), set to zero or other start value
107  ! - if dimensioned array, set to zero
108  ! - if allocatable array, deallocate
109  ! - private variables: ignore, since we would need to
110  ! make them public in order to clear them, which may do more
111  ! harm than good.
112  !
113  ! 7. Parameters, Variables and types :
114  !
115  ! 8. Called by :
116  !
117  ! Subroutine SCRIP_interface
118  !
119  ! 9. Subroutines and functions used :
120  !
121  ! None
122  !
123  ! 10. Error messages:
124  !
125  ! 11. Remarks :
126  !
127  ! We "clear" all variables with "save" attribute,
128  ! both "module variables" and "subroutine variables"
129  ! including all variables that are initialized with a value
130  ! in the type declaration, e.g. "real :: x=5.0"
131  !
132  ! 12. Structure :
133  !
134  ! 13. Switches :
135  !
136  ! 14. Source code :
137 
138  use scrip_kindsmod
139  use scrip_timers
140  use scrip_remap_vars
142  use scrip_iounitsmod
143  use scrip_grids
144 
145  implicit none
146 
147  call timers_init ! takes care of all variables in timers.f
148 
149  !.....scrip_remap_vars.f :
154  num_maps=0
155  num_wts=0
156  map_type=0
157  norm_opt=0
159  if(allocated(grid1_add_map1))deallocate(grid1_add_map1)
160  if(allocated(grid2_add_map1))deallocate(grid2_add_map1)
161  if(allocated(grid1_add_map2))deallocate(grid1_add_map2)
162  if(allocated(grid2_add_map2))deallocate(grid2_add_map2)
163  if(allocated(wts_map1))deallocate(wts_map1)
164  if(allocated(wts_map2))deallocate(wts_map2)
165 
166  !.....remap_conserv.f :
167  !.....scalars:
173  avoid_pole_count = 0
194  !.....arrays :
195  if(allocated(link_add1))deallocate(link_add1)
196  if(allocated(link_add2))deallocate(link_add2)
197 
203 
204  if(allocated(srch_add_find_adj_cell))deallocate(srch_add_find_adj_cell)
209 
210  if(allocated(srch_add_locate_segstart))deallocate(srch_add_locate_segstart)
215 
216  if(allocated(srch_add_locate_point))deallocate(srch_add_locate_point)
221 
222  !.....scrip_grids.f :
223  grid1_size=0
224  grid2_size=0
225  grid1_rank=0
226  grid2_rank=0
227  grid1_corners=0
228  grid2_corners=0
229  grid1_name=''
230  grid2_name=''
231  grid1_units=''
232  grid2_units=''
233  luse_grid_centers=.false.
234  luse_grid1_area=.false.
235  luse_grid2_area=.false.
236  restrict_type=''
237  num_srch_bins=0
238  if(allocated(bin_addr1))deallocate(bin_addr1)
239  if(allocated(bin_addr2))deallocate(bin_addr2)
240  if(allocated(bin_lats))deallocate(bin_lats)
241  if(allocated(bin_lons))deallocate(bin_lons)
242  if(allocated(grid1_dims))deallocate(grid1_dims)
243  if(allocated(grid2_dims))deallocate(grid2_dims)
244  if(allocated(grid1_mask))deallocate(grid1_mask)
245  if(allocated(grid2_mask))deallocate(grid2_mask)
246  if(allocated(grid1_center_lat))deallocate(grid1_center_lat)
247  if(allocated(grid1_center_lon))deallocate(grid1_center_lon)
248  if(allocated(grid2_center_lat))deallocate(grid2_center_lat)
249  if(allocated(grid2_center_lon))deallocate(grid2_center_lon)
250  if(allocated(grid1_area))deallocate(grid1_area)
251  if(allocated(grid2_area))deallocate(grid2_area)
252  if(allocated(grid1_area_in))deallocate(grid1_area_in)
253  if(allocated(grid2_area_in))deallocate(grid2_area_in)
254  if(allocated(grid1_frac))deallocate(grid1_frac)
255  if(allocated(grid2_frac))deallocate(grid2_frac)
256  if(allocated(grid1_corner_lat))deallocate(grid1_corner_lat)
257  if(allocated(grid1_corner_lon))deallocate(grid1_corner_lon)
258  if(allocated(grid2_corner_lat))deallocate(grid2_corner_lat)
259  if(allocated(grid2_corner_lon))deallocate(grid2_corner_lon)
260  if(allocated(grid1_bound_box))deallocate(grid1_bound_box)
261  if(allocated(grid2_bound_box))deallocate(grid2_bound_box)
262  if(allocated(special_polar_cell1))deallocate(special_polar_cell1)
263  if(allocated(special_polar_cell2))deallocate(special_polar_cell2)
264  if(allocated(grid1_centroid_lat))deallocate(grid1_centroid_lat)
265  if(allocated(grid1_centroid_lon))deallocate(grid1_centroid_lon)
266  if(allocated(grid2_centroid_lat))deallocate(grid2_centroid_lat)
267  if(allocated(grid2_centroid_lon))deallocate(grid2_centroid_lon)
268 
269  !#######################################################################

References scrip_remap_conservative::avoid_pole_count, scrip_remap_conservative::avoid_pole_offset, scrip_grids::bin_addr1, scrip_grids::bin_addr2, scrip_grids::bin_lats, scrip_grids::bin_lons, scrip_remap_conservative::first_call_find_adj_cell, scrip_remap_conservative::first_call_get_srch_cells, scrip_remap_conservative::first_call_locate_point, scrip_remap_conservative::first_call_locate_segstart, scrip_remap_conservative::first_call_store_link_cnsrv, scrip_remap_vars::grid1_add_map1, scrip_remap_vars::grid1_add_map2, scrip_grids::grid1_area, scrip_grids::grid1_area_in, scrip_grids::grid1_bound_box, scrip_grids::grid1_center_lat, scrip_grids::grid1_center_lon, scrip_grids::grid1_centroid_lat, scrip_grids::grid1_centroid_lon, scrip_grids::grid1_corner_lat, scrip_grids::grid1_corner_lon, scrip_grids::grid1_corners, scrip_grids::grid1_dims, scrip_grids::grid1_frac, scrip_grids::grid1_mask, scrip_grids::grid1_name, scrip_grids::grid1_rank, scrip_grids::grid1_size, scrip_grids::grid1_units, scrip_remap_vars::grid2_add_map1, scrip_remap_vars::grid2_add_map2, scrip_grids::grid2_area, scrip_grids::grid2_area_in, scrip_grids::grid2_bound_box, scrip_grids::grid2_center_lat, scrip_grids::grid2_center_lon, scrip_grids::grid2_centroid_lat, scrip_grids::grid2_centroid_lon, scrip_grids::grid2_corner_lat, scrip_grids::grid2_corner_lon, scrip_grids::grid2_corners, scrip_grids::grid2_dims, scrip_grids::grid2_frac, scrip_grids::grid2_mask, scrip_grids::grid2_name, scrip_grids::grid2_rank, scrip_grids::grid2_size, scrip_grids::grid2_units, scrip_remap_conservative::last_cell_add_get_srch_cells, scrip_remap_conservative::last_cell_find_adj_cell, scrip_remap_conservative::last_cell_grid_num_find_adj_cell, scrip_remap_conservative::last_cell_grid_num_get_srch_cells, scrip_remap_conservative::last_cell_grid_num_locate_point, scrip_remap_conservative::last_cell_grid_num_locate_segstart, scrip_remap_conservative::last_cell_locate_point, scrip_remap_conservative::last_cell_locate_segstart, scrip_remap_conservative::last_srch_grid_num_get_srch_cells, scrip_remap_conservative::last_srch_grid_num_locate_point, scrip_remap_conservative::last_srch_grid_num_locate_segstart, scrip_remap_conservative::link_add1, scrip_remap_conservative::link_add2, scrip_grids::luse_grid1_area, scrip_grids::luse_grid2_area, scrip_grids::luse_grid_centers, scrip_remap_vars::map_type, scrip_remap_vars::max_links_map1, scrip_remap_vars::max_links_map2, scrip_remap_vars::norm_opt, scrip_remap_vars::num_links_map1, scrip_remap_vars::num_links_map2, scrip_remap_vars::num_maps, scrip_grids::num_srch_bins, scrip_remap_conservative::num_srch_cell_locate_points, scrip_remap_conservative::num_srch_cells_find_adj_cell, scrip_remap_conservative::num_srch_cells_loc_get_srch_cells, scrip_remap_conservative::num_srch_cells_locate_segstart, scrip_remap_vars::num_wts, scrip_remap_vars::resize_increment, scrip_grids::restrict_type, scrip_grids::special_polar_cell1, scrip_grids::special_polar_cell2, scrip_remap_conservative::srch_add_find_adj_cell, scrip_remap_conservative::srch_add_loc_get_srch_cells, scrip_remap_conservative::srch_add_locate_point, scrip_remap_conservative::srch_add_locate_segstart, scrip_remap_conservative::srch_center_lat_find_adj_cell, scrip_remap_conservative::srch_center_lat_loc_get_srch_cells, scrip_remap_conservative::srch_center_lat_locate_point, scrip_remap_conservative::srch_center_lat_locate_segstart, scrip_remap_conservative::srch_center_lon_find_adj_cell, scrip_remap_conservative::srch_center_lon_loc_get_srch_cells, scrip_remap_conservative::srch_center_lon_locate_point, scrip_remap_conservative::srch_center_lon_locate_segstart, scrip_remap_conservative::srch_corner_lat_find_adj_cell, scrip_remap_conservative::srch_corner_lat_loc_get_srch_cells, scrip_remap_conservative::srch_corner_lat_locate_point, scrip_remap_conservative::srch_corner_lat_locate_segstart, scrip_remap_conservative::srch_corner_lon_find_adj_cell, scrip_remap_conservative::srch_corner_lon_loc_get_srch_cells, scrip_remap_conservative::srch_corner_lon_locate_point, scrip_remap_conservative::srch_corner_lon_locate_segstart, scrip_remap_conservative::srch_corners_find_adj_cell, scrip_remap_conservative::srch_corners_loc_get_srch_cells, scrip_remap_conservative::srch_corners_locate_point, scrip_remap_conservative::srch_corners_locate_segstart, scrip_timers::timers_init(), scrip_remap_vars::wts_map1, and scrip_remap_vars::wts_map2.

Referenced by wmscrpmd::scrip_wrapper().

◆ scrip_driverexit()

subroutine scrip_interface::scrip_driverexit ( integer (scrip_i4), intent(in)  errorCode,
character*(*), intent(in)  errormsg 
)

Definition at line 626 of file scrip_interface.F90.

626 
627  ! !DESCRIPTION:
628  ! This routine exits the SCRIP driver program. It first calls the
629  ! SCRIP error print function to print any errors encountered and then
630  ! exits the message environment before stopping.
631  !
632  ! !USES:
633 
634  use scrip_kindsmod
635 
636  ! !INPUT PARAMETERS:
637 
638  integer (SCRIP_i4), intent(in) :: &
639  errorCode ! error flag to detect any errors encountered
640 
641  CHARACTER*(*), INTENT(IN) :: errormsg
642  !-----------------------------------------------------------------------
643  !
644  ! call SCRIP error print function to output any logged errors that
645  ! were encountered during execution. Then stop.
646  !
647  !-----------------------------------------------------------------------
648 
649  write(*,*)'error encountered : ',errorcode
650  write(*,*)errormsg
651 
652  stop
653 
654  !-----------------------------------------------------------------------
655 

Referenced by scrip().

◆ sort_add_v2()

subroutine scrip_interface::sort_add_v2 ( integer (scrip_i4), dimension(:), intent(inout)  add1,
integer (scrip_i4), dimension(:), intent(inout)  add2,
real (scrip_r8), dimension(:,:), intent(inout)  weights 
)

Definition at line 660 of file scrip_interface.F90.

660  !#######################################################################
661 
662  !-----------------------------------------------------------------------
663  !
664  ! this routine sorts address and weight arrays based on the
665  ! destination address with the source address as a secondary
666  ! sorting criterion. the method is a standard heap sort.
667  !
668  ! sort_add_v2 is identical to subroutine sort_add, but is moved into
669  ! scrip_interface.ftn and converted to f90 format (line continuations)
670  !
671  !-----------------------------------------------------------------------
672 
673  implicit none
674 
675  !-----------------------------------------------------------------------
676  !
677  ! Input and Output arrays
678  !
679  !-----------------------------------------------------------------------
680 
681  integer (SCRIP_i4), intent(inout), dimension(:) :: &
682  add1, & ! destination address array (num_links)
683  add2 ! source address array
684 
685  real (SCRIP_r8), intent(inout), dimension(:,:) :: &
686  weights ! remapping weights (num_wts, num_links)
687 
688  !-----------------------------------------------------------------------
689  !
690  ! local variables
691  !
692  !-----------------------------------------------------------------------
693 
694  integer (SCRIP_i4) :: &
695  num_links, & ! num of links for this mapping
696  num_wts, & ! num of weights for this mapping
697  add1_tmp, add2_tmp,& ! temp for addresses during swap
698  lvl, final_lvl, & ! level indexes for heap sort levels
699  chk_lvl1, chk_lvl2, max_lvl
700 
701  real (SCRIP_r8), dimension(SIZE(weights,DIM=1)) :: &
702  wgttmp ! temp for holding wts during swap
703 
704  !-----------------------------------------------------------------------
705  !
706  ! determine total number of links to sort and number of weights
707  !
708  !-----------------------------------------------------------------------
709 
710  num_links = SIZE(add1)
711  num_wts = SIZE(weights, dim=1)
712 
713  !-----------------------------------------------------------------------
714  !
715  ! start at the lowest level (N/2) of the tree and sift lower
716  ! values to the bottom of the tree, promoting the larger numbers
717  !
718  !-----------------------------------------------------------------------
719 
720  do lvl=num_links/2,1,-1
721 
722  final_lvl = lvl
723  add1_tmp = add1(lvl)
724  add2_tmp = add2(lvl)
725  wgttmp(:) = weights(:,lvl)
726 
727  !***
728  !*** loop until proper level is found for this link, or reach
729  !*** bottom
730  !***
731 
732  sift_loop1: do
733 
734  !***
735  !*** find the largest of the two daughters
736  !***
737 
738  chk_lvl1 = 2*final_lvl
739  chk_lvl2 = 2*final_lvl+1
740  if (chk_lvl1 .EQ. num_links) chk_lvl2 = chk_lvl1
741 
742  if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR. &
743  ((add1(chk_lvl1) == add1(chk_lvl2)) .AND. &
744  (add2(chk_lvl1) > add2(chk_lvl2)))) then
745  max_lvl = chk_lvl1
746  else
747  max_lvl = chk_lvl2
748  endif
749 
750  !***
751  !*** if the parent is greater than both daughters,
752  !*** the correct level has been found
753  !***
754 
755  if ((add1_tmp .GT. add1(max_lvl)) .OR. &
756  ((add1_tmp .EQ. add1(max_lvl)) .AND. &
757  (add2_tmp .GT. add2(max_lvl)))) then
758  add1(final_lvl) = add1_tmp
759  add2(final_lvl) = add2_tmp
760  weights(:,final_lvl) = wgttmp(:)
761  exit sift_loop1
762 
763  !***
764  !*** otherwise, promote the largest daughter and push
765  !*** down one level in the tree. if haven't reached
766  !*** the end of the tree, repeat the process. otherwise
767  !*** store last values and exit the loop
768  !***
769 
770  else
771  add1(final_lvl) = add1(max_lvl)
772  add2(final_lvl) = add2(max_lvl)
773  weights(:,final_lvl) = weights(:,max_lvl)
774 
775  final_lvl = max_lvl
776  if (2*final_lvl > num_links) then
777  add1(final_lvl) = add1_tmp
778  add2(final_lvl) = add2_tmp
779  weights(:,final_lvl) = wgttmp(:)
780  exit sift_loop1
781  endif
782  endif
783  end do sift_loop1
784  end do
785 
786  !-----------------------------------------------------------------------
787  !
788  ! now that the heap has been sorted, strip off the top (largest)
789  ! value and promote the values below
790  !
791  !-----------------------------------------------------------------------
792 
793  do lvl=num_links,3,-1
794 
795  !***
796  !*** move the top value and insert it into the correct place
797  !***
798 
799  add1_tmp = add1(lvl)
800  add1(lvl) = add1(1)
801 
802  add2_tmp = add2(lvl)
803  add2(lvl) = add2(1)
804 
805  wgttmp(:) = weights(:,lvl)
806  weights(:,lvl) = weights(:,1)
807 
808  !***
809  !*** as above this loop sifts the tmp values down until proper
810  !*** level is reached
811  !***
812 
813  final_lvl = 1
814 
815  sift_loop2: do
816 
817  !***
818  !*** find the largest of the two daughters
819  !***
820 
821  chk_lvl1 = 2*final_lvl
822  chk_lvl2 = 2*final_lvl+1
823  if (chk_lvl2 >= lvl) chk_lvl2 = chk_lvl1
824 
825  if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR. &
826  ((add1(chk_lvl1) == add1(chk_lvl2)) .AND. &
827  (add2(chk_lvl1) > add2(chk_lvl2)))) then
828  max_lvl = chk_lvl1
829  else
830  max_lvl = chk_lvl2
831  endif
832 
833  !***
834  !*** if the parent is greater than both daughters,
835  !*** the correct level has been found
836  !***
837 
838  if ((add1_tmp > add1(max_lvl)) .OR. &
839  ((add1_tmp == add1(max_lvl)) .AND. &
840  (add2_tmp > add2(max_lvl)))) then
841  add1(final_lvl) = add1_tmp
842  add2(final_lvl) = add2_tmp
843  weights(:,final_lvl) = wgttmp(:)
844  exit sift_loop2
845 
846  !***
847  !*** otherwise, promote the largest daughter and push
848  !*** down one level in the tree. if haven't reached
849  !*** the end of the tree, repeat the process. otherwise
850  !*** store last values and exit the loop
851  !***
852 
853  else
854  add1(final_lvl) = add1(max_lvl)
855  add2(final_lvl) = add2(max_lvl)
856  weights(:,final_lvl) = weights(:,max_lvl)
857 
858  final_lvl = max_lvl
859  if (2*final_lvl >= lvl) then
860  add1(final_lvl) = add1_tmp
861  add2(final_lvl) = add2_tmp
862  weights(:,final_lvl) = wgttmp(:)
863  exit sift_loop2
864  endif
865  endif
866  end do sift_loop2
867  end do
868 
869  !***
870  !*** swap the last two entries
871  !***
872 
873 
874  add1_tmp = add1(2)
875  add1(2) = add1(1)
876  add1(1) = add1_tmp
877 
878  add2_tmp = add2(2)
879  add2(2) = add2(1)
880  add2(1) = add2_tmp
881 
882  wgttmp(:) = weights(:,2)
883  weights(:,2) = weights(:,1)
884  weights(:,1) = wgttmp(:)
885 
886  !#######################################################################

Referenced by scrip().

Variable Documentation

◆ wgtdata

type(weight_data), dimension(:), allocatable scrip_interface::wgtdata

Definition at line 73 of file scrip_interface.F90.

73  type(weight_data), allocatable :: wgtdata(:)

Referenced by wmscrpmd::scrip_wrapper(), and wmgridmd::wmghgh().

scrip_remap_conservative::srch_corners_locate_point
integer(scrip_i4), save srch_corners_locate_point
Definition: scrip_remap_conservative.f:130
scrip_remap_conservative::avoid_pole_offset
real(scrip_r8), save avoid_pole_offset
Definition: scrip_remap_conservative.f:83
scrip_grids::grid2_centroid_lat
real(scrip_r8), dimension(:), allocatable, target, save grid2_centroid_lat
Definition: scrip_grids.f:103
scrip_remap_vars::norm_opt_none
integer(scrip_i4), parameter norm_opt_none
Definition: scrip_remap_vars.f:54
scrip_grids::grid1_corners
integer(scrip_i4), save grid1_corners
Definition: scrip_grids.f:68
scrip_remap_vars::grid1_add_map1
integer(scrip_i4), dimension(:), allocatable, save grid1_add_map1
Definition: scrip_remap_vars.f:77
scrip_grids::grid1_mask
logical(scrip_logical), dimension(:), allocatable, target, save grid1_mask
Definition: scrip_grids.f:93
scrip_remap_conservative::srch_corners_locate_segstart
integer(scrip_i4), save srch_corners_locate_segstart
Definition: scrip_remap_conservative.f:102
scrip_remap_conservative::srch_center_lat_locate_point
real(scrip_r8), dimension(:), allocatable, save srch_center_lat_locate_point
Definition: scrip_remap_conservative.f:144
scrip_remap_read
Definition: scrip_remap_read.f:40
scrip_remap_conservative::srch_corner_lon_locate_segstart
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lon_locate_segstart
Definition: scrip_remap_conservative.f:111
scrip_remap_conservative::srch_center_lon_locate_segstart
real(scrip_r8), dimension(:), allocatable, save srch_center_lon_locate_segstart
Definition: scrip_remap_conservative.f:117
scrip_remap_conservative::last_cell_grid_num_get_srch_cells
integer(scrip_i4), save last_cell_grid_num_get_srch_cells
Definition: scrip_remap_conservative.f:166
scrip_remap_vars::map_type_bilinear
integer(scrip_i4), parameter map_type_bilinear
Definition: scrip_remap_vars.f:59
scrip_grids::grid1_centroid_lat
real(scrip_r8), dimension(:), allocatable, target, save grid1_centroid_lat
Definition: scrip_grids.f:103
scrip_remap_vars::num_maps
integer(scrip_i4), save num_maps
Definition: scrip_remap_vars.f:66
scrip_remap_conservative::last_cell_locate_segstart
integer(scrip_i4), save last_cell_locate_segstart
Definition: scrip_remap_conservative.f:96
scrip_grids::grid2_units
character(scrip_charlength), save grid2_units
Definition: scrip_grids.f:80
scrip_grids::luse_grid2_area
logical(scrip_logical), save luse_grid2_area
Definition: scrip_grids.f:126
scrip_remap_conservative::first_call_locate_point
logical(scrip_logical), save first_call_locate_point
Definition: scrip_remap_conservative.f:121
scrip_remap_conservative::link_add2
integer(scrip_i4), dimension(:,:), allocatable, save link_add2
Definition: scrip_remap_conservative.f:86
scrip_grids::grid2_centroid_lon
real(scrip_r8), dimension(:), allocatable, target, save grid2_centroid_lon
Definition: scrip_grids.f:103
scrip_remap_conservative::srch_corner_lon_loc_get_srch_cells
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lon_loc_get_srch_cells
Definition: scrip_remap_conservative.f:158
scrip_grids::num_srch_bins
integer(scrip_i4), save num_srch_bins
Definition: scrip_grids.f:152
scrip_remap_conservative::last_srch_grid_num_locate_segstart
integer(scrip_i4), save last_srch_grid_num_locate_segstart
Definition: scrip_remap_conservative.f:96
scrip_remap_conservative::last_cell_grid_num_locate_point
integer(scrip_i4), save last_cell_grid_num_locate_point
Definition: scrip_remap_conservative.f:124
scrip_remap_vars::num_links_map1
integer(scrip_i4), save num_links_map1
Definition: scrip_remap_vars.f:66
scrip_timers::max_timers
integer(scrip_i4), parameter max_timers
Definition: scrip_timers.f:47
scrip_grids::grid2_rank
integer(scrip_i4), save grid2_rank
Definition: scrip_grids.f:68
scrip_remap_conservative::remap_conserv
subroutine remap_conserv(l_master, l_test)
Definition: scrip_remap_conservative.f:246
scrip_grids::grid2_mask
logical(scrip_logical), dimension(:), allocatable, target, save grid2_mask
Definition: scrip_grids.f:93
scrip_remap_vars::norm_opt
integer(scrip_i4), save norm_opt
Definition: scrip_remap_vars.f:66
scrip_grids::grid2_center_lat
real(scrip_r8), dimension(:), allocatable, target, save grid2_center_lat
Definition: scrip_grids.f:103
scrip_remap_conservative::srch_add_find_adj_cell
integer(scrip_i4), dimension(:), allocatable, save srch_add_find_adj_cell
Definition: scrip_remap_conservative.f:186
scrip_grids::bin_lats
real(scrip_r8), dimension(:,:), allocatable, save bin_lats
Definition: scrip_grids.f:159
scrip_remap_vars::max_links_map2
integer(scrip_i4), save max_links_map2
Definition: scrip_remap_vars.f:66
scrip_remap_conservative::last_cell_locate_point
integer(scrip_i4), save last_cell_locate_point
Definition: scrip_remap_conservative.f:124
scrip_remap_conservative::first_call_locate_segstart
logical(scrip_logical), save first_call_locate_segstart
Definition: scrip_remap_conservative.f:93
scrip_remap_conservative::srch_center_lon_locate_point
real(scrip_r8), dimension(:), allocatable, save srch_center_lon_locate_point
Definition: scrip_remap_conservative.f:144
scrip_grids::grid1_rank
integer(scrip_i4), save grid1_rank
Definition: scrip_grids.f:68
scrip_remap_vars::num_links_map2
integer(scrip_i4), save num_links_map2
Definition: scrip_remap_vars.f:66
scrip_grids::grid1_area
real(scrip_r8), dimension(:), allocatable, target, save grid1_area
Definition: scrip_grids.f:103
scrip_grids::north_thresh
real(scrip_r8), save north_thresh
Definition: scrip_grids.f:169
scrip_remap_vars::map_type
integer(scrip_i4), save map_type
Definition: scrip_remap_vars.f:66
scrip_grids::luse_grid_centers
logical(scrip_logical), save luse_grid_centers
Definition: scrip_grids.f:126
scrip_remap_conservative::srch_corner_lat_locate_point
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lat_locate_point
Definition: scrip_remap_conservative.f:138
scrip_grids
Definition: scrip_grids.f:49
scrip_constants::tiny
real(kind=scrip_r8), parameter tiny
Definition: scrip_constants.f:50
scrip_remap_conservative::srch_corner_lon_find_adj_cell
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lon_find_adj_cell
Definition: scrip_remap_conservative.f:188
scrip_grids::special_polar_cell2
logical(scrip_logical), dimension(:), allocatable, target, save special_polar_cell2
Definition: scrip_grids.f:93
scrip_grids::restrict_type
character(scrip_charlength), save restrict_type
Definition: scrip_grids.f:149
scrip_remap_conservative::last_srch_grid_num_get_srch_cells
integer(scrip_i4), save last_srch_grid_num_get_srch_cells
Definition: scrip_remap_conservative.f:166
scrip_grids::grid2_corners
integer(scrip_i4), save grid2_corners
Definition: scrip_grids.f:68
scrip_grids::grid1_area_in
real(scrip_r8), dimension(:), allocatable, target, save grid1_area_in
Definition: scrip_grids.f:103
scrip_remap_conservative::num_srch_cell_locate_points
integer(scrip_i4), save num_srch_cell_locate_points
Definition: scrip_remap_conservative.f:130
scrip_remap_vars::map_type_particle
integer(scrip_i4), parameter map_type_particle
Definition: scrip_remap_vars.f:59
scrip_errormod::scrip_errorset
subroutine, public scrip_errorset(errorCode, rtnName, errorMsg)
Definition: scrip_errormod.f90:81
scrip_remap_conservative::last_cell_find_adj_cell
integer(scrip_i4), save last_cell_find_adj_cell
Definition: scrip_remap_conservative.f:180
scrip_grids::grid1_centroid_lon
real(scrip_r8), dimension(:), allocatable, target, save grid1_centroid_lon
Definition: scrip_grids.f:103
scrip_grids::grid1_bound_box
real(scrip_r8), dimension(:,:), allocatable, target, save grid1_bound_box
Definition: scrip_grids.f:131
scrip_grids::bin_lons
real(scrip_r8), dimension(:,:), allocatable, save bin_lons
Definition: scrip_grids.f:159
scrip_remap_vars::grid2_add_map2
integer(scrip_i4), dimension(:), allocatable, save grid2_add_map2
Definition: scrip_remap_vars.f:77
scrip_grids::grid2_corner_lat
real(scrip_r8), dimension(:,:), allocatable, target, save grid2_corner_lat
Definition: scrip_grids.f:120
scrip_remap_vars::norm_opt_frcarea
integer(scrip_i4), parameter norm_opt_frcarea
Definition: scrip_remap_vars.f:54
scrip_grids::grid1_dims
integer(scrip_i4), dimension(:), allocatable, save grid1_dims
Definition: scrip_grids.f:74
scrip_remap_vars::max_links_map1
integer(scrip_i4), save max_links_map1
Definition: scrip_remap_vars.f:66
scrip_remap_conservative::num_srch_cells_locate_segstart
integer(scrip_i4), save num_srch_cells_locate_segstart
Definition: scrip_remap_conservative.f:102
scrip_remap_conservative::srch_center_lat_find_adj_cell
real(scrip_r8), dimension(:), allocatable, save srch_center_lat_find_adj_cell
Definition: scrip_remap_conservative.f:191
scrip_remap_write
Definition: scrip_remap_write.f:40
scrip_grids::luse_grid1_area
logical(scrip_logical), save luse_grid1_area
Definition: scrip_grids.f:126
scrip_grids::grid1_center_lon
real(scrip_r8), dimension(:), allocatable, target, save grid1_center_lon
Definition: scrip_grids.f:103
scrip_grids::grid2_area
real(scrip_r8), dimension(:), allocatable, target, save grid2_area
Definition: scrip_grids.f:103
scrip_timers::timer_clear
subroutine timer_clear(timer)
Definition: scrip_timers.f:105
scrip_errormod::scrip_success
integer(scrip_i4), parameter, public scrip_success
Definition: scrip_errormod.f90:42
scrip_grids::grid1_corner_lat
real(scrip_r8), dimension(:,:), allocatable, target, save grid1_corner_lat
Definition: scrip_grids.f:120
scrip_remap_conservative::srch_center_lon_loc_get_srch_cells
real(scrip_r8), dimension(:), allocatable, save srch_center_lon_loc_get_srch_cells
Definition: scrip_remap_conservative.f:162
scrip_remap_conservative::srch_corner_lat_find_adj_cell
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lat_find_adj_cell
Definition: scrip_remap_conservative.f:188
scrip_grids::grid2_size
integer(scrip_i4), save grid2_size
Definition: scrip_grids.f:68
scrip_grids::grid1_name
character(scrip_charlength), save grid1_name
Definition: scrip_grids.f:77
scrip_grids::bin_addr2
integer(scrip_i4), dimension(:,:), allocatable, save bin_addr2
Definition: scrip_grids.f:155
scrip_grids::grid2_name
character(scrip_charlength), save grid2_name
Definition: scrip_grids.f:77
scrip_remap_conservative::last_cell_grid_num_locate_segstart
integer(scrip_i4), save last_cell_grid_num_locate_segstart
Definition: scrip_remap_conservative.f:96
scrip_remap_conservative::nthreads
integer(scrip_i4) nthreads
Definition: scrip_remap_conservative.f:73
scrip_remap_conservative::last_cell_add_get_srch_cells
integer(scrip_i4), save last_cell_add_get_srch_cells
Definition: scrip_remap_conservative.f:166
scrip_remap_write::write_remap_ww3
subroutine write_remap_ww3(map1_name, interp_file1, output_opt, l_master, errorCode)
Definition: scrip_remap_write.f:953
scrip_grids::grid_init
subroutine grid_init(errorCode, l_master, l_test)
Definition: scrip_grids.f:188
scrip_remap_vars::wts_map1
real(scrip_r8), dimension(:,:), allocatable, save wts_map1
Definition: scrip_remap_vars.f:83
scrip_remap_conservative::link_add1
integer(scrip_i4), dimension(:,:), allocatable, save link_add1
Definition: scrip_remap_conservative.f:86
scrip_remap_conservative::srch_add_locate_point
integer(scrip_i4), dimension(:), allocatable, save srch_add_locate_point
Definition: scrip_remap_conservative.f:134
scrip_remap_vars::resize_remap_vars
subroutine resize_remap_vars(nmap, increment)
Definition: scrip_remap_vars.f:173
scrip_constants
Definition: scrip_constants.f:38
scrip_remap_conservative::last_cell_grid_num_find_adj_cell
integer(scrip_i4), save last_cell_grid_num_find_adj_cell
Definition: scrip_remap_conservative.f:180
scrip_grids::grid2_area_in
real(scrip_r8), dimension(:), allocatable, target, save grid2_area_in
Definition: scrip_grids.f:103
scrip_remap_vars::map_type_distwgt
integer(scrip_i4), parameter map_type_distwgt
Definition: scrip_remap_vars.f:59
scrip_remap_conservative
Definition: scrip_remap_conservative.f:60
scrip_remap_conservative::num_srch_cells_loc_get_srch_cells
integer(scrip_i4), save num_srch_cells_loc_get_srch_cells
Definition: scrip_remap_conservative.f:148
scrip_grids::grid2_dims
integer(scrip_i4), dimension(:), allocatable, save grid2_dims
Definition: scrip_grids.f:74
scrip_remap_read::read_remap_ww3
subroutine read_remap_ww3(map_name, interp_file, errorCode)
Definition: scrip_remap_read.f:115
scrip_grids::grid1_center_lat
real(scrip_r8), dimension(:), allocatable, target, save grid1_center_lat
Definition: scrip_grids.f:103
scrip_remap_vars::init_remap_vars
subroutine init_remap_vars
Definition: scrip_remap_vars.f:97
scrip_remap_conservative::srch_corner_lat_loc_get_srch_cells
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lat_loc_get_srch_cells
Definition: scrip_remap_conservative.f:158
scrip_kindsmod
Definition: scrip_kindsmod.f90:3
scrip_remap_conservative::num_srch_cells_find_adj_cell
integer(scrip_i4), save num_srch_cells_find_adj_cell
Definition: scrip_remap_conservative.f:180
scrip_grids::grid1_corner_lon
real(scrip_r8), dimension(:,:), allocatable, target, save grid1_corner_lon
Definition: scrip_grids.f:120
scrip_remap_conservative::srch_corner_lon_locate_point
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lon_locate_point
Definition: scrip_remap_conservative.f:138
scrip_grids::npseg
integer(scrip_i4), save npseg
Definition: scrip_grids.f:178
scrip_remap_vars::map_type_conserv
integer(scrip_i4), parameter map_type_conserv
Definition: scrip_remap_vars.f:59
scrip_errormod
Definition: scrip_errormod.f90:3
scrip_remap_vars::grid1_add_map2
integer(scrip_i4), dimension(:), allocatable, save grid1_add_map2
Definition: scrip_remap_vars.f:77
scrip_remap_vars::grid2_add_map1
integer(scrip_i4), dimension(:), allocatable, save grid2_add_map1
Definition: scrip_remap_vars.f:77
scrip_grids::south_thresh
real(scrip_r8), save south_thresh
Definition: scrip_grids.f:169
scrip_remap_conservative::srch_corner_lat_locate_segstart
real(scrip_r8), dimension(:,:), allocatable, save srch_corner_lat_locate_segstart
Definition: scrip_remap_conservative.f:111
scrip_grids::grid1_size
integer(scrip_i4), save grid1_size
Definition: scrip_grids.f:68
scrip_remap_conservative::first_call_find_adj_cell
logical(scrip_logical), save first_call_find_adj_cell
Definition: scrip_remap_conservative.f:174
scrip_remap_conservative::srch_center_lon_find_adj_cell
real(scrip_r8), dimension(:), allocatable, save srch_center_lon_find_adj_cell
Definition: scrip_remap_conservative.f:191
scrip_remap_vars::map_type_bicubic
integer(scrip_i4), parameter map_type_bicubic
Definition: scrip_remap_vars.f:59
scrip_remap_vars::wts_map2
real(scrip_r8), dimension(:,:), allocatable, save wts_map2
Definition: scrip_remap_vars.f:83
scrip_grids::grid2_bound_box
real(scrip_r8), dimension(:,:), allocatable, target, save grid2_bound_box
Definition: scrip_grids.f:131
scrip_timers
Definition: scrip_timers.f:39
scrip_grids::grid2_corner_lon
real(scrip_r8), dimension(:,:), allocatable, target, save grid2_corner_lon
Definition: scrip_grids.f:120
scrip_grids::grid1_frac
real(scrip_r8), dimension(:), allocatable, target, save grid1_frac
Definition: scrip_grids.f:103
scrip_remap_conservative::srch_corners_find_adj_cell
integer(scrip_i4), save srch_corners_find_adj_cell
Definition: scrip_remap_conservative.f:180
scrip_remap_vars
Definition: scrip_remap_vars.f:40
scrip_remap_conservative::srch_center_lat_loc_get_srch_cells
real(scrip_r8), dimension(:), allocatable, save srch_center_lat_loc_get_srch_cells
Definition: scrip_remap_conservative.f:162
scrip_remap_write::write_remap
subroutine write_remap(map1_name, map2_name, interp_file1, interp_file2, output_opt, l_master, errorCode)
Definition: scrip_remap_write.f:122
scrip_remap_conservative::srch_corners_loc_get_srch_cells
integer(scrip_i4), save srch_corners_loc_get_srch_cells
Definition: scrip_remap_conservative.f:148
scrip_remap_conservative::srch_add_loc_get_srch_cells
integer(scrip_i4), dimension(:), allocatable, save srch_add_loc_get_srch_cells
Definition: scrip_remap_conservative.f:154
scrip_grids::grid1_units
character(scrip_charlength), save grid1_units
Definition: scrip_grids.f:80
scrip_remap_vars::norm_opt_dstarea
integer(scrip_i4), parameter norm_opt_dstarea
Definition: scrip_remap_vars.f:54
scrip_remap_conservative::srch_add_locate_segstart
integer(scrip_i4), dimension(:), allocatable, save srch_add_locate_segstart
Definition: scrip_remap_conservative.f:107
scrip_remap_conservative::avoid_pole_count
integer(scrip_i4), save avoid_pole_count
Definition: scrip_remap_conservative.f:80
scrip_iounitsmod
Definition: scrip_iounitsmod.f90:3
scrip_remap_conservative::last_srch_grid_num_locate_point
integer(scrip_i4), save last_srch_grid_num_locate_point
Definition: scrip_remap_conservative.f:124
scrip_iounitsmod::scrip_stdout
integer(scrip_i4), parameter, public scrip_stdout
Definition: scrip_iounitsmod.f90:47
scrip_remap_conservative::first_call_store_link_cnsrv
logical(scrip_logical), save first_call_store_link_cnsrv
Definition: scrip_remap_conservative.f:90
scrip_grids::special_polar_cell1
logical(scrip_logical), dimension(:), allocatable, target, save special_polar_cell1
Definition: scrip_grids.f:93
scrip_remap_vars::num_wts
integer(scrip_i4), save num_wts
Definition: scrip_remap_vars.f:66
scrip_remap_conservative::first_call_get_srch_cells
logical(scrip_logical), save first_call_get_srch_cells
Definition: scrip_remap_conservative.f:171
scrip_grids::grid2_frac
real(scrip_r8), dimension(:), allocatable, target, save grid2_frac
Definition: scrip_grids.f:103
scrip_grids::bin_addr1
integer(scrip_i4), dimension(:,:), allocatable, save bin_addr1
Definition: scrip_grids.f:155
scrip_grids::grid2_center_lon
real(scrip_r8), dimension(:), allocatable, target, save grid2_center_lon
Definition: scrip_grids.f:103
scrip_timers::timers_init
subroutine timers_init
Definition: scrip_timers.f:303
scrip_remap_conservative::srch_center_lat_locate_segstart
real(scrip_r8), dimension(:), allocatable, save srch_center_lat_locate_segstart
Definition: scrip_remap_conservative.f:117
scrip_remap_vars::resize_increment
integer(scrip_i4), save resize_increment
Definition: scrip_remap_vars.f:66