1"""DatetimeIndex analog for cftime.datetime objects"""
2# The pandas.Index subclass defined here was copied and adapted for
3# use with cftime.datetime objects based on the source code defining
4# pandas.DatetimeIndex.
5
6# For reference, here is a copy of the pandas copyright notice:
7
8# (c) 2011-2012, Lambda Foundry, Inc. and PyData Development Team
9# All rights reserved.
10
11# Copyright (c) 2008-2011 AQR Capital Management, LLC
12# All rights reserved.
13
14# Redistribution and use in source and binary forms, with or without
15# modification, are permitted provided that the following conditions are
16# met:
17
18#     * Redistributions of source code must retain the above copyright
19#        notice, this list of conditions and the following disclaimer.
20
21#     * Redistributions in binary form must reproduce the above
22#        copyright notice, this list of conditions and the following
23#        disclaimer in the documentation and/or other materials provided
24#        with the distribution.
25
26#     * Neither the name of the copyright holder nor the names of any
27#        contributors may be used to endorse or promote products derived
28#        from this software without specific prior written permission.
29
30# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDER AND CONTRIBUTORS
31# "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
32# LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
33# A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
34# OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
35# SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
36# LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
37# DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
38# THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
39# (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
40# OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
41
42import re
43import warnings
44from datetime import timedelta
45from distutils.version import LooseVersion
46from typing import Tuple, Type
47
48import numpy as np
49import pandas as pd
50
51from xarray.core.utils import is_scalar
52
53from ..core.common import _contains_cftime_datetimes
54from ..core.options import OPTIONS
55from .times import _STANDARD_CALENDARS, cftime_to_nptime, infer_calendar_name
56
57try:
58    import cftime
59except ImportError:
60    cftime = None
61
62
63# constants for cftimeindex.repr
64CFTIME_REPR_LENGTH = 19
65ITEMS_IN_REPR_MAX_ELSE_ELLIPSIS = 100
66REPR_ELLIPSIS_SHOW_ITEMS_FRONT_END = 10
67
68
69OUT_OF_BOUNDS_TIMEDELTA_ERRORS: Tuple[Type[Exception], ...]
70try:
71    OUT_OF_BOUNDS_TIMEDELTA_ERRORS = (pd.errors.OutOfBoundsTimedelta, OverflowError)
72except AttributeError:
73    OUT_OF_BOUNDS_TIMEDELTA_ERRORS = (OverflowError,)
74
75
76def named(name, pattern):
77    return "(?P<" + name + ">" + pattern + ")"
78
79
80def optional(x):
81    return "(?:" + x + ")?"
82
83
84def trailing_optional(xs):
85    if not xs:
86        return ""
87    return xs[0] + optional(trailing_optional(xs[1:]))
88
89
90def build_pattern(date_sep=r"\-", datetime_sep=r"T", time_sep=r"\:"):
91    pieces = [
92        (None, "year", r"\d{4}"),
93        (date_sep, "month", r"\d{2}"),
94        (date_sep, "day", r"\d{2}"),
95        (datetime_sep, "hour", r"\d{2}"),
96        (time_sep, "minute", r"\d{2}"),
97        (time_sep, "second", r"\d{2}"),
98    ]
99    pattern_list = []
100    for sep, name, sub_pattern in pieces:
101        pattern_list.append((sep if sep else "") + named(name, sub_pattern))
102        # TODO: allow timezone offsets?
103    return "^" + trailing_optional(pattern_list) + "$"
104
105
106_BASIC_PATTERN = build_pattern(date_sep="", time_sep="")
107_EXTENDED_PATTERN = build_pattern()
108_CFTIME_PATTERN = build_pattern(datetime_sep=" ")
109_PATTERNS = [_BASIC_PATTERN, _EXTENDED_PATTERN, _CFTIME_PATTERN]
110
111
112def parse_iso8601_like(datetime_string):
113    for pattern in _PATTERNS:
114        match = re.match(pattern, datetime_string)
115        if match:
116            return match.groupdict()
117    raise ValueError(
118        f"no ISO-8601 or cftime-string-like match for string: {datetime_string}"
119    )
120
121
122def _parse_iso8601_with_reso(date_type, timestr):
123    if cftime is None:
124        raise ModuleNotFoundError("No module named 'cftime'")
125
126    default = date_type(1, 1, 1)
127    result = parse_iso8601_like(timestr)
128    replace = {}
129
130    for attr in ["year", "month", "day", "hour", "minute", "second"]:
131        value = result.get(attr, None)
132        if value is not None:
133            # Note ISO8601 conventions allow for fractional seconds.
134            # TODO: Consider adding support for sub-second resolution?
135            replace[attr] = int(value)
136            resolution = attr
137    return default.replace(**replace), resolution
138
139
140def _parsed_string_to_bounds(date_type, resolution, parsed):
141    """Generalization of
142    pandas.tseries.index.DatetimeIndex._parsed_string_to_bounds
143    for use with non-standard calendars and cftime.datetime
144    objects.
145    """
146    if resolution == "year":
147        return (
148            date_type(parsed.year, 1, 1),
149            date_type(parsed.year + 1, 1, 1) - timedelta(microseconds=1),
150        )
151    elif resolution == "month":
152        if parsed.month == 12:
153            end = date_type(parsed.year + 1, 1, 1) - timedelta(microseconds=1)
154        else:
155            end = date_type(parsed.year, parsed.month + 1, 1) - timedelta(
156                microseconds=1
157            )
158        return date_type(parsed.year, parsed.month, 1), end
159    elif resolution == "day":
160        start = date_type(parsed.year, parsed.month, parsed.day)
161        return start, start + timedelta(days=1, microseconds=-1)
162    elif resolution == "hour":
163        start = date_type(parsed.year, parsed.month, parsed.day, parsed.hour)
164        return start, start + timedelta(hours=1, microseconds=-1)
165    elif resolution == "minute":
166        start = date_type(
167            parsed.year, parsed.month, parsed.day, parsed.hour, parsed.minute
168        )
169        return start, start + timedelta(minutes=1, microseconds=-1)
170    elif resolution == "second":
171        start = date_type(
172            parsed.year,
173            parsed.month,
174            parsed.day,
175            parsed.hour,
176            parsed.minute,
177            parsed.second,
178        )
179        return start, start + timedelta(seconds=1, microseconds=-1)
180    else:
181        raise KeyError
182
183
184def get_date_field(datetimes, field):
185    """Adapted from pandas.tslib.get_date_field"""
186    return np.array([getattr(date, field) for date in datetimes])
187
188
189def _field_accessor(name, docstring=None, min_cftime_version="0.0"):
190    """Adapted from pandas.tseries.index._field_accessor"""
191
192    def f(self, min_cftime_version=min_cftime_version):
193        if cftime is None:
194            raise ModuleNotFoundError("No module named 'cftime'")
195
196        version = cftime.__version__
197
198        if LooseVersion(version) >= LooseVersion(min_cftime_version):
199            return get_date_field(self._data, name)
200        else:
201            raise ImportError(
202                "The {!r} accessor requires a minimum "
203                "version of cftime of {}. Found an "
204                "installed version of {}.".format(name, min_cftime_version, version)
205            )
206
207    f.__name__ = name
208    f.__doc__ = docstring
209    return property(f)
210
211
212def get_date_type(self):
213    if self._data.size:
214        return type(self._data[0])
215    else:
216        return None
217
218
219def assert_all_valid_date_type(data):
220    if cftime is None:
221        raise ModuleNotFoundError("No module named 'cftime'")
222
223    if len(data) > 0:
224        sample = data[0]
225        date_type = type(sample)
226        if not isinstance(sample, cftime.datetime):
227            raise TypeError(
228                "CFTimeIndex requires cftime.datetime "
229                "objects. Got object of {}.".format(date_type)
230            )
231        if not all(isinstance(value, date_type) for value in data):
232            raise TypeError(
233                "CFTimeIndex requires using datetime "
234                "objects of all the same type.  Got\n{}.".format(data)
235            )
236
237
238def format_row(times, indent=0, separator=", ", row_end=",\n"):
239    """Format a single row from format_times."""
240    return indent * " " + separator.join(map(str, times)) + row_end
241
242
243def format_times(
244    index,
245    max_width,
246    offset,
247    separator=", ",
248    first_row_offset=0,
249    intermediate_row_end=",\n",
250    last_row_end="",
251):
252    """Format values of cftimeindex as pd.Index."""
253    n_per_row = max(max_width // (CFTIME_REPR_LENGTH + len(separator)), 1)
254    n_rows = int(np.ceil(len(index) / n_per_row))
255
256    representation = ""
257    for row in range(n_rows):
258        indent = first_row_offset if row == 0 else offset
259        row_end = last_row_end if row == n_rows - 1 else intermediate_row_end
260        times_for_row = index[row * n_per_row : (row + 1) * n_per_row]
261        representation += format_row(
262            times_for_row, indent=indent, separator=separator, row_end=row_end
263        )
264
265    return representation
266
267
268def format_attrs(index, separator=", "):
269    """Format attributes of CFTimeIndex for __repr__."""
270    attrs = {
271        "dtype": f"'{index.dtype}'",
272        "length": f"{len(index)}",
273        "calendar": f"'{index.calendar}'",
274        "freq": f"'{index.freq}'" if len(index) >= 3 else None,
275    }
276
277    attrs_str = [f"{k}={v}" for k, v in attrs.items()]
278    attrs_str = f"{separator}".join(attrs_str)
279    return attrs_str
280
281
282class CFTimeIndex(pd.Index):
283    """Custom Index for working with CF calendars and dates
284
285    All elements of a CFTimeIndex must be cftime.datetime objects.
286
287    Parameters
288    ----------
289    data : array or CFTimeIndex
290        Sequence of cftime.datetime objects to use in index
291    name : str, default: None
292        Name of the resulting index
293
294    See Also
295    --------
296    cftime_range
297    """
298
299    year = _field_accessor("year", "The year of the datetime")
300    month = _field_accessor("month", "The month of the datetime")
301    day = _field_accessor("day", "The days of the datetime")
302    hour = _field_accessor("hour", "The hours of the datetime")
303    minute = _field_accessor("minute", "The minutes of the datetime")
304    second = _field_accessor("second", "The seconds of the datetime")
305    microsecond = _field_accessor("microsecond", "The microseconds of the datetime")
306    dayofyear = _field_accessor(
307        "dayofyr", "The ordinal day of year of the datetime", "1.0.2.1"
308    )
309    dayofweek = _field_accessor("dayofwk", "The day of week of the datetime", "1.0.2.1")
310    days_in_month = _field_accessor(
311        "daysinmonth", "The number of days in the month of the datetime", "1.1.0.0"
312    )
313    date_type = property(get_date_type)
314
315    def __new__(cls, data, name=None):
316        assert_all_valid_date_type(data)
317        if name is None and hasattr(data, "name"):
318            name = data.name
319
320        result = object.__new__(cls)
321        result._data = np.array(data, dtype="O")
322        result.name = name
323        result._cache = {}
324        return result
325
326    def __repr__(self):
327        """
328        Return a string representation for this object.
329        """
330        klass_name = type(self).__name__
331        display_width = OPTIONS["display_width"]
332        offset = len(klass_name) + 2
333
334        if len(self) <= ITEMS_IN_REPR_MAX_ELSE_ELLIPSIS:
335            datastr = format_times(
336                self.values, display_width, offset=offset, first_row_offset=0
337            )
338        else:
339            front_str = format_times(
340                self.values[:REPR_ELLIPSIS_SHOW_ITEMS_FRONT_END],
341                display_width,
342                offset=offset,
343                first_row_offset=0,
344                last_row_end=",",
345            )
346            end_str = format_times(
347                self.values[-REPR_ELLIPSIS_SHOW_ITEMS_FRONT_END:],
348                display_width,
349                offset=offset,
350                first_row_offset=offset,
351            )
352            datastr = "\n".join([front_str, f"{' '*offset}...", end_str])
353
354        attrs_str = format_attrs(self)
355        # oneliner only if smaller than display_width
356        full_repr_str = f"{klass_name}([{datastr}], {attrs_str})"
357        if len(full_repr_str) > display_width:
358            # if attrs_str too long, one per line
359            if len(attrs_str) >= display_width - offset:
360                attrs_str = attrs_str.replace(",", f",\n{' '*(offset-2)}")
361            full_repr_str = f"{klass_name}([{datastr}],\n{' '*(offset-1)}{attrs_str})"
362
363        return full_repr_str
364
365    def _partial_date_slice(self, resolution, parsed):
366        """Adapted from
367        pandas.tseries.index.DatetimeIndex._partial_date_slice
368
369        Note that when using a CFTimeIndex, if a partial-date selection
370        returns a single element, it will never be converted to a scalar
371        coordinate; this is in slight contrast to the behavior when using
372        a DatetimeIndex, which sometimes will return a DataArray with a scalar
373        coordinate depending on the resolution of the datetimes used in
374        defining the index.  For example:
375
376        >>> from cftime import DatetimeNoLeap
377        >>> da = xr.DataArray(
378        ...     [1, 2],
379        ...     coords=[[DatetimeNoLeap(2001, 1, 1), DatetimeNoLeap(2001, 2, 1)]],
380        ...     dims=["time"],
381        ... )
382        >>> da.sel(time="2001-01-01")
383        <xarray.DataArray (time: 1)>
384        array([1])
385        Coordinates:
386          * time     (time) object 2001-01-01 00:00:00
387        >>> da = xr.DataArray(
388        ...     [1, 2],
389        ...     coords=[[pd.Timestamp(2001, 1, 1), pd.Timestamp(2001, 2, 1)]],
390        ...     dims=["time"],
391        ... )
392        >>> da.sel(time="2001-01-01")
393        <xarray.DataArray ()>
394        array(1)
395        Coordinates:
396            time     datetime64[ns] 2001-01-01
397        >>> da = xr.DataArray(
398        ...     [1, 2],
399        ...     coords=[[pd.Timestamp(2001, 1, 1, 1), pd.Timestamp(2001, 2, 1)]],
400        ...     dims=["time"],
401        ... )
402        >>> da.sel(time="2001-01-01")
403        <xarray.DataArray (time: 1)>
404        array([1])
405        Coordinates:
406          * time     (time) datetime64[ns] 2001-01-01T01:00:00
407        """
408        start, end = _parsed_string_to_bounds(self.date_type, resolution, parsed)
409
410        times = self._data
411
412        if self.is_monotonic:
413            if len(times) and (
414                (start < times[0] and end < times[0])
415                or (start > times[-1] and end > times[-1])
416            ):
417                # we are out of range
418                raise KeyError
419
420            # a monotonic (sorted) series can be sliced
421            left = times.searchsorted(start, side="left")
422            right = times.searchsorted(end, side="right")
423            return slice(left, right)
424
425        lhs_mask = times >= start
426        rhs_mask = times <= end
427        return np.flatnonzero(lhs_mask & rhs_mask)
428
429    def _get_string_slice(self, key):
430        """Adapted from pandas.tseries.index.DatetimeIndex._get_string_slice"""
431        parsed, resolution = _parse_iso8601_with_reso(self.date_type, key)
432        try:
433            loc = self._partial_date_slice(resolution, parsed)
434        except KeyError:
435            raise KeyError(key)
436        return loc
437
438    def _get_nearest_indexer(self, target, limit, tolerance):
439        """Adapted from pandas.Index._get_nearest_indexer"""
440        left_indexer = self.get_indexer(target, "pad", limit=limit)
441        right_indexer = self.get_indexer(target, "backfill", limit=limit)
442        left_distances = abs(self.values[left_indexer] - target.values)
443        right_distances = abs(self.values[right_indexer] - target.values)
444
445        if self.is_monotonic_increasing:
446            condition = (left_distances < right_distances) | (right_indexer == -1)
447        else:
448            condition = (left_distances <= right_distances) | (right_indexer == -1)
449        indexer = np.where(condition, left_indexer, right_indexer)
450
451        if tolerance is not None:
452            indexer = self._filter_indexer_tolerance(target, indexer, tolerance)
453        return indexer
454
455    def _filter_indexer_tolerance(self, target, indexer, tolerance):
456        """Adapted from pandas.Index._filter_indexer_tolerance"""
457        if isinstance(target, pd.Index):
458            distance = abs(self.values[indexer] - target.values)
459        else:
460            distance = abs(self.values[indexer] - target)
461        indexer = np.where(distance <= tolerance, indexer, -1)
462        return indexer
463
464    def get_loc(self, key, method=None, tolerance=None):
465        """Adapted from pandas.tseries.index.DatetimeIndex.get_loc"""
466        if isinstance(key, str):
467            return self._get_string_slice(key)
468        else:
469            return pd.Index.get_loc(self, key, method=method, tolerance=tolerance)
470
471    def _maybe_cast_slice_bound(self, label, side, kind=None):
472        """Adapted from
473        pandas.tseries.index.DatetimeIndex._maybe_cast_slice_bound
474
475        Note that we have never used the kind argument in CFTimeIndex and it is
476        deprecated as of pandas version 1.3.0.  It exists only for compatibility
477        reasons.  We can remove it when our minimum version of pandas is 1.3.0.
478        """
479        if not isinstance(label, str):
480            return label
481
482        parsed, resolution = _parse_iso8601_with_reso(self.date_type, label)
483        start, end = _parsed_string_to_bounds(self.date_type, resolution, parsed)
484        if self.is_monotonic_decreasing and len(self) > 1:
485            return end if side == "left" else start
486        return start if side == "left" else end
487
488    # TODO: Add ability to use integer range outside of iloc?
489    # e.g. series[1:5].
490    def get_value(self, series, key):
491        """Adapted from pandas.tseries.index.DatetimeIndex.get_value"""
492        if np.asarray(key).dtype == np.dtype(bool):
493            return series.iloc[key]
494        elif isinstance(key, slice):
495            return series.iloc[self.slice_indexer(key.start, key.stop, key.step)]
496        else:
497            return series.iloc[self.get_loc(key)]
498
499    def __contains__(self, key):
500        """Adapted from
501        pandas.tseries.base.DatetimeIndexOpsMixin.__contains__"""
502        try:
503            result = self.get_loc(key)
504            return (
505                is_scalar(result)
506                or type(result) == slice
507                or (isinstance(result, np.ndarray) and result.size)
508            )
509        except (KeyError, TypeError, ValueError):
510            return False
511
512    def contains(self, key):
513        """Needed for .loc based partial-string indexing"""
514        return self.__contains__(key)
515
516    def shift(self, n, freq):
517        """Shift the CFTimeIndex a multiple of the given frequency.
518
519        See the documentation for :py:func:`~xarray.cftime_range` for a
520        complete listing of valid frequency strings.
521
522        Parameters
523        ----------
524        n : int
525            Periods to shift by
526        freq : str or datetime.timedelta
527            A frequency string or datetime.timedelta object to shift by
528
529        Returns
530        -------
531        CFTimeIndex
532
533        See Also
534        --------
535        pandas.DatetimeIndex.shift
536
537        Examples
538        --------
539        >>> index = xr.cftime_range("2000", periods=1, freq="M")
540        >>> index
541        CFTimeIndex([2000-01-31 00:00:00],
542                    dtype='object', length=1, calendar='gregorian', freq=None)
543        >>> index.shift(1, "M")
544        CFTimeIndex([2000-02-29 00:00:00],
545                    dtype='object', length=1, calendar='gregorian', freq=None)
546        """
547        from .cftime_offsets import to_offset
548
549        if not isinstance(n, int):
550            raise TypeError(f"'n' must be an int, got {n}.")
551        if isinstance(freq, timedelta):
552            return self + n * freq
553        elif isinstance(freq, str):
554            return self + n * to_offset(freq)
555        else:
556            raise TypeError(
557                "'freq' must be of type "
558                "str or datetime.timedelta, got {}.".format(freq)
559            )
560
561    def __add__(self, other):
562        if isinstance(other, pd.TimedeltaIndex):
563            other = other.to_pytimedelta()
564        return CFTimeIndex(np.array(self) + other)
565
566    def __radd__(self, other):
567        if isinstance(other, pd.TimedeltaIndex):
568            other = other.to_pytimedelta()
569        return CFTimeIndex(other + np.array(self))
570
571    def __sub__(self, other):
572        if _contains_datetime_timedeltas(other):
573            return CFTimeIndex(np.array(self) - other)
574        elif isinstance(other, pd.TimedeltaIndex):
575            return CFTimeIndex(np.array(self) - other.to_pytimedelta())
576        elif _contains_cftime_datetimes(np.array(other)):
577            try:
578                return pd.TimedeltaIndex(np.array(self) - np.array(other))
579            except OUT_OF_BOUNDS_TIMEDELTA_ERRORS:
580                raise ValueError(
581                    "The time difference exceeds the range of values "
582                    "that can be expressed at the nanosecond resolution."
583                )
584        else:
585            return NotImplemented
586
587    def __rsub__(self, other):
588        try:
589            return pd.TimedeltaIndex(other - np.array(self))
590        except OUT_OF_BOUNDS_TIMEDELTA_ERRORS:
591            raise ValueError(
592                "The time difference exceeds the range of values "
593                "that can be expressed at the nanosecond resolution."
594            )
595
596    def to_datetimeindex(self, unsafe=False):
597        """If possible, convert this index to a pandas.DatetimeIndex.
598
599        Parameters
600        ----------
601        unsafe : bool
602            Flag to turn off warning when converting from a CFTimeIndex with
603            a non-standard calendar to a DatetimeIndex (default ``False``).
604
605        Returns
606        -------
607        pandas.DatetimeIndex
608
609        Raises
610        ------
611        ValueError
612            If the CFTimeIndex contains dates that are not possible in the
613            standard calendar or outside the pandas.Timestamp-valid range.
614
615        Warns
616        -----
617        RuntimeWarning
618            If converting from a non-standard calendar to a DatetimeIndex.
619
620        Warnings
621        --------
622        Note that for non-standard calendars, this will change the calendar
623        type of the index.  In that case the result of this method should be
624        used with caution.
625
626        Examples
627        --------
628        >>> times = xr.cftime_range("2000", periods=2, calendar="gregorian")
629        >>> times
630        CFTimeIndex([2000-01-01 00:00:00, 2000-01-02 00:00:00],
631                    dtype='object', length=2, calendar='gregorian', freq=None)
632        >>> times.to_datetimeindex()
633        DatetimeIndex(['2000-01-01', '2000-01-02'], dtype='datetime64[ns]', freq=None)
634        """
635        nptimes = cftime_to_nptime(self)
636        calendar = infer_calendar_name(self)
637        if calendar not in _STANDARD_CALENDARS and not unsafe:
638            warnings.warn(
639                "Converting a CFTimeIndex with dates from a non-standard "
640                "calendar, {!r}, to a pandas.DatetimeIndex, which uses dates "
641                "from the standard calendar.  This may lead to subtle errors "
642                "in operations that depend on the length of time between "
643                "dates.".format(calendar),
644                RuntimeWarning,
645                stacklevel=2,
646            )
647        return pd.DatetimeIndex(nptimes)
648
649    def strftime(self, date_format):
650        """
651        Return an Index of formatted strings specified by date_format, which
652        supports the same string format as the python standard library. Details
653        of the string format can be found in `python string format doc
654        <https://docs.python.org/3/library/datetime.html#strftime-strptime-behavior>`__
655
656        Parameters
657        ----------
658        date_format : str
659            Date format string (e.g. "%Y-%m-%d")
660
661        Returns
662        -------
663        pandas.Index
664            Index of formatted strings
665
666        Examples
667        --------
668        >>> rng = xr.cftime_range(
669        ...     start="2000", periods=5, freq="2MS", calendar="noleap"
670        ... )
671        >>> rng.strftime("%B %d, %Y, %r")
672        Index(['January 01, 2000, 12:00:00 AM', 'March 01, 2000, 12:00:00 AM',
673               'May 01, 2000, 12:00:00 AM', 'July 01, 2000, 12:00:00 AM',
674               'September 01, 2000, 12:00:00 AM'],
675              dtype='object')
676        """
677        return pd.Index([date.strftime(date_format) for date in self._data])
678
679    @property
680    def asi8(self):
681        """Convert to integers with units of microseconds since 1970-01-01."""
682        from ..core.resample_cftime import exact_cftime_datetime_difference
683
684        epoch = self.date_type(1970, 1, 1)
685        return np.array(
686            [
687                _total_microseconds(exact_cftime_datetime_difference(epoch, date))
688                for date in self.values
689            ],
690            dtype=np.int64,
691        )
692
693    @property
694    def calendar(self):
695        """The calendar used by the datetimes in the index."""
696        from .times import infer_calendar_name
697
698        return infer_calendar_name(self)
699
700    @property
701    def freq(self):
702        """The frequency used by the dates in the index."""
703        from .frequencies import infer_freq
704
705        return infer_freq(self)
706
707    def _round_via_method(self, freq, method):
708        """Round dates using a specified method."""
709        from .cftime_offsets import CFTIME_TICKS, to_offset
710
711        offset = to_offset(freq)
712        if not isinstance(offset, CFTIME_TICKS):
713            raise ValueError(f"{offset} is a non-fixed frequency")
714
715        unit = _total_microseconds(offset.as_timedelta())
716        values = self.asi8
717        rounded = method(values, unit)
718        return _cftimeindex_from_i8(rounded, self.date_type, self.name)
719
720    def floor(self, freq):
721        """Round dates down to fixed frequency.
722
723        Parameters
724        ----------
725        freq : str
726            The frequency level to round the index to.  Must be a fixed
727            frequency like 'S' (second) not 'ME' (month end).  See `frequency
728            aliases <https://pandas.pydata.org/pandas-docs/stable/user_guide/timeseries.html#offset-aliases>`_
729            for a list of possible values.
730
731        Returns
732        -------
733        CFTimeIndex
734        """
735        return self._round_via_method(freq, _floor_int)
736
737    def ceil(self, freq):
738        """Round dates up to fixed frequency.
739
740        Parameters
741        ----------
742        freq : str
743            The frequency level to round the index to.  Must be a fixed
744            frequency like 'S' (second) not 'ME' (month end).  See `frequency
745            aliases <https://pandas.pydata.org/pandas-docs/stable/user_guide/timeseries.html#offset-aliases>`_
746            for a list of possible values.
747
748        Returns
749        -------
750        CFTimeIndex
751        """
752        return self._round_via_method(freq, _ceil_int)
753
754    def round(self, freq):
755        """Round dates to a fixed frequency.
756
757        Parameters
758        ----------
759        freq : str
760            The frequency level to round the index to.  Must be a fixed
761            frequency like 'S' (second) not 'ME' (month end).  See `frequency
762            aliases <https://pandas.pydata.org/pandas-docs/stable/user_guide/timeseries.html#offset-aliases>`_
763            for a list of possible values.
764
765        Returns
766        -------
767        CFTimeIndex
768        """
769        return self._round_via_method(freq, _round_to_nearest_half_even)
770
771
772def _parse_iso8601_without_reso(date_type, datetime_str):
773    date, _ = _parse_iso8601_with_reso(date_type, datetime_str)
774    return date
775
776
777def _parse_array_of_cftime_strings(strings, date_type):
778    """Create a numpy array from an array of strings.
779
780    For use in generating dates from strings for use with interp.  Assumes the
781    array is either 0-dimensional or 1-dimensional.
782
783    Parameters
784    ----------
785    strings : array of strings
786        Strings to convert to dates
787    date_type : cftime.datetime type
788        Calendar type to use for dates
789
790    Returns
791    -------
792    np.array
793    """
794    return np.array(
795        [_parse_iso8601_without_reso(date_type, s) for s in strings.ravel()]
796    ).reshape(strings.shape)
797
798
799def _contains_datetime_timedeltas(array):
800    """Check if an input array contains datetime.timedelta objects."""
801    array = np.atleast_1d(array)
802    return isinstance(array[0], timedelta)
803
804
805def _cftimeindex_from_i8(values, date_type, name):
806    """Construct a CFTimeIndex from an array of integers.
807
808    Parameters
809    ----------
810    values : np.array
811        Integers representing microseconds since 1970-01-01.
812    date_type : cftime.datetime
813        Type of date for the index.
814    name : str
815        Name of the index.
816
817    Returns
818    -------
819    CFTimeIndex
820    """
821    epoch = date_type(1970, 1, 1)
822    dates = np.array([epoch + timedelta(microseconds=int(value)) for value in values])
823    return CFTimeIndex(dates, name=name)
824
825
826def _total_microseconds(delta):
827    """Compute the total number of microseconds of a datetime.timedelta.
828
829    Parameters
830    ----------
831    delta : datetime.timedelta
832        Input timedelta.
833
834    Returns
835    -------
836    int
837    """
838    return delta / timedelta(microseconds=1)
839
840
841def _floor_int(values, unit):
842    """Copied from pandas."""
843    return values - np.remainder(values, unit)
844
845
846def _ceil_int(values, unit):
847    """Copied from pandas."""
848    return values + np.remainder(-values, unit)
849
850
851def _round_to_nearest_half_even(values, unit):
852    """Copied from pandas."""
853    if unit % 2:
854        return _ceil_int(values - unit // 2, unit)
855    quotient, remainder = np.divmod(values, unit)
856    mask = np.logical_or(
857        remainder > (unit // 2), np.logical_and(remainder == (unit // 2), quotient % 2)
858    )
859    quotient[mask] += 1
860    return quotient * unit
861