FMS  2024.03
Flexible Modeling System
get_grid_version.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup get_grid_version_mod get_grid_version_mod
20 !> @ingroup data_override
21 !> @brief get_grid implementations and helper routines for @ref data_override_mod
22 
23 !> @addtogroup get_grid_version_mod
24 !> @{
25 module get_grid_version_mod
26 use constants_mod, only: deg_to_rad
27 use platform_mod, only: r4_kind, r8_kind, fms_path_len
28 use mpp_mod, only : mpp_error,fatal,note, mpp_min, mpp_max
29 use mpp_domains_mod, only : domain2d, operator(.NE.),operator(.EQ.)
30 use mpp_domains_mod, only : mpp_get_global_domain, mpp_get_data_domain
31 use fms2_io_mod, only : fmsnetcdfdomainfile_t, fmsnetcdffile_t, open_file, close_file, &
32  variable_exists, read_data, get_variable_size, get_variable_num_dimensions
33 use mosaic2_mod, only : get_mosaic_tile_grid
34 
35 implicit none
36 
38  module procedure get_grid_version_1_r4
39  module procedure get_grid_version_1_r8
40 end interface get_grid_version_1
41 
43  module procedure get_grid_version_2_r4
44  module procedure get_grid_version_2_r8
45 end interface get_grid_version_2
46 
47 contains
48 
49 !> Get lon and lat of three model (target) grids from grid_spec.nc
50 subroutine check_grid_sizes(domain_name, Domain, nlon, nlat)
51 character(len=12), intent(in) :: domain_name
52 type (domain2d), intent(in) :: Domain
53 integer, intent(in) :: nlon, nlat
54 
55 character(len=184) :: error_message
56 integer :: xsize, ysize
57 
58 call mpp_get_global_domain(domain, xsize=xsize, ysize=ysize)
59 if(nlon .NE. xsize .OR. nlat .NE. ysize) then
60  error_message = 'Error in data_override_init. Size of grid as specified by '// &
61  ' does not conform to that specified by grid_spec.nc.'// &
62  ' From : by From grid_spec.nc: by '
63  error_message( 59: 70) = domain_name
64  error_message(130:141) = domain_name
65  write(error_message(143:146),'(i4)') xsize
66  write(error_message(150:153),'(i4)') ysize
67  write(error_message(174:177),'(i4)') nlon
68  write(error_message(181:184),'(i4)') nlat
69  call mpp_error(fatal,error_message)
70 endif
71 end subroutine check_grid_sizes
72 
73 #include "get_grid_version_r4.fh"
74 #include "get_grid_version_r8.fh"
75 
76 end module get_grid_version_mod
77 !> @}
78 ! close documentation grouping
Close a netcdf or domain file opened with open_file or open_virtual_file.
Definition: fms2_io.F90:166
Opens a given netcdf or domain file.
Definition: fms2_io.F90:122
Read data from a defined field in a file.
Definition: fms2_io.F90:292
subroutine check_grid_sizes(domain_name, Domain, nlon, nlat)
Get lon and lat of three model (target) grids from grid_spec.nc.
subroutine, public get_mosaic_tile_grid(grid_file, fileobj, domain, tile_count)
Gets the name of a mosaic tile grid file.
Definition: mosaic2.F90:429
These routines retrieve the axis specifications associated with the data domains. The domain is a der...
These routines retrieve the axis specifications associated with the global domains....
The domain2D type contains all the necessary information to define the global, compute and data domai...
Error handler.
Definition: mpp.F90:382
Reduction operations. Find the max of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:538
Reduction operations. Find the min of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:560