source: XIOS/dev/branch_yushan_merged/extern/remap/src/polyg.cpp @ 1176

Last change on this file since 1176 was 1176, checked in by yushan, 7 years ago

save modif

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