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

Last change on this file since 6331 was 6331, checked in by aclsce, 17 months ago

Moved oasis-mct_5.0 in oasis3-mct/branches directory.

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