source: XIOS/trunk/extern/remap/src/polyg.cpp @ 815

Last change on this file since 815 was 688, checked in by mhnguyen, 9 years ago

Integrating remap library into XIOS

+) Change name of some files of remap library to be compatible with XIOS
+) Implement function to fill in automatically boundary longitude and latitude

Test
+) On Curie
+) test_remap correct

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