FV3 Bundle
gradient_c2l.c
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  **********************************************************************/
19 #include <math.h>
20 #include <stdlib.h>
21 #include "constant.h"
22 #include "mosaic_util.h"
23 #include "gradient_c2l.h"
24 #include <stdio.h>
25 
26 /*------------------------------------------------------------------------------
27  Routine to compute gradient terms for SCRIP:
28  SJL: Oct 5, 2007
29  NOTe: pin has halo size = 1.
30  the size of pin will be (nx+2,ny+2), T-cell center, with halo = 1
31  the size of dx will be (nx, ny+1), N-cell center
32  the size of dy will be (nx+1, ny), E-cell center
33  the size of area will be (nx, ny), T-cell center.
34  The size of edge_w will be (ny+1), C-cell center
35  The size of edge_e will be (ny+1), C-cell center
36  The size of edge_s will be (nx+1), C-cell center
37  The size of edge_n will be (nx+1), C-cell center
38  The size of en_n will be (nx, ny+1,3),N-cell center
39  The size of en_e will be (nx+1,ny,3), E-cell center
40  The size of vlon will be (nx, ny, 3) T-cell center
41  The size of vlat will be (nx, ny, 3), T-cell center
42  ----------------------------------------------------------------------------*/
43 void grad_c2l_(const int *nlon, const int *nlat, const double *pin, const double *dx, const double *dy, const double *area,
44  const double *edge_w, const double *edge_e, const double *edge_s, const double *edge_n,
45  const double *en_n, const double *en_e, const double *vlon, const double *vlat,
46  double *grad_x, double *grad_y, const int *on_west_edge, const int *on_east_edge,
47  const int *on_south_edge, const int *on_north_edge)
48 {
49  grad_c2l(nlon, nlat, pin, dx, dy, area, edge_w, edge_e, edge_s, edge_n, en_n, en_e, vlon, vlat, grad_x, grad_y,
50  on_west_edge, on_east_edge, on_south_edge, on_north_edge);
51 }
52 
53 void grad_c2l(const int *nlon, const int *nlat, const double *pin, const double *dx, const double *dy, const double *area,
54  const double *edge_w, const double *edge_e, const double *edge_s, const double *edge_n,
55  const double *en_n, const double *en_e, const double *vlon, const double *vlat,
56  double *grad_x, double *grad_y, const int *on_west_edge, const int *on_east_edge,
57  const int *on_south_edge, const int *on_north_edge)
58 {
59 
60  double *pb, *pdx, *pdy, *grad3;
61  int nx, ny, nxp, nyp, i, j, m0, m1, n;
62 
63  nx = *nlon;
64  ny = *nlat;
65  nxp = nx+1;
66  nyp = ny+1;
67  pb = (double *)malloc(nxp*nyp*sizeof(double));
68  pdx = (double *)malloc(3*nx*(ny+1)*sizeof(double));
69  pdy = (double *)malloc(3*(nx+1)*ny*sizeof(double));
70  grad3 = (double *)malloc(3*nx*ny*sizeof(double));
71  a2b_ord2(nx, ny, pin, edge_w, edge_e, edge_s, edge_n, pb, *on_west_edge, *on_east_edge,*on_south_edge, *on_north_edge);
72 
73  for(j=0; j<nyp; j++) for(i=0; i<nx; i++) {
74  m0 = j*nx+i;
75  m1 = j*nxp+i;
76  for(n=0; n<3; n++) {
77  pdx[3*m0+n] = 0.5*(pb[m1]+pb[m1+1])*dx[m0]*en_n[3*m0+n];
78  }
79  }
80 
81  for(j=0; j<ny; j++) for(i=0; i<nxp; i++) {
82  m0 = j*nxp+i;
83  for(n=0; n<3; n++) {
84  pdy[3*m0+n] = 0.5*(pb[m0]+pb[m0+nxp])*dy[m0]*en_e[3*m0+n];
85  }
86  }
87 
88  /* Compute 3D grad of the input scalar field by Green's theorem */
89  for(j=0; j<ny; j++) for(i=0; i<nx; i++) {
90  m0 = 3*(j*nx+i);
91  for(n=0; n<3; n++) {
92  grad3[m0+n] = pdx[3*((j+1)*nx+i)+n]-pdx[m0+n]-pdy[3*(j*nxp+i)+n]+pdy[3*(j*nxp+i+1)+n];
93  }
94  }
95 
96  /* Compute inner product: V3 * grad (pe) */
97  for(j=0; j<ny; j++) for(i=0; i<nx; i++) {
98  m0 = j*nx+i;
99  m1 = 3*m0;
100  /* dq / d(Lamda)*/
101  grad_x[m0] = (vlon[m1]*grad3[m1] + vlon[m1+1]*grad3[m1+1] + vlon[m1+2]*grad3[m1+2])/area[m0];
102  grad_x[m0] *= RADIUS;
103  /* dq / d(theta) */
104  grad_y[m0] = (vlat[m1]*grad3[m1] + vlat[m1+1]*grad3[m1+1] + vlat[m1+2]*grad3[m1+2] )/area[m0];
105  grad_y[m0] *= RADIUS;
106  }
107 
108  free(pb);
109  free(pdx);
110  free(pdy);
111  free(grad3);
112 
113 } /* grad_c2l */
114 
115 /*------------------------------------------------------------------------------
116  qin: A-grid field, size (nx+2, ny+2)
117  qout: B-grid field, size (nx+1, ny+1)
118  ----------------------------------------------------------------------------*/
119 void a2b_ord2(int nx, int ny, const double *qin, const double *edge_w, const double *edge_e,
120  const double *edge_s, const double *edge_n, double *qout,
121  int on_west_edge, int on_east_edge, int on_south_edge, int on_north_edge)
122 {
123  int nxp, nyp, i, j;
124  int istart, iend, jstart, jend;
125  double *q1, *q2;
126  const double r3 = 1./3.;
127 
128  nxp = nx+1;
129  nyp = ny+1;
130  q1 = (double *)malloc((nx+2)*sizeof(double));
131  q2 = (double *)malloc((ny+2)*sizeof(double));
132 
133 
134  if(on_west_edge)
135  istart = 1;
136  else
137  istart = 0;
138  if(on_east_edge)
139  iend = nx;
140  else
141  iend = nxp;
142  if(on_south_edge)
143  jstart = 1;
144  else
145  jstart = 0;
146  if(on_north_edge)
147  jend = ny;
148  else
149  jend = nyp;
150 
151  /* internal region ( 1: nx-1, 1:ny-1) */
152  for(j=jstart; j<jend; j++) for(i=istart; i<iend; i++) {
153  qout[j*nxp+i] = 0.25*(qin[j*(nx+2)+i] + qin[j*(nx+2)+i+1] +
154  qin[(j+1)*(nx+2)+i] + qin[(j+1)*(nx+2)+i+1] );
155  }
156 
157  /* Fix the 4 Corners */
158  if(on_west_edge && on_south_edge)qout[ 0] = r3*(qin[1* (nx+2)+1 ]+qin[1* (nx+2) ]+qin[ 1]); /* sw_corner */
159  if(on_east_edge && on_south_edge)qout[ nx] = r3*(qin[1* (nx+2)+nx]+qin[ nx ]+qin[1* (nx+2)+nxp]); /* se_corner */
160  if(on_east_edge && on_north_edge)qout[ny*nxp+nx] = r3*(qin[ny*(nx+2)+nx]+qin[ny*(nx+2)+nxp]+qin[nyp*(nx+2)+nx ]); /* ne_corner */
161  if(on_west_edge && on_north_edge)qout[ny*nxp ] = r3*(qin[ny*(nx+2)+1 ]+qin[ny*(nx+2) ]+qin[nyp*(nx+2)+1 ]); /* nw_corner */
162 
163  /* West Edges */
164  if(on_west_edge) {
165  for(j=jstart; j<=jend; j++){
166  q2[j] = 0.5*(qin[j*(nx+2)] + qin[j*(nx+2)+1]);
167  }
168  for(j=jstart; j<jend; j++){
169  qout[j*nxp] = edge_w[j]*q2[j] + (1-edge_w[j])*q2[j+1];
170  }
171  }
172 
173  /* East Edges */
174  if(on_east_edge) {
175  for(j=jstart; j<=jend; j++){
176  q2[j] = 0.5*(qin[j*(nx+2)+nx] + qin[j*(nx+2)+nxp]);
177  }
178  for(j=jstart; j<jend; j++){
179  qout[j*nxp+nx] = edge_e[j]*q2[j] + (1-edge_e[j])*q2[j+1];
180  }
181  }
182 
183  /* south edge */
184  if(on_south_edge) {
185  for(i=istart; i<=iend; i++){
186  q1[i] = 0.5*(qin[i] + qin[(nx+2)+i]);
187  }
188  for(i=istart; i<iend; i++){
189  qout[i] = edge_s[i]*q1[i] + (1 - edge_s[i])*q1[i+1];
190  }
191  }
192 
193  /* north edge */
194  if(on_north_edge) {
195  for(i=istart; i<=iend; i++){
196  q1[i] = 0.5*(qin[ny*(nx+2)+i] + qin[nyp*(nx+2)+i]);
197  }
198  for(i=istart; i<iend; i++){
199  qout[ny*nxp+i] = edge_n[i]*q1[i] + (1 - edge_n[i])*q1[i+1];
200  }
201  }
202 
203  free(q1);
204  free(q2);
205 
206 } /* a2b_ord2 */
207 
208 
209 void get_edge(int nx, int ny, const double *lont, const double *latt,
210  const double *lonc, const double *latc, double *edge_w, double *edge_e, double *edge_s, double *edge_n,
211  int on_west_edge, int on_east_edge, int on_south_edge, int on_north_edge)
212 {
213  int i, j, nxp, nyp;
214  int istart, iend, jstart, jend;
215  double p1[2], p2[2];
216  double *py, *px;
217  double d1, d2;
218 
219  nxp = nx + 1;
220  nyp = ny + 1;
221 
222  for(i=0; i<nxp; i++) {
223  edge_s[i] = 0.5; /* dummy value */
224  edge_n[i] = 0.5; /* dummy value */
225  }
226  for(j=0; j<nyp; j++) {
227  edge_w[j] = 0.5; /* dummy value */
228  edge_e[j] = 0.5; /* dummy value */
229  }
230 
231  px = (double *)malloc(2*(nx+2)*sizeof(double));
232  py = (double *)malloc(2*(ny+2)*sizeof(double));
233 
234  if(on_west_edge)
235  istart = 1;
236  else
237  istart = 0;
238  if(on_east_edge)
239  iend = nx;
240  else
241  iend = nxp;
242  if(on_south_edge)
243  jstart = 1;
244  else
245  jstart = 0;
246  if(on_north_edge)
247  jend = ny;
248  else
249  jend = nyp;
250  /* west edge */
251 
252  if(on_west_edge) {
253  i=0;
254  for(j=jstart; j<=jend; j++) {
255  /* get mid point sphere */
256  p1[0] = lont[j*(nx+2)+i ]; p1[1] = latt[j*(nx+2)+i ];
257  p2[0] = lont[j*(nx+2)+i+1]; p2[1] = latt[j*(nx+2)+i+1];
258  mid_pt_sphere(p1, p2, &(py[2*j]));
259  }
260 
261  for(j=jstart; j<jend; j++) {
262  p1[0] = lonc[j*nxp+i];
263  p1[1] = latc[j*nxp+i];
264  d1 = great_circle_distance(py+2*j, p1);
265  d2 = great_circle_distance(py+2*(j+1), p1);
266  edge_w[j] = d2/(d1+d2);
267  }
268  }
269  /* east edge */
270  if(on_east_edge) {
271  i=nx;
272  for(j=jstart; j<=jend; j++) {
273  /* get mid point sphere */
274  p1[0] = lont[j*(nx+2)+i ]; p1[1] = latt[j*(nx+2)+i ];
275  p2[0] = lont[j*(nx+2)+i+1]; p2[1] = latt[j*(nx+2)+i+1];
276  mid_pt_sphere(p1, p2, &(py[2*j]));
277  }
278 
279  for(j=jstart; j<jend; j++) {
280  p1[0] = lonc[j*nxp+i];
281  p1[1] = latc[j*nxp+i];
282  d1 = great_circle_distance(&(py[2*j]), p1);
283  d2 = great_circle_distance(&(py[2*(j+1)]), p1);
284  edge_e[j] = d2/(d1+d2);
285  }
286  }
287 
288  /* south edge */
289  if(on_south_edge) {
290  j=0;
291  for(i=istart; i<=iend; i++) {
292  p1[0] = lont[j *(nx+2)+i]; p1[1] = latt[j *(nx+2)+i];
293  p2[0] = lont[(j+1)*(nx+2)+i]; p2[1] = latt[(j+1)*(nx+2)+i];
294  mid_pt_sphere(p1, p2, &(px[2*i]));
295  }
296  for(i=istart; i<iend; i++) {
297  p1[0] = lonc[j*nxp+i];
298  p1[1] = latc[j*nxp+i];
299  d1 = great_circle_distance(&(px[2*i]), p1);
300  d2 = great_circle_distance(&(px[2*(i+1)]), p1);
301  edge_s[i] = d2/(d1+d2);
302  }
303  }
304  /* north edge */
305  if(on_north_edge) {
306  j=ny;
307  for(i=istart; i<=iend; i++) {
308  p1[0] = lont[j *(nx+2)+i]; p1[1] = latt[j *(nx+2)+i];
309  p2[0] = lont[(j+1)*(nx+2)+i]; p2[1] = latt[(j+1)*(nx+2)+i];
310  mid_pt_sphere(p1, p2, &(px[2*i]));
311  }
312  for(i=istart; i<iend; i++) {
313  p1[0] = lonc[j*nxp+i];
314  p1[1] = latc[j*nxp+i];
315  d1 = great_circle_distance(&(px[2*i]), p1);
316  d2 = great_circle_distance(&(px[2*(i+1)]), p1);
317  edge_n[i] = d2/(d1+d2);
318  }
319  }
320 
321  free(px);
322  free(py);
323 
324 } /* get_edge */
325 
326 void mid_pt_sphere(const double *p1, const double *p2, double *pm)
327 {
328  double e1[3], e2[3], e3[3];
329 
330  latlon2xyz(1, p1, p1+1, e1, e1+1, e1+2);
331  latlon2xyz(1, p2, p2+1, e2, e2+1, e2+2);
332  mid_pt3_cart(e1, e2, e3);
333  xyz2latlon(1, e3, e3+1, e3+2, pm, pm+1);
334 
335 }
336 
337 void mid_pt3_cart(const double *p1, const double *p2, double *e)
338 {
339  double dd;
340 
341  e[0] = p1[0] + p2[0];
342  e[1] = p1[1] + p2[1];
343  e[2] = p1[2] + p2[2];
344  dd = sqrt( e[0]*e[0] + e[1]*e[1] + e[2]*e[2] );
345  e[0] /= dd;
346  e[1] /= dd;
347  e[2] /= dd;
348 }
349 
350 /**********************************************************************************************
351  This routine is used to calculate grid information for second order conservative interpolation
352  from cubic grid to other grid
353  the size of xt will be (nx+2,ny+2), T-cell center, with halo = 1
354  the size of yt will be (nx+2,ny+2), T-cell center, with halo = 1
355  the size of xc will be (nx+1,ny+1), C-cell center
356  the size of yc will be (nx+1,ny+1), C-cell center
357  the size of dx will be (nx, ny+1), N-cell center
358  the size of dy will be (nx+1, ny), E-cell center
359  the size of area will be (nx, ny), T-cell center.
360  The size of edge_w will be (ny-1), C-cell center, without two end point
361  The size of edge_e will be (ny-1), C-cell center, without two end point
362  The size of edge_s will be (nx-1), C-cell center, without two end point
363  The size of edge_n will be (nx-1), C-cell center, without two end point
364  The size of en_n will be (nx, ny+1,3),N-cell center
365  The size of en_e will be (nx+1,ny,3), E-cell center
366  The size of vlon will be (nx, ny) T-cell center
367  The size of vlat will be (nx, ny), T-cell center
368 **********************************************************************************************/
369 void calc_c2l_grid_info_(int *nx_pt, int *ny_pt, const double *xt, const double *yt, const double *xc, const double *yc,
370  double *dx, double *dy, double *area, double *edge_w, double *edge_e, double *edge_s,
371  double *edge_n, double *en_n, double *en_e, double *vlon, double *vlat,
372  int *on_west_edge, int *on_east_edge, int *on_south_edge, int *on_north_edge)
373 {
374  calc_c2l_grid_info(nx_pt, ny_pt, xt, yt, xc, yc, dx, dy, area, edge_w, edge_e, edge_s, edge_n,
375  en_n, en_e, vlon, vlat, on_west_edge, on_east_edge, on_south_edge, on_north_edge);
376 
377 }
378 
379 void calc_c2l_grid_info(int *nx_pt, int *ny_pt, const double *xt, const double *yt, const double *xc, const double *yc,
380  double *dx, double *dy, double *area, double *edge_w, double *edge_e, double *edge_s,
381  double *edge_n, double *en_n, double *en_e, double *vlon, double *vlat,
382  int *on_west_edge, int *on_east_edge, int *on_south_edge, int *on_north_edge)
383 {
384  double *x, *y, *z, *xt_tmp, *yt_tmp;
385  int nx, ny, nxp, nyp, i, j;
386  double p1[3], p2[3], p3[3], p4[3];
387 
388 
389  nx = *nx_pt;
390  ny = *ny_pt;
391  nxp = nx+1;
392  nyp = ny+1;
393 
394  for(j=0; j<nyp; j++) for(i=0; i<nx; i++) {
395  p1[0] = xc[j*nxp+i];
396  p1[1] = yc[j*nxp+i];
397  p2[0] = xc[j*nxp+i+1];
398  p2[1] = yc[j*nxp+i+1];
399  dx[j*nx+i] = great_circle_distance(p1, p2);
400  }
401 
402  for(j=0; j<ny; j++) for(i=0; i<nxp; i++) {
403  p1[0] = xc[j*nxp+i];
404  p1[1] = yc[j*nxp+i];
405  p2[0] = xc[(j+1)*nxp+i];
406  p2[1] = yc[(j+1)*nxp+i];
407  dy[j*nxp+i] = great_circle_distance(p1, p2);
408  }
409 
410  for(j=0; j<ny; j++) for(i=0; i<nx; i++) {
411  p1[0] = xc[j*nxp+i]; /* ll lon */
412  p1[1] = yc[j*nxp+i]; /* ll lat */
413  p2[0] = xc[(j+1)*nxp+i]; /* ul lon */
414  p2[1] = yc[(j+1)*nxp+i]; /* ul lat */
415  p3[0] = xc[j*nxp+i+1]; /* lr lon */
416  p3[1] = yc[j*nxp+i+1]; /* lr lat */
417  p4[0] = xc[(j+1)*nxp+i+1]; /* ur lon */
418  p4[1] = yc[(j+1)*nxp+i+1]; /* ur lat */
419  area[j*nx+i] = spherical_excess_area(p1, p2, p3, p4, RADIUS);
420  }
421 
422  x = (double *)malloc(nxp*nyp*sizeof(double));
423  y = (double *)malloc(nxp*nyp*sizeof(double));
424  z = (double *)malloc(nxp*nyp*sizeof(double));
425 
426  latlon2xyz(nxp*nyp, xc, yc, x, y, z);
427  for(j=0; j<nyp; j++) for(i=0; i<nx; i++) {
428  p1[0] = x[j*nxp+i];
429  p1[1] = y[j*nxp+i];
430  p1[2] = z[j*nxp+i];
431  p2[0] = x[j*nxp+i+1];
432  p2[1] = y[j*nxp+i+1];
433  p2[2] = z[j*nxp+i+1];
434  vect_cross(p1, p2, en_n+3*(j*nx+i) );
435  normalize_vect(en_n+3*(j*nx+i));
436  }
437 
438  for(j=0; j<ny; j++) for(i=0; i<nxp; i++) {
439  p2[0] = x[j*nxp+i];
440  p2[1] = y[j*nxp+i];
441  p2[2] = z[j*nxp+i];
442  p1[0] = x[(j+1)*nxp+i];
443  p1[1] = y[(j+1)*nxp+i];
444  p1[2] = z[(j+1)*nxp+i];
445  vect_cross(p1, p2, en_e+3*(j*nxp+i) );
446  normalize_vect(en_e+3*(j*nxp+i));
447  }
448 
449  xt_tmp = (double *)malloc(nx*ny*sizeof(double));
450  yt_tmp = (double *)malloc(nx*ny*sizeof(double));
451  for(j=0; j<ny; j++)for(i=0; i<nx; i++) {
452  xt_tmp[j*nx+i] = xt[(j+1)*(nx+2)+i+1];
453  yt_tmp[j*nx+i] = yt[(j+1)*(nx+2)+i+1];
454  }
455  unit_vect_latlon(nx*ny, xt_tmp, yt_tmp, vlon, vlat);
456  get_edge(nx, ny, xt, yt, xc, yc, edge_w, edge_e, edge_s, edge_n, *on_west_edge, *on_east_edge,
457  *on_south_edge, *on_north_edge);
458 
459  free(x);
460  free(y);
461  free(z);
462  free(xt_tmp);
463  free(yt_tmp);
464 
465 }
l_size ! loop over number of fields ke do je do i
real, parameter r3
void mid_pt3_cart(const double *p1, const double *p2, double *e)
Definition: gradient_c2l.c:337
integer, parameter nx
void mid_pt_sphere(const double *p1, const double *p2, double *pm)
Definition: gradient_c2l.c:326
void vect_cross(const double *p1, const double *p2, double *e)
Definition: mosaic_util.c:648
integer, parameter, public double
Definition: Type_Kinds.f90:106
void xyz2latlon(int np, const double *x, const double *y, const double *z, double *lon, double *lat)
Definition: mosaic_util.c:227
integer, parameter ny
l_size ! loop over number of fields ke do j
void a2b_ord2(int nx, int ny, const double *qin, const double *edge_w, const double *edge_e, const double *edge_s, const double *edge_n, double *qout, int on_west_edge, int on_east_edge, int on_south_edge, int on_north_edge)
Definition: gradient_c2l.c:119
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:43
integer nlon
No description.
pure elemental type(point) function latlon2xyz(lat, lon)
Definition: diag_grid.F90:1196
double spherical_excess_area(const double *p_ll, const double *p_ul, const double *p_lr, const double *p_ur, double radius)
Definition: mosaic_util.c:606
real(fp), parameter, public e
real, dimension(:,:), allocatable area
double great_circle_distance(double *p1, double *p2)
Definition: mosaic_util.c:488
void get_edge(int nx, int ny, const double *lont, const double *latt, const double *lonc, const double *latc, double *edge_w, double *edge_e, double *edge_s, double *edge_n, int on_west_edge, int on_east_edge, int on_south_edge, int on_north_edge)
Definition: gradient_c2l.c:209
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:369
#define RADIUS
Definition: constant.h:19
real, dimension(9) pm
void normalize_vect(double *e)
Definition: mosaic_util.c:679
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
real, parameter p2
Definition: sw_core_nlm.F90:47
real, parameter p1
Definition: sw_core_nlm.F90:46
integer nlat
No description.
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
void unit_vect_latlon(int size, const double *lon, const double *lat, double *vlon, double *vlat)
Definition: mosaic_util.c:697