source: CPL/oasis3-mct/branches/OASIS3-MCT_2.0_branch/lib/scrip/src/remap_write.f @ 4775

Last change on this file since 4775 was 4775, checked in by aclsce, 4 years ago
  • Imported oasis3-mct from Cerfacs svn server (not suppotred anymore).

The version has been extracted from https://oasis3mct.cerfacs.fr/svn/branches/OASIS3-MCT_2.0_branch/oasis3-mct@1818

File size: 53.3 KB
Line 
1C****
2C                    ************************
3C                    *     OASIS MODULE     *
4C                    *     ------------     *
5C                    ************************
6C**** 
7C***********************************************************************
8C     This module belongs to the SCRIP library. It is modified to run
9C     within OASIS. 
10C     Modifications:
11C       - Added CASE for GAUSWGT
12C
13C     Modified by            D. Declat,  CERFACS              27.06.2002
14C***********************************************************************
15!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16!
17!     This module contains routines for writing the remapping data to
18!     a file.  Before writing the data for each mapping, the links are
19!     sorted by destination grid address.
20!
21!-----------------------------------------------------------------------
22!
23!     CVS:$Id: remap_write.f 2826 2010-12-10 11:14:21Z valcke $
24!
25!     Copyright (c) 1997, 1998 the Regents of the University of
26!       California.
27!
28!     This software and ancillary information (herein called software)
29!     called SCRIP is made available under the terms described here. 
30!     The software has been approved for release with associated
31!     LA-CC Number 98-45.
32!
33!     Unless otherwise indicated, this software has been authored
34!     by an employee or employees of the University of California,
35!     operator of the Los Alamos National Laboratory under Contract
36!     No. W-7405-ENG-36 with the U.S. Department of Energy.  The U.S.
37!     Government has rights to use, reproduce, and distribute this
38!     software.  The public may copy and use this software without
39!     charge, provided that this Notice and any statement of authorship
40!     are reproduced on all copies.  Neither the Government nor the
41!     University makes any warranty, express or implied, or assumes
42!     any liability or responsibility for the use of this software.
43!
44!     If software is modified to produce derivative works, such modified
45!     software should be clearly marked, so as not to confuse it with
46!     the version available from Los Alamos National Laboratory.
47!
48!***********************************************************************
49
50      module remap_write
51
52!-----------------------------------------------------------------------
53
54      use kinds_mod     ! defines common data types
55      use constants     ! defines common scalar constants
56      use grids         ! module containing grid information
57      use remap_vars    ! module containing remap information
58      use netcdf_mod    ! module with netCDF stuff
59      USE mod_oasis_flush
60
61      implicit none
62
63!-----------------------------------------------------------------------
64!
65!     module variables
66!
67!-----------------------------------------------------------------------
68
69      character(char_len), private :: 
70     &   map_method       ! character string for map_type
71     &,  normalize_opt    ! character string for normalization option
72     &,  history          ! character string for history information
73     &,  convention       ! character string for output convention
74
75      character(8), private :: 
76     &   cdate            ! character date string
77
78      integer (kind=int_kind), dimension(:), allocatable, private ::
79     &   src_mask_int     ! integer masks to determine
80     &,  dst_mask_int     ! cells that participate in map
81
82!-----------------------------------------------------------------------
83!
84!     various netCDF identifiers used by output routines
85!
86!-----------------------------------------------------------------------
87
88      integer (kind=int_kind), private ::
89     &   ncstat               ! error flag for netCDF calls
90     &,  nc_file_id           ! id for netCDF file
91     &,  nc_srcgrdsize_id     ! id for source grid size
92     &,  nc_dstgrdsize_id     ! id for destination grid size
93     &,  nc_srcgrdcorn_id     ! id for number of source grid corners
94     &,  nc_dstgrdcorn_id     ! id for number of dest grid corners
95     &,  nc_srcgrdrank_id     ! id for source grid rank
96     &,  nc_dstgrdrank_id     ! id for dest grid rank
97     &,  nc_numlinks_id       ! id for number of links in mapping
98     &,  nc_numwgts_id        ! id for number of weights for mapping
99     &,  nc_srcgrddims_id     ! id for source grid dimensions
100     &,  nc_dstgrddims_id     ! id for dest grid dimensions
101     &,  nc_srcgrdcntrlat_id  ! id for source grid center latitude
102     &,  nc_dstgrdcntrlat_id  ! id for dest grid center latitude
103     &,  nc_srcgrdcntrlon_id  ! id for source grid center longitude
104     &,  nc_dstgrdcntrlon_id  ! id for dest grid center longitude
105     &,  nc_srcgrdimask_id    ! id for source grid mask
106      integer (kind=int_kind), private ::
107     &   nc_dstgrdimask_id    ! id for dest grid mask
108     &,  nc_srcgrdcrnrlat_id  ! id for latitude of source grid corners
109     &,  nc_srcgrdcrnrlon_id  ! id for longitude of source grid corners
110     &,  nc_dstgrdcrnrlat_id  ! id for latitude of dest grid corners
111     &,  nc_dstgrdcrnrlon_id  ! id for longitude of dest grid corners
112     &,  nc_srcgrdarea_id     ! id for area of source grid cells
113     &,  nc_dstgrdarea_id     ! id for area of dest grid cells
114     &,  nc_srcgrdfrac_id     ! id for area fraction on source grid
115     &,  nc_dstgrdfrac_id     ! id for area fraction on dest grid
116     &,  nc_srcadd_id         ! id for map source address
117     &,  nc_dstadd_id         ! id for map destination address
118     &,  nc_rmpmatrix_id      ! id for remapping matrix
119
120      integer (kind=int_kind), dimension(2), private ::
121     &   nc_dims2_id  ! netCDF ids for 2d array dims
122
123!***********************************************************************
124
125      contains
126
127!***********************************************************************
128
129      subroutine write_remap(map1_name, map2_name, 
130     &                       interp_file1, interp_file2, output_opt)
131
132!-----------------------------------------------------------------------
133!
134!     calls correct output routine based on output format choice
135!
136!-----------------------------------------------------------------------
137
138!-----------------------------------------------------------------------
139!
140!     input variables
141!
142!-----------------------------------------------------------------------
143
144      character(char_len), intent(in) ::
145     &            map1_name,    ! name for mapping grid1 to grid2
146     &            map2_name,    ! name for mapping grid2 to grid1
147     &            interp_file1, ! filename for map1 remap data
148     &            interp_file2, ! filename for map2 remap data
149     &            output_opt    ! option for output conventions
150
151!-----------------------------------------------------------------------
152!
153!     local variables
154!
155!-----------------------------------------------------------------------
156!-----------------------------------------------------------------------
157!
158!     define some common variables to be used in all routines
159!
160!-----------------------------------------------------------------------
161      IF (nlogprt .GE. 2) THEN
162         WRITE (UNIT = nulou,FMT = *)' '
163         WRITE (UNIT = nulou,FMT = *)'Entering routine write_remap '
164         CALL OASIS_FLUSH_SCRIP(nulou)
165      ENDIF
166      select case(norm_opt)
167      case (norm_opt_none)
168        normalize_opt = 'none'
169      case (norm_opt_frcarea)
170        normalize_opt = 'fracarea'
171      case (norm_opt_dstarea)
172        normalize_opt = 'destarea'
173      case (norm_opt_nonorm)
174        normalize_opt = 'no norm'
175      end select
176      select case(map_type)
177      case(map_type_conserv)
178        map_method = 'Conservative remapping'
179      case(map_type_bilinear)
180        map_method = 'Bilinear remapping'
181      case(map_type_distwgt)
182        map_method = 'Distance weighted avg of nearest neighbors'
183      case(map_type_gauswgt)
184        map_method = 'Gaussian weighted avg of nearest neighbors'
185      case(map_type_bicubic)
186        map_method = 'Bicubic remapping'
187      case default
188        stop 'Invalid Map Type'
189      end select
190
191      call date_and_time(date=cdate)
192      write (history,1000) cdate(5:6),cdate(7:8),cdate(1:4)
193 1000 format('Created: ',a2,'-',a2,'-',a4)
194
195!-----------------------------------------------------------------------
196!
197!     call appropriate output routine
198!
199!-----------------------------------------------------------------------
200
201      select case(output_opt)
202      case ('scrip')
203        call write_remap_scrip(map1_name, interp_file1, 1)
204      case ('ncar-csm')
205        call write_remap_csm  (map1_name, interp_file1, 1)
206      case default
207        stop 'unknown output file convention'
208      end select
209
210!-----------------------------------------------------------------------
211!
212!     call appropriate output routine for second mapping if required
213!
214!-----------------------------------------------------------------------
215
216      if (num_maps > 1) then
217        select case(output_opt)
218        case ('scrip')
219          call write_remap_scrip(map2_name, interp_file2, 2)
220        case ('ncar-csm')
221          call write_remap_csm  (map2_name, interp_file2, 2)
222        case default
223          stop 'unknown output file convention'
224        end select
225      endif
226
227!-----------------------------------------------------------------------
228      IF (nlogprt .GE. 2) THEN
229         WRITE (UNIT = nulou,FMT = *)' '
230         WRITE (UNIT = nulou,FMT = *)'Leaving routine write_remap '
231         CALL OASIS_FLUSH_SCRIP(nulou)
232      ENDIF
233
234      end subroutine write_remap
235
236!***********************************************************************
237
238      subroutine write_remap_scrip(map_name, interp_file, direction)
239
240!-----------------------------------------------------------------------
241!
242!     writes remap data to a netCDF file using SCRIP conventions
243!
244!-----------------------------------------------------------------------
245
246!-----------------------------------------------------------------------
247!
248!     input variables
249!
250!-----------------------------------------------------------------------
251
252      character(char_len), intent(in) ::
253     &            map_name     ! name for mapping
254     &,           interp_file  ! filename for remap data
255
256      integer (kind=int_kind), intent(in) ::
257     &  direction              ! direction of map (1=grid1 to grid2
258                               !                   2=grid2 to grid1)
259
260!-----------------------------------------------------------------------
261!
262!     local variables
263!
264!-----------------------------------------------------------------------
265
266      character(char_len) ::
267     &  grid1_ctmp        ! character temp for grid1 names
268     &, grid2_ctmp        ! character temp for grid2 names
269
270      integer (kind=int_kind) ::
271     &  itmp1             ! integer temp
272     &, itmp2             ! integer temp
273     &, itmp3             ! integer temp
274     &, itmp4             ! integer temp
275
276      integer (kind=int_kind) :: il_nftype
277!-----------------------------------------------------------------------
278      IF (nlogprt .GE. 2) THEN
279       WRITE (UNIT = nulou,FMT = *)' '
280       WRITE (UNIT = nulou,FMT = *)'Entering routine write_remap_scrip '
281       CALL OASIS_FLUSH_SCRIP(nulou)
282      ENDIF
283!-----------------------------------------------------------------------
284!
285!     create netCDF file for mapping and define some global attributes
286!
287!-----------------------------------------------------------------------
288      ncstat = nf_create (interp_file, NF_CLOBBER, nc_file_id)
289      call netcdf_error_handler(ncstat)
290      !***
291      !*** map name
292      !***
293      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'title',
294     &                          len_trim(map_name), map_name)
295      call netcdf_error_handler(ncstat)
296
297      !***
298      !*** normalization option
299      !***
300      ncstat = nf_put_att_text(nc_file_id, NF_GLOBAL, 'normalization',
301     &                         len_trim(normalize_opt), normalize_opt)
302      call netcdf_error_handler(ncstat)
303
304      !***
305      !*** map method
306      !***
307      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'map_method',
308     &                          len_trim(map_method), map_method)
309      call netcdf_error_handler(ncstat)
310
311      !***
312      !*** history
313      !***
314      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'history',
315     &                          len_trim(history), history)
316      call netcdf_error_handler(ncstat)
317
318      !***
319      !*** file convention
320      !***
321      convention = 'SCRIP'
322      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'conventions',
323     &                          len_trim(convention), convention)
324      call netcdf_error_handler(ncstat)
325
326      !***
327      !*** source and destination grid names
328      !***
329
330      if (direction == 1) then
331        grid1_ctmp = 'source_grid'
332        grid2_ctmp = 'dest_grid'
333      else
334        grid1_ctmp = 'dest_grid'
335        grid2_ctmp = 'source_grid'
336      endif
337
338      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid1_ctmp),
339     &                          len_trim(grid1_name), grid1_name)
340      call netcdf_error_handler(ncstat)
341      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid2_ctmp),
342     &                          len_trim(grid2_name), grid2_name)
343      call netcdf_error_handler(ncstat)
344
345!-----------------------------------------------------------------------
346!
347!     prepare netCDF dimension info
348!
349!-----------------------------------------------------------------------
350
351      !***
352      !*** define grid size dimensions
353      !***
354
355      if (direction == 1) then
356        itmp1 = grid1_size
357        itmp2 = grid2_size
358      else
359        itmp1 = grid2_size
360        itmp2 = grid1_size
361      endif
362
363      ncstat = nf_def_dim (nc_file_id, 'src_grid_size', itmp1, 
364     &                     nc_srcgrdsize_id)
365      call netcdf_error_handler(ncstat)
366
367      ncstat = nf_def_dim (nc_file_id, 'dst_grid_size', itmp2, 
368     &                     nc_dstgrdsize_id)
369      call netcdf_error_handler(ncstat)
370
371      !***
372      !*** define grid corner dimension
373      !***
374
375      if (direction == 1) then
376        itmp1 = grid1_corners
377        itmp2 = grid2_corners
378      else
379        itmp1 = grid2_corners
380        itmp2 = grid1_corners
381      endif
382
383      ncstat = nf_def_dim (nc_file_id, 'src_grid_corners', 
384     &                     itmp1, nc_srcgrdcorn_id)
385      call netcdf_error_handler(ncstat)
386
387      ncstat = nf_def_dim (nc_file_id, 'dst_grid_corners', 
388     &                     itmp2, nc_dstgrdcorn_id)
389      call netcdf_error_handler(ncstat)
390      !***
391      !*** define grid rank dimension
392      !***
393
394      if (direction == 1) then
395        itmp1 = grid1_rank
396        itmp2 = grid2_rank
397      else
398        itmp1 = grid2_rank
399        itmp2 = grid1_rank
400      endif
401
402      ncstat = nf_def_dim (nc_file_id, 'src_grid_rank', 
403     &                     itmp1, nc_srcgrdrank_id)
404      call netcdf_error_handler(ncstat)
405
406      ncstat = nf_def_dim (nc_file_id, 'dst_grid_rank', 
407     &                     itmp2, nc_dstgrdrank_id)
408      call netcdf_error_handler(ncstat)
409      !***
410      !*** define map size dimensions
411      !***
412
413      if (direction == 1) then
414        itmp1 = num_links_map1
415      else
416        itmp1 = num_links_map2
417      endif
418
419      ncstat = nf_def_dim (nc_file_id, 'num_links', 
420     &                     itmp1, nc_numlinks_id)
421      call netcdf_error_handler(ncstat)
422      ncstat = nf_def_dim (nc_file_id, 'num_wgts', 
423     &                     num_wts, nc_numwgts_id)
424      call netcdf_error_handler(ncstat)
425      !***
426      !*** define grid dimensions
427      !***
428
429      ncstat = nf_def_var (nc_file_id, 'src_grid_dims', NF_INT,
430     &                     1, nc_srcgrdrank_id, nc_srcgrddims_id)
431      call netcdf_error_handler(ncstat)
432      ncstat = nf_def_var (nc_file_id, 'dst_grid_dims', NF_INT,
433     &                     1, nc_dstgrdrank_id, nc_dstgrddims_id)
434      call netcdf_error_handler(ncstat)
435!-----------------------------------------------------------------------
436!
437!     define all arrays for netCDF descriptors
438!
439!-----------------------------------------------------------------------
440
441      !***
442      !*** define grid center latitude array
443      !***
444
445      il_nftype = NF_DOUBLE
446      IF (ll_single) il_nftype = NF_REAL
447
448      ncstat = nf_def_var (nc_file_id, 'src_grid_center_lat', 
449     &                     il_nftype, 1, nc_srcgrdsize_id, 
450     &                     nc_srcgrdcntrlat_id)
451      call netcdf_error_handler(ncstat)
452      ncstat = nf_def_var (nc_file_id, 'dst_grid_center_lat', 
453     &                     il_nftype, 1, nc_dstgrdsize_id, 
454     &                     nc_dstgrdcntrlat_id)
455      call netcdf_error_handler(ncstat)
456      !***
457      !*** define grid center longitude array
458      !***
459
460      ncstat = nf_def_var (nc_file_id, 'src_grid_center_lon', 
461     &                     il_nftype, 1, nc_srcgrdsize_id, 
462     &                     nc_srcgrdcntrlon_id)
463      call netcdf_error_handler(ncstat)
464      ncstat = nf_def_var (nc_file_id, 'dst_grid_center_lon', 
465     &                     il_nftype, 1, nc_dstgrdsize_id, 
466     &                     nc_dstgrdcntrlon_id)
467      call netcdf_error_handler(ncstat)
468      !***
469      !*** define grid corner lat/lon arrays
470      !***
471
472      nc_dims2_id(1) = nc_srcgrdcorn_id
473      nc_dims2_id(2) = nc_srcgrdsize_id
474
475      ncstat = nf_def_var (nc_file_id, 'src_grid_corner_lat', 
476     &                     il_nftype, 2, nc_dims2_id, 
477     &                     nc_srcgrdcrnrlat_id)
478      call netcdf_error_handler(ncstat)
479      ncstat = nf_def_var (nc_file_id, 'src_grid_corner_lon', 
480     &                     il_nftype, 2, nc_dims2_id, 
481     &                     nc_srcgrdcrnrlon_id)
482      call netcdf_error_handler(ncstat)
483
484      nc_dims2_id(1) = nc_dstgrdcorn_id
485      nc_dims2_id(2) = nc_dstgrdsize_id
486
487      ncstat = nf_def_var (nc_file_id, 'dst_grid_corner_lat', 
488     &                     il_nftype, 2, nc_dims2_id, 
489     &                     nc_dstgrdcrnrlat_id)
490      call netcdf_error_handler(ncstat)
491      ncstat = nf_def_var (nc_file_id, 'dst_grid_corner_lon', 
492     &                     il_nftype, 2, nc_dims2_id, 
493     &                     nc_dstgrdcrnrlon_id)
494      call netcdf_error_handler(ncstat)
495
496      !***
497      !*** define units for all coordinate arrays
498      !***
499
500      if (direction == 1) then
501        grid1_ctmp = grid1_units
502        grid2_ctmp = grid2_units
503      else
504        grid1_ctmp = grid2_units
505        grid2_ctmp = grid1_units
506      endif
507
508      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlat_id, 
509     &                          'units', 7, grid1_ctmp)
510      call netcdf_error_handler(ncstat)
511      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlat_id, 
512     &                          'units', 7, grid2_ctmp)
513      call netcdf_error_handler(ncstat)
514
515      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlon_id, 
516     &                          'units', 7, grid1_ctmp)
517      call netcdf_error_handler(ncstat)
518
519      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlon_id, 
520     &                          'units', 7, grid2_ctmp)
521      call netcdf_error_handler(ncstat)
522
523      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlat_id, 
524     &                          'units', 7, grid1_ctmp)
525      call netcdf_error_handler(ncstat)
526      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlon_id, 
527     &                          'units', 7, grid1_ctmp)
528      call netcdf_error_handler(ncstat)
529
530      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlat_id, 
531     &                          'units', 7, grid2_ctmp)
532      call netcdf_error_handler(ncstat)
533
534      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlon_id, 
535     &                          'units', 7, grid2_ctmp)
536      call netcdf_error_handler(ncstat)
537
538      !***
539      !*** define grid mask
540      !***
541
542      ncstat = nf_def_var (nc_file_id, 'src_grid_imask', NF_INT,
543     &                     1, nc_srcgrdsize_id, nc_srcgrdimask_id)
544      call netcdf_error_handler(ncstat)
545
546      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdimask_id, 
547     &                          'units', 8, 'unitless')
548      call netcdf_error_handler(ncstat)
549      ncstat = nf_def_var (nc_file_id, 'dst_grid_imask', NF_INT,
550     &                     1, nc_dstgrdsize_id, nc_dstgrdimask_id)
551      call netcdf_error_handler(ncstat)
552
553      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdimask_id, 
554     &                          'units', 8, 'unitless')
555      call netcdf_error_handler(ncstat)
556      !***
557      !*** define grid area arrays
558      !***
559
560      ncstat = nf_def_var (nc_file_id, 'src_grid_area', 
561     &                     il_nftype, 1, nc_srcgrdsize_id, 
562     &                     nc_srcgrdarea_id)
563      call netcdf_error_handler(ncstat)
564      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdarea_id, 
565     &                          'units', 14, 'square radians')
566      call netcdf_error_handler(ncstat)
567
568      ncstat = nf_def_var (nc_file_id, 'dst_grid_area', 
569     &                     il_nftype, 1, nc_dstgrdsize_id, 
570     &                     nc_dstgrdarea_id)
571      call netcdf_error_handler(ncstat)
572
573      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdarea_id, 
574     &                          'units', 14, 'square radians')
575      call netcdf_error_handler(ncstat)
576
577      !***
578      !*** define grid fraction arrays
579      !***
580
581      ncstat = nf_def_var (nc_file_id, 'src_grid_frac', 
582     &                     il_nftype, 1, nc_srcgrdsize_id, 
583     &                     nc_srcgrdfrac_id)
584      call netcdf_error_handler(ncstat)
585      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdfrac_id, 
586     &                          'units', 8, 'unitless')
587      call netcdf_error_handler(ncstat)
588
589      ncstat = nf_def_var (nc_file_id, 'dst_grid_frac', 
590     &                     il_nftype, 1, nc_dstgrdsize_id, 
591     &                     nc_dstgrdfrac_id)
592      call netcdf_error_handler(ncstat)
593
594      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdfrac_id, 
595     &                          'units', 8, 'unitless')
596      call netcdf_error_handler(ncstat)
597      !***
598      !*** define mapping arrays
599      !***
600
601      ncstat = nf_def_var (nc_file_id, 'src_address', 
602     &                     NF_INT, 1, nc_numlinks_id, 
603     &                     nc_srcadd_id)
604      call netcdf_error_handler(ncstat)
605
606      ncstat = nf_def_var (nc_file_id, 'dst_address', 
607     &                     NF_INT, 1, nc_numlinks_id, 
608     &                     nc_dstadd_id)
609      call netcdf_error_handler(ncstat)
610
611      nc_dims2_id(1) = nc_numwgts_id
612      nc_dims2_id(2) = nc_numlinks_id
613
614      ncstat = nf_def_var (nc_file_id, 'remap_matrix', 
615     &                     il_nftype, 2, nc_dims2_id, 
616     &                     nc_rmpmatrix_id)
617      call netcdf_error_handler(ncstat)
618
619      !***
620      !*** end definition stage
621      !***
622
623      ncstat = nf_enddef(nc_file_id)
624      call netcdf_error_handler(ncstat)
625!-----------------------------------------------------------------------
626!
627!     compute integer masks
628!
629!-----------------------------------------------------------------------
630
631      if (direction == 1) then
632        allocate (src_mask_int(grid1_size),
633     &            dst_mask_int(grid2_size))
634
635        where (grid2_mask)
636          dst_mask_int = 1
637        elsewhere
638          dst_mask_int = 0
639        endwhere
640
641        where (grid1_mask)
642          src_mask_int = 1
643        elsewhere
644          src_mask_int = 0
645        endwhere
646      else
647        allocate (src_mask_int(grid2_size),
648     &            dst_mask_int(grid1_size))
649
650        where (grid1_mask)
651          dst_mask_int = 1
652        elsewhere
653          dst_mask_int = 0
654        endwhere
655
656        where (grid2_mask)
657          src_mask_int = 1
658        elsewhere
659          src_mask_int = 0
660        endwhere
661      endif
662
663!-----------------------------------------------------------------------
664!
665!     change units of lat/lon coordinates if input units different
666!     from radians
667!
668!-----------------------------------------------------------------------
669
670      if (grid1_units(1:7) == 'degrees' .and. direction == 1) then
671        grid1_center_lat = grid1_center_lat/deg2rad
672        grid1_center_lon = grid1_center_lon/deg2rad
673        grid1_corner_lat = grid1_corner_lat/deg2rad
674        grid1_corner_lon = grid1_corner_lon/deg2rad
675      endif
676
677      if (grid2_units(1:7) == 'degrees' .and. direction == 1) then
678        grid2_center_lat = grid2_center_lat/deg2rad
679        grid2_center_lon = grid2_center_lon/deg2rad
680        grid2_corner_lat = grid2_corner_lat/deg2rad
681        grid2_corner_lon = grid2_corner_lon/deg2rad
682      endif
683
684!-----------------------------------------------------------------------
685!
686!     write mapping data
687!
688!-----------------------------------------------------------------------
689
690      if (direction == 1) then
691        itmp1 = nc_srcgrddims_id
692        itmp2 = nc_dstgrddims_id
693      else
694        itmp2 = nc_srcgrddims_id
695        itmp1 = nc_dstgrddims_id
696      endif
697
698      ncstat = nf_put_var_int(nc_file_id, itmp1, grid1_dims)
699      call netcdf_error_handler(ncstat)
700
701      ncstat = nf_put_var_int(nc_file_id, itmp2, grid2_dims)
702      call netcdf_error_handler(ncstat)
703
704      ncstat = nf_put_var_int(nc_file_id, nc_srcgrdimask_id, 
705     &                        src_mask_int)
706      call netcdf_error_handler(ncstat)
707      ncstat = nf_put_var_int(nc_file_id, nc_dstgrdimask_id,
708     &                        dst_mask_int)
709      call netcdf_error_handler(ncstat)
710
711      deallocate(src_mask_int, dst_mask_int)
712
713      if (direction == 1) then
714        itmp1 = nc_srcgrdcntrlat_id
715        itmp2 = nc_srcgrdcntrlon_id
716        itmp3 = nc_srcgrdcrnrlat_id
717        itmp4 = nc_srcgrdcrnrlon_id
718      else
719        itmp1 = nc_dstgrdcntrlat_id
720        itmp2 = nc_dstgrdcntrlon_id
721        itmp3 = nc_dstgrdcrnrlat_id
722        itmp4 = nc_dstgrdcrnrlon_id
723      endif
724      IF (ll_single) THEN
725          ncstat=nf_put_var_real(nc_file_id, itmp1, grid1_center_lat)
726      ELSE
727          ncstat=nf_put_var_double(nc_file_id, itmp1, grid1_center_lat)
728      ENDIF
729      call netcdf_error_handler(ncstat)
730
731      IF (ll_single) THEN
732          ncstat=nf_put_var_real(nc_file_id, itmp2, grid1_center_lon)
733      ELSE
734          ncstat=nf_put_var_double(nc_file_id, itmp2, grid1_center_lon)
735      ENDIF
736      call netcdf_error_handler(ncstat)
737
738      if (.not. luse_grid_centers) THEN
739          IF (ll_single) THEN
740              ncstat = nf_put_var_real(nc_file_id, itmp3, 
741     $            grid1_corner_lat)
742          ELSE
743              ncstat = nf_put_var_double(nc_file_id, itmp3, 
744     $            grid1_corner_lat)
745          ENDIF
746          call netcdf_error_handler(ncstat)
747          IF (ll_single) THEN
748              ncstat = nf_put_var_real(nc_file_id, itmp4, 
749     $            grid1_corner_lon)
750          ELSE
751              ncstat = nf_put_var_double(nc_file_id, itmp4, 
752     $            grid1_corner_lon)
753          ENDIF
754          call netcdf_error_handler(ncstat)
755      endif
756
757      if (direction == 1) then
758        itmp1 = nc_dstgrdcntrlat_id
759        itmp2 = nc_dstgrdcntrlon_id
760        itmp3 = nc_dstgrdcrnrlat_id
761        itmp4 = nc_dstgrdcrnrlon_id
762      else
763        itmp1 = nc_srcgrdcntrlat_id
764        itmp2 = nc_srcgrdcntrlon_id
765        itmp3 = nc_srcgrdcrnrlat_id
766        itmp4 = nc_srcgrdcrnrlon_id
767      endif
768     
769      IF (ll_single) THEN
770          ncstat=nf_put_var_real(nc_file_id, itmp1, grid2_center_lat)
771      ELSE
772          ncstat=nf_put_var_double(nc_file_id, itmp1, grid2_center_lat)
773      ENDIF
774      call netcdf_error_handler(ncstat)
775
776      IF (ll_single) THEN
777          ncstat=nf_put_var_real(nc_file_id, itmp2, grid2_center_lon)
778      ELSE
779          ncstat=nf_put_var_double(nc_file_id, itmp2, grid2_center_lon)
780      ENDIF
781      call netcdf_error_handler(ncstat)
782
783      if (.not. luse_grid_centers) THEN
784          IF (ll_single) THEN
785              ncstat = nf_put_var_real(nc_file_id, itmp3, 
786     $            grid2_corner_lat)
787          ELSE
788              ncstat = nf_put_var_double(nc_file_id, itmp3, 
789     $            grid2_corner_lat)
790          ENDIF
791          call netcdf_error_handler(ncstat)
792
793          IF (ll_single) THEN
794              ncstat = nf_put_var_real(nc_file_id, itmp4, 
795     $            grid2_corner_lon)
796          ELSE
797              ncstat = nf_put_var_double(nc_file_id, itmp4, 
798     $            grid2_corner_lon)
799          ENDIF
800          call netcdf_error_handler(ncstat)
801
802      endif
803      if (direction == 1) then
804        itmp1 = nc_srcgrdarea_id
805        itmp2 = nc_srcgrdfrac_id
806        itmp3 = nc_dstgrdarea_id
807        itmp4 = nc_dstgrdfrac_id
808      else
809        itmp1 = nc_dstgrdarea_id
810        itmp2 = nc_dstgrdfrac_id
811        itmp3 = nc_srcgrdarea_id
812        itmp4 = nc_srcgrdfrac_id
813      endif
814
815      if (luse_grid1_area) THEN
816          IF (ll_single) THEN
817              ncstat=nf_put_var_real(nc_file_id, itmp1, grid1_area_in)
818          ELSE
819              ncstat=nf_put_var_double(nc_file_id, itmp1, grid1_area_in)
820          ENDIF
821      ELSE
822          IF (ll_single) THEN
823              ncstat = nf_put_var_real(nc_file_id, itmp1, grid1_area)
824          ELSE
825              ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_area)
826          ENDIF
827      endif
828      call netcdf_error_handler(ncstat)
829
830      IF (ll_single) THEN
831          ncstat = nf_put_var_real(nc_file_id, itmp2, grid1_frac)
832      ELSE
833          ncstat = nf_put_var_double(nc_file_id, itmp2, grid1_frac)
834      ENDIF
835      call netcdf_error_handler(ncstat)
836
837      if (luse_grid2_area) THEN
838          IF (ll_single) THEN
839              ncstat = nf_put_var_real
840     $            (nc_file_id, itmp3, grid2_area_in)
841          ELSE
842              ncstat = nf_put_var_double
843     $            (nc_file_id, itmp3, grid2_area_in)
844          ENDIF
845      ELSE
846          IF (ll_single) THEN
847              ncstat = nf_put_var_real(nc_file_id, itmp3, grid2_area)
848          ELSE
849              ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_area)
850          ENDIF
851      endif
852      call netcdf_error_handler(ncstat)
853
854      IF (ll_single) THEN
855          ncstat = nf_put_var_real(nc_file_id, itmp4, grid2_frac)
856      ELSE
857          ncstat = nf_put_var_double(nc_file_id, itmp4, grid2_frac)
858      ENDIF
859      call netcdf_error_handler(ncstat)
860      if (direction == 1) then
861        ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id, 
862     &                          grid1_add_map1)
863        call netcdf_error_handler(ncstat)
864
865        ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id, 
866     &                          grid2_add_map1)
867        call netcdf_error_handler(ncstat)
868        IF (ll_single) THEN
869            ncstat = nf_put_var_real(nc_file_id, nc_rmpmatrix_id, 
870     &          wts_map1)
871        ELSE
872            ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id, 
873     &          wts_map1)
874        ENDIF
875        call netcdf_error_handler(ncstat)
876      else
877        ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id, 
878     &                          grid2_add_map2)
879        call netcdf_error_handler(ncstat)
880
881        ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id, 
882     &                          grid1_add_map2)
883        call netcdf_error_handler(ncstat)
884        IF (ll_single) THEN
885            ncstat = nf_put_var_real(nc_file_id, nc_rmpmatrix_id, 
886     &          wts_map2)
887        ELSE
888            ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id, 
889     &          wts_map2)
890        ENDIF
891        call netcdf_error_handler(ncstat)
892      endif
893
894      ncstat = nf_close(nc_file_id)
895      call netcdf_error_handler(ncstat)
896!
897      IF (nlogprt .GE. 2) THEN
898       WRITE (UNIT = nulou,FMT = *)' '
899       WRITE (UNIT = nulou,FMT = *)'Leaving routine write_remap_scrip '
900       CALL OASIS_FLUSH_SCRIP(nulou)
901      ENDIF
902!-----------------------------------------------------------------------
903
904      end subroutine write_remap_scrip
905
906!***********************************************************************
907
908      subroutine write_remap_csm(map_name, interp_file, direction)
909
910!-----------------------------------------------------------------------
911!
912!     writes remap data to a netCDF file using NCAR-CSM conventions
913!
914!-----------------------------------------------------------------------
915
916!-----------------------------------------------------------------------
917!
918!     input variables
919!
920!-----------------------------------------------------------------------
921
922      character(char_len), intent(in) ::
923     &            map_name     ! name for mapping
924     &,           interp_file  ! filename for remap data
925
926      integer (kind=int_kind), intent(in) ::
927     &  direction              ! direction of map (1=grid1 to grid2
928                               !                   2=grid2 to grid1)
929
930!-----------------------------------------------------------------------
931!
932!     local variables
933!
934!-----------------------------------------------------------------------
935
936      character(char_len) ::
937     &  grid1_ctmp        ! character temp for grid1 names
938     &, grid2_ctmp        ! character temp for grid2 names
939
940      integer (kind=int_kind) ::
941     &  itmp1             ! integer temp
942     &, itmp2             ! integer temp
943     &, itmp3             ! integer temp
944     &, itmp4             ! integer temp
945     &, nc_numwgts1_id    ! extra netCDF id for additional weights
946     &, nc_src_isize_id   ! extra netCDF id for ni_a
947     &, nc_src_jsize_id   ! extra netCDF id for nj_a
948     &, nc_dst_isize_id   ! extra netCDF id for ni_b
949     &, nc_dst_jsize_id   ! extra netCDF id for nj_b
950     &, nc_rmpmatrix2_id  ! extra netCDF id for high-order remap matrix
951
952      integer (kind=int_kind) :: il_nftype
953
954      real (kind=dbl_kind), dimension(:),allocatable ::
955     &  wts1              ! CSM wants single array for 1st-order wts
956
957      real (kind=dbl_kind), dimension(:,:),allocatable ::
958     &  wts2              ! write remaining weights in different array
959
960!-----------------------------------------------------------------------
961!
962!     create netCDF file for mapping and define some global attributes
963!
964!-----------------------------------------------------------------------
965
966      ncstat = nf_create (interp_file, NF_CLOBBER, nc_file_id)
967      call netcdf_error_handler(ncstat)
968
969      !***
970      !*** map name
971      !***
972      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'title',
973     &                          len_trim(map_name), map_name)
974      call netcdf_error_handler(ncstat)
975
976      !***
977      !*** normalization option
978      !***
979      ncstat = nf_put_att_text(nc_file_id, NF_GLOBAL, 'normalization',
980     &                         len_trim(normalize_opt), normalize_opt)
981      call netcdf_error_handler(ncstat)
982
983      !***
984      !*** map method
985      !***
986      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'map_method',
987     &                          len_trim(map_method), map_method)
988      call netcdf_error_handler(ncstat)
989
990      !***
991      !*** history
992      !***
993      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'history',
994     &                          len_trim(history), history)
995      call netcdf_error_handler(ncstat)
996
997      !***
998      !*** file convention
999      !***
1000      convention = 'NCAR-CSM'
1001      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'conventions',
1002     &                          len_trim(convention), convention)
1003      call netcdf_error_handler(ncstat)
1004
1005      !***
1006      !*** source and destination grid names
1007      !***
1008
1009      if (direction == 1) then
1010        grid1_ctmp = 'domain_a'
1011        grid2_ctmp = 'domain_b'
1012      else
1013        grid1_ctmp = 'domain_b'
1014        grid2_ctmp = 'domain_a'
1015      endif
1016
1017      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid1_ctmp),
1018     &                          len_trim(grid1_name), grid1_name)
1019      call netcdf_error_handler(ncstat)
1020
1021      ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid2_ctmp),
1022     &                          len_trim(grid2_name), grid2_name)
1023      call netcdf_error_handler(ncstat)
1024
1025!-----------------------------------------------------------------------
1026!
1027!     prepare netCDF dimension info
1028!
1029!-----------------------------------------------------------------------
1030
1031      !***
1032      !*** define grid size dimensions
1033      !***
1034
1035      if (direction == 1) then
1036        itmp1 = grid1_size
1037        itmp2 = grid2_size
1038      else
1039        itmp1 = grid2_size
1040        itmp2 = grid1_size
1041      endif
1042
1043      ncstat = nf_def_dim (nc_file_id, 'n_a', itmp1, nc_srcgrdsize_id)
1044      call netcdf_error_handler(ncstat)
1045
1046      ncstat = nf_def_dim (nc_file_id, 'n_b', itmp2, nc_dstgrdsize_id)
1047      call netcdf_error_handler(ncstat)
1048
1049      !***
1050      !*** define grid corner dimension
1051      !***
1052
1053      if (direction == 1) then
1054        itmp1 = grid1_corners
1055        itmp2 = grid2_corners
1056      else
1057        itmp1 = grid2_corners
1058        itmp2 = grid1_corners
1059      endif
1060
1061      ncstat = nf_def_dim (nc_file_id, 'nv_a', itmp1, nc_srcgrdcorn_id)
1062      call netcdf_error_handler(ncstat)
1063
1064      ncstat = nf_def_dim (nc_file_id, 'nv_b', itmp2, nc_dstgrdcorn_id)
1065      call netcdf_error_handler(ncstat)
1066
1067      !***
1068      !*** define grid rank dimension
1069      !***
1070
1071      if (direction == 1) then
1072        itmp1 = grid1_rank
1073        itmp2 = grid2_rank
1074      else
1075        itmp1 = grid2_rank
1076        itmp2 = grid1_rank
1077      endif
1078
1079      ncstat = nf_def_dim (nc_file_id, 'src_grid_rank', 
1080     &                     itmp1, nc_srcgrdrank_id)
1081      call netcdf_error_handler(ncstat)
1082
1083      ncstat = nf_def_dim (nc_file_id, 'dst_grid_rank', 
1084     &                     itmp2, nc_dstgrdrank_id)
1085      call netcdf_error_handler(ncstat)
1086
1087      !***
1088      !*** define first two dims as if 2-d cartesian domain
1089      !***
1090
1091      if (direction == 1) then
1092        itmp1 = grid1_dims(1)
1093        if (grid1_rank > 1) then
1094          itmp2 = grid1_dims(2)
1095        else
1096          itmp2 = 0
1097        endif
1098        itmp3 = grid2_dims(1)
1099        if (grid2_rank > 1) then
1100          itmp4 = grid2_dims(2)
1101        else
1102          itmp4 = 0
1103        endif
1104      else
1105        itmp1 = grid2_dims(1)
1106        if (grid2_rank > 1) then
1107          itmp2 = grid2_dims(2)
1108        else
1109          itmp2 = 0
1110        endif
1111        itmp3 = grid1_dims(1)
1112        if (grid1_rank > 1) then
1113          itmp4 = grid1_dims(2)
1114        else
1115          itmp4 = 0
1116        endif
1117      endif
1118
1119      ncstat = nf_def_dim (nc_file_id, 'ni_a', itmp1, nc_src_isize_id)
1120      call netcdf_error_handler(ncstat)
1121
1122      ncstat = nf_def_dim (nc_file_id, 'nj_a', itmp2, nc_src_jsize_id)
1123      call netcdf_error_handler(ncstat)
1124
1125      ncstat = nf_def_dim (nc_file_id, 'ni_b', itmp3, nc_dst_isize_id)
1126      call netcdf_error_handler(ncstat)
1127
1128      ncstat = nf_def_dim (nc_file_id, 'nj_b', itmp4, nc_dst_jsize_id)
1129      call netcdf_error_handler(ncstat)
1130
1131      !***
1132      !*** define map size dimensions
1133      !***
1134
1135      if (direction == 1) then
1136        itmp1 = num_links_map1
1137      else
1138        itmp1 = num_links_map2
1139      endif
1140
1141      ncstat = nf_def_dim (nc_file_id, 'n_s', itmp1, nc_numlinks_id)
1142      call netcdf_error_handler(ncstat)
1143
1144      ncstat = nf_def_dim (nc_file_id, 'num_wgts', 
1145     &                     num_wts, nc_numwgts_id)
1146      call netcdf_error_handler(ncstat)
1147
1148      if (num_wts > 1) then
1149        ncstat = nf_def_dim (nc_file_id, 'num_wgts1', 
1150     &                       num_wts-1, nc_numwgts1_id)
1151        call netcdf_error_handler(ncstat)
1152      endif
1153
1154      !***
1155      !*** define grid dimensions
1156      !***
1157
1158      ncstat = nf_def_var (nc_file_id, 'src_grid_dims', NF_INT,
1159     &                     1, nc_srcgrdrank_id, nc_srcgrddims_id)
1160      call netcdf_error_handler(ncstat)
1161
1162      ncstat = nf_def_var (nc_file_id, 'dst_grid_dims', NF_INT,
1163     &                     1, nc_dstgrdrank_id, nc_dstgrddims_id)
1164      call netcdf_error_handler(ncstat)
1165
1166!-----------------------------------------------------------------------
1167!
1168!     define all arrays for netCDF descriptors
1169!
1170!-----------------------------------------------------------------------
1171
1172      !***
1173      !*** define grid center latitude array
1174      !***
1175
1176      il_nftype = NF_DOUBLE
1177      IF (ll_single) il_nftype = NF_REAL
1178
1179      ncstat = nf_def_var (nc_file_id, 'yc_a',
1180     &                     il_nftype, 1, nc_srcgrdsize_id, 
1181     &                     nc_srcgrdcntrlat_id)
1182      call netcdf_error_handler(ncstat)
1183
1184      ncstat = nf_def_var (nc_file_id, 'yc_b', 
1185     &                     il_nftype, 1, nc_dstgrdsize_id, 
1186     &                     nc_dstgrdcntrlat_id)
1187      call netcdf_error_handler(ncstat)
1188
1189      !***
1190      !*** define grid center longitude array
1191      !***
1192
1193      ncstat = nf_def_var (nc_file_id, 'xc_a', 
1194     &                     il_nftype, 1, nc_srcgrdsize_id, 
1195     &                     nc_srcgrdcntrlon_id)
1196      call netcdf_error_handler(ncstat)
1197
1198      ncstat = nf_def_var (nc_file_id, 'xc_b', 
1199     &                     il_nftype, 1, nc_dstgrdsize_id, 
1200     &                     nc_dstgrdcntrlon_id)
1201      call netcdf_error_handler(ncstat)
1202
1203      !***
1204      !*** define grid corner lat/lon arrays
1205      !***
1206
1207      nc_dims2_id(1) = nc_srcgrdcorn_id
1208      nc_dims2_id(2) = nc_srcgrdsize_id
1209
1210      ncstat = nf_def_var (nc_file_id, 'yv_a', 
1211     &                     il_nftype, 2, nc_dims2_id, 
1212     &                     nc_srcgrdcrnrlat_id)
1213      call netcdf_error_handler(ncstat)
1214
1215      ncstat = nf_def_var (nc_file_id, 'xv_a', 
1216     &                     il_nftype, 2, nc_dims2_id, 
1217     &                     nc_srcgrdcrnrlon_id)
1218      call netcdf_error_handler(ncstat)
1219
1220      nc_dims2_id(1) = nc_dstgrdcorn_id
1221      nc_dims2_id(2) = nc_dstgrdsize_id
1222
1223      ncstat = nf_def_var (nc_file_id, 'yv_b', 
1224     &                     il_nftype, 2, nc_dims2_id, 
1225     &                     nc_dstgrdcrnrlat_id)
1226      call netcdf_error_handler(ncstat)
1227
1228      ncstat = nf_def_var (nc_file_id, 'xv_b', 
1229     &                     il_nftype, 2, nc_dims2_id, 
1230     &                     nc_dstgrdcrnrlon_id)
1231      call netcdf_error_handler(ncstat)
1232
1233      !***
1234      !*** CSM wants all in degrees
1235      !***
1236
1237      grid1_units = 'degrees'
1238      grid2_units = 'degrees'
1239
1240      if (direction == 1) then
1241        grid1_ctmp = grid1_units
1242        grid2_ctmp = grid2_units
1243      else
1244        grid1_ctmp = grid2_units
1245        grid2_ctmp = grid1_units
1246      endif
1247
1248      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlat_id, 
1249     &                          'units', 7, grid1_ctmp)
1250      call netcdf_error_handler(ncstat)
1251
1252      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlat_id, 
1253     &                          'units', 7, grid2_ctmp)
1254      call netcdf_error_handler(ncstat)
1255
1256      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlon_id, 
1257     &                          'units', 7, grid1_ctmp)
1258      call netcdf_error_handler(ncstat)
1259
1260      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlon_id, 
1261     &                          'units', 7, grid2_ctmp)
1262      call netcdf_error_handler(ncstat)
1263
1264      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlat_id, 
1265     &                          'units', 7, grid1_ctmp)
1266      call netcdf_error_handler(ncstat)
1267
1268      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlon_id, 
1269     &                          'units', 7, grid1_ctmp)
1270      call netcdf_error_handler(ncstat)
1271
1272      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlat_id, 
1273     &                          'units', 7, grid2_ctmp)
1274      call netcdf_error_handler(ncstat)
1275
1276      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlon_id, 
1277     &                          'units', 7, grid2_ctmp)
1278      call netcdf_error_handler(ncstat)
1279
1280      !***
1281      !*** define grid mask
1282      !***
1283
1284      ncstat = nf_def_var (nc_file_id, 'mask_a', NF_INT,
1285     &                     1, nc_srcgrdsize_id, nc_srcgrdimask_id)
1286      call netcdf_error_handler(ncstat)
1287
1288      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdimask_id, 
1289     &                          'units', 8, 'unitless')
1290      call netcdf_error_handler(ncstat)
1291
1292      ncstat = nf_def_var (nc_file_id, 'mask_b', NF_INT,
1293     &                     1, nc_dstgrdsize_id, nc_dstgrdimask_id)
1294      call netcdf_error_handler(ncstat)
1295
1296      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdimask_id, 
1297     &                          'units', 8, 'unitless')
1298      call netcdf_error_handler(ncstat)
1299
1300      !***
1301      !*** define grid area arrays
1302      !***
1303
1304      ncstat = nf_def_var (nc_file_id, 'area_a', 
1305     &                     il_nftype, 1, nc_srcgrdsize_id, 
1306     &                     nc_srcgrdarea_id)
1307      call netcdf_error_handler(ncstat)
1308
1309      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdarea_id, 
1310     &                          'units', 14, 'square radians')
1311      call netcdf_error_handler(ncstat)
1312
1313      ncstat = nf_def_var (nc_file_id, 'area_b', 
1314     &                     il_nftype, 1, nc_dstgrdsize_id, 
1315     &                     nc_dstgrdarea_id)
1316      call netcdf_error_handler(ncstat)
1317
1318      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdarea_id, 
1319     &                          'units', 14, 'square radians')
1320      call netcdf_error_handler(ncstat)
1321
1322      !***
1323      !*** define grid fraction arrays
1324      !***
1325
1326      ncstat = nf_def_var (nc_file_id, 'frac_a', 
1327     &                     il_nftype, 1, nc_srcgrdsize_id, 
1328     &                     nc_srcgrdfrac_id)
1329      call netcdf_error_handler(ncstat)
1330
1331      ncstat = nf_put_att_text (nc_file_id, nc_srcgrdfrac_id, 
1332     &                          'units', 8, 'unitless')
1333      call netcdf_error_handler(ncstat)
1334
1335      ncstat = nf_def_var (nc_file_id, 'frac_b', 
1336     &                     il_nftype, 1, nc_dstgrdsize_id, 
1337     &                     nc_dstgrdfrac_id)
1338      call netcdf_error_handler(ncstat)
1339
1340      ncstat = nf_put_att_text (nc_file_id, nc_dstgrdfrac_id, 
1341     &                          'units', 8, 'unitless')
1342      call netcdf_error_handler(ncstat)
1343
1344      !***
1345      !*** define mapping arrays
1346      !***
1347
1348      ncstat = nf_def_var (nc_file_id, 'col', 
1349     &                     NF_INT, 1, nc_numlinks_id, 
1350     &                     nc_srcadd_id)
1351      call netcdf_error_handler(ncstat)
1352
1353      ncstat = nf_def_var (nc_file_id, 'row', 
1354     &                     NF_INT, 1, nc_numlinks_id, 
1355     &                     nc_dstadd_id)
1356      call netcdf_error_handler(ncstat)
1357
1358      ncstat = nf_def_var (nc_file_id, 'S', 
1359     &                     il_nftype, 1, nc_numlinks_id, 
1360     &                     nc_rmpmatrix_id)
1361      call netcdf_error_handler(ncstat)
1362
1363      if (num_wts > 1) then
1364        nc_dims2_id(1) = nc_numwgts1_id
1365        nc_dims2_id(2) = nc_numlinks_id
1366
1367        ncstat = nf_def_var (nc_file_id, 'S2', 
1368     &                     il_nftype, 2, nc_dims2_id, 
1369     &                     nc_rmpmatrix2_id)
1370        call netcdf_error_handler(ncstat)
1371      endif
1372
1373      !***
1374      !*** end definition stage
1375      !***
1376
1377      ncstat = nf_enddef(nc_file_id)
1378      call netcdf_error_handler(ncstat)
1379
1380!-----------------------------------------------------------------------
1381!
1382!     compute integer masks
1383!
1384!-----------------------------------------------------------------------
1385
1386      if (direction == 1) then
1387        allocate (src_mask_int(grid1_size),
1388     &            dst_mask_int(grid2_size))
1389
1390        where (grid2_mask)
1391          dst_mask_int = 1
1392        elsewhere
1393          dst_mask_int = 0
1394        endwhere
1395
1396        where (grid1_mask)
1397          src_mask_int = 1
1398        elsewhere
1399          src_mask_int = 0
1400        endwhere
1401      else
1402        allocate (src_mask_int(grid2_size),
1403     &            dst_mask_int(grid1_size))
1404
1405        where (grid1_mask)
1406          dst_mask_int = 1
1407        elsewhere
1408          dst_mask_int = 0
1409        endwhere
1410
1411        where (grid2_mask)
1412          src_mask_int = 1
1413        elsewhere
1414          src_mask_int = 0
1415        endwhere
1416      endif
1417
1418!-----------------------------------------------------------------------
1419!
1420!     change units of lat/lon coordinates if input units different
1421!     from radians. if this is the second mapping, the conversion has
1422!     alread been done.
1423!
1424!-----------------------------------------------------------------------
1425
1426      if (grid1_units(1:7) == 'degrees' .and. direction == 1) then
1427        grid1_center_lat = grid1_center_lat/deg2rad
1428        grid1_center_lon = grid1_center_lon/deg2rad
1429        grid1_corner_lat = grid1_corner_lat/deg2rad
1430        grid1_corner_lon = grid1_corner_lon/deg2rad
1431      endif
1432
1433      if (grid2_units(1:7) == 'degrees' .and. direction == 1) then
1434        grid2_center_lat = grid2_center_lat/deg2rad
1435        grid2_center_lon = grid2_center_lon/deg2rad
1436        grid2_corner_lat = grid2_corner_lat/deg2rad
1437        grid2_corner_lon = grid2_corner_lon/deg2rad
1438      endif
1439
1440!-----------------------------------------------------------------------
1441!
1442!     write mapping data
1443!
1444!-----------------------------------------------------------------------
1445
1446      if (direction == 1) then
1447        itmp1 = nc_srcgrddims_id
1448        itmp2 = nc_dstgrddims_id
1449      else
1450        itmp2 = nc_srcgrddims_id
1451        itmp1 = nc_dstgrddims_id
1452      endif
1453
1454      ncstat = nf_put_var_int(nc_file_id, itmp1, grid1_dims)
1455      call netcdf_error_handler(ncstat)
1456
1457      ncstat = nf_put_var_int(nc_file_id, itmp2, grid2_dims)
1458      call netcdf_error_handler(ncstat)
1459
1460      ncstat = nf_put_var_int(nc_file_id, nc_srcgrdimask_id, 
1461     &                        src_mask_int)
1462      call netcdf_error_handler(ncstat)
1463
1464      ncstat = nf_put_var_int(nc_file_id, nc_dstgrdimask_id,
1465     &                        dst_mask_int)
1466      call netcdf_error_handler(ncstat)
1467
1468      deallocate(src_mask_int, dst_mask_int)
1469
1470      if (direction == 1) then
1471        itmp1 = nc_srcgrdcntrlat_id
1472        itmp2 = nc_srcgrdcntrlon_id
1473        itmp3 = nc_srcgrdcrnrlat_id
1474        itmp4 = nc_srcgrdcrnrlon_id
1475      else
1476        itmp1 = nc_dstgrdcntrlat_id
1477        itmp2 = nc_dstgrdcntrlon_id
1478        itmp3 = nc_dstgrdcrnrlat_id
1479        itmp4 = nc_dstgrdcrnrlon_id
1480      endif
1481
1482      ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_center_lat)
1483      call netcdf_error_handler(ncstat)
1484
1485      ncstat = nf_put_var_double(nc_file_id, itmp2, grid1_center_lon)
1486      call netcdf_error_handler(ncstat)
1487
1488      ncstat = nf_put_var_double(nc_file_id, itmp3, grid1_corner_lat)
1489      call netcdf_error_handler(ncstat)
1490
1491      ncstat = nf_put_var_double(nc_file_id, itmp4, grid1_corner_lon)
1492      call netcdf_error_handler(ncstat)
1493
1494      if (direction == 1) then
1495        itmp1 = nc_dstgrdcntrlat_id
1496        itmp2 = nc_dstgrdcntrlon_id
1497        itmp3 = nc_dstgrdcrnrlat_id
1498        itmp4 = nc_dstgrdcrnrlon_id
1499      else
1500        itmp1 = nc_srcgrdcntrlat_id
1501        itmp2 = nc_srcgrdcntrlon_id
1502        itmp3 = nc_srcgrdcrnrlat_id
1503        itmp4 = nc_srcgrdcrnrlon_id
1504      endif
1505
1506      ncstat = nf_put_var_double(nc_file_id, itmp1, grid2_center_lat)
1507      call netcdf_error_handler(ncstat)
1508
1509      ncstat = nf_put_var_double(nc_file_id, itmp2, grid2_center_lon)
1510      call netcdf_error_handler(ncstat)
1511
1512      ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_corner_lat)
1513      call netcdf_error_handler(ncstat)
1514
1515      ncstat = nf_put_var_double(nc_file_id, itmp4, grid2_corner_lon)
1516      call netcdf_error_handler(ncstat)
1517
1518      if (direction == 1) then
1519        itmp1 = nc_srcgrdarea_id
1520        itmp2 = nc_srcgrdfrac_id
1521        itmp3 = nc_dstgrdarea_id
1522        itmp4 = nc_dstgrdfrac_id
1523      else
1524        itmp1 = nc_dstgrdarea_id
1525        itmp2 = nc_dstgrdfrac_id
1526        itmp3 = nc_srcgrdarea_id
1527        itmp4 = nc_srcgrdfrac_id
1528      endif
1529
1530      if (luse_grid1_area) then
1531        ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_area_in)
1532      else
1533        ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_area)
1534      endif
1535      call netcdf_error_handler(ncstat)
1536
1537      ncstat = nf_put_var_double(nc_file_id, itmp2, grid1_frac)
1538      call netcdf_error_handler(ncstat)
1539
1540      if (luse_grid2_area) then
1541        ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_area)
1542      else
1543        ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_area)
1544      endif
1545      call netcdf_error_handler(ncstat)
1546
1547      ncstat = nf_put_var_double(nc_file_id, itmp4, grid2_frac)
1548      call netcdf_error_handler(ncstat)
1549
1550      if (direction == 1) then
1551        ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id, 
1552     &                          grid1_add_map1)
1553        call netcdf_error_handler(ncstat)
1554
1555        ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id, 
1556     &                          grid2_add_map1)
1557        call netcdf_error_handler(ncstat)
1558
1559        if (num_wts == 1) then
1560          ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id, 
1561     &                               wts_map1)
1562          call netcdf_error_handler(ncstat)
1563        else
1564          allocate(wts1(num_links_map1),wts2(num_wts-1,num_links_map1))
1565
1566          wts1 = wts_map1(1,:)
1567          wts2 = wts_map1(2:,:)
1568
1569          ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id, 
1570     &                               wts1)
1571          call netcdf_error_handler(ncstat)
1572          ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix2_id, 
1573     &                               wts2)
1574          call netcdf_error_handler(ncstat)
1575          deallocate(wts1,wts2)
1576        endif
1577      else
1578        ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id, 
1579     &                          grid2_add_map2)
1580        call netcdf_error_handler(ncstat)
1581
1582        ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id, 
1583     &                          grid1_add_map2)
1584        call netcdf_error_handler(ncstat)
1585
1586        if (num_wts == 1) then
1587          ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id, 
1588     &                               wts_map2)
1589          call netcdf_error_handler(ncstat)
1590        else
1591          allocate(wts1(num_links_map2),wts2(num_wts-1,num_links_map2))
1592
1593          wts1 = wts_map2(1,:)
1594          wts2 = wts_map2(2:,:)
1595
1596          ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id, 
1597     &                               wts1)
1598          call netcdf_error_handler(ncstat)
1599          ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix2_id, 
1600     &                               wts2)
1601          call netcdf_error_handler(ncstat)
1602          deallocate(wts1,wts2)
1603        endif
1604      endif
1605
1606      ncstat = nf_close(nc_file_id)
1607      call netcdf_error_handler(ncstat)
1608
1609!-----------------------------------------------------------------------
1610
1611      end subroutine write_remap_csm
1612
1613!***********************************************************************
1614
1615      end module remap_write
1616
1617!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1618
Note: See TracBrowser for help on using the repository browser.