FV3 Bundle
mosaic_util.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 <stdlib.h>
20 #include <stdio.h>
21 #include <math.h>
22 #include <string.h>
23 #ifdef use_libMPI
24 #include <mpi.h>
25 #endif
26 #include "mosaic_util.h"
27 #include "constant.h"
28 
29 #define HPI (0.5*M_PI)
30 #define TPI (2.0*M_PI)
31 #define TOLORENCE (1.e-6)
32 #define EPSLN8 (1.e-8)
33 #define EPSLN10 (1.e-10)
34 #define EPSLN15 (1.e-15)
35 #define EPSLN30 (1.e-30)
36 /***********************************************************
37  void error_handler(char *str)
38  error handler: will print out error message and then abort
39 ***********************************************************/
41 
43 {
44  reproduce_siena = 1;
45 }
46 
47 #ifndef __AIX
49 {
50  reproduce_siena = 1;
51 }
52 #endif
53 
54 
55 void error_handler(const char *msg)
56 {
57  fprintf(stderr, "FATAL Error: %s\n", msg );
58 #ifdef use_libMPI
59  MPI_Abort(MPI_COMM_WORLD, -1);
60 #else
61  exit(1);
62 #endif
63 } /* error_handler */
64 
65 /*********************************************************************
66 
67  int nearest_index(double value, const double *array, int ia)
68 
69  return index of nearest data point within "array" corresponding to "value".
70  if "value" is outside the domain of "array" then nearest_index = 0
71  or = size(array)-1 depending on whether array(0) or array(ia-1) is
72  closest to "value"
73 
74  Arguments:
75  value: arbitrary data...same units as elements in "array"
76  array: array of data points (must be monotonically increasing)
77  ia : size of array.
78 
79 ********************************************************************/
80 int nearest_index(double value, const double *array, int ia)
81 {
82  int index, i;
83  int keep_going;
84 
85  for(i=1; i<ia; i++){
86  if (array[i] < array[i-1])
87  error_handler("nearest_index: array must be monotonically increasing");
88  }
89  if (value < array[0] )
90  index = 0;
91  else if ( value > array[ia-1])
92  index = ia-1;
93  else
94  {
95  i=0;
96  keep_going = 1;
97  while (i < ia && keep_going) {
98  i = i+1;
99  if (value <= array[i]) {
100  index = i;
101  if (array[i]-value > value-array[i-1]) index = i-1;
102  keep_going = 0;
103  }
104  }
105  }
106  return index;
107 
108 }
109 
110 /******************************************************************/
111 
112 void tokenize(const char * const string, const char *tokens, unsigned int varlen,
113  unsigned int maxvar, char * pstring, unsigned int * const nstr)
114 {
115  size_t i, j, nvar, len, ntoken;
116  int found, n;
117 
118  nvar = 0; j = 0;
119  len = strlen(string);
120  ntoken = strlen(tokens);
121  /* here we use the fact that C array [][] is contiguous in memory */
122  if(string[0] == 0)error_handler("Error from tokenize: to-be-parsed string is empty");
123 
124  for(i = 0; i < len; i ++){
125  if(string[i] != ' ' && string[i] != '\t'){
126  found = 0;
127  for(n=0; n<ntoken; n++) {
128  if(string[i] == tokens[n] ) {
129  found = 1;
130  break;
131  }
132  }
133  if(found) {
134  if( j != 0) { /* remove :: */
135  *(pstring + (nvar++)*varlen + j) = 0;
136  j = 0;
137  if(nvar >= maxvar) error_handler("Error from tokenize: number of variables exceeds limit");
138  }
139  }
140  else {
141  *(pstring + nvar*varlen + j++) = string[i];
142  if(j >= varlen ) error_handler("error from tokenize: variable name length exceeds limit during tokenization");
143  }
144  }
145  }
146  *(pstring + nvar*varlen + j) = 0;
147 
148  *nstr = ++nvar;
149 
150 }
151 
152 /*******************************************************************************
153  double maxval_double(int size, double *data)
154  get the maximum value of double array
155 *******************************************************************************/
156 double maxval_double(int size, const double *data)
157 {
158  int n;
159  double maxval;
160 
161  maxval = data[0];
162  for(n=1; n<size; n++){
163  if( data[n] > maxval ) maxval = data[n];
164  }
165 
166  return maxval;
167 
168 } /* maxval_double */
169 
170 
171 /*******************************************************************************
172  double minval_double(int size, double *data)
173  get the minimum value of double array
174 *******************************************************************************/
175 double minval_double(int size, const double *data)
176 {
177  int n;
178  double minval;
179 
180  minval = data[0];
181  for(n=1; n<size; n++){
182  if( data[n] < minval ) minval = data[n];
183  }
184 
185  return minval;
186 
187 } /* minval_double */
188 
189 /*******************************************************************************
190  double avgval_double(int size, double *data)
191  get the average value of double array
192 *******************************************************************************/
193 double avgval_double(int size, const double *data)
194 {
195  int n;
196  double avgval;
197 
198  avgval = 0;
199  for(n=0; n<size; n++) avgval += data[n];
200  avgval /= size;
201 
202  return avgval;
203 
204 } /* avgval_double */
205 
206 
207 /*******************************************************************************
208  void latlon2xyz
209  Routine to map (lon, lat) to (x,y,z)
210 ******************************************************************************/
211 void latlon2xyz(int size, const double *lon, const double *lat, double *x, double *y, double *z)
212 {
213  int n;
214 
215  for(n=0; n<size; n++) {
216  x[n] = cos(lat[n])*cos(lon[n]);
217  y[n] = cos(lat[n])*sin(lon[n]);
218  z[n] = sin(lat[n]);
219  }
220 
221 } /* latlon2xyz */
222 
223 /*------------------------------------------------------------
224  void xyz2laton(np, p, xs, ys)
225  Transfer cartesian coordinates to spherical coordinates
226  ----------------------------------------------------------*/
227 void xyz2latlon( int np, const double *x, const double *y, const double *z, double *lon, double *lat)
228 {
229 
230  double xx, yy, zz;
231  double dist;
232  int i;
233 
234  for(i=0; i<np; i++) {
235  xx = x[i];
236  yy = y[i];
237  zz = z[i];
238  dist = sqrt(xx*xx+yy*yy+zz*zz);
239  xx /= dist;
240  yy /= dist;
241  zz /= dist;
242 
243  if ( fabs(xx)+fabs(yy) < EPSLN10 )
244  lon[i] = 0;
245  else
246  lon[i] = atan2(yy, xx);
247  lat[i] = asin(zz);
248 
249  if ( lon[i] < 0.) lon[i] = 2.*M_PI + lon[i];
250  }
251 
252 } /* xyz2latlon */
253 
254 /*------------------------------------------------------------------------------
255  double box_area(double ll_lon, double ll_lat, double ur_lon, double ur_lat)
256  return the area of a lat-lon grid box. grid is in radians.
257  ----------------------------------------------------------------------------*/
258 double box_area(double ll_lon, double ll_lat, double ur_lon, double ur_lat)
259 {
260  double dx = ur_lon-ll_lon;
261 
262  if(dx > M_PI) dx = dx - 2.0*M_PI;
263  if(dx < -M_PI) dx = dx + 2.0*M_PI;
264 
265  return (dx*(sin(ur_lat)-sin(ll_lat))*RADIUS*RADIUS ) ;
266 
267 } /* box_area */
268 
269 
270 /*------------------------------------------------------------------------------
271  double poly_area(const x[], const y[], int n)
272  obtains area of input polygon by line integrating -sin(lat)d(lon)
273  Vertex coordinates must be in degrees.
274  Vertices must be listed counter-clockwise around polygon.
275  grid is in radians.
276  ----------------------------------------------------------------------------*/
277 double poly_area_dimensionless(const double x[], const double y[], int n)
278 {
279  double area = 0.0;
280  int i;
281 
282  for (i=0;i<n;i++) {
283  int ip = (i+1) % n;
284  double dx = (x[ip]-x[i]);
285  double lat1, lat2;
286  double dy, dat;
287 
288  lat1 = y[ip];
289  lat2 = y[i];
290  if(dx > M_PI) dx = dx - 2.0*M_PI;
291  if(dx < -M_PI) dx = dx + 2.0*M_PI;
292  if (dx==0.0) continue;
293 
294  if ( fabs(lat1-lat2) < SMALL_VALUE) /* cheap area calculation along latitude */
295  area -= dx*sin(0.5*(lat1+lat2));
296  else {
297  if(reproduce_siena) {
298  area += dx*(cos(lat1)-cos(lat2))/(lat1-lat2);
299  }
300  else {
301  dy = 0.5*(lat1-lat2);
302  dat = sin(dy)/dy;
303  area -= dx*sin(0.5*(lat1+lat2))*dat;
304  }
305  }
306  }
307  if(area < 0)
308  return (-area/(4*M_PI));
309  else
310  return (area/(4*M_PI));
311 
312 } /* poly_area */
313 
314 double poly_area(const double x[], const double y[], int n)
315 {
316  double area = 0.0;
317  int i;
318 
319  for (i=0;i<n;i++) {
320  int ip = (i+1) % n;
321  double dx = (x[ip]-x[i]);
322  double lat1, lat2;
323  double dy, dat;
324 
325  lat1 = y[ip];
326  lat2 = y[i];
327  if(dx > M_PI) dx = dx - 2.0*M_PI;
328  if(dx < -M_PI) dx = dx + 2.0*M_PI;
329  if (dx==0.0) continue;
330 
331  if ( fabs(lat1-lat2) < SMALL_VALUE) /* cheap area calculation along latitude */
332  area -= dx*sin(0.5*(lat1+lat2));
333  else {
334  if(reproduce_siena) {
335  area += dx*(cos(lat1)-cos(lat2))/(lat1-lat2);
336  }
337  else {
338  dy = 0.5*(lat1-lat2);
339  dat = sin(dy)/dy;
340  area -= dx*sin(0.5*(lat1+lat2))*dat;
341  }
342  }
343  }
344  if(area < 0)
345  return -area*RADIUS*RADIUS;
346  else
347  return area*RADIUS*RADIUS;
348 
349 } /* poly_area */
350 
351 double poly_area_no_adjust(const double x[], const double y[], int n)
352 {
353  double area = 0.0;
354  int i;
355 
356  for (i=0;i<n;i++) {
357  int ip = (i+1) % n;
358  double dx = (x[ip]-x[i]);
359  double lat1, lat2;
360 
361  lat1 = y[ip];
362  lat2 = y[i];
363  if (dx==0.0) continue;
364 
365  if ( fabs(lat1-lat2) < SMALL_VALUE) /* cheap area calculation along latitude */
366  area -= dx*sin(0.5*(lat1+lat2));
367  else
368  area += dx*(cos(lat1)-cos(lat2))/(lat1-lat2);
369  }
370  if(area < 0)
371  return area*RADIUS*RADIUS;
372  else
373  return area*RADIUS*RADIUS;
374 } /* poly_area_no_adjust */
375 
376 int delete_vtx(double x[], double y[], int n, int n_del)
377 {
378  for (;n_del<n-1;n_del++) {
379  x[n_del] = x[n_del+1];
380  y[n_del] = y[n_del+1];
381  }
382 
383  return (n-1);
384 } /* delete_vtx */
385 
386 int insert_vtx(double x[], double y[], int n, int n_ins, double lon_in, double lat_in)
387 {
388  int i;
389 
390  for (i=n-1;i>=n_ins;i--) {
391  x[i+1] = x[i];
392  y[i+1] = y[i];
393  }
394 
395  x[n_ins] = lon_in;
396  y[n_ins] = lat_in;
397  return (n+1);
398 } /* insert_vtx */
399 
400 void v_print(double x[], double y[], int n)
401 {
402  int i;
403 
404  for (i=0;i<n;i++){
405  printf(" %20g %20g\n", x[i], y[i]);
406  }
407 } /* v_print */
408 
409 int fix_lon(double x[], double y[], int n, double tlon)
410 {
411  double x_sum, dx, dx_;
412  int i, nn = n, pole = 0;
413 
414  for (i=0;i<nn;i++) if (fabs(y[i])>=HPI-TOLORENCE) pole = 1;
415  if (0&&pole) {
416  printf("fixing pole cell\n");
417  v_print(x, y, nn);
418  printf("---------");
419  }
420 
421  /* all pole points must be paired */
422  for (i=0;i<nn;i++) if (fabs(y[i])>=HPI-TOLORENCE) {
423  int im=(i+nn-1)%nn, ip=(i+1)%nn;
424 
425  if (y[im]==y[i] && y[ip]==y[i]) {
426  nn = delete_vtx(x, y, nn, i);
427  i--;
428  } else if (y[im]!=y[i] && y[ip]!=y[i]) {
429  nn = insert_vtx(x, y, nn, i, x[i], y[i]);
430  i++;
431  }
432  }
433  /* first of pole pair has longitude of previous vertex */
434  /* second of pole pair has longitude of subsequent vertex */
435  for (i=0;i<nn;i++) if (fabs(y[i])>=HPI-TOLORENCE) {
436  int im=(i+nn-1)%nn, ip=(i+1)%nn;
437 
438  if (y[im]!=y[i]){
439  x[i] = x[im];
440  }
441  if (y[ip]!=y[i]){
442  x[i] = x[ip];
443  }
444  }
445 
446  if (nn){
447  x_sum = x[0];
448  }
449  else{
450  return(0);
451  }
452  for (i=1;i<nn;i++) {
453  dx_ = x[i]-x[i-1];
454 
455  if (dx_ < -M_PI) dx_ = dx_ + TPI;
456  else if (dx_ > M_PI) dx_ = dx_ - TPI;
457  x_sum += (x[i] = x[i-1] + dx_);
458  }
459 
460  dx = (x_sum/nn)-tlon;
461  if (dx < -M_PI){
462  for (i=0;i<nn;i++){
463  x[i] += TPI;
464  }
465  }
466  else if (dx > M_PI){
467  for (i=0;i<nn;i++){
468  x[i] -= TPI;
469  }
470  }
471 
472  if (0&&pole) {
473  printf("area=%g\n", poly_area(x, y,nn));
474  v_print(x, y, nn);
475  printf("---------");
476  }
477 
478  return (nn);
479 } /* fix_lon */
480 
481 
482 /*------------------------------------------------------------------------------
483  double great_circle_distance()
484  computes distance between two points along a great circle
485  (the shortest distance between 2 points on a sphere)
486  returned in units of meter
487  ----------------------------------------------------------------------------*/
488 double great_circle_distance(double *p1, double *p2)
489 {
490  double dist, beta;
491 
492  /* This algorithm is not accurate for small distance
493  dist = RADIUS*ACOS(SIN(p1[1])*SIN(p2[1]) + COS(p1[1])*COS(p2[1])*COS(p1[0]-p2[0]));
494  */
495  beta = 2.*asin( sqrt( sin((p1[1]-p2[1])/2.)*sin((p1[1]-p2[1])/2.) +
496  cos(p1[1])*cos(p2[1])*(sin((p1[0]-p2[0])/2.)*sin((p1[0]-p2[0])/2.)) ) );
497  dist = RADIUS*beta;
498  return dist;
499 
500 } /* great_circle_distance */
501 
502 
503 /* Compute the great circle area of a polygon on a sphere */
504 double great_circle_area(int n, const double *x, const double *y, const double *z) {
505  int i;
506  double pnt0[3], pnt1[3], pnt2[3];
507  double sum, area;
508 
509  /* sum angles around polygon */
510  sum=0.0;
511  for ( i=0; i<n; i++) {
512  /* points that make up a side of polygon */
513  pnt0[0] = x[i];
514  pnt0[1] = y[i];
515  pnt0[2] = z[i];
516  pnt1[0] = x[(i+1)%n];
517  pnt1[1] = y[(i+1)%n];
518  pnt1[2] = z[(i+1)%n];
519  pnt2[0] = x[(i+2)%n];
520  pnt2[1] = y[(i+2)%n];
521  pnt2[2] = z[(i+2)%n];
522 
523  /* compute angle for pnt1 */
524  sum += spherical_angle(pnt1, pnt2, pnt0);
525 
526  }
527  area = (sum - (n-2.)*M_PI) * RADIUS* RADIUS;
528  return area;
529 }
530 
531 /*------------------------------------------------------------------------------
532  double spherical_angle(const double *p1, const double *p2, const double *p3)
533  p3
534  /
535  /
536  p1 ---> angle
537  \
538  \
539  p2
540  -----------------------------------------------------------------------------*/
541 double spherical_angle(const double *v1, const double *v2, const double *v3)
542 {
543  double angle;
544 #ifdef NO_QUAD_PRECISION
545  double px, py, pz, qx, qy, qz, ddd;
546 #ifndef SQRT_
547 #define SQRT_ sqrt
548 #else
549 #error "SQRT_ Previously Defined"
550 #endif /* SQRT_ */
551 #ifndef ABS_
552 #define ABS_ fabsl
553 #else
554 #error "ABS_ Previously Defined"
555 #endif /* ABS_ */
556 #else
557  long double px, py, pz, qx, qy, qz, ddd;
558 #ifndef SQRT_
559 #define SQRT_ sqrtl
560 #else
561 #error "SQRT_ Previously Defined"
562 #endif /* SQRT_ */
563 #ifndef ABS_
564 #define ABS_ fabs
565 #else
566 #error "ABS_ Previously Defined"
567 #endif /* ABS_ */
568 #endif
569 
570  /* vector product between v1 and v2 */
571  px = v1[1]*v2[2] - v1[2]*v2[1];
572  py = v1[2]*v2[0] - v1[0]*v2[2];
573  pz = v1[0]*v2[1] - v1[1]*v2[0];
574  /* vector product between v1 and v3 */
575  qx = v1[1]*v3[2] - v1[2]*v3[1];
576  qy = v1[2]*v3[0] - v1[0]*v3[2];
577  qz = v1[0]*v3[1] - v1[1]*v3[0];
578 
579  ddd = (px*px+py*py+pz*pz)*(qx*qx+qy*qy+qz*qz);
580  if ( ddd <= 0.0 )
581  angle = 0. ;
582  else {
583  ddd = (px*qx+py*qy+pz*qz) / SQRT_(ddd);
584  if( ABS_(ddd-1) < EPSLN30 ) ddd = 1;
585  if( ABS_(ddd+1) < EPSLN30 ) ddd = -1;
586  if ( ddd>1. || ddd<-1. ) {
587  /*FIX (lmh) to correctly handle co-linear points (angle near pi or 0) */
588  if (ddd < 0.)
589  angle = M_PI;
590  else
591  angle = 0.;
592  }
593  else
594  angle = ((double)acosl( ddd ));
595  }
596 
597  return angle;
598 } /* spherical_angle */
599 
600 /*------------------------------------------------------------------------------
601  double spherical_excess_area(p_lL, p_uL, p_lR, p_uR)
602  get the surface area of a cell defined as a quadrilateral
603  on the sphere. Area is computed as the spherical excess
604  [area units are m^2]
605  ----------------------------------------------------------------------------*/
606 double spherical_excess_area(const double* p_ll, const double* p_ul,
607  const double* p_lr, const double* p_ur, double radius)
608 {
609  double area, ang1, ang2, ang3, ang4;
610  double v1[3], v2[3], v3[3];
611 
612  /* S-W: 1 */
613  latlon2xyz(1, p_ll, p_ll+1, v1, v1+1, v1+2);
614  latlon2xyz(1, p_lr, p_lr+1, v2, v2+1, v2+2);
615  latlon2xyz(1, p_ul, p_ul+1, v3, v3+1, v3+2);
616  ang1 = spherical_angle(v1, v2, v3);
617 
618  /* S-E: 2 */
619  latlon2xyz(1, p_lr, p_lr+1, v1, v1+1, v1+2);
620  latlon2xyz(1, p_ur, p_ur+1, v2, v2+1, v2+2);
621  latlon2xyz(1, p_ll, p_ll+1, v3, v3+1, v3+2);
622  ang2 = spherical_angle(v1, v2, v3);
623 
624  /* N-E: 3 */
625  latlon2xyz(1, p_ur, p_ur+1, v1, v1+1, v1+2);
626  latlon2xyz(1, p_ul, p_ul+1, v2, v2+1, v2+2);
627  latlon2xyz(1, p_lr, p_lr+1, v3, v3+1, v3+2);
628  ang3 = spherical_angle(v1, v2, v3);
629 
630  /* N-W: 4 */
631  latlon2xyz(1, p_ul, p_ul+1, v1, v1+1, v1+2);
632  latlon2xyz(1, p_ur, p_ur+1, v2, v2+1, v2+2);
633  latlon2xyz(1, p_ll, p_ll+1, v3, v3+1, v3+2);
634  ang4 = spherical_angle(v1, v2, v3);
635 
636  area = (ang1 + ang2 + ang3 + ang4 - 2.*M_PI) * radius* radius;
637 
638  return area;
639 
640 } /* spherical_excess_area */
641 
642 
643 /*----------------------------------------------------------------------
644  void vect_cross(e, p1, p2)
645  Perform cross products of 3D vectors: e = P1 X P2
646  -------------------------------------------------------------------*/
647 
648 void vect_cross(const double *p1, const double *p2, double *e )
649 {
650 
651  e[0] = p1[1]*p2[2] - p1[2]*p2[1];
652  e[1] = p1[2]*p2[0] - p1[0]*p2[2];
653  e[2] = p1[0]*p2[1] - p1[1]*p2[0];
654 
655 } /* vect_cross */
656 
657 
658 /*----------------------------------------------------------------------
659  double* vect_cross(p1, p2)
660  return cross products of 3D vectors: = P1 X P2
661  -------------------------------------------------------------------*/
662 
663 double dot(const double *p1, const double *p2)
664 {
665 
666  return( p1[0]*p2[0] + p1[1]*p2[1] + p1[2]*p2[2] );
667 
668 }
669 
670 
671 double metric(const double *p) {
672  return (sqrt(p[0]*p[0] + p[1]*p[1]+p[2]*p[2]) );
673 }
674 
675 
676 /* ----------------------------------------------------------------
677  make a unit vector
678  --------------------------------------------------------------*/
679 void normalize_vect(double *e)
680 {
681  double pdot;
682  int k;
683 
684  pdot = e[0]*e[0] + e[1] * e[1] + e[2] * e[2];
685  pdot = sqrt( pdot );
686 
687  for(k=0; k<3; k++) e[k] /= pdot;
688 }
689 
690 
691 /*------------------------------------------------------------------
692  void unit_vect_latlon(int size, lon, lat, vlon, vlat)
693 
694  calculate unit vector for latlon in cartesian coordinates
695 
696  ---------------------------------------------------------------------*/
697 void unit_vect_latlon(int size, const double *lon, const double *lat, double *vlon, double *vlat)
698 {
699  double sin_lon, cos_lon, sin_lat, cos_lat;
700  int n;
701 
702  for(n=0; n<size; n++) {
703  sin_lon = sin(lon[n]);
704  cos_lon = cos(lon[n]);
705  sin_lat = sin(lat[n]);
706  cos_lat = cos(lat[n]);
707 
708  vlon[3*n] = -sin_lon;
709  vlon[3*n+1] = cos_lon;
710  vlon[3*n+2] = 0.;
711 
712  vlat[3*n] = -sin_lat*cos_lon;
713  vlat[3*n+1] = -sin_lat*sin_lon;
714  vlat[3*n+2] = cos_lat;
715  }
716 } /* unit_vect_latlon */
717 
718 
719 /* Intersect a line and a plane
720  Intersects between the plane ( three points ) (entries in counterclockwise order)
721  and the line determined by the endpoints l1 and l2 (t=0.0 at l1 and t=1.0 at l2)
722  returns true if the two intersect and the output variables are valid
723  outputs p containing the coordinates in the tri and t the coordinate in the line
724  of the intersection.
725  NOTE: the intersection doesn't have to be inside the tri or line for this to return true
726 */
727 int intersect_tri_with_line(const double *plane, const double *l1, const double *l2, double *p,
728  double *t) {
729 
730  long double M[3*3], inv_M[3*3];
731  long double V[3];
732  long double X[3];
733  int is_invert=0;
734 
735  const double *pnt0=plane;
736  const double *pnt1=plane+3;
737  const double *pnt2=plane+6;
738 
739  /* To do intersection just solve the set of linear equations for both
740  Setup M
741  */
742  M[0]=l1[0]-l2[0]; M[1]=pnt1[0]-pnt0[0]; M[2]=pnt2[0]-pnt0[0];
743  M[3]=l1[1]-l2[1]; M[4]=pnt1[1]-pnt0[1]; M[5]=pnt2[1]-pnt0[1];
744  M[6]=l1[2]-l2[2]; M[7]=pnt1[2]-pnt0[2]; M[8]=pnt2[2]-pnt0[2];
745 
746 
747  /* Invert M */
748  is_invert = invert_matrix_3x3(M,inv_M);
749  if (!is_invert) return 0;
750 
751  /* Set variable holding vector */
752  V[0]=l1[0]-pnt0[0];
753  V[1]=l1[1]-pnt0[1];
754  V[2]=l1[2]-pnt0[2];
755 
756  /* Calculate solution */
757  mult(inv_M, V, X);
758 
759  /* Get answer out */
760  *t=((double)X[0]);
761  p[0]=((double)X[1]);
762  p[1]=((double)X[2]);
763 
764  return 1;
765 }
766 
767 
768 void mult(long double m[], long double v[], long double out_v[]) {
769 
770  out_v[0]=m[0]*v[0]+m[1]*v[1]+m[2]*v[2];
771  out_v[1]=m[3]*v[0]+m[4]*v[1]+m[5]*v[2];
772  out_v[2]=m[6]*v[0]+m[7]*v[1]+m[8]*v[2];
773 
774 }
775 
776 
777 /* returns 1 if matrix is inverted, 0 otherwise */
778 int invert_matrix_3x3(long double m[], long double m_inv[]) {
779 
780 
781  const long double det = m[0] * (m[4]*m[8] - m[5]*m[7])
782  -m[1] * (m[3]*m[8] - m[5]*m[6])
783  +m[2] * (m[3]*m[7] - m[4]*m[6]);
784 #ifdef test_invert_matrix_3x3
785  printf("det = %Lf\n", det);
786 #endif
787  if (fabsl(det) < EPSLN15 ) return 0;
788 
789  const long double deti = 1.0/det;
790 
791  m_inv[0] = (m[4]*m[8] - m[5]*m[7]) * deti;
792  m_inv[1] = (m[2]*m[7] - m[1]*m[8]) * deti;
793  m_inv[2] = (m[1]*m[5] - m[2]*m[4]) * deti;
794 
795  m_inv[3] = (m[5]*m[6] - m[3]*m[8]) * deti;
796  m_inv[4] = (m[0]*m[8] - m[2]*m[6]) * deti;
797  m_inv[5] = (m[2]*m[3] - m[0]*m[5]) * deti;
798 
799  m_inv[6] = (m[3]*m[7] - m[4]*m[6]) * deti;
800  m_inv[7] = (m[1]*m[6] - m[0]*m[7]) * deti;
801  m_inv[8] = (m[0]*m[4] - m[1]*m[3]) * deti;
802 
803  return 1;
804 }
805 
806 #ifndef MAXNODELIST
807 #define MAXNODELIST 100
808 #endif
809 
812 
813 void rewindList(void)
814 {
815  int n;
816 
817  curListPos = 0;
818  if(!nodeList) nodeList = (struct Node *)malloc(MAXNODELIST*sizeof(struct Node));
819  for(n=0; n<MAXNODELIST; n++) initNode(nodeList+n);
820 
821 }
822 
823 struct Node *getNext()
824 {
825  struct Node *temp=NULL;
826  int n;
827 
828  if(!nodeList) {
829  nodeList = (struct Node *)malloc(MAXNODELIST*sizeof(struct Node));
830  for(n=0; n<MAXNODELIST; n++) initNode(nodeList+n);
831  }
832 
833  temp = nodeList+curListPos;
834  curListPos++;
835  if(curListPos > MAXNODELIST) error_handler("getNext: curListPos >= MAXNODELIST");
836 
837  return (temp);
838 }
839 
840 
841 void initNode(struct Node *node)
842 {
843  node->x = 0;
844  node->y = 0;
845  node->z = 0;
846  node->u = 0;
847  node->intersect = 0;
848  node->inbound = 0;
849  node->isInside = 0;
850  node->Next = NULL;
851  node->initialized=0;
852 
853 }
854 
855 void addEnd(struct Node *list, double x, double y, double z, int intersect, double u, int inbound, int inside)
856 {
857 
858  struct Node *temp=NULL;
859 
860  if(list == NULL) error_handler("addEnd: list is NULL");
861 
862  if(list->initialized) {
863 
864  /* (x,y,z) might already in the list when intersect is true and u=0 or 1 */
865  temp = list;
866  while (temp) {
867  if(samePoint(temp->x, temp->y, temp->z, x, y, z)) return;
868  temp=temp->Next;
869  }
870  temp = list;
871  while(temp->Next)
872  temp=temp->Next;
873 
874  /* Append at the end of the list. */
875  temp->Next = getNext();
876  temp = temp->Next;
877  }
878  else {
879  temp = list;
880  }
881 
882  temp->x = x;
883  temp->y = y;
884  temp->z = z;
885  temp->u = u;
886  temp->intersect = intersect;
887  temp->inbound = inbound;
888  temp->initialized=1;
889  temp->isInside = inside;
890 }
891 
892 /* return 1 if the point (x,y,z) is added in the list, return 0 if it is already in the list */
893 
894 int addIntersect(struct Node *list, double x, double y, double z, int intersect, double u1, double u2, int inbound,
895  int is1, int ie1, int is2, int ie2)
896 {
897 
898  double u1_cur, u2_cur;
899  int i1_cur, i2_cur;
900  struct Node *temp=NULL;
901 
902  if(list == NULL) error_handler("addEnd: list is NULL");
903 
904  /* first check to make sure this point is not in the list */
905  u1_cur = u1;
906  i1_cur = is1;
907  u2_cur = u2;
908  i2_cur = is2;
909  if(u1_cur == 1) {
910  u1_cur = 0;
911  i1_cur = ie1;
912  }
913  if(u2_cur == 1) {
914  u2_cur = 0;
915  i2_cur = ie2;
916  }
917 
918  if(list->initialized) {
919  temp = list;
920  while(temp) {
921  if( temp->u == u1_cur && temp->subj_index == i1_cur) return 0;
922  if( temp->u_clip == u2_cur && temp->clip_index == i2_cur) return 0;
923  if( !temp->Next ) break;
924  temp=temp->Next;
925  }
926 
927  /* Append at the end of the list. */
928  temp->Next = getNext();
929  temp = temp->Next;
930  }
931  else {
932  temp = list;
933  }
934 
935  temp->x = x;
936  temp->y = y;
937  temp->z = z;
938  temp->intersect = intersect;
939  temp->inbound = inbound;
940  temp->initialized=1;
941  temp->isInside = 0;
942  temp->u = u1_cur;
943  temp->subj_index = i1_cur;
944  temp->u_clip = u2_cur;
945  temp->clip_index = i2_cur;
946 
947  return 1;
948 }
949 
950 
951 int length(struct Node *list)
952 {
953  struct Node *cur_ptr=NULL;
954  int count=0;
955 
956  cur_ptr=list;
957 
958  while(cur_ptr)
959  {
960  if(cur_ptr->initialized ==0) break;
961  cur_ptr=cur_ptr->Next;
962  count++;
963  }
964  return(count);
965 }
966 
967 /* two points are the same if there are close enough */
968 int samePoint(double x1, double y1, double z1, double x2, double y2, double z2)
969 {
970  if( fabs(x1-x2) > EPSLN10 || fabs(y1-y2) > EPSLN10 || fabs(z1-z2) > EPSLN10 )
971  return 0;
972  else
973  return 1;
974 }
975 
976 
977 
978 int sameNode(struct Node node1, struct Node node2)
979 {
980  if( node1.x == node2.x && node1.y == node2.y && node1.z==node2.z )
981  return 1;
982  else
983  return 0;
984 }
985 
986 
987 void addNode(struct Node *list, struct Node inNode)
988 {
989 
990  addEnd(list, inNode.x, inNode.y, inNode.z, inNode.intersect, inNode.u, inNode.inbound, inNode.isInside);
991 
992 }
993 
994 struct Node *getNode(struct Node *list, struct Node inNode)
995 {
996  struct Node *thisNode=NULL;
997  struct Node *temp=NULL;
998 
999  temp = list;
1000  while( temp ) {
1001  if( sameNode( *temp, inNode ) ) {
1002  thisNode = temp;
1003  temp = NULL;
1004  break;
1005  }
1006  temp = temp->Next;
1007  }
1008 
1009  return thisNode;
1010 }
1011 
1012 struct Node *getNextNode(struct Node *list)
1013 {
1014  return list->Next;
1015 }
1016 
1017 void copyNode(struct Node *node_out, struct Node node_in)
1018 {
1019 
1020  node_out->x = node_in.x;
1021  node_out->y = node_in.y;
1022  node_out->z = node_in.z;
1023  node_out->u = node_in.u;
1024  node_out->intersect = node_in.intersect;
1025  node_out->inbound = node_in.inbound;
1026  node_out->Next = NULL;
1027  node_out->initialized = node_in.initialized;
1028  node_out->isInside = node_in.isInside;
1029 }
1030 
1031 void printNode(struct Node *list, char *str)
1032 {
1033  struct Node *temp;
1034 
1035  if(list == NULL) error_handler("printNode: list is NULL");
1036  if(str) printf(" %s \n", str);
1037  temp = list;
1038  while(temp) {
1039  if(temp->initialized ==0) break;
1040  printf(" (x, y, z, interset, inbound, isInside) = (%19.15f,%19.15f,%19.15f,%d,%d,%d)\n",
1041  temp->x, temp->y, temp->z, temp->intersect, temp->inbound, temp->isInside);
1042  temp = temp->Next;
1043  }
1044  printf("\n");
1045 }
1046 
1047 int intersectInList(struct Node *list, double x, double y, double z)
1048 {
1049  struct Node *temp;
1050  int found=0;
1051 
1052  temp = list;
1053  found = 0;
1054  while ( temp ) {
1055  if( temp->x == x && temp->y == y && temp->z == z ) {
1056  found = 1;
1057  break;
1058  }
1059  temp=temp->Next;
1060  }
1061  if (!found) error_handler("intersectInList: point (x,y,z) is not found in the list");
1062  if( temp->intersect == 2 )
1063  return 1;
1064  else
1065  return 0;
1066 
1067 }
1068 
1069 
1070 /* The following insert a intersection after non-intersect point (x2,y2,z2), if the point
1071  after (x2,y2,z2) is an intersection, if u is greater than the u value of the intersection,
1072  insert after, otherwise insert before
1073 */
1074 void insertIntersect(struct Node *list, double x, double y, double z, double u1, double u2, int inbound,
1075  double x2, double y2, double z2)
1076 {
1077  struct Node *temp1=NULL, *temp2=NULL;
1078  struct Node *temp;
1079  double u_cur;
1080  int found=0;
1081 
1082  temp1 = list;
1083  found = 0;
1084  while ( temp1 ) {
1085  if( temp1->x == x2 && temp1->y == y2 && temp1->z == z2 ) {
1086  found = 1;
1087  break;
1088  }
1089  temp1=temp1->Next;
1090  }
1091  if (!found) error_handler("inserAfter: point (x,y,z) is not found in the list");
1092 
1093  /* when u = 0 or u = 1, set the grid point to be the intersection point to solve truncation error isuse */
1094  u_cur = u1;
1095  if(u1 == 1) {
1096  u_cur = 0;
1097  temp1 = temp1->Next;
1098  if(!temp1) temp1 = list;
1099  }
1100  if(u_cur==0) {
1101  temp1->intersect = 2;
1102  temp1->isInside = 1;
1103  temp1->u = u_cur;
1104  temp1->x = x;
1105  temp1->y = y;
1106  temp1->z = z;
1107  return;
1108  }
1109 
1110  /* when u2 != 0 and u2 !=1, can decide if one end of the point is outside depending on inbound value */
1111  if(u2 != 0 && u2 != 1) {
1112  if(inbound == 1) { /* goes outside, then temp1->Next is an outside point */
1113  /* find the next non-intersect point */
1114  temp2 = temp1->Next;
1115  if(!temp2) temp2 = list;
1116  while(temp2->intersect) {
1117  temp2=temp2->Next;
1118  if(!temp2) temp2 = list;
1119  }
1120 
1121  temp2->isInside = 0;
1122  }
1123  else if(inbound ==2) { /* goes inside, then temp1 is an outside point */
1124  temp1->isInside = 0;
1125  }
1126  }
1127 
1128  temp2 = temp1->Next;
1129  while ( temp2 ) {
1130  if( temp2->intersect == 1 ) {
1131  if( temp2->u > u_cur ) {
1132  break;
1133  }
1134  }
1135  else
1136  break;
1137  temp1 = temp2;
1138  temp2 = temp2->Next;
1139  }
1140 
1141  /* assign value */
1142  temp = getNext();
1143  temp->x = x;
1144  temp->y = y;
1145  temp->z = z;
1146  temp->u = u_cur;
1147  temp->intersect = 1;
1148  temp->inbound = inbound;
1149  temp->isInside = 1;
1150  temp->initialized = 1;
1151  temp1->Next = temp;
1152  temp->Next = temp2;
1153 
1154 }
1155 
1156 double gridArea(struct Node *grid) {
1157  double x[20], y[20], z[20];
1158  struct Node *temp=NULL;
1159  double area;
1160  int n;
1161 
1162  temp = grid;
1163  n = 0;
1164  while( temp ) {
1165  x[n] = temp->x;
1166  y[n] = temp->y;
1167  z[n] = temp->z;
1168  n++;
1169  temp = temp->Next;
1170  }
1171 
1172  area = great_circle_area(n, x, y, z);
1173 
1174  return area;
1175 
1176 }
1177 
1178 int isIntersect(struct Node node) {
1179 
1180  return node.intersect;
1181 
1182 }
1183 
1184 
1185 int getInbound( struct Node node )
1186 {
1187  return node.inbound;
1188 }
1189 
1190 struct Node *getLast(struct Node *list)
1191 {
1192  struct Node *temp1;
1193 
1194  temp1 = list;
1195  if( temp1 ) {
1196  while( temp1->Next ) {
1197  temp1 = temp1->Next;
1198  }
1199  }
1200 
1201  return temp1;
1202 }
1203 
1204 
1205 int getFirstInbound( struct Node *list, struct Node *nodeOut)
1206 {
1207  struct Node *temp=NULL;
1208 
1209  temp=list;
1210 
1211  while(temp) {
1212  if( temp->inbound == 2 ) {
1213  copyNode(nodeOut, *temp);
1214  return 1;
1215  }
1216  temp=temp->Next;
1217  }
1218 
1219  return 0;
1220 }
1221 
1222 void getCoordinate(struct Node node, double *x, double *y, double *z)
1223 {
1224 
1225 
1226  *x = node.x;
1227  *y = node.y;
1228  *z = node.z;
1229 
1230 }
1231 
1232 void getCoordinates(struct Node *node, double *p)
1233 {
1234 
1235 
1236  p[0] = node->x;
1237  p[1] = node->y;
1238  p[2] = node->z;
1239 
1240 }
1241 
1242 void setCoordinate(struct Node *node, double x, double y, double z)
1243 {
1244 
1245 
1246  node->x = x;
1247  node->y = y;
1248  node->z = z;
1249 
1250 }
1251 
1252 /* set inbound value for the points in interList that has inbound =0,
1253  this will also set some inbound value of the points in list1
1254 */
1255 
1256 void setInbound(struct Node *interList, struct Node *list)
1257 {
1258 
1259  struct Node *temp1=NULL, *temp=NULL;
1260  struct Node *temp1_prev=NULL, *temp1_next=NULL;
1261  int prev_is_inside, next_is_inside;
1262 
1263  /* for each point in interList, search through list to decide the inbound value the interList point */
1264  /* For each inbound point, the prev node should be outside and the next is inside. */
1265  if(length(interList) == 0) return;
1266 
1267  temp = interList;
1268 
1269  while(temp) {
1270  if( !temp->inbound) {
1271  /* search in grid1 to find the prev and next point of temp, when prev point is outside and next point is inside
1272  inbound = 2, else inbound = 1*/
1273  temp1 = list;
1274  temp1_prev = NULL;
1275  temp1_next = NULL;
1276  while(temp1) {
1277  if(sameNode(*temp1, *temp)) {
1278  if(!temp1_prev) temp1_prev = getLast(list);
1279  temp1_next = temp1->Next;
1280  if(!temp1_next) temp1_next = list;
1281  break;
1282  }
1283  temp1_prev = temp1;
1284  temp1 = temp1->Next;
1285  }
1286  if(!temp1_next) error_handler("Error from create_xgrid.c: temp is not in list1");
1287  if( temp1_prev->isInside == 0 && temp1_next->isInside == 1)
1288  temp->inbound = 2; /* go inside */
1289  else
1290  temp->inbound = 1;
1291  }
1292  temp=temp->Next;
1293  }
1294 }
1295 
1296 int isInside(struct Node *node) {
1297 
1298  if(node->isInside == -1) error_handler("Error from mosaic_util.c: node->isInside is not set");
1299  return(node->isInside);
1300 
1301 }
1302 
1303 /* #define debug_test_create_xgrid */
1304 
1305 /* check if node is inside polygon list or not */
1306 int insidePolygon( struct Node *node, struct Node *list)
1307 {
1308  int is_inside;
1309  double pnt0[3], pnt1[3], pnt2[3];
1310  double anglesum;
1311  struct Node *p1=NULL, *p2=NULL;
1312 
1313  anglesum = 0;
1314 
1315  pnt0[0] = node->x;
1316  pnt0[1] = node->y;
1317  pnt0[2] = node->z;
1318 
1319  p1 = list;
1320  p2 = list->Next;
1321  is_inside = 0;
1322 
1323 
1324  while(p1) {
1325  pnt1[0] = p1->x;
1326  pnt1[1] = p1->y;
1327  pnt1[2] = p1->z;
1328  pnt2[0] = p2->x;
1329  pnt2[1] = p2->y;
1330  pnt2[2] = p2->z;
1331  if( samePoint(pnt0[0], pnt0[1], pnt0[2], pnt1[0], pnt1[1], pnt1[2]) ){
1332  return 1;
1333  }
1334  anglesum += spherical_angle(pnt0, pnt2, pnt1);
1335  p1 = p1->Next;
1336  p2 = p2->Next;
1337  if(p2==NULL){
1338  p2 = list;
1339  }
1340  }
1341 
1342  if( fabs(anglesum - 2*M_PI) < EPSLN8 ){
1343  is_inside = 1;
1344  }
1345  else{
1346  is_inside = 0;
1347  }
1348 
1349 #ifdef debug_test_create_xgrid
1350  printf("anglesum-2PI is %19.15f, is_inside = %d\n", anglesum- 2*M_PI, is_inside);
1351 #endif
1352 
1353  return is_inside;
1354 
1355 }
1356 
1357 int inside_a_polygon(double *lon1, double *lat1, int *npts, double *lon2, double *lat2)
1358 {
1359 
1360  double x2[20], y2[20], z2[20];
1361  double x1, y1, z1;
1362  double min_x2, max_x2, min_y2, max_y2, min_z2, max_z2;
1363  int isinside, i;
1364 
1365  struct Node *grid1=NULL, *grid2=NULL;
1366 
1367  /* first convert to cartesian grid */
1368  latlon2xyz(*npts, lon2, lat2, x2, y2, z2);
1369  latlon2xyz(1, lon1, lat1, &x1, &y1, &z1);
1370 
1371  max_x2 = maxval_double(*npts, x2);
1372  if(x1 >= max_x2+RANGE_CHECK_CRITERIA) return 0;
1373  min_x2 = minval_double(*npts, x2);
1374  if(min_x2 >= x1+RANGE_CHECK_CRITERIA) return 0;
1375 
1376  max_y2 = maxval_double(*npts, y2);
1377  if(y1 >= max_y2+RANGE_CHECK_CRITERIA) return 0;
1378  min_y2 = minval_double(*npts, y2);
1379  if(min_y2 >= y1+RANGE_CHECK_CRITERIA) return 0;
1380 
1381  max_z2 = maxval_double(*npts, z2);
1382  if(z1 >= max_z2+RANGE_CHECK_CRITERIA) return 0;
1383  min_z2 = minval_double(*npts, z2);
1384  if(min_z2 >= z1+RANGE_CHECK_CRITERIA) return 0;
1385 
1386 
1387  /* add x2,y2,z2 to a Node */
1388  rewindList();
1389  grid1 = getNext();
1390  grid2 = getNext();
1391 
1392  addEnd(grid1, x1, y1, z1, 0, 0, 0, -1);
1393  for(i=0; i<*npts; i++) addEnd(grid2, x2[i], y2[i], z2[i], 0, 0, 0, -1);
1394 
1395  isinside = insidePolygon(grid1, grid2);
1396 
1397  return isinside;
1398 
1399 }
1400 
1401 #ifndef __AIX
1402 int inside_a_polygon_(double *lon1, double *lat1, int *npts, double *lon2, double *lat2)
1403 {
1404 
1405  int isinside;
1406 
1407  isinside = inside_a_polygon(lon1, lat1, npts, lon2, lat2);
1408 
1409  return isinside;
1410 
1411 }
1412 #endif
int getFirstInbound(struct Node *list, struct Node *nodeOut)
Definition: mosaic_util.c:1205
int inbound
Definition: mosaic_util.h:39
real, parameter, public radius
Radius of the Earth [m].
Definition: constants.F90:72
l_size ! loop over number of fields ke do je do i
#define MAXNODELIST
Definition: mosaic_util.c:807
void error_handler(const char *msg)
Definition: mosaic_util.c:55
integer, parameter, public strlen
double gridArea(struct Node *grid)
Definition: mosaic_util.c:1156
real function intersect(x1, y1, x2, y2, x)
double metric(const double *p)
Definition: mosaic_util.c:671
void setInbound(struct Node *interList, struct Node *list)
Definition: mosaic_util.c:1256
real(fvprc), dimension(:), allocatable lon
#define SQRT_
#define RANGE_CHECK_CRITERIA
Definition: mosaic_util.h:29
integer, parameter, public ip
Definition: Type_Kinds.f90:97
double y
Definition: mosaic_util.h:37
#define EPSLN8
Definition: mosaic_util.c:32
int nearest_index(double value, const double *array, int ia)
Definition: mosaic_util.c:80
struct Node * getLast(struct Node *list)
Definition: mosaic_util.c:1190
int initialized
Definition: mosaic_util.h:40
#define EPSLN30
Definition: mosaic_util.c:35
void rewindList(void)
Definition: mosaic_util.c:813
int intersect_tri_with_line(const double *plane, const double *l1, const double *l2, double *p, double *t)
Definition: mosaic_util.c:727
void tokenize(const char *const string, const char *tokens, unsigned int varlen, unsigned int maxvar, char *pstring, unsigned int *const nstr)
Definition: mosaic_util.c:112
int isInside
Definition: mosaic_util.h:41
int clip_index
Definition: mosaic_util.h:43
*f90 *************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************! this routine is used to retrieve scalar boundary data for symmetric domain. subroutine MPP_GET_BOUNDARY_2D_(field, domain, ebuffer, sbuffer, wbuffer, nbuffer, flags, &position, complete, tile_count) type(domain2D), intent(in) ::domain MPP_TYPE_, intent(in) ::field(:,:) MPP_TYPE_, intent(inout), optional ::ebuffer(:), sbuffer(:), wbuffer(:), nbuffer(:) integer, intent(in), optional ::flags, position, tile_count logical, intent(in), optional ::complete MPP_TYPE_ ::field3D(size(field, 1), size(field, 2), 1) MPP_TYPE_, allocatable, dimension(:,:) ::ebuffer2D, sbuffer2D, wbuffer2D, nbuffer2D integer ::xcount, ycount integer ::ntile logical ::need_ebuffer, need_sbuffer, need_wbuffer, need_nbuffer integer(LONG_KIND), dimension(MAX_DOMAIN_FIELDS, MAX_TILES), save ::f_addrs=-9999 integer(LONG_KIND), dimension(4, MAX_DOMAIN_FIELDS, MAX_TILES), save ::b_addrs=-9999 integer, save ::bsize(4)=0, isize=0, jsize=0, ksize=0, pos, list=0, l_size=0, upflags integer ::buffer_size(4) integer ::max_ntile, tile, update_position, ishift, jshift logical ::do_update, is_complete, set_mismatch character(len=3) ::text MPP_TYPE_ ::d_type type(overlapSpec), pointer ::bound=> NULL() ntile
void vect_cross(const double *p1, const double *p2, double *e)
Definition: mosaic_util.c:648
real, dimension(:,:), allocatable temp1
integer, parameter, public double
Definition: Type_Kinds.f90:106
void addEnd(struct Node *list, double x, double y, double z, int intersect, double u, int inbound, int inside)
Definition: mosaic_util.c:855
void xyz2latlon(int np, const double *x, const double *y, const double *z, double *lon, double *lat)
Definition: mosaic_util.c:227
double spherical_angle(const double *v1, const double *v2, const double *v3)
Definition: mosaic_util.c:541
void copyNode(struct Node *node_out, struct Node node_in)
Definition: mosaic_util.c:1017
int inside_a_polygon_(double *lon1, double *lat1, int *npts, double *lon2, double *lat2)
Definition: mosaic_util.c:1402
str
Definition: c2f.py:15
real(r8), dimension(cast_m, cast_n) p
subroutine det(x1, y1, z1, x2, y2, z2, x0, y0, z0, output)
l_size ! loop over number of fields ke do j
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
struct Node * getNextNode(struct Node *list)
Definition: mosaic_util.c:1012
void set_reproduce_siena_true_(void)
Definition: mosaic_util.c:48
integer, parameter m
int insidePolygon(struct Node *node, struct Node *list)
Definition: mosaic_util.c:1306
double great_circle_area(int n, const double *x, const double *y, const double *z)
Definition: mosaic_util.c:504
double avgval_double(int size, const double *data)
Definition: mosaic_util.c:193
struct Node * getNode(struct Node *list, struct Node inNode)
Definition: mosaic_util.c:994
real(fp), parameter, public e
int subj_index
Definition: mosaic_util.h:42
#define HPI
Definition: mosaic_util.c:29
integer nvar
No description.
int insert_vtx(double x[], double y[], int n, int n_ins, double lon_in, double lat_in)
Definition: mosaic_util.c:386
double z
Definition: mosaic_util.h:37
#define TPI
Definition: mosaic_util.c:30
real, dimension(:,:), allocatable area
double u
Definition: mosaic_util.h:37
#define SMALL_VALUE
Definition: mosaic_util.h:34
real(kind=kind_real), parameter u1
double great_circle_distance(double *p1, double *p2)
Definition: mosaic_util.c:488
struct Node * getNext()
Definition: mosaic_util.c:823
double poly_area_dimensionless(const double x[], const double y[], int n)
Definition: mosaic_util.c:277
real, dimension(:,:), allocatable temp2
void v_print(double x[], double y[], int n)
Definition: mosaic_util.c:400
struct Node * nodeList
Definition: mosaic_util.c:810
double box_area(double ll_lon, double ll_lat, double ur_lon, double ur_lat)
Definition: mosaic_util.c:258
#define TOLORENCE
Definition: mosaic_util.c:31
real(fvprc), dimension(:), allocatable lat
int intersect
Definition: mosaic_util.h:38
int isIntersect(struct Node node)
Definition: mosaic_util.c:1178
#define RADIUS
Definition: constant.h:19
int sameNode(struct Node node1, struct Node node2)
Definition: mosaic_util.c:978
void insertIntersect(struct Node *list, double x, double y, double z, double u1, double u2, int inbound, double x2, double y2, double z2)
Definition: mosaic_util.c:1074
void addNode(struct Node *list, struct Node inNode)
Definition: mosaic_util.c:987
int length(struct Node *list)
Definition: mosaic_util.c:951
integer, parameter, public npts
int samePoint(double x1, double y1, double z1, double x2, double y2, double z2)
Definition: mosaic_util.c:968
int inside_a_polygon(double *lon1, double *lat1, int *npts, double *lon2, double *lat2)
Definition: mosaic_util.c:1357
void normalize_vect(double *e)
Definition: mosaic_util.c:679
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************subroutine READ_RECORD_CORE_(unit, field, nwords, data, start, axsiz) integer, intent(in) ::unit type(fieldtype), intent(in) ::field integer, intent(in) ::nwords MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in) ::start(:), axsiz(:) integer(SHORT_KIND) ::i2vals(nwords)!rab used in conjunction with transfer intrinsic to determine size of a variable integer(KIND=1) ::one_byte(8) integer ::word_sz!#ifdef __sgi integer(INT_KIND) ::ivals(nwords) real(FLOAT_KIND) ::rvals(nwords)!#else! integer ::ivals(nwords)! real ::rvals(nwords)!#endif real(DOUBLE_KIND) ::r8vals(nwords) pointer(ptr1, i2vals) pointer(ptr2, ivals) pointer(ptr3, rvals) pointer(ptr4, r8vals) if(mpp_io_stack_size< nwords) call mpp_io_set_stack_size(nwords) call mpp_error(FATAL, 'MPP_READ currently requires use_netCDF option') end subroutine READ_RECORD_CORE_ subroutine READ_RECORD_(unit, field, nwords, data, time_level, domain, position, tile_count, start_in, axsiz_in)!routine that is finally called by all mpp_read routines to perform the read!a non-netCDF record contains:! field ID! a set of 4 coordinates(is:ie, js:je) giving the data subdomain! a timelevel and a timestamp(=NULLTIME if field is static)! 3D real data(stored as 1D)!if you are using direct access I/O, the RECL argument to OPEN must be large enough for the above!in a global direct access file, record position on PE is given by %record.!Treatment of timestamp:! We assume that static fields have been passed without a timestamp.! Here that is converted into a timestamp of NULLTIME.! For non-netCDF fields, field is treated no differently, but is written! with a timestamp of NULLTIME. There is no check in the code to prevent! the user from repeatedly writing a static field. integer, intent(in) ::unit, nwords type(fieldtype), intent(in) ::field MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in), optional ::time_level type(domain2D), intent(in), optional ::domain integer, intent(in), optional ::position, tile_count integer, intent(in), optional ::start_in(:), axsiz_in(:) integer, dimension(size(field%axes(:))) ::start, axsiz integer ::tlevel !, subdomain(4) integer ::i, error, is, ie, js, je, isg, ieg, jsg, jeg type(domain2d), pointer ::io_domain=> tlevel if(PRESENT(start_in) .AND. PRESENT(axsiz_in)) then if(size(start(! the data domain and compute domain must refer to the subdomain being passed ! In this ! since that attempts to gather all data on PE size(field%axes(:)) axsiz(i)
int fix_lon(double x[], double y[], int n, double tlon)
Definition: mosaic_util.c:409
int intersectInList(struct Node *list, double x, double y, double z)
Definition: mosaic_util.c:1047
void set_reproduce_siena_true(void)
Definition: mosaic_util.c:42
#define ABS_
real(r8), dimension(cast_m, cast_n) t
double poly_area(const double x[], const double y[], int n)
Definition: mosaic_util.c:314
double x
Definition: mosaic_util.h:37
void getCoordinate(struct Node node, double *x, double *y, double *z)
Definition: mosaic_util.c:1222
double u_clip
Definition: mosaic_util.h:37
real, parameter p2
Definition: sw_core_nlm.F90:47
double minval_double(int size, const double *data)
Definition: mosaic_util.c:175
void setCoordinate(struct Node *node, double x, double y, double z)
Definition: mosaic_util.c:1242
double dot(const double *p1, const double *p2)
Definition: mosaic_util.c:663
int getInbound(struct Node node)
Definition: mosaic_util.c:1185
real function, private angle(t)
angle determines the position within the earth&#39;s orbit at time t in the year (t = 0 at NH autumnal eq...
Definition: astronomy.F90:2094
real, parameter p1
Definition: sw_core_nlm.F90:46
struct Node * Next
Definition: mosaic_util.h:44
void latlon2xyz(int size, const double *lon, const double *lat, double *x, double *y, double *z)
Definition: mosaic_util.c:211
int isInside(struct Node *node)
Definition: mosaic_util.c:1296
void mult(long double m[], long double v[], long double out_v[])
Definition: mosaic_util.c:768
int invert_matrix_3x3(long double m[], long double m_inv[])
Definition: mosaic_util.c:778
int curListPos
Definition: mosaic_util.c:811
double maxval_double(int size, const double *data)
Definition: mosaic_util.c:156
int addIntersect(struct Node *list, double x, double y, double z, int intersect, double u1, double u2, int inbound, int is1, int ie1, int is2, int ie2)
Definition: mosaic_util.c:894
logical function, public inside(p, lv, xv, yv, zv, nv, listv, ier)
double poly_area_no_adjust(const double x[], const double y[], int n)
Definition: mosaic_util.c:351
void initNode(struct Node *node)
Definition: mosaic_util.c:841
int delete_vtx(double x[], double y[], int n, int n_del)
Definition: mosaic_util.c:376
real(kind=kind_real), parameter u2
#define EPSLN15
Definition: mosaic_util.c:34
void getCoordinates(struct Node *node, double *p)
Definition: mosaic_util.c:1232
#define EPSLN10
Definition: mosaic_util.c:33
int reproduce_siena
Definition: mosaic_util.c:40
void printNode(struct Node *list, char *str)
Definition: mosaic_util.c:1031
void unit_vect_latlon(int size, const double *lon, const double *lat, double *vlon, double *vlat)
Definition: mosaic_util.c:697