1 /* zzhullax.f -- translated by f2c (version 19980913).
2    You must link the resulting object file with the libraries:
3 	-lf2c -lm   (in that order)
4 */
5 
6 #include "f2c.h"
7 
8 /* Table of constant values */
9 
10 static doublereal c_b20 = -1.;
11 static doublereal c_b36 = .5;
12 
13 /* $Procedure   ZZHULLAX ( Pyramidal FOV convex hull to FOV axis ) */
zzhullax_(char * inst,integer * n,doublereal * bounds,doublereal * axis,ftnlen inst_len)14 /* Subroutine */ int zzhullax_(char *inst, integer *n, doublereal *bounds,
15 	doublereal *axis, ftnlen inst_len)
16 {
17     /* System generated locals */
18     integer bounds_dim2, i__1, i__2;
19     doublereal d__1;
20 
21     /* Builtin functions */
22     integer s_rnge(char *, integer, char *, integer);
23 
24     /* Local variables */
25     extern /* Subroutine */ int vhat_(doublereal *, doublereal *);
26     doublereal xvec[3], yvec[3], zvec[3];
27     integer xidx;
28     extern doublereal vsep_(doublereal *, doublereal *);
29     integer next;
30     logical pass1;
31     integer i__, m;
32     doublereal r__, v[3], delta;
33     extern /* Subroutine */ int chkin_(char *, ftnlen), errch_(char *, char *,
34 	     ftnlen, ftnlen);
35     logical found;
36     extern /* Subroutine */ int errdp_(char *, doublereal *, ftnlen), vlcom_(
37 	    doublereal *, doublereal *, doublereal *, doublereal *,
38 	    doublereal *);
39     integer minix, maxix;
40     doublereal trans[9]	/* was [3][3] */;
41     extern /* Subroutine */ int ucrss_(doublereal *, doublereal *, doublereal
42 	    *), vcrss_(doublereal *, doublereal *, doublereal *);
43     extern logical vzero_(doublereal *);
44     extern /* Subroutine */ int vrotv_(doublereal *, doublereal *, doublereal
45 	    *, doublereal *);
46     doublereal cp[3];
47     extern doublereal pi_(void);
48     logical ok;
49     extern doublereal halfpi_(void);
50     extern /* Subroutine */ int reclat_(doublereal *, doublereal *,
51 	    doublereal *, doublereal *), sigerr_(char *, ftnlen);
52     doublereal minlon;
53     extern /* Subroutine */ int chkout_(char *, ftnlen);
54     doublereal maxlon;
55     extern /* Subroutine */ int vhatip_(doublereal *), vsclip_(doublereal *,
56 	    doublereal *), setmsg_(char *, ftnlen), errint_(char *, integer *,
57 	     ftnlen);
58     extern logical return_(void);
59     doublereal lat, sep, lon;
60     extern /* Subroutine */ int mxv_(doublereal *, doublereal *, doublereal *)
61 	    ;
62     doublereal ray1[3], ray2[3];
63 
64 /* $ Abstract */
65 
66 /*     SPICE Private routine intended solely for the support of SPICE */
67 /*     routines.  Users should not call this routine directly due */
68 /*     to the volatile nature of this routine. */
69 
70 /*     Identify a face of the convex hull of an instrument's */
71 /*     polygonal FOV, and use this face to generate an axis of the */
72 /*     FOV. */
73 
74 /* $ Disclaimer */
75 
76 /*     THIS SOFTWARE AND ANY RELATED MATERIALS WERE CREATED BY THE */
77 /*     CALIFORNIA INSTITUTE OF TECHNOLOGY (CALTECH) UNDER A U.S. */
78 /*     GOVERNMENT CONTRACT WITH THE NATIONAL AERONAUTICS AND SPACE */
79 /*     ADMINISTRATION (NASA). THE SOFTWARE IS TECHNOLOGY AND SOFTWARE */
80 /*     PUBLICLY AVAILABLE UNDER U.S. EXPORT LAWS AND IS PROVIDED "AS-IS" */
81 /*     TO THE RECIPIENT WITHOUT WARRANTY OF ANY KIND, INCLUDING ANY */
82 /*     WARRANTIES OF PERFORMANCE OR MERCHANTABILITY OR FITNESS FOR A */
83 /*     PARTICULAR USE OR PURPOSE (AS SET FORTH IN UNITED STATES UCC */
84 /*     SECTIONS 2312-2313) OR FOR ANY PURPOSE WHATSOEVER, FOR THE */
85 /*     SOFTWARE AND RELATED MATERIALS, HOWEVER USED. */
86 
87 /*     IN NO EVENT SHALL CALTECH, ITS JET PROPULSION LABORATORY, OR NASA */
88 /*     BE LIABLE FOR ANY DAMAGES AND/OR COSTS, INCLUDING, BUT NOT */
89 /*     LIMITED TO, INCIDENTAL OR CONSEQUENTIAL DAMAGES OF ANY KIND, */
90 /*     INCLUDING ECONOMIC DAMAGE OR INJURY TO PROPERTY AND LOST PROFITS, */
91 /*     REGARDLESS OF WHETHER CALTECH, JPL, OR NASA BE ADVISED, HAVE */
92 /*     REASON TO KNOW, OR, IN FACT, SHALL KNOW OF THE POSSIBILITY. */
93 
94 /*     RECIPIENT BEARS ALL RISK RELATING TO QUALITY AND PERFORMANCE OF */
95 /*     THE SOFTWARE AND ANY RELATED MATERIALS, AND AGREES TO INDEMNIFY */
96 /*     CALTECH AND NASA FOR ALL THIRD-PARTY CLAIMS RESULTING FROM THE */
97 /*     ACTIONS OF RECIPIENT IN THE USE OF THE SOFTWARE. */
98 
99 /* $ Required_Reading */
100 
101 /*     CK */
102 /*     FRAMES */
103 /*     GF */
104 /*     IK */
105 /*     KERNEL */
106 
107 /* $ Keywords */
108 
109 /*     FOV */
110 /*     GEOMETRY */
111 /*     INSTRUMENT */
112 
113 /* $ Declarations */
114 /* $ Brief_I/O */
115 
116 /*     VARIABLE  I/O  DESCRIPTION */
117 /*     --------  ---  -------------------------------------------------- */
118 /*     MARGIN     P   Minimum complement of FOV cone angle. */
119 /*     INST       I   Instrument name. */
120 /*     N          I   Number of FOV boundary vectors. */
121 /*     BOUNDS     I   FOV boundary vectors. */
122 /*     AXIS       O   Instrument FOV axis vector. */
123 
124 /* $ Detailed_Input */
125 
126 /*     INST       is the name of an instrument with which the field of */
127 /*                view (FOV) of interest is associated. This name is */
128 /*                used only to generate long error messages. */
129 
130 /*     N          is the number of boundary vectors in the array */
131 /*                BOUNDS. */
132 
133 /*     BOUNDS     is an array of N vectors emanating from a common */
134 /*                vertex and defining the edges of a pyramidal region in */
135 /*                three-dimensional space: this the region within the */
136 /*                FOV of the instrument designated by INST. The Ith */
137 /*                vector of BOUNDS resides in elements (1:3,I) of this */
138 /*                array. */
139 
140 /*                The vectors contained in BOUNDS are called the */
141 /*                "boundary vectors" of the FOV. */
142 
143 /*                The boundary vectors  must satisfy the constraints: */
144 
145 /*                   1)  The boundary vectors  must be contained within */
146 /*                       a right circular cone of angular radius less */
147 /*                       than than (pi/2) - MARGIN radians; in other */
148 /*                       words, there must be a vector A such that all */
149 /*                       boundary vectors have angular separation from */
150 /*                       A of less than (pi/2)-MARGIN radians. */
151 
152 /*                   2)  There must be a pair of vectors U, V in BOUNDS */
153 /*                       such that all other boundary vectors lie in */
154 /*                       the same half space bounded by the plane */
155 /*                       containing U and V. Furthermore, all other */
156 /*                       boundary vectors must have orthogonal */
157 /*                       projections onto a plane normal to this plane */
158 /*                       such that the projections have angular */
159 /*                       separation of at least 2*MARGIN radians from */
160 /*                       the plane spanned by U and V. */
161 
162 /*                Given the first constraint above, there is plane PL */
163 /*                such that each of the set of rays extending the */
164 /*                boundary vectors intersects PL. (In fact, there is an */
165 /*                infinite set of such planes.) The boundary vectors */
166 /*                must be ordered so that the set of line segments */
167 /*                connecting the intercept on PL of the ray extending */
168 /*                the Ith vector to that of the (I+1)st, with the Nth */
169 /*                intercept connected to the first, form a polygon (the */
170 /*                "FOV polygon") constituting the intersection of the */
171 /*                FOV pyramid with PL. This polygon may wrap in either */
172 /*                the positive or negative sense about a ray emanating */
173 /*                from the FOV vertex and passing through the plane */
174 /*                region bounded by the FOV polygon. */
175 
176 /*                The FOV polygon need not be convex; it may be */
177 /*                self-intersecting as well. */
178 
179 /*                No pair of consecutive vectors in BOUNDS may be */
180 /*                linearly dependent. */
181 
182 /*                The boundary vectors need not have unit length. */
183 
184 
185 /* $ Detailed_Output */
186 
187 /*     AXIS       is a unit vector normal to a plane containing the */
188 /*                FOV polygon. All boundary vectors have angular */
189 /*                separation from AXIS of not more than */
190 
191 /*                   ( pi/2 ) - MARGIN */
192 
193 /*                radians. */
194 
195 /*                This routine signals an error if it cannot find */
196 /*                a satisfactory value of AXIS. */
197 
198 /* $ Parameters */
199 
200 /*     MARGIN     is a small positive number used to constrain the */
201 /*                orientation of the boundary vectors. See the two */
202 /*                constraints described in the Detailed_Input section */
203 /*                above for specifics. */
204 
205 /* $ Exceptions */
206 
207 /*     1)  In the input vector count N is not at least 3, the error */
208 /*         SPICE(INVALIDCOUNT) is signaled. */
209 
210 /*     2)  If any pair of consecutive boundary vectors has cross */
211 /*         product zero, the error SPICE(DEGENERATECASE) is signaled. */
212 /*         For this test, the first vector is considered the successor */
213 /*         of the Nth. */
214 
215 /*     3)  If this routine can't find a face of the convex hull of */
216 /*         the set of boundary vectors such that this face satisfies */
217 /*         constraint (2) of the Detailed_Input section above, the */
218 /*         error SPICE(FACENOTFOUND) is signaled. */
219 
220 /*     4)  If any boundary vectors have longitude too close to 0 */
221 /*         or too close to pi radians in the face frame (see discussion */
222 /*         of the search algorithm's steps 3 and 4 in Particulars */
223 /*         below), the respective errors SPICE(NOTSUPPORTED) or */
224 /*         SPICE(FOVTOOWIDE) are signaled. */
225 
226 /*     5)  If any boundary vectors have angular separation of more than */
227 /*         (pi/2)-MARGIN radians from the candidate FOV axis, the */
228 /*         error SPICE(FOVTOOWIDE) is signaled. */
229 
230 /* $ Files */
231 
232 /*     The boundary vectors input to this routine are typically */
233 /*     obtained from an IK file. */
234 
235 /* $ Particulars */
236 
237 /*     Normally implementation is not discussed in SPICE headers, but we */
238 /*     make an exception here because this routine's implementation and */
239 /*     specification are deeply intertwined. */
240 
241 /*     This routine produces an "axis" for a polygonal FOV using the */
242 /*     following approach: */
243 
244 /*        1)  Test pairs of consecutive FOV boundary vectors to see */
245 /*            whether there's a pair such that the plane region bounded */
246 /*            by these vectors is */
247 
248 /*            a)  part of the convex hull of the set of boundary vectors */
249 
250 /*            b)  such that all other boundary vectors have angular */
251 /*                separation of at least MARGIN from the plane */
252 /*                containing these vectors */
253 
254 /*            This search has O(N**2) run time dependency on N. */
255 
256 /*            If this test produces a candidate face of the convex hull, */
257 /*            proceed to step 3. */
258 
259 
260 /*        2)  If step (1) fails, repeat the search for a candidate */
261 /*            convex hull face, but this time search over every pair of */
262 /*            distinct boundary vectors. */
263 
264 /*            This search has O(N**3) run time dependency on N. */
265 
266 /*            If this search fails, signal an error. */
267 
268 
269 /*        3)  Produce a set of basis vectors for a reference frame, */
270 /*            which we'll call the "face frame," using as the +X axis */
271 /*            the angle bisector of the vectors bounding the candidate */
272 /*            face, the +Y axis the inward normal vector to this face, */
273 /*            and the +Z axis completing a right-handed basis. */
274 
275 
276 /*        4)  Transform each boundary vector, other than the two vectors */
277 /*            defining the selected convex hull face, to the face frame */
278 /*            and compute the vector's longitude in that frame. Find the */
279 /*            maximum and minimum longitudes of the vectors in the face */
280 /*            frame. */
281 
282 /*            If any vector's longitude is less than 2*MARGIN or greater */
283 /*            than pi - 2*MARGIN radians, signal an error. */
284 
285 
286 /*        5)  Let DELTA be the difference between pi and the maximum */
287 /*            longitude found in step (4). Rotate the +Y axis (which */
288 /*            points in the inward normal direction relative to the */
289 /*            selected face) by -DELTA/2 radians about the +Z axis of */
290 /*            the face frame. This rotation aligns the +Y axis with the */
291 /*            central longitude of the set of boundary vectors. The */
292 /*            resulting vector is our candidate FOV axis. */
293 
294 
295 /*        6)  Check the angular separation of the candidate FOV axis */
296 /*            against each boundary vector. If any vector has angular */
297 /*            separation of more than (pi/2)-MARGIN radians from the */
298 /*            axis, signal an error. */
299 
300 
301 /*     Note that there are reasonable FOVs that cannot be handled by the */
302 /*     algorithm described here. For example, any FOV whose cross */
303 /*     section is a regular convex polygon can be made unusable by */
304 /*     adding boundary vectors aligned with the angle bisectors of each */
305 /*     face of the pyramid defined by the FOV's boundary vectors. The */
306 /*     resulting set of boundary vectors has no face in its convex hull */
307 /*     such that all other boundary vectors have positive angular */
308 /*     separation from that face. */
309 
310 /*     Because of this limitation, this algorithm should be used only */
311 /*     after a simple FOV axis-finding approach, such as using as the */
312 /*     FOV axis the average of the boundary vectors, has been tried */
313 /*     unsuccessfully. */
314 
315 /*     Note that it's easy to construct FOVs where the average of the */
316 /*     boundary vectors doesn't yield a viable axis: a FOV of angular */
317 /*     width nearly equal to pi radians, with a sufficiently large */
318 /*     number of boundary vectors on one side and few boundary vectors */
319 /*     on the other, is one such example. This routine can find an */
320 /*     axis for many such intractable FOVs---that's why this routine */
321 /*     should be called after the simple approach fails. */
322 
323 /* $ Examples */
324 
325 /*     See SPICELIB private routine ZZFOVAXI. */
326 
327 /* $ Restrictions */
328 
329 /*     1) This is a SPICE private routine. User applications should not */
330 /*        call this routine. */
331 
332 /*     2) There are "reasonable" polygonal FOVs that cannot be handled */
333 /*        by this routine. See the discussion in Particulars above. */
334 
335 /* $ Literature_References */
336 
337 /*     None. */
338 
339 /* $ Author_and_Institution */
340 
341 /*     N.J. Bachman    (JPL) */
342 
343 /* $ Version */
344 
345 /* -    SPICELIB 1.0.0, 05-MAR-2009 (NJB) */
346 
347 /* -& */
348 /* $ Index_Entries */
349 
350 /*     Create axis vector for polygonal FOV */
351 
352 /* -& */
353 
354 /*     SPICELIB functions */
355 
356 
357 /*     Local parameters */
358 
359 
360 /*     Local variables */
361 
362     /* Parameter adjustments */
363     bounds_dim2 = *n;
364 
365     /* Function Body */
366     if (return_()) {
367 	return 0;
368     }
369     chkin_("ZZHULLAX", (ftnlen)8);
370 
371 /*     Nothing found yet. */
372 
373     found = FALSE_;
374     xidx = 0;
375 
376 /*     We must have at least 3 boundary vectors. */
377 
378     if (*n < 3) {
379 	setmsg_("Polygonal FOV requires at least 3 boundary vectors but numb"
380 		"er supplied for # was #.", (ftnlen)83);
381 	errch_("#", inst, (ftnlen)1, inst_len);
382 	errint_("#", n, (ftnlen)1);
383 	sigerr_("SPICE(INVALIDCOUNT)", (ftnlen)19);
384 	chkout_("ZZHULLAX", (ftnlen)8);
385 	return 0;
386     }
387 
388 /*     Find an exterior face of the pyramid defined by the */
389 /*     input boundary vectors. Since most polygonal FOVs will have */
390 /*     an exterior face bounded by two consecutive rays, we'll */
391 /*     try pairs of consecutive rays first. If this fails, we'll */
392 /*     try each pair of rays. */
393 
394     i__ = 1;
395     while(i__ <= *n && ! found) {
396 
397 /*        Set the index of the next ray. When we get to the */
398 /*        last boundary vector, the next ray is the first. */
399 
400 	if (i__ == *n) {
401 	    next = 1;
402 	} else {
403 	    next = i__ + 1;
404 	}
405 
406 /*        Find the cross product of the first ray with the */
407 /*        second. Depending on the ordering of the boundary */
408 /*        vectors, this could be an inward or outward normal, */
409 /*        in the case the current face is exterior. */
410 
411 	vcrss_(&bounds[(i__1 = i__ * 3 - 3) < bounds_dim2 * 3 && 0 <= i__1 ?
412 		i__1 : s_rnge("bounds", i__1, "zzhullax_", (ftnlen)408)], &
413 		bounds[(i__2 = next * 3 - 3) < bounds_dim2 * 3 && 0 <= i__2 ?
414 		i__2 : s_rnge("bounds", i__2, "zzhullax_", (ftnlen)408)], cp);
415 
416 /*        We insist on consecutive boundary vectors being */
417 /*        linearly independent. */
418 
419 	if (vzero_(cp)) {
420 	    setmsg_("Polygonal FOV must have linearly independent consecutiv"
421 		    "e boundary but vectors at indices # and # have cross pro"
422 		    "duct equal to the zero vector. Instrument is #.", (ftnlen)
423 		    158);
424 	    errint_("#", &i__, (ftnlen)1);
425 	    errint_("#", &next, (ftnlen)1);
426 	    errch_("#", inst, (ftnlen)1, inst_len);
427 	    sigerr_("SPICE(DEGENERATECASE)", (ftnlen)21);
428 	    chkout_("ZZHULLAX", (ftnlen)8);
429 	    return 0;
430 	}
431 
432 /*        See whether the other boundary vectors have angular */
433 /*        separation of at least MARGIN from the plane containing */
434 /*        the current face. */
435 
436 	pass1 = TRUE_;
437 	ok = TRUE_;
438 	m = 1;
439 	while(m <= *n && ok) {
440 
441 /*           Find the angular separation of CP and the Mth vector if the */
442 /*           latter is not an edge of the current face. */
443 
444 	    if (m != i__ && m != next) {
445 		sep = vsep_(cp, &bounds[(i__1 = m * 3 - 3) < bounds_dim2 * 3
446 			&& 0 <= i__1 ? i__1 : s_rnge("bounds", i__1, "zzhull"
447 			"ax_", (ftnlen)446)]);
448 		if (pass1) {
449 
450 /*                 Adjust CP if necessary so that it points */
451 /*                 toward the interior of the pyramid. */
452 
453 		    if (sep > halfpi_()) {
454 
455 /*                    Invert the cross product vector and adjust SEP */
456 /*                    accordingly. Within this "M" loop, all other */
457 /*                    angular separations will be computed using the new */
458 /*                    value of CP. */
459 
460 			vsclip_(&c_b20, cp);
461 			sep = pi_() - sep;
462 		    }
463 		    pass1 = FALSE_;
464 		}
465 		ok = sep < halfpi_() - 1e-12;
466 	    }
467 	    if (ok) {
468 
469 /*              Consider the next boundary vector. */
470 
471 		++m;
472 	    }
473 	}
474 
475 /*        We've tested each boundary vector against the current face, or */
476 /*        else the loop terminated early because a vector with */
477 /*        insufficient angular separation from the plane containing the */
478 /*        face was found. */
479 
480 	if (ok) {
481 
482 /*           The current face is exterior. It's bounded by rays I and */
483 /*           NEXT. */
484 
485 	    xidx = i__;
486 	    found = TRUE_;
487 	} else {
488 
489 /*           Look at the next face of the pyramid. */
490 
491 	    ++i__;
492 	}
493     }
494 
495 /*     If we didn't find an exterior face, we'll have to look at each */
496 /*     face bounded by a pair of rays, even if those rays are not */
497 /*     adjacent. (This can be a very slow process is N is large.) */
498 
499     if (! found) {
500 	i__ = 1;
501 	while(i__ <= *n && ! found) {
502 
503 /*           Consider all ray pairs (I,NEXT) where NEXT > I. */
504 
505 	    next = i__ + 1;
506 	    while(next <= *n && ! found) {
507 
508 /*              Find the cross product of the first ray with the second. */
509 /*              If the current face is exterior, CP could be an inward */
510 /*              or outward normal, depending on the ordering of the */
511 /*              boundary vectors. */
512 
513 		vcrss_(&bounds[(i__1 = i__ * 3 - 3) < bounds_dim2 * 3 && 0 <=
514 			i__1 ? i__1 : s_rnge("bounds", i__1, "zzhullax_", (
515 			ftnlen)530)], &bounds[(i__2 = next * 3 - 3) <
516 			bounds_dim2 * 3 && 0 <= i__2 ? i__2 : s_rnge("bounds",
517 			 i__2, "zzhullax_", (ftnlen)530)], cp);
518 
519 /*              It's allowable for non-consecutive boundary vectors to */
520 /*              be linearly dependent, but if we have such a pair, */
521 /*              it doesn't define an exterior face. */
522 
523 		if (! vzero_(cp)) {
524 
525 /*                 The rays having direction vectors indexed I and NEXT */
526 /*                 define a semi-infinite sector of a plane that might */
527 /*                 be of interest. */
528 
529 /*                 Check whether all of the boundary vectors that are */
530 /*                 not edges of the current face have angular separation */
531 /*                 of at least MARGIN from the plane containing the */
532 /*                 current face. */
533 
534 		    pass1 = TRUE_;
535 		    ok = TRUE_;
536 		    m = 1;
537 		    while(m <= *n && ok) {
538 
539 /*                    Find the angular separation of CP and the Mth */
540 /*                    vector if the latter is not an edge of the current */
541 /*                    face. */
542 
543 			if (m != i__ && m != next) {
544 			    sep = vsep_(cp, &bounds[(i__1 = m * 3 - 3) <
545 				    bounds_dim2 * 3 && 0 <= i__1 ? i__1 :
546 				    s_rnge("bounds", i__1, "zzhullax_", (
547 				    ftnlen)560)]);
548 			    if (pass1) {
549 
550 /*                          Adjust CP if necessary so that it points */
551 /*                          toward the interior of the pyramid. */
552 
553 				if (sep > halfpi_()) {
554 
555 /*                             Invert the cross product vector and */
556 /*                             adjust SEP accordingly. Within this "M" */
557 /*                             loop, all other angular separations will */
558 /*                             be computed using the new value of CP. */
559 
560 				    vsclip_(&c_b20, cp);
561 				    sep = pi_() - sep;
562 				}
563 				pass1 = FALSE_;
564 			    }
565 			    ok = sep < halfpi_() - 1e-12;
566 			}
567 			if (ok) {
568 
569 /*                       Consider the next boundary vector. */
570 
571 			    ++m;
572 			}
573 		    }
574 
575 /*                 We've tested each boundary vector against the current */
576 /*                 face, or else the loop terminated early because a */
577 /*                 vector with insufficient angular separation from the */
578 /*                 plane containing the face was found. */
579 
580 		    if (ok) {
581 
582 /*                    The current face is exterior. It's bounded by rays */
583 /*                    I and NEXT. */
584 			xidx = i__;
585 			found = TRUE_;
586 		    }
587 
588 /*                 End of angular separation test block. */
589 
590 		}
591 
592 /*              End of non-zero cross product block. */
593 
594 		if (! found) {
595 
596 /*                 Look at the face bounded by the rays */
597 /*                 at indices I and NEXT+1. */
598 
599 		    ++next;
600 		}
601 	    }
602 
603 /*           End of NEXT loop. */
604 
605 	    if (! found) {
606 
607 /*              Look at the face bounded by the pairs of rays */
608 /*              including the ray at index I+1. */
609 
610 		++i__;
611 	    }
612 	}
613 
614 /*        End of I loop. */
615 
616     }
617 
618 /*     End of search for exterior face using each pair of rays. */
619 
620 /*     If we still haven't found an exterior face, we can't continue. */
621 
622     if (! found) {
623 	setmsg_("Unable to find face of convex hull of FOV of instrument #.",
624 		(ftnlen)58);
625 	errch_("#", inst, (ftnlen)1, inst_len);
626 	sigerr_("SPICE(FACENOTFOUND)", (ftnlen)19);
627 	chkout_("ZZHULLAX", (ftnlen)8);
628 	return 0;
629     }
630 
631 /*     Arrival at this point means that the rays at indices */
632 /*     XIDX and NEXT define a plane such that all boundary */
633 /*     vectors lie in a half-space bounded by that plane. */
634 
635 /*     We're now going to define a set of orthonormal basis vectors: */
636 
637 /*        +X  points along the angle bisector of the bounding vectors */
638 /*            of the exterior face. */
639 
640 /*        +Y  points along CP. */
641 
642 /*        +Z  is the cross product of +X and +Y. */
643 
644 /*     We'll call the reference frame having these basis vectors */
645 /*     the "face frame." */
646 
647 
648     vhat_(&bounds[(i__1 = i__ * 3 - 3) < bounds_dim2 * 3 && 0 <= i__1 ? i__1 :
649 	     s_rnge("bounds", i__1, "zzhullax_", (ftnlen)683)], ray1);
650     vhat_(&bounds[(i__1 = next * 3 - 3) < bounds_dim2 * 3 && 0 <= i__1 ? i__1
651 	    : s_rnge("bounds", i__1, "zzhullax_", (ftnlen)684)], ray2);
652     vlcom_(&c_b36, ray1, &c_b36, ray2, xvec);
653     vhatip_(xvec);
654     vhat_(cp, yvec);
655     ucrss_(xvec, yvec, zvec);
656 
657 /*     Create a transformation matrix to map the input boundary */
658 /*     vectors into the face frame. */
659 
660     for (i__ = 1; i__ <= 3; ++i__) {
661 	trans[(i__1 = i__ * 3 - 3) < 9 && 0 <= i__1 ? i__1 : s_rnge("trans",
662 		i__1, "zzhullax_", (ftnlen)698)] = xvec[(i__2 = i__ - 1) < 3
663 		&& 0 <= i__2 ? i__2 : s_rnge("xvec", i__2, "zzhullax_", (
664 		ftnlen)698)];
665 	trans[(i__1 = i__ * 3 - 2) < 9 && 0 <= i__1 ? i__1 : s_rnge("trans",
666 		i__1, "zzhullax_", (ftnlen)699)] = yvec[(i__2 = i__ - 1) < 3
667 		&& 0 <= i__2 ? i__2 : s_rnge("yvec", i__2, "zzhullax_", (
668 		ftnlen)699)];
669 	trans[(i__1 = i__ * 3 - 1) < 9 && 0 <= i__1 ? i__1 : s_rnge("trans",
670 		i__1, "zzhullax_", (ftnlen)700)] = zvec[(i__2 = i__ - 1) < 3
671 		&& 0 <= i__2 ? i__2 : s_rnge("zvec", i__2, "zzhullax_", (
672 		ftnlen)700)];
673     }
674 
675 /*     Now we're going to compute the longitude of each boundary in the */
676 /*     face frame. The vectors with indices XIDX and NEXT are excluded. */
677 /*     We expect all longitudes to be between MARGIN and pi - MARGIN. */
678 
679     minlon = pi_();
680     maxlon = 0.;
681     minix = 1;
682     maxix = 1;
683     i__1 = *n;
684     for (i__ = 1; i__ <= i__1; ++i__) {
685 	if (i__ != xidx && i__ != next) {
686 
687 /*           The current vector is not a boundary of our edge, */
688 /*           so find its longitude. */
689 
690 	    mxv_(trans, &bounds[(i__2 = i__ * 3 - 3) < bounds_dim2 * 3 && 0 <=
691 		     i__2 ? i__2 : s_rnge("bounds", i__2, "zzhullax_", (
692 		    ftnlen)720)], v);
693 	    reclat_(v, &r__, &lon, &lat);
694 
695 /*           Update the longitude bounds. */
696 
697 	    if (lon < minlon) {
698 		minix = i__;
699 		minlon = lon;
700 	    }
701 	    if (lon > maxlon) {
702 		maxix = i__;
703 		maxlon = lon;
704 	    }
705 	}
706     }
707 
708 /*     If the longitude bounds are not as expected, don't try */
709 /*     to continue. */
710 
711     if (minlon < 2e-12) {
712 	setmsg_("Minimum boundary vector longitude in exterior face frame is"
713 		" # radians. Minimum occurs at index #. This FOV does not con"
714 		"form to the requirements of this routine. Instrument is #.", (
715 		ftnlen)177);
716 	errdp_("#", &minlon, (ftnlen)1);
717 	errint_("#", &minix, (ftnlen)1);
718 	errch_("#", inst, (ftnlen)1, inst_len);
719 	sigerr_("SPICE(NOTSUPPORTED)", (ftnlen)19);
720 	chkout_("ZZHULLAX", (ftnlen)8);
721 	return 0;
722     } else if (maxlon > pi_() - 2e-12) {
723 	setmsg_("Maximum boundary vector longitude in exterior face frame is"
724 		" # radians. Maximum occurs at index #. This FOV does not con"
725 		"form to the requirements of this routine. Instrument is #.", (
726 		ftnlen)177);
727 	errdp_("#", &maxlon, (ftnlen)1);
728 	errint_("#", &maxix, (ftnlen)1);
729 	errch_("#", inst, (ftnlen)1, inst_len);
730 	sigerr_("SPICE(FOVTOOWIDE)", (ftnlen)17);
731 	chkout_("ZZHULLAX", (ftnlen)8);
732 	return 0;
733     }
734 
735 /*     Let delta represent the amount we can rotate the exterior */
736 /*     face clockwise about +Z without contacting another boundary */
737 /*     vector. */
738 
739     delta = pi_() - maxlon;
740 
741 /*     Rotate +Y by -DELTA/2 about +Z. The result is our candidate */
742 /*     FOV axis. Make the axis vector unit length. */
743 
744     d__1 = -delta / 2;
745     vrotv_(yvec, zvec, &d__1, axis);
746     vhatip_(axis);
747 
748 /*     If we have a viable result, ALL boundary vectors have */
749 /*     angular separation less than HALFPI-MARGIN from AXIS. */
750 
751     i__1 = *n;
752     for (i__ = 1; i__ <= i__1; ++i__) {
753 	sep = vsep_(&bounds[(i__2 = i__ * 3 - 3) < bounds_dim2 * 3 && 0 <=
754 		i__2 ? i__2 : s_rnge("bounds", i__2, "zzhullax_", (ftnlen)794)
755 		], axis);
756 	if (sep > halfpi_() - 1e-12) {
757 	    setmsg_("Boundary vector at index # has angular separation of # "
758 		    "radians from candidate FOV axis. This FOV does not confo"
759 		    "rm to the requirements of this routine. Instrument is #.",
760 		     (ftnlen)167);
761 	    errint_("#", &i__, (ftnlen)1);
762 	    errdp_("#", &sep, (ftnlen)1);
763 	    errch_("#", inst, (ftnlen)1, inst_len);
764 	    sigerr_("SPICE(FOVTOOWIDE)", (ftnlen)17);
765 	    chkout_("ZZHULLAX", (ftnlen)8);
766 	    return 0;
767 	}
768     }
769     chkout_("ZZHULLAX", (ftnlen)8);
770     return 0;
771 } /* zzhullax_ */
772 
773