FV3 Bundle
gradient.F90
Go to the documentation of this file.
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 !***********************************************************************
20 ! <CONTACT EMAIL="Zhi.Liang@noaa.gov">
21 ! Zhi Liang
22 ! </CONTACT>
23 
24 ! <HISTORY SRC="http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/"/>
25 
26 ! <OVERVIEW>
27 ! <TT>gradient_mod</TT> implements some utility routines to calculate gradient.
28 ! </OVERVIEW>
29 
30 ! <DESCRIPTION>
31 ! <TT>gradient_mod</TT> implements some utility routines to calculate gradient.
32 ! Currently only gradient on cubic grid is implemented. Also a public interface
33 ! is provided to calculate grid information needed to calculate gradient.
34 
35 use mpp_mod, only : mpp_error, fatal
36 use constants_mod, only : radius
37 
38 implicit none
39 private
40 
41 
42 public :: gradient_cubic
43 public :: calc_cubic_grid_info
44 
45 ! Include variable "version" to be written to log file.
46 #include<file_version.h>
47 
48 contains
49 
50 
51 !#####################################################################
52 ! NOTe: pin has halo size = 1.
53 ! the size of pin will be (nx+2,ny+2), T-cell center, with halo = 1
54 ! the size of dx will be (nx, ny+1), N-cell center
55 ! the size of dy will be (nx+1, ny), E-cell center
56 ! the size of area will be (nx, ny), T-cell center.
57 ! The size of edge_w will be (ny+1), C-cell center
58 ! The size of edge_e will be (ny+1), C-cell center
59 ! The size of edge_s will be (nx+1), C-cell center
60 ! The size of edge_n will be (nx+1), C-cell center
61 ! The size of en_n will be (3,nx,ny+1), N-cell center
62 ! The size of en_e will be (3,nx+1,ny), E-cell center
63 ! The size of vlon will be (3,nx, ny) T-cell center
64 ! The size of vlat will be (3,nx, ny), T-cell center
65 
66 subroutine gradient_cubic(pin, dx, dy, area, edge_w, edge_e, edge_s, edge_n, &
67  en_n, en_e, vlon, vlat, grad_x, grad_y, on_west_edge, &
68  on_east_edge, on_south_edge, on_north_edge)
69 
70  real, dimension(:,: ), intent(in ) :: pin, dx, dy, area
71  real, dimension(: ), intent(in ) :: edge_w, edge_e, edge_s, edge_n
72  real, dimension(:,:,:), intent(in ) :: en_n, en_e
73  real, dimension(:,:,:), intent(in ) :: vlon, vlat
74  real, dimension(:,: ), intent(out) :: grad_x, grad_y
75  logical, intent(in ) :: on_west_edge, on_east_edge, on_south_edge, on_north_edge
76  integer :: nx, ny
77 
78 
79  nx = size(grad_x,1)
80  ny = size(grad_x,2)
81 
82  if(size(pin,1) .NE. nx+2 .OR. size(pin,2) .NE. ny+2)call mpp_error(fatal, "gradient_mod:size of pin should be (nx+2, ny+2)")
83  if(size(dx,1) .NE. nx .OR. size(dx,2) .NE. ny+1 ) call mpp_error(fatal, "gradient_mod: size of dx should be (nx,ny+1)")
84  if(size(dy,1) .NE. nx+1 .OR. size(dy,2) .NE. ny ) call mpp_error(fatal, "gradient_mod: size of dy should be (nx+1,ny)")
85  if(size(area,1) .NE. nx .OR. size(area,2) .NE. ny ) call mpp_error(fatal, "gradient_mod: size of area should be (nx,ny)")
86  if(size(vlon,1) .NE. 3 .OR. size(vlon,2) .NE. nx .OR. size(vlon,3) .NE. ny) &
87  call mpp_error(fatal, "gradient_mod: size of vlon should be (3,nx,ny)")
88  if(size(vlat,1) .NE. 3 .OR. size(vlat,2) .NE. nx .OR. size(vlat,3) .NE. ny) &
89  call mpp_error(fatal, "gradient_mod: size of vlat should be (3,nx,ny)")
90  if(size(edge_w) .NE. ny+1) call mpp_error(fatal, "gradient_mod: size of edge_w should be (ny+1)")
91  if(size(edge_e) .NE. ny+1) call mpp_error(fatal, "gradient_mod: size of edge_e should be (ny+1)")
92  if(size(edge_s) .NE. nx+1) call mpp_error(fatal, "gradient_mod: size of edge_s should be (nx+1)")
93  if(size(edge_n) .NE. nx+1) call mpp_error(fatal, "gradient_mod: size of edge_n should be (nx+1)")
94  if(size(en_n,1) .NE. 3 .OR. size(en_n,2) .NE. nx .OR. size(en_n,3) .NE. ny+1 ) &
95  call mpp_error(fatal, "gradient_mod:size of en_n should be (3, nx, ny+1)")
96  if(size(en_e,1) .NE. 3 .OR. size(en_e,2) .NE. nx+1 .OR. size(en_e,3) .NE. ny ) &
97  call mpp_error(fatal, "gradient_mod:size of en_e should be (3, nx+1, ny)")
98 
99  call grad_c2l(nx, ny, pin, dx, dy, area, edge_w, edge_e, edge_s, edge_n, en_n, en_e, vlon, vlat, &
100  grad_x, grad_y, on_west_edge, on_east_edge, on_south_edge, on_north_edge)
101 
102  return
103 
104 end subroutine gradient_cubic
105 
106 
107 subroutine calc_cubic_grid_info(xt, yt, xc, yc, dx, dy, area, edge_w, edge_e, edge_s, edge_n, &
108  en_n, en_e, vlon, vlat, on_west_edge, on_east_edge, on_south_edge, on_north_edge )
109  real, dimension(:,: ), intent(in ) :: xt, yt, xc, yc
110  real, dimension(:,: ), intent(out) :: dx, dy, area
111  real, dimension(: ), intent(out) :: edge_w, edge_e, edge_s, edge_n
112  real, dimension(:,:,:), intent(out) :: en_n, en_e
113  real, dimension(:,:,:), intent(out) :: vlon, vlat
114  logical, intent(in ) :: on_west_edge, on_east_edge, on_south_edge, on_north_edge
115  integer :: nx, ny, nxp, nyp
116 
117 
118  nx = size(area,1)
119  ny = size(area,2)
120  nxp = nx+1
121  nyp = ny+1
122 
123  if(size(xt,1) .NE. nx+2 .OR. size(xt,2) .NE. ny+2 ) call mpp_error(fatal, "gradient_mod: size of xt should be (nx+2,ny+2)")
124  if(size(yt,1) .NE. nx+2 .OR. size(yt,2) .NE. ny+2 ) call mpp_error(fatal, "gradient_mod: size of yt should be (nx+2,ny+2)")
125  if(size(xc,1) .NE. nxp .OR. size(xc,2) .NE. nyp ) call mpp_error(fatal, "gradient_mod: size of xc should be (nx+1,ny+1)")
126  if(size(yc,1) .NE. nxp .OR. size(yc,2) .NE. nyp ) call mpp_error(fatal, "gradient_mod: size of yc should be (nx+1,ny+1)")
127  if(size(dx,1) .NE. nx .OR. size(dx,2) .NE. nyp ) call mpp_error(fatal, "gradient_mod: size of dx should be (nx,ny+1)")
128  if(size(dy,1) .NE. nxp .OR. size(dy,2) .NE. ny ) call mpp_error(fatal, "gradient_mod: size of dy should be (nx+1,ny)")
129  if(size(area,1) .NE. nx .OR. size(area,2) .NE. ny ) call mpp_error(fatal, "gradient_mod: size of area should be (nx,ny)")
130  if(size(vlon,1) .NE. 3 .OR. size(vlon,2) .NE. nx .OR. size(vlon,3) .NE. ny) &
131  call mpp_error(fatal, "gradient_mod: size of vlon should be (3,nx,ny)")
132  if(size(vlat,1) .NE. 3 .OR. size(vlat,2) .NE. nx .OR. size(vlat,3) .NE. ny) &
133  call mpp_error(fatal, "gradient_mod: size of vlat should be (3,nx,ny)")
134  if(size(edge_w) .NE. ny+1) call mpp_error(fatal, "gradient_mod: size of edge_w should be (ny-1)")
135  if(size(edge_e) .NE. ny+1) call mpp_error(fatal, "gradient_mod: size of edge_e should be (ny-1)")
136  if(size(edge_s) .NE. nx+1) call mpp_error(fatal, "gradient_mod: size of edge_s should be (nx-1)")
137  if(size(edge_n) .NE. nx+1) call mpp_error(fatal, "gradient_mod: size of edge_n should be (nx-1)")
138  if(size(en_n,1) .NE. 3 .OR. size(en_n,2) .NE. nx .OR. size(en_n,3) .NE. nyp ) &
139  call mpp_error(fatal, "gradient_mod:size of en_n should be (3, nx, ny+1)")
140  if(size(en_e,1) .NE. 3 .OR. size(en_e,2) .NE. nxp .OR. size(en_e,3) .NE. ny ) &
141  call mpp_error(fatal, "gradient_mod:size of en_e should be (3, nx+1, ny)")
142 
143 
144  call calc_c2l_grid_info(nx, ny, xt, yt, xc, yc, dx, dy, area, edge_w, edge_e, edge_s, edge_n, &
145  en_n, en_e, vlon, vlat, on_west_edge, on_east_edge, on_south_edge, on_north_edge )
146 
147 
148  return
149 
150 end subroutine calc_cubic_grid_info
151 
152 end module gradient_mod
real, parameter, public radius
Radius of the Earth [m].
Definition: constants.F90:72
Definition: mpp.F90:39
void grad_c2l(const int *nlon, const int *nlat, const double *pin, const double *dx, const double *dy, const double *area, const double *edge_w, const double *edge_e, const double *edge_s, const double *edge_n, const double *en_n, const double *en_e, const double *vlon, const double *vlat, double *grad_x, double *grad_y, const int *on_west_edge, const int *on_east_edge, const int *on_south_edge, const int *on_north_edge)
Definition: gradient_c2l.c:53
subroutine, public calc_cubic_grid_info(xt, yt, xc, yc, dx, dy, area, edge_w, edge_e, edge_s, edge_n, en_n, en_e, vlon, vlat, on_west_edge, on_east_edge, on_south_edge, on_north_edge)
Definition: gradient.F90:109
subroutine, public gradient_cubic(pin, dx, dy, area, edge_w, edge_e, edge_s, edge_n, en_n, en_e, vlon, vlat, grad_x, grad_y, on_west_edge, on_east_edge, on_south_edge, on_north_edge)
Definition: gradient.F90:69
void calc_c2l_grid_info(int *nx_pt, int *ny_pt, const double *xt, const double *yt, const double *xc, const double *yc, double *dx, double *dy, double *area, double *edge_w, double *edge_e, double *edge_s, double *edge_n, double *en_n, double *en_e, double *vlon, double *vlat, int *on_west_edge, int *on_east_edge, int *on_south_edge, int *on_north_edge)
Definition: gradient_c2l.c:379