source: XIOS2/dev/dev_ym/XIOS_COUPLING/extern/remap/src/polyg.cpp

Last change on this file was 2282, checked in by jderouillat, 2 years ago

Add missing return statements detected using -fsanitize=return. Return errors at runtime if reached.

File size: 9.6 KB
Line 
1/* utilities related to polygons */
2
3#include <cassert>
4#include <iostream>
5#include "elt.hpp"
6#include "errhandle.hpp"
7
8#include "polyg.hpp"
9
10namespace sphereRemap {
11
12using namespace std;
13
14/* given `N` `vertex`es, `N` `edge`s and `N` `d`s (d for small circles)
15   and `g` the barycentre,
16   this function reverses the order of arrays of `vertex`es, `edge`s and `d`s
17   but only if it is required!
18  (because computing intersection area requires both polygons to have same orientation)
19*/
20void orient(int N, Coord *vertex, Coord *edge, double *d, const Coord &g)
21{
22  Coord ga = vertex[0] - g;
23  Coord gb = vertex[1] - g;
24  Coord vertical = crossprod(ga, gb);
25  if (N > 2 && scalarprod(g, vertical) < 0)  // (GAxGB).G
26  {
27    for (int i = 0; i < N/2; i++)
28      swap(vertex[i], vertex[N-1-i]);
29
30    for (int i = 0; i < (N-1)/2; i++)
31    {
32      swap(edge[N-2-i], edge[i]);
33      swap(d[i], d[N-2-i]);
34    }
35  }
36}
37
38
39void normals(Coord *x, int n, Coord *a)
40{
41  for (int i = 0; i<n; i++)
42    a[i] = crossprod(x[(i+1)%n], x[i]);
43}
44
45Coord barycentre(const Coord *x, int n)
46{
47  if (n == 0) return ORIGIN;
48  Coord bc = ORIGIN;
49  for (int i = 0; i < n; i++)
50    bc = bc + x[i];
51  /* both distances can be equal down to roundoff when norm(bc) < mashineepsilon
52     which can occur when weighted with tiny area */
53
54  assert(squaredist(bc, proj(bc)) <= squaredist(bc, proj(bc * (-1.0))));
55  //if (squaredist(bc, proj(bc)) > squaredist(bc, proj(bc * (-1.0)))) return proj(bc * (-1.0));
56
57  return proj(bc);
58}
59
60/** computes the barycentre of the area which is the difference
61  of a side ABO of the spherical tetrahedron and of the straight tetrahedron */
62static Coord tetrah_side_diff_centre(Coord a, Coord b)
63{
64  Coord n = crossprod(a,b);
65        double sinc2 = n.x*n.x + n.y*n.y + n.z*n.z;
66  assert(sinc2 < 1.0 + EPS);
67
68        // exact: u = asin(sinc)/sinc - 1; asin(sinc) = geodesic length of arc ab
69  // approx:
70        // double u = sinc2/6. + (3./40.)*sinc2*sinc2;
71
72  // exact 
73  if (sinc2 > 1.0 - EPS) /* if round-off leads to sinc > 1 asin produces NaN */
74    return n * (M_PI_2 - 1);
75  double sinc = sqrt(sinc2);
76        double u = asin(sinc)/sinc - 1;
77
78        return n*u;
79}
80
81/* compute the barycentre as the negative sum of all the barycentres of sides of the tetrahedron */
82Coord gc_normalintegral(const Coord *x, int n)
83{
84  Coord m = barycentre(x, n);
85  Coord bc = crossprod(x[n-1]-m, x[0]-m) + tetrah_side_diff_centre(x[n-1], x[0]);
86  for (int i = 1; i < n; i++)
87    bc = bc + crossprod(x[i-1]-m, x[i]-m) + tetrah_side_diff_centre(x[i-1], x[i]);
88  return bc*0.5;
89}
90
91Coord exact_barycentre(const Coord *x, int n)
92{
93  if (n >= 3)
94  {
95    return  proj(gc_normalintegral(x, n));
96  }
97  else if (n == 0) return ORIGIN;
98  else if (n == 2) return midpoint(x[0], x[1]);
99  else if (n == 1) return x[0];
100
101  error_exit( "Missing return in : Coord exact_barycentre(const Coord *x, int n)" );
102  return ORIGIN;
103}
104
105Coord sc_gc_moon_normalintegral(Coord a, Coord b, Coord pole)
106{
107  double hemisphere = (a.z > 0) ? 1: -1;
108
109  double lat = hemisphere * (M_PI_2 - acos(a.z));
110  double lon1 = atan2(a.y, a.x);
111  double lon2 = atan2(b.y, b.x);
112  double lon_diff = lon2 - lon1;
113
114  // wraparound at lon=-pi=pi
115  if (lon_diff < -M_PI) lon_diff += 2.0*M_PI;
116  else if (lon_diff > M_PI) lon_diff -= 2.0*M_PI;
117
118  // integral of the normal over the surface bound by great arcs a-pole and b-pole and small arc a-b
119  Coord sc_normalintegral = Coord(0.5*(sin(lon2)-sin(lon1))*(M_PI_2 - lat - 0.5*sin(2.0*lat)),
120                                  0.5*(cos(lon1)-cos(lon2))*(M_PI_2 - lat - 0.5*sin(2.0*lat)),
121                                  hemisphere * lon_diff * 0.25 * (cos(2.0*lat) + 1.0));
122  Coord p = Coord(0,0,hemisphere); // TODO assumes north pole is (0,0,1)
123  Coord t[] = {a, b, p};
124  if (hemisphere < 0) swap(t[0], t[1]);
125  return (sc_normalintegral - gc_normalintegral(t, 3)) * hemisphere;
126}
127
128
129double triarea(const Coord& A, const Coord& B, const Coord& C)
130{
131  double a = ds(B, C);
132  double b = ds(C, A);
133  double c = ds(A, B);
134  double tmp ;
135
136  if (a<b) { tmp=a ; a=b ; b=tmp ; }
137  if (c > a) { tmp=a ; a=c ; c=b, b=tmp;  }
138  else if (c > b) { tmp=c ; c=b ; b=tmp ; }
139 
140  double s = 0.5 * (a + b + c);
141  double t = tan(0.25*(a+(b+c))) * tan(0.25*(c-(a-b))) * tan(0.25*( c + (a-b) )) * tan(0.25*( a + (b - c)));
142  if (t>0) return 4 * atan(sqrt(t));
143  else
144  {
145    std::cout<<"double triarea(const Coord& A, const Coord& B, const Coord& C) : t < 0 !!! t="<<t<<endl ;
146    return 0 ;
147  }
148}
149
150/** Computes area of two two-sided polygon
151   needs to have one small and one great circle, otherwise zero
152   (name origin: lun is moon in french)
153*/
154double alun(double b, double d)
155{
156  double a  = acos(d);
157  assert(b <= 2 * a);
158  double s  = a + 0.5 * b;
159  double t = tan(0.5 * s) * tan(0.5 * (s - a)) * tan(0.5 * (s - a)) * tan(0.5 * (s - b));
160  double r  = sqrt(1 - d*d);
161  double p  = 2 * asin(sin(0.5*b) / r);
162  return p*(1 - d) - 4*atan(sqrt(t));
163}
164
165/**
166  This function returns the area of a spherical element
167that can be composed of great and small circle arcs.
168 The caller must ensure this function is not called when `alun` should be called instaed.
169  This function also sets `gg` to the barycentre of the element.
170  "air" stands for area and  "bar" for barycentre.
171*/
172double airbar(int N, const Coord *x, const Coord *c, double *d, const Coord& pole, Coord& gg)
173{
174  if (N < 3)
175    return 0; /* polygons with less then three vertices have zero area */
176  Coord t[3];
177  t[0] = barycentre(x, N);
178  Coord *g = new Coord[N];
179  double area = 0;
180  Coord gg_exact = gc_normalintegral(x, N);
181  for (int i = 0; i < N; i++)
182  {
183    /* for "spherical circle segment" sum triangular part and "small moon" and => account for small circle */
184    int ii = (i + 1) % N;
185    t[1] = x[i];
186    t[2] = x[ii];
187                double sc=scalarprod(crossprod(t[1] - t[0], t[2] - t[0]), t[0]) ;
188    assert(sc >= -1e-10); // Error: tri a l'env (wrong orientation)
189    double area_gc = triarea(t[0], t[1], t[2]);
190    double area_sc_gc_moon = 0;
191    if (d[i]) /* handle small circle case */
192    {
193      Coord m = midpoint(t[1], t[2]);
194      double mext = scalarprod(m, c[i]) - d[i];
195      char sgl = (mext > 0) ? -1 : 1;
196      area_sc_gc_moon = sgl * alun(arcdist(t[1], t[2]), fabs(scalarprod(t[1], pole)));
197      gg_exact = gg_exact + sc_gc_moon_normalintegral(t[1], t[2], pole);
198    }
199    area += area_gc + area_sc_gc_moon; /* for "spherical circle segment" sum triangular part (at) and "small moon" and => account for small circle */
200    g[i] = barycentre(t, 3) * (area_gc + area_sc_gc_moon);
201  }
202  gg = barycentre(g, N);
203  gg_exact = proj(gg_exact);
204  delete[] g;
205  gg = gg_exact;
206  return area;
207}
208
209double polygonarea(Coord *vertices, int N)
210{
211  assert(N >= 3);
212
213  /* compute polygon area as sum of triangles */
214  Coord centre = barycentre(vertices, N);
215  double area = 0;
216  for (int i = 0; i < N; i++)
217    area += triarea(centre, vertices[i], vertices[(i+1)%N]);
218  return area;
219}
220
221int packedPolygonSize(const Elt& e)
222{
223  return sizeof(e.id) + sizeof(e.src_id) + sizeof(e.x) + sizeof(e.val)  + sizeof(e.given_area)+
224         sizeof(e.n) + e.n*(sizeof(double)+sizeof(Coord));
225}
226
227void packPolygon(const Elt& e, char *buffer, int& pos) 
228{
229  *((GloId *) &(buffer[pos])) = e.id;
230  pos += sizeof(e.id);
231  *((GloId *) &(buffer[pos])) = e.src_id;
232  pos += sizeof(e.src_id);
233
234  *((Coord *) &(buffer[pos])) = e.x;
235  pos += sizeof(e.x);
236
237  *((double*) &(buffer[pos])) = e.val;
238  pos += sizeof(e.val);
239
240  *((double*) &(buffer[pos])) = e.given_area;
241  pos += sizeof(e.val);
242
243  *((int *) & (buffer[pos])) = e.n;
244  pos += sizeof(e.n);
245
246  for (int i = 0; i < e.n; i++)
247  {
248    *((double *) & (buffer[pos])) = e.d[i];
249    pos += sizeof(e.d[i]);
250
251    *((Coord *) &(buffer[pos])) = e.vertex[i];
252    pos += sizeof(e.vertex[i]);
253  } 
254
255}
256
257void unpackPolygon(Elt& e, const char *buffer, int& pos) 
258{
259  e.id = *((GloId *) &(buffer[pos]));
260  pos += sizeof(e.id);
261  e.src_id = *((GloId *) &(buffer[pos]));
262  pos += sizeof(e.src_id);
263
264  e.x = *((Coord *) & (buffer[pos]) );
265  pos += sizeof(e.x);
266
267  e.val = *((double *) & (buffer[pos]));
268  pos += sizeof(double);
269
270  e.given_area = *((double *) & (buffer[pos]));
271  pos += sizeof(double);
272
273  e.n = *((int *) & (buffer[pos]));
274  pos += sizeof(int);
275
276  for (int i = 0; i < e.n; i++)
277  {
278    e.d[i] = *((double *) & (buffer[pos]));
279    pos += sizeof(double);
280
281    e.vertex[i] = *((Coord *) & (buffer[pos]));
282    pos += sizeof(Coord);
283  }
284}
285
286int packIntersectionSize(const Elt& elt) 
287{
288  return elt.is.size() * (2*sizeof(int)+ sizeof(GloId) + 5*sizeof(double));
289}
290
291void packIntersection(const Elt& e, char* buffer,int& pos) 
292{
293  for (list<Polyg *>::const_iterator it = e.is.begin(); it != e.is.end(); ++it)
294  {
295    *((int *) &(buffer[0])) += 1;
296
297    *((int *) &(buffer[pos])) = e.id.ind;
298    pos += sizeof(int);
299
300    *((double *) &(buffer[pos])) = e.area;
301    pos += sizeof(double);
302
303
304    *((GloId *) &(buffer[pos])) = (*it)->id;
305    pos += sizeof(GloId);
306 
307    *((int *) &(buffer[pos])) = (*it)->n;
308    pos += sizeof(int);
309    *((double *) &(buffer[pos])) = (*it)->area;
310    pos += sizeof(double);
311
312    *((Coord *) &(buffer[pos])) = (*it)->x;
313    pos += sizeof(Coord);
314  }
315}
316
317void unpackIntersection(Elt* e, const char* buffer) 
318{
319  int ind;
320  int pos = 0;
321 
322  int n = *((int *) & (buffer[pos]));
323  pos += sizeof(int);
324  for (int i = 0; i < n; i++)
325  {
326    ind = *((int*) &(buffer[pos]));
327    pos+=sizeof(int);
328
329    Elt& elt= e[ind];
330
331    elt.area=*((double *) & (buffer[pos]));
332    pos += sizeof(double);
333
334
335    Polyg *polygon = new Polyg;
336
337    polygon->id =  *((GloId *) & (buffer[pos]));
338    pos += sizeof(GloId);
339
340    polygon->n =  *((int *) & (buffer[pos]));
341    pos += sizeof(int);
342
343    polygon->area =  *((double *) & (buffer[pos]));
344    pos += sizeof(double);
345
346    polygon->x = *((Coord *) & (buffer[pos]));
347    pos += sizeof(Coord);
348
349    elt.is.push_back(polygon);
350  }
351}
352
353}
Note: See TracBrowser for help on using the repository browser.