source: XIOS3/branches/xios-3.0-beta/extern/remap/src/polyg.cpp

Last change on this file was 2536, checked in by jderouillat, 12 months ago

Fix some minor merge error

File size: 10.4 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
131
132
133Coord sc_gc_moon_normalintegral(Coord a, Coord b, Coord pole)
134{
135  double hemisphere = (a.z > 0) ? 1: -1;
136
137  double lat = hemisphere * (M_PI_2 - acos(a.z));
138  double lon1 = atan2(a.y, a.x);
139  double lon2 = atan2(b.y, b.x);
140  double lon_diff = lon2 - lon1;
141
142  // wraparound at lon=-pi=pi
143  if (lon_diff < -M_PI) lon_diff += 2.0*M_PI;
144  else if (lon_diff > M_PI) lon_diff -= 2.0*M_PI;
145
146  // integral of the normal over the surface bound by great arcs a-pole and b-pole and small arc a-b
147  Coord sc_normalintegral = Coord(0.5*(sin(lon2)-sin(lon1))*(M_PI_2 - lat - 0.5*sin(2.0*lat)),
148                                  0.5*(cos(lon1)-cos(lon2))*(M_PI_2 - lat - 0.5*sin(2.0*lat)),
149                                  hemisphere * lon_diff * 0.25 * (cos(2.0*lat) + 1.0));
150  Coord p = Coord(0,0,hemisphere); // TODO assumes north pole is (0,0,1)
151  Coord t[] = {a, b, p};
152  if (hemisphere < 0) swap(t[0], t[1]);
153  return (sc_normalintegral - gc_normalintegral(t, 3)) * hemisphere;
154}
155
156
157double triarea(const Coord& A, const Coord& B, const Coord& C)
158{
159  double a = ds(B, C);
160  double b = ds(C, A);
161  double c = ds(A, B);
162  double tmp ;
163
164  if (a<b) { tmp=a ; a=b ; b=tmp ; }
165  if (c > a) { tmp=a ; a=c ; c=b, b=tmp;  }
166  else if (c > b) { tmp=c ; c=b ; b=tmp ; }
167 
168  double s = 0.5 * (a + b + c);
169  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)));
170  if (t>0) return 4 * atan(sqrt(t));
171  else
172  {
173    std::cout<<"double triarea(const Coord& A, const Coord& B, const Coord& C) : t < 0 !!! t="<<t<<endl ;
174    return 0 ;
175  }
176}
177
178/** Computes area of two two-sided polygon
179   needs to have one small and one great circle, otherwise zero
180   (name origin: lun is moon in french)
181*/
182double alun(double b, double d)
183{
184  double a  = acos(d);
185  assert(b <= 2 * a);
186  double s  = a + 0.5 * b;
187  double t = tan(0.5 * s) * tan(0.5 * (s - a)) * tan(0.5 * (s - a)) * tan(0.5 * (s - b));
188  double r  = sqrt(1 - d*d);
189  double p  = 2 * asin(sin(0.5*b) / r);
190  return p*(1 - d) - 4*atan(sqrt(t));
191}
192
193/**
194  This function returns the area of a spherical element
195that can be composed of great and small circle arcs.
196 The caller must ensure this function is not called when `alun` should be called instaed.
197  This function also sets `gg` to the barycentre of the element.
198  "air" stands for area and  "bar" for barycentre.
199*/
200double airbar(int N, const Coord *x, const Coord *c, double *d, const Coord& pole, Coord& gg)
201{
202  if (N < 3)
203    return 0; /* polygons with less then three vertices have zero area */
204  Coord t[3];
205  t[0] = barycentre(x, N);
206  Coord *g = new Coord[N];
207  double area = 0;
208  Coord gg_exact = gc_normalintegral(x, N);
209  for (int i = 0; i < N; i++)
210  {
211    /* for "spherical circle segment" sum triangular part and "small moon" and => account for small circle */
212    int ii = (i + 1) % N;
213    t[1] = x[i];
214    t[2] = x[ii];
215                double sc=scalarprod(crossprod(t[1] - t[0], t[2] - t[0]), t[0]) ;
216    assert(sc >= -1e-10); // Error: tri a l'env (wrong orientation)
217    double area_gc = triarea(t[0], t[1], t[2]);
218    double area_sc_gc_moon = 0;
219    if (d[i]) /* handle small circle case */
220    {
221      Coord m = midpoint(t[1], t[2]);
222      double mext = scalarprod(m, c[i]) - d[i];
223      char sgl = (mext > 0) ? -1 : 1;
224      area_sc_gc_moon = sgl * alun(arcdist(t[1], t[2]), fabs(scalarprod(t[1], pole)));
225      gg_exact = gg_exact + sc_gc_moon_normalintegral(t[1], t[2], pole);
226    }
227    area += area_gc + area_sc_gc_moon; /* for "spherical circle segment" sum triangular part (at) and "small moon" and => account for small circle */
228    g[i] = barycentre(t, 3) * (area_gc + area_sc_gc_moon);
229  }
230  gg = barycentre(g, N);
231  gg_exact = proj(gg_exact);
232  delete[] g;
233  gg = gg_exact;
234  return area;
235}
236
237double polygonarea(Coord *vertices, int N)
238{
239  assert(N >= 3);
240
241  /* compute polygon area as sum of triangles */
242  Coord centre = barycentre(vertices, N);
243  double area = 0;
244  for (int i = 0; i < N; i++)
245    area += triarea(centre, vertices[i], vertices[(i+1)%N]);
246  return area;
247}
248
249int packedPolygonSize(const Elt& e)
250{
251  return sizeof(e.id) + sizeof(e.src_id) + sizeof(e.x) + sizeof(e.val)  + sizeof(e.given_area)+
252         sizeof(e.n) + e.n*(sizeof(double)+sizeof(Coord));
253}
254
255void packPolygon(const Elt& e, char *buffer, int& pos) 
256{
257  *((GloId *) &(buffer[pos])) = e.id;
258  pos += sizeof(e.id);
259  *((GloId *) &(buffer[pos])) = e.src_id;
260  pos += sizeof(e.src_id);
261
262  *((Coord *) &(buffer[pos])) = e.x;
263  pos += sizeof(e.x);
264
265  *((double*) &(buffer[pos])) = e.val;
266  pos += sizeof(e.val);
267
268  *((double*) &(buffer[pos])) = e.given_area;
269  pos += sizeof(e.val);
270
271  *((int *) & (buffer[pos])) = e.n;
272  pos += sizeof(e.n);
273
274  for (int i = 0; i < e.n; i++)
275  {
276    *((double *) & (buffer[pos])) = e.d[i];
277    pos += sizeof(e.d[i]);
278
279    *((Coord *) &(buffer[pos])) = e.vertex[i];
280    pos += sizeof(e.vertex[i]);
281  } 
282
283}
284
285void unpackPolygon(Elt& e, const char *buffer, int& pos) 
286{
287  e.id = *((GloId *) &(buffer[pos]));
288  pos += sizeof(e.id);
289  e.src_id = *((GloId *) &(buffer[pos]));
290  pos += sizeof(e.src_id);
291
292  e.x = *((Coord *) & (buffer[pos]) );
293  pos += sizeof(e.x);
294
295  e.val = *((double *) & (buffer[pos]));
296  pos += sizeof(double);
297
298  e.given_area = *((double *) & (buffer[pos]));
299  pos += sizeof(double);
300
301  e.n = *((int *) & (buffer[pos]));
302  pos += sizeof(int);
303
304  for (int i = 0; i < e.n; i++)
305  {
306    e.d[i] = *((double *) & (buffer[pos]));
307    pos += sizeof(double);
308
309    e.vertex[i] = *((Coord *) & (buffer[pos]));
310    pos += sizeof(Coord);
311  }
312}
313
314int packIntersectionSize(const Elt& elt) 
315{
316  return elt.is.size() * (2*sizeof(int)+ sizeof(GloId) + 5*sizeof(double));
317}
318
319void packIntersection(const Elt& e, char* buffer,int& pos) 
320{
321  for (list<Polyg *>::const_iterator it = e.is.begin(); it != e.is.end(); ++it)
322  {
323    *((int *) &(buffer[0])) += 1;
324
325    *((int *) &(buffer[pos])) = e.id.ind;
326    pos += sizeof(int);
327
328    *((double *) &(buffer[pos])) = e.area;
329    pos += sizeof(double);
330
331
332    *((GloId *) &(buffer[pos])) = (*it)->id;
333    pos += sizeof(GloId);
334 
335    *((int *) &(buffer[pos])) = (*it)->n;
336    pos += sizeof(int);
337    *((double *) &(buffer[pos])) = (*it)->area;
338    pos += sizeof(double);
339
340    *((Coord *) &(buffer[pos])) = (*it)->x;
341    pos += sizeof(Coord);
342  }
343}
344
345void unpackIntersection(Elt* e, const char* buffer) 
346{
347  int ind;
348  int pos = 0;
349 
350  int n = *((int *) & (buffer[pos]));
351  pos += sizeof(int);
352  for (int i = 0; i < n; i++)
353  {
354    ind = *((int*) &(buffer[pos]));
355    pos+=sizeof(int);
356
357    Elt& elt= e[ind];
358
359    elt.area=*((double *) & (buffer[pos]));
360    pos += sizeof(double);
361
362
363    Polyg *polygon = new Polyg;
364
365    polygon->id =  *((GloId *) & (buffer[pos]));
366    pos += sizeof(GloId);
367
368    polygon->n =  *((int *) & (buffer[pos]));
369    pos += sizeof(int);
370
371    polygon->area =  *((double *) & (buffer[pos]));
372    pos += sizeof(double);
373
374    polygon->x = *((Coord *) & (buffer[pos]));
375    pos += sizeof(Coord);
376
377    elt.is.push_back(polygon);
378  }
379}
380
381}
Note: See TracBrowser for help on using the repository browser.