source: XIOS3/trunk/extern/remap/src/polyg.cpp @ 2537

Last change on this file since 2537 was 2537, checked in by jderouillat, 11 months ago

Add a missing return statment in the remapper

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