lm_space.py 32.8 KB
Newer Older
csongor's avatar
csongor committed
1 2 3 4 5 6 7 8 9 10 11

from __future__ import division

import os
import numpy as np
import pylab as pl
from matplotlib.colors import LogNorm as ln
from matplotlib.ticker import LogFormatter as lf

from d2o import STRATEGIES as DISTRIBUTION_STRATEGIES

12
from nifty.spaces.space import Space
theos's avatar
theos committed
13

csongor's avatar
csongor committed
14 15 16
from nifty.config import about,\
                         nifty_configuration as gc,\
                         dependency_injector as gdi
theos's avatar
theos committed
17

18
LMSpaceParadict = None
Jait Dixit's avatar
Jait Dixit committed
19
# from nifty.nifty_power_indices import lm_power_indices
csongor's avatar
csongor committed
20 21 22 23 24 25 26 27
from nifty.nifty_random import random

gl = gdi.get('libsharp_wrapper_gl')
hp = gdi.get('healpy')

LM_DISTRIBUTION_STRATEGIES = DISTRIBUTION_STRATEGIES['global']


28
class LMSpace(Space):
csongor's avatar
csongor committed
29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122
    """
        ..       __
        ..     /  /
        ..    /  /    __ ____ ___
        ..   /  /   /   _    _   |
        ..  /  /_  /  / /  / /  /
        ..  \___/ /__/ /__/ /__/  space class

        NIFTY subclass for spherical harmonics components, for representations
        of fields on the two-sphere.

        Parameters
        ----------
        lmax : int
            Maximum :math:`\ell`-value up to which the spherical harmonics
            coefficients are to be used.
        mmax : int, *optional*
            Maximum :math:`m`-value up to which the spherical harmonics
            coefficients are to be used (default: `lmax`).
        dtype : numpy.dtype, *optional*
            Data type of the field values (default: numpy.complex128).

        See Also
        --------
        hp_space : A class for the HEALPix discretization of the sphere [#]_.
        gl_space : A class for the Gauss-Legendre discretization of the
            sphere [#]_.

        Notes
        -----
        Hermitian symmetry, i.e. :math:`a_{\ell -m} = \overline{a}_{\ell m}` is
        always assumed for the spherical harmonics components, i.e. only fields
        on the two-sphere with real-valued representations in position space
        can be handled.

        References
        ----------
        .. [#] K.M. Gorski et al., 2005, "HEALPix: A Framework for
               High-Resolution Discretization and Fast Analysis of Data
               Distributed on the Sphere", *ApJ* 622..759G.
        .. [#] M. Reinecke and D. Sverre Seljebotn, 2013, "Libsharp - spherical
               harmonic transforms revisited";
               `arXiv:1303.4945 <http://www.arxiv.org/abs/1303.4945>`_

        Attributes
        ----------
        para : numpy.ndarray
            One-dimensional array containing the two numbers `lmax` and
            `mmax`.
        dtype : numpy.dtype
            Data type of the field values.
        discrete : bool
            Parameter captioning the fact that an :py:class:`lm_space` is
            always discrete.
        vol : numpy.ndarray
            Pixel volume of the :py:class:`lm_space`, which is always 1.
    """

    def __init__(self, lmax, mmax=None, dtype=np.dtype('complex128')):
        """
            Sets the attributes for an lm_space class instance.

            Parameters
            ----------
            lmax : int
                Maximum :math:`\ell`-value up to which the spherical harmonics
                coefficients are to be used.
            mmax : int, *optional*
                Maximum :math:`m`-value up to which the spherical harmonics
                coefficients are to be used (default: `lmax`).
            dtype : numpy.dtype, *optional*
                Data type of the field values (default: numpy.complex128).

            Returns
            -------
            None.

            Raises
            ------
            ImportError
                If neither the libsharp_wrapper_gl nor the healpy module are
                available.
            ValueError
                If input `nside` is invaild.

        """

        # check imports
        if not gc['use_libsharp'] and not gc['use_healpy']:
            raise ImportError(about._errors.cstring(
                "ERROR: neither libsharp_wrapper_gl nor healpy activated."))

        self._cache_dict = {'check_codomain': {}}

theos's avatar
theos committed
123
        self.paradict = LMSpaceParadict(lmax=lmax, mmax=mmax)
csongor's avatar
csongor committed
124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169

        # check data type
        dtype = np.dtype(dtype)
        if dtype not in [np.dtype('complex64'), np.dtype('complex128')]:
            about.warnings.cprint("WARNING: data type set to complex128.")
            dtype = np.dtype('complex128')
        self.dtype = dtype

        self.harmonic = True
        self.distances = (np.float(1),)

        self.power_indices = lm_power_indices(
                    lmax=self.paradict['lmax'],
                    dim=self.dim,
                    allowed_distribution_strategies=LM_DISTRIBUTION_STRATEGIES)

    @property
    def para(self):
        temp = np.array([self.paradict['lmax'],
                         self.paradict['mmax']], dtype=int)
        return temp

    @para.setter
    def para(self, x):
        self.paradict['lmax'] = x[0]
        self.paradict['mmax'] = x[1]

    def __hash__(self):
        result_hash = 0
        for (key, item) in vars(self).items():
            if key in ['_cache_dict', 'power_indices']:
                continue
            result_hash ^= item.__hash__() * hash(key)
        return result_hash

    def _identifier(self):
        # Extract the identifying parts from the vars(self) dict.
        temp = [(ii[0],
                 ((lambda x: tuple(x) if
                  isinstance(x, np.ndarray) else x)(ii[1])))
                for ii in vars(self).iteritems()
                if ii[0] not in ['_cache_dict', 'power_indices']]
        # Return the sorted identifiers as a tuple.
        return tuple(sorted(temp))

    def copy(self):
170 171 172
        return LMSpace(lmax=self.paradict['lmax'],
                       mmax=self.paradict['mmax'],
                       dtype=self.dtype)
csongor's avatar
csongor committed
173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266

    @property
    def shape(self):
        lmax = self.paradict['lmax']
        mmax = self.paradict['mmax']
        return (np.int((mmax + 1) * (lmax + 1) - ((mmax + 1) * mmax) // 2),)


    @property
    def meta_volume(self):
        """
            Calculates the meta volumes.

            The meta volumes are the volumes associated with each component of
            a field, taking into account field components that are not
            explicitly included in the array of field values but are determined
            by symmetry conditions.

            Parameters
            ----------
            total : bool, *optional*
                Whether to return the total meta volume of the space or the
                individual ones of each field component (default: False).

            Returns
            -------
            mol : {numpy.ndarray, float}
                Meta volume of the field components or the complete space.

            Notes
            -----
            The spherical harmonics components with :math:`m=0` have meta
            volume 1, the ones with :math:`m>0` have meta volume 2, sinnce they
            each determine another component with negative :math:`m`.
        """
        return np.float(self.dof())

    @property
    def meta_volume_split(self):
        mol = self.cast(1, dtype=np.float)
        mol[self.paradict['lmax'] + 1:] = 2  # redundant: (l,m) and (l,-m)
        return mol

    def complement_cast(self, x, axis=None, **kwargs):
        if axis is None:
            lmax = self.paradict['lmax']
            complexity_mask = x[:lmax+1].iscomplex()
            if complexity_mask.any():
                about.warnings.cprint("WARNING: Taking the absolute values for " +
                                      "all complex entries where lmax==0")
                x[:lmax+1] = abs(x[:lmax+1])
        else:
            # TODO hermitianize only on specific axis
            lmax = self.paradict['lmax']
            complexity_mask = x[:lmax+1].iscomplex()
            if complexity_mask.any():
                about.warnings.cprint("WARNING: Taking the absolute values for " +
                                      "all complex entries where lmax==0")
                x[:lmax+1] = abs(x[:lmax+1])
        return x

    # TODO: Extend to binning/log
    def enforce_power(self, spec, size=None, kindex=None):
        if kindex is None:
            kindex_size = self.paradict['lmax'] + 1
            kindex = np.arange(kindex_size,
                               dtype=np.array(self.distances).dtype)
        return self._enforce_power_helper(spec=spec,
                                          size=size,
                                          kindex=kindex)

    def _check_codomain(self, codomain):
        """
            Checks whether a given codomain is compatible to the
            :py:class:`lm_space` or not.

            Parameters
            ----------
            codomain : nifty.space
                Space to be checked for compatibility.

            Returns
            -------
            check : bool
                Whether or not the given codomain is compatible to the space.

            Notes
            -----
            Compatible codomains are instances of :py:class:`lm_space`,
            :py:class:`gl_space`, and :py:class:`hp_space`.
        """
        if codomain is None:
            return False

267 268
        from nifty.spaces.hp_space import HPSpace
        from nifty.spaces.gl_space import GLSpace
csongor's avatar
csongor committed
269 270 271 272
        if not isinstance(codomain, Space):
            raise TypeError(about._errors.cstring(
                "ERROR: The given codomain must be a nifty lm_space."))

273
        elif isinstance(codomain, GLSpace):
csongor's avatar
csongor committed
274 275 276 277 278 279 280 281
            # lmax==mmax
            # nlat==lmax+1
            # nlon==2*lmax+1
            if ((self.paradict['lmax'] == self.paradict['mmax']) and
                    (codomain.paradict['nlat'] == self.paradict['lmax']+1) and
                    (codomain.paradict['nlon'] == 2*self.paradict['lmax']+1)):
                return True

282
        elif isinstance(codomain, HPSpace):
csongor's avatar
csongor committed
283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322
            # lmax==mmax
            # 3*nside-1==lmax
            if ((self.paradict['lmax'] == self.paradict['mmax']) and
                    (3*codomain.paradict['nside']-1 == self.paradict['lmax'])):
                return True

        return False

    def get_codomain(self, coname=None, **kwargs):
        """
            Generates a compatible codomain to which transformations are
            reasonable, i.e.\  a pixelization of the two-sphere.

            Parameters
            ----------
            coname : string, *optional*
                String specifying a desired codomain (default: None).

            Returns
            -------
            codomain : nifty.space
                A compatible codomain.

            Notes
            -----
            Possible arguments for `coname` are ``'gl'`` in which case a Gauss-
            Legendre pixelization [#]_ of the sphere is generated, and ``'hp'``
            in which case a HEALPix pixelization [#]_ is generated.

            References
            ----------
            .. [#] K.M. Gorski et al., 2005, "HEALPix: A Framework for
                   High-Resolution Discretization and Fast Analysis of Data
                   Distributed on the Sphere", *ApJ* 622..759G.
            .. [#] M. Reinecke and D. Sverre Seljebotn, 2013,
                   "Libsharp - spherical
                   harmonic transforms revisited";
                   `arXiv:1303.4945 <http://www.arxiv.org/abs/1303.4945>`_

        """
323 324
        from hp_space import HPSpace
        from gl_space import GLSpace
csongor's avatar
csongor committed
325 326 327 328 329 330 331 332 333
        if coname == 'gl' or (coname is None and gc['lm2gl']):
            if self.dtype == np.dtype('complex64'):
                new_dtype = np.float32
            elif self.dtype == np.dtype('complex128'):
                new_dtype = np.float64
            else:
                raise NotImplementedError
            nlat = self.paradict['lmax'] + 1
            nlon = self.paradict['lmax'] * 2 + 1
334
            return GLSpace(nlat=nlat, nlon=nlon, dtype=new_dtype)
csongor's avatar
csongor committed
335 336 337

        elif coname == 'hp' or (coname is None and not gc['lm2gl']):
            nside = (self.paradict['lmax']+1) // 3
338
            return HPSpace(nside=nside)
csongor's avatar
csongor committed
339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414

        else:
            raise ValueError(about._errors.cstring(
                "ERROR: unsupported or incompatible codomain '"+coname+"'."))

    def get_random_values(self, **kwargs):
        """
            Generates random field values according to the specifications given
            by the parameters, taking into account complex-valuedness and
            hermitian symmetry.

            Returns
            -------
            x : numpy.ndarray
                Valid field values.

            Other parameters
            ----------------
            random : string, *optional*
                Specifies the probability distribution from which the random
                numbers are to be drawn.
                Supported distributions are:

                - "pm1" (uniform distribution over {+1,-1} or {+1,+i,-1,-i}
                - "gau" (normal distribution with zero-mean and a given
                    standard
                    deviation or variance)
                - "syn" (synthesizes from a given power spectrum)
                - "uni" (uniform distribution over [vmin,vmax[)

                (default: None).
            dev : float, *optional*
                Standard deviation (default: 1).
            var : float, *optional*
                Variance, overriding `dev` if both are specified
                (default: 1).
            spec : {scalar, list, numpy.array, nifty.field, function},
                *optional*
                Power spectrum (default: 1).
            vmin : float, *optional*
                Lower limit for a uniform distribution (default: 0).
            vmax : float, *optional*
                Upper limit for a uniform distribution (default: 1).
        """
        arg = random.parse_arguments(self, **kwargs)

#        if arg is None:
#            x = 0
#
#        elif arg['random'] == "pm1":
#            x = random.pm1(dtype=self.dtype, shape=self.shape)
#
#        elif arg['random'] == "gau":
#            x = random.gau(dtype=self.dtype,
#                           shape=self.shape,
#                           mean=arg['mean'],
#                           std=arg['std'])

        if arg['random'] == "syn":
            lmax = self.paradict['lmax']
            mmax = self.paradict['mmax']
            if self.dtype == np.dtype('complex64'):
                if gc['use_libsharp']:
                    sample = gl.synalm_f(arg['spec'], lmax=lmax, mmax=mmax)
                else:
                    sample = hp.synalm(
                                arg['spec'].astype(np.complex128),
                                lmax=lmax, mmax=mmax).astype(np.complex64,
                                                             copy=False)
            else:
                if gc['use_healpy']:
                    sample = hp.synalm(arg['spec'], lmax=lmax, mmax=mmax)
                else:
                    sample = gl.synalm(arg['spec'], lmax=lmax, mmax=mmax)

        else:
415
            sample = super(LMSpace, self).get_random_values(**arg)
csongor's avatar
csongor committed
416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501

#        elif arg['random'] == "uni":
#            x = random.uni(dtype=self.dtype,
#                           shape=self.shape,
#                           vmin=arg['vmin'],
#                           vmax=arg['vmax'])
#
#        else:
#            raise KeyError(about._errors.cstring(
#                "ERROR: unsupported random key '" + str(arg['random']) + "'."))
        sample = self.cast(sample)
        return sample

#    def calc_dot(self, x, y):
#        """
#            Computes the discrete inner product of two given arrays of field
#            values.
#
#            Parameters
#            ----------
#            x : numpy.ndarray
#                First array
#            y : numpy.ndarray
#                Second array
#
#            Returns
#            -------
#            dot : scalar
#                Inner product of the two arrays.
#        """
#        x = self.cast(x)
#        y = self.cast(y)
#
#        lmax = self.paradict['lmax']
#
#        x_low = x[:lmax + 1]
#        x_high = x[lmax + 1:]
#        y_low = y[:lmax + 1]
#        y_high = y[lmax + 1:]
#
#        dot = (x_low.real * y_low.real).sum()
#        dot += 2 * (x_high.real * y_high.real).sum()
#        dot += 2 * (x_high.imag * y_high.imag).sum()
#        return dot

    def dot_contraction(self, x, axes):
        assert len(axes) == 1
        axis = axes[0]
        lmax = self.paradict['lmax']

        # extract the low and high parts of x
        extractor = ()
        extractor += (slice(None),)*axis
        low_extractor = extractor + (slice(None, lmax+1), )
        high_extractor = extractor + (slice(lmax+1), )

        result = x[low_extractor].sum(axes) + 2 * x[high_extractor].sum(axes)
        return result

    def calc_transform(self, x, codomain=None, **kwargs):
        """
            Computes the transform of a given array of field values.

            Parameters
            ----------
            x : numpy.ndarray
                Array to be transformed.
            codomain : nifty.space, *optional*
                codomain space to which the transformation shall map
                (default: self).

            Returns
            -------
            Tx : numpy.ndarray
                Transformed array
        """
        x = self.cast(x)

        if codomain is None:
            codomain = self.get_codomain()

        # Check if the given codomain is suitable for the transformation
        if not self.check_codomain(codomain):
            raise ValueError(about._errors.cstring(
                "ERROR: unsupported codomain."))

502
        # if self.distribution_strategy != 'not':
csongor's avatar
csongor committed
503 504 505 506 507
        #     about.warnings.cprint(
        #         "WARNING: Field data is consolidated to all nodes for "
        #         "external alm2map method!")

        np_x = x.get_full_data()
508 509
        from hp_space import HPSpace
        from gl_space import GLSpace
510
        if isinstance(codomain, GLSpace):
csongor's avatar
csongor committed
511 512 513 514 515 516 517 518 519 520 521 522 523 524
            nlat = codomain.paradict['nlat']
            nlon = codomain.paradict['nlon']
            lmax = self.paradict['lmax']
            mmax = self.paradict['mmax']

            # transform
            if self.dtype == np.dtype('complex64'):
                np_Tx = gl.alm2map_f(np_x, nlat=nlat, nlon=nlon,
                                     lmax=lmax, mmax=mmax, cl=False)
            else:
                np_Tx = gl.alm2map(np_x, nlat=nlat, nlon=nlon,
                                   lmax=lmax, mmax=mmax, cl=False)
            Tx = codomain.cast(np_Tx)

525
        elif isinstance(codomain, HPSpace):
csongor's avatar
csongor committed
526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571
            nside = codomain.paradict['nside']
            lmax = self.paradict['lmax']
            mmax = self.paradict['mmax']

            # transform
            np_x = np_x.astype(np.complex128, copy=False)
            np_Tx = hp.alm2map(np_x, nside, lmax=lmax,
                               mmax=mmax, pixwin=False, fwhm=0.0, sigma=None,
                               pol=True, inplace=False)
            Tx = codomain.cast(np_Tx)

        else:
            raise ValueError(about._errors.cstring(
                "ERROR: unsupported transformation."))

        return codomain.cast(Tx)

    def calc_smooth(self, x, sigma=0, **kwargs):
        """
            Smoothes an array of field values by convolution with a Gaussian
            kernel in position space.

            Parameters
            ----------
            x : numpy.ndarray
                Array of field values to be smoothed.
            sigma : float, *optional*
                Standard deviation of the Gaussian kernel, specified in units
                of length in position space; for testing: a sigma of -1 will be
                reset to a reasonable value (default: 0).

            Returns
            -------
            Gx : numpy.ndarray
                Smoothed array.
        """
        x = self.cast(x)
        # check sigma
        if sigma == 0:
            return self.unary_operation(x, op='copy')
        elif sigma == -1:
            about.infos.cprint("INFO: invalid sigma reset.")
            sigma = np.sqrt(2) * np.pi / (self.paradict['lmax'] + 1)
        elif sigma < 0:
            raise ValueError(about._errors.cstring("ERROR: invalid sigma."))

572
        # if self.distribution_strategy != 'not':
csongor's avatar
csongor committed
573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614
        #     about.warnings.cprint(
        #         "WARNING: Field data is consolidated to all nodes for "
        #         "external smoothalm method!")

        np_x = x.get_full_data()

        if gc['use_healpy']:
            np_smoothed_x = hp.smoothalm(np_x,
                                         fwhm=0.0,
                                         sigma=sigma,
                                         pol=True,
                                         mmax=self.paradict['mmax'],
                                         verbose=False,
                                         inplace=False)
        else:
            np_smoothed_x = gl.smoothalm(np_x,
                                         lmax=self.paradict['lmax'],
                                         mmax=self.paradict['mmax'],
                                         fwhm=0.0,
                                         sigma=sigma,
                                         overwrite=False)
        return self.cast(np_smoothed_x)

    def calc_power(self, x, **kwargs):
        """
            Computes the power of an array of field values.

            Parameters
            ----------
            x : numpy.ndarray
                Array containing the field values of which the power is to be
                calculated.

            Returns
            -------
            spec : numpy.ndarray
                Power contained in the input array.
        """
        x = self.cast(x)
        lmax = self.paradict['lmax']
        mmax = self.paradict['mmax']

615
        # if self.distribution_strategy != 'not':
csongor's avatar
csongor committed
616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702
        #     about.warnings.cprint(
        #         "WARNING: Field data is consolidated to all nodes for "
        #         "external anaalm/alm2cl method!")

        np_x = x.get_full_data()

        # power spectrum
        if self.dtype == np.dtype('complex64'):
            if gc['use_libsharp']:
                result = gl.anaalm_f(np_x, lmax=lmax, mmax=mmax)
            else:
                np_x = np_x.astype(np.complex128, copy=False)
                result = hp.alm2cl(np_x,
                                   alms2=None,
                                   lmax=lmax,
                                   mmax=mmax,
                                   lmax_out=lmax,
                                   nspec=None)
        else:
            if gc['use_healpy']:
                result = hp.alm2cl(np_x,
                                   alms2=None,
                                   lmax=lmax,
                                   mmax=mmax,
                                   lmax_out=lmax,
                                   nspec=None)
            else:
                result = gl.anaalm(np_x,
                                   lmax=lmax,
                                   mmax=mmax)

        if self.dtype == np.dtype('complex64'):
            result = result.astype(np.float32, copy=False)
        elif self.dtype == np.dtype('complex128'):
            result = result.astype(np.float64, copy=False)
        else:
            raise NotImplementedError(about._errors.cstring(
                "ERROR: dtype %s not known to calc_power method." %
                str(self.dtype)))

    def get_plot(self, x, title="", vmin=None, vmax=None, power=True,
                 norm=None, cmap=None, cbar=True, other=None, legend=False,
                 mono=True, **kwargs):
        """
            Creates a plot of field values according to the specifications
            given by the parameters.

            Parameters
            ----------
            x : numpy.ndarray
                Array containing the field values.

            Returns
            -------
            None

            Other parameters
            ----------------
            title : string, *optional*
                Title of the plot (default: "").
            vmin : float, *optional*
                Minimum value to be displayed (default: ``min(x)``).
            vmax : float, *optional*
                Maximum value to be displayed (default: ``max(x)``).
            power : bool, *optional*
                Whether to plot the power contained in the field or the field
                values themselves (default: True).
            unit : string, *optional*
                Unit of the field values (default: "").
            norm : string, *optional*
                Scaling of the field values before plotting (default: None).
            cmap : matplotlib.colors.LinearSegmentedColormap, *optional*
                Color map to be used for two-dimensional plots (default: None).
            cbar : bool, *optional*
                Whether to show the color bar or not (default: True).
            other : {single object, tuple of objects}, *optional*
                Object or tuple of objects to be added, where objects can be
                scalars, arrays, or fields (default: None).
            legend : bool, *optional*
                Whether to show the legend or not (default: False).
            mono : bool, *optional*
                Whether to plot the monopole or not (default: True).
            save : string, *optional*
                Valid file name where the figure is to be stored, by default
                the figure is not saved (default: False).

        """
703 704
        from nifty.field import Field

csongor's avatar
csongor committed
705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736
        try:
            x = x.get_full_data()
        except AttributeError:
            pass

        if(not pl.isinteractive())and(not bool(kwargs.get("save", False))):
            about.warnings.cprint("WARNING: interactive mode off.")

        if(power):
            x = self.calc_power(x)

            fig = pl.figure(num=None, figsize=(6.4, 4.8), dpi=None, facecolor="none",
                            edgecolor="none", frameon=False, FigureClass=pl.Figure)
            ax0 = fig.add_axes([0.12, 0.12, 0.82, 0.76])

            xaxes = np.arange(self.para[0] + 1, dtype=np.int)
            if(vmin is None):
                vmin = np.min(x[:mono].tolist(
                ) + (xaxes * (2 * xaxes + 1) * x)[1:].tolist(), axis=None, out=None)
            if(vmax is None):
                vmax = np.max(x[:mono].tolist(
                ) + (xaxes * (2 * xaxes + 1) * x)[1:].tolist(), axis=None, out=None)
            ax0.loglog(xaxes[1:], (xaxes * (2 * xaxes + 1) * x)[1:], color=[0.0,
                                                                            0.5, 0.0], label="graph 0", linestyle='-', linewidth=2.0, zorder=1)
            if(mono):
                ax0.scatter(0.5 * (xaxes[1] + xaxes[2]), x[0], s=20, color=[0.0, 0.5, 0.0], marker='o',
                            cmap=None, norm=None, vmin=None, vmax=None, alpha=None, linewidths=None, verts=None, zorder=1)

            if(other is not None):
                if(isinstance(other, tuple)):
                    other = list(other)
                    for ii in xrange(len(other)):
737
                        if(isinstance(other[ii], Field)):
csongor's avatar
csongor committed
738 739 740
                            other[ii] = other[ii].power(**kwargs)
                        else:
                            other[ii] = self.enforce_power(other[ii])
741
                elif(isinstance(other, Field)):
csongor's avatar
csongor committed
742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851
                    other = [other.power(**kwargs)]
                else:
                    other = [self.enforce_power(other)]
                imax = max(1, len(other) - 1)
                for ii in xrange(len(other)):
                    ax0.loglog(xaxes[1:], (xaxes * (2 * xaxes + 1) * other[ii])[1:], color=[max(0.0, 1.0 - (2 * ii / imax)**2), 0.5 * ((2 * ii - imax) / imax)
                                                                                            ** 2, max(0.0, 1.0 - (2 * (ii - imax) / imax)**2)], label="graph " + str(ii + 1), linestyle='-', linewidth=1.0, zorder=-ii)
                    if(mono):
                        ax0.scatter(0.5 * (xaxes[1] + xaxes[2]), other[ii][0], s=20, color=[max(0.0, 1.0 - (2 * ii / imax)**2), 0.5 * ((2 * ii - imax) / imax)**2, max(
                            0.0, 1.0 - (2 * (ii - imax) / imax)**2)], marker='o', cmap=None, norm=None, vmin=None, vmax=None, alpha=None, linewidths=None, verts=None, zorder=-ii)
                if(legend):
                    ax0.legend()

            ax0.set_xlim(xaxes[1], xaxes[-1])
            ax0.set_xlabel(r"$\ell$")
            ax0.set_ylim(vmin, vmax)
            ax0.set_ylabel(r"$\ell(2\ell+1) C_\ell$")
            ax0.set_title(title)

        else:
            if(np.iscomplexobj(x)):
                if(title):
                    title += " "
                if(bool(kwargs.get("save", False))):
                    save_ = os.path.splitext(
                        os.path.basename(str(kwargs.get("save"))))
                    kwargs.update(save=save_[0] + "_absolute" + save_[1])
                self.get_plot(np.absolute(x), title=title + "(absolute)", vmin=vmin, vmax=vmax,
                              power=False, norm=norm, cmap=cmap, cbar=cbar, other=None, legend=False, **kwargs)
#                self.get_plot(np.real(x),title=title+"(real part)",vmin=vmin,vmax=vmax,power=False,norm=norm,cmap=cmap,cbar=cbar,other=None,legend=False,**kwargs)
#                self.get_plot(np.imag(x),title=title+"(imaginary part)",vmin=vmin,vmax=vmax,power=False,norm=norm,cmap=cmap,cbar=cbar,other=None,legend=False,**kwargs)
                if(cmap is None):
                    cmap = pl.cm.hsv_r
                if(bool(kwargs.get("save", False))):
                    kwargs.update(save=save_[0] + "_phase" + save_[1])
                self.get_plot(np.angle(x, deg=False), title=title + "(phase)", vmin=-3.1416, vmax=3.1416, power=False,
                              norm=None, cmap=cmap, cbar=cbar, other=None, legend=False, **kwargs)  # values in [-pi,pi]
                return None  # leave method
            else:
                if(vmin is None):
                    vmin = np.min(x, axis=None, out=None)
                if(vmax is None):
                    vmax = np.max(x, axis=None, out=None)
                if(norm == "log")and(vmin <= 0):
                    raise ValueError(about._errors.cstring(
                        "ERROR: nonpositive value(s)."))

                # not a number
                xmesh = np.nan * \
                    np.empty(self.para[::-1] + 1, dtype=np.float16, order='C')
                xmesh[4, 1] = None
                xmesh[1, 4] = None
                lm = 0
                for mm in xrange(self.para[1] + 1):
                    xmesh[mm][mm:] = x[lm:lm + self.para[0] + 1 - mm]
                    lm += self.para[0] + 1 - mm

                s_ = np.array([1, self.para[1] / self.para[0]
                               * (1.0 + 0.159 * bool(cbar))])
                fig = pl.figure(num=None, figsize=(
                    6.4 * s_[0], 6.4 * s_[1]), dpi=None, facecolor="none", edgecolor="none", frameon=False, FigureClass=pl.Figure)
                ax0 = fig.add_axes(
                    [0.06 / s_[0], 0.06 / s_[1], 1.0 - 0.12 / s_[0], 1.0 - 0.12 / s_[1]])
                ax0.set_axis_bgcolor([0.0, 0.0, 0.0, 0.0])

                xaxes = np.arange(self.para[0] + 2, dtype=np.int) - 0.5
                yaxes = np.arange(self.para[1] + 2, dtype=np.int) - 0.5
                if(norm == "log"):
                    n_ = ln(vmin=vmin, vmax=vmax)
                else:
                    n_ = None
                sub = ax0.pcolormesh(xaxes, yaxes, np.ma.masked_where(np.isnan(
                    xmesh), xmesh), cmap=cmap, norm=n_, vmin=vmin, vmax=vmax, clim=(vmin, vmax))
                ax0.set_xlim(xaxes[0], xaxes[-1])
                ax0.set_xticks([0], minor=False)
                ax0.set_xlabel(r"$\ell$")
                ax0.set_ylim(yaxes[0], yaxes[-1])
                ax0.set_yticks([0], minor=False)
                ax0.set_ylabel(r"$m$")
                ax0.set_aspect("equal")
                if(cbar):
                    if(norm == "log"):
                        f_ = lf(10, labelOnlyBase=False)
                        b_ = sub.norm.inverse(
                            np.linspace(0, 1, sub.cmap.N + 1))
                        v_ = np.linspace(
                            sub.norm.vmin, sub.norm.vmax, sub.cmap.N)
                    else:
                        f_ = None
                        b_ = None
                        v_ = None
                    fig.colorbar(sub, ax=ax0, orientation="horizontal", fraction=0.1, pad=0.05, shrink=0.75, aspect=20, ticks=[
                                 vmin, vmax], format=f_, drawedges=False, boundaries=b_, values=v_)
                ax0.set_title(title)

        if(bool(kwargs.get("save", False))):
            fig.savefig(str(kwargs.get("save")), dpi=None, facecolor="none", edgecolor="none", orientation="portrait",
                        papertype=None, format=None, transparent=False, bbox_inches=None, pad_inches=0.1)
            pl.close(fig)
        else:
            fig.canvas.draw()

    def getlm(self):  # > compute all (l,m)
        index = np.arange(self.dim)
        n = 2 * self.paradict['lmax'] + 1
        m = np.ceil(
            (n - np.sqrt(n**2 - 8 * (index - self.paradict['lmax']))) / 2
                    ).astype(np.int)
        l = index - self.paradict['lmax'] * m + m * (m - 1) // 2
        return l, m