spice/core/
raw.rs

1/*!
2A Rust idiomatic CSPICE wrapper built with [procedural macros][`spice_derive`].
3*/
4
5use crate::{c, cstr, fcstr, malloc, mallocstr};
6use spice_derive::{cspice_proc, return_output};
7use std::ops::{Deref, DerefMut};
8
9#[cfg(any(feature = "lock", doc))]
10use {crate::core::lock::SpiceLock, spice_derive::impl_for};
11
12#[allow(clippy::upper_case_acronyms)]
13pub type DLADSC = c::SpiceDLADescr;
14#[allow(clippy::upper_case_acronyms)]
15pub type DSKDSC = c::SpiceDSKDescr;
16#[allow(clippy::upper_case_acronyms)]
17pub type CELL = c::SpiceCell;
18pub const CELL_MAXID: usize = 10_000;
19
20/**
21A cell is a data structure intended to provide safe array access within the applications.
22
23See the [C documentation](https://naif.jpl.nasa.gov/pub/naif/toolkit_docs/C/req/cells.html).
24*/
25#[derive(Debug)]
26pub struct Cell(c::SpiceCell);
27
28impl Cell {
29    /**
30    Declare a cell from integer.
31    */
32    pub fn new_int() -> Self {
33        let base = malloc!(i32, CELL_MAXID + c::SPICE_CELL_CTRLSZ as usize);
34        Self(CELL {
35            dtype: c::_SpiceDataType_SPICE_INT,
36            length: 0i32,
37            size: CELL_MAXID as i32,
38            card: 0i32,
39            isSet: 1i32,
40            adjust: 0i32,
41            init: 0i32,
42            base: base as *mut libc::c_void,
43            data: base.wrapping_add(c::SPICE_CELL_CTRLSZ as usize) as *mut libc::c_void,
44        })
45    }
46
47    /**
48    Declare data from a cell at index.
49    */
50    pub fn get_data_int(&self, index: usize) -> i32 {
51        unsafe { *(self.data as *mut i32).wrapping_add(index) }
52    }
53}
54
55impl Deref for Cell {
56    type Target = c::SpiceCell;
57
58    fn deref(&self) -> &Self::Target {
59        &self.0
60    }
61}
62
63impl DerefMut for Cell {
64    fn deref_mut(&mut self) -> &mut Self::Target {
65        &mut self.0
66    }
67}
68
69cspice_proc! {
70    /**
71    Translate the SPICE integer code of a body into a common name for that body.
72
73    This function has a [neat version][crate::neat::bodc2n].
74    */
75    pub fn bodc2n(code: i32, lenout: i32) -> (String, bool) {}
76}
77
78cspice_proc! {
79    /**
80    Determine whether values exist for some item for any body in the kernel pool.
81    */
82    #[return_output]
83    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
84    pub fn bodfnd(body: i32, item: &str) -> bool {}
85}
86
87cspice_proc! {
88    /**
89    Translate the name of a body or object to the corresponding SPICE integer ID code.
90    */
91    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
92    pub fn bodn2c(name: &str) -> (i32, bool) {}
93}
94
95/**
96Fetch from the kernel pool the double precision values of an item associated with a body.
97*/
98pub fn bodvrd(bodynm: &str, item: &str, maxn: usize) -> Vec<f64> {
99    let bodynm = cstr!(bodynm);
100    let item = cstr!(item);
101    let mut dim = 0;
102    let mut values = vec![0.0; maxn];
103    unsafe { crate::c::bodvrd_c(bodynm, item, maxn as _, &mut dim, values.as_mut_ptr()) };
104    values.truncate(dim as _);
105    values
106}
107
108cspice_proc! {
109    /**
110    close a das file.
111    */
112    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
113    pub fn dascls(handle: i32) {}
114}
115
116cspice_proc! {
117    /**
118    Open a DAS file for reading.
119    */
120    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
121    pub fn dasopr(fname: &str) -> i32 {}
122}
123
124/**
125Return the value of Delta ET (ET-UTC) for an input epoch.
126*/
127pub fn deltet(epoch: f64, eptype: &str) -> f64 {
128    let eptype = cstr!(eptype);
129    let mut delta = 0.0;
130    unsafe {
131        crate::c::deltet_c(epoch, eptype, &mut delta);
132    }
133    delta
134}
135
136cspice_proc! {
137    /**
138    Begin a forward segment search in a DLA file.
139    */
140    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
141    pub fn dlabfs(handle: i32) -> (DLADSC, bool) {}
142}
143
144cspice_proc! {
145    /**
146    Return the DSK descriptor from a DSK segment identified  by a DAS handle and DLA descriptor.
147    */
148    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
149    pub fn dskgd(handle: i32, dladsc: DLADSC) -> DSKDSC {}
150}
151
152cspice_proc! {
153    /**
154    Compute the unit normal vector for a specified plate from a type 2 DSK segment.
155    */
156    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
157    pub fn dskn02(handle: i32, dladsc: DLADSC, plid: i32) -> [f64; 3] {}
158}
159
160cspice_proc! {
161    /**
162    Find the set of body ID codes of all objects for which topographic data are provided in a
163    specified DSK file.
164    */
165    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
166    pub fn dskobj(dsk: &str) -> Cell {}
167}
168
169/**
170Fetch triangular plates from a type 2 DSK segment.
171
172This function has a [neat version][crate::neat::dskp02].
173*/
174pub fn dskp02(handle: i32, mut dladsc: DLADSC, start: usize, room: usize) -> Vec<[i32; 3]> {
175    let mut n = 0;
176    let mut plates = vec![[0; 3]; room];
177
178    unsafe {
179        crate::c::dskp02_c(
180            handle,
181            &mut dladsc,
182            start as _,
183            room as _,
184            &mut n,
185            plates.as_mut_ptr(),
186        );
187    }
188
189    plates.truncate(n as _);
190    plates
191}
192
193/**
194Fetch vertices from a type 2 DSK segment.
195
196This function has a [neat version][crate::neat::dskv02].
197*/
198pub fn dskv02(handle: i32, mut dladsc: DLADSC, start: usize, room: usize) -> Vec<[f64; 3]> {
199    let mut n = 0;
200    let mut vrtces = vec![[0.0; 3]; room];
201
202    unsafe {
203        crate::c::dskv02_c(
204            handle,
205            &mut dladsc,
206            start as _,
207            room as _,
208            &mut n,
209            vrtces.as_mut_ptr(),
210        );
211    }
212
213    vrtces.truncate(n as _);
214    vrtces
215}
216
217cspice_proc! {
218    /**
219    Determine the plate ID and body-fixed coordinates of the intersection of a specified ray with
220    the surface defined by a type 2 DSK plate model.
221    */
222    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
223    pub fn dskx02(
224        handle: i32,
225        dladsc: DLADSC,
226        vertex: [f64; 3],
227        raydir: [f64; 3]
228    ) -> (i32, [f64; 3], bool) {
229    }
230}
231
232cspice_proc! {
233    /**
234    Return plate model size parameters---plate count and
235    vertex count---for a type 2 DSK segment.
236
237    Plate first, vertices second.
238
239
240    See [`neat::dskp02`] for the raw interface.
241    */
242    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
243    pub fn dskz02(handle: i32, dladsc: DLADSC) -> (i32, i32) {}
244}
245
246/**
247Return the d.p. value of a kernel variable from the kernel pool.
248*/
249pub fn gdpool(name: &str, start: usize, room: usize) -> Vec<f64> {
250    let name = cstr!(name);
251    let start = start as _;
252    let mut n = 0;
253    let mut values = vec![0.0; room];
254    let mut found = 0;
255    unsafe {
256        crate::c::gdpool_c(
257            name,
258            start,
259            room as _,
260            &mut n,
261            values.as_mut_ptr(),
262            &mut found,
263        )
264    }
265    // let found = found != 0;
266    values.truncate(n as _);
267    values
268}
269
270cspice_proc! {
271    /**
272    Convert geodetic coordinates to rectangular coordinates.
273     */
274    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
275    pub fn georec(lon: f64, lat: f64, alt: f64, re: f64, f: f64) -> [f64; 3] {}
276}
277
278/**
279Return the field-of-view (FOV) parameters for a specified
280instrument. The instrument is specified by its NAIF ID code.
281*/
282pub fn getfov(
283    instid: isize,
284    room: usize,
285    shapelen: usize,
286    framelen: usize,
287) -> (String, String, [f64; 3], Vec<[f64; 3]>) {
288    let shape = mallocstr!(shapelen);
289    let frame = mallocstr!(framelen);
290
291    let mut bsight = [0.0; 3];
292    let mut n = 0;
293    let mut bounds = vec![[0.0; 3]; room];
294
295    unsafe {
296        crate::c::getfov_c(
297            instid as _,
298            room as _,
299            shapelen as _,
300            framelen as _,
301            shape,
302            frame,
303            bsight.as_mut_ptr(),
304            &mut n,
305            bounds.as_mut_ptr(),
306        )
307    };
308
309    bounds.truncate(n as _);
310    (fcstr!(shape), fcstr!(frame), bsight, bounds)
311}
312
313cspice_proc! {
314    /**
315    Compute the illumination angles---phase, incidence, and emission---at a specified point on a
316    target body. Return logical flags indicating whether the surface point is visible from the
317    observer's position and whether the surface point is illuminated.
318
319    The target body's surface is represented using topographic data provided by DSK files, or by a
320    reference ellipsoid.
321
322    The illumination source is a specified ephemeris object.
323    */
324    #[allow(clippy::too_many_arguments)]
325    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
326    pub fn illumf(
327        method: &str,
328        target: &str,
329        ilusrc: &str,
330        et: f64,
331        fixref: &str,
332        abcorr: &str,
333        obsrvr: &str,
334        spoint: [f64; 3]
335    ) -> (f64, [f64; 3], f64, f64, f64, bool, bool) {}
336}
337
338cspice_proc! {
339    /**
340    Load one or more SPICE kernels into a program.
341    */
342    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
343    pub fn furnsh(name: &str) {}
344}
345
346cspice_proc! {
347    /**
348    Clear the KEEPER subsystem: unload all kernels, clear the kernel pool, and re-initialize the
349    subsystem. Existing watches on kernel variables are retained.
350    */
351    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
352    pub fn kclear() {}
353}
354
355/**
356Fetch vertices from a type 2 DSK segment.
357
358This function has a [neat version][crate::neat::kdata].
359*/
360pub fn kdata(
361    which: i32,
362    kind: &str,
363    fillen: i32,
364    typlen: i32,
365    srclen: i32,
366) -> (String, String, String, i32, bool) {
367    #[allow(unused_unsafe)]
368    unsafe {
369        let varout_0 = mallocstr!(fillen);
370        let varout_1 = mallocstr!(typlen);
371        let varout_2 = mallocstr!(srclen);
372        let mut varout_3 = 0i32;
373        let mut varout_4 = 0i32;
374        crate::c::kdata_c(
375            which,
376            cstr!(kind),
377            fillen,
378            typlen,
379            srclen,
380            varout_0,
381            varout_1,
382            varout_2,
383            &mut varout_3,
384            &mut varout_4,
385        );
386        (
387            fcstr!(varout_0),
388            fcstr!(varout_1),
389            fcstr!(varout_2),
390            varout_3,
391            varout_4 != 0,
392        )
393    }
394}
395
396cspice_proc! {
397    /**
398    Return the current number of kernels that have been loaded via the KEEPER interface that are of
399    a specified type.
400    */
401    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
402    pub fn ktotal(kind: &str) -> i32 {}
403}
404
405cspice_proc! {
406    /**
407    Convert from latitudinal coordinates to rectangular coordinates.
408    */
409    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
410    pub fn latrec(radius: f64, longitude: f64, latitude: f64) -> [f64; 3] {}
411}
412
413cspice_proc! {
414    /**
415       Multiply a 3x3 double precision matrix with a 3-dimensional double precision vector.
416    */
417    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
418    pub fn mxv(m1: [[f64; 3]; 3], vin: [f64; 3]) -> [f64; 3] {}
419}
420
421cspice_proc! {
422    /**
423    Determines the occultation condition (not occulted, partially, etc.) of one target relative to
424    another target as seen by an observer at a given time, with targets modeled as points,
425    ellipsoids, or digital shapes (DSK)
426    */
427    #[allow(clippy::too_many_arguments)]
428    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
429    pub fn occult(
430        targ1: &str,
431        shape1: &str,
432        frame1: &str,
433        targ2: &str,
434        shape2: &str,
435        frame2: &str,
436        abcorr: &str,
437        obsrvr: &str,
438        et: f64,
439    ) -> i32 {}
440}
441
442cspice_proc! {
443    /**
444    Return the matrix that transforms position vectors from one specified frame to another at a
445    specified epoch.
446    */
447    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
448    pub fn pxform(from: &str, to: &str, et: f64) -> [[f64; 3]; 3] {}
449}
450
451cspice_proc! {
452    /**
453    Return the 3x3 matrix that transforms position vectors from one specified frame at a specified
454    epoch to another specified frame at another specified epoch.
455    */
456    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
457    pub fn pxfrm2(from: &str, to: &str, etfrom: f64, etto: f64) -> [[f64; 3]; 3] {}
458}
459
460cspice_proc! {
461    /**
462    Convert range, right ascension, and declination to rectangular coordinates
463     */
464    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
465    pub fn radrec(range: f64, ra: f64, dec: f64) -> [f64; 3] {}
466}
467
468/**
469Convert rectangular coordinates to planetographic coordinates.
470*/
471pub fn recpgr(body: &str, rectan: [f64; 3], re: f64, f: f64) -> [f64; 3] {
472    let body = cstr!(body);
473    let mut rectan: [f64; 3] = rectan;
474    let mut lon = 0.0;
475    let mut lat = 0.0;
476    let mut alt = 0.0;
477    unsafe { crate::c::recpgr_c(body, &mut rectan as _, re, f, &mut lon, &mut lat, &mut alt) };
478    [lon, lat, alt]
479}
480
481cspice_proc! {
482    /**
483    Convert rectangular coordinates to range, right ascension, and declination.
484    */
485    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
486    pub fn recrad(rectan: [f64; 3]) -> (f64, f64, f64) {}
487}
488
489cspice_proc! {
490    /**
491    Compute, for a given observer and a ray emanating from the
492    observer, the surface intercept of the ray on a target body at
493    a specified epoch, optionally corrected for light time and
494    stellar aberration.
495
496    The surface of the target body may be represented by a triaxial
497    ellipsoid or by topographic data provided by DSK files.
498
499    This routine supersedes srfxpt.
500    */
501    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
502    pub fn sincpt(
503        method:&str,
504        target: &str,
505        et: f64,
506        fixred: &str,
507        abcorr: &str,
508        obsrvr: &str,
509        dref: &str,
510        dvec: [f64; 3]) -> ([f64; 3], f64, [f64; 3], bool) {}
511}
512
513cspice_proc! {
514    /**
515    Close a SPK file opened for read or write.
516    */
517    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
518    pub fn spkcls(handle: i32) {}
519}
520
521cspice_proc! {
522    /**
523    Create a new SPK file, returning the handle of the opened file
524     */
525    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
526    pub fn spkopn(fname: &str, ifname: &str, ncomch: i32) -> i32 {}
527}
528
529cspice_proc! {
530    /**
531    Write a type 9 segment to an SPK file.
532    */
533    #[allow(clippy::too_many_arguments)]
534    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
535    pub fn spkw09(handle: i32, body: i32, center: i32, frame: &str, first: f64, last: f64, segid: &str, degree: i32, n: i32, states: &mut [[f64; 6]], epochs: &mut [f64]) {}
536}
537
538cspice_proc! {
539    /**
540    Return the position of a target body relative to an observing body, optionally corrected for
541    light time (planetary aberration) and stellar aberration.
542    */
543    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
544    pub fn spkpos(targ: &str, et: f64, frame: &str, abcorr: &str, obs: &str) -> ([f64; 3], f64) {}
545}
546
547cspice_proc! {
548    /**
549    Return the state (position and velocity) of a target body
550    relative to an observing body, optionally corrected for light
551    time (planetary aberration) and stellar aberration.
552    */
553    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
554    pub fn spkezr(targ: &str, et: f64, frame: &str, abcorr: &str, obs: &str) -> ([f64; 6], f64) {}
555}
556
557cspice_proc! {
558    /**
559    Convert a string representing an epoch to a double precision value representing the number of
560    TDB seconds past the J2000 epoch corresponding to the input epoch.
561    */
562    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
563    pub fn str2et(targ: &str) -> f64 {}
564}
565
566/**
567Compute the rectangular coordinates of the sub-observer point on
568a target body at a specified epoch, optionally corrected for
569light time and stellar aberration.
570
571The surface of the target body may be represented by a triaxial
572ellipsoid or by topographic data provided by DSK files.
573*/
574pub fn subpnt(
575    method: &str,
576    target: &str,
577    et: f64,
578    fixref: &str,
579    abcorr: &str,
580    obsrvr: &str,
581) -> ([f64; 3], f64, [f64; 3]) {
582    let method = cstr!(method);
583    let target = cstr!(target);
584    let fixref = cstr!(fixref);
585    let abcorr = cstr!(abcorr);
586    let obsrvr = cstr!(obsrvr);
587    let mut sp = [0.0; 3];
588    let mut et_sp = 0.0;
589    let mut vec_sp = [0.0; 3];
590    unsafe {
591        crate::c::subpnt_c(
592            method,
593            target,
594            et,
595            fixref,
596            abcorr,
597            obsrvr,
598            &mut sp as _,
599            &mut et_sp,
600            &mut vec_sp as _,
601        )
602    };
603    (sp, et_sp, vec_sp)
604}
605
606cspice_proc! {
607    /**
608    Determine the intersection of a line-of-sight vector with the surface of an ellipsoid.
609    */
610    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
611    pub fn surfpt(positn: [f64; 3], u: [f64; 3], a: f64, b: f64, c: f64) -> ([f64; 3], bool) {}
612}
613
614/**
615Convert an input epoch represented in TDB seconds past the TDB epoch of J2000 to a character string formatted to the
616specifications of a user's format picture.
617
618This function has a [neat version][crate::neat::timout].
619*/
620pub fn timout(et: f64, pictur: &str, lenout: usize) -> String {
621    let varout_0 = mallocstr!(lenout);
622    unsafe {
623        crate::c::timout_c(et, cstr!(pictur), lenout as i32, varout_0);
624    }
625    fcstr!(varout_0)
626}
627
628/**
629Transform time from one uniform scale to another. The uniform time scales are
630TAI, GPS, TT, TDT, TDB, ET, JED, JDTDB, JDTDT.
631*/
632pub fn unitim(epoch: f64, insys: &str, outsys: &str) -> f64 {
633    let insys = cstr!(insys);
634    let outsys = cstr!(outsys);
635    unsafe { crate::c::unitim_c(epoch, insys, outsys) }
636}
637
638cspice_proc! {
639    /**
640    Unload a SPICE kernel.
641    */
642    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
643    pub fn unload(name: &str) {}
644}
645
646cspice_proc! {
647    /**
648    Find the separation angle in radians between two double precision, 3-dimensional vectors. This
649    angle is defined as zero if either vector is zero.
650    */
651    #[return_output]
652    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
653    pub fn vsep(v1: [f64; 3], v2: [f64; 3]) -> f64 {}
654}
655
656cspice_proc! {
657    /**
658    Compute the dot product of two double precision, 3-dimensional vectors.
659     */
660    #[return_output]
661    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
662    pub fn vdot(v1: [f64; 3], v2: [f64; 3]) -> f64 {}
663}
664
665cspice_proc! {
666    /**
667    Compute the cross product of two 3-dimensional vectors.
668     */
669    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
670    pub fn vcrss(v1: [f64; 3], v2: [f64; 3]) -> [f64; 3] {}
671}
672
673cspice_proc! {
674    /**
675    Transpose a 3x3 matrix.
676     */
677    #[cfg_attr(any(feature = "lock", doc), impl_for(SpiceLock))]
678    pub fn xpose(m1: [[f64; 3]; 3]) -> [[f64; 3]; 3] {}
679}