source: trunk/SRC/Interpolation/inquad.pro @ 232

Last change on this file since 232 was 232, checked in by pinsard, 17 years ago

improvements/corrections of some *.pro headers

  • Property svn:keywords set to Id
File size: 9.8 KB
Line 
1;+
2;
3; @file_comments
4; to find if an (x,y) point is in a quadrilateral (x1,x2,x3,x4)
5;
6; @categories
7; Grid
8;
9; @param x {in}{required}
10; @param y {in}{required}
11; the coordinates of the point we want to know where it is.
12; Must be a scalar if /ONSPHERE activated else can be scalar or array.
13;
14; @param x1 {in}{required}
15; @param y1 {in}{required}
16; @param x2 {in}{required}
17; @param y2 {in}{required}
18; @param x3 {in}{required}
19; @param y3 {in}{required}
20; @param x4 {in}{required}
21; @param y4 {in}{required}
22; the coordinates of the quadrilateral given in the CLOCKWISE order.
23; Scalar or array.
24;
25; @keyword DOUBLE
26; use double precision to perform the computation
27;
28; @keyword ONSPHERE
29; to specify that the quadrilateral are on a sphere and
30; that their coordinates are longitude-latitude coordinates. In this
31; case, est-west periodicity, poles singularity and other pbs
32; related to longitude-latitude coordinates are managed
33; automatically.
34;
35; @keyword ZOOMRADIUS {default=4}
36; the zoom (circle centered on the (x,y) with a radius of
37; zoomradius degree where we look for the the quadrilateral which
38; contains the (x,y) point) used for the satellite projection
39; when /ONSPHERE is activated.
40; 4 seems to be the minimum which can be used.
41; Can be increase if the cell size is larger than 5 degrees.
42;
43; @keyword NOPRINT
44; to suppress the print messages.
45;
46; @keyword NEWCOORD
47;
48; @returns
49; a n element vector. Where n is the number of elements of
50; x. res[i]=j means that the point number i is located in the
51; quadrilateral number j with (0 <= j <= n_elements(x0)-1)
52;
53; @restrictions
54; I think degenerated quadrilateral (e.g. flat of twisted) is not work.
55; This has to be tested.
56;
57; @examples
58;
59; IDL> x = 1.*[1, 2, 6, 7, 3]
60; IDL> y = 1.*[1, 3, 3, 4, 7]
61; IDL> x1 = 1.*[0,4,2]
62; IDL> y1 = 1.*[1,4,8]
63; IDL> x2 = 1.*[1,6,4]
64; IDL> y2 = 1.*[5,6,8]
65; IDL> x3 = 1.*[3,8,4]
66; IDL> y3 = 1.*[4,4,6]
67; IDL> x4 = 1.*[2,6,2]
68; IDL> y4 = 1.*[0,2,6]
69; IDL> splot, [0,10], [0,10], xstyle = 1, ystyle = 1,/nodata
70; IDL> for i=0,2 do oplot, [x4[i],x1[i],x2[i],x3[i],x4[i]],[y4[i],y1[i],y2[i],y3[i],y4[i]]
71; IDL> oplot, x, y, color = 20, psym = 1, thick = 2
72; IDL> print, inquad(x, y, x1, y1, x2, y2, x3, y3, x4, y4)
73;
74; On a sphere see
75; <pro>clickincell</pro> ...
76;
77; @history
78;      Sebastien Masson (smasson\@lodyc.jussieu.fr)
79;      August 2003
80;      Based on Convert_clic_ij.pro written by Gurvan Madec
81;
82; @version
83; $Id$
84;
85;-
86;
87FUNCTION inquad, x, y, x1, y1, x2, y2, x3, y3, x4, y4, ONSPHERE = onsphere,  DOUBLE = double, ZOOMRADIUS = zoomradius, NOPRINT = noprint, NEWCOORD = newcoord
88;
89  compile_opt idl2, strictarrsubs
90;
91  ntofind = n_elements(x)
92  nquad = n_elements(x2)
93;
94  IF keyword_set(onsphere) THEN BEGIN
95; save the inputs parameters
96    xin = x
97    yin = y
98    x1in = x1
99    y1in = y1
100    x2in = x2
101    y2in = y2
102    x3in = x3
103    y3in = y3
104    x4in = x4
105    y4in = y4
106; for map_set
107    x = x MOD 360
108    x1 = x1 MOD 360
109    x2 = x2 MOD 360
110    x3 = x3 MOD 360
111    x4 = x4 MOD 360
112; save !map
113    save = {map:!map, x:!x, y:!y, z:!z, p:!p}
114; do a satellite projection
115    IF NOT keyword_set(zoomradius) THEN zoomradius = 4
116    map_set, y[0], x[0], 0, /satellite, sat_p = [1+zoomradius*20/6371.229, 0, 0], /noerase, /iso, /noborder
117; use normal coordinates to reject cells which are out of the projection.
118    tmp  = convert_coord(x, y, /DATA, /TO_NORMAL, DOUBLE = double)
119    tmp1 = convert_coord(x1, y1, /DATA, /TO_NORMAL, DOUBLE = double)
120    tmp2 = convert_coord(x2, y2, /DATA, /TO_NORMAL, DOUBLE = double)
121    tmp3 = convert_coord(x3, y3, /DATA, /TO_NORMAL, DOUBLE = double)
122    tmp4 = convert_coord(x4, y4, /DATA, /TO_NORMAL, DOUBLE = double)
123; remove cell which have one corner with coordinates equal to NaN
124    test = finite(tmp1[0, *]+tmp1[1, *]+tmp2[0, *]+tmp2[1, *] $
125                  +tmp3[0, *]+tmp3[1, *]+tmp4[0, *]+tmp4[1, *])
126    good = where(temporary(test) EQ 1)
127;
128    IF good[0] EQ -1 THEN BEGIN
129      IF NOT keyword_set(noprint) THEN print, 'The point is out of the cells'
130; restore the input parameters
131      x = temporary(xin)
132      y = temporary(yin)
133      x1 = temporary(x1in)
134      y1 = temporary(y1in)
135      x2 = temporary(x2in)
136      y2 = temporary(y2in)
137      x3 = temporary(x3in)
138      y3 = temporary(y3in)
139      x4 = temporary(x4in)
140      y4 = temporary(y4in)
141; restore old !map...
142      !map = save.map
143      !x = save.x
144      !y = save.y
145      !z = save.z
146      !p = save.p
147      RETURN,  -1
148    ENDIF
149;
150    x  = tmp[0]
151    y  = tmp[1]
152    x1 = tmp1[0, good]
153    y1 = tmp1[1, good]
154    x2 = tmp2[0, good]
155    y2 = tmp2[1, good]
156    x3 = tmp3[0, good]
157    y3 = tmp3[1, good]
158    x4 = tmp4[0, good]
159    y4 = tmp4[1, good]
160;
161    tmp1 = -1 & tmp2 = -1 & tmp3 = -1 & tmp4 = -1
162; remove cells which are obviously bad
163    test = (x1 GT x AND x2 GT x AND x3 GT x AND x4 GT x) $
164      OR (x1 LT x AND x2 LT x AND x3 LT x AND x4 LT x) $
165      OR (y1 GT y AND y2 GT y AND y3 GT y AND y4 GT y) $
166      OR (y1 LT y AND y2 LT y AND y3 LT y AND y4 LT y)
167    good2 = where(temporary(test) EQ 0)
168;
169    IF good2[0] EQ -1 THEN BEGIN
170      IF NOT keyword_set(noprint) THEN print, 'The point is out of the cells'
171; restore the input parameters
172      x = temporary(xin)
173      y = temporary(yin)
174      x1 = temporary(x1in)
175      y1 = temporary(y1in)
176      x2 = temporary(x2in)
177      y2 = temporary(y2in)
178      x3 = temporary(x3in)
179      y3 = temporary(y3in)
180      x4 = temporary(x4in)
181      y4 = temporary(y4in)
182; restore old !map...
183      !map = save.map
184      !x = save.x
185      !y = save.y
186      !z = save.z
187      !p = save.p
188      RETURN,  -1
189    ENDIF
190;
191    nquad = n_elements(good2)
192    x1 = x1[good2]
193    y1 = y1[good2]
194    x2 = x2[good2]
195    y2 = y2[good2]
196    x3 = x3[good2]
197    y3 = y3[good2]
198    x4 = x4[good2]
199    y4 = y4[good2]
200  ENDIF
201;
202;
203; the point is inside the quadrilateral if test eq 1
204; with test equal to:
205;     test = ((x-x1)*(y2-y1) GE (x2-x1)*(y-y1)) $
206;       *((x-x2)*(y3-y2) GT (x3-x2)*(y-y2)) $
207;       *((x-x3)*(y4-y3) GT (x4-x3)*(y-y3)) $
208;       *((x-x4)*(y1-y4) GE (x1-x4)*(y-y4))
209;
210; computation of test without any do loop for ntofind points (x,y) and
211; nquad quadrilateral((x1,x2,x3,x4),(y1,y2,y3,y4))
212; test dimensions are (ntofind, nquad)
213; column i of test corresponds to the intersection of point i with all
214; quadrilateral.
215; row j of test corresponds to all the points localized in cell j
216  test = $
217; (x-x1)
218  ((x[*]#replicate(1, nquad)-replicate(1, ntofind)#x1[*]) $
219; *(y2-y1)
220  *(replicate(1, ntofind)#(y2-y1)[*]) $
221; GE (x2-x1)
222  GE ((replicate(1, ntofind)#(x2-x1)[*]) $
223; *(y-y1)
224  *(y[*]#replicate(1, nquad)-replicate(1, ntofind)#y1[*])))
225;-------
226  test = temporary(test) $
227; *(x-x2)
228  *((x[*]#replicate(1, nquad)-replicate(1, ntofind)#x2[*]) $
229; *(y3-y2)
230  *(replicate(1, ntofind)#(y3-y2)[*]) $
231; GE (x3-x2)
232  GE ((replicate(1, ntofind)#(x3-x2)[*]) $
233; *(y-y2)
234  *(y[*]#replicate(1, nquad)-replicate(1, ntofind)#y2[*])))
235;-------
236  test = temporary(test) $
237; *(x-x3)
238  *((x[*]#replicate(1, nquad)-replicate(1, ntofind)#x3[*]) $
239; *(y4-y3)
240  *(replicate(1, ntofind)#(y4-y3)[*]) $
241; GE (x4-x3)
242  GE ((replicate(1, ntofind)#(x4-x3)[*]) $
243; *(y-y3)
244  *(y[*]#replicate(1, nquad)-replicate(1, ntofind)#y3[*])))
245;-------
246  test = temporary(test) $
247; *(x-x4)
248  *((x[*]#replicate(1, nquad)-replicate(1, ntofind)#x4[*]) $
249; *(y1-y4)
250  *(replicate(1, ntofind)#(y1-y4)[*]) $
251; GE (x1-x4)
252  GE ((replicate(1, ntofind)#(x1-x4)[*]) $
253; *(y-y4)
254  *(y[*]#replicate(1, nquad)-replicate(1, ntofind)#y4[*])))
255;
256; check test if ntofind gt 1
257; if ntofind gt 1, each point must be localised in one uniq cell.
258  IF ntofind GT 1 THEN BEGIN
259; each column of test must have only 1 position equal to one
260    chtest = total(test, 2)
261; points out of the cells
262    IF (where(chtest EQ 0))[0] NE -1 THEN BEGIN
263      IF NOT keyword_set(noprint) THEN print, 'Points number '+strjoin(strtrim(where(chtest EQ 0), 1), ', ')+' are out of the grid'
264      stop
265    ENDIF
266; points in more than one cell
267    IF (where(chtest GT 1))[0] NE -1 THEN BEGIN
268      IF NOT keyword_set(noprint) THEN print, 'Points number '+strjoin(strtrim(where(chtest GT 1), 1), ', ')+' are in more than one cell'
269      stop
270    ENDIF
271  ENDIF
272; find the points for which test eq 1
273  found = where(temporary(test) EQ 1)
274; if ntofind eq 1, the point may be localised in more than one grid
275; cell ou may also be out of the cells
276  IF ntofind EQ 1 THEN BEGIN
277    CASE 1 OF
278      found[0] EQ -1:BEGIN
279        IF NOT keyword_set(noprint) THEN print, 'The point is out of the cells'
280        IF keyword_set(onsphere) THEN BEGIN
281; restore old !map...
282          !map = save.map
283          !x = save.x
284          !y = save.y
285          !z = save.z
286          !p = save.p
287        ENDIF
288        return,  -1
289      END
290      n_elements(found) GT ntofind:BEGIN
291        IF NOT keyword_set(noprint) THEN print, 'The point is in more than one cell'
292      END
293      ELSE:
294    ENDCASE
295  ENDIF ELSE BEGIN
296; if ntofind GT 1, found must be sorted
297; i position of found. this corresponds to one x,y point
298    forsort = found MOD ntofind
299; j position of found. this corresponds to cell in which is one x,y
300; point
301    found = temporary(found)/ntofind
302; found must be sorted according to forsort
303    found = found[sort(forsort)]
304  ENDELSE
305;
306  IF keyword_set(onsphere) THEN BEGIN
307    IF arg_present(newcoord) THEN BEGIN
308      found = found[0]
309      newcoord = [[x1[found], y1[found]] $
310                  , [x2[found], y2[found]] $
311                  , [x3[found], y3[found]] $
312                  , [x4[found], y4[found]] $
313                  , [x, y]]
314    ENDIF
315;
316    found = good[good2[found]]
317; restore the input parameters
318    x = temporary(xin)
319    y = temporary(yin)
320    x1 = temporary(x1in)
321    y1 = temporary(y1in)
322    x2 = temporary(x2in)
323    y2 = temporary(y2in)
324    x3 = temporary(x3in)
325    y3 = temporary(y3in)
326    x4 = temporary(x4in)
327    y4 = temporary(y4in)
328; restore old !map...
329    !map = save.map
330    !x = save.x
331    !y = save.y
332    !z = save.z
333    !p = save.p
334  ENDIF
335;;
336  RETURN, found
337END
Note: See TracBrowser for help on using the repository browser.