Mark5BHeader¶
-
class
baseband.mark5b.
Mark5BHeader
(words, kday=None, ref_time=None, verify=True)[source] [edit on github]¶ Bases:
baseband.base.header.VLBIHeaderBase
Decoder/encoder of a Mark5B Frame Header.
See page 15 of https://www.haystack.mit.edu/tech/vlbi/mark5/docs/Mark%205B%20users%20manual.pdf
- Parameters
- wordstuple of int, or None
Four 32-bit unsigned int header words. If
None
, set to a tuple of zeros for later initialisation.- kdayint or None
Explicit thousands of MJD of the observation time (needed to remove ambiguity in the Mark 5B time stamp). Can instead pass an approximate
ref_time
.- ref_time
Time
or None Reference time within 500 days of the observation time, used to infer the full MJD. Used only if
kday
is not given.- verifybool, optional
Whether to do basic verification of integrity. Default:
True
.
- Returns
- header
Mark5BHeader
- header
Attributes Summary
Whether the data are complex (always False for Mark5B).
Fractional seconds (decoded from ‘bcd_fraction’).
Size of the frame in bytes (always 10016 for Mark5B).
Last three digits of MJD (decoded from ‘bcd_jday’).
Thousands of MJD, to complement
jday
from header.Whether the header can be modified.
Size of the header in bytes.
Size of the payload in bytes (always 10000 for Mark5B).
Integer seconds on day (decoded from ‘bcd_seconds’).
Convert year, BCD time code to Time object.
Methods Summary
copy
(**kwargs)Create a mutable and independent copy of the header.
fromfile
(fh, *args, **kwargs)Read VLBI Header from file.
fromkeys
(*args, **kwargs)Initialise a header from parsed values.
fromvalues
(*[, verify])Initialise a header from parsed values.
get_time
([frame_rate])Convert year, BCD time code to Time object.
infer_kday
(ref_time)Uses a reference time to set a header’s
kday
.invariant_pattern
([invariants])Pattern and mask shared between headers of a type or stream.
Set of keys of invariant header parts.
keys
()All keys defined for this header.
set_time
(time[, frame_rate])Convert Time object to BCD timestamp elements and ‘frame_nr’.
tofile
(fh)Write VLBI frame header to filehandle.
update
(*[, time, frame_rate, crc, verify])Update the header by setting keywords or properties.
verify
()Verify header integrity.
Attributes Documentation
-
complex_data
¶ Whether the data are complex (always False for Mark5B).
-
fraction
¶ Fractional seconds (decoded from ‘bcd_fraction’).
The fraction is stored to 0.1 ms accuracy. Following mark5access, this is “unrounded” to give the exact time of the start of the frame for any total bit rate below 512 Mbps. For rates above this value, it is no longer guaranteed that subsequent frames have unique rates.
Note to the above: since a Mark5B frame contains 80000 bits, the total bit rate for which times can be unique would in principle be 800 Mbps. However, standard VLBI only uses bit rates that are powers of 2 in MHz.
-
frame_nbytes
¶ Size of the frame in bytes (always 10016 for Mark5B).
-
jday
¶ Last three digits of MJD (decoded from ‘bcd_jday’).
-
kday
= None¶ Thousands of MJD, to complement
jday
from header.
-
mutable
¶ Whether the header can be modified.
-
nbytes
¶ Size of the header in bytes.
-
payload_nbytes
¶ Size of the payload in bytes (always 10000 for Mark5B).
-
seconds
¶ Integer seconds on day (decoded from ‘bcd_seconds’).
-
time
¶ Convert year, BCD time code to Time object.
Calculate time using
jday
,seconds
, andfraction
properties (which reflect the bcd-encoded ‘bcd_jday’, ‘bcd_seconds’ and ‘bcd_fraction’ header items), pluskday
from the initialisation. See https://www.haystack.mit.edu/tech/vlbi/mark5/docs/Mark%205B%20users%20manual.pdfNote that some non-compliant files do not have ‘bcd_fraction’ set. For those, the time can still be calculated using the header’s ‘frame_nr’ by passing in a frame rate.
Furthermore, fractional seconds are stored only to 0.1 ms accuracy. In the code, this is “unrounded” to give the exact time of the start of the frame for any total bit rate below 512 Mbps. For higher rates, it is no longer guaranteed that subsequent frames have unique
fraction
, and one should pass in an explicit frame rate instead.
Methods Documentation
-
copy
(**kwargs)[source] [edit on github]¶ Create a mutable and independent copy of the header.
Keyword arguments can be passed on as needed by possible subclasses.
-
classmethod
fromfile
(fh, *args, **kwargs) [edit on github]¶ Read VLBI Header from file.
Arguments are the same as for class initialisation. The header constructed will be immutable.
-
classmethod
fromkeys
(*args, **kwargs) [edit on github]¶ Initialise a header from parsed values.
Like fromvalues, but without any interpretation of keywords.
- Raises
- KeyErrorif not all keys required are present in
kwargs
- KeyErrorif not all keys required are present in
-
classmethod
fromvalues
(*, verify=True, **kwargs)[source] [edit on github]¶ Initialise a header from parsed values.
Here, the parsed values must be given as keyword arguments, i.e., for any
header = cls(<data>)
,cls.fromvalues(**header) == header
.However, unlike for the
Mark5BHeader.fromkeys()
class method, data can also be set using arguments named after methods, such asjday
andseconds
.Given defaults:
sync_pattern : 0xABADDEED
Values set by other keyword arguments (if present):
bcd_jday : from
jday
ortime
bcd_seconds : fromseconds
ortime
bcd_fraction : fromfraction
ortime
(may needframe_rate
) frame_nr : fromtime
(may needframe_rate
)
-
get_time
(frame_rate=None)[source] [edit on github]¶ Convert year, BCD time code to Time object.
Calculate time using
jday
,seconds
, andfraction
properties (which reflect the bcd-encoded ‘bcd_jday’, ‘bcd_seconds’ and ‘bcd_fraction’ header items), pluskday
from the initialisation. See https://www.haystack.mit.edu/tech/vlbi/mark5/docs/Mark%205B%20users%20manual.pdfNote that some non-compliant files do not have ‘bcd_fraction’ set. For those, the time can still be calculated using the header’s ‘frame_nr’ by passing in a frame rate.
Furthermore, fractional seconds are stored only to 0.1 ms accuracy. In the code, this is “unrounded” to give the exact time of the start of the frame for any total bit rate below 512 Mbps. For higher rates, it is no longer guaranteed that subsequent frames have unique
fraction
, and one should pass in an explicit frame rate instead.
-
infer_kday
(ref_time)[source] [edit on github]¶ Uses a reference time to set a header’s
kday
.- Parameters
- ref_time
Time
Reference time within 500 days of the observation time.
- ref_time
-
classmethod
invariant_pattern
(invariants=None, **kwargs) [edit on github]¶ Pattern and mask shared between headers of a type or stream.
This is mostly for use inside
locate_frames()
.- Parameters
- invariantsset of str, optional
Set of keys to header parts that are shared between all headers of a given type or within a given stream/file. Default: from
invariants()
.- **kwargs
Keyword arguments needed to instantiate an empty header. (Mostly for Mark 4).
- Returns
- patternlist of int
The pattern that is shared between headers. If called on an instance, just the header words; if called on a class, words with defaults for the relevant parts set.
- masklist of int
For each entry in
pattern
a bit mask with bits set for the parts that are invariant.
-
classmethod
invariants
() [edit on github]¶ Set of keys of invariant header parts.
On the class, this returns keys of parts that are shared by all headers for the type, on an instance, those that are shared with other headers in the same file.
If neither are defined, returns ‘sync_pattern’ if the header containts that key.
-
keys
() [edit on github]¶ All keys defined for this header.
-
set_time
(time, frame_rate=None)[source] [edit on github]¶ Convert Time object to BCD timestamp elements and ‘frame_nr’.
For non-integer seconds, the frame number will be calculated if not given explicitly. Doing so requires the frame rate.
-
tofile
(fh) [edit on github]¶ Write VLBI frame header to filehandle.
-
update
(*, time=None, frame_rate=None, crc=None, verify=True, **kwargs)[source] [edit on github]¶ Update the header by setting keywords or properties.
Here, any keywords matching header keys are applied first, and any remaining ones are used to set header properties, in the order set by the class (in
_properties
).- Parameters
- time
Time
, optional A possible new time. This is updated last.
- frame_rate
Quantity
, optional The frame rate to use in calculating the frame number from the time. Needed for times at non-integer seconds.
- crcint or None, optional
If
None
(default), recalculate the CRC after updating.- verifybool, optional
If
True
(default), verify integrity after updating.- **kwargs
Arguments used to set keywords and properties.
- time
-
verify
()[source] [edit on github]¶ Verify header integrity.