c-lamp_rewrite.f90
changeset 1 2bbb71ef2f67
equal deleted inserted replaced
-1:000000000000 0:ed26b01c1523
       
     1   program clamp_rewrite
       
     2 
       
     3 ! Program to use CMOR routines to rewrite model results from the simulations
       
     4 ! for the Carbon-Land Model Intercomparison Project (C-LAMP).
       
     5 !
       
     6 ! Forrest M. Hoffman <forrest@climatemodeling.org>
       
     7 ! Created: Fri Jun  8 11:50:01 EDT 2007
       
     8 
       
     9 
       
    10     use cmor_users_functions
       
    11     use netcdf
       
    12     use kind_mod
       
    13     use clm_mod, only: namelen, varnamelen, coord, var_data, &
       
    14       clm_read_coord, clm_copy_grid_data, clm_read_data, clm_convert_data, &
       
    15       clm_free_data
       
    16     implicit none
       
    17 
       
    18     integer, parameter :: maxf = 2500
       
    19     integer, parameter :: maxv = 200
       
    20     integer :: i, numf, numv
       
    21 
       
    22     integer, allocatable :: iaxis(:)          ! CMOR handle for axes (time, [hour,] [z-level,] latitude, longitude)
       
    23     integer              :: ivar              ! CMOR handle for variable
       
    24     integer              :: ierr              ! error flag
       
    25     integer              :: axis_num          ! axis counter
       
    26     integer              :: realization       ! model run realization
       
    27     character (len=128)  :: input_path        ! input path
       
    28     character (len=128)  :: fnames(maxf) = '' ! input file names
       
    29     integer              :: tshift(maxf) = 0  ! time shift for files
       
    30     character (len=128)  :: output_path       ! output path
       
    31     character (len=128)  :: input_table       ! CMOR table path
       
    32     character (len=32)   :: vnames(maxv) = '' ! input variable names
       
    33     character (len=256)  :: source = ''       ! source global attribute
       
    34     character (len=128)  :: experiment = ''   ! experiment_id global attribute
       
    35     logical              :: casa_flux_bug     ! scale fluxes due to bug in CASA'
       
    36 
       
    37     real(r4) :: imissing = 1.0e+36 ! missing-data flag on input
       
    38     real(r4) :: omissing = 1.0e+36 ! missing-data flag on output
       
    39 
       
    40     casa_flux_bug = .false.
       
    41 
       
    42     namelist /inparm/ input_table, output_path, experiment, source, &
       
    43       realization, input_path, fnames, tshift, vnames, casa_flux_bug
       
    44 
       
    45     open(10, file='namelist', status='old')
       
    46     ierr = 1
       
    47     do while(ierr /= 0)
       
    48       read(10, inparm, iostat=ierr)
       
    49       if (ierr < 0) then
       
    50         stop 'End of file on namelist read'
       
    51       end if
       
    52     end do
       
    53     close(10)
       
    54 
       
    55     ! Count number of input data files
       
    56     numf = 0
       
    57     do i = 1, maxf
       
    58        if (trim(fnames(i)) /= '') numf = numf + 1
       
    59     end do
       
    60     ! Count number of variables
       
    61     numv = 0
       
    62     do i = 1, maxv
       
    63        if (trim(vnames(i)) /= '') numv = numv + 1
       
    64     end do
       
    65 
       
    66     print *, 'Reading input files for coordinates and counts...'
       
    67     call clm_read_coord(numf, input_path, fnames, .true.)
       
    68 
       
    69     do i = 1, numv
       
    70       select case (vnames(i))
       
    71       ! gca
       
    72       case ('area')
       
    73          call clm_copy_grid_data(vnames(i))
       
    74       ! lbm
       
    75       case ('landmask')
       
    76          call clm_copy_grid_data(vnames(i))
       
    77       ! orog
       
    78       case ('topo')
       
    79          call clm_copy_grid_data(vnames(i))
       
    80       ! sftlf
       
    81       case ('landfrac')
       
    82          call clm_copy_grid_data(vnames(i))
       
    83       case default
       
    84          print *, 'Reading input files for variable ', vnames(i)
       
    85          call clm_read_data(numf, input_path, fnames, tshift, vnames(i))
       
    86 
       
    87          !print *,'var_data%varname = ',var_data%varname
       
    88          !print *,'var_data%ndims = ',var_data%ndims
       
    89          !print *,'var_data%time = ',var_data%time
       
    90          !print *,'var_data%long_name = ',trim(var_data%long_name)
       
    91          !print *,'var_data%missing_value = ',var_data%missing_value
       
    92          !print *,'var_data%units = ',trim(var_data%units)
       
    93       end select
       
    94 
       
    95       print *, 'Performing unit conversion on data...'
       
    96       call clm_convert_data(casa_flux_bug)
       
    97 
       
    98       print *, 'Initializing CMOR...'
       
    99       ierr = cmor_setup(inpath='./',netcdf_file_action='preserve', &
       
   100         set_verbosity=2, exit_control=2)
       
   101 
       
   102       print *, 'Identifying output data sets for CMOR...'
       
   103       ierr = cmor_dataset(                                              &
       
   104         outpath       = output_path,                                    &
       
   105         experiment_id = experiment,                                     &
       
   106         institution   = 'ORNL (Oak Ridge National Laboratory, Oak Ridge, Tennessee, USA)', &
       
   107         source        = trim(source),                                   &
       
   108         realization   = realization,                                    &
       
   109         calendar      = 'noleap',                                       &
       
   110         history       = 'Extracted from case '//trim(coord%case_id),    &
       
   111         comment       = 'Initial dataset: '//trim(coord%inidat)//       &
       
   112           '; surface dataset: '//trim(coord%surdat)//'; pft dataset: '//&
       
   113           trim(coord%pftdat)//'; rtm dataset: '//trim(coord%rtmdat),    &
       
   114         references    = 'http://www.climatemodeling.org/c-lamp',        &
       
   115         contact       = 'Forrest M. Hoffman <forrest@climatemodeling.org>')
       
   116 
       
   117       allocate(iaxis(var_data%ndims), stat=ierr)
       
   118       if (ierr /= 0) then
       
   119          print *, 'Cannot allocate iaxis'
       
   120          stop
       
   121       end if
       
   122 
       
   123       print *, 'Defining coordinates for CMOR output data...'
       
   124       axis_num = 1
       
   125 
       
   126       iaxis(axis_num) = cmor_axis(                                      &
       
   127         table       = input_table,                                      &
       
   128         table_entry = 'longitude',                                      &
       
   129         units       = 'degrees_east',                                   &
       
   130         length      = coord%xsize,                                      &
       
   131         coord_vals  = coord%x)
       
   132       axis_num = axis_num + 1
       
   133 
       
   134       iaxis(axis_num) = cmor_axis(                                      &
       
   135         table       = input_table,                                      &
       
   136         table_entry = 'latitude',                                       &
       
   137         units       = 'degrees_north',                                  &
       
   138         length      = coord%ysize,                                      &
       
   139         coord_vals  = coord%y)
       
   140       axis_num = axis_num + 1
       
   141 
       
   142       if (var_data%soil_layer_flag) then
       
   143          iaxis(axis_num) = cmor_axis(                                   &
       
   144            table       = input_table,                                   &
       
   145            table_entry = 'depth_soil',                                  &
       
   146            units       = 'm',                                           &
       
   147            length      = coord%zsoi_size,                               &
       
   148            coord_vals  = coord%zsoi)
       
   149          axis_num = axis_num + 1
       
   150       end if
       
   151 
       
   152       if (var_data%hour_flag) then
       
   153          iaxis(axis_num) = cmor_axis(                                   &
       
   154            table       = input_table,                                   &
       
   155            table_entry = 'hour',                                        &
       
   156            units       = 'hours',                                       &
       
   157            length      = coord%hr_size,                                 &
       
   158            coord_vals  = coord%hr)
       
   159          axis_num = axis_num + 1
       
   160       end if
       
   161 
       
   162       if (var_data%ndims > 2) then
       
   163         iaxis(axis_num) = cmor_axis(                                    &
       
   164           table       = input_table,                                    &
       
   165           table_entry = 'time',                                         &
       
   166           units       = 'days since 1798-01-01 00:00:00')
       
   167         axis_num = axis_num + 1
       
   168       end if
       
   169 
       
   170       print *, 'Defining CMOR output data variables...'
       
   171       if (var_data%positive == '') then
       
   172          ivar = cmor_variable(                                          &
       
   173            table         = input_table,                                 &
       
   174            table_entry   = var_data%out_varname,                        &
       
   175            original_name = var_data%varname,                            &
       
   176            units         = var_data%units,                              &
       
   177            missing_value = var_data%missing_value,                      &
       
   178            axis_ids      = iaxis)
       
   179       else
       
   180          ivar = cmor_variable(                                          &
       
   181            table         = input_table,                                 &
       
   182            table_entry   = var_data%out_varname,                        &
       
   183            original_name = var_data%varname,                            &
       
   184            units         = var_data%units,                              &
       
   185            positive      = var_data%positive,                           &
       
   186            missing_value = var_data%missing_value,                      &
       
   187            axis_ids      = iaxis)
       
   188       end if
       
   189 
       
   190       print *, 'Writing CMOR output...'
       
   191       select case (var_data%ndims)
       
   192          case (2)
       
   193             if (var_data%int_type) then
       
   194                ierr = cmor_write(                                       &
       
   195                  var_id    = ivar,                                      &
       
   196                  data      = var_data%int2d)
       
   197             else
       
   198                ierr = cmor_write(                                       &
       
   199                  var_id    = ivar,                                      &
       
   200                  data      = var_data%var2d)
       
   201             end if
       
   202          case (3)
       
   203             ierr = cmor_write(                                          &
       
   204               var_id    = ivar,                                         &
       
   205               data      = var_data%var3d,                               &
       
   206               time_vals = var_data%time,                                &
       
   207               time_bnds = var_data%time_bounds)
       
   208          case (4)
       
   209             ierr = cmor_write(                                          &
       
   210               var_id    = ivar,                                         &
       
   211               data      = var_data%var4d,                               &
       
   212               time_vals = var_data%time,                                &
       
   213               time_bnds = var_data%time_bounds)
       
   214          case default
       
   215             print *, 'Unable to handle data with ', var_data%ndims, 'dimensions'
       
   216             stop
       
   217       end select
       
   218 
       
   219       print *, 'Closing CMOR file(s)...'
       
   220       ierr = cmor_close()
       
   221 
       
   222       deallocate(iaxis)
       
   223 
       
   224       print *, 'Freeing data...'
       
   225       call clm_free_data()
       
   226 
       
   227     end do ! loop over variables
       
   228 
       
   229   end program clamp_rewrite