1# Licensed under a 3-clause BSD style license - see LICENSE.rst
2# TODO: Test FITS parsing
3
4# STDLIB
5import io
6import re
7import gzip
8import base64
9import codecs
10import urllib.request
11import warnings
12
13# THIRD-PARTY
14import numpy as np
15from numpy import ma
16
17# LOCAL
18from astropy.io import fits
19from astropy import __version__ as astropy_version
20from astropy.utils.collections import HomogeneousList
21from astropy.utils.xml.writer import XMLWriter
22from astropy.utils.exceptions import AstropyDeprecationWarning
23
24from . import converters
25from .exceptions import (warn_or_raise, vo_warn, vo_raise, vo_reraise,
26                         warn_unknown_attrs, W06, W07, W08, W09, W10, W11, W12,
27                         W13, W15, W17, W18, W19, W20, W21, W22, W26, W27, W28,
28                         W29, W32, W33, W35, W36, W37, W38, W40, W41, W42, W43,
29                         W44, W45, W50, W52, W53, W54, E06, E08, E09, E10, E11,
30                         E12, E13, E15, E16, E17, E18, E19, E20, E21, E22, E23,
31                         E25)
32from . import ucd as ucd_mod
33from . import util
34from . import xmlutil
35
36try:
37    from . import tablewriter
38    _has_c_tabledata_writer = True
39except ImportError:
40    _has_c_tabledata_writer = False
41
42
43__all__ = [
44    'Link', 'Info', 'Values', 'Field', 'Param', 'CooSys', 'TimeSys',
45    'FieldRef', 'ParamRef', 'Group', 'Table', 'Resource',
46    'VOTableFile', 'Element'
47    ]
48
49
50# The default number of rows to read in each chunk before converting
51# to an array.
52DEFAULT_CHUNK_SIZE = 256
53RESIZE_AMOUNT = 1.5
54
55######################################################################
56# FACTORY FUNCTIONS
57
58
59def _resize(masked, new_size):
60    """
61    Masked arrays can not be resized inplace, and `np.resize` and
62    `ma.resize` are both incompatible with structured arrays.
63    Therefore, we do all this.
64    """
65    new_array = ma.zeros((new_size,), dtype=masked.dtype)
66    length = min(len(masked), new_size)
67    new_array[:length] = masked[:length]
68
69    return new_array
70
71
72def _lookup_by_attr_factory(attr, unique, iterator, element_name, doc):
73    """
74    Creates a function useful for looking up an element by a given
75    attribute.
76
77    Parameters
78    ----------
79    attr : str
80        The attribute name
81
82    unique : bool
83        Should be `True` if the attribute is unique and therefore this
84        should return only one value.  Otherwise, returns a list of
85        values.
86
87    iterator : generator
88        A generator that iterates over some arbitrary set of elements
89
90    element_name : str
91        The XML element name of the elements being iterated over (used
92        for error messages only).
93
94    doc : str
95        A docstring to apply to the generated function.
96
97    Returns
98    -------
99    factory : function
100        A function that looks up an element by the given attribute.
101    """
102
103    def lookup_by_attr(self, ref, before=None):
104        """
105        Given a string *ref*, finds the first element in the iterator
106        where the given attribute == *ref*.  If *before* is provided,
107        will stop searching at the object *before*.  This is
108        important, since "forward references" are not allowed in the
109        VOTABLE format.
110        """
111        for element in getattr(self, iterator)():
112            if element is before:
113                if getattr(element, attr, None) == ref:
114                    vo_raise(
115                        f"{element_name} references itself",
116                        element._config, element._pos, KeyError)
117                break
118            if getattr(element, attr, None) == ref:
119                yield element
120
121    def lookup_by_attr_unique(self, ref, before=None):
122        for element in lookup_by_attr(self, ref, before=before):
123            return element
124        raise KeyError(
125            "No {} with {} '{}' found before the referencing {}".format(
126                element_name, attr, ref, element_name))
127
128    if unique:
129        lookup_by_attr_unique.__doc__ = doc
130        return lookup_by_attr_unique
131    else:
132        lookup_by_attr.__doc__ = doc
133        return lookup_by_attr
134
135
136def _lookup_by_id_or_name_factory(iterator, element_name, doc):
137    """
138    Like `_lookup_by_attr_factory`, but looks in both the "ID" and
139    "name" attributes.
140    """
141
142    def lookup_by_id_or_name(self, ref, before=None):
143        """
144        Given an key *ref*, finds the first element in the iterator
145        with the attribute ID == *ref* or name == *ref*.  If *before*
146        is provided, will stop searching at the object *before*.  This
147        is important, since "forward references" are not allowed in
148        the VOTABLE format.
149        """
150        for element in getattr(self, iterator)():
151            if element is before:
152                if ref in (element.ID, element.name):
153                    vo_raise(
154                        f"{element_name} references itself",
155                        element._config, element._pos, KeyError)
156                break
157            if ref in (element.ID, element.name):
158                return element
159        raise KeyError(
160            "No {} with ID or name '{}' found before the referencing {}".format(
161                element_name, ref, element_name))
162
163    lookup_by_id_or_name.__doc__ = doc
164    return lookup_by_id_or_name
165
166
167def _get_default_unit_format(config):
168    """
169    Get the default unit format as specified in the VOTable spec.
170    """
171    # The unit format changed between VOTable versions 1.3 and 1.4,
172    # see issue #10791.
173    if config['version_1_4_or_later']:
174        return 'vounit'
175    else:
176        return 'cds'
177
178
179def _get_unit_format(config):
180    """
181    Get the unit format based on the configuration.
182    """
183    if config.get('unit_format') is None:
184        format = _get_default_unit_format(config)
185    else:
186        format = config['unit_format']
187    return format
188
189
190######################################################################
191# ATTRIBUTE CHECKERS
192def check_astroyear(year, field, config=None, pos=None):
193    """
194    Raises a `~astropy.io.votable.exceptions.VOTableSpecError` if
195    *year* is not a valid astronomical year as defined by the VOTABLE
196    standard.
197
198    Parameters
199    ----------
200    year : str
201        An astronomical year string
202
203    field : str
204        The name of the field this year was found in (used for error
205        message)
206
207    config, pos : optional
208        Information about the source of the value
209    """
210    if (year is not None and
211        re.match(r"^[JB]?[0-9]+([.][0-9]*)?$", year) is None):
212        warn_or_raise(W07, W07, (field, year), config, pos)
213        return False
214    return True
215
216
217def check_string(string, attr_name, config=None, pos=None):
218    """
219    Raises a `~astropy.io.votable.exceptions.VOTableSpecError` if
220    *string* is not a string or Unicode string.
221
222    Parameters
223    ----------
224    string : str
225        An astronomical year string
226
227    attr_name : str
228        The name of the field this year was found in (used for error
229        message)
230
231    config, pos : optional
232        Information about the source of the value
233    """
234    if string is not None and not isinstance(string, str):
235        warn_or_raise(W08, W08, attr_name, config, pos)
236        return False
237    return True
238
239
240def resolve_id(ID, id, config=None, pos=None):
241    if ID is None and id is not None:
242        warn_or_raise(W09, W09, (), config, pos)
243        return id
244    return ID
245
246
247def check_ucd(ucd, config=None, pos=None):
248    """
249    Warns or raises a
250    `~astropy.io.votable.exceptions.VOTableSpecError` if *ucd* is not
251    a valid `unified content descriptor`_ string as defined by the
252    VOTABLE standard.
253
254    Parameters
255    ----------
256    ucd : str
257        A UCD string.
258
259    config, pos : optional
260        Information about the source of the value
261    """
262    if config is None:
263        config = {}
264    if config.get('version_1_1_or_later'):
265        try:
266            ucd_mod.parse_ucd(
267                ucd,
268                check_controlled_vocabulary=config.get(
269                    'version_1_2_or_later', False),
270                has_colon=config.get('version_1_2_or_later', False))
271        except ValueError as e:
272            # This weird construction is for Python 3 compatibility
273            if config.get('verify', 'ignore') == 'exception':
274                vo_raise(W06, (ucd, str(e)), config, pos)
275            elif config.get('verify', 'ignore') == 'warn':
276                vo_warn(W06, (ucd, str(e)), config, pos)
277                return False
278            else:
279                return False
280    return True
281
282
283######################################################################
284# PROPERTY MIXINS
285class _IDProperty:
286    @property
287    def ID(self):
288        """
289        The XML ID_ of the element.  May be `None` or a string
290        conforming to XML ID_ syntax.
291        """
292        return self._ID
293
294    @ID.setter
295    def ID(self, ID):
296        xmlutil.check_id(ID, 'ID', self._config, self._pos)
297        self._ID = ID
298
299    @ID.deleter
300    def ID(self):
301        self._ID = None
302
303
304class _NameProperty:
305    @property
306    def name(self):
307        """An optional name for the element."""
308        return self._name
309
310    @name.setter
311    def name(self, name):
312        xmlutil.check_token(name, 'name', self._config, self._pos)
313        self._name = name
314
315    @name.deleter
316    def name(self):
317        self._name = None
318
319
320class _XtypeProperty:
321    @property
322    def xtype(self):
323        """Extended data type information."""
324        return self._xtype
325
326    @xtype.setter
327    def xtype(self, xtype):
328        if xtype is not None and not self._config.get('version_1_2_or_later'):
329            warn_or_raise(
330                W28, W28, ('xtype', self._element_name, '1.2'),
331                self._config, self._pos)
332        check_string(xtype, 'xtype', self._config, self._pos)
333        self._xtype = xtype
334
335    @xtype.deleter
336    def xtype(self):
337        self._xtype = None
338
339
340class _UtypeProperty:
341    _utype_in_v1_2 = False
342
343    @property
344    def utype(self):
345        """The usage-specific or `unique type`_ of the element."""
346        return self._utype
347
348    @utype.setter
349    def utype(self, utype):
350        if (self._utype_in_v1_2 and
351            utype is not None and
352            not self._config.get('version_1_2_or_later')):
353            warn_or_raise(
354                W28, W28, ('utype', self._element_name, '1.2'),
355                self._config, self._pos)
356        check_string(utype, 'utype', self._config, self._pos)
357        self._utype = utype
358
359    @utype.deleter
360    def utype(self):
361        self._utype = None
362
363
364class _UcdProperty:
365    _ucd_in_v1_2 = False
366
367    @property
368    def ucd(self):
369        """The `unified content descriptor`_ for the element."""
370        return self._ucd
371
372    @ucd.setter
373    def ucd(self, ucd):
374        if ucd is not None and ucd.strip() == '':
375            ucd = None
376        if ucd is not None:
377            if (self._ucd_in_v1_2 and
378                not self._config.get('version_1_2_or_later')):
379                warn_or_raise(
380                    W28, W28, ('ucd', self._element_name, '1.2'),
381                    self._config, self._pos)
382            check_ucd(ucd, self._config, self._pos)
383        self._ucd = ucd
384
385    @ucd.deleter
386    def ucd(self):
387        self._ucd = None
388
389
390class _DescriptionProperty:
391    @property
392    def description(self):
393        """
394        An optional string describing the element.  Corresponds to the
395        DESCRIPTION_ element.
396        """
397        return self._description
398
399    @description.setter
400    def description(self, description):
401        self._description = description
402
403    @description.deleter
404    def description(self):
405        self._description = None
406
407
408######################################################################
409# ELEMENT CLASSES
410class Element:
411    """
412    A base class for all classes that represent XML elements in the
413    VOTABLE file.
414    """
415    _element_name = ''
416    _attr_list = []
417
418    def _add_unknown_tag(self, iterator, tag, data, config, pos):
419        warn_or_raise(W10, W10, tag, config, pos)
420
421    def _ignore_add(self, iterator, tag, data, config, pos):
422        warn_unknown_attrs(tag, data.keys(), config, pos)
423
424    def _add_definitions(self, iterator, tag, data, config, pos):
425        if config.get('version_1_1_or_later'):
426            warn_or_raise(W22, W22, (), config, pos)
427        warn_unknown_attrs(tag, data.keys(), config, pos)
428
429    def parse(self, iterator, config):
430        """
431        For internal use. Parse the XML content of the children of the
432        element.
433
434        Parameters
435        ----------
436        iterator : xml iterable
437            An iterator over XML elements as returned by
438            `~astropy.utils.xml.iterparser.get_xml_iterator`.
439
440        config : dict
441            The configuration dictionary that affects how certain
442            elements are read.
443
444        Returns
445        -------
446        self : `~astropy.io.votable.tree.Element`
447            Returns self as a convenience.
448        """
449        raise NotImplementedError()
450
451    def to_xml(self, w, **kwargs):
452        """
453        For internal use. Output the element to XML.
454
455        Parameters
456        ----------
457        w : astropy.utils.xml.writer.XMLWriter object
458            An XML writer to write to.
459
460        kwargs : dict
461            Any configuration parameters to control the output.
462        """
463        raise NotImplementedError()
464
465
466class SimpleElement(Element):
467    """
468    A base class for simple elements, such as FIELD, PARAM and INFO
469    that don't require any special parsing or outputting machinery.
470    """
471
472    def __init__(self):
473        Element.__init__(self)
474
475    def __repr__(self):
476        buff = io.StringIO()
477        SimpleElement.to_xml(self, XMLWriter(buff))
478        return buff.getvalue().strip()
479
480    def parse(self, iterator, config):
481        for start, tag, data, pos in iterator:
482            if start and tag != self._element_name:
483                self._add_unknown_tag(iterator, tag, data, config, pos)
484            elif tag == self._element_name:
485                break
486
487        return self
488
489    def to_xml(self, w, **kwargs):
490        w.element(self._element_name,
491                  attrib=w.object_attrs(self, self._attr_list))
492
493
494class SimpleElementWithContent(SimpleElement):
495    """
496    A base class for simple elements, such as FIELD, PARAM and INFO
497    that don't require any special parsing or outputting machinery.
498    """
499
500    def __init__(self):
501        SimpleElement.__init__(self)
502
503        self._content = None
504
505    def parse(self, iterator, config):
506        for start, tag, data, pos in iterator:
507            if start and tag != self._element_name:
508                self._add_unknown_tag(iterator, tag, data, config, pos)
509            elif tag == self._element_name:
510                if data:
511                    self.content = data
512                break
513
514        return self
515
516    def to_xml(self, w, **kwargs):
517        w.element(self._element_name, self._content,
518                  attrib=w.object_attrs(self, self._attr_list))
519
520    @property
521    def content(self):
522        """The content of the element."""
523        return self._content
524
525    @content.setter
526    def content(self, content):
527        check_string(content, 'content', self._config, self._pos)
528        self._content = content
529
530    @content.deleter
531    def content(self):
532        self._content = None
533
534
535class Link(SimpleElement, _IDProperty):
536    """
537    LINK_ elements: used to reference external documents and servers through a URI.
538
539    The keyword arguments correspond to setting members of the same
540    name, documented below.
541    """
542    _attr_list = ['ID', 'content_role', 'content_type', 'title', 'value',
543                  'href', 'action']
544    _element_name = 'LINK'
545
546    def __init__(self, ID=None, title=None, value=None, href=None, action=None,
547                 id=None, config=None, pos=None, **kwargs):
548        if config is None:
549            config = {}
550        self._config = config
551        self._pos = pos
552
553        SimpleElement.__init__(self)
554
555        content_role = kwargs.get('content-role') or kwargs.get('content_role')
556        content_type = kwargs.get('content-type') or kwargs.get('content_type')
557
558        if 'gref' in kwargs:
559            warn_or_raise(W11, W11, (), config, pos)
560
561        self.ID = resolve_id(ID, id, config, pos)
562        self.content_role = content_role
563        self.content_type = content_type
564        self.title = title
565        self.value = value
566        self.href = href
567        self.action = action
568
569        warn_unknown_attrs(
570            'LINK', kwargs.keys(), config, pos,
571            ['content-role', 'content_role', 'content-type', 'content_type',
572             'gref'])
573
574    @property
575    def content_role(self):
576        """
577        Defines the MIME role of the referenced object.  Must be one of:
578
579          None, 'query', 'hints', 'doc', 'location' or 'type'
580        """
581        return self._content_role
582
583    @content_role.setter
584    def content_role(self, content_role):
585        if ((content_role == 'type' and
586             not self._config['version_1_3_or_later']) or
587             content_role not in
588             (None, 'query', 'hints', 'doc', 'location')):
589            vo_warn(W45, (content_role,), self._config, self._pos)
590        self._content_role = content_role
591
592    @content_role.deleter
593    def content_role(self):
594        self._content_role = None
595
596    @property
597    def content_type(self):
598        """Defines the MIME content type of the referenced object."""
599        return self._content_type
600
601    @content_type.setter
602    def content_type(self, content_type):
603        xmlutil.check_mime_content_type(content_type, self._config, self._pos)
604        self._content_type = content_type
605
606    @content_type.deleter
607    def content_type(self):
608        self._content_type = None
609
610    @property
611    def href(self):
612        """
613        A URI to an arbitrary protocol.  The vo package only supports
614        http and anonymous ftp.
615        """
616        return self._href
617
618    @href.setter
619    def href(self, href):
620        xmlutil.check_anyuri(href, self._config, self._pos)
621        self._href = href
622
623    @href.deleter
624    def href(self):
625        self._href = None
626
627    def to_table_column(self, column):
628        meta = {}
629        for key in self._attr_list:
630            val = getattr(self, key, None)
631            if val is not None:
632                meta[key] = val
633
634        column.meta.setdefault('links', [])
635        column.meta['links'].append(meta)
636
637    @classmethod
638    def from_table_column(cls, d):
639        return cls(**d)
640
641
642class Info(SimpleElementWithContent, _IDProperty, _XtypeProperty,
643           _UtypeProperty):
644    """
645    INFO_ elements: arbitrary key-value pairs for extensions to the standard.
646
647    The keyword arguments correspond to setting members of the same
648    name, documented below.
649    """
650    _element_name = 'INFO'
651    _attr_list_11 = ['ID', 'name', 'value']
652    _attr_list_12 = _attr_list_11 + ['xtype', 'ref', 'unit', 'ucd', 'utype']
653    _utype_in_v1_2 = True
654
655    def __init__(self, ID=None, name=None, value=None, id=None, xtype=None,
656                 ref=None, unit=None, ucd=None, utype=None,
657                 config=None, pos=None, **extra):
658        if config is None:
659            config = {}
660        self._config = config
661        self._pos = pos
662
663        SimpleElementWithContent.__init__(self)
664
665        self.ID = (resolve_id(ID, id, config, pos) or
666                        xmlutil.fix_id(name, config, pos))
667        self.name = name
668        self.value = value
669        self.xtype = xtype
670        self.ref = ref
671        self.unit = unit
672        self.ucd = ucd
673        self.utype = utype
674
675        if config.get('version_1_2_or_later'):
676            self._attr_list = self._attr_list_12
677        else:
678            self._attr_list = self._attr_list_11
679            if xtype is not None:
680                warn_unknown_attrs('INFO', ['xtype'], config, pos)
681            if ref is not None:
682                warn_unknown_attrs('INFO', ['ref'], config, pos)
683            if unit is not None:
684                warn_unknown_attrs('INFO', ['unit'], config, pos)
685            if ucd is not None:
686                warn_unknown_attrs('INFO', ['ucd'], config, pos)
687            if utype is not None:
688                warn_unknown_attrs('INFO', ['utype'], config, pos)
689
690        warn_unknown_attrs('INFO', extra.keys(), config, pos)
691
692    @property
693    def name(self):
694        """[*required*] The key of the key-value pair."""
695        return self._name
696
697    @name.setter
698    def name(self, name):
699        if name is None:
700            warn_or_raise(W35, W35, ('name'), self._config, self._pos)
701        xmlutil.check_token(name, 'name', self._config, self._pos)
702        self._name = name
703
704    @property
705    def value(self):
706        """
707        [*required*] The value of the key-value pair.  (Always stored
708        as a string or unicode string).
709        """
710        return self._value
711
712    @value.setter
713    def value(self, value):
714        if value is None:
715            warn_or_raise(W35, W35, ('value'), self._config, self._pos)
716        check_string(value, 'value', self._config, self._pos)
717        self._value = value
718
719    @property
720    def content(self):
721        """The content inside the INFO element."""
722        return self._content
723
724    @content.setter
725    def content(self, content):
726        check_string(content, 'content', self._config, self._pos)
727        self._content = content
728
729    @content.deleter
730    def content(self):
731        self._content = None
732
733    @property
734    def ref(self):
735        """
736        Refer to another INFO_ element by ID_, defined previously in
737        the document.
738        """
739        return self._ref
740
741    @ref.setter
742    def ref(self, ref):
743        if ref is not None and not self._config.get('version_1_2_or_later'):
744            warn_or_raise(W28, W28, ('ref', 'INFO', '1.2'),
745                          self._config, self._pos)
746        xmlutil.check_id(ref, 'ref', self._config, self._pos)
747        # TODO: actually apply the reference
748        # if ref is not None:
749        #     try:
750        #         other = self._votable.get_values_by_id(ref, before=self)
751        #     except KeyError:
752        #         vo_raise(
753        #             "VALUES ref='%s', which has not already been defined." %
754        #             self.ref, self._config, self._pos, KeyError)
755        #     self.null = other.null
756        #     self.type = other.type
757        #     self.min = other.min
758        #     self.min_inclusive = other.min_inclusive
759        #     self.max = other.max
760        #     self.max_inclusive = other.max_inclusive
761        #     self._options[:] = other.options
762        self._ref = ref
763
764    @ref.deleter
765    def ref(self):
766        self._ref = None
767
768    @property
769    def unit(self):
770        """A string specifying the units_ for the INFO_."""
771        return self._unit
772
773    @unit.setter
774    def unit(self, unit):
775        if unit is None:
776            self._unit = None
777            return
778
779        from astropy import units as u
780
781        if not self._config.get('version_1_2_or_later'):
782            warn_or_raise(W28, W28, ('unit', 'INFO', '1.2'),
783                          self._config, self._pos)
784
785        # First, parse the unit in the default way, so that we can
786        # still emit a warning if the unit is not to spec.
787        default_format = _get_default_unit_format(self._config)
788        unit_obj = u.Unit(
789            unit, format=default_format, parse_strict='silent')
790        if isinstance(unit_obj, u.UnrecognizedUnit):
791            warn_or_raise(W50, W50, (unit,),
792                          self._config, self._pos)
793
794        format = _get_unit_format(self._config)
795        if format != default_format:
796            unit_obj = u.Unit(
797                unit, format=format, parse_strict='silent')
798
799        self._unit = unit_obj
800
801    @unit.deleter
802    def unit(self):
803        self._unit = None
804
805    def to_xml(self, w, **kwargs):
806        attrib = w.object_attrs(self, self._attr_list)
807        if 'unit' in attrib:
808            attrib['unit'] = self.unit.to_string('cds')
809        w.element(self._element_name, self._content,
810                  attrib=attrib)
811
812
813class Values(Element, _IDProperty):
814    """
815    VALUES_ element: used within FIELD_ and PARAM_ elements to define the domain of values.
816
817    The keyword arguments correspond to setting members of the same
818    name, documented below.
819    """
820
821    def __init__(self, votable, field, ID=None, null=None, ref=None,
822                 type="legal", id=None, config=None, pos=None, **extras):
823        if config is None:
824            config = {}
825        self._config = config
826        self._pos = pos
827
828        Element.__init__(self)
829
830        self._votable = votable
831        self._field = field
832        self.ID = resolve_id(ID, id, config, pos)
833        self.null = null
834        self._ref = ref
835        self.type = type
836
837        self.min = None
838        self.max = None
839        self.min_inclusive = True
840        self.max_inclusive = True
841        self._options = []
842
843        warn_unknown_attrs('VALUES', extras.keys(), config, pos)
844
845    def __repr__(self):
846        buff = io.StringIO()
847        self.to_xml(XMLWriter(buff))
848        return buff.getvalue().strip()
849
850    @property
851    def null(self):
852        """
853        For integral datatypes, *null* is used to define the value
854        used for missing values.
855        """
856        return self._null
857
858    @null.setter
859    def null(self, null):
860        if null is not None and isinstance(null, str):
861            try:
862                null_val = self._field.converter.parse_scalar(
863                    null, self._config, self._pos)[0]
864            except Exception:
865                warn_or_raise(W36, W36, null, self._config, self._pos)
866                null_val = self._field.converter.parse_scalar(
867                    '0', self._config, self._pos)[0]
868        else:
869            null_val = null
870        self._null = null_val
871
872    @null.deleter
873    def null(self):
874        self._null = None
875
876    @property
877    def type(self):
878        """
879        [*required*] Defines the applicability of the domain defined
880        by this VALUES_ element.  Must be one of the following
881        strings:
882
883          - 'legal': The domain of this column applies in general to
884            this datatype. (default)
885
886          - 'actual': The domain of this column applies only to the
887            data enclosed in the parent table.
888        """
889        return self._type
890
891    @type.setter
892    def type(self, type):
893        if type not in ('legal', 'actual'):
894            vo_raise(E08, type, self._config, self._pos)
895        self._type = type
896
897    @property
898    def ref(self):
899        """
900        Refer to another VALUES_ element by ID_, defined previously in
901        the document, for MIN/MAX/OPTION information.
902        """
903        return self._ref
904
905    @ref.setter
906    def ref(self, ref):
907        xmlutil.check_id(ref, 'ref', self._config, self._pos)
908        if ref is not None:
909            try:
910                other = self._votable.get_values_by_id(ref, before=self)
911            except KeyError:
912                warn_or_raise(W43, W43, ('VALUES', self.ref), self._config,
913                              self._pos)
914                ref = None
915            else:
916                self.null = other.null
917                self.type = other.type
918                self.min = other.min
919                self.min_inclusive = other.min_inclusive
920                self.max = other.max
921                self.max_inclusive = other.max_inclusive
922                self._options[:] = other.options
923        self._ref = ref
924
925    @ref.deleter
926    def ref(self):
927        self._ref = None
928
929    @property
930    def min(self):
931        """
932        The minimum value of the domain.  See :attr:`min_inclusive`.
933        """
934        return self._min
935
936    @min.setter
937    def min(self, min):
938        if hasattr(self._field, 'converter') and min is not None:
939            self._min = self._field.converter.parse(min)[0]
940        else:
941            self._min = min
942
943    @min.deleter
944    def min(self):
945        self._min = None
946
947    @property
948    def min_inclusive(self):
949        """When `True`, the domain includes the minimum value."""
950        return self._min_inclusive
951
952    @min_inclusive.setter
953    def min_inclusive(self, inclusive):
954        if inclusive == 'yes':
955            self._min_inclusive = True
956        elif inclusive == 'no':
957            self._min_inclusive = False
958        else:
959            self._min_inclusive = bool(inclusive)
960
961    @min_inclusive.deleter
962    def min_inclusive(self):
963        self._min_inclusive = True
964
965    @property
966    def max(self):
967        """
968        The maximum value of the domain.  See :attr:`max_inclusive`.
969        """
970        return self._max
971
972    @max.setter
973    def max(self, max):
974        if hasattr(self._field, 'converter') and max is not None:
975            self._max = self._field.converter.parse(max)[0]
976        else:
977            self._max = max
978
979    @max.deleter
980    def max(self):
981        self._max = None
982
983    @property
984    def max_inclusive(self):
985        """When `True`, the domain includes the maximum value."""
986        return self._max_inclusive
987
988    @max_inclusive.setter
989    def max_inclusive(self, inclusive):
990        if inclusive == 'yes':
991            self._max_inclusive = True
992        elif inclusive == 'no':
993            self._max_inclusive = False
994        else:
995            self._max_inclusive = bool(inclusive)
996
997    @max_inclusive.deleter
998    def max_inclusive(self):
999        self._max_inclusive = True
1000
1001    @property
1002    def options(self):
1003        """
1004        A list of string key-value tuples defining other OPTION
1005        elements for the domain.  All options are ignored -- they are
1006        stored for round-tripping purposes only.
1007        """
1008        return self._options
1009
1010    def parse(self, iterator, config):
1011        if self.ref is not None:
1012            for start, tag, data, pos in iterator:
1013                if start:
1014                    warn_or_raise(W44, W44, tag, config, pos)
1015                else:
1016                    if tag != 'VALUES':
1017                        warn_or_raise(W44, W44, tag, config, pos)
1018                    break
1019        else:
1020            for start, tag, data, pos in iterator:
1021                if start:
1022                    if tag == 'MIN':
1023                        if 'value' not in data:
1024                            vo_raise(E09, 'MIN', config, pos)
1025                        self.min = data['value']
1026                        self.min_inclusive = data.get('inclusive', 'yes')
1027                        warn_unknown_attrs(
1028                            'MIN', data.keys(), config, pos,
1029                            ['value', 'inclusive'])
1030                    elif tag == 'MAX':
1031                        if 'value' not in data:
1032                            vo_raise(E09, 'MAX', config, pos)
1033                        self.max = data['value']
1034                        self.max_inclusive = data.get('inclusive', 'yes')
1035                        warn_unknown_attrs(
1036                            'MAX', data.keys(), config, pos,
1037                            ['value', 'inclusive'])
1038                    elif tag == 'OPTION':
1039                        if 'value' not in data:
1040                            vo_raise(E09, 'OPTION', config, pos)
1041                        xmlutil.check_token(
1042                            data.get('name'), 'name', config, pos)
1043                        self.options.append(
1044                            (data.get('name'), data.get('value')))
1045                        warn_unknown_attrs(
1046                            'OPTION', data.keys(), config, pos,
1047                            ['value', 'name'])
1048                elif tag == 'VALUES':
1049                    break
1050
1051        return self
1052
1053    def is_defaults(self):
1054        """
1055        Are the settings on this ``VALUE`` element all the same as the
1056        XML defaults?
1057        """
1058        # If there's nothing meaningful or non-default to write,
1059        # don't write anything.
1060        return (self.ref is None and self.null is None and self.ID is None and
1061                self.max is None and self.min is None and self.options == [])
1062
1063    def to_xml(self, w, **kwargs):
1064        def yes_no(value):
1065            if value:
1066                return 'yes'
1067            return 'no'
1068
1069        if self.is_defaults():
1070            return
1071
1072        if self.ref is not None:
1073            w.element('VALUES', attrib=w.object_attrs(self, ['ref']))
1074        else:
1075            with w.tag('VALUES',
1076                       attrib=w.object_attrs(
1077                           self, ['ID', 'null', 'ref'])):
1078                if self.min is not None:
1079                    w.element(
1080                        'MIN',
1081                        value=self._field.converter.output(self.min, False),
1082                        inclusive=yes_no(self.min_inclusive))
1083                if self.max is not None:
1084                    w.element(
1085                        'MAX',
1086                        value=self._field.converter.output(self.max, False),
1087                        inclusive=yes_no(self.max_inclusive))
1088                for name, value in self.options:
1089                    w.element(
1090                        'OPTION',
1091                        name=name,
1092                        value=value)
1093
1094    def to_table_column(self, column):
1095        # Have the ref filled in here
1096        meta = {}
1097        for key in ['ID', 'null']:
1098            val = getattr(self, key, None)
1099            if val is not None:
1100                meta[key] = val
1101        if self.min is not None:
1102            meta['min'] = {
1103                'value': self.min,
1104                'inclusive': self.min_inclusive}
1105        if self.max is not None:
1106            meta['max'] = {
1107                'value': self.max,
1108                'inclusive': self.max_inclusive}
1109        if len(self.options):
1110            meta['options'] = dict(self.options)
1111
1112        column.meta['values'] = meta
1113
1114    def from_table_column(self, column):
1115        if column.info.meta is None or 'values' not in column.info.meta:
1116            return
1117
1118        meta = column.info.meta['values']
1119        for key in ['ID', 'null']:
1120            val = meta.get(key, None)
1121            if val is not None:
1122                setattr(self, key, val)
1123        if 'min' in meta:
1124            self.min = meta['min']['value']
1125            self.min_inclusive = meta['min']['inclusive']
1126        if 'max' in meta:
1127            self.max = meta['max']['value']
1128            self.max_inclusive = meta['max']['inclusive']
1129        if 'options' in meta:
1130            self._options = list(meta['options'].items())
1131
1132
1133class Field(SimpleElement, _IDProperty, _NameProperty, _XtypeProperty,
1134            _UtypeProperty, _UcdProperty):
1135    """
1136    FIELD_ element: describes the datatype of a particular column of data.
1137
1138    The keyword arguments correspond to setting members of the same
1139    name, documented below.
1140
1141    If *ID* is provided, it is used for the column name in the
1142    resulting recarray of the table.  If no *ID* is provided, *name*
1143    is used instead.  If neither is provided, an exception will be
1144    raised.
1145    """
1146    _attr_list_11 = ['ID', 'name', 'datatype', 'arraysize', 'ucd',
1147                     'unit', 'width', 'precision', 'utype', 'ref']
1148    _attr_list_12 = _attr_list_11 + ['xtype']
1149    _element_name = 'FIELD'
1150
1151    def __init__(self, votable, ID=None, name=None, datatype=None,
1152                 arraysize=None, ucd=None, unit=None, width=None,
1153                 precision=None, utype=None, ref=None, type=None, id=None,
1154                 xtype=None,
1155                 config=None, pos=None, **extra):
1156        if config is None:
1157            if hasattr(votable, '_get_version_checks'):
1158                config = votable._get_version_checks()
1159            else:
1160                config = {}
1161        self._config = config
1162        self._pos = pos
1163
1164        SimpleElement.__init__(self)
1165
1166        if config.get('version_1_2_or_later'):
1167            self._attr_list = self._attr_list_12
1168        else:
1169            self._attr_list = self._attr_list_11
1170            if xtype is not None:
1171                warn_unknown_attrs(self._element_name, ['xtype'], config, pos)
1172
1173        # TODO: REMOVE ME ----------------------------------------
1174        # This is a terrible hack to support Simple Image Access
1175        # Protocol results from archive.noao.edu.  It creates a field
1176        # for the coordinate projection type of type "double", which
1177        # actually contains character data.  We have to hack the field
1178        # to store character data, or we can't read it in.  A warning
1179        # will be raised when this happens.
1180        if (config.get('verify', 'ignore') != 'exception' and name == 'cprojection' and
1181            ID == 'cprojection' and ucd == 'VOX:WCS_CoordProjection' and
1182            datatype == 'double'):
1183            datatype = 'char'
1184            arraysize = '3'
1185            vo_warn(W40, (), config, pos)
1186        # ----------------------------------------
1187
1188        self.description = None
1189        self._votable = votable
1190
1191        self.ID = (resolve_id(ID, id, config, pos) or
1192                   xmlutil.fix_id(name, config, pos))
1193        self.name = name
1194        if name is None:
1195            if (self._element_name == 'PARAM' and
1196                not config.get('version_1_1_or_later')):
1197                pass
1198            else:
1199                warn_or_raise(W15, W15, self._element_name, config, pos)
1200            self.name = self.ID
1201
1202        if self._ID is None and name is None:
1203            vo_raise(W12, self._element_name, config, pos)
1204
1205        datatype_mapping = {
1206            'string': 'char',
1207            'unicodeString': 'unicodeChar',
1208            'int16': 'short',
1209            'int32': 'int',
1210            'int64': 'long',
1211            'float32': 'float',
1212            'float64': 'double',
1213            # The following appear in some Vizier tables
1214            'unsignedInt': 'long',
1215            'unsignedShort': 'int'
1216        }
1217
1218        datatype_mapping.update(config.get('datatype_mapping', {}))
1219
1220        if datatype in datatype_mapping:
1221            warn_or_raise(W13, W13, (datatype, datatype_mapping[datatype]),
1222                          config, pos)
1223            datatype = datatype_mapping[datatype]
1224
1225        self.ref = ref
1226        self.datatype = datatype
1227        self.arraysize = arraysize
1228        self.ucd = ucd
1229        self.unit = unit
1230        self.width = width
1231        self.precision = precision
1232        self.utype = utype
1233        self.type = type
1234        self._links = HomogeneousList(Link)
1235        self.title = self.name
1236        self.values = Values(self._votable, self)
1237        self.xtype = xtype
1238
1239        self._setup(config, pos)
1240
1241        warn_unknown_attrs(self._element_name, extra.keys(), config, pos)
1242
1243    @classmethod
1244    def uniqify_names(cls, fields):
1245        """
1246        Make sure that all names and titles in a list of fields are
1247        unique, by appending numbers if necessary.
1248        """
1249        unique = {}
1250        for field in fields:
1251            i = 2
1252            new_id = field.ID
1253            while new_id in unique:
1254                new_id = field.ID + f"_{i:d}"
1255                i += 1
1256            if new_id != field.ID:
1257                vo_warn(W32, (field.ID, new_id), field._config, field._pos)
1258            field.ID = new_id
1259            unique[new_id] = field.ID
1260
1261        for field in fields:
1262            i = 2
1263            if field.name is None:
1264                new_name = field.ID
1265                implicit = True
1266            else:
1267                new_name = field.name
1268                implicit = False
1269            if new_name != field.ID:
1270                while new_name in unique:
1271                    new_name = field.name + f" {i:d}"
1272                    i += 1
1273
1274            if (not implicit and
1275                new_name != field.name):
1276                vo_warn(W33, (field.name, new_name), field._config, field._pos)
1277            field._unique_name = new_name
1278            unique[new_name] = field.name
1279
1280    def _setup(self, config, pos):
1281        if self.values._ref is not None:
1282            self.values.ref = self.values._ref
1283        self.converter = converters.get_converter(self, config, pos)
1284
1285    @property
1286    def datatype(self):
1287        """
1288        [*required*] The datatype of the column.  Valid values (as
1289        defined by the spec) are:
1290
1291          'boolean', 'bit', 'unsignedByte', 'short', 'int', 'long',
1292          'char', 'unicodeChar', 'float', 'double', 'floatComplex', or
1293          'doubleComplex'
1294
1295        Many VOTABLE files in the wild use 'string' instead of 'char',
1296        so that is also a valid option, though 'string' will always be
1297        converted to 'char' when writing the file back out.
1298        """
1299        return self._datatype
1300
1301    @datatype.setter
1302    def datatype(self, datatype):
1303        if datatype is None:
1304            if self._config.get('version_1_1_or_later'):
1305                warn_or_raise(E10, E10, self._element_name, self._config,
1306                              self._pos)
1307            datatype = 'char'
1308        if datatype not in converters.converter_mapping:
1309            vo_raise(E06, (datatype, self.ID), self._config, self._pos)
1310        self._datatype = datatype
1311
1312    @property
1313    def precision(self):
1314        """
1315        Along with :attr:`width`, defines the `numerical accuracy`_
1316        associated with the data.  These values are used to limit the
1317        precision when writing floating point values back to the XML
1318        file.  Otherwise, it is purely informational -- the Numpy
1319        recarray containing the data itself does not use this
1320        information.
1321        """
1322        return self._precision
1323
1324    @precision.setter
1325    def precision(self, precision):
1326        if precision is not None and not re.match(r"^[FE]?[0-9]+$", precision):
1327            vo_raise(E11, precision, self._config, self._pos)
1328        self._precision = precision
1329
1330    @precision.deleter
1331    def precision(self):
1332        self._precision = None
1333
1334    @property
1335    def width(self):
1336        """
1337        Along with :attr:`precision`, defines the `numerical
1338        accuracy`_ associated with the data.  These values are used to
1339        limit the precision when writing floating point values back to
1340        the XML file.  Otherwise, it is purely informational -- the
1341        Numpy recarray containing the data itself does not use this
1342        information.
1343        """
1344        return self._width
1345
1346    @width.setter
1347    def width(self, width):
1348        if width is not None:
1349            width = int(width)
1350            if width <= 0:
1351                vo_raise(E12, width, self._config, self._pos)
1352        self._width = width
1353
1354    @width.deleter
1355    def width(self):
1356        self._width = None
1357
1358    # ref on FIELD and PARAM behave differently than elsewhere -- here
1359    # they're just informational, such as to refer to a coordinate
1360    # system.
1361    @property
1362    def ref(self):
1363        """
1364        On FIELD_ elements, ref is used only for informational
1365        purposes, for example to refer to a COOSYS_ or TIMESYS_ element.
1366        """
1367        return self._ref
1368
1369    @ref.setter
1370    def ref(self, ref):
1371        xmlutil.check_id(ref, 'ref', self._config, self._pos)
1372        self._ref = ref
1373
1374    @ref.deleter
1375    def ref(self):
1376        self._ref = None
1377
1378    @property
1379    def unit(self):
1380        """A string specifying the units_ for the FIELD_."""
1381        return self._unit
1382
1383    @unit.setter
1384    def unit(self, unit):
1385        if unit is None:
1386            self._unit = None
1387            return
1388
1389        from astropy import units as u
1390
1391        # First, parse the unit in the default way, so that we can
1392        # still emit a warning if the unit is not to spec.
1393        default_format = _get_default_unit_format(self._config)
1394        unit_obj = u.Unit(
1395            unit, format=default_format, parse_strict='silent')
1396        if isinstance(unit_obj, u.UnrecognizedUnit):
1397            warn_or_raise(W50, W50, (unit,),
1398                          self._config, self._pos)
1399
1400        format = _get_unit_format(self._config)
1401        if format != default_format:
1402            unit_obj = u.Unit(
1403                unit, format=format, parse_strict='silent')
1404
1405        self._unit = unit_obj
1406
1407    @unit.deleter
1408    def unit(self):
1409        self._unit = None
1410
1411    @property
1412    def arraysize(self):
1413        """
1414        Specifies the size of the multidimensional array if this
1415        FIELD_ contains more than a single value.
1416
1417        See `multidimensional arrays`_.
1418        """
1419        return self._arraysize
1420
1421    @arraysize.setter
1422    def arraysize(self, arraysize):
1423        if (arraysize is not None and
1424            not re.match(r"^([0-9]+x)*[0-9]*[*]?(s\W)?$", arraysize)):
1425            vo_raise(E13, arraysize, self._config, self._pos)
1426        self._arraysize = arraysize
1427
1428    @arraysize.deleter
1429    def arraysize(self):
1430        self._arraysize = None
1431
1432    @property
1433    def type(self):
1434        """
1435        The type attribute on FIELD_ elements is reserved for future
1436        extensions.
1437        """
1438        return self._type
1439
1440    @type.setter
1441    def type(self, type):
1442        self._type = type
1443
1444    @type.deleter
1445    def type(self):
1446        self._type = None
1447
1448    @property
1449    def values(self):
1450        """
1451        A :class:`Values` instance (or `None`) defining the domain
1452        of the column.
1453        """
1454        return self._values
1455
1456    @values.setter
1457    def values(self, values):
1458        assert values is None or isinstance(values, Values)
1459        self._values = values
1460
1461    @values.deleter
1462    def values(self):
1463        self._values = None
1464
1465    @property
1466    def links(self):
1467        """
1468        A list of :class:`Link` instances used to reference more
1469        details about the meaning of the FIELD_.  This is purely
1470        informational and is not used by the `astropy.io.votable`
1471        package.
1472        """
1473        return self._links
1474
1475    def parse(self, iterator, config):
1476        for start, tag, data, pos in iterator:
1477            if start:
1478                if tag == 'VALUES':
1479                    self.values.__init__(
1480                        self._votable, self, config=config, pos=pos, **data)
1481                    self.values.parse(iterator, config)
1482                elif tag == 'LINK':
1483                    link = Link(config=config, pos=pos, **data)
1484                    self.links.append(link)
1485                    link.parse(iterator, config)
1486                elif tag == 'DESCRIPTION':
1487                    warn_unknown_attrs(
1488                        'DESCRIPTION', data.keys(), config, pos)
1489                elif tag != self._element_name:
1490                    self._add_unknown_tag(iterator, tag, data, config, pos)
1491            else:
1492                if tag == 'DESCRIPTION':
1493                    if self.description is not None:
1494                        warn_or_raise(
1495                            W17, W17, self._element_name, config, pos)
1496                    self.description = data or None
1497                elif tag == self._element_name:
1498                    break
1499
1500        if self.description is not None:
1501            self.title = " ".join(x.strip() for x in
1502                                  self.description.splitlines())
1503        else:
1504            self.title = self.name
1505
1506        self._setup(config, pos)
1507
1508        return self
1509
1510    def to_xml(self, w, **kwargs):
1511        attrib = w.object_attrs(self, self._attr_list)
1512        if 'unit' in attrib:
1513            attrib['unit'] = self.unit.to_string('cds')
1514        with w.tag(self._element_name, attrib=attrib):
1515            if self.description is not None:
1516                w.element('DESCRIPTION', self.description, wrap=True)
1517            if not self.values.is_defaults():
1518                self.values.to_xml(w, **kwargs)
1519            for link in self.links:
1520                link.to_xml(w, **kwargs)
1521
1522    def to_table_column(self, column):
1523        """
1524        Sets the attributes of a given `astropy.table.Column` instance
1525        to match the information in this `Field`.
1526        """
1527        for key in ['ucd', 'width', 'precision', 'utype', 'xtype']:
1528            val = getattr(self, key, None)
1529            if val is not None:
1530                column.meta[key] = val
1531        if not self.values.is_defaults():
1532            self.values.to_table_column(column)
1533        for link in self.links:
1534            link.to_table_column(column)
1535        if self.description is not None:
1536            column.description = self.description
1537        if self.unit is not None:
1538            # TODO: Use units framework when it's available
1539            column.unit = self.unit
1540        if (isinstance(self.converter, converters.FloatingPoint) and
1541                self.converter.output_format != '{!r:>}'):
1542            column.format = self.converter.output_format
1543        elif isinstance(self.converter, converters.Char):
1544            column.info.meta['_votable_string_dtype'] = 'char'
1545        elif isinstance(self.converter, converters.UnicodeChar):
1546            column.info.meta['_votable_string_dtype'] = 'unicodeChar'
1547
1548    @classmethod
1549    def from_table_column(cls, votable, column):
1550        """
1551        Restores a `Field` instance from a given
1552        `astropy.table.Column` instance.
1553        """
1554        kwargs = {}
1555        meta = column.info.meta
1556        if meta:
1557            for key in ['ucd', 'width', 'precision', 'utype', 'xtype']:
1558                val = meta.get(key, None)
1559                if val is not None:
1560                    kwargs[key] = val
1561        # TODO: Use the unit framework when available
1562        if column.info.unit is not None:
1563            kwargs['unit'] = column.info.unit
1564        kwargs['name'] = column.info.name
1565        result = converters.table_column_to_votable_datatype(column)
1566        kwargs.update(result)
1567
1568        field = cls(votable, **kwargs)
1569
1570        if column.info.description is not None:
1571            field.description = column.info.description
1572        field.values.from_table_column(column)
1573        if meta and 'links' in meta:
1574            for link in meta['links']:
1575                field.links.append(Link.from_table_column(link))
1576
1577        # TODO: Parse format into precision and width
1578        return field
1579
1580
1581class Param(Field):
1582    """
1583    PARAM_ element: constant-valued columns in the data.
1584
1585    :class:`Param` objects are a subclass of :class:`Field`, and have
1586    all of its methods and members.  Additionally, it defines :attr:`value`.
1587    """
1588    _attr_list_11 = Field._attr_list_11 + ['value']
1589    _attr_list_12 = Field._attr_list_12 + ['value']
1590    _element_name = 'PARAM'
1591
1592    def __init__(self, votable, ID=None, name=None, value=None, datatype=None,
1593                 arraysize=None, ucd=None, unit=None, width=None,
1594                 precision=None, utype=None, type=None, id=None, config=None,
1595                 pos=None, **extra):
1596        self._value = value
1597        Field.__init__(self, votable, ID=ID, name=name, datatype=datatype,
1598                       arraysize=arraysize, ucd=ucd, unit=unit,
1599                       precision=precision, utype=utype, type=type,
1600                       id=id, config=config, pos=pos, **extra)
1601
1602    @property
1603    def value(self):
1604        """
1605        [*required*] The constant value of the parameter.  Its type is
1606        determined by the :attr:`~Field.datatype` member.
1607        """
1608        return self._value
1609
1610    @value.setter
1611    def value(self, value):
1612        if value is None:
1613            value = ""
1614        if isinstance(value, str):
1615            self._value = self.converter.parse(
1616                value, self._config, self._pos)[0]
1617        else:
1618            self._value = value
1619
1620    def _setup(self, config, pos):
1621        Field._setup(self, config, pos)
1622        self.value = self._value
1623
1624    def to_xml(self, w, **kwargs):
1625        tmp_value = self._value
1626        self._value = self.converter.output(tmp_value, False)
1627        # We must always have a value
1628        if self._value is None:
1629            self._value = ""
1630        Field.to_xml(self, w, **kwargs)
1631        self._value = tmp_value
1632
1633
1634class CooSys(SimpleElement):
1635    """
1636    COOSYS_ element: defines a coordinate system.
1637
1638    The keyword arguments correspond to setting members of the same
1639    name, documented below.
1640    """
1641    _attr_list = ['ID', 'equinox', 'epoch', 'system']
1642    _element_name = 'COOSYS'
1643
1644    def __init__(self, ID=None, equinox=None, epoch=None, system=None, id=None,
1645                 config=None, pos=None, **extra):
1646        if config is None:
1647            config = {}
1648        self._config = config
1649        self._pos = pos
1650
1651        # COOSYS was deprecated in 1.2 but then re-instated in 1.3
1652        if (config.get('version_1_2_or_later') and
1653                not config.get('version_1_3_or_later')):
1654            warn_or_raise(W27, W27, (), config, pos)
1655
1656        SimpleElement.__init__(self)
1657
1658        self.ID = resolve_id(ID, id, config, pos)
1659        self.equinox = equinox
1660        self.epoch = epoch
1661        self.system = system
1662
1663        warn_unknown_attrs('COOSYS', extra.keys(), config, pos)
1664
1665    @property
1666    def ID(self):
1667        """
1668        [*required*] The XML ID of the COOSYS_ element, used for
1669        cross-referencing.  May be `None` or a string conforming to
1670        XML ID_ syntax.
1671        """
1672        return self._ID
1673
1674    @ID.setter
1675    def ID(self, ID):
1676        if self._config.get('version_1_1_or_later'):
1677            if ID is None:
1678                vo_raise(E15, (), self._config, self._pos)
1679        xmlutil.check_id(ID, 'ID', self._config, self._pos)
1680        self._ID = ID
1681
1682    @property
1683    def system(self):
1684        """
1685        Specifies the type of coordinate system.  Valid choices are:
1686
1687          'eq_FK4', 'eq_FK5', 'ICRS', 'ecl_FK4', 'ecl_FK5', 'galactic',
1688          'supergalactic', 'xy', 'barycentric', or 'geo_app'
1689        """
1690        return self._system
1691
1692    @system.setter
1693    def system(self, system):
1694        if system not in ('eq_FK4', 'eq_FK5', 'ICRS', 'ecl_FK4', 'ecl_FK5',
1695                          'galactic', 'supergalactic', 'xy', 'barycentric',
1696                          'geo_app'):
1697            warn_or_raise(E16, E16, system, self._config, self._pos)
1698        self._system = system
1699
1700    @system.deleter
1701    def system(self):
1702        self._system = None
1703
1704    @property
1705    def equinox(self):
1706        """
1707        A parameter required to fix the equatorial or ecliptic systems
1708        (as e.g. "J2000" as the default "eq_FK5" or "B1950" as the
1709        default "eq_FK4").
1710        """
1711        return self._equinox
1712
1713    @equinox.setter
1714    def equinox(self, equinox):
1715        check_astroyear(equinox, 'equinox', self._config, self._pos)
1716        self._equinox = equinox
1717
1718    @equinox.deleter
1719    def equinox(self):
1720        self._equinox = None
1721
1722    @property
1723    def epoch(self):
1724        """
1725        Specifies the epoch of the positions.  It must be a string
1726        specifying an astronomical year.
1727        """
1728        return self._epoch
1729
1730    @epoch.setter
1731    def epoch(self, epoch):
1732        check_astroyear(epoch, 'epoch', self._config, self._pos)
1733        self._epoch = epoch
1734
1735    @epoch.deleter
1736    def epoch(self):
1737        self._epoch = None
1738
1739
1740class TimeSys(SimpleElement):
1741    """
1742    TIMESYS_ element: defines a time system.
1743
1744    The keyword arguments correspond to setting members of the same
1745    name, documented below.
1746    """
1747    _attr_list = ['ID', 'timeorigin', 'timescale', 'refposition']
1748    _element_name = 'TIMESYS'
1749
1750    def __init__(self, ID=None, timeorigin=None, timescale=None, refposition=None, id=None,
1751                 config=None, pos=None, **extra):
1752        if config is None:
1753            config = {}
1754        self._config = config
1755        self._pos = pos
1756
1757        # TIMESYS is supported starting in version 1.4
1758        if not config['version_1_4_or_later']:
1759            warn_or_raise(
1760                W54, W54, config['version'], config, pos)
1761
1762        SimpleElement.__init__(self)
1763
1764        self.ID = resolve_id(ID, id, config, pos)
1765        self.timeorigin = timeorigin
1766        self.timescale = timescale
1767        self.refposition = refposition
1768
1769        warn_unknown_attrs('TIMESYS', extra.keys(), config, pos,
1770                           ['ID', 'timeorigin', 'timescale', 'refposition'])
1771
1772    @property
1773    def ID(self):
1774        """
1775        [*required*] The XML ID of the TIMESYS_ element, used for
1776        cross-referencing.  Must be a string conforming to
1777        XML ID_ syntax.
1778        """
1779        return self._ID
1780
1781    @ID.setter
1782    def ID(self, ID):
1783        if ID is None:
1784            vo_raise(E22, (), self._config, self._pos)
1785        xmlutil.check_id(ID, 'ID', self._config, self._pos)
1786        self._ID = ID
1787
1788    @property
1789    def timeorigin(self):
1790        """
1791        Specifies the time origin of the time coordinate,
1792        given as a Julian Date for the the time scale and
1793        reference point defined. It is usually given as a
1794        floating point literal; for convenience, the magic
1795        strings "MJD-origin" (standing for 2400000.5) and
1796        "JD-origin" (standing for 0) are also allowed.
1797
1798        The timeorigin attribute MUST be given unless the
1799        time’s representation contains a year of a calendar
1800        era, in which case it MUST NOT be present. In VOTables,
1801        these representations currently are Gregorian calendar
1802        years with xtype="timestamp", or years in the Julian
1803        or Besselian calendar when a column has yr, a, or Ba as
1804        its unit and no time origin is given.
1805        """
1806        return self._timeorigin
1807
1808    @timeorigin.setter
1809    def timeorigin(self, timeorigin):
1810        if (timeorigin is not None and
1811                timeorigin != 'MJD-origin' and timeorigin != 'JD-origin'):
1812            try:
1813                timeorigin = float(timeorigin)
1814            except ValueError:
1815                warn_or_raise(E23, E23, timeorigin, self._config, self._pos)
1816        self._timeorigin = timeorigin
1817
1818    @timeorigin.deleter
1819    def timeorigin(self):
1820        self._timeorigin = None
1821
1822    @property
1823    def timescale(self):
1824        """
1825        [*required*] String specifying the time scale used. Values
1826        should be taken from the IVOA timescale vocabulary (documented
1827        at http://www.ivoa.net/rdf/timescale).
1828        """
1829        return self._timescale
1830
1831    @timescale.setter
1832    def timescale(self, timescale):
1833        self._timescale = timescale
1834
1835    @timescale.deleter
1836    def timescale(self):
1837        self._timescale = None
1838
1839    @property
1840    def refposition(self):
1841        """
1842        [*required*] String specifying the reference position. Values
1843        should be taken from the IVOA refposition vocabulary (documented
1844        at http://www.ivoa.net/rdf/refposition).
1845        """
1846        return self._refposition
1847
1848    @refposition.setter
1849    def refposition(self, refposition):
1850        self._refposition = refposition
1851
1852    @refposition.deleter
1853    def refposition(self):
1854        self._refposition = None
1855
1856
1857class FieldRef(SimpleElement, _UtypeProperty, _UcdProperty):
1858    """
1859    FIELDref_ element: used inside of GROUP_ elements to refer to remote FIELD_ elements.
1860    """
1861    _attr_list_11 = ['ref']
1862    _attr_list_12 = _attr_list_11 + ['ucd', 'utype']
1863    _element_name = "FIELDref"
1864    _utype_in_v1_2 = True
1865    _ucd_in_v1_2 = True
1866
1867    def __init__(self, table, ref, ucd=None, utype=None, config=None, pos=None,
1868                 **extra):
1869        """
1870        *table* is the :class:`Table` object that this :class:`FieldRef`
1871        is a member of.
1872
1873        *ref* is the ID to reference a :class:`Field` object defined
1874        elsewhere.
1875        """
1876        if config is None:
1877            config = {}
1878        self._config = config
1879        self._pos = pos
1880
1881        SimpleElement.__init__(self)
1882        self._table = table
1883        self.ref = ref
1884        self.ucd = ucd
1885        self.utype = utype
1886
1887        if config.get('version_1_2_or_later'):
1888            self._attr_list = self._attr_list_12
1889        else:
1890            self._attr_list = self._attr_list_11
1891            if ucd is not None:
1892                warn_unknown_attrs(self._element_name, ['ucd'], config, pos)
1893            if utype is not None:
1894                warn_unknown_attrs(self._element_name, ['utype'], config, pos)
1895
1896    @property
1897    def ref(self):
1898        """The ID_ of the FIELD_ that this FIELDref_ references."""
1899        return self._ref
1900
1901    @ref.setter
1902    def ref(self, ref):
1903        xmlutil.check_id(ref, 'ref', self._config, self._pos)
1904        self._ref = ref
1905
1906    @ref.deleter
1907    def ref(self):
1908        self._ref = None
1909
1910    def get_ref(self):
1911        """
1912        Lookup the :class:`Field` instance that this :class:`FieldRef`
1913        references.
1914        """
1915        for field in self._table._votable.iter_fields_and_params():
1916            if isinstance(field, Field) and field.ID == self.ref:
1917                return field
1918        vo_raise(
1919            f"No field named '{self.ref}'",
1920            self._config, self._pos, KeyError)
1921
1922
1923class ParamRef(SimpleElement, _UtypeProperty, _UcdProperty):
1924    """
1925    PARAMref_ element: used inside of GROUP_ elements to refer to remote PARAM_ elements.
1926
1927    The keyword arguments correspond to setting members of the same
1928    name, documented below.
1929
1930    It contains the following publicly-accessible members:
1931
1932      *ref*: An XML ID referring to a <PARAM> element.
1933    """
1934    _attr_list_11 = ['ref']
1935    _attr_list_12 = _attr_list_11 + ['ucd', 'utype']
1936    _element_name = "PARAMref"
1937    _utype_in_v1_2 = True
1938    _ucd_in_v1_2 = True
1939
1940    def __init__(self, table, ref, ucd=None, utype=None, config=None, pos=None):
1941        if config is None:
1942            config = {}
1943
1944        self._config = config
1945        self._pos = pos
1946
1947        Element.__init__(self)
1948        self._table = table
1949        self.ref = ref
1950        self.ucd = ucd
1951        self.utype = utype
1952
1953        if config.get('version_1_2_or_later'):
1954            self._attr_list = self._attr_list_12
1955        else:
1956            self._attr_list = self._attr_list_11
1957            if ucd is not None:
1958                warn_unknown_attrs(self._element_name, ['ucd'], config, pos)
1959            if utype is not None:
1960                warn_unknown_attrs(self._element_name, ['utype'], config, pos)
1961
1962    @property
1963    def ref(self):
1964        """The ID_ of the PARAM_ that this PARAMref_ references."""
1965        return self._ref
1966
1967    @ref.setter
1968    def ref(self, ref):
1969        xmlutil.check_id(ref, 'ref', self._config, self._pos)
1970        self._ref = ref
1971
1972    @ref.deleter
1973    def ref(self):
1974        self._ref = None
1975
1976    def get_ref(self):
1977        """
1978        Lookup the :class:`Param` instance that this :class:``PARAMref``
1979        references.
1980        """
1981        for param in self._table._votable.iter_fields_and_params():
1982            if isinstance(param, Param) and param.ID == self.ref:
1983                return param
1984        vo_raise(
1985            f"No params named '{self.ref}'",
1986            self._config, self._pos, KeyError)
1987
1988
1989class Group(Element, _IDProperty, _NameProperty, _UtypeProperty,
1990            _UcdProperty, _DescriptionProperty):
1991    """
1992    GROUP_ element: groups FIELD_ and PARAM_ elements.
1993
1994    This information is currently ignored by the vo package---that is
1995    the columns in the recarray are always flat---but the grouping
1996    information is stored so that it can be written out again to the
1997    XML file.
1998
1999    The keyword arguments correspond to setting members of the same
2000    name, documented below.
2001    """
2002
2003    def __init__(self, table, ID=None, name=None, ref=None, ucd=None,
2004                 utype=None, id=None, config=None, pos=None, **extra):
2005        if config is None:
2006            config = {}
2007        self._config = config
2008        self._pos = pos
2009
2010        Element.__init__(self)
2011        self._table = table
2012
2013        self.ID = (resolve_id(ID, id, config, pos)
2014                            or xmlutil.fix_id(name, config, pos))
2015        self.name = name
2016        self.ref = ref
2017        self.ucd = ucd
2018        self.utype = utype
2019        self.description = None
2020
2021        self._entries = HomogeneousList(
2022            (FieldRef, ParamRef, Group, Param))
2023
2024        warn_unknown_attrs('GROUP', extra.keys(), config, pos)
2025
2026    def __repr__(self):
2027        return f'<GROUP>... {len(self._entries)} entries ...</GROUP>'
2028
2029    @property
2030    def ref(self):
2031        """
2032        Currently ignored, as it's not clear from the spec how this is
2033        meant to work.
2034        """
2035        return self._ref
2036
2037    @ref.setter
2038    def ref(self, ref):
2039        xmlutil.check_id(ref, 'ref', self._config, self._pos)
2040        self._ref = ref
2041
2042    @ref.deleter
2043    def ref(self):
2044        self._ref = None
2045
2046    @property
2047    def entries(self):
2048        """
2049        [read-only] A list of members of the GROUP_.  This list may
2050        only contain objects of type :class:`Param`, :class:`Group`,
2051        :class:`ParamRef` and :class:`FieldRef`.
2052        """
2053        return self._entries
2054
2055    def _add_fieldref(self, iterator, tag, data, config, pos):
2056        fieldref = FieldRef(self._table, config=config, pos=pos, **data)
2057        self.entries.append(fieldref)
2058
2059    def _add_paramref(self, iterator, tag, data, config, pos):
2060        paramref = ParamRef(self._table, config=config, pos=pos, **data)
2061        self.entries.append(paramref)
2062
2063    def _add_param(self, iterator, tag, data, config, pos):
2064        if isinstance(self._table, VOTableFile):
2065            votable = self._table
2066        else:
2067            votable = self._table._votable
2068        param = Param(votable, config=config, pos=pos, **data)
2069        self.entries.append(param)
2070        param.parse(iterator, config)
2071
2072    def _add_group(self, iterator, tag, data, config, pos):
2073        group = Group(self._table, config=config, pos=pos, **data)
2074        self.entries.append(group)
2075        group.parse(iterator, config)
2076
2077    def parse(self, iterator, config):
2078        tag_mapping = {
2079            'FIELDref': self._add_fieldref,
2080            'PARAMref': self._add_paramref,
2081            'PARAM': self._add_param,
2082            'GROUP': self._add_group,
2083            'DESCRIPTION': self._ignore_add}
2084
2085        for start, tag, data, pos in iterator:
2086            if start:
2087                tag_mapping.get(tag, self._add_unknown_tag)(
2088                    iterator, tag, data, config, pos)
2089            else:
2090                if tag == 'DESCRIPTION':
2091                    if self.description is not None:
2092                        warn_or_raise(W17, W17, 'GROUP', config, pos)
2093                    self.description = data or None
2094                elif tag == 'GROUP':
2095                    break
2096        return self
2097
2098    def to_xml(self, w, **kwargs):
2099        with w.tag(
2100            'GROUP',
2101            attrib=w.object_attrs(
2102                self, ['ID', 'name', 'ref', 'ucd', 'utype'])):
2103            if self.description is not None:
2104                w.element("DESCRIPTION", self.description, wrap=True)
2105            for entry in self.entries:
2106                entry.to_xml(w, **kwargs)
2107
2108    def iter_fields_and_params(self):
2109        """
2110        Recursively iterate over all :class:`Param` elements in this
2111        :class:`Group`.
2112        """
2113        for entry in self.entries:
2114            if isinstance(entry, Param):
2115                yield entry
2116            elif isinstance(entry, Group):
2117                for field in entry.iter_fields_and_params():
2118                    yield field
2119
2120    def iter_groups(self):
2121        """
2122        Recursively iterate over all sub-:class:`Group` instances in
2123        this :class:`Group`.
2124        """
2125        for entry in self.entries:
2126            if isinstance(entry, Group):
2127                yield entry
2128                for group in entry.iter_groups():
2129                    yield group
2130
2131
2132class Table(Element, _IDProperty, _NameProperty, _UcdProperty,
2133            _DescriptionProperty):
2134    """
2135    TABLE_ element: optionally contains data.
2136
2137    It contains the following publicly-accessible and mutable
2138    attribute:
2139
2140        *array*: A Numpy masked array of the data itself, where each
2141        row is a row of votable data, and columns are named and typed
2142        based on the <FIELD> elements of the table.  The mask is
2143        parallel to the data array, except for variable-length fields.
2144        For those fields, the numpy array's column type is "object"
2145        (``"O"``), and another masked array is stored there.
2146
2147    If the Table contains no data, (for example, its enclosing
2148    :class:`Resource` has :attr:`~Resource.type` == 'meta') *array*
2149    will have zero-length.
2150
2151    The keyword arguments correspond to setting members of the same
2152    name, documented below.
2153    """
2154
2155    def __init__(self, votable, ID=None, name=None, ref=None, ucd=None,
2156                 utype=None, nrows=None, id=None, config=None, pos=None,
2157                 **extra):
2158        if config is None:
2159            config = {}
2160        self._config = config
2161        self._pos = pos
2162        self._empty = False
2163
2164        Element.__init__(self)
2165        self._votable = votable
2166
2167        self.ID = (resolve_id(ID, id, config, pos)
2168                   or xmlutil.fix_id(name, config, pos))
2169        self.name = name
2170        xmlutil.check_id(ref, 'ref', config, pos)
2171        self._ref = ref
2172        self.ucd = ucd
2173        self.utype = utype
2174        if nrows is not None:
2175            nrows = int(nrows)
2176            if nrows < 0:
2177                raise ValueError("'nrows' cannot be negative.")
2178        self._nrows = nrows
2179        self.description = None
2180        self.format = 'tabledata'
2181
2182        self._fields = HomogeneousList(Field)
2183        self._params = HomogeneousList(Param)
2184        self._groups = HomogeneousList(Group)
2185        self._links = HomogeneousList(Link)
2186        self._infos = HomogeneousList(Info)
2187
2188        self.array = ma.array([])
2189
2190        warn_unknown_attrs('TABLE', extra.keys(), config, pos)
2191
2192    def __repr__(self):
2193        return repr(self.to_table())
2194
2195    def __bytes__(self):
2196        return bytes(self.to_table())
2197
2198    def __str__(self):
2199        return str(self.to_table())
2200
2201    @property
2202    def ref(self):
2203        return self._ref
2204
2205    @ref.setter
2206    def ref(self, ref):
2207        """
2208        Refer to another TABLE, previously defined, by the *ref* ID_
2209        for all metadata (FIELD_, PARAM_ etc.) information.
2210        """
2211        # When the ref changes, we want to verify that it will work
2212        # by actually going and looking for the referenced table.
2213        # If found, set a bunch of properties in this table based
2214        # on the other one.
2215        xmlutil.check_id(ref, 'ref', self._config, self._pos)
2216        if ref is not None:
2217            try:
2218                table = self._votable.get_table_by_id(ref, before=self)
2219            except KeyError:
2220                warn_or_raise(
2221                    W43, W43, ('TABLE', self.ref), self._config, self._pos)
2222                ref = None
2223            else:
2224                self._fields = table.fields
2225                self._params = table.params
2226                self._groups = table.groups
2227                self._links = table.links
2228        else:
2229            del self._fields[:]
2230            del self._params[:]
2231            del self._groups[:]
2232            del self._links[:]
2233        self._ref = ref
2234
2235    @ref.deleter
2236    def ref(self):
2237        self._ref = None
2238
2239    @property
2240    def format(self):
2241        """
2242        [*required*] The serialization format of the table.  Must be
2243        one of:
2244
2245          'tabledata' (TABLEDATA_), 'binary' (BINARY_), 'binary2' (BINARY2_)
2246          'fits' (FITS_).
2247
2248        Note that the 'fits' format, since it requires an external
2249        file, can not be written out.  Any file read in with 'fits'
2250        format will be read out, by default, in 'tabledata' format.
2251
2252        See :ref:`astropy:votable-serialization`.
2253        """
2254        return self._format
2255
2256    @format.setter
2257    def format(self, format):
2258        format = format.lower()
2259        if format == 'fits':
2260            vo_raise("fits format can not be written out, only read.",
2261                     self._config, self._pos, NotImplementedError)
2262        if format == 'binary2':
2263            if not self._config['version_1_3_or_later']:
2264                vo_raise(
2265                    "binary2 only supported in votable 1.3 or later",
2266                    self._config, self._pos)
2267        elif format not in ('tabledata', 'binary'):
2268            vo_raise(f"Invalid format '{format}'",
2269                     self._config, self._pos)
2270        self._format = format
2271
2272    @property
2273    def nrows(self):
2274        """
2275        [*immutable*] The number of rows in the table, as specified in
2276        the XML file.
2277        """
2278        return self._nrows
2279
2280    @property
2281    def fields(self):
2282        """
2283        A list of :class:`Field` objects describing the types of each
2284        of the data columns.
2285        """
2286        return self._fields
2287
2288    @property
2289    def params(self):
2290        """
2291        A list of parameters (constant-valued columns) for the
2292        table.  Must contain only :class:`Param` objects.
2293        """
2294        return self._params
2295
2296    @property
2297    def groups(self):
2298        """
2299        A list of :class:`Group` objects describing how the columns
2300        and parameters are grouped.  Currently this information is
2301        only kept around for round-tripping and informational
2302        purposes.
2303        """
2304        return self._groups
2305
2306    @property
2307    def links(self):
2308        """
2309        A list of :class:`Link` objects (pointers to other documents
2310        or servers through a URI) for the table.
2311        """
2312        return self._links
2313
2314    @property
2315    def infos(self):
2316        """
2317        A list of :class:`Info` objects for the table.  Allows for
2318        post-operational diagnostics.
2319        """
2320        return self._infos
2321
2322    def is_empty(self):
2323        """
2324        Returns True if this table doesn't contain any real data
2325        because it was skipped over by the parser (through use of the
2326        ``table_number`` kwarg).
2327        """
2328        return self._empty
2329
2330    def create_arrays(self, nrows=0, config=None):
2331        """
2332        Create a new array to hold the data based on the current set
2333        of fields, and store them in the *array* and member variable.
2334        Any data in the existing array will be lost.
2335
2336        *nrows*, if provided, is the number of rows to allocate.
2337        """
2338        if nrows is None:
2339            nrows = 0
2340
2341        fields = self.fields
2342
2343        if len(fields) == 0:
2344            array = np.recarray((nrows,), dtype='O')
2345            mask = np.zeros((nrows,), dtype='b')
2346        else:
2347            # for field in fields: field._setup(config)
2348            Field.uniqify_names(fields)
2349
2350            dtype = []
2351            for x in fields:
2352                if x._unique_name == x.ID:
2353                    id = x.ID
2354                else:
2355                    id = (x._unique_name, x.ID)
2356                dtype.append((id, x.converter.format))
2357
2358            array = np.recarray((nrows,), dtype=np.dtype(dtype))
2359            descr_mask = []
2360            for d in array.dtype.descr:
2361                new_type = (d[1][1] == 'O' and 'O') or 'bool'
2362                if len(d) == 2:
2363                    descr_mask.append((d[0], new_type))
2364                elif len(d) == 3:
2365                    descr_mask.append((d[0], new_type, d[2]))
2366            mask = np.zeros((nrows,), dtype=descr_mask)
2367
2368        self.array = ma.array(array, mask=mask)
2369
2370    def _resize_strategy(self, size):
2371        """
2372        Return a new (larger) size based on size, used for
2373        reallocating an array when it fills up.  This is in its own
2374        function so the resizing strategy can be easily replaced.
2375        """
2376        # Once we go beyond 0, make a big step -- after that use a
2377        # factor of 1.5 to help keep memory usage compact
2378        if size == 0:
2379            return 512
2380        return int(np.ceil(size * RESIZE_AMOUNT))
2381
2382    def _add_field(self, iterator, tag, data, config, pos):
2383        field = Field(self._votable, config=config, pos=pos, **data)
2384        self.fields.append(field)
2385        field.parse(iterator, config)
2386
2387    def _add_param(self, iterator, tag, data, config, pos):
2388        param = Param(self._votable, config=config, pos=pos, **data)
2389        self.params.append(param)
2390        param.parse(iterator, config)
2391
2392    def _add_group(self, iterator, tag, data, config, pos):
2393        group = Group(self, config=config, pos=pos, **data)
2394        self.groups.append(group)
2395        group.parse(iterator, config)
2396
2397    def _add_link(self, iterator, tag, data, config, pos):
2398        link = Link(config=config, pos=pos, **data)
2399        self.links.append(link)
2400        link.parse(iterator, config)
2401
2402    def _add_info(self, iterator, tag, data, config, pos):
2403        if not config.get('version_1_2_or_later'):
2404            warn_or_raise(W26, W26, ('INFO', 'TABLE', '1.2'), config, pos)
2405        info = Info(config=config, pos=pos, **data)
2406        self.infos.append(info)
2407        info.parse(iterator, config)
2408
2409    def parse(self, iterator, config):
2410        columns = config.get('columns')
2411
2412        # If we've requested to read in only a specific table, skip
2413        # all others
2414        table_number = config.get('table_number')
2415        current_table_number = config.get('_current_table_number')
2416        skip_table = False
2417        if current_table_number is not None:
2418            config['_current_table_number'] += 1
2419            if (table_number is not None and
2420                table_number != current_table_number):
2421                skip_table = True
2422                self._empty = True
2423
2424        table_id = config.get('table_id')
2425        if table_id is not None:
2426            if table_id != self.ID:
2427                skip_table = True
2428                self._empty = True
2429
2430        if self.ref is not None:
2431            # This table doesn't have its own datatype descriptors, it
2432            # just references those from another table.
2433
2434            # This is to call the property setter to go and get the
2435            # referenced information
2436            self.ref = self.ref
2437
2438            for start, tag, data, pos in iterator:
2439                if start:
2440                    if tag == 'DATA':
2441                        warn_unknown_attrs(
2442                            'DATA', data.keys(), config, pos)
2443                        break
2444                else:
2445                    if tag == 'TABLE':
2446                        return self
2447                    elif tag == 'DESCRIPTION':
2448                        if self.description is not None:
2449                            warn_or_raise(W17, W17, 'RESOURCE', config, pos)
2450                        self.description = data or None
2451        else:
2452            tag_mapping = {
2453                'FIELD': self._add_field,
2454                'PARAM': self._add_param,
2455                'GROUP': self._add_group,
2456                'LINK': self._add_link,
2457                'INFO': self._add_info,
2458                'DESCRIPTION': self._ignore_add}
2459
2460            for start, tag, data, pos in iterator:
2461                if start:
2462                    if tag == 'DATA':
2463                        if len(self.fields) == 0:
2464                            warn_or_raise(E25, E25, None, config, pos)
2465                        warn_unknown_attrs(
2466                            'DATA', data.keys(), config, pos)
2467                        break
2468
2469                    tag_mapping.get(tag, self._add_unknown_tag)(
2470                        iterator, tag, data, config, pos)
2471                else:
2472                    if tag == 'DESCRIPTION':
2473                        if self.description is not None:
2474                            warn_or_raise(W17, W17, 'RESOURCE', config, pos)
2475                        self.description = data or None
2476                    elif tag == 'TABLE':
2477                        # For error checking purposes
2478                        Field.uniqify_names(self.fields)
2479                        # We still need to create arrays, even if the file
2480                        # contains no DATA section
2481                        self.create_arrays(nrows=0, config=config)
2482                        return self
2483
2484        self.create_arrays(nrows=self._nrows, config=config)
2485        fields = self.fields
2486        names = [x.ID for x in fields]
2487        # Deal with a subset of the columns, if requested.
2488        if not columns:
2489            colnumbers = list(range(len(fields)))
2490        else:
2491            if isinstance(columns, str):
2492                columns = [columns]
2493            columns = np.asarray(columns)
2494            if issubclass(columns.dtype.type, np.integer):
2495                if np.any(columns < 0) or np.any(columns > len(fields)):
2496                    raise ValueError(
2497                        "Some specified column numbers out of range")
2498                colnumbers = columns
2499            elif issubclass(columns.dtype.type, np.character):
2500                try:
2501                    colnumbers = [names.index(x) for x in columns]
2502                except ValueError:
2503                    raise ValueError(
2504                        f"Columns '{columns}' not found in fields list")
2505            else:
2506                raise TypeError("Invalid columns list")
2507
2508        if (not skip_table) and (len(fields) > 0):
2509            for start, tag, data, pos in iterator:
2510                if start:
2511                    if tag == 'TABLEDATA':
2512                        warn_unknown_attrs(
2513                            'TABLEDATA', data.keys(), config, pos)
2514                        self.array = self._parse_tabledata(
2515                            iterator, colnumbers, fields, config)
2516                        break
2517                    elif tag == 'BINARY':
2518                        warn_unknown_attrs(
2519                            'BINARY', data.keys(), config, pos)
2520                        self.array = self._parse_binary(
2521                            1, iterator, colnumbers, fields, config, pos)
2522                        break
2523                    elif tag == 'BINARY2':
2524                        if not config['version_1_3_or_later']:
2525                            warn_or_raise(
2526                                W52, W52, config['version'], config, pos)
2527                        self.array = self._parse_binary(
2528                            2, iterator, colnumbers, fields, config, pos)
2529                        break
2530                    elif tag == 'FITS':
2531                        warn_unknown_attrs(
2532                            'FITS', data.keys(), config, pos, ['extnum'])
2533                        try:
2534                            extnum = int(data.get('extnum', 0))
2535                            if extnum < 0:
2536                                raise ValueError("'extnum' cannot be negative.")
2537                        except ValueError:
2538                            vo_raise(E17, (), config, pos)
2539                        self.array = self._parse_fits(
2540                            iterator, extnum, config)
2541                        break
2542                    else:
2543                        warn_or_raise(W37, W37, tag, config, pos)
2544                        break
2545
2546        for start, tag, data, pos in iterator:
2547            if not start and tag == 'DATA':
2548                break
2549
2550        for start, tag, data, pos in iterator:
2551            if start and tag == 'INFO':
2552                if not config.get('version_1_2_or_later'):
2553                    warn_or_raise(
2554                        W26, W26, ('INFO', 'TABLE', '1.2'), config, pos)
2555                info = Info(config=config, pos=pos, **data)
2556                self.infos.append(info)
2557                info.parse(iterator, config)
2558            elif not start and tag == 'TABLE':
2559                break
2560
2561        return self
2562
2563    def _parse_tabledata(self, iterator, colnumbers, fields, config):
2564        # Since we don't know the number of rows up front, we'll
2565        # reallocate the record array to make room as we go.  This
2566        # prevents the need to scan through the XML twice.  The
2567        # allocation is by factors of 1.5.
2568        invalid = config.get('invalid', 'exception')
2569
2570        # Need to have only one reference so that we can resize the
2571        # array
2572        array = self.array
2573        del self.array
2574
2575        parsers = [field.converter.parse for field in fields]
2576        binparsers = [field.converter.binparse for field in fields]
2577
2578        numrows = 0
2579        alloc_rows = len(array)
2580        colnumbers_bits = [i in colnumbers for i in range(len(fields))]
2581        row_default = [x.converter.default for x in fields]
2582        mask_default = [True] * len(fields)
2583        array_chunk = []
2584        mask_chunk = []
2585        chunk_size = config.get('chunk_size', DEFAULT_CHUNK_SIZE)
2586        for start, tag, data, pos in iterator:
2587            if tag == 'TR':
2588                # Now parse one row
2589                row = row_default[:]
2590                row_mask = mask_default[:]
2591                i = 0
2592                for start, tag, data, pos in iterator:
2593                    if start:
2594                        binary = (data.get('encoding', None) == 'base64')
2595                        warn_unknown_attrs(
2596                            tag, data.keys(), config, pos, ['encoding'])
2597                    else:
2598                        if tag == 'TD':
2599                            if i >= len(fields):
2600                                vo_raise(E20, len(fields), config, pos)
2601
2602                            if colnumbers_bits[i]:
2603                                try:
2604                                    if binary:
2605                                        rawdata = base64.b64decode(
2606                                            data.encode('ascii'))
2607                                        buf = io.BytesIO(rawdata)
2608                                        buf.seek(0)
2609                                        try:
2610                                            value, mask_value = binparsers[i](
2611                                                buf.read)
2612                                        except Exception as e:
2613                                            vo_reraise(
2614                                                e, config, pos,
2615                                                "(in row {:d}, col '{}')".format(
2616                                                    len(array_chunk),
2617                                                    fields[i].ID))
2618                                    else:
2619                                        try:
2620                                            value, mask_value = parsers[i](
2621                                                data, config, pos)
2622                                        except Exception as e:
2623                                            vo_reraise(
2624                                                e, config, pos,
2625                                                "(in row {:d}, col '{}')".format(
2626                                                    len(array_chunk),
2627                                                    fields[i].ID))
2628                                except Exception as e:
2629                                    if invalid == 'exception':
2630                                        vo_reraise(e, config, pos)
2631                                else:
2632                                    row[i] = value
2633                                    row_mask[i] = mask_value
2634                        elif tag == 'TR':
2635                            break
2636                        else:
2637                            self._add_unknown_tag(
2638                                iterator, tag, data, config, pos)
2639                        i += 1
2640
2641                if i < len(fields):
2642                    vo_raise(E21, (i, len(fields)), config, pos)
2643
2644                array_chunk.append(tuple(row))
2645                mask_chunk.append(tuple(row_mask))
2646
2647                if len(array_chunk) == chunk_size:
2648                    while numrows + chunk_size > alloc_rows:
2649                        alloc_rows = self._resize_strategy(alloc_rows)
2650                    if alloc_rows != len(array):
2651                        array = _resize(array, alloc_rows)
2652                    array[numrows:numrows + chunk_size] = array_chunk
2653                    array.mask[numrows:numrows + chunk_size] = mask_chunk
2654                    numrows += chunk_size
2655                    array_chunk = []
2656                    mask_chunk = []
2657
2658            elif not start and tag == 'TABLEDATA':
2659                break
2660
2661        # Now, resize the array to the exact number of rows we need and
2662        # put the last chunk values in there.
2663        alloc_rows = numrows + len(array_chunk)
2664
2665        array = _resize(array, alloc_rows)
2666        array[numrows:] = array_chunk
2667        if alloc_rows != 0:
2668            array.mask[numrows:] = mask_chunk
2669        numrows += len(array_chunk)
2670
2671        if (self.nrows is not None and
2672            self.nrows >= 0 and
2673            self.nrows != numrows):
2674            warn_or_raise(W18, W18, (self.nrows, numrows), config, pos)
2675        self._nrows = numrows
2676
2677        return array
2678
2679    def _get_binary_data_stream(self, iterator, config):
2680        have_local_stream = False
2681        for start, tag, data, pos in iterator:
2682            if tag == 'STREAM':
2683                if start:
2684                    warn_unknown_attrs(
2685                        'STREAM', data.keys(), config, pos,
2686                        ['type', 'href', 'actuate', 'encoding', 'expires',
2687                         'rights'])
2688                    if 'href' not in data:
2689                        have_local_stream = True
2690                        if data.get('encoding', None) != 'base64':
2691                            warn_or_raise(
2692                                W38, W38, data.get('encoding', None),
2693                                config, pos)
2694                    else:
2695                        href = data['href']
2696                        xmlutil.check_anyuri(href, config, pos)
2697                        encoding = data.get('encoding', None)
2698                else:
2699                    buffer = data
2700                    break
2701
2702        if have_local_stream:
2703            buffer = base64.b64decode(buffer.encode('ascii'))
2704            string_io = io.BytesIO(buffer)
2705            string_io.seek(0)
2706            read = string_io.read
2707        else:
2708            if not href.startswith(('http', 'ftp', 'file')):
2709                vo_raise(
2710                    "The vo package only supports remote data through http, " +
2711                    "ftp or file",
2712                    self._config, self._pos, NotImplementedError)
2713            fd = urllib.request.urlopen(href)
2714            if encoding is not None:
2715                if encoding == 'gzip':
2716                    fd = gzip.GzipFile(href, 'rb', fileobj=fd)
2717                elif encoding == 'base64':
2718                    fd = codecs.EncodedFile(fd, 'base64')
2719                else:
2720                    vo_raise(
2721                        f"Unknown encoding type '{encoding}'",
2722                        self._config, self._pos, NotImplementedError)
2723            read = fd.read
2724
2725        def careful_read(length):
2726            result = read(length)
2727            if len(result) != length:
2728                raise EOFError
2729            return result
2730
2731        return careful_read
2732
2733    def _parse_binary(self, mode, iterator, colnumbers, fields, config, pos):
2734        fields = self.fields
2735
2736        careful_read = self._get_binary_data_stream(iterator, config)
2737
2738        # Need to have only one reference so that we can resize the
2739        # array
2740        array = self.array
2741        del self.array
2742
2743        binparsers = [field.converter.binparse for field in fields]
2744
2745        numrows = 0
2746        alloc_rows = len(array)
2747        while True:
2748            # Resize result arrays if necessary
2749            if numrows >= alloc_rows:
2750                alloc_rows = self._resize_strategy(alloc_rows)
2751                array = _resize(array, alloc_rows)
2752
2753            row_data = []
2754            row_mask_data = []
2755
2756            try:
2757                if mode == 2:
2758                    mask_bits = careful_read(int((len(fields) + 7) / 8))
2759                    row_mask_data = list(converters.bitarray_to_bool(
2760                        mask_bits, len(fields)))
2761
2762                    # Ignore the mask for string columns (see issue 8995)
2763                    for i, f in enumerate(fields):
2764                        if row_mask_data[i] and (f.datatype == 'char' or f.datatype == 'unicodeChar'):
2765                            row_mask_data[i] = False
2766
2767                for i, binparse in enumerate(binparsers):
2768                    try:
2769                        value, value_mask = binparse(careful_read)
2770                    except EOFError:
2771                        raise
2772                    except Exception as e:
2773                        vo_reraise(
2774                            e, config, pos, "(in row {:d}, col '{}')".format(
2775                                numrows, fields[i].ID))
2776                    row_data.append(value)
2777                    if mode == 1:
2778                        row_mask_data.append(value_mask)
2779                    else:
2780                        row_mask_data[i] = row_mask_data[i] or value_mask
2781            except EOFError:
2782                break
2783
2784            row = [x.converter.default for x in fields]
2785            row_mask = [False] * len(fields)
2786            for i in colnumbers:
2787                row[i] = row_data[i]
2788                row_mask[i] = row_mask_data[i]
2789
2790            array[numrows] = tuple(row)
2791            array.mask[numrows] = tuple(row_mask)
2792            numrows += 1
2793
2794        array = _resize(array, numrows)
2795
2796        return array
2797
2798    def _parse_fits(self, iterator, extnum, config):
2799        for start, tag, data, pos in iterator:
2800            if tag == 'STREAM':
2801                if start:
2802                    warn_unknown_attrs(
2803                        'STREAM', data.keys(), config, pos,
2804                        ['type', 'href', 'actuate', 'encoding', 'expires',
2805                         'rights'])
2806                    href = data['href']
2807                    encoding = data.get('encoding', None)
2808                else:
2809                    break
2810
2811        if not href.startswith(('http', 'ftp', 'file')):
2812            vo_raise(
2813                "The vo package only supports remote data through http, "
2814                "ftp or file",
2815                self._config, self._pos, NotImplementedError)
2816
2817        fd = urllib.request.urlopen(href)
2818        if encoding is not None:
2819            if encoding == 'gzip':
2820                fd = gzip.GzipFile(href, 'r', fileobj=fd)
2821            elif encoding == 'base64':
2822                fd = codecs.EncodedFile(fd, 'base64')
2823            else:
2824                vo_raise(
2825                    f"Unknown encoding type '{encoding}'",
2826                    self._config, self._pos, NotImplementedError)
2827
2828        hdulist = fits.open(fd)
2829
2830        array = hdulist[int(extnum)].data
2831        if array.dtype != self.array.dtype:
2832            warn_or_raise(W19, W19, (), self._config, self._pos)
2833
2834        return array
2835
2836    def to_xml(self, w, **kwargs):
2837        specified_format = kwargs.get('tabledata_format')
2838        if specified_format is not None:
2839            format = specified_format
2840        else:
2841            format = self.format
2842        if format == 'fits':
2843            format = 'tabledata'
2844
2845        with w.tag(
2846            'TABLE',
2847            attrib=w.object_attrs(
2848                self,
2849                ('ID', 'name', 'ref', 'ucd', 'utype', 'nrows'))):
2850
2851            if self.description is not None:
2852                w.element("DESCRIPTION", self.description, wrap=True)
2853
2854            for element_set in (self.fields, self.params):
2855                for element in element_set:
2856                    element._setup({}, None)
2857
2858            if self.ref is None:
2859                for element_set in (self.fields, self.params, self.groups,
2860                                    self.links):
2861                    for element in element_set:
2862                        element.to_xml(w, **kwargs)
2863            elif kwargs['version_1_2_or_later']:
2864                index = list(self._votable.iter_tables()).index(self)
2865                group = Group(self, ID=f"_g{index}")
2866                group.to_xml(w, **kwargs)
2867
2868            if len(self.array):
2869                with w.tag('DATA'):
2870                    if format == 'tabledata':
2871                        self._write_tabledata(w, **kwargs)
2872                    elif format == 'binary':
2873                        self._write_binary(1, w, **kwargs)
2874                    elif format == 'binary2':
2875                        self._write_binary(2, w, **kwargs)
2876
2877            if kwargs['version_1_2_or_later']:
2878                for element in self._infos:
2879                    element.to_xml(w, **kwargs)
2880
2881    def _write_tabledata(self, w, **kwargs):
2882        fields = self.fields
2883        array = self.array
2884
2885        with w.tag('TABLEDATA'):
2886            w._flush()
2887            if (_has_c_tabledata_writer and
2888                not kwargs.get('_debug_python_based_parser')):
2889                supports_empty_values = [
2890                    field.converter.supports_empty_values(kwargs)
2891                    for field in fields]
2892                fields = [field.converter.output for field in fields]
2893                indent = len(w._tags) - 1
2894                tablewriter.write_tabledata(
2895                    w.write, array.data, array.mask, fields,
2896                    supports_empty_values, indent, 1 << 8)
2897            else:
2898                write = w.write
2899                indent_spaces = w.get_indentation_spaces()
2900                tr_start = indent_spaces + "<TR>\n"
2901                tr_end = indent_spaces + "</TR>\n"
2902                td = indent_spaces + " <TD>{}</TD>\n"
2903                td_empty = indent_spaces + " <TD/>\n"
2904                fields = [(i, field.converter.output,
2905                           field.converter.supports_empty_values(kwargs))
2906                          for i, field in enumerate(fields)]
2907                for row in range(len(array)):
2908                    write(tr_start)
2909                    array_row = array.data[row]
2910                    mask_row = array.mask[row]
2911                    for i, output, supports_empty_values in fields:
2912                        data = array_row[i]
2913                        masked = mask_row[i]
2914                        if supports_empty_values and np.all(masked):
2915                            write(td_empty)
2916                        else:
2917                            try:
2918                                val = output(data, masked)
2919                            except Exception as e:
2920                                vo_reraise(
2921                                    e,
2922                                    additional="(in row {:d}, col '{}')".format(
2923                                        row, self.fields[i].ID))
2924                            if len(val):
2925                                write(td.format(val))
2926                            else:
2927                                write(td_empty)
2928                    write(tr_end)
2929
2930    def _write_binary(self, mode, w, **kwargs):
2931        fields = self.fields
2932        array = self.array
2933        if mode == 1:
2934            tag_name = 'BINARY'
2935        else:
2936            tag_name = 'BINARY2'
2937
2938        with w.tag(tag_name):
2939            with w.tag('STREAM', encoding='base64'):
2940                fields_basic = [(i, field.converter.binoutput)
2941                                for (i, field) in enumerate(fields)]
2942
2943                data = io.BytesIO()
2944                for row in range(len(array)):
2945                    array_row = array.data[row]
2946                    array_mask = array.mask[row]
2947
2948                    if mode == 2:
2949                        flattened = np.array([np.all(x) for x in array_mask])
2950                        data.write(converters.bool_to_bitarray(flattened))
2951
2952                    for i, converter in fields_basic:
2953                        try:
2954                            chunk = converter(array_row[i], array_mask[i])
2955                            assert type(chunk) == bytes
2956                        except Exception as e:
2957                            vo_reraise(
2958                                e, additional=f"(in row {row:d}, col '{fields[i].ID}')")
2959                        data.write(chunk)
2960
2961                w._flush()
2962                w.write(base64.b64encode(data.getvalue()).decode('ascii'))
2963
2964    def to_table(self, use_names_over_ids=False):
2965        """
2966        Convert this VO Table to an `astropy.table.Table` instance.
2967
2968        Parameters
2969        ----------
2970        use_names_over_ids : bool, optional
2971           When `True` use the ``name`` attributes of columns as the
2972           names of columns in the `astropy.table.Table` instance.
2973           Since names are not guaranteed to be unique, this may cause
2974           some columns to be renamed by appending numbers to the end.
2975           Otherwise (default), use the ID attributes as the column
2976           names.
2977
2978        .. warning::
2979           Variable-length array fields may not be restored
2980           identically when round-tripping through the
2981           `astropy.table.Table` instance.
2982        """
2983        from astropy.table import Table
2984
2985        meta = {}
2986        for key in ['ID', 'name', 'ref', 'ucd', 'utype', 'description']:
2987            val = getattr(self, key, None)
2988            if val is not None:
2989                meta[key] = val
2990
2991        if use_names_over_ids:
2992            names = [field.name for field in self.fields]
2993            unique_names = []
2994            for i, name in enumerate(names):
2995                new_name = name
2996                i = 2
2997                while new_name in unique_names:
2998                    new_name = f'{name}{i}'
2999                    i += 1
3000                unique_names.append(new_name)
3001            names = unique_names
3002        else:
3003            names = [field.ID for field in self.fields]
3004
3005        table = Table(self.array, names=names, meta=meta)
3006
3007        for name, field in zip(names, self.fields):
3008            column = table[name]
3009            field.to_table_column(column)
3010
3011        return table
3012
3013    @classmethod
3014    def from_table(cls, votable, table):
3015        """
3016        Create a `Table` instance from a given `astropy.table.Table`
3017        instance.
3018        """
3019        kwargs = {}
3020        for key in ['ID', 'name', 'ref', 'ucd', 'utype']:
3021            val = table.meta.get(key)
3022            if val is not None:
3023                kwargs[key] = val
3024        new_table = cls(votable, **kwargs)
3025        if 'description' in table.meta:
3026            new_table.description = table.meta['description']
3027
3028        for colname in table.colnames:
3029            column = table[colname]
3030            new_table.fields.append(Field.from_table_column(votable, column))
3031
3032        if table.mask is None:
3033            new_table.array = ma.array(np.asarray(table))
3034        else:
3035            new_table.array = ma.array(np.asarray(table),
3036                                       mask=np.asarray(table.mask))
3037
3038        return new_table
3039
3040    def iter_fields_and_params(self):
3041        """
3042        Recursively iterate over all FIELD and PARAM elements in the
3043        TABLE.
3044        """
3045        for param in self.params:
3046            yield param
3047        for field in self.fields:
3048            yield field
3049        for group in self.groups:
3050            for field in group.iter_fields_and_params():
3051                yield field
3052
3053    get_field_by_id = _lookup_by_attr_factory(
3054        'ID', True, 'iter_fields_and_params', 'FIELD or PARAM',
3055        """
3056        Looks up a FIELD or PARAM element by the given ID.
3057        """)
3058
3059    get_field_by_id_or_name = _lookup_by_id_or_name_factory(
3060        'iter_fields_and_params', 'FIELD or PARAM',
3061        """
3062        Looks up a FIELD or PARAM element by the given ID or name.
3063        """)
3064
3065    get_fields_by_utype = _lookup_by_attr_factory(
3066        'utype', False, 'iter_fields_and_params', 'FIELD or PARAM',
3067        """
3068        Looks up a FIELD or PARAM element by the given utype and
3069        returns an iterator emitting all matches.
3070        """)
3071
3072    def iter_groups(self):
3073        """
3074        Recursively iterate over all GROUP elements in the TABLE.
3075        """
3076        for group in self.groups:
3077            yield group
3078            for g in group.iter_groups():
3079                yield g
3080
3081    get_group_by_id = _lookup_by_attr_factory(
3082        'ID', True, 'iter_groups', 'GROUP',
3083        """
3084        Looks up a GROUP element by the given ID.  Used by the group's
3085        "ref" attribute
3086        """)
3087
3088    get_groups_by_utype = _lookup_by_attr_factory(
3089        'utype', False, 'iter_groups', 'GROUP',
3090        """
3091        Looks up a GROUP element by the given utype and returns an
3092        iterator emitting all matches.
3093        """)
3094
3095    def iter_info(self):
3096        for info in self.infos:
3097            yield info
3098
3099
3100class Resource(Element, _IDProperty, _NameProperty, _UtypeProperty,
3101               _DescriptionProperty):
3102    """
3103    RESOURCE_ element: Groups TABLE_ and RESOURCE_ elements.
3104
3105    The keyword arguments correspond to setting members of the same
3106    name, documented below.
3107    """
3108
3109    def __init__(self, name=None, ID=None, utype=None, type='results',
3110                 id=None, config=None, pos=None, **kwargs):
3111        if config is None:
3112            config = {}
3113        self._config = config
3114        self._pos = pos
3115
3116        Element.__init__(self)
3117        self.name = name
3118        self.ID = resolve_id(ID, id, config, pos)
3119        self.utype = utype
3120        self.type = type
3121        self._extra_attributes = kwargs
3122        self.description = None
3123
3124        self._coordinate_systems = HomogeneousList(CooSys)
3125        self._time_systems = HomogeneousList(TimeSys)
3126        self._groups = HomogeneousList(Group)
3127        self._params = HomogeneousList(Param)
3128        self._infos = HomogeneousList(Info)
3129        self._links = HomogeneousList(Link)
3130        self._tables = HomogeneousList(Table)
3131        self._resources = HomogeneousList(Resource)
3132
3133        warn_unknown_attrs('RESOURCE', kwargs.keys(), config, pos)
3134
3135    def __repr__(self):
3136        buff = io.StringIO()
3137        w = XMLWriter(buff)
3138        w.element(
3139            self._element_name,
3140            attrib=w.object_attrs(self, self._attr_list))
3141        return buff.getvalue().strip()
3142
3143    @property
3144    def type(self):
3145        """
3146        [*required*] The type of the resource.  Must be either:
3147
3148          - 'results': This resource contains actual result values
3149            (default)
3150
3151          - 'meta': This resource contains only datatype descriptions
3152            (FIELD_ elements), but no actual data.
3153        """
3154        return self._type
3155
3156    @type.setter
3157    def type(self, type):
3158        if type not in ('results', 'meta'):
3159            vo_raise(E18, type, self._config, self._pos)
3160        self._type = type
3161
3162    @property
3163    def extra_attributes(self):
3164        """
3165        A dictionary of string keys to string values containing any
3166        extra attributes of the RESOURCE_ element that are not defined
3167        in the specification.  (The specification explicitly allows
3168        for extra attributes here, but nowhere else.)
3169        """
3170        return self._extra_attributes
3171
3172    @property
3173    def coordinate_systems(self):
3174        """
3175        A list of coordinate system definitions (COOSYS_ elements) for
3176        the RESOURCE_.  Must contain only `CooSys` objects.
3177        """
3178        return self._coordinate_systems
3179
3180    @property
3181    def time_systems(self):
3182        """
3183        A list of time system definitions (TIMESYS_ elements) for
3184        the RESOURCE_.  Must contain only `TimeSys` objects.
3185        """
3186        return self._time_systems
3187
3188    @property
3189    def infos(self):
3190        """
3191        A list of informational parameters (key-value pairs) for the
3192        resource.  Must only contain `Info` objects.
3193        """
3194        return self._infos
3195
3196    @property
3197    def groups(self):
3198        """
3199        A list of groups
3200        """
3201        return self._groups
3202
3203    @property
3204    def params(self):
3205        """
3206        A list of parameters (constant-valued columns) for the
3207        resource.  Must contain only `Param` objects.
3208        """
3209        return self._params
3210
3211    @property
3212    def links(self):
3213        """
3214        A list of links (pointers to other documents or servers
3215        through a URI) for the resource.  Must contain only `Link`
3216        objects.
3217        """
3218        return self._links
3219
3220    @property
3221    def tables(self):
3222        """
3223        A list of tables in the resource.  Must contain only
3224        `Table` objects.
3225        """
3226        return self._tables
3227
3228    @property
3229    def resources(self):
3230        """
3231        A list of nested resources inside this resource.  Must contain
3232        only `Resource` objects.
3233        """
3234        return self._resources
3235
3236    def _add_table(self, iterator, tag, data, config, pos):
3237        table = Table(self._votable, config=config, pos=pos, **data)
3238        self.tables.append(table)
3239        table.parse(iterator, config)
3240
3241    def _add_info(self, iterator, tag, data, config, pos):
3242        info = Info(config=config, pos=pos, **data)
3243        self.infos.append(info)
3244        info.parse(iterator, config)
3245
3246    def _add_group(self, iterator, tag, data, config, pos):
3247        group = Group(self, config=config, pos=pos, **data)
3248        self.groups.append(group)
3249        group.parse(iterator, config)
3250
3251    def _add_param(self, iterator, tag, data, config, pos):
3252        param = Param(self._votable, config=config, pos=pos, **data)
3253        self.params.append(param)
3254        param.parse(iterator, config)
3255
3256    def _add_coosys(self, iterator, tag, data, config, pos):
3257        coosys = CooSys(config=config, pos=pos, **data)
3258        self.coordinate_systems.append(coosys)
3259        coosys.parse(iterator, config)
3260
3261    def _add_timesys(self, iterator, tag, data, config, pos):
3262        timesys = TimeSys(config=config, pos=pos, **data)
3263        self.time_systems.append(timesys)
3264        timesys.parse(iterator, config)
3265
3266    def _add_resource(self, iterator, tag, data, config, pos):
3267        resource = Resource(config=config, pos=pos, **data)
3268        self.resources.append(resource)
3269        resource.parse(self._votable, iterator, config)
3270
3271    def _add_link(self, iterator, tag, data, config, pos):
3272        link = Link(config=config, pos=pos, **data)
3273        self.links.append(link)
3274        link.parse(iterator, config)
3275
3276    def parse(self, votable, iterator, config):
3277        self._votable = votable
3278
3279        tag_mapping = {
3280            'TABLE': self._add_table,
3281            'INFO': self._add_info,
3282            'PARAM': self._add_param,
3283            'GROUP': self._add_group,
3284            'COOSYS': self._add_coosys,
3285            'TIMESYS': self._add_timesys,
3286            'RESOURCE': self._add_resource,
3287            'LINK': self._add_link,
3288            'DESCRIPTION': self._ignore_add
3289            }
3290
3291        for start, tag, data, pos in iterator:
3292            if start:
3293                tag_mapping.get(tag, self._add_unknown_tag)(
3294                    iterator, tag, data, config, pos)
3295            elif tag == 'DESCRIPTION':
3296                if self.description is not None:
3297                    warn_or_raise(W17, W17, 'RESOURCE', config, pos)
3298                self.description = data or None
3299            elif tag == 'RESOURCE':
3300                break
3301
3302        del self._votable
3303
3304        return self
3305
3306    def to_xml(self, w, **kwargs):
3307        attrs = w.object_attrs(self, ('ID', 'type', 'utype'))
3308        attrs.update(self.extra_attributes)
3309        with w.tag('RESOURCE', attrib=attrs):
3310            if self.description is not None:
3311                w.element("DESCRIPTION", self.description, wrap=True)
3312            for element_set in (self.coordinate_systems, self.time_systems,
3313                                self.params, self.infos, self.links,
3314                                self.tables, self.resources):
3315                for element in element_set:
3316                    element.to_xml(w, **kwargs)
3317
3318    def iter_tables(self):
3319        """
3320        Recursively iterates over all tables in the resource and
3321        nested resources.
3322        """
3323        for table in self.tables:
3324            yield table
3325        for resource in self.resources:
3326            for table in resource.iter_tables():
3327                yield table
3328
3329    def iter_fields_and_params(self):
3330        """
3331        Recursively iterates over all FIELD_ and PARAM_ elements in
3332        the resource, its tables and nested resources.
3333        """
3334        for param in self.params:
3335            yield param
3336        for table in self.tables:
3337            for param in table.iter_fields_and_params():
3338                yield param
3339        for resource in self.resources:
3340            for param in resource.iter_fields_and_params():
3341                yield param
3342
3343    def iter_coosys(self):
3344        """
3345        Recursively iterates over all the COOSYS_ elements in the
3346        resource and nested resources.
3347        """
3348        for coosys in self.coordinate_systems:
3349            yield coosys
3350        for resource in self.resources:
3351            for coosys in resource.iter_coosys():
3352                yield coosys
3353
3354    def iter_timesys(self):
3355        """
3356        Recursively iterates over all the TIMESYS_ elements in the
3357        resource and nested resources.
3358        """
3359        for timesys in self.time_systems:
3360            yield timesys
3361        for resource in self.resources:
3362            for timesys in resource.iter_timesys():
3363                yield timesys
3364
3365    def iter_info(self):
3366        """
3367        Recursively iterates over all the INFO_ elements in the
3368        resource and nested resources.
3369        """
3370        for info in self.infos:
3371            yield info
3372        for table in self.tables:
3373            for info in table.iter_info():
3374                yield info
3375        for resource in self.resources:
3376            for info in resource.iter_info():
3377                yield info
3378
3379
3380class VOTableFile(Element, _IDProperty, _DescriptionProperty):
3381    """
3382    VOTABLE_ element: represents an entire file.
3383
3384    The keyword arguments correspond to setting members of the same
3385    name, documented below.
3386
3387    *version* is settable at construction time only, since conformance
3388    tests for building the rest of the structure depend on it.
3389    """
3390
3391    def __init__(self, ID=None, id=None, config=None, pos=None, version="1.4"):
3392        if config is None:
3393            config = {}
3394        self._config = config
3395        self._pos = pos
3396
3397        Element.__init__(self)
3398        self.ID = resolve_id(ID, id, config, pos)
3399        self.description = None
3400
3401        self._coordinate_systems = HomogeneousList(CooSys)
3402        self._time_systems = HomogeneousList(TimeSys)
3403        self._params = HomogeneousList(Param)
3404        self._infos = HomogeneousList(Info)
3405        self._resources = HomogeneousList(Resource)
3406        self._groups = HomogeneousList(Group)
3407
3408        version = str(version)
3409        if version == '1.0':
3410            warnings.warn('VOTable 1.0 support is deprecated in astropy 4.3 and will be '
3411                          'removed in a future release', AstropyDeprecationWarning)
3412        elif (version != '1.0') and (version not in self._version_namespace_map):
3413            allowed_from_map = "', '".join(self._version_namespace_map)
3414            raise ValueError(f"'version' should be in ('1.0', '{allowed_from_map}').")
3415
3416        self._version = version
3417
3418    def __repr__(self):
3419        n_tables = len(list(self.iter_tables()))
3420        return f'<VOTABLE>... {n_tables} tables ...</VOTABLE>'
3421
3422    @property
3423    def version(self):
3424        """
3425        The version of the VOTable specification that the file uses.
3426        """
3427        return self._version
3428
3429    @version.setter
3430    def version(self, version):
3431        version = str(version)
3432        if version not in self._version_namespace_map:
3433            allowed_from_map = "', '".join(self._version_namespace_map)
3434            raise ValueError(
3435                f"astropy.io.votable only supports VOTable versions '{allowed_from_map}'")
3436        self._version = version
3437
3438    @property
3439    def coordinate_systems(self):
3440        """
3441        A list of coordinate system descriptions for the file.  Must
3442        contain only `CooSys` objects.
3443        """
3444        return self._coordinate_systems
3445
3446    @property
3447    def time_systems(self):
3448        """
3449        A list of time system descriptions for the file.  Must
3450        contain only `TimeSys` objects.
3451        """
3452        return self._time_systems
3453
3454    @property
3455    def params(self):
3456        """
3457        A list of parameters (constant-valued columns) that apply to
3458        the entire file.  Must contain only `Param` objects.
3459        """
3460        return self._params
3461
3462    @property
3463    def infos(self):
3464        """
3465        A list of informational parameters (key-value pairs) for the
3466        entire file.  Must only contain `Info` objects.
3467        """
3468        return self._infos
3469
3470    @property
3471    def resources(self):
3472        """
3473        A list of resources, in the order they appear in the file.
3474        Must only contain `Resource` objects.
3475        """
3476        return self._resources
3477
3478    @property
3479    def groups(self):
3480        """
3481        A list of groups, in the order they appear in the file.  Only
3482        supported as a child of the VOTABLE element in VOTable 1.2 or
3483        later.
3484        """
3485        return self._groups
3486
3487    def _add_param(self, iterator, tag, data, config, pos):
3488        param = Param(self, config=config, pos=pos, **data)
3489        self.params.append(param)
3490        param.parse(iterator, config)
3491
3492    def _add_resource(self, iterator, tag, data, config, pos):
3493        resource = Resource(config=config, pos=pos, **data)
3494        self.resources.append(resource)
3495        resource.parse(self, iterator, config)
3496
3497    def _add_coosys(self, iterator, tag, data, config, pos):
3498        coosys = CooSys(config=config, pos=pos, **data)
3499        self.coordinate_systems.append(coosys)
3500        coosys.parse(iterator, config)
3501
3502    def _add_timesys(self, iterator, tag, data, config, pos):
3503        timesys = TimeSys(config=config, pos=pos, **data)
3504        self.time_systems.append(timesys)
3505        timesys.parse(iterator, config)
3506
3507    def _add_info(self, iterator, tag, data, config, pos):
3508        info = Info(config=config, pos=pos, **data)
3509        self.infos.append(info)
3510        info.parse(iterator, config)
3511
3512    def _add_group(self, iterator, tag, data, config, pos):
3513        if not config.get('version_1_2_or_later'):
3514            warn_or_raise(W26, W26, ('GROUP', 'VOTABLE', '1.2'), config, pos)
3515        group = Group(self, config=config, pos=pos, **data)
3516        self.groups.append(group)
3517        group.parse(iterator, config)
3518
3519    def _get_version_checks(self):
3520        config = {}
3521        config['version_1_1_or_later'] = \
3522            util.version_compare(self.version, '1.1') >= 0
3523        config['version_1_2_or_later'] = \
3524            util.version_compare(self.version, '1.2') >= 0
3525        config['version_1_3_or_later'] = \
3526            util.version_compare(self.version, '1.3') >= 0
3527        config['version_1_4_or_later'] = \
3528            util.version_compare(self.version, '1.4') >= 0
3529        return config
3530
3531    # Map VOTable version numbers to namespace URIs and schema information.
3532    _version_namespace_map = {
3533        # Version 1.0 isn't well-supported, but is allowed on parse (with a warning).
3534        # It used DTD rather than schema, so this information would not be useful.
3535        # By omitting 1.0 from this dict we can use the keys as the list of versions
3536        # that are allowed in various other checks.
3537        "1.1": {
3538            "namespace_uri": "http://www.ivoa.net/xml/VOTable/v1.1",
3539            "schema_location_attr": "xsi:noNamespaceSchemaLocation",
3540            "schema_location_value": "http://www.ivoa.net/xml/VOTable/v1.1"
3541        },
3542        "1.2": {
3543            "namespace_uri": "http://www.ivoa.net/xml/VOTable/v1.2",
3544            "schema_location_attr": "xsi:noNamespaceSchemaLocation",
3545            "schema_location_value": "http://www.ivoa.net/xml/VOTable/v1.2"
3546        },
3547        # With 1.3 we'll be more explicit with the schema location.
3548        # - xsi:schemaLocation uses the namespace name along with the URL
3549        #   to reference it.
3550        # - For convenience, but somewhat confusingly, the namespace URIs
3551        #   are also usable URLs for accessing an applicable schema.
3552        #   However to avoid confusion, we'll use the explicit schema URL.
3553        "1.3": {
3554            "namespace_uri": "http://www.ivoa.net/xml/VOTable/v1.3",
3555            "schema_location_attr": "xsi:schemaLocation",
3556            "schema_location_value":
3557            "http://www.ivoa.net/xml/VOTable/v1.3 http://www.ivoa.net/xml/VOTable/VOTable-1.3.xsd"
3558        },
3559        # With 1.4 namespace URIs stopped incrementing with minor version changes
3560        # so we use the same URI as with 1.3.  See this IVOA note for more info:
3561        # http://www.ivoa.net/documents/Notes/XMLVers/20180529/
3562        "1.4": {
3563            "namespace_uri": "http://www.ivoa.net/xml/VOTable/v1.3",
3564            "schema_location_attr": "xsi:schemaLocation",
3565            "schema_location_value":
3566            "http://www.ivoa.net/xml/VOTable/v1.3 http://www.ivoa.net/xml/VOTable/VOTable-1.4.xsd"
3567        }
3568    }
3569
3570    def parse(self, iterator, config):
3571        config['_current_table_number'] = 0
3572
3573        for start, tag, data, pos in iterator:
3574            if start:
3575                if tag == 'xml':
3576                    pass
3577                elif tag == 'VOTABLE':
3578                    if 'version' not in data:
3579                        warn_or_raise(W20, W20, self.version, config, pos)
3580                        config['version'] = self.version
3581                    else:
3582                        config['version'] = self._version = data['version']
3583                        if config['version'].lower().startswith('v'):
3584                            warn_or_raise(
3585                                W29, W29, config['version'], config, pos)
3586                            self._version = config['version'] = \
3587                                            config['version'][1:]
3588                        if config['version'] not in self._version_namespace_map:
3589                            vo_warn(W21, config['version'], config, pos)
3590
3591                    if 'xmlns' in data:
3592                        ns_info = self._version_namespace_map.get(config['version'], {})
3593                        correct_ns = ns_info.get('namespace_uri')
3594                        if data['xmlns'] != correct_ns:
3595                            vo_warn(W41, (correct_ns, data['xmlns']), config, pos)
3596                    else:
3597                        vo_warn(W42, (), config, pos)
3598
3599                    break
3600                else:
3601                    vo_raise(E19, (), config, pos)
3602        config.update(self._get_version_checks())
3603
3604        tag_mapping = {
3605            'PARAM': self._add_param,
3606            'RESOURCE': self._add_resource,
3607            'COOSYS': self._add_coosys,
3608            'TIMESYS': self._add_timesys,
3609            'INFO': self._add_info,
3610            'DEFINITIONS': self._add_definitions,
3611            'DESCRIPTION': self._ignore_add,
3612            'GROUP': self._add_group}
3613
3614        for start, tag, data, pos in iterator:
3615            if start:
3616                tag_mapping.get(tag, self._add_unknown_tag)(
3617                    iterator, tag, data, config, pos)
3618            elif tag == 'DESCRIPTION':
3619                if self.description is not None:
3620                    warn_or_raise(W17, W17, 'VOTABLE', config, pos)
3621                self.description = data or None
3622
3623        if not len(self.resources) and config['version_1_2_or_later']:
3624            warn_or_raise(W53, W53, (), config, pos)
3625
3626        return self
3627
3628    def to_xml(self, fd, compressed=False, tabledata_format=None,
3629               _debug_python_based_parser=False, _astropy_version=None):
3630        """
3631        Write to an XML file.
3632
3633        Parameters
3634        ----------
3635        fd : str or file-like
3636            Where to write the file. If a file-like object, must be writable.
3637
3638        compressed : bool, optional
3639            When `True`, write to a gzip-compressed file.  (Default:
3640            `False`)
3641
3642        tabledata_format : str, optional
3643            Override the format of the table(s) data to write.  Must
3644            be one of ``tabledata`` (text representation), ``binary`` or
3645            ``binary2``.  By default, use the format that was specified
3646            in each `Table` object as it was created or read in.  See
3647            :ref:`astropy:votable-serialization`.
3648        """
3649        if tabledata_format is not None:
3650            if tabledata_format.lower() not in (
3651                    'tabledata', 'binary', 'binary2'):
3652                raise ValueError(f"Unknown format type '{format}'")
3653
3654        kwargs = {
3655            'version': self.version,
3656            'tabledata_format':
3657                tabledata_format,
3658            '_debug_python_based_parser': _debug_python_based_parser,
3659            '_group_number': 1}
3660        kwargs.update(self._get_version_checks())
3661
3662        with util.convert_to_writable_filelike(
3663            fd, compressed=compressed) as fd:
3664            w = XMLWriter(fd)
3665            version = self.version
3666            if _astropy_version is None:
3667                lib_version = astropy_version
3668            else:
3669                lib_version = _astropy_version
3670
3671            xml_header = """
3672<?xml version="1.0" encoding="utf-8"?>
3673<!-- Produced with astropy.io.votable version {lib_version}
3674     http://www.astropy.org/ -->\n"""
3675            w.write(xml_header.lstrip().format(**locals()))
3676
3677            # Build the VOTABLE tag attributes.
3678            votable_attr = {
3679                'version': version,
3680                'xmlns:xsi': "http://www.w3.org/2001/XMLSchema-instance"
3681            }
3682            ns_info = self._version_namespace_map.get(version, {})
3683            namespace_uri = ns_info.get('namespace_uri')
3684            if namespace_uri:
3685                votable_attr['xmlns'] = namespace_uri
3686            schema_location_attr = ns_info.get('schema_location_attr')
3687            schema_location_value = ns_info.get('schema_location_value')
3688            if schema_location_attr and schema_location_value:
3689                votable_attr[schema_location_attr] = schema_location_value
3690
3691            with w.tag('VOTABLE', votable_attr):
3692                if self.description is not None:
3693                    w.element("DESCRIPTION", self.description, wrap=True)
3694                element_sets = [self.coordinate_systems, self.time_systems,
3695                                self.params, self.infos, self.resources]
3696                if kwargs['version_1_2_or_later']:
3697                    element_sets[0] = self.groups
3698                for element_set in element_sets:
3699                    for element in element_set:
3700                        element.to_xml(w, **kwargs)
3701
3702    def iter_tables(self):
3703        """
3704        Iterates over all tables in the VOTable file in a "flat" way,
3705        ignoring the nesting of resources etc.
3706        """
3707        for resource in self.resources:
3708            for table in resource.iter_tables():
3709                yield table
3710
3711    def get_first_table(self):
3712        """
3713        Often, you know there is only one table in the file, and
3714        that's all you need.  This method returns that first table.
3715        """
3716        for table in self.iter_tables():
3717            if not table.is_empty():
3718                return table
3719        raise IndexError("No table found in VOTABLE file.")
3720
3721    get_table_by_id = _lookup_by_attr_factory(
3722        'ID', True, 'iter_tables', 'TABLE',
3723        """
3724        Looks up a TABLE_ element by the given ID.  Used by the table
3725        "ref" attribute.
3726        """)
3727
3728    get_tables_by_utype = _lookup_by_attr_factory(
3729        'utype', False, 'iter_tables', 'TABLE',
3730        """
3731        Looks up a TABLE_ element by the given utype, and returns an
3732        iterator emitting all matches.
3733        """)
3734
3735    def get_table_by_index(self, idx):
3736        """
3737        Get a table by its ordinal position in the file.
3738        """
3739        for i, table in enumerate(self.iter_tables()):
3740            if i == idx:
3741                return table
3742        raise IndexError(
3743            f"No table at index {idx:d} found in VOTABLE file.")
3744
3745    def iter_fields_and_params(self):
3746        """
3747        Recursively iterate over all FIELD_ and PARAM_ elements in the
3748        VOTABLE_ file.
3749        """
3750        for resource in self.resources:
3751            for field in resource.iter_fields_and_params():
3752                yield field
3753
3754    get_field_by_id = _lookup_by_attr_factory(
3755        'ID', True, 'iter_fields_and_params', 'FIELD',
3756        """
3757        Looks up a FIELD_ element by the given ID_.  Used by the field's
3758        "ref" attribute.
3759        """)
3760
3761    get_fields_by_utype = _lookup_by_attr_factory(
3762        'utype', False, 'iter_fields_and_params', 'FIELD',
3763        """
3764        Looks up a FIELD_ element by the given utype and returns an
3765        iterator emitting all matches.
3766        """)
3767
3768    get_field_by_id_or_name = _lookup_by_id_or_name_factory(
3769        'iter_fields_and_params', 'FIELD',
3770        """
3771        Looks up a FIELD_ element by the given ID_ or name.
3772        """)
3773
3774    def iter_values(self):
3775        """
3776        Recursively iterate over all VALUES_ elements in the VOTABLE_
3777        file.
3778        """
3779        for field in self.iter_fields_and_params():
3780            yield field.values
3781
3782    get_values_by_id = _lookup_by_attr_factory(
3783        'ID', True, 'iter_values', 'VALUES',
3784        """
3785        Looks up a VALUES_ element by the given ID.  Used by the values
3786        "ref" attribute.
3787        """)
3788
3789    def iter_groups(self):
3790        """
3791        Recursively iterate over all GROUP_ elements in the VOTABLE_
3792        file.
3793        """
3794        for table in self.iter_tables():
3795            for group in table.iter_groups():
3796                yield group
3797
3798    get_group_by_id = _lookup_by_attr_factory(
3799        'ID', True, 'iter_groups', 'GROUP',
3800        """
3801        Looks up a GROUP_ element by the given ID.  Used by the group's
3802        "ref" attribute
3803        """)
3804
3805    get_groups_by_utype = _lookup_by_attr_factory(
3806        'utype', False, 'iter_groups', 'GROUP',
3807        """
3808        Looks up a GROUP_ element by the given utype and returns an
3809        iterator emitting all matches.
3810        """)
3811
3812    def iter_coosys(self):
3813        """
3814        Recursively iterate over all COOSYS_ elements in the VOTABLE_
3815        file.
3816        """
3817        for coosys in self.coordinate_systems:
3818            yield coosys
3819        for resource in self.resources:
3820            for coosys in resource.iter_coosys():
3821                yield coosys
3822
3823    get_coosys_by_id = _lookup_by_attr_factory(
3824        'ID', True, 'iter_coosys', 'COOSYS',
3825        """Looks up a COOSYS_ element by the given ID.""")
3826
3827    def iter_timesys(self):
3828        """
3829        Recursively iterate over all TIMESYS_ elements in the VOTABLE_
3830        file.
3831        """
3832        for timesys in self.time_systems:
3833            yield timesys
3834        for resource in self.resources:
3835            for timesys in resource.iter_timesys():
3836                yield timesys
3837
3838    get_timesys_by_id = _lookup_by_attr_factory(
3839        'ID', True, 'iter_timesys', 'TIMESYS',
3840        """Looks up a TIMESYS_ element by the given ID.""")
3841
3842    def iter_info(self):
3843        """
3844        Recursively iterate over all INFO_ elements in the VOTABLE_
3845        file.
3846        """
3847        for info in self.infos:
3848            yield info
3849        for resource in self.resources:
3850            for info in resource.iter_info():
3851                yield info
3852
3853    get_info_by_id = _lookup_by_attr_factory(
3854        'ID', True, 'iter_info', 'INFO',
3855        """Looks up a INFO element by the given ID.""")
3856
3857    def set_all_tables_format(self, format):
3858        """
3859        Set the output storage format of all tables in the file.
3860        """
3861        for table in self.iter_tables():
3862            table.format = format
3863
3864    @classmethod
3865    def from_table(cls, table, table_id=None):
3866        """
3867        Create a `VOTableFile` instance from a given
3868        `astropy.table.Table` instance.
3869
3870        Parameters
3871        ----------
3872        table_id : str, optional
3873            Set the given ID attribute on the returned Table instance.
3874        """
3875        votable_file = cls()
3876        resource = Resource()
3877        votable = Table.from_table(votable_file, table)
3878        if table_id is not None:
3879            votable.ID = table_id
3880        resource.tables.append(votable)
3881        votable_file.resources.append(resource)
3882        return votable_file
3883