WAVEWATCH III  beta 0.0.1
scrip_remap_write.f
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 !
3 ! This module contains routines for writing the remapping data to
4 ! a file. Before writing the data for each mapping, the links are
5 ! sorted by destination grid address.
6 !
7 !-----------------------------------------------------------------------
8 !
9 ! CVS:$Id: remap_write.f,v 1.7 2001/08/21 21:06:42 pwjones Exp $
10 !
11 ! Copyright (c) 1997, 1998 the Regents of the University of
12 ! California.
13 !
14 ! This software and ancillary information (herein called software)
15 ! called SCRIP is made available under the terms described here.
16 ! The software has been approved for release with associated
17 ! LA-CC Number 98-45.
18 !
19 ! Unless otherwise indicated, this software has been authored
20 ! by an employee or employees of the University of California,
21 ! operator of the Los Alamos National Laboratory under Contract
22 ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
23 ! Government has rights to use, reproduce, and distribute this
24 ! software. The public may copy and use this software without
25 ! charge, provided that this Notice and any statement of authorship
26 ! are reproduced on all copies. Neither the Government nor the
27 ! University makes any warranty, express or implied, or assumes
28 ! any liability or responsibility for the use of this software.
29 !
30 ! If software is modified to produce derivative works, such modified
31 ! software should be clearly marked, so as not to confuse it with
32 ! the version available from Los Alamos National Laboratory.
33 !
34 ! This code has been modified from the version available from
35 ! Los Alamos National Laboratory, for the purpose of running it
36 ! within WW3.
37 !
38 !***********************************************************************
39 
41 
42 !-----------------------------------------------------------------------
43 
44  use scrip_kindsmod ! defines common data types
45  use scrip_errormod ! SCRIP error handler
46  use scrip_netcdfmod ! netCDF error handler
47  use netcdf ! netCDF library
48 
49  use scrip_constants ! defines common scalar constants
50  use scrip_grids ! module containing grid information
51  use scrip_remap_vars ! module containing remap information
52 
53  implicit none
54 
55 !-----------------------------------------------------------------------
56 !
57 ! module variables
58 !
59 !-----------------------------------------------------------------------
60 
61  character(SCRIP_charLength), private ::
62  & map_method ! character string for map_type
63  &, normalize_opt ! character string for normalization option
64  &, history ! character string for history information
65  &, convention ! character string for output convention
66 
67  character(8), private ::
68  & cdate ! character date string
69 
70  integer (SCRIP_i4), dimension(:), allocatable, private ::
71  & src_mask_int ! integer masks to determine
72  &, dst_mask_int ! cells that participate in map
73 
74 !-----------------------------------------------------------------------
75 !
76 ! various netCDF identifiers used by output routines
77 !
78 !-----------------------------------------------------------------------
79 
80  integer (SCRIP_i4), private ::
81  & ncstat ! error flag for netCDF calls
82  &, nc_file_id ! id for netCDF file
83  &, nc_srcgrdsize_id ! id for source grid size
84  &, nc_dstgrdsize_id ! id for destination grid size
85  &, nc_srcgrdcorn_id ! id for number of source grid corners
86  &, nc_dstgrdcorn_id ! id for number of dest grid corners
87  &, nc_srcgrdrank_id ! id for source grid rank
88  &, nc_dstgrdrank_id ! id for dest grid rank
89  &, nc_numlinks_id ! id for number of links in mapping
90  &, nc_numwgts_id ! id for number of weights for mapping
91  &, nc_srcgrddims_id ! id for source grid dimensions
92  &, nc_dstgrddims_id ! id for dest grid dimensions
93  &, nc_srcgrdcntrlat_id ! id for source grid center latitude
94  &, nc_dstgrdcntrlat_id ! id for dest grid center latitude
95  &, nc_srcgrdcntrlon_id ! id for source grid center longitude
96  &, nc_dstgrdcntrlon_id ! id for dest grid center longitude
97  &, nc_srcgrdimask_id ! id for source grid mask
98  &, nc_dstgrdimask_id ! id for dest grid mask
99  &, nc_srcgrdcrnrlat_id ! id for latitude of source grid corners
100  &, nc_srcgrdcrnrlon_id ! id for longitude of source grid corners
101  &, nc_dstgrdcrnrlat_id ! id for latitude of dest grid corners
102  &, nc_dstgrdcrnrlon_id ! id for longitude of dest grid corners
103  &, nc_srcgrdarea_id ! id for area of source grid cells
104  &, nc_dstgrdarea_id ! id for area of dest grid cells
105  &, nc_srcgrdfrac_id ! id for area fraction on source grid
106  &, nc_dstgrdfrac_id ! id for area fraction on dest grid
107  &, nc_srcadd_id ! id for map source address
108  &, nc_dstadd_id ! id for map destination address
109  &, nc_rmpmatrix_id ! id for remapping matrix
110 
111  integer (SCRIP_i4), dimension(2), private ::
112  & nc_dims2_id ! netCDF ids for 2d array dims
113 
114 !***********************************************************************
115 
116  contains
117 
118 !***********************************************************************
119 
120  subroutine write_remap(map1_name, map2_name, interp_file1,
121  & interp_file2, output_opt, l_master, errorCode)
122 
123 !-----------------------------------------------------------------------
124 !
125 ! calls correct output routine based on output format choice
126 !
127 !-----------------------------------------------------------------------
128 
129 !-----------------------------------------------------------------------
130 !
131 ! input variables
132 !
133 !-----------------------------------------------------------------------
134 
135  character(SCRIP_charLength), intent(in) ::
136  & map1_name, ! name for mapping grid1 to grid2
137  & map2_name, ! name for mapping grid2 to grid1
138  & interp_file1, ! filename for map1 remap data
139  & interp_file2, ! filename for map2 remap data
140  & output_opt ! option for output conventions
141 
142  logical, intent(in) ::
143  & l_master ! Am I the master processor?
144 
145 !-----------------------------------------------------------------------
146 !
147 ! output variables
148 !
149 !-----------------------------------------------------------------------
150 
151  integer (SCRIP_i4), intent(out) ::
152  & errorcode ! returned error code
153 
154 !-----------------------------------------------------------------------
155 !
156 ! local variables
157 !
158 !-----------------------------------------------------------------------
159 
160  character (11), parameter :: rtnName = 'write_remap'
161 
162 !-----------------------------------------------------------------------
163 !
164 ! define some common variables to be used in all routines
165 !
166 !-----------------------------------------------------------------------
167 
168  errorcode = scrip_success
169 
170  select case(norm_opt)
171  case (norm_opt_none)
172  normalize_opt = 'none'
173  case (norm_opt_frcarea)
174  normalize_opt = 'fracarea'
175  case (norm_opt_dstarea)
176  normalize_opt = 'destarea'
177  end select
178 
179  select case(map_type)
180  case(map_type_conserv)
181  map_method = 'Conservative remapping'
182  case(map_type_bilinear)
183  map_method = 'Bilinear remapping'
184  case(map_type_distwgt)
185  map_method = 'Distance weighted avg of nearest neighbors'
186  case(map_type_bicubic)
187  map_method = 'Bicubic remapping'
188  case(map_type_particle)
189  map_method = 'Particle remapping'
190  case default
191  call scrip_errorset(errorcode, rtnname, 'Invalid Map Type')
192  return
193  end select
194 
195  call date_and_time(date=cdate)
196  write (history,1000) cdate(5:6),cdate(7:8),cdate(1:4)
197  1000 format('Created: ',a2,'-',a2,'-',a4)
198 
199 !-----------------------------------------------------------------------
200 !
201 ! sort address and weight arrays
202 !
203 !-----------------------------------------------------------------------
204 
205 ! New Apr4 2013: sort_add is instead called from scrip_interface.ftn
206 ! prior to entering this routine...no need to call it again here.
207 
208 ! call sort_add(grid2_add_map1, grid1_add_map1, wts_map1)
209  if (num_maps > 1) then
211  endif
212 
213 !-----------------------------------------------------------------------
214 !
215 ! call appropriate output routine
216 !
217 !-----------------------------------------------------------------------
218 
219  select case(output_opt)
220  case ('scrip')
221  if (l_master)
222  & call write_remap_scrip(map1_name, interp_file1, 1, errorcode)
223  if (scrip_errorcheck(errorcode, rtnname,
224  & 'error in write_remap_scrip')) return
225  case ('ncar-csm')
226  call write_remap_csm (map1_name, interp_file1, 1, errorcode)
227  if (scrip_errorcheck(errorcode, rtnname,
228  & 'error in write_remap_csm')) return
229  case default
230  call scrip_errorset(errorcode, rtnname,
231  & 'unknown output file convention')
232  return
233  end select
234 
235 !-----------------------------------------------------------------------
236 !
237 ! call appropriate output routine for second mapping if required
238 !
239 !-----------------------------------------------------------------------
240 
241  if (num_maps > 1) then
242  select case(output_opt)
243  case ('scrip')
244  if (l_master)
245  & call write_remap_scrip(map2_name, interp_file2, 2, errorcode)
246  if (scrip_errorcheck(errorcode, rtnname,
247  & 'error in write_remap_scrip')) return
248  case ('ncar-csm')
249  call write_remap_csm (map2_name, interp_file2, 2, errorcode)
250  if (scrip_errorcheck(errorcode, rtnname,
251  & 'error in write_remap_csm')) return
252  case default
253  call scrip_errorset(errorcode, rtnname,
254  & 'unknown output file convention')
255  return
256  end select
257  endif
258 
259 !-----------------------------------------------------------------------
260 
261  end subroutine write_remap
262 
263 !***********************************************************************
264 
265  subroutine write_remap_scrip(map_name, interp_file, direction,
266  & errorCode)
267 
268 !-----------------------------------------------------------------------
269 !
270 ! writes remap data to a netCDF file using SCRIP conventions
271 !
272 !-----------------------------------------------------------------------
273 
274 !-----------------------------------------------------------------------
275 !
276 ! input variables
277 !
278 !-----------------------------------------------------------------------
279 
280  character(SCRIP_charLength), intent(in) ::
281  & map_name ! name for mapping
282  &, interp_file ! filename for remap data
283 
284  integer (SCRIP_i4), intent(in) ::
285  & direction ! direction of map (1=grid1 to grid2
286  ! 2=grid2 to grid1)
287 
288 !-----------------------------------------------------------------------
289 !
290 ! output variables
291 !
292 !-----------------------------------------------------------------------
293 
294  integer (SCRIP_i4), intent(out) ::
295  & errorCode ! returned error code
296 
297 !-----------------------------------------------------------------------
298 !
299 ! local variables
300 !
301 !-----------------------------------------------------------------------
302 
303  character(SCRIP_charLength) ::
304  & grid1_ctmp ! character temp for grid1 names
305  &, grid2_ctmp ! character temp for grid2 names
306 
307  integer (SCRIP_i4) ::
308  & itmp1 ! integer temp
309  &, itmp2 ! integer temp
310  &, itmp3 ! integer temp
311  &, itmp4 ! integer temp
312 
313  character (17), parameter :: rtnName = 'write_remap_scrip'
314 
315 !-----------------------------------------------------------------------
316 !
317 ! create netCDF file for mapping and define some global attributes
318 !
319 !-----------------------------------------------------------------------
320 
321  errorcode = scrip_success
322 
323  ncstat = nf90_create(interp_file, nf90_clobber, nc_file_id)
324  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
325  & 'error creating remap file')) return
326 
327  !***
328  !*** map name
329  !***
330  ncstat = nf90_put_att(nc_file_id, nf90_global, 'title', map_name)
331  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
332  & 'error writing remap name')) return
333 
334  !***
335  !*** normalization option
336  !***
337  ncstat = nf90_put_att(nc_file_id, nf90_global, 'normalization',
338  & normalize_opt)
339  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
340  & 'error writing normalize option')) return
341 
342  !***
343  !*** map method
344  !***
345  ncstat = nf90_put_att(nc_file_id, nf90_global, 'map_method',
346  & map_method)
347  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
348  & 'error writing remap method')) return
349 
350  !***
351  !*** history
352  !***
353  ncstat = nf90_put_att(nc_file_id, nf90_global, 'history',
354  & history)
355  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
356  & 'error writing history')) return
357 
358  !***
359  !*** file convention
360  !***
361  convention = 'SCRIP'
362  ncstat = nf90_put_att(nc_file_id, nf90_global, 'conventions',
363  & convention)
364  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
365  & 'error writing output convention')) return
366 
367  !***
368  !*** source and destination grid names
369  !***
370 
371  if (direction == 1) then
372  grid1_ctmp = 'source_grid'
373  grid2_ctmp = 'dest_grid'
374  else
375  grid1_ctmp = 'dest_grid'
376  grid2_ctmp = 'source_grid'
377  endif
378 
379  ncstat = nf90_put_att(nc_file_id, nf90_global, trim(grid1_ctmp),
380  & grid1_name)
381  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
382  & 'error writing source grid name')) return
383 
384  ncstat = nf90_put_att(nc_file_id, nf90_global, trim(grid2_ctmp),
385  & grid2_name)
386  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
387  & 'error writing destination grid name')) return
388 
389 !-----------------------------------------------------------------------
390 !
391 ! prepare netCDF dimension info
392 !
393 !-----------------------------------------------------------------------
394 
395  !***
396  !*** define grid size dimensions
397  !***
398 
399  if (direction == 1) then
400  itmp1 = grid1_size
401  itmp2 = grid2_size
402  else
403  itmp1 = grid2_size
404  itmp2 = grid1_size
405  endif
406 
407  ncstat = nf90_def_dim(nc_file_id, 'src_grid_size', itmp1,
408  & nc_srcgrdsize_id)
409  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
410  & 'error defining source grid size')) return
411 
412  ncstat = nf90_def_dim(nc_file_id, 'dst_grid_size', itmp2,
413  & nc_dstgrdsize_id)
414  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
415  & 'error defining destination grid size')) return
416 
417  !***
418  !*** define grid corner dimension
419  !***
420 
421  if (direction == 1) then
422  itmp1 = grid1_corners
423  itmp2 = grid2_corners
424  else
425  itmp1 = grid2_corners
426  itmp2 = grid1_corners
427  endif
428 
429  ncstat = nf90_def_dim(nc_file_id, 'src_grid_corners',
430  & itmp1, nc_srcgrdcorn_id)
431  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
432  & 'error defining num corners on source grid')) return
433 
434  ncstat = nf90_def_dim(nc_file_id, 'dst_grid_corners',
435  & itmp2, nc_dstgrdcorn_id)
436  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
437  & 'error defining num corners on destination grid')) return
438 
439  !***
440  !*** define grid rank dimension
441  !***
442 
443  if (direction == 1) then
444  itmp1 = grid1_rank
445  itmp2 = grid2_rank
446  else
447  itmp1 = grid2_rank
448  itmp2 = grid1_rank
449  endif
450 
451  ncstat = nf90_def_dim(nc_file_id, 'src_grid_rank',
452  & itmp1, nc_srcgrdrank_id)
453  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
454  & 'error defining source grid rank')) return
455 
456  ncstat = nf90_def_dim(nc_file_id, 'dst_grid_rank',
457  & itmp2, nc_dstgrdrank_id)
458  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
459  & 'error defining destination grid rank')) return
460 
461  !***
462  !*** define map size dimensions
463  !***
464 
465  if (direction == 1) then
466  itmp1 = num_links_map1
467  else
468  itmp1 = num_links_map2
469  endif
470 
471  ncstat = nf90_def_dim(nc_file_id, 'num_links',
472  & itmp1, nc_numlinks_id)
473  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
474  & 'error defining remap size')) return
475 
476  ncstat = nf90_def_dim(nc_file_id, 'num_wgts',
477  & num_wts, nc_numwgts_id)
478  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
479  & 'error defining number of weights')) return
480 
481  !***
482  !*** define grid dimensions
483  !***
484 
485  ncstat = nf90_def_var(nc_file_id, 'src_grid_dims', nf90_int,
486  & nc_srcgrdrank_id, nc_srcgrddims_id)
487  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
488  & 'error defining source grid dims')) return
489 
490  ncstat = nf90_def_var(nc_file_id, 'dst_grid_dims', nf90_int,
491  & nc_dstgrdrank_id, nc_dstgrddims_id)
492  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
493  & 'error defining destination grid dims')) return
494 
495 !-----------------------------------------------------------------------
496 !
497 ! define all arrays for netCDF descriptors
498 !
499 !-----------------------------------------------------------------------
500 
501  !***
502  !*** define grid center latitude array
503  !***
504 
505  ncstat = nf90_def_var(nc_file_id, 'src_grid_center_lat',
506  & nf90_double, nc_srcgrdsize_id,
507  & nc_srcgrdcntrlat_id)
508  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
509  & 'error defining source grid center lat')) return
510 
511  ncstat = nf90_def_var(nc_file_id, 'dst_grid_center_lat',
512  & nf90_double, nc_dstgrdsize_id,
513  & nc_dstgrdcntrlat_id)
514  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
515  & 'error defining destination grid center lat')) return
516 
517  !***
518  !*** define grid center longitude array
519  !***
520 
521  ncstat = nf90_def_var(nc_file_id, 'src_grid_center_lon',
522  & nf90_double, nc_srcgrdsize_id,
523  & nc_srcgrdcntrlon_id)
524  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
525  & 'error defining source grid center lon')) return
526 
527  ncstat = nf90_def_var(nc_file_id, 'dst_grid_center_lon',
528  & nf90_double, nc_dstgrdsize_id,
529  & nc_dstgrdcntrlon_id)
530  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
531  & 'error defining destination grid center lon')) return
532 
533  !***
534  !*** define grid corner lat/lon arrays
535  !***
536 
537  nc_dims2_id(1) = nc_srcgrdcorn_id
538  nc_dims2_id(2) = nc_srcgrdsize_id
539 
540  ncstat = nf90_def_var(nc_file_id, 'src_grid_corner_lat',
541  & nf90_double, nc_dims2_id,
542  & nc_srcgrdcrnrlat_id)
543  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
544  & 'error defining source grid corner lats')) return
545 
546  ncstat = nf90_def_var(nc_file_id, 'src_grid_corner_lon',
547  & nf90_double, nc_dims2_id,
548  & nc_srcgrdcrnrlon_id)
549  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
550  & 'error defining source grid corner lons')) return
551 
552  nc_dims2_id(1) = nc_dstgrdcorn_id
553  nc_dims2_id(2) = nc_dstgrdsize_id
554 
555  ncstat = nf90_def_var(nc_file_id, 'dst_grid_corner_lat',
556  & nf90_double, nc_dims2_id,
557  & nc_dstgrdcrnrlat_id)
558  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
559  & 'error defining destination grid corner lats')) return
560 
561  ncstat = nf90_def_var(nc_file_id, 'dst_grid_corner_lon',
562  & nf90_double, nc_dims2_id,
563  & nc_dstgrdcrnrlon_id)
564  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
565  & 'error defining destination grid corner lons')) return
566 
567  !***
568  !*** define units for all coordinate arrays
569  !***
570 
571  if (direction == 1) then
572  grid1_ctmp = grid1_units
573  grid2_ctmp = grid2_units
574  else
575  grid1_ctmp = grid2_units
576  grid2_ctmp = grid1_units
577  endif
578 
579  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcntrlat_id,
580  & 'units', grid1_ctmp)
581  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
582  & 'error writing source grid units')) return
583 
584  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcntrlat_id,
585  & 'units', grid2_ctmp)
586  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
587  & 'error writing destination grid units')) return
588 
589  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcntrlon_id,
590  & 'units', grid1_ctmp)
591  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
592  & 'error writing source grid units')) return
593 
594  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcntrlon_id,
595  & 'units', grid2_ctmp)
596  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
597  & 'error writing destination grid units')) return
598 
599  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcrnrlat_id,
600  & 'units', grid1_ctmp)
601  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
602  & 'error writing source grid units')) return
603 
604  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcrnrlon_id,
605  & 'units', grid1_ctmp)
606  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
607  & 'error writing source grid units')) return
608 
609  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcrnrlat_id,
610  & 'units', grid2_ctmp)
611  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
612  & 'error writing destination grid units')) return
613 
614  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcrnrlon_id,
615  & 'units', grid2_ctmp)
616  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
617  & 'error writing destination grid units')) return
618 
619  !***
620  !*** define grid mask
621  !***
622 
623  ncstat = nf90_def_var(nc_file_id, 'src_grid_imask', nf90_int,
624  & nc_srcgrdsize_id, nc_srcgrdimask_id)
625  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
626  & 'error defining source grid mask')) return
627 
628  ncstat = nf90_put_att(nc_file_id, nc_srcgrdimask_id,
629  & 'units', 'unitless')
630  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
631  & 'error writing source grid mask units')) return
632 
633  ncstat = nf90_def_var(nc_file_id, 'dst_grid_imask', nf90_int,
634  & nc_dstgrdsize_id, nc_dstgrdimask_id)
635  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
636  & 'error defining destination grid mask')) return
637 
638  ncstat = nf90_put_att(nc_file_id, nc_dstgrdimask_id,
639  & 'units', 'unitless')
640  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
641  & 'error writing destination grid mask units')) return
642 
643  !***
644  !*** define grid area arrays
645  !***
646 
647  ncstat = nf90_def_var(nc_file_id, 'src_grid_area',
648  & nf90_double, nc_srcgrdsize_id,
649  & nc_srcgrdarea_id)
650  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
651  & 'error defining source grid area')) return
652 
653  ncstat = nf90_put_att(nc_file_id, nc_srcgrdarea_id,
654  & 'units', 'square radians')
655  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
656  & 'error writing source area units')) return
657 
658  ncstat = nf90_def_var(nc_file_id, 'dst_grid_area',
659  & nf90_double, nc_dstgrdsize_id,
660  & nc_dstgrdarea_id)
661  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
662  & 'error defining destination grid area')) return
663 
664  ncstat = nf90_put_att(nc_file_id, nc_dstgrdarea_id,
665  & 'units', 'square radians')
666  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
667  & 'error writing destination area units')) return
668 
669  !***
670  !*** define grid fraction arrays
671  !***
672 
673  ncstat = nf90_def_var(nc_file_id, 'src_grid_frac',
674  & nf90_double, nc_srcgrdsize_id,
675  & nc_srcgrdfrac_id)
676  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
677  & 'error defining source grid fraction')) return
678 
679  ncstat = nf90_put_att(nc_file_id, nc_srcgrdfrac_id,
680  & 'units', 'unitless')
681  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
682  & 'error writing source fraction units')) return
683 
684  ncstat = nf90_def_var(nc_file_id, 'dst_grid_frac',
685  & nf90_double, nc_dstgrdsize_id,
686  & nc_dstgrdfrac_id)
687  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
688  & 'error defining destination fraction')) return
689 
690  ncstat = nf90_put_att(nc_file_id, nc_dstgrdfrac_id,
691  & 'units', 'unitless')
692  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
693  & 'error writing destination frac units')) return
694 
695  !***
696  !*** define mapping arrays
697  !***
698 
699  ncstat = nf90_def_var(nc_file_id, 'src_address',
700  & nf90_int, nc_numlinks_id,
701  & nc_srcadd_id)
702  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
703  & 'error defining source addresses')) return
704 
705  ncstat = nf90_def_var(nc_file_id, 'dst_address',
706  & nf90_int, nc_numlinks_id,
707  & nc_dstadd_id)
708  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
709  & 'error defining destination addresses')) return
710 
711  nc_dims2_id(1) = nc_numwgts_id
712  nc_dims2_id(2) = nc_numlinks_id
713 
714  ncstat = nf90_def_var(nc_file_id, 'remap_matrix',
715  & nf90_double, nc_dims2_id,
716  & nc_rmpmatrix_id)
717  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
718  & 'error defining remapping weights')) return
719 
720  !***
721  !*** end definition stage
722  !***
723 
724  ncstat = nf90_enddef(nc_file_id)
725  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
726  & 'error ending definition phase')) return
727 
728 !-----------------------------------------------------------------------
729 !
730 ! compute integer masks
731 !
732 !-----------------------------------------------------------------------
733 
734  if (direction == 1) then
735  allocate (src_mask_int(grid1_size),
736  & dst_mask_int(grid2_size))
737 
738  where (grid2_mask)
739  dst_mask_int = 1
740  elsewhere
741  dst_mask_int = 0
742  endwhere
743 
744  where (grid1_mask)
745  src_mask_int = 1
746  elsewhere
747  src_mask_int = 0
748  endwhere
749  else
750  allocate (src_mask_int(grid2_size),
751  & dst_mask_int(grid1_size))
752 
753  where (grid1_mask)
754  dst_mask_int = 1
755  elsewhere
756  dst_mask_int = 0
757  endwhere
758 
759  where (grid2_mask)
760  src_mask_int = 1
761  elsewhere
762  src_mask_int = 0
763  endwhere
764  endif
765 
766 !-----------------------------------------------------------------------
767 !
768 ! change units of lat/lon coordinates if input units different
769 ! from radians
770 !
771 !-----------------------------------------------------------------------
772 
773  if (grid1_units(1:7) == 'degrees' .and. direction == 1) then
778  endif
779 
780  if (grid2_units(1:7) == 'degrees' .and. direction == 1) then
785  endif
786 
787 !-----------------------------------------------------------------------
788 !
789 ! write mapping data
790 !
791 !-----------------------------------------------------------------------
792 
793  if (direction == 1) then
794  itmp1 = nc_srcgrddims_id
795  itmp2 = nc_dstgrddims_id
796  else
797  itmp2 = nc_srcgrddims_id
798  itmp1 = nc_dstgrddims_id
799  endif
800 
801  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_dims)
802  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
803  & 'error writing source grid dims')) return
804 
805  ncstat = nf90_put_var(nc_file_id, itmp2, grid2_dims)
806  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
807  & 'error writing destination grid dims')) return
808 
809  ncstat = nf90_put_var(nc_file_id, nc_srcgrdimask_id, src_mask_int)
810  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
811  & 'error writing source grid mask')) return
812 
813  ncstat = nf90_put_var(nc_file_id, nc_dstgrdimask_id, dst_mask_int)
814  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
815  & 'error writing destination grid mask')) return
816 
817  deallocate(src_mask_int, dst_mask_int)
818 
819  if (direction == 1) then
820  itmp1 = nc_srcgrdcntrlat_id
821  itmp2 = nc_srcgrdcntrlon_id
822  itmp3 = nc_srcgrdcrnrlat_id
823  itmp4 = nc_srcgrdcrnrlon_id
824  else
825  itmp1 = nc_dstgrdcntrlat_id
826  itmp2 = nc_dstgrdcntrlon_id
827  itmp3 = nc_dstgrdcrnrlat_id
828  itmp4 = nc_dstgrdcrnrlon_id
829  endif
830 
831  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_center_lat)
832  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
833  & 'error writing source grid center lats')) return
834 
835  ncstat = nf90_put_var(nc_file_id, itmp2, grid1_center_lon)
836  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
837  & 'error writing source grid center lons')) return
838 
839  ncstat = nf90_put_var(nc_file_id, itmp3, grid1_corner_lat)
840  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
841  & 'error writing source grid corner lats')) return
842 
843  ncstat = nf90_put_var(nc_file_id, itmp4, grid1_corner_lon)
844  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
845  & 'error writing source grid corner lons')) return
846 
847  if (direction == 1) then
848  itmp1 = nc_dstgrdcntrlat_id
849  itmp2 = nc_dstgrdcntrlon_id
850  itmp3 = nc_dstgrdcrnrlat_id
851  itmp4 = nc_dstgrdcrnrlon_id
852  else
853  itmp1 = nc_srcgrdcntrlat_id
854  itmp2 = nc_srcgrdcntrlon_id
855  itmp3 = nc_srcgrdcrnrlat_id
856  itmp4 = nc_srcgrdcrnrlon_id
857  endif
858 
859  ncstat = nf90_put_var(nc_file_id, itmp1, grid2_center_lat)
860  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
861  & 'error writing destination grid center lats')) return
862 
863  ncstat = nf90_put_var(nc_file_id, itmp2, grid2_center_lon)
864  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
865  & 'error writing destination grid center lons')) return
866 
867  ncstat = nf90_put_var(nc_file_id, itmp3, grid2_corner_lat)
868  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
869  & 'error writing destination grid corner lats')) return
870 
871  ncstat = nf90_put_var(nc_file_id, itmp4, grid2_corner_lon)
872  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
873  & 'error writing destination grid corner lons')) return
874 
875  if (direction == 1) then
876  itmp1 = nc_srcgrdarea_id
877  itmp2 = nc_srcgrdfrac_id
878  itmp3 = nc_dstgrdarea_id
879  itmp4 = nc_dstgrdfrac_id
880  else
881  itmp1 = nc_dstgrdarea_id
882  itmp2 = nc_dstgrdfrac_id
883  itmp3 = nc_srcgrdarea_id
884  itmp4 = nc_srcgrdfrac_id
885  endif
886 
887  if (luse_grid1_area) then
888  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_area_in)
889  else
890  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_area)
891  endif
892  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
893  & 'error writing grid1 area')) return
894 
895  ncstat = nf90_put_var(nc_file_id, itmp2, grid1_frac)
896  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
897  & 'error writing grid1 frac')) return
898 
899  if (luse_grid2_area) then
900  ncstat = nf90_put_var(nc_file_id, itmp3, grid2_area_in)
901  else
902  ncstat = nf90_put_var(nc_file_id, itmp3, grid2_area)
903  endif
904  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
905  & 'error writing grid2 area')) return
906 
907  ncstat = nf90_put_var(nc_file_id, itmp4, grid2_frac)
908  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
909  & 'error writing grid2 frac')) return
910 
911  if (direction == 1) then
912  ncstat = nf90_put_var(nc_file_id, nc_srcadd_id,
913  & grid1_add_map1)
914  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
915  & 'error writing source addresses')) return
916 
917  ncstat = nf90_put_var(nc_file_id, nc_dstadd_id,
918  & grid2_add_map1)
919  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
920  & 'error writing destination addresses')) return
921 
922  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix_id, wts_map1)
923  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
924  & 'error writing weights')) return
925  else
926  ncstat = nf90_put_var(nc_file_id, nc_srcadd_id,
927  & grid2_add_map2)
928  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
929  & 'error writing source addresses')) return
930 
931  ncstat = nf90_put_var(nc_file_id, nc_dstadd_id,
932  & grid1_add_map2)
933  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
934  & 'error writing destination addresses')) return
935 
936  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix_id, wts_map2)
937  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
938  & 'error writing weights')) return
939  endif
940 
941  ncstat = nf90_close(nc_file_id)
942  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
943  & 'error closing file')) return
944 
945 !-----------------------------------------------------------------------
946 
947  end subroutine write_remap_scrip
948 
949 !***********************************************************************
950 
951  subroutine write_remap_ww3(map1_name, interp_file1,
952  & output_opt, l_master, errorCode)
953 
954 !-----------------------------------------------------------------------
955 !
956 ! calls correct output routine for WW3 based on output format choice
957 !
958 ! only output variables needed for WW3 remapping
959 !
960 !-----------------------------------------------------------------------
961  use netcdf
962 
963 !-----------------------------------------------------------------------
964 !
965 ! input variables
966 !
967 !-----------------------------------------------------------------------
968 
969  character(SCRIP_charLength), intent(in) ::
970  & map1_name, ! name for mapping grid1 to grid2
971  & interp_file1, ! filename for map1 remap data
972  & output_opt ! option for output conventions
973  logical, intent(in) ::
974  & l_master ! Am I the master processor?
975 
976 !-----------------------------------------------------------------------
977 !
978 ! output variables
979 !
980 !-----------------------------------------------------------------------
981 
982  integer (SCRIP_i4), intent(out) ::
983  & errorCode ! returned error code
984 
985 !-----------------------------------------------------------------------
986 !
987 ! local variables
988 !
989 !-----------------------------------------------------------------------
990 
991  character (15), parameter :: rtnName = 'write_remap_ww3'
992  character(scrip_charlength) ::
993  & map1_name_pass, ! name for mapping grid1 to grid2
994  & interp_file1_pass ! filename for map1 remap data
995 
996 !-----------------------------------------------------------------------
997 !
998 ! define some common variables to be used in all routines
999 !
1000 !-----------------------------------------------------------------------
1001 
1002  errorcode = scrip_success
1003  map1_name_pass = map1_name
1004  interp_file1_pass = interp_file1
1005 
1006  select case(norm_opt)
1007  case (norm_opt_none)
1008  normalize_opt = 'none'
1009  case (norm_opt_frcarea)
1010  normalize_opt = 'fracarea'
1011  case (norm_opt_dstarea)
1012  normalize_opt = 'destarea'
1013  end select
1014 
1015  select case(map_type)
1016  case(map_type_conserv)
1017  map_method = 'Conservative remapping'
1018  case(map_type_bilinear)
1019  map_method = 'Bilinear remapping'
1020  case(map_type_distwgt)
1021  map_method = 'Distance weighted avg of nearest neighbors'
1022  case(map_type_bicubic)
1023  map_method = 'Bicubic remapping'
1024  case(map_type_particle)
1025  map_method = 'Particle remapping'
1026  case default
1027  call scrip_errorset(errorcode, rtnname, 'Invalid Map Type')
1028  return
1029  end select
1030 
1031  call date_and_time(date=cdate)
1032  write (history,1000) cdate(5:6),cdate(7:8),cdate(1:4)
1033  1000 format('Created: ',a2,'-',a2,'-',a4)
1034 
1035 !-----------------------------------------------------------------------
1036 !
1037 ! sort address and weight arrays
1038 !
1039 !-----------------------------------------------------------------------
1040 
1042  if (num_maps > 1) then
1044  endif
1045 
1046 !-----------------------------------------------------------------------
1047 !
1048 ! call appropriate output routine
1049 !
1050 !-----------------------------------------------------------------------
1051 
1052  select case(output_opt)
1053  case ('scrip')
1054  if (l_master) then
1055  call write_remap_scrip_ww3(map1_name_pass,
1056  & interp_file1_pass, errorcode)
1057  endif
1058  if (scrip_errorcheck(errorcode, rtnname,
1059  & 'error in write_remap_scrip')) return
1060  case default
1061  call scrip_errorset(errorcode, rtnname,
1062  & 'unknown output file convention')
1063  return
1064  end select
1065 
1066 !-----------------------------------------------------------------------
1067 
1068  end subroutine write_remap_ww3
1069 
1070 !***********************************************************************
1071 
1072  subroutine write_remap_scrip_ww3(map_name, interp_file, errorCode)
1074 !-----------------------------------------------------------------------
1075 !
1076 ! writes remap data to a netCDF file using SCRIP conventions
1077 ! only writes variables needed for WW3 remap
1078 !
1079 !-----------------------------------------------------------------------
1080 
1081 !-----------------------------------------------------------------------
1082 !
1083 ! input variables
1084 !
1085 !-----------------------------------------------------------------------
1086 
1087  character(SCRIP_charLength), intent(in) ::
1088  & map_name ! name for mapping
1089  &, interp_file ! filename for remap data
1090 
1091 !-----------------------------------------------------------------------
1092 !
1093 ! output variables
1094 !
1095 !-----------------------------------------------------------------------
1096 
1097  integer (SCRIP_i4), intent(out) ::
1098  & errorCode ! returned error code
1099 
1100 !-----------------------------------------------------------------------
1101 !
1102 ! local variables
1103 !
1104 !-----------------------------------------------------------------------
1105 
1106  character(SCRIP_charLength) ::
1107  & grid1_ctmp ! character temp for grid1 names
1108  &, grid2_ctmp ! character temp for grid2 names
1109  &, map_name_ctmp ! character temp for name for mapping
1110  &, interp_file_ctmp ! character temp filename for remap data
1111 
1112  integer (SCRIP_i4) ::
1113  & nc_file_id ! netCDF file id
1114  &, itmp1 ! integer temp
1115  &, itmp2 ! integer temp
1116  &, itmp3 ! integer temp
1117  &, itmp4 ! integer temp
1118 
1119  character (21), parameter :: rtnName = 'write_remap_scrip_ww3'
1120 
1121 !-----------------------------------------------------------------------
1122 !
1123 ! create netCDF file for mapping and define some global attributes
1124 !
1125 !-----------------------------------------------------------------------
1126 
1127  errorcode = scrip_success
1128  map_name_ctmp = map_name
1129  interp_file_ctmp = trim(interp_file)
1130 
1131  itmp1 = len(interp_file_ctmp)
1132  itmp2 = len_trim(interp_file_ctmp)
1133  ncstat = nf90_create(interp_file_ctmp, nf90_clobber, nc_file_id)
1134  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1135  & 'error creating remap file')) return
1136 
1137  !***
1138  !*** map name
1139  !***
1140  ncstat = nf90_put_att(nc_file_id, nf90_global, 'title', map_name)
1141  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1142  & 'error writing remap name')) return
1143 
1144  !***
1145  !*** normalization option
1146  !***
1147  ncstat = nf90_put_att(nc_file_id, nf90_global, 'normalization',
1148  & normalize_opt)
1149  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1150  & 'error writing normalize option')) return
1151 
1152  !***
1153  !*** map method
1154  !***
1155  ncstat = nf90_put_att(nc_file_id, nf90_global, 'map_method',
1156  & map_method)
1157  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1158  & 'error writing remap method')) return
1159 
1160  !***
1161  !*** history
1162  !***
1163  ncstat = nf90_put_att(nc_file_id, nf90_global, 'history',
1164  & history)
1165  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1166  & 'error writing history')) return
1167 
1168  !***
1169  !*** file convention
1170  !***
1171  convention = 'SCRIP'
1172  ncstat = nf90_put_att(nc_file_id, nf90_global, 'conventions',
1173  & convention)
1174  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1175  & 'error writing output convention')) return
1176 
1177  !***
1178  !*** source and destination grid names
1179  !***
1180 
1181  grid1_ctmp = 'source_grid'
1182  grid2_ctmp = 'dest_grid'
1183 
1184  ncstat = nf90_put_att(nc_file_id, nf90_global, trim(grid1_ctmp),
1185  & grid1_name)
1186  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1187  & 'error writing source grid name')) return
1188 
1189  ncstat = nf90_put_att(nc_file_id, nf90_global, trim(grid2_ctmp),
1190  & grid2_name)
1191  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1192  & 'error writing destination grid name')) return
1193 
1194 !-----------------------------------------------------------------------
1195 !
1196 ! prepare netCDF dimension info
1197 !
1198 !-----------------------------------------------------------------------
1199 
1200  !***
1201  !*** define grid size dimensions
1202  !***
1203 
1204  itmp2 = grid2_size
1205 
1206  ncstat = nf90_def_dim(nc_file_id, 'dst_grid_size', itmp2,
1207  & nc_dstgrdsize_id)
1208  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1209  & 'error defining destination grid size')) return
1210 
1211  !***
1212  !*** define map size dimensions
1213  !***
1214 
1215  itmp1 = num_links_map1
1216 
1217  ncstat = nf90_def_dim(nc_file_id, 'num_links',
1218  & itmp1, nc_numlinks_id)
1219  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1220  & 'error defining remap size')) return
1221 
1222  ncstat = nf90_def_dim(nc_file_id, 'num_wgts',
1223  & num_wts, nc_numwgts_id)
1224  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1225  & 'error defining number of weights')) return
1226 
1227 !-----------------------------------------------------------------------
1228 !
1229 ! define all arrays for netCDF descriptors
1230 !
1231 !-----------------------------------------------------------------------
1232 
1233 
1234  !***
1235  !*** define grid fraction arrays
1236  !***
1237 
1238 
1239  ncstat = nf90_def_var(nc_file_id, 'dst_grid_frac',
1240  & nf90_double, nc_dstgrdsize_id,
1241  & nc_dstgrdfrac_id)
1242  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1243  & 'error defining destination fraction')) return
1244 
1245  ncstat = nf90_put_att(nc_file_id, nc_dstgrdfrac_id,
1246  & 'units', 'unitless')
1247  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1248  & 'error writing destination frac units')) return
1249 
1250  !***
1251  !*** define mapping arrays
1252  !***
1253 
1254  ncstat = nf90_def_var(nc_file_id, 'src_address',
1255  & nf90_int, nc_numlinks_id,
1256  & nc_srcadd_id)
1257  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1258  & 'error defining source addresses')) return
1259 
1260  ncstat = nf90_def_var(nc_file_id, 'dst_address',
1261  & nf90_int, nc_numlinks_id,
1262  & nc_dstadd_id)
1263  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1264  & 'error defining destination addresses')) return
1265 
1266  nc_dims2_id(1) = nc_numwgts_id
1267  nc_dims2_id(2) = nc_numlinks_id
1268 
1269  ncstat = nf90_def_var(nc_file_id, 'remap_matrix',
1270  & nf90_double, nc_dims2_id,
1271  & nc_rmpmatrix_id)
1272  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1273  & 'error defining remapping weights')) return
1274 
1275  !***
1276  !*** end definition stage
1277  !***
1278 
1279  ncstat = nf90_enddef(nc_file_id)
1280  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1281  & 'error ending definition phase')) return
1282 
1283 
1284  !***
1285  !*** write variable arrays
1286  !***
1287 
1288  itmp4 = nc_dstgrdfrac_id
1289 
1290  ncstat = nf90_put_var(nc_file_id, itmp4, grid2_frac)
1291  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1292  & 'error writing grid2 frac')) return
1293 
1294  ncstat = nf90_put_var(nc_file_id, nc_srcadd_id,
1295  & grid1_add_map1)
1296  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1297  & 'error writing source addresses')) return
1298 
1299  ncstat = nf90_put_var(nc_file_id, nc_dstadd_id,
1300  & grid2_add_map1)
1301  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1302  & 'error writing destination addresses')) return
1303 
1304  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix_id, wts_map1)
1305  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1306  & 'error writing weights')) return
1307 
1308 
1309  ncstat = nf90_close(nc_file_id)
1310  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1311  & 'error closing file')) return
1312 
1313 
1314 !-----------------------------------------------------------------------
1315 
1316  end subroutine write_remap_scrip_ww3
1317 
1318 !***********************************************************************
1319 
1320  subroutine write_remap_csm(map_name, interp_file, direction,
1321  & errorCode)
1323 !-----------------------------------------------------------------------
1324 !
1325 ! writes remap data to a netCDF file using NCAR-CSM conventions
1326 !
1327 !-----------------------------------------------------------------------
1328 
1329 !-----------------------------------------------------------------------
1330 !
1331 ! input variables
1332 !
1333 !-----------------------------------------------------------------------
1334 
1335  character(SCRIP_charLength), intent(in) ::
1336  & map_name ! name for mapping
1337  &, interp_file ! filename for remap data
1338 
1339  integer (SCRIP_i4), intent(in) ::
1340  & direction ! direction of map (1=grid1 to grid2
1341  ! 2=grid2 to grid1)
1342 
1343 !-----------------------------------------------------------------------
1344 !
1345 ! output variables
1346 !
1347 !-----------------------------------------------------------------------
1348 
1349  integer (SCRIP_i4), intent(out) ::
1350  & errorCode ! returned error code
1351 
1352 !-----------------------------------------------------------------------
1353 !
1354 ! local variables
1355 !
1356 !-----------------------------------------------------------------------
1357 
1358  character(SCRIP_charLength) ::
1359  & grid1_ctmp ! character temp for grid1 names
1360  &, grid2_ctmp ! character temp for grid2 names
1361 
1362  integer (SCRIP_i4) ::
1363  & itmp1 ! integer temp
1364  &, itmp2 ! integer temp
1365  &, itmp3 ! integer temp
1366  &, itmp4 ! integer temp
1367  &, nc_numwgts1_id ! extra netCDF id for additional weights
1368  &, nc_src_isize_id ! extra netCDF id for ni_a
1369  &, nc_src_jsize_id ! extra netCDF id for nj_a
1370  &, nc_dst_isize_id ! extra netCDF id for ni_b
1371  &, nc_dst_jsize_id ! extra netCDF id for nj_b
1372  &, nc_rmpmatrix2_id ! extra netCDF id for high-order remap matrix
1373 
1374  real (SCRIP_r8), dimension(:),allocatable ::
1375  & wts1 ! CSM wants single array for 1st-order wts
1376 
1377  real (SCRIP_r8), dimension(:,:),allocatable ::
1378  & wts2 ! write remaining weights in different array
1379 
1380  character (15), parameter :: rtnName = 'write_remap_csm'
1381 
1382 !-----------------------------------------------------------------------
1383 !
1384 ! create netCDF file for mapping and define some global attributes
1385 !
1386 !-----------------------------------------------------------------------
1387 
1388  errorcode = scrip_success
1389 
1390  ncstat = nf90_create(interp_file, nf90_clobber, nc_file_id)
1391  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1392  & 'error opening file')) return
1393 
1394  !***
1395  !*** map name
1396  !***
1397  ncstat = nf90_put_att(nc_file_id, nf90_global, 'title', map_name)
1398  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1399  & 'error writing remap name')) return
1400 
1401  !***
1402  !*** normalization option
1403  !***
1404  ncstat = nf90_put_att(nc_file_id, nf90_global, 'normalization',
1405  & normalize_opt)
1406  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1407  & 'error writing normalization option')) return
1408 
1409  !***
1410  !*** map method
1411  !***
1412  ncstat = nf90_put_att(nc_file_id, nf90_global, 'map_method',
1413  & map_method)
1414  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1415  & 'error writing remap method')) return
1416 
1417  !***
1418  !*** history
1419  !***
1420  ncstat = nf90_put_att(nc_file_id, nf90_global, 'history',
1421  & history)
1422  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1423  & 'error writing history')) return
1424 
1425  !***
1426  !*** file convention
1427  !***
1428  convention = 'NCAR-CSM'
1429  ncstat = nf90_put_att(nc_file_id, nf90_global, 'conventions',
1430  & convention)
1431  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1432  & 'error writing output convention')) return
1433 
1434  !***
1435  !*** source and destination grid names
1436  !***
1437 
1438  if (direction == 1) then
1439  grid1_ctmp = 'domain_a'
1440  grid2_ctmp = 'domain_b'
1441  else
1442  grid1_ctmp = 'domain_b'
1443  grid2_ctmp = 'domain_a'
1444  endif
1445 
1446  ncstat = nf90_put_att(nc_file_id, nf90_global, trim(grid1_ctmp),
1447  & grid1_name)
1448  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1449  & 'error writing grid1 name')) return
1450 
1451  ncstat = nf90_put_att(nc_file_id, nf90_global, trim(grid2_ctmp),
1452  & grid2_name)
1453  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1454  & 'error writing grid2 name')) return
1455 
1456 !-----------------------------------------------------------------------
1457 !
1458 ! prepare netCDF dimension info
1459 !
1460 !-----------------------------------------------------------------------
1461 
1462  !***
1463  !*** define grid size dimensions
1464  !***
1465 
1466  if (direction == 1) then
1467  itmp1 = grid1_size
1468  itmp2 = grid2_size
1469  else
1470  itmp1 = grid2_size
1471  itmp2 = grid1_size
1472  endif
1473 
1474  ncstat = nf90_def_dim(nc_file_id, 'n_a', itmp1, nc_srcgrdsize_id)
1475  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1476  & 'error defining source grid size')) return
1477 
1478  ncstat = nf90_def_dim(nc_file_id, 'n_b', itmp2, nc_dstgrdsize_id)
1479  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1480  & 'error defining destination grid size')) return
1481 
1482  !***
1483  !*** define grid corner dimension
1484  !***
1485 
1486  if (direction == 1) then
1487  itmp1 = grid1_corners
1488  itmp2 = grid2_corners
1489  else
1490  itmp1 = grid2_corners
1491  itmp2 = grid1_corners
1492  endif
1493 
1494  ncstat = nf90_def_dim(nc_file_id, 'nv_a', itmp1, nc_srcgrdcorn_id)
1495  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1496  & 'error defining number of corners on source grid')) return
1497 
1498  ncstat = nf90_def_dim(nc_file_id, 'nv_b', itmp2, nc_dstgrdcorn_id)
1499  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1500  & 'error defining number of corners on destination grid')) return
1501 
1502  !***
1503  !*** define grid rank dimension
1504  !***
1505 
1506  if (direction == 1) then
1507  itmp1 = grid1_rank
1508  itmp2 = grid2_rank
1509  else
1510  itmp1 = grid2_rank
1511  itmp2 = grid1_rank
1512  endif
1513 
1514  ncstat = nf90_def_dim(nc_file_id, 'src_grid_rank',
1515  & itmp1, nc_srcgrdrank_id)
1516  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1517  & 'error defining source grid rank')) return
1518 
1519  ncstat = nf90_def_dim(nc_file_id, 'dst_grid_rank',
1520  & itmp2, nc_dstgrdrank_id)
1521  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1522  & 'error defining destination grid rank')) return
1523 
1524  !***
1525  !*** define first two dims as if 2-d cartesian domain
1526  !***
1527 
1528  if (direction == 1) then
1529  itmp1 = grid1_dims(1)
1530  if (grid1_rank > 1) then
1531  itmp2 = grid1_dims(2)
1532  else
1533  itmp2 = 0
1534  endif
1535  itmp3 = grid2_dims(1)
1536  if (grid2_rank > 1) then
1537  itmp4 = grid2_dims(2)
1538  else
1539  itmp4 = 0
1540  endif
1541  else
1542  itmp1 = grid2_dims(1)
1543  if (grid2_rank > 1) then
1544  itmp2 = grid2_dims(2)
1545  else
1546  itmp2 = 0
1547  endif
1548  itmp3 = grid1_dims(1)
1549  if (grid1_rank > 1) then
1550  itmp4 = grid1_dims(2)
1551  else
1552  itmp4 = 0
1553  endif
1554  endif
1555 
1556  ncstat = nf90_def_dim(nc_file_id, 'ni_a', itmp1, nc_src_isize_id)
1557  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1558  & 'error defining source isize')) return
1559 
1560  ncstat = nf90_def_dim(nc_file_id, 'nj_a', itmp2, nc_src_jsize_id)
1561  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1562  & 'error defining source jsize')) return
1563 
1564  ncstat = nf90_def_dim(nc_file_id, 'ni_b', itmp3, nc_dst_isize_id)
1565  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1566  & 'error defining destination isize')) return
1567 
1568  ncstat = nf90_def_dim(nc_file_id, 'nj_b', itmp4, nc_dst_jsize_id)
1569  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1570  & 'error defining destination jsize')) return
1571 
1572  !***
1573  !*** define map size dimensions
1574  !***
1575 
1576  if (direction == 1) then
1577  itmp1 = num_links_map1
1578  else
1579  itmp1 = num_links_map2
1580  endif
1581 
1582  ncstat = nf90_def_dim(nc_file_id, 'n_s', itmp1, nc_numlinks_id)
1583  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1584  & 'error defining remap size')) return
1585 
1586  ncstat = nf90_def_dim(nc_file_id, 'num_wgts',
1587  & num_wts, nc_numwgts_id)
1588  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1589  & 'error defining number of weights')) return
1590 
1591  if (num_wts > 1) then
1592  ncstat = nf90_def_dim(nc_file_id, 'num_wgts1',
1593  & num_wts-1, nc_numwgts1_id)
1594  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1595  & 'error defining number of weights1')) return
1596  endif
1597 
1598  !***
1599  !*** define grid dimensions
1600  !***
1601 
1602  ncstat = nf90_def_var(nc_file_id, 'src_grid_dims', nf90_int,
1603  & nc_srcgrdrank_id, nc_srcgrddims_id)
1604  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1605  & 'error defining source grid dims')) return
1606 
1607  ncstat = nf90_def_var(nc_file_id, 'dst_grid_dims', nf90_int,
1608  & nc_dstgrdrank_id, nc_dstgrddims_id)
1609  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1610  & 'error defining destination grid dims')) return
1611 
1612 !-----------------------------------------------------------------------
1613 !
1614 ! define all arrays for netCDF descriptors
1615 !
1616 !-----------------------------------------------------------------------
1617 
1618  !***
1619  !*** define grid center latitude array
1620  !***
1621 
1622  ncstat = nf90_def_var(nc_file_id, 'yc_a',
1623  & nf90_double, nc_srcgrdsize_id,
1624  & nc_srcgrdcntrlat_id)
1625  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1626  & 'error defining source grid center lats')) return
1627 
1628  ncstat = nf90_def_var(nc_file_id, 'yc_b',
1629  & nf90_double, nc_dstgrdsize_id,
1630  & nc_dstgrdcntrlat_id)
1631  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1632  & 'error defining destination grid center lats')) return
1633 
1634  !***
1635  !*** define grid center longitude array
1636  !***
1637 
1638  ncstat = nf90_def_var(nc_file_id, 'xc_a',
1639  & nf90_double, nc_srcgrdsize_id,
1640  & nc_srcgrdcntrlon_id)
1641  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1642  & 'error defining source grid center lons')) return
1643 
1644  ncstat = nf90_def_var(nc_file_id, 'xc_b',
1645  & nf90_double, nc_dstgrdsize_id,
1646  & nc_dstgrdcntrlon_id)
1647  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1648  & 'error defining destination grid center lons')) return
1649 
1650  !***
1651  !*** define grid corner lat/lon arrays
1652  !***
1653 
1654  nc_dims2_id(1) = nc_srcgrdcorn_id
1655  nc_dims2_id(2) = nc_srcgrdsize_id
1656 
1657  ncstat = nf90_def_var(nc_file_id, 'yv_a',
1658  & nf90_double, nc_dims2_id,
1659  & nc_srcgrdcrnrlat_id)
1660  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1661  & 'error defining source grid corner lats')) return
1662 
1663  ncstat = nf90_def_var(nc_file_id, 'xv_a',
1664  & nf90_double, nc_dims2_id,
1665  & nc_srcgrdcrnrlon_id)
1666  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1667  & 'error defining source grid corner lons')) return
1668 
1669  nc_dims2_id(1) = nc_dstgrdcorn_id
1670  nc_dims2_id(2) = nc_dstgrdsize_id
1671 
1672  ncstat = nf90_def_var(nc_file_id, 'yv_b',
1673  & nf90_double, nc_dims2_id,
1674  & nc_dstgrdcrnrlat_id)
1675  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1676  & 'error defining destination grid corner lats')) return
1677 
1678  ncstat = nf90_def_var(nc_file_id, 'xv_b',
1679  & nf90_double, nc_dims2_id,
1680  & nc_dstgrdcrnrlon_id)
1681  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1682  & 'error defining destination grid corner lons')) return
1683 
1684  !***
1685  !*** CSM wants all in degrees
1686  !***
1687 
1688  grid1_units = 'degrees'
1689  grid2_units = 'degrees'
1690 
1691  if (direction == 1) then
1692  grid1_ctmp = grid1_units
1693  grid2_ctmp = grid2_units
1694  else
1695  grid1_ctmp = grid2_units
1696  grid2_ctmp = grid1_units
1697  endif
1698 
1699  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcntrlat_id,
1700  & 'units', grid1_ctmp)
1701  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1702  & 'error writing grid units')) return
1703 
1704  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcntrlat_id,
1705  & 'units', grid2_ctmp)
1706  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1707  & 'error writing grid units')) return
1708 
1709  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcntrlon_id,
1710  & 'units', grid1_ctmp)
1711  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1712  & 'error writing grid units')) return
1713 
1714  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcntrlon_id,
1715  & 'units', grid2_ctmp)
1716  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1717  & 'error writing grid units')) return
1718 
1719  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcrnrlat_id,
1720  & 'units', grid1_ctmp)
1721  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1722  & 'error writing grid units')) return
1723 
1724  ncstat = nf90_put_att(nc_file_id, nc_srcgrdcrnrlon_id,
1725  & 'units', grid1_ctmp)
1726  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1727  & 'error writing grid units')) return
1728 
1729  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcrnrlat_id,
1730  & 'units', grid2_ctmp)
1731  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1732  & 'error writing grid units')) return
1733 
1734  ncstat = nf90_put_att(nc_file_id, nc_dstgrdcrnrlon_id,
1735  & 'units', grid2_ctmp)
1736  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1737  & 'error writing grid units')) return
1738 
1739  !***
1740  !*** define grid mask
1741  !***
1742 
1743  ncstat = nf90_def_var(nc_file_id, 'mask_a', nf90_int,
1744  & nc_srcgrdsize_id, nc_srcgrdimask_id)
1745  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1746  & 'error defining source grid mask')) return
1747 
1748  ncstat = nf90_put_att(nc_file_id, nc_srcgrdimask_id,
1749  & 'units', 'unitless')
1750  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1751  & 'error writing source mask units')) return
1752 
1753  ncstat = nf90_def_var(nc_file_id, 'mask_b', nf90_int,
1754  & nc_dstgrdsize_id, nc_dstgrdimask_id)
1755  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1756  & 'error defining destination grid mask')) return
1757 
1758  ncstat = nf90_put_att(nc_file_id, nc_dstgrdimask_id,
1759  & 'units', 'unitless')
1760  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1761  & 'error writing destination mask units')) return
1762 
1763  !***
1764  !*** define grid area arrays
1765  !***
1766 
1767  ncstat = nf90_def_var(nc_file_id, 'area_a',
1768  & nf90_double, nc_srcgrdsize_id,
1769  & nc_srcgrdarea_id)
1770  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1771  & 'error defining source grid area')) return
1772 
1773  ncstat = nf90_put_att(nc_file_id, nc_srcgrdarea_id,
1774  & 'units', 'square radians')
1775  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1776  & 'error defining source area units')) return
1777 
1778  ncstat = nf90_def_var(nc_file_id, 'area_b',
1779  & nf90_double, nc_dstgrdsize_id,
1780  & nc_dstgrdarea_id)
1781  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1782  & 'error defining destination grid area')) return
1783 
1784  ncstat = nf90_put_att(nc_file_id, nc_dstgrdarea_id,
1785  & 'units', 'square radians')
1786  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1787  & 'error defining destination area units')) return
1788 
1789  !***
1790  !*** define grid fraction arrays
1791  !***
1792 
1793  ncstat = nf90_def_var(nc_file_id, 'frac_a',
1794  & nf90_double, nc_srcgrdsize_id,
1795  & nc_srcgrdfrac_id)
1796  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1797  & 'error defining source grid frac')) return
1798 
1799  ncstat = nf90_put_att(nc_file_id, nc_srcgrdfrac_id,
1800  & 'units', 'unitless')
1801  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1802  & 'error defining source frac units')) return
1803 
1804  ncstat = nf90_def_var(nc_file_id, 'frac_b',
1805  & nf90_double, nc_dstgrdsize_id,
1806  & nc_dstgrdfrac_id)
1807  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1808  & 'error defining destination grid frac')) return
1809 
1810  ncstat = nf90_put_att(nc_file_id, nc_dstgrdfrac_id,
1811  & 'units', 'unitless')
1812  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1813  & 'error defining destination frac units')) return
1814 
1815  !***
1816  !*** define mapping arrays
1817  !***
1818 
1819  ncstat = nf90_def_var(nc_file_id, 'col',
1820  & nf90_int, nc_numlinks_id,
1821  & nc_srcadd_id)
1822  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1823  & 'error defining source addresses')) return
1824 
1825  ncstat = nf90_def_var(nc_file_id, 'row',
1826  & nf90_int, nc_numlinks_id,
1827  & nc_dstadd_id)
1828  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1829  & 'error defining destination addresses')) return
1830 
1831  ncstat = nf90_def_var(nc_file_id, 'S',
1832  & nf90_double, nc_numlinks_id,
1833  & nc_rmpmatrix_id)
1834  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1835  & 'error defining weights')) return
1836 
1837  if (num_wts > 1) then
1838  nc_dims2_id(1) = nc_numwgts1_id
1839  nc_dims2_id(2) = nc_numlinks_id
1840 
1841  ncstat = nf90_def_var(nc_file_id, 'S2',
1842  & nf90_double, nc_dims2_id,
1843  & nc_rmpmatrix2_id)
1844  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1845  & 'error defining weights2')) return
1846  endif
1847 
1848  !***
1849  !*** end definition stage
1850  !***
1851 
1852  ncstat = nf90_enddef(nc_file_id)
1853  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1854  & 'error ending definition phase')) return
1855 
1856 !-----------------------------------------------------------------------
1857 !
1858 ! compute integer masks
1859 !
1860 !-----------------------------------------------------------------------
1861 
1862  if (direction == 1) then
1863  allocate (src_mask_int(grid1_size),
1864  & dst_mask_int(grid2_size))
1865 
1866  where (grid2_mask)
1867  dst_mask_int = 1
1868  elsewhere
1869  dst_mask_int = 0
1870  endwhere
1871 
1872  where (grid1_mask)
1873  src_mask_int = 1
1874  elsewhere
1875  src_mask_int = 0
1876  endwhere
1877  else
1878  allocate (src_mask_int(grid2_size),
1879  & dst_mask_int(grid1_size))
1880 
1881  where (grid1_mask)
1882  dst_mask_int = 1
1883  elsewhere
1884  dst_mask_int = 0
1885  endwhere
1886 
1887  where (grid2_mask)
1888  src_mask_int = 1
1889  elsewhere
1890  src_mask_int = 0
1891  endwhere
1892  endif
1893 
1894 !-----------------------------------------------------------------------
1895 !
1896 ! change units of lat/lon coordinates if input units different
1897 ! from radians. if this is the second mapping, the conversion has
1898 ! alread been done.
1899 !
1900 !-----------------------------------------------------------------------
1901 
1902  if (grid1_units(1:7) == 'degrees' .and. direction == 1) then
1907  endif
1908 
1909  if (grid2_units(1:7) == 'degrees' .and. direction == 1) then
1914  endif
1915 
1916 !-----------------------------------------------------------------------
1917 !
1918 ! write mapping data
1919 !
1920 !-----------------------------------------------------------------------
1921 
1922  if (direction == 1) then
1923  itmp1 = nc_srcgrddims_id
1924  itmp2 = nc_dstgrddims_id
1925  else
1926  itmp2 = nc_srcgrddims_id
1927  itmp1 = nc_dstgrddims_id
1928  endif
1929 
1930  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_dims)
1931  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1932  & 'error writing grid1 dims')) return
1933 
1934  ncstat = nf90_put_var(nc_file_id, itmp2, grid2_dims)
1935  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1936  & 'error writing grid2 dims')) return
1937 
1938  ncstat = nf90_put_var(nc_file_id, nc_srcgrdimask_id,
1939  & src_mask_int)
1940  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1941  & 'error writing source grid mask')) return
1942 
1943  ncstat = nf90_put_var(nc_file_id, nc_dstgrdimask_id,
1944  & dst_mask_int)
1945  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1946  & 'error writing destination grid mask')) return
1947 
1948  deallocate(src_mask_int, dst_mask_int)
1949 
1950  if (direction == 1) then
1951  itmp1 = nc_srcgrdcntrlat_id
1952  itmp2 = nc_srcgrdcntrlon_id
1953  itmp3 = nc_srcgrdcrnrlat_id
1954  itmp4 = nc_srcgrdcrnrlon_id
1955  else
1956  itmp1 = nc_dstgrdcntrlat_id
1957  itmp2 = nc_dstgrdcntrlon_id
1958  itmp3 = nc_dstgrdcrnrlat_id
1959  itmp4 = nc_dstgrdcrnrlon_id
1960  endif
1961 
1962  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_center_lat)
1963  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1964  & 'error writing grid1 center lats')) return
1965 
1966  ncstat = nf90_put_var(nc_file_id, itmp2, grid1_center_lon)
1967  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1968  & 'error writing grid1 center lons')) return
1969 
1970  ncstat = nf90_put_var(nc_file_id, itmp3, grid1_corner_lat)
1971  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1972  & 'error writing grid1 corner lats')) return
1973 
1974  ncstat = nf90_put_var(nc_file_id, itmp4, grid1_corner_lon)
1975  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1976  & 'error writing grid1 corner lons')) return
1977 
1978  if (direction == 1) then
1979  itmp1 = nc_dstgrdcntrlat_id
1980  itmp2 = nc_dstgrdcntrlon_id
1981  itmp3 = nc_dstgrdcrnrlat_id
1982  itmp4 = nc_dstgrdcrnrlon_id
1983  else
1984  itmp1 = nc_srcgrdcntrlat_id
1985  itmp2 = nc_srcgrdcntrlon_id
1986  itmp3 = nc_srcgrdcrnrlat_id
1987  itmp4 = nc_srcgrdcrnrlon_id
1988  endif
1989 
1990  ncstat = nf90_put_var(nc_file_id, itmp1, grid2_center_lat)
1991  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1992  & 'error writing grid2 center lats')) return
1993 
1994  ncstat = nf90_put_var(nc_file_id, itmp2, grid2_center_lon)
1995  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
1996  & 'error writing grid2 center lons')) return
1997 
1998  ncstat = nf90_put_var(nc_file_id, itmp3, grid2_corner_lat)
1999  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2000  & 'error writing grid2 corner lats')) return
2001 
2002  ncstat = nf90_put_var(nc_file_id, itmp4, grid2_corner_lon)
2003  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2004  & 'error writing grid2 corner lons')) return
2005 
2006  if (direction == 1) then
2007  itmp1 = nc_srcgrdarea_id
2008  itmp2 = nc_srcgrdfrac_id
2009  itmp3 = nc_dstgrdarea_id
2010  itmp4 = nc_dstgrdfrac_id
2011  else
2012  itmp1 = nc_dstgrdarea_id
2013  itmp2 = nc_dstgrdfrac_id
2014  itmp3 = nc_srcgrdarea_id
2015  itmp4 = nc_srcgrdfrac_id
2016  endif
2017 
2018  if (luse_grid1_area) then
2019  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_area_in)
2020  else
2021  ncstat = nf90_put_var(nc_file_id, itmp1, grid1_area)
2022  endif
2023  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2024  & 'error writing grid1 area')) return
2025 
2026  ncstat = nf90_put_var(nc_file_id, itmp2, grid1_frac)
2027  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2028  & 'error writing grid1 frac')) return
2029 
2030  if (luse_grid2_area) then
2031  ncstat = nf90_put_var(nc_file_id, itmp3, grid2_area)
2032  else
2033  ncstat = nf90_put_var(nc_file_id, itmp3, grid2_area)
2034  endif
2035  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2036  & 'error writing grid2 area')) return
2037 
2038  ncstat = nf90_put_var(nc_file_id, itmp4, grid2_frac)
2039  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2040  & 'error writing grid2 frac')) return
2041 
2042  if (direction == 1) then
2043  ncstat = nf90_put_var(nc_file_id, nc_srcadd_id,
2044  & grid1_add_map1)
2045  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2046  & 'error writing source addresses')) return
2047 
2048  ncstat = nf90_put_var(nc_file_id, nc_dstadd_id,
2049  & grid2_add_map1)
2050  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2051  & 'error writing destination addresses')) return
2052 
2053  if (num_wts == 1) then
2054  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix_id,
2055  & wts_map1(1,:))
2056  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2057  & 'error writing weights')) return
2058  else
2059  allocate(wts1(num_links_map1),wts2(num_wts-1,num_links_map1))
2060 
2061  wts1 = wts_map1(1,:)
2062  wts2 = wts_map1(2:,:)
2063 
2064  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix_id, wts1)
2065  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2066  & 'error writing weights1')) return
2067  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix2_id, wts2)
2068  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2069  & 'error writing weights2')) return
2070  deallocate(wts1,wts2)
2071  endif
2072  else
2073  ncstat = nf90_put_var(nc_file_id, nc_srcadd_id,
2074  & grid2_add_map2)
2075  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2076  & 'error writing source addresses')) return
2077 
2078  ncstat = nf90_put_var(nc_file_id, nc_dstadd_id,
2079  & grid1_add_map2)
2080  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2081  & 'error writing destination addresses')) return
2082 
2083  if (num_wts == 1) then
2084  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix_id,
2085  & wts_map2(1,:))
2086  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2087  & 'error writing weights')) return
2088  else
2089  allocate(wts1(num_links_map2),wts2(num_wts-1,num_links_map2))
2090 
2091  wts1 = wts_map2(1,:)
2092  wts2 = wts_map2(2:,:)
2093 
2094  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix_id, wts1)
2095  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2096  & 'error writing weights1')) return
2097  ncstat = nf90_put_var(nc_file_id, nc_rmpmatrix2_id, wts2)
2098  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2099  & 'error writing weights2')) return
2100  deallocate(wts1,wts2)
2101  endif
2102  endif
2103 
2104  ncstat = nf90_close(nc_file_id)
2105  if (scrip_netcdferrorcheck(ncstat, errorcode, rtnname,
2106  & 'error closing file')) return
2107 
2108 !-----------------------------------------------------------------------
2109 
2110  end subroutine write_remap_csm
2111 
2112 !***********************************************************************
2113 
2114  subroutine sort_add(add1, add2, weights)
2116 !-----------------------------------------------------------------------
2117 !
2118 ! this routine sorts address and weight arrays based on the
2119 ! destination address with the source address as a secondary
2120 ! sorting criterion. the method is a standard heap sort.
2121 !
2122 !-----------------------------------------------------------------------
2123 
2124  implicit none
2125 
2126 !-----------------------------------------------------------------------
2127 !
2128 ! Input and Output arrays
2129 !
2130 !-----------------------------------------------------------------------
2131 
2132  integer (SCRIP_i4), intent(inout), dimension(:) ::
2133  & add1, ! destination address array (num_links)
2134  & add2 ! source address array
2135 
2136  real (SCRIP_r8), intent(inout), dimension(:,:) ::
2137  & weights ! remapping weights (num_wts, num_links)
2138 
2139 !-----------------------------------------------------------------------
2140 !
2141 ! local variables
2142 !
2143 !-----------------------------------------------------------------------
2144 
2145  integer (SCRIP_i4) ::
2146  & num_links, ! num of links for this mapping
2147  & num_wts, ! num of weights for this mapping
2148  & add1_tmp, add2_tmp, ! temp for addresses during swap
2149  & lvl, final_lvl, ! level indexes for heap sort levels
2150  & chk_lvl1, chk_lvl2, max_lvl
2151 
2152  real (SCRIP_r8), dimension(SIZE(weights,DIM=1)) ::
2153  & wgttmp ! temp for holding wts during swap
2154 
2155 !-----------------------------------------------------------------------
2156 !
2157 ! determine total number of links to sort and number of weights
2158 !
2159 !-----------------------------------------------------------------------
2160 
2161  num_links = SIZE(add1)
2162  num_wts = SIZE(weights, dim=1)
2163 
2164 !-----------------------------------------------------------------------
2165 !
2166 ! start at the lowest level (N/2) of the tree and sift lower
2167 ! values to the bottom of the tree, promoting the larger numbers
2168 !
2169 !-----------------------------------------------------------------------
2170 
2171  do lvl=num_links/2,1,-1
2172 
2173  final_lvl = lvl
2174  add1_tmp = add1(lvl)
2175  add2_tmp = add2(lvl)
2176  wgttmp(:) = weights(:,lvl)
2177 
2178  !***
2179  !*** loop until proper level is found for this link, or reach
2180  !*** bottom
2181  !***
2182 
2183  sift_loop1: do
2184 
2185  !***
2186  !*** find the largest of the two daughters
2187  !***
2188 
2189  chk_lvl1 = 2*final_lvl
2190  chk_lvl2 = 2*final_lvl+1
2191  if (chk_lvl1 .EQ. num_links) chk_lvl2 = chk_lvl1
2192 
2193  if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR.
2194  & ((add1(chk_lvl1) == add1(chk_lvl2)) .AND.
2195  & (add2(chk_lvl1) > add2(chk_lvl2)))) then
2196  max_lvl = chk_lvl1
2197  else
2198  max_lvl = chk_lvl2
2199  endif
2200 
2201  !***
2202  !*** if the parent is greater than both daughters,
2203  !*** the correct level has been found
2204  !***
2205 
2206  if ((add1_tmp .GT. add1(max_lvl)) .OR.
2207  & ((add1_tmp .EQ. add1(max_lvl)) .AND.
2208  & (add2_tmp .GT. add2(max_lvl)))) then
2209  add1(final_lvl) = add1_tmp
2210  add2(final_lvl) = add2_tmp
2211  weights(:,final_lvl) = wgttmp(:)
2212  exit sift_loop1
2213 
2214  !***
2215  !*** otherwise, promote the largest daughter and push
2216  !*** down one level in the tree. if haven't reached
2217  !*** the end of the tree, repeat the process. otherwise
2218  !*** store last values and exit the loop
2219  !***
2220 
2221  else
2222  add1(final_lvl) = add1(max_lvl)
2223  add2(final_lvl) = add2(max_lvl)
2224  weights(:,final_lvl) = weights(:,max_lvl)
2225 
2226  final_lvl = max_lvl
2227  if (2*final_lvl > num_links) then
2228  add1(final_lvl) = add1_tmp
2229  add2(final_lvl) = add2_tmp
2230  weights(:,final_lvl) = wgttmp(:)
2231  exit sift_loop1
2232  endif
2233  endif
2234  end do sift_loop1
2235  end do
2236 
2237 !-----------------------------------------------------------------------
2238 !
2239 ! now that the heap has been sorted, strip off the top (largest)
2240 ! value and promote the values below
2241 !
2242 !-----------------------------------------------------------------------
2243 
2244  do lvl=num_links,3,-1
2245 
2246  !***
2247  !*** move the top value and insert it into the correct place
2248  !***
2249 
2250  add1_tmp = add1(lvl)
2251  add1(lvl) = add1(1)
2252 
2253  add2_tmp = add2(lvl)
2254  add2(lvl) = add2(1)
2255 
2256  wgttmp(:) = weights(:,lvl)
2257  weights(:,lvl) = weights(:,1)
2258 
2259  !***
2260  !*** as above this loop sifts the tmp values down until proper
2261  !*** level is reached
2262  !***
2263 
2264  final_lvl = 1
2265 
2266  sift_loop2: do
2267 
2268  !***
2269  !*** find the largest of the two daughters
2270  !***
2271 
2272  chk_lvl1 = 2*final_lvl
2273  chk_lvl2 = 2*final_lvl+1
2274  if (chk_lvl2 >= lvl) chk_lvl2 = chk_lvl1
2275 
2276  if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR.
2277  & ((add1(chk_lvl1) == add1(chk_lvl2)) .AND.
2278  & (add2(chk_lvl1) > add2(chk_lvl2)))) then
2279  max_lvl = chk_lvl1
2280  else
2281  max_lvl = chk_lvl2
2282  endif
2283 
2284  !***
2285  !*** if the parent is greater than both daughters,
2286  !*** the correct level has been found
2287  !***
2288 
2289  if ((add1_tmp > add1(max_lvl)) .OR.
2290  & ((add1_tmp == add1(max_lvl)) .AND.
2291  & (add2_tmp > add2(max_lvl)))) then
2292  add1(final_lvl) = add1_tmp
2293  add2(final_lvl) = add2_tmp
2294  weights(:,final_lvl) = wgttmp(:)
2295  exit sift_loop2
2296 
2297  !***
2298  !*** otherwise, promote the largest daughter and push
2299  !*** down one level in the tree. if haven't reached
2300  !*** the end of the tree, repeat the process. otherwise
2301  !*** store last values and exit the loop
2302  !***
2303 
2304  else
2305  add1(final_lvl) = add1(max_lvl)
2306  add2(final_lvl) = add2(max_lvl)
2307  weights(:,final_lvl) = weights(:,max_lvl)
2308 
2309  final_lvl = max_lvl
2310  if (2*final_lvl >= lvl) then
2311  add1(final_lvl) = add1_tmp
2312  add2(final_lvl) = add2_tmp
2313  weights(:,final_lvl) = wgttmp(:)
2314  exit sift_loop2
2315  endif
2316  endif
2317  end do sift_loop2
2318  end do
2319 
2320  !***
2321  !*** swap the last two entries
2322  !***
2323 
2324 
2325  add1_tmp = add1(2)
2326  add1(2) = add1(1)
2327  add1(1) = add1_tmp
2328 
2329  add2_tmp = add2(2)
2330  add2(2) = add2(1)
2331  add2(1) = add2_tmp
2332 
2333  wgttmp(:) = weights(:,2)
2334  weights(:,2) = weights(:,1)
2335  weights(:,1) = wgttmp(:)
2336 
2337 !-----------------------------------------------------------------------
2338 
2339  end subroutine sort_add
2340 
2341 !***********************************************************************
2342 
2343  end module scrip_remap_write
2344 
2345 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
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_vars::map_type_bilinear
integer(scrip_i4), parameter map_type_bilinear
Definition: scrip_remap_vars.f:59
scrip_remap_vars::num_maps
integer(scrip_i4), save num_maps
Definition: scrip_remap_vars.f:66
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_netcdfmod
Definition: scrip_netcdfmod.f90:3
scrip_remap_vars::num_links_map1
integer(scrip_i4), save num_links_map1
Definition: scrip_remap_vars.f:66
scrip_grids::grid2_rank
integer(scrip_i4), save grid2_rank
Definition: scrip_grids.f:68
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_kindsmod::scrip_charlength
integer, parameter, public scrip_charlength
Definition: scrip_kindsmod.f90:38
scrip_grids::deg2rad
real(scrip_r8), parameter deg2rad
Definition: scrip_grids.f:84
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_remap_vars::map_type
integer(scrip_i4), save map_type
Definition: scrip_remap_vars.f:66
scrip_grids
Definition: scrip_grids.f:49
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_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_write::write_remap_csm
subroutine write_remap_csm(map_name, interp_file, direction, errorCode)
Definition: scrip_remap_write.f:1322
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_write::write_remap_scrip_ww3
subroutine write_remap_scrip_ww3(map_name, interp_file, errorCode)
Definition: scrip_remap_write.f:1073
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_errormod::scrip_success
integer(scrip_i4), parameter, public scrip_success
Definition: scrip_errormod.f90:42
scrip_remap_write::write_remap_scrip
subroutine write_remap_scrip(map_name, interp_file, direction, errorCode)
Definition: scrip_remap_write.f:267
scrip_grids::grid1_corner_lat
real(scrip_r8), dimension(:,:), allocatable, target, save grid1_corner_lat
Definition: scrip_grids.f:120
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::grid2_name
character(scrip_charlength), save grid2_name
Definition: scrip_grids.f:77
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_remap_vars::wts_map1
real(scrip_r8), dimension(:,:), allocatable, save wts_map1
Definition: scrip_remap_vars.f:83
scrip_constants
Definition: scrip_constants.f:38
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_grids::grid2_dims
integer(scrip_i4), dimension(:), allocatable, save grid2_dims
Definition: scrip_grids.f:74
scrip_grids::grid1_center_lat
real(scrip_r8), dimension(:), allocatable, target, save grid1_center_lat
Definition: scrip_grids.f:103
scrip_kindsmod
Definition: scrip_kindsmod.f90:3
scrip_grids::grid1_corner_lon
real(scrip_r8), dimension(:,:), allocatable, target, save grid1_corner_lon
Definition: scrip_grids.f:120
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::grid1_size
integer(scrip_i4), save grid1_size
Definition: scrip_grids.f:68
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_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_vars
Definition: scrip_remap_vars.f:40
scrip_netcdfmod::scrip_netcdferrorcheck
logical(scrip_logical) function, public scrip_netcdferrorcheck(netcdfStat, errorCode, rtnName, errorMsg)
Definition: scrip_netcdfmod.f90:45
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_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_write::sort_add
subroutine sort_add(add1, add2, weights)
Definition: scrip_remap_write.f:2115
scrip_remap_vars::num_wts
integer(scrip_i4), save num_wts
Definition: scrip_remap_vars.f:66
scrip_grids::grid2_frac
real(scrip_r8), dimension(:), allocatable, target, save grid2_frac
Definition: scrip_grids.f:103
scrip_errormod::scrip_errorcheck
logical(scrip_logical) function, public scrip_errorcheck(errorCode, rtnName, errorMsg)
Definition: scrip_errormod.f90:143
scrip_grids::grid2_center_lon
real(scrip_r8), dimension(:), allocatable, target, save grid2_center_lon
Definition: scrip_grids.f:103