Support for G_#.TAG files

class magic.MagicGraph(ivar=None, datadir='.', quiet=True, ave=False, tag=None, precision=<class 'numpy.float32'>)[source]

This class allows to read the 3-D graphic outputs of the MagIC code (G_#.TAG and G_ave.TAG) files. Those are binary unformatted outputs, there are therefore two ways to load them:

  • If buildLib=True in magic.cfg and the fortran libraries were correctly built, then the reader uses a fortran program that is expected to be much faster than the pure python routine.

  • If buildLib=False, then a pure python program is used to read the G files.

>>> # Regular G files
>>> gr = MagicGraph(ivar=1, tag='N0m2a')
>>> print(gr.vr.shape) # shape of vr
>>> print(gr.ek) # print ekman number
>>> print(gr.minc) # azimuthal symmetry
>>> # Averaged G file with double precision
>>> gr = MagicGraph(ave=True, tag='N0m2', precision=np.float64)
__init__(ivar=None, datadir='.', quiet=True, ave=False, tag=None, precision=<class 'numpy.float32'>)[source]
Parameters:
  • ave (bool) – when set to True, it tries to find an average G file (G_ave.TAG)

  • ivar (int) – the number of the G file

  • tag (str) – extension TAG of the G file. If not specified, the most recent G_#.TAG file found in the directory will be selected.

  • quiet (bool) – when set to True, makes the output silent

  • datadir (str) – directory of the G file (default is . )

  • precision (str) – single or double precision (default np.float32)

read_record_marker(filename, endian, quiet=True)[source]

This function is used to read a Graphic file that contains record markers.

Parameters:
  • filename (str) – name of the graphic file

  • endian (str) – endianness of the file

  • quiet (bool) – when set to True, makes the output silent

read_stream(filename, endian)[source]

This function is used to read a Graphic file that has no record marker.

Parameters:
  • filename (str) – name of the graphic file

  • endian (str) – endianness of the file

rearangeLat(field)[source]

This function is used to unfold the colatitudes

Parameters:

field (numpy.ndarray) – input array with MagIC ordering of colatitudes (i.e. successively Northern Hemisphere and Southern Hemisphere)

Returns:

an array with the regular ordering of the colatitudes

Return type:

numpy.ndarray

class magic.Surf(ivar=None, datadir='.', vort=False, ave=False, tag=None, precision=<class 'numpy.float32'>)[source]

This class allows to display the content of a graphic file (G_#.TAG or G_ave.TAG). It allows to plot radial, azimuthal and equatorial cuts as well as phi-averages.

>>> # To read G_1.test
>>> s = Surf(ivar=1, ave=False, tag='test')
>>> # To read the latest G file in the working directory (double precision)
>>> s = Surf(precision=np.float64)
>>> # Possible plots
>>> s.equat(field='vr')
>>> s.avg(field='vp')
>>> s.surf(field='entropy', r=0.8)
>>> s.slice(field='Br', lon_0=[0, 30])
__init__(ivar=None, datadir='.', vort=False, ave=False, tag=None, precision=<class 'numpy.float32'>)[source]
Parameters:
  • ivar (int) – index of the graphic file

  • ave (bool) – when set to True, it tries to read a time-averaged graphic file

  • tag (str) – TAG suffix extension of the graphic file

  • vort (bool) – a boolean to specify whether one wants to compute the 3-D vorticiy components (take care of the memory imprint)

  • datadir (str) – the working directory

  • precision (str) – the storage precision of the graphic file (single or double precision). Default is np.float32 (single)

__weakref__

list of weak references to the object (if defined)

avg(field='vphi', levels=65, cm=None, normed=None, vmax=None, vmin=None, cbar=True, tit=True, pol=False, tor=False, mer=False, merLevels=16, polLevels=16, ic=False, lines=False)[source]

Plot the azimutal average of a given field.

>>> s = Surf()
>>> # Axisymmetric zonal flows, 65 contour levels
>>> s.avg(field='vp', levels=65, cm='seismic')
>>> # Minimal plot (no cbar, not title)
>>> s.avg(field='Br', tit=False, cbar=False)
>>> # Axisymmetric Bphi + poloidal field lines
>>> s.avg(field='Bp', pol=True, polLevels=8)
>>> # Omega-effect, contours truncated from -1e3 to 1e3
>>> s.avg(field='omeffect', vmax=1e3, vmin=-1e3)
Parameters:
  • field (str) – the field you want to display

  • levels (int) – the number of levels in the contourf plot

  • cm (str) – name of the colormap (‘jet’, ‘seismic’, ‘RdYlBu_r’, etc.)

  • tit (bool) – display the title of the figure when set to True

  • cbar (bool) – display the colorbar when set to True

  • vmax (float) – maximum value of the contour levels

  • vmin (float) – minimum value of the contour levels

  • normed (bool) – when set to True, the colormap is centered around zero. Default is True, except for entropy/temperature plots.

  • pol (bool) – diplay the poloidal field lines contours when set to True

  • tor (bool) – diplay the toroidal axisymmetric field contours when set to True

  • mer (bool) – display the meridional circulation contours when set to True

  • merLevels (int) – number of contour levels to display meridional circulation

  • polLevels (int) – number of contour levels to display poloidal field lines

  • ic (bool) – when set to True, also display the contour levels in the inner core

  • lines (bool) – when set to True, over-plot solid lines to highlight the limits between two adjacent contour levels

equat(field='vr', levels=65, cm=None, normed=None, vmax=None, vmin=None, cbar=True, tit=True, avg=False, normRad=False, ic=False)[source]

Plot the equatorial cut of a given field

>>> s = Surf()
>>> # Equatorial cut of the z-vorticity, 65 contour levels
>>> s.equat(field='vortz', levels=65, cm='seismic')
>>> # Minimal plot (no cbar, not title)
>>> s.equat(field='bphi', tit=False, cbar=False)
>>> # Control the limit of the colormap from -1e3 to 1e3
>>> s.equat(field='vr', vmin=-1e3, vmax=1e3, levels=33)
>>> # Normalise the contour levels radius by radius
>>> s.equat(field='jphi', normRad=True)
Parameters:
  • field (str) – the name of the input physical quantity you want to display

  • avg (bool) – when set to True, an additional figure which shows the radial profile of the input physical quantity (azimuthal average) is also displayed

  • normRad (bool) – when set to True, the contour levels are normalised radius by radius (default is False)

  • levels (int) – the number of levels in the contour

  • cm (str) – name of the colormap (‘jet’, ‘seismic’, ‘RdYlBu_r’, etc.)

  • tit (bool) – display the title of the figure when set to True

  • cbar (bool) – display the colorbar when set to True

  • vmax (float) – maximum value of the contour levels

  • vmin (float) – minimum value of the contour levels

  • normed (bool) – when set to True, the colormap is centered around zero. Default is True, except for entropy/temperature plots.

  • ic (bool) – when set to True, also display the contour levels in the inner core

slice(field='Bphi', lon_0=0.0, levels=65, cm=None, normed=None, vmin=None, vmax=None, cbar=True, tit=True, grid=False, nGridLevs=16, normRad=False, ic=False)[source]

Plot an azimuthal slice of a given field.

>>> s = Surf()
>>> # vphi at 0, 30, 60 degrees in longitude
>>> s.slice(field='vp', lon_0=[0, 30, 60], levels=65, cm='seismic')
>>> # Minimal plot (no cbar, not title)
>>> s.avg(field='vp', lon_0=32, tit=False, cbar=False)
>>> # Axisymmetric Bphi + poloidal field lines
>>> s.avg(field='Bp', pol=True, polLevels=8)
>>> # Omega-effect, contours truncated from -1e3 to 1e3
>>> s.avg(field='omeffect', vmax=1e3, vmin=-1e3)
Parameters:
  • field (str) – the field you want to display

  • lon_0 (float or list) – the longitude of the slice in degrees, or a list of longitudes

  • levels (int) – the number of levels in the contourf plot

  • cm (str) – name of the colormap (‘jet’, ‘seismic’, ‘RdYlBu_r’, etc.)

  • tit (bool) – display the title of the figure when set to True

  • cbar (bool) – display the colorbar when set to True

  • vmax (float) – maximum value of the contour levels

  • vmin (float) – minimum value of the contour levels

  • grid (bool) – display or hide the grid

  • nGridLevs (int) – number of grid levels

  • normRad (bool) – when set to True, the contour levels are normalised radius by radius (default is False)

  • ic (bool) – when set to True, also display the contour levels in the inner core

surf(field='Bphi', proj='hammer', lon_0=0.0, r=0.85, vmax=None, vmin=None, lat_0=30.0, levels=65, cm=None, ic=False, lon_shift=0, normed=None, cbar=True, tit=True, lines=False)[source]

Plot the surface distribution of an input field at a given input radius (normalised by the outer boundary radius).

>>> s = Surf()
>>> # Radial flow component at ``r=0.95 r_o``, 65 contour levels
>>> s.surf(field='vr', r=0.95, levels=65, cm='seismic')
>>> # Minimal plot (no cbar, not title)
>>> s.surf(field='entropyfluct', r=0.6, tit=False, cbar=False)
>>> # Control the limit of the colormap from -1e3 to 1e3
>>> s.surf(field='vp', r=1., vmin=-1e3, vmax=1e3, levels=33)
>>> # If basemap is installed, additional projections are available
>>> s.surf(field='Br', r=0.95, proj='ortho', lat_0=45, lon_0=45)
Parameters:
  • field (str) – the name of the field you want to display

  • proj (str) – the type of projection. Default is Hammer, in case you want to use ‘ortho’ or ‘moll’, then Basemap is required.

  • r (float) – the radius at which you want to display the input data (in normalised units with the radius of the outer boundary)

  • levels (int) – the number of levels in the contour

  • cm (str) – name of the colormap (‘jet’, ‘seismic’, ‘RdYlBu_r’, etc.)

  • lon_shift (int) – translate map in azimuth (in degrees)

  • lon_0 (float) – central azimuth (only used with Basemap)

  • lat_0 (float) – central latitude (only used with Basemap)

  • tit (bool) – display the title of the figure when set to True

  • tit – display the title of the figure when set to True

  • cbar (bool) – display the colorbar when set to True

  • lines (bool) – when set to True, over-plot solid lines to highlight the limits between two adjacent contour levels

  • vmax (float) – maximum value of the contour levels

  • vmin (float) – minimum value of the contour levels

  • normed (bool) – when set to True, the colormap is centered around zero. Default is True, except for entropy/temperature plots.

  • lines – when set to True, over-plot solid lines to highlight the limits between two adjacent contour levels

Support for checkpoint_#.TAG files

magic.checkpoint.Graph2Rst(gr, filename='checkpoint_ave')[source]

This function allows to transform an input Graphic file into a checkpoint file format that can be read by MagIC to restart a simulation.

>>> # Load a Graphic File
>>> gr = MagicGraph()
>>> # Produce the file checkpoint_ave.from_G
>>> Graph2Rst(gr, filename='checkpoint_ave.from_G')
Parameters:
  • gr (magic.MagicGraph) – the input graphic file one wants to convert into a restart file

  • filename (str) – name of the checkpoint file

class magic.checkpoint.MagicCheckpoint(l_read=True, filename=None, endian='l')[source]

This class allows to manipulate checkpoint files produced by MagIC. It can read it as

>>> chk = MagicCheckpoint(filename='checkpoint_end.test')
>>> print(chk.wpol.shape, chk.l_max)

This class can also be used to intepolate from FD to Cheb or the opposite >>> chk.cheb2fd(96) >>> chk.write(‘checkpoint_fd.test’)

One can also transform a Graphic file into a checkpoint >>> gr = MagicGraph() >>> chk = MagicCheckpoint(l_read=False) >>> chk.graph2rst(gr)

Finally one can convert checkpoints from XSHELLS >>> chk = MagicCheckpoint(l_read=False) >>> chk.xshells2magic(‘st0’, 161, rscheme=’cheb’, cond_state=’deltaT’)

__init__(l_read=True, filename=None, endian='l')[source]
Parameters:
  • l_read (bool) – a boolean to decide whether one reads a checkpoint or not

  • filename (str) – name of the checkpoint file to be read

__weakref__

list of weak references to the object (if defined)

cheb2fd(n_r_max, fd_stretch=0.3, fd_ratio=0.1)[source]

This routine is used to convert a checkpoint that has a Gauss-Lobatto grid into a finite-difference grid.

Parameters:
  • n_r_max (int) – number of radial grid points of the finite difference grid

  • fd_stretch (float) – stretching of the radial grid

  • fd_ratio (float) – ratio of smallest to largest grid spacing

fd2cheb(n_r_max)[source]

This routine is used to convert a checkpoint that has finite differences in radius into a Gauss-Lobatto grid.

Parameters:

n_r_max (int) – number of radial grid points of the Gauss-Lobatto grid

graph2rst(gr, filename='checkpoint_ave.from_chk')[source]
Parameters:
  • gr (magic.MagicGraph) – the input graphic file one wants to convert into a restart file

  • filename (str) – name of the checkpoint file

read(filename, endian='l')[source]

This routine is used to read a checkpoint file.

Parameters:

filename (str) – name of the checkpoint file

write(filename)[source]

This routine is used to store a checkpoint file. It only stores the state vector not the past quantities required to restart a multistep scheme.

Parameters:

filename (str) – name of the checkpoint file

xshells2magic(xsh_trailing, n_r_max, rscheme='cheb', cond_state='deltaT', scale_b=1.0, filename='checkpoint_end.from_xhells')[source]

This routine is used to convert XSHELLS field[U,B,T].xsh_trailing files into a MagIC checkpoint file.

>>> chk = MagicCheckPoint()
>>> # Convert field[U,T,B].st1ns_hr2 into a MagIC checkpoint file
>>> chk.xshells2magic('st1ns_hr2', 512, rscheme='fd', cond_state='mixed',
                      scale_b=4.472136e-4)
Parameters:
  • xsh_trailing (str) – trailing of the field[U,B,T].xsh_trailing files

  • n_r_max (int) – number of radial grid points to be used

  • rscheme (str) – the type of radial scheme (‘cheb’ or ‘fd’)

  • cond_state (str) – the type of conducting state: - ‘deltaT’: fixed temperature contrast - ‘mixed’: hybrid forcing (STEP1-2 like)

  • scale_b (float) – a rescaling factor for the magnetic field

magic.checkpoint.get_map(lm_max, lmax, mmin, mmax, minc)[source]

This routine determines the look-up tables to convert the indices (l, m) to the single index lm.

Parameters:
  • lm_max (int) – total number of lm combinations.

  • lmax (int) – maximum spherical harmonic degree

  • mmin (int) – minimum spherical harmonic order

  • mmax (int) – maximum spherical harmonic order

  • minc (int) – azimuthal symmetry

Returns:

returns a list of three look-up tables: idx, lm2l, lm2m

Return type:

list

magic.checkpoint.get_truncation(n_theta_max, nalias, minc)[source]

This routine determines l_max, m_max and lm_max from the values of n_theta_max, minc and nalias.

Parameters:
  • n_theta_max (int) – number of points along the colatitude

  • nalias (int) – dealiasing paramete (20 is fully dealiased)

  • minc (int) – azimuthal symmetry

Returns:

returns a list of three integers: l_max, m_max and lm_max

Return type:

list

magic.checkpoint.interp_one_field(field, rold, rnew, rfac=None)[source]

This routine interpolates a complex input field from an old radial grid to a new one.

Parameters:
  • field (numpy.ndarray) – the field to be interpolated

  • rold (numpy.ndarray) – the old radial grid points

  • rnew (numpy.ndarray) – the new radial grid points

  • rfac (numpy.ndarray) – a rescaling function that depends on the radius

Returns:

the field interpolated on the new radial grid

Return type:

numpy.ndarray

Support for movie files *_mov.TAG

class magic.Movie(file=None, iplot=True, step=1, png=False, lastvar=None, nvar='all', levels=12, cm='RdYlBu_r', cut=0.5, bgcolor=None, fluct=False, normed=False, avg=False, std=False, dpi=80, normRad=False, precision=<class 'numpy.float32'>, deminc=True, ifield=0, centeredCm=True, datadir='.')[source]

This class allows to read the movie files generated by the MagIC code.

>>> m = Movie()
>>> # This returns a list of the available movies in the directory
>>> # and lets you decide which one you want to read
>>> # Reads and display AV_mov.test
>>> m = Movie(file='AV_mov.test')
>>> print(m.data) # access to the data
>>> # Read three movie files (no display)
>>> m1 = Movie(file='AV_mov.testa', iplot=False)
>>> m2 = Movie(file='AV_mov.testb', iplot=False)
>>> m3 = Movie(file='AV_mov.testc', iplot=False)
>>> # Stack them together
>>> m = m1+m2+m3
>>> # Display
>>> m.plot(levels=33, cm='seismic', cut=0.5)
>>> # Store the outputs in movie/img_#.png
>>> # Only from the timesteps 280 to 380
>>> m = Movie(file='AB_mov.test', png=True, nvar=100, lastvar=380)
__add__(new)[source]

Built-in function to sum two movies. In case, the spatial grid have been changed an interpolation onto the new grid is used.

Parameters:

new (magic.Movie) – the new movie file to be added

__init__(file=None, iplot=True, step=1, png=False, lastvar=None, nvar='all', levels=12, cm='RdYlBu_r', cut=0.5, bgcolor=None, fluct=False, normed=False, avg=False, std=False, dpi=80, normRad=False, precision=<class 'numpy.float32'>, deminc=True, ifield=0, centeredCm=True, datadir='.')[source]
Parameters:
  • nvar (int) – the number of timesteps of the movie file we want to plot starting from the last line

  • png (bool) – if png=True, write the png files instead of display

  • iplot (bool) – if iplot=True, display otherwise just read

  • lastvar (int) – the number of the last timesteps to be read

  • step (int) – the stepping between two timesteps

  • levels (int) – the number of contour levels

  • cm (str) – the name of the color map

  • fluct (bool) – if fluct=True, substract the axisymmetric part

  • normed (bool) – the colormap is rescaled every timestep when set to True, otherwise it is calculated from the global extrema

  • avg (bool) – if avg=True, time-average is displayed

  • centeredCm (bool) – when set to True, the colormap is centered between -vmax and vmax

  • std (bool) – if std=True, standard deviation is displayed

  • dpi (int) – dot per inch when saving PNGs

  • normRad (bool) – if normRad=True, then we normalise for each radial level

  • precision (str) – precision of the input file, np.float32 for single precision, np.float64 for double precision

  • cut (float) – adjust the contour extrema to max(abs(data))*cut

  • bgcolor (str) – background color of the figure

  • deminc (bool) – a logical to indicate if one wants do get rid of the possible azimuthal symmetry

  • ifield (int) – in case of a multiple-field movie file, you can change the default field displayed using the parameter ifield

  • datadir (str) – working directory

__weakref__

list of weak references to the object (if defined)

avgStd(ifield=0, std=False, cut=0.5, centeredCm=True, levels=12, cmap='RdYlBu_r', ic=False)[source]

Plot time-average or standard deviation

Parameters:
  • ifield (int) – in case of a multiple-field movie file, you can change the default field displayed using the parameter ifield

  • std (bool) – the standard deviation is computed instead the average when std is True

  • levels (int) – number of contour levels

  • cmap (str) – name of the colormap

  • cut (float) – adjust the contour extrema to max(abs(data))*cut

  • centeredCm (bool) – when set to True, the colormap is centered between -vmax and vmax

plot(ifield=0, cut=0.5, centeredCm=True, levels=12, cmap='RdYlBu_r', png=False, step=1, normed=False, dpi=80, bgcolor=None, deminc=True, ic=False)[source]

Plotting function (it can also write the png files)

Parameters:
  • ifield (int) – in case of a multiple-field movie file, you can change the default field displayed using the parameter ifield

  • levels (int) – number of contour levels

  • cmap (str) – name of the colormap

  • cut (float) – adjust the contour extrema to max(abs(data))*cut

  • png (bool) – save the movie as a series of png files when set to True

  • dpi (int) – dot per inch when saving PNGs

  • bgcolor (str) – background color of the figure

  • normed (bool) – the colormap is rescaled every timestep when set to True, otherwise it is calculated from the global extrema

  • step (int) – the stepping between two timesteps

  • deminc (bool) – a logical to indicate if one wants do get rid of the possible azimuthal symmetry

  • centeredCm (bool) – when set to True, the colormap is centered between -vmax and vmax

timeLongitude(ifield=0, removeMean=True, lat0=0.0, levels=12, cm='RdYlBu_r', deminc=True)[source]

Plot the time-longitude diagram (input latitude can be chosen)

Parameters:
  • ifield (int) – in case of a multiple-field movie file, you can change the default field displayed using the parameter ifield

  • lat0 (float) – value of the latitude

  • levels (int) – number of contour levels

  • cm (str) – name of the colormap

  • deminc (bool) – a logical to indicate if one wants do get rid of the possible azimuthal symmetry

  • removeMean (bool) – remove the time-averaged part when set to True

class magic.Movie3D(file=None, step=1, lastvar=None, nvar='all', nrout=48, ratio_out=2.0, potExtra=False, precision=<class 'numpy.float32'>)[source]

This class allows to read the 3D movie files (B|V)_3D_.TAG and transform them into a series of VTS files ./vtsFiles/B3D_#.TAG that can be further read using paraview.

>>> Movie3D(file='B_3D.TAG')
__init__(file=None, step=1, lastvar=None, nvar='all', nrout=48, ratio_out=2.0, potExtra=False, precision=<class 'numpy.float32'>)[source]
Parameters:
  • file (str) – file name

  • nvar (int) – the number of timesteps of the movie file we want to plot starting from the last line

  • lastvar (int) – the number of the last timestep to be read

  • step (int) – the stepping between two timesteps

  • precision (str) – precision of the input file, np.float32 for single precision, np.float64 for double precision

  • potExtra (bool) – when set to True, potential extrapolation of the magnetic field outside the fluid domain is also computed

  • ratio_out (float) – ratio of desired external radius to the CMB radius. This is is only used when potExtra=True

  • nrout (int) – number of additional radial grid points to compute the potential extrapolation. This is only used when potExtra=True

__weakref__

list of weak references to the object (if defined)

Support for B_cmb_coeff.TAG and (V|B)_coeff_r#.TAG files

class magic.coeff.MagicCoeffCmb(tag=None, datadir='.', ratio_cmb_surface=1, scale_b=1, iplot=True, lCut=None, precision=<class 'numpy.float64'>, ave=False, sv=False, quiet=False)[source]

This class allows to read the B_coeff_cmb.TAG files. It first read the poloidal potential at the CMB and then transform it to the Gauss coefficients \(g_{\ell m}\) and \(h_{\ell m}\) using the getGauss function.

>>> # Reads the files B_coeff_cmb.testa, B_coeff_cmb.testb
>>> # and B_coeff_cmb.testc and stack them in one single time series
>>> cmb = MagicCoeffCmb(tag='test[a-c]')
>>> print(cmb.ell, cmb.glm) # print \ell and g_{\ell m}
>>> print(cmb.glm[:, cmb.idx[1, 0]]) # time-series of the axisymmetric dipole
>>> plot(cmb.time, cmb.dglmdt[:, cmb.idx[2, 0]]) # Secular variation of the quadrupole
>>> # Display the time-evolution of the CMB field
>>> cmb.movieCmb(levels=12, cm='seismic')
>>> # Save the time-evolution of the CMB field
>>> cmb.movieCmb(levels=12, cm='seismic', png=True)
__add__(new)[source]

Built-in function to sum two cmb files

Note

So far this function only works for two cmb files with the same grid sizes. At some point, we might introduce grid extrapolation to allow any summation/

__init__(tag=None, datadir='.', ratio_cmb_surface=1, scale_b=1, iplot=True, lCut=None, precision=<class 'numpy.float64'>, ave=False, sv=False, quiet=False)[source]

A class to read the B_coeff_cmb files

Parameters:
  • tag (str) – if you specify a pattern, it tries to read the corresponding files

  • ratio_cmb_surface (float) – ratio of CMB to surface radius (default is 1)

  • scale_b (float) – magnetic field unit (default is 1)

  • iplot (int) – a logical to toggle the plot (default is True)

  • precision (char) – single or double precision

  • ave (bool) – load a time-averaged CMB file when set to True

  • sv (bool) – load a dt_b CMB file when set to True

  • quiet (bool) – verbose when toggled to True (default is True)

  • lCut (int) – reduce the spherical harmonic truncation to l <= lCut

  • datadir (str) – working directory

movieCmb(cut=0.5, levels=12, cm='RdYlBu_r', png=False, step=1, normed=False, dpi=80, bgcolor=None, deminc=True, removeMean=False, precision=<class 'numpy.float64'>, contour=False, mer=False)[source]

Plotting function (it can also write the png files)

Parameters:
  • levels (int) – number of contour levels

  • cm (str) – name of the colormap

  • cut (float) – adjust the contour extrema to max(abs(data))*cut

  • png (bool) – save the movie as a series of png files when set to True

  • dpi (int) – dot per inch when saving PNGs

  • bgcolor (str) – background color of the figure

  • normed (bool) – the colormap is rescaled every timestep when set to True, otherwise it is calculated from the global extrema

  • step (int) – the stepping between two timesteps

  • deminc (bool) – a logical to indicate if one wants do get rid of the possible azimuthal symmetry

  • precision (char) – single or double precision

  • contour (bool) – also display the solid contour levels when set to True

  • mer (bool) – display meridians and circles when set to True

  • removeMean (bool) – remove the time-averaged part when set to True

plot()[source]

Display some results when iplot is set to True

timeLongitude(removeMean=True, lat0=0.0, levels=12, cm='RdYlBu_r', deminc=True, shtns_lib='shtns')[source]

Plot the time-longitude diagram of Br (input latitude can be chosen)

Warning

the python bindings of SHTns are mandatory to use this plotting function!

Parameters:
  • lat0 (float) – value of the latitude

  • levels (int) – number of contour levels

  • cm (str) – name of the colormap

  • deminc (bool) – a logical to indicate if one wants do get rid of the possible azimuthal symmetry

  • shtns_lib (char) – version of shtns library used: can be either ‘shtns’ or ‘shtns-magic’

  • removeMean (bool) – remove the time-averaged part when set to True

truncate(lCut)[source]
Parameters:

lCut (int) – truncate to spherical harmonic degree lCut

class magic.coeff.MagicCoeffR(tag, datadir='.', ratio_cmb_surface=1, scale_b=1, iplot=True, field='B', r=1, precision=<class 'numpy.float64'>, lCut=None, quiet=False, step=1)[source]

This class allows to read the B_coeff_r#.TAG and V_coeff_r#.TAG files. It reads the poloidal and toroidal potentials and reconstruct the time series (or the energy) contained in any given mode.

>>> # Reads the files V_coeff_r2.test*
>>> cr = MagicCoeffR(tag='test*', field='V', r=2)
>>> print(cr.ell, cr.wlm) # print \ell and w_{\ell m}
>>> # Time-evolution of the poloidal energy in the (\ell=10, m=10) mode
>>> plot(cr.time, cr.epolLM[:, cr.idx[10, 10]])
__init__(tag, datadir='.', ratio_cmb_surface=1, scale_b=1, iplot=True, field='B', r=1, precision=<class 'numpy.float64'>, lCut=None, quiet=False, step=1)[source]
Parameters:
  • tag (str) – if you specify a pattern, it tries to read the corresponding files

  • ratio_cmb_surface (float) – ratio of surface ratio to CMB radius (default is 1)

  • scale_b (float) – magnetic field unit (default is 1)

  • iplot (bool) – a logical to toggle the plot (default is True)

  • field (str) – ‘B’, ‘V’, ‘T’ or ‘Xi’ (magnetic field, velocity field, temperature or composition)

  • r (int) – an integer to characterise which file we want to plot

  • precision (str) – single or double precision

  • lCut (int) – reduce the spherical harmonic truncation to l <= lCut

  • quiet (bool) – verbose when toggled to True (default is True)

  • datadir (str) – working directory

  • step (int) – step>1 allows to down sample the data

cwt(ell, w0=20, nfreq=256, fmin_fac=8)[source]

Build a time-frequency spectrum at a given degree \(\ell\) using a continuous wavelet transform with morlet wavelets.

Parameters:
  • w0 (float) – a parameter to normalize the width of the wavelet

  • fmin_fac (float) – a factor to adjust the minimum frequency considered in the time-frequency domain. Minimum frequency is given by fmin=1/(time[-1]-time[0]), such that the minimum frequency retained is fmin_fac*fmin

  • ell (int) – spherical harmonic degree at which ones want to build the time frequency diagram

  • nfreq (int) – number of frequency bins

fft()[source]

Fourier transform of the poloidal potential

movieRad(cut=0.5, levels=12, cm='RdYlBu_r', png=False, step=1, normed=False, dpi=80, bgcolor=None, deminc=True, removeMean=False, precision=<class 'numpy.float64'>, contour=False, mer=False)[source]

Plotting function (it can also write the png files)

Parameters:
  • levels (int) – number of contour levels

  • cm (str) – name of the colormap

  • cut (float) – adjust the contour extrema to max(abs(data))*cut

  • png (bool) – save the movie as a series of png files when set to True

  • dpi (int) – dot per inch when saving PNGs

  • bgcolor (str) – background color of the figure

  • normed (bool) – the colormap is rescaled every timestep when set to True, otherwise it is calculated from the global extrema

  • step (int) – the stepping between two timesteps

  • deminc (bool) – a logical to indicate if one wants do get rid of the possible azimuthal symmetry

  • precision (char) – single or double precision

  • contour (bool) – also display the solid contour levels when set to True

  • mer (bool) – display meridians and circles when set to True

  • removeMean (bool) – remove the time-averaged part when set to True

truncate(lCut, field='B')[source]
Parameters:
  • lCut (int) – truncate to spherical harmonic degree lCut

  • field (char) – name of the field (‘V’, ‘B’ or ‘T’)

magic.coeff.deriv(x, y, axis=0)[source]

This function is a simple second order derivative

Parameters:
  • x (numpy.ndarray) – input x-axis

  • y (numpy.ndarray) – input array

Returns:

an array that contains the derivatives

Return type:

numpy.ndarray

magic.coeff.getGauss(alm, blm, ell, m, scale_b, ratio_cmb_surface, rcmb)[source]

Get the Gauss coefficients from the real and imaginary parts of the poloidal potential

Parameters:
  • alm (numpy.ndarray) – real part of the poloidal potential

  • blm (numpy.ndarray) – imaginary part of the poloidal potential

  • ell (numpy.ndarray) – spherical harmonic degree ell

  • scale_b (float) – magnetic field unit (default is 1)

  • ratio_cmb_surface (float) – ratio of CMB to surface radius (default is 1)

  • rcmb (float) – radius of the outer boundary

magic.coeff.rearangeLat(field)[source]

This function is used to unfold the colatitudes

Parameters:

field (numpy.ndarray) – input array with MagIC ordering of colatitudes (i.e. successively Northern Hemisphere and Southern Hemisphere)

Returns:

an array with the regular ordering of the colatitudes

Return type:

numpy.ndarray

Support for B[rp]Spec.TAG

class magic.MagicRSpec(tag, field='Br', precision=<class 'numpy.float32'>, avg=False)[source]

This class allows to read the rB[r|p]Spec.TAG files. Those files contain the time-evolution of the poloidal/toroidal magnetic energy for all radii and for spherical harmonic degrees from 1 to 6. This is an unformatted fortran file.

>>> # Read all the `BrSpec.test*` files in the current working directory and
>>> # stack them.
>>> rsp = MagicRSpec(tag='test*', field='Br')
__init__(tag, field='Br', precision=<class 'numpy.float32'>, avg=False)[source]
Parameters:
  • tag (str) – if you specify a pattern, it tries to read the corresponding files and stack them.

  • field (str) – nature of the radial spectra. Possible choices are ‘Bt’ or ‘Bp’

  • precision (str) – single or double precision (default single, i.e. np.float32)

  • avg (bool) – when set to True, display time averaged quantities

plotAvg()[source]

Plotting function for time-averaged profiles

Support for [V|B|T]_lmr_#.TAG

class magic.MagicPotential(field='V', datadir='.', tag=None, ave=False, ipot=None, precision=<class 'numpy.float32'>, verbose=True, ic=False)[source]

This class allows to load and display the content of the potential files: V_lmr.TAG, B_lmr.TAG and T_lmr.TAG. This class allows to transform the poloidal/toroidal potential in spectral space to the physical quantities in the physical space. It allows to plot radial and equatorial cuts as well as phi-averages.

>>> # To read T_lmr.test
>>> p = MagicPotential(field='T', ipot=1, tag='test')
>>> # To read the latest V_lmr file in the working directory
>>> p = MagicPotential(field='V')
>>> # Get the poloidal potential (lm, nR)
>>> wlm = p.pol
>>> # Obtain the value of w(l=12, m=12, nR=33)
>>> print( p.pol[p.idx[12,12], 32] )
>>> # Possible plots
>>> p.equat(field='vr')
>>> p.avg(field='vp')
>>> p.surf(field='vt', r=0.8)
__init__(field='V', datadir='.', tag=None, ave=False, ipot=None, precision=<class 'numpy.float32'>, verbose=True, ic=False)[source]
Parameters:
  • field (str) – ‘B’, ‘V’, ‘T’ or ‘Xi’ (magnetic field, velocity field, temperature or chemical composition)

  • datadir (str) – the working directory

  • tag (str) – if you specify a pattern, it tries to read the corresponding files

  • ave (bool) – plot a time-averaged spectrum when set to True

  • ipot (int) – the number of the lmr file you want to plot

  • precision (str) – single or double precision

  • verbose (bool) – some info about the SHT layout

  • ic (bool) – read or don’t read the inner core

avg(field='vphi', levels=65, cm='seismic', normed=True, vmax=None, vmin=None, cbar=True, tit=True)[source]

Plot the azimutal average of a given field.

>>> p = MagicPotential(field='V')
>>> # Axisymmetric zonal flows, 65 contour levels
>>> p.avg(field='vp', levels=65, cm='seismic')
>>> # Minimal plot (no cbar, not title)
>>> p.avg(field='vr', tit=False, cbar=False)
Parameters:
  • field (str) – the field you want to display

  • levels (int) – the number of levels in the contourf plot

  • cm (str) – name of the colormap (‘jet’, ‘seismic’, ‘RdYlBu_r’, etc.)

  • tit (bool) – display the title of the figure when set to True

  • cbar (bool) – display the colorbar when set to True

  • vmax (float) – maximum value of the contour levels

  • vmin (float) – minimum value of the contour levels

  • normed (bool) – when set to True, the colormap is centered around zero. Default is True, except for entropy/temperature plots.

equat(field='vr', levels=65, cm='seismic', normed=True, vmax=None, vmin=None, cbar=True, tit=True, normRad=False)[source]

Plot the equatorial cut of a given field

>>> p = MagicPotential(field='B')
>>> # Equatorial cut of the Br
>>> p.equat(field='Br')
>>> # Normalise the contour levels radius by radius
>>> p.equat(field='Bphi', normRad=True)
Parameters:
  • field (str) – the name of the input physical quantity you want to display

  • normRad (bool) – when set to True, the contour levels are normalised radius by radius (default is False)

  • levels (int) – the number of levels in the contour

  • cm (str) – name of the colormap (‘jet’, ‘seismic’, ‘RdYlBu_r’, etc.)

  • tit (bool) – display the title of the figure when set to True

  • cbar (bool) – display the colorbar when set to True

  • vmax (float) – maximum value of the contour levels

  • vmin (float) – minimum value of the contour levels

  • normed (bool) – when set to True, the colormap is centered around zero. Default is True, except for entropy/temperature plots.

read(filename, field, endian, record_marker, ic=False, precision=<class 'numpy.float32'>)[source]

This routine defines a reader for the various versions of the lmr files.

Parameters:
  • filename (str) – name of the input lmr file

  • field (str) – ‘B’, ‘V’, ‘T’ or ‘Xi’ (magnetic field, velocity field, temperature or chemical composition)

  • endian (str) – a character string that specifies the endianness of the input file (‘B’ for big endian or ‘l’ for little endian)

  • record_marker (bool) – a boolean to specify whether the file contains record marker

  • ic (bool) – read or don’t read the inner core

  • precision (str) – single or double precision

surf(field='vr', proj='hammer', lon_0=0.0, r=0.85, vmax=None, vmin=None, lat_0=30.0, levels=65, cm='seismic', lon_shift=0, normed=True, cbar=True, tit=True, lines=False)[source]

Plot the surface distribution of an input field at a given input radius (normalised by the outer boundary radius).

>>> p = MagicPotential(field='V')
>>> # Radial flow component at ``r=0.95 r_o``, 65 contour levels
>>> p.surf(field='vr', r=0.95, levels=65, cm='seismic')
>>> # Control the limit of the colormap from -1e3 to 1e3
>>> p.surf(field='vp', r=1., vmin=-1e3, vmax=1e3, levels=33)
Parameters:
  • field (str) – the name of the field you want to display

  • proj (str) – the type of projection. Default is Hammer, in case you want to use ‘ortho’ or ‘moll’, then Basemap is required.

  • lon_0 (float) – central azimuth (only used with Basemap)

  • lat_0 (float) – central latitude (only used with Basemap)

  • r (float) – the radius at which you want to display the input data (in normalised units with the radius of the outer boundary)

  • levels (int) – the number of levels in the contour

  • cm (str) – name of the colormap (‘jet’, ‘seismic’, ‘RdYlBu_r’, etc.)

  • tit (bool) – display the title of the figure when set to True

  • cbar (bool) – display the colorbar when set to True

  • lines (bool) – when set to True, over-plot solid lines to highlight the limits between two adjacent contour levels

  • vmax (float) – maximum value of the contour levels

  • vmin (float) – minimum value of the contour levels

  • normed (bool) – when set to True, the colormap is centered around zero.