Additional tools

class pyPLUTO.Tools

This Class has all the functions doing basic mathematical operations to the vector or scalar fields. It is called after pyPLUTO.pload object is defined.

Div(u1, u2, x1, x2, dx1, dx2, geometry=None)

This method calculates the divergence of the 2D vector fields u1 and u2.

Inputs:

u1 – 2D vector along x1 whose divergence is to be determined.

u2 – 2D vector along x2 whose divergence is to be determined.

x1 – The ‘x’ array

x2 – The ‘y’ array

dx1 – The grid spacing in ‘x’ direction.

dx2 – The grid spacing in ‘y’ direction.

geometry – The keyword geometry is by default set to ‘cartesian’. It can be set to either one of the following : cartesian, cylindrical, spherical or polar. To calculate the divergence of the vector fields, respective geometric corrections are taken into account based on the value of this keyword.

Outputs:

A 2D array with same shape as u1(or u2) having the values of divergence.
Grad(phi, x1, x2, dx1, dx2, polar=False)

This method calculates the gradient of the 2D scalar phi.

Inputs:

phi – 2D scalar whose gradient is to be determined.

x1 – The ‘x’ array

x2 – The ‘y’ array

dx1 – The grid spacing in ‘x’ direction.

dx2 – The grid spacing in ‘y’ direction.

polar – The keyword should be set to True inorder to estimate the Gradient in polar co-ordinates. By default it is set to False.

Outputs:

This routine outputs a 3D array with shape = (len(x1),len(x2),2), such that [:,:,0] element corresponds to the gradient values of phi wrt to x1 and [:,:,1] are the gradient values of phi wrt to x2.
RTh2Cyl(R, Th, X1, X2)

This method does the transformation from spherical coordinates to cylindrical ones.

Inputs:

R - 2D array of spherical radius coordinates.

Th - 2D array of spherical theta-angle coordinates.

X1 - 2D array of radial component of given vector

X2 - 2D array of thetoidal component of given vector

Outputs:

This routine outputs two 2D arrays after transformation.

Usage:

import pyPLUTO as pp

import numpy as np

D = pp.pload(0)

ppt=pp.Tools()

TH,R=np.meshgrid(D.x2,D.x1)

Br,Bz=ppt.RTh2Cyl(R,TH,D.bx1,D.bx2)

D.bx1 and D.bx2 should be vectors in spherical coordinates. After transformation (Br,Bz) corresponds to vector in cilindrical coordinates.

congrid(a, newdims, method='linear', centre=False, minusone=False)

Arbitrary resampling of source array to new dimension sizes. Currently only supports maintaining the same number of dimensions. To use 1-D arrays, first promote them to shape (x,1).

Uses the same parameters and creates the same co-ordinate lookup points as IDL’’s congrid routine, which apparently originally came from a VAX/VMS routine of the same name.

Inputs:

a – The 2D array that needs resampling into new dimensions.

newdims – A tuple which represents the shape of the resampled data

method – This keyword decides the method used for interpolation.

neighbour - closest value from original data

nearest and linear - uses n x 1-D interpolations using scipy.interpolate.interp1d (see Numerical Recipes for validity of use of n 1-D interpolations)

spline - uses ndimage.map_coordinates

centre – This keyword decides the positions of interpolation points.

True - interpolation points are at the centres of the bins

False - points are at the front edge of the bin

minusone – This prevents extrapolation one element beyond bounds of input array

For example- inarray.shape = (i,j) & new dimensions = (x,y)

False - inarray is resampled by factors of (i/x) * (j/y)

True - inarray is resampled by(i-1)/(x-1) * (j-1)/(y-1)

Outputs:

A 2D array with resampled data having a shape corresponding to newdims.
deriv(Y, X=None)

Calculates the derivative of Y with respect to X.

Inputs:

Y : 1-D array to be differentiated.

X : 1-D array with len(X) = len(Y).

If X is not specified then by default X is chosen to be an equally spaced array having same number of elements as Y.

Outputs:

This returns an 1-D array having the same no. of elements as Y (or X) and contains the values of dY/dX.
getUniformGrid(r, th, rho, Nr, Nth)

This method transforms data with non-uniform grid (stretched) to uniform. Useful for stretched grid calculations.

Inputs:

r - 1D vector of X1 coordinate (could be any, e.g D.x1).

th - 1D vector of X2 coordinate (could be any, e.g D.x2).

rho- 2D array of data.

Nr - new size of X1 vector.

Nth- new size of X2 vector.

Outputs:

This routine outputs 2D uniform array Nr x Nth dimension

Usage:

import pyPLUTO as pp

import numpy as np

D = pp.pload(0)

ppt=pp.Tools()

X1new, X2new, res = ppt.getUniformGrid(D.x1,D.x2,D.rho,20,30)

X1new - X1 interpolated grid len(X1new)=20 X2new - X2 interpolated grid len(X2new)=30 res - 2D array of interpolated variable

myInterpol(RR, N)

This method interpolates (linear interpolation) vector 1D vector RR to 1D N-length vector. Useful for stretched grid calculations.

Inputs:

RR - 1D array to interpolate.

N - Number of grids to interpolate to.

Outputs:

This routine outputs interpolated 1D array to the new grid (len=N).

Usage:

import pyPLUTO as pp

import numpy as np

D = pp.pload(0)

ppt=pp.Tools()

x=linspace(0,1,10) #len(x)=10

y=x*x

Ri,Ni=ppt.myInterpol(y,100) #len(Ri)=100

Ri - interpolated numbers; Ni - grid for Ri

sph2cyl(D, Dx, rphi=None, theta0=None)

This method transforms spherical data into cylindrical applying interpolation. Works for stretched grid as well, transforms poloidal (R-Theta) data by default. Fix theta and set rphi=True to get (R-Phi) transformation.

Inputs:

D - structure from ‘pload’ method.

Dx - variable to be transformed (D.rho for example).

Outputs:

This routine outputs transformed (sph->cyl) variable and grid.

Usage:

import pyPLUTO as pp

import numpy as np

D = pp.pload(0)

ppt=pp.Tools()

R,Z,res = ppt.sph2cyl(D,D.rho.transpose())

R - 2D array with cylindrical radius values Z - 2D array with cylindrical Z values res - 2D array of transformed variable

Previous topic

Viewing the Data

This Page