PLUTO
pyPLUTO.pyPLUTO.Tools Class Reference
Inheritance diagram for pyPLUTO.pyPLUTO.Tools:
Collaboration diagram for pyPLUTO.pyPLUTO.Tools:

Public Member Functions

def deriv
 
def Grad
 
def Div
 
def RTh2Cyl (self, R, Th, X1, X2)
 
def myInterpol (self, RR, N)
 
def getUniformGrid (self, r, th, rho, Nr, Nth)
 
def congrid
 

Detailed Description

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.

Definition at line 742 of file pyPLUTO.py.

Member Function Documentation

def pyPLUTO.pyPLUTO.Tools.congrid (   self,
  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.\n
  newdims -- A tuple which represents the shape of the resampled data\n
  method -- This keyword decides the method used for interpolation.\n
    neighbour - closest value from original data\n
    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)\n
    spline - uses ndimage.map_coordinates\n
  centre -- This keyword decides the positions of interpolation points.\n
    True - interpolation points are at the centres of the bins\n
    False - points are at the front edge of the bin\n
  minusone -- This prevents extrapolation one element beyond bounds of input array\n
    For example- inarray.shape = (i,j) & new dimensions = (x,y)\n
    False - inarray is resampled by factors of (i/x) * (j/y)\n
    True - inarray is resampled by(i-1)/(x-1) * (j-1)/(y-1)\n

    **Outputs:**

      A 2D array with resampled data having a shape corresponding to newdims. 

Definition at line 1129 of file pyPLUTO.py.

1129  def congrid(self, a, newdims, method='linear', centre=False, minusone=False):
1130  """
1131  Arbitrary resampling of source array to new dimension sizes.
1132  Currently only supports maintaining the same number of dimensions.
1133  To use 1-D arrays, first promote them to shape (x,1).
1134 
1135  Uses the same parameters and creates the same co-ordinate lookup points
1136  as IDL''s congrid routine, which apparently originally came from a VAX/VMS
1137  routine of the same name.
1138 
1139  **Inputs:**
1140 
1141  a -- The 2D array that needs resampling into new dimensions.\n
1142  newdims -- A tuple which represents the shape of the resampled data\n
1143  method -- This keyword decides the method used for interpolation.\n
1144  neighbour - closest value from original data\n
1145  nearest and linear - uses n x 1-D interpolations using scipy.interpolate.interp1d
1146  (see Numerical Recipes for validity of use of n 1-D interpolations)\n
1147  spline - uses ndimage.map_coordinates\n
1148  centre -- This keyword decides the positions of interpolation points.\n
1149  True - interpolation points are at the centres of the bins\n
1150  False - points are at the front edge of the bin\n
1151  minusone -- This prevents extrapolation one element beyond bounds of input array\n
1152  For example- inarray.shape = (i,j) & new dimensions = (x,y)\n
1153  False - inarray is resampled by factors of (i/x) * (j/y)\n
1154  True - inarray is resampled by(i-1)/(x-1) * (j-1)/(y-1)\n
1155 
1156  **Outputs:**
1157 
1158  A 2D array with resampled data having a shape corresponding to newdims.
1159 
1160  """
1161  if not a.dtype in [np.float64, np.float32]:
1162  a = np.cast[float](a)
1163 
1164  m1 = np.cast[int](minusone)
1165  ofs = np.cast[int](centre) * 0.5
1166  old = np.array( a.shape )
1167  ndims = len( a.shape )
1168  if len( newdims ) != ndims:
1169  print "[congrid] dimensions error. " \
1170  "This routine currently only support " \
1171  "rebinning to the same number of dimensions."
1172  return None
1173  newdims = np.asarray( newdims, dtype=float )
1174  dimlist = []
1175 
1176  if method == 'neighbour':
1177  for i in range( ndims ):
1178  base = np.indices(newdims)[i]
1179  dimlist.append( (old[i] - m1) / (newdims[i] - m1) \
1180  * (base + ofs) - ofs )
1181  cd = np.array( dimlist ).round().astype(int)
1182  newa = a[list( cd )]
1183  return newa
1184 
1185  elif method in ['nearest','linear']:
1186  # calculate new dims
1187  for i in range( ndims ):
1188  base = np.arange( newdims[i] )
1189  dimlist.append( (old[i] - m1) / (newdims[i] - m1) \
1190  * (base + ofs) - ofs )
1191  # specify old dims
1192  olddims = [np.arange(i, dtype = np.float) for i in list( a.shape )]
1193 
1194  # first interpolation - for ndims = any
1195  mint = scipy.interpolate.interp1d( olddims[-1], a, kind=method )
1196  newa = mint( dimlist[-1] )
1197 
1198  trorder = [ndims - 1] + range( ndims - 1 )
1199  for i in range( ndims - 2, -1, -1 ):
1200  newa = newa.transpose( trorder )
1201 
1202  mint = scipy.interpolate.interp1d( olddims[i], newa, kind=method )
1203  newa = mint( dimlist[i] )
1204 
1205  if ndims > 1:
1206  # need one more transpose to return to original dimensions
1207  newa = newa.transpose( trorder )
1208 
1209  return newa
1210  elif method in ['spline']:
1211  oslices = [ slice(0,j) for j in old ]
1212  oldcoords = np.ogrid[oslices]
1213  nslices = [ slice(0,j) for j in list(newdims) ]
1214  newcoords = np.mgrid[nslices]
1215 
1216  newcoords_dims = range(n.rank(newcoords))
1217  #make first index last
1218  newcoords_dims.append(newcoords_dims.pop(0))
1219  newcoords_tr = newcoords.transpose(newcoords_dims)
1220  # makes a view that affects newcoords
1221 
1222  newcoords_tr += ofs
1223 
1224  deltas = (np.asarray(old) - m1) / (newdims - m1)
1225  newcoords_tr *= deltas
1226 
1227  newcoords_tr -= ofs
1228 
1229  newa = scipy.ndimage.map_coordinates(a, newcoords)
1230  return newa
1231  else:
1232  print "Congrid error: Unrecognized interpolation type.\n", \
1233  "Currently only \'neighbour\', \'nearest\',\'linear\',", \
1234  "and \'spline\' are supported."
1235  return None
1236 
1237 
1238 
1239 
1240 
1241 
1242 
1243 
1244 
1245 
1246 
def pyPLUTO.pyPLUTO.Tools.deriv (   self,
  Y,
  X = None 
)
Calculates the derivative of Y with respect to X.

**Inputs:**

  Y : 1-D array to be differentiated.\n
  X : 1-D array with len(X) = len(Y).\n

  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.

Definition at line 751 of file pyPLUTO.py.

751  def deriv(self,Y,X=None):
752  """
753  Calculates the derivative of Y with respect to X.
754 
755  **Inputs:**
756 
757  Y : 1-D array to be differentiated.\n
758  X : 1-D array with len(X) = len(Y).\n
759 
760  If X is not specified then by default X is chosen to be an equally spaced array having same number of elements
761  as Y.
762 
763  **Outputs:**
764 
765  This returns an 1-D array having the same no. of elements as Y (or X) and contains the values of dY/dX.
766 
767  """
768  n = len(Y)
769  n2 = n-2
770  if X==None : X = np.arange(n)
771  Xarr = np.asarray(X,dtype='float')
772  Yarr = np.asarray(Y,dtype='float')
773  x12 = Xarr - np.roll(Xarr,-1) #x1 - x2
774  x01 = np.roll(Xarr,1) - Xarr #x0 - x1
775  x02 = np.roll(Xarr,1) - np.roll(Xarr,-1) #x0 - x2
776  DfDx = np.roll(Yarr,1) * (x12 / (x01*x02)) + Yarr * (1./x12 - 1./x01) - np.roll(Yarr,-1) * (x01 / (x02 * x12))
777  # Formulae for the first and last points:
778 
779  DfDx[0] = Yarr[0] * (x01[1]+x02[1])/(x01[1]*x02[1]) - Yarr[1] * x02[1]/(x01[1]*x12[1]) + Yarr[2] * x01[1]/(x02[1]*x12[1])
780  DfDx[n-1] = -Yarr[n-3] * x12[n2]/(x01[n2]*x02[n2]) + Yarr[n-2]*x02[n2]/(x01[n2]*x12[n2]) - Yarr[n-1]*(x02[n2]+x12[n2])/(x02[n2]*x12[n2])
781 
782  return DfDx
783 

Here is the caller graph for this function:

def pyPLUTO.pyPLUTO.Tools.Div (   self,
  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.\n
  u2 -- 2D vector along x2 whose divergence is to be determined.\n
  x1 -- The 'x' array\n
  x2 -- The 'y' array\n
  dx1 -- The grid spacing in 'x' direction.\n
  dx2 -- The grid spacing in 'y' direction.\n
  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.

Definition at line 817 of file pyPLUTO.py.

817  def Div(self,u1,u2,x1,x2,dx1,dx2,geometry=None):
818  """ This method calculates the divergence of the 2D vector fields u1 and u2.
819 
820  **Inputs:**
821 
822  u1 -- 2D vector along x1 whose divergence is to be determined.\n
823  u2 -- 2D vector along x2 whose divergence is to be determined.\n
824  x1 -- The 'x' array\n
825  x2 -- The 'y' array\n
826  dx1 -- The grid spacing in 'x' direction.\n
827  dx2 -- The grid spacing in 'y' direction.\n
828  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.
829 
830  **Outputs:**
831 
832  A 2D array with same shape as u1(or u2) having the values of divergence.
833 
834  """
835  (n1, n2) = u1.shape
836  Divergence = np.zeros(shape=(n1,n2))
837  du1 = np.zeros(shape=(n1,n2))
838  du2 = np.zeros(shape=(n1,n2))
839 
840  A1 = np.zeros(shape=n1)
841  A2 = np.zeros(shape=n2)
842 
843  dV1 = np.zeros(shape=(n1,n2))
844  dV2 = np.zeros(shape=(n1,n2))
845 
846  if geometry == None : geometry = 'cartesian'
847 
848  #------------------------------------------------
849  # define area and volume elements for the
850  # different coordinate systems
851  #------------------------------------------------
852 
853  if geometry == 'cartesian' :
854  A1[:] = 1.0
855  A2[:] = 1.0
856  dV1 = np.outer(dx1,A2)
857  dV2 = np.outer(A1,dx2)
858 
859  if geometry == 'cylindrical' :
860  A1 = x1
861  A2[:] = 1.0
862  dV1 = np.meshgrid(x1*dx1,A2)[0].T*np.meshgrid(x1*dx1,A2)[1].T
863  for i in range(n1) : dV2[i,:] = dx2[:]
864 
865  if geometry == 'polar' :
866  A1 = x1
867  A2[:] = 1.0
868  dV1 = np.meshgrid(x1,A2)[0].T*np.meshgrid(x1,A2)[1].T
869  dV2 = np.meshgrid(x1,dx2)[0].T*np.meshgrid(x1,dx2)[1].T
870 
871  if geometry == 'spherical' :
872  A1 = x1*x1
873  A2 = np.sin(x2)
874  for j in range(n2): dV1[:,j] = A1*dx1
875  dV2 = np.meshgrid(x1,np.sin(x2)*dx2)[0].T*np.meshgrid(x1,np.sin(x2)*dx2)[1].T
876 
877  # ------------------------------------------------
878  # Make divergence
879  # ------------------------------------------------
880 
881 
882  for i in range(1,n1-1):
883  du1[i,:] = 0.5*(A1[i+1]*u1[i+1,:] - A1[i-1]*u1[i-1,:])/dV1[i,:]
884  for j in range(1,n2-1):
885  du2[:,j] = 0.5*(A2[j+1]*u2[:,j+1] - A2[j-1]*u2[:,j-1])/dV2[:,j]
886 
887  Divergence = du1 + du2
888  return Divergence
889 
890 
891 
def pyPLUTO.pyPLUTO.Tools.getUniformGrid (   self,
  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).\n
  th - 1D vector of X2 coordinate (could be any, e.g D.x2).\n
  rho- 2D array of data.\n
  Nr - new size of X1 vector.\n
  Nth- new size of X2 vector.\n
  
**Outputs:**

  This routine outputs 2D uniform array Nr x Nth dimension
  
**Usage:**

  ``import pyPLUTO as pp``\n
  ``import numpy as np``\n
  ``D = pp.pload(0)``\n
  ``ppt=pp.Tools()``\n
  ``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

Definition at line 962 of file pyPLUTO.py.

962  def getUniformGrid(self,r,th,rho,Nr,Nth):
963  """ This method transforms data with non-uniform grid (stretched) to uniform. Useful for stretched grid calculations.
964 
965  **Inputs:**
966 
967  r - 1D vector of X1 coordinate (could be any, e.g D.x1).\n
968  th - 1D vector of X2 coordinate (could be any, e.g D.x2).\n
969  rho- 2D array of data.\n
970  Nr - new size of X1 vector.\n
971  Nth- new size of X2 vector.\n
972 
973  **Outputs:**
974 
975  This routine outputs 2D uniform array Nr x Nth dimension
976 
977  **Usage:**
978 
979  ``import pyPLUTO as pp``\n
980  ``import numpy as np``\n
981  ``D = pp.pload(0)``\n
982  ``ppt=pp.Tools()``\n
983  ``X1new, X2new, res = ppt.getUniformGrid(D.x1,D.x2,D.rho,20,30)``
984 
985  X1new - X1 interpolated grid len(X1new)=20
986  X2new - X2 interpolated grid len(X2new)=30
987  res - 2D array of interpolated variable
988 
989  """
990 
991  Ri,NRi=self.myInterpol(r,Nr)
992  Ra=np.int32(NRi);Wr=NRi-Ra
993 
994  YY=np.ones([Nr,len(th)])
995  for i in range(len(th)):
996  YY[:,i]=(1-Wr)*rho[Ra,i] + Wr*rho[Ra+1,i]
997 
998  THi,NTHi=self.myInterpol(th,Nth)
999  THa=np.int32(NTHi);Wth=NTHi-THa
1000 
1001  ZZ=np.ones([Nr,Nth])
1002  for i in range(Nr):
1003  ZZ[i,:]=(1-Wth)*YY[i,THa] + Wth*YY[i,THa+1]
1004 
1005  return Ri,THi,ZZ
1006 
1007  def sph2cyl(self,D,Dx,rphi=None,theta0=None):
1008  """ 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.
1009 
1010  **Inputs:**
1011 
1012  D - structure from 'pload' method.\n
1013  Dx - variable to be transformed (D.rho for example).\n
1014 
1015  **Outputs:**
1016 
1017  This routine outputs transformed (sph->cyl) variable and grid.
1018 
1019  **Usage:**
1020 
1021  ``import pyPLUTO as pp``\n
1022  ``import numpy as np``\n
1023  ``D = pp.pload(0)``\n
1024  ``ppt=pp.Tools()``\n
1025  ``R,Z,res = ppt.sph2cyl(D,D.rho.transpose())``
1026 
1027  R - 2D array with cylindrical radius values
1028  Z - 2D array with cylindrical Z values
1029  res - 2D array of transformed variable
1030 
1031  """
1032 
1033  if rphi is None or rphi == False:
1034  rx=D.x1
1035  th=D.x2
1036  else:
1037  rx=D.x1*np.sin(theta0)
1038  th=D.x3
1039 
1040  rx,th,Dx=self.getUniformGrid(rx,th,Dx.T,200,200)
1041  Dx=Dx.T
1042 
1043  if rphi is None or rphi == False:
1044 
1045  r0=np.min(np.sin(th)*rx[0])
1046  rN=rx[-1]
1047  dr=rN-r0
1048  z0=np.min(np.cos(th)*rN)
1049  zN=np.max(np.cos(th)*rN)
1050  dz=zN-z0
1051  dth=th[-1]-th[0]
1052  rl=np.int32(len(rx)*dr/(rx[-1]-rx[0]))
1053  zl=np.int32(rl* dz/dr)
1054  thl=len(th)
1055  r=np.linspace(r0, rN, rl)
1056  z=np.linspace(z0, zN, zl)
1057  else:
1058  r0=np.min([np.sin(th)*rx[0] , np.sin(th)*rx[-1]])
1059  rN=np.max([np.sin(th)*rx[0] , np.sin(th)*rx[-1]])
1060  dr=rN-r0
1061  z0=np.min(np.cos(th)*rN)
1062  zN=np.max(np.cos(th)*rN)
1063  dz=zN-z0
1064  dth=th[-1]-th[0]
1065  rl=np.int32(len(rx)*dr/(rx[-1]-rx[0]))
1066  zl=np.int32(rl* dz/dr)
1067  thl=len(th)
1068  r=np.linspace(r0, rN, rl)
1069  z=np.linspace(z0, zN, zl)
1070 
1071  R,Z = np.meshgrid(r, z)
1072  Rs = np.sqrt(R*R + Z*Z)
1073 
1074 
1075  Th = np.arccos(Z/Rs)
1076  kv_34=find(R<0)
1077  Th.flat[kv_34]=2*np.pi - Th.flat[kv_34]
1078 
1079 
1080  ddr=rx[1]-rx[0]
1081  ddth=th[1]-th[0]
1082 
1083  Rs_copy=Rs.copy()
1084  Th_copy=Th.copy()
1085 
1086  nR1=find(Rs<rx[0])
1087  Rs.flat[nR1]=rx[0]
1088  nR2=find(Rs>rN)
1089  Rs.flat[nR2]=rN
1090 
1091  nTh1=find(Th>th[-1])
1092  Th.flat[nTh1]=th[-1]
1093  nTh2=find(Th<th[0])
1094  Th.flat[nTh2]=th[0]
1095 
1096 
1097  ra = ((len(rx)-1.001)/(np.max(Rs.flat)-np.min(Rs.flat)) *(Rs-np.min(Rs.flat)))
1098  tha = ((thl-1.001)/dth *(Th-th[0]))
1099 
1100  rn = np.int32(ra)
1101  thn = np.int32(tha)
1102  dra=ra-rn
1103  dtha=tha-thn
1104  w1=1-dra
1105  w2=dra
1106  w3=1-dtha
1107  w4=dtha
1108  lrx=len(rx)
1109  NN1=np.int32(rn+thn*lrx)
1110  NN2=np.int32((rn+1)+thn*lrx)
1111  NN3=np.int32(rn+(thn+1)*lrx)
1112  NN4=np.int32((rn+1)+(thn+1)*lrx)
1113  n=np.transpose(np.arange(0,np.product(np.shape(R))))
1114  DD=Dx.copy()
1115  F=R.copy()
1116  F.flat[n]=w1.flat[n]*(w3.flat[n]*Dx.flat[NN1.flat[n]] + w4.flat[n]*Dx.flat[NN3.flat[n]] )+\
1117  w2.flat[n]*(w3.flat[n]*Dx.flat[NN2.flat[n]] + w4.flat[n]*Dx.flat[NN4.flat[n]] )
1118 
1119  nR1=find(Rs_copy<rx[0]-ddr/1.5)
1120  nR2=find(Rs_copy>rN+ddr/1.5)
1121  nTh1=find(Th_copy>th[-1]+ddth/1.5)
1122  nTh2=find(Th_copy<th[0]-ddth/1.5)
1123 
1124  nmask=np.concatenate((nR1,nR2,nTh1,nTh2))
1125  F.flat[nmask]=np.nan;
1126  return R,Z,F
1127 
1128 
def getUniformGrid(self, r, th, rho, Nr, Nth)
Definition: pyPLUTO.py:962
def find(fname, str, action, want)
Definition: file.py:150
def myInterpol(self, RR, N)
Definition: pyPLUTO.py:926

Here is the call graph for this function:

Here is the caller graph for this function:

def pyPLUTO.pyPLUTO.Tools.Grad (   self,
  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.\n
  x1 -- The 'x' array\n
  x2 -- The 'y' array\n
  dx1 -- The grid spacing in 'x' direction.\n
  dx2 -- The grid spacing in 'y' direction.\n
  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.

Definition at line 784 of file pyPLUTO.py.

784  def Grad(self,phi,x1,x2,dx1,dx2,polar=False):
785  """ This method calculates the gradient of the 2D scalar phi.
786 
787  **Inputs:**
788 
789  phi -- 2D scalar whose gradient is to be determined.\n
790  x1 -- The 'x' array\n
791  x2 -- The 'y' array\n
792  dx1 -- The grid spacing in 'x' direction.\n
793  dx2 -- The grid spacing in 'y' direction.\n
794  polar -- The keyword should be set to True inorder to estimate the Gradient in polar co-ordinates. By default it is set to False.
795 
796  **Outputs:**
797 
798  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.
799 
800  """
801  (n1, n2) = phi.shape
802  grad_phi = np.zeros(shape=(n1,n2,2))
803  h2 = np.ones(shape=(n1,n2))
804  if polar == True:
805  for j in range(n2):
806  h2[:,j] = x1
807 
808  for i in range(n1):
809  scrh1 = phi[i,:]
810  grad_phi[i,:,1] = self.deriv(scrh1,x2)/h2[i,:]
811  for j in range(n2):
812  scrh2 = phi[:,j]
813  grad_phi[:,j,0] = self.deriv(scrh2,x1)
814 
815  return grad_phi
816 

Here is the call graph for this function:

def pyPLUTO.pyPLUTO.Tools.myInterpol (   self,
  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
  N  - Number of grids to interpolate to.\n
  
**Outputs:**

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

  ``import pyPLUTO as pp``\n
  ``import numpy as np``\n
  ``D = pp.pload(0)``\n
  ``ppt=pp.Tools()``\n
  ``x=linspace(0,1,10) #len(x)=10``\n
  ``y=x*x``\n
  ``Ri,Ni=ppt.myInterpol(y,100) #len(Ri)=100``

  Ri - interpolated numbers;
  Ni - grid for Ri

Definition at line 926 of file pyPLUTO.py.

926  def myInterpol(self,RR,N):
927  """ This method interpolates (linear interpolation) vector 1D vector RR to 1D N-length vector. Useful for stretched grid calculations.
928 
929  **Inputs:**
930 
931  RR - 1D array to interpolate.\n
932  N - Number of grids to interpolate to.\n
933 
934  **Outputs:**
935 
936  This routine outputs interpolated 1D array to the new grid (len=N).
937 
938  **Usage:**
939 
940  ``import pyPLUTO as pp``\n
941  ``import numpy as np``\n
942  ``D = pp.pload(0)``\n
943  ``ppt=pp.Tools()``\n
944  ``x=linspace(0,1,10) #len(x)=10``\n
945  ``y=x*x``\n
946  ``Ri,Ni=ppt.myInterpol(y,100) #len(Ri)=100``
947 
948  Ri - interpolated numbers;
949  Ni - grid for Ri
950 
951  """
952 
953  NN=np.linspace(0,len(RR)-1,len(RR))
954  spline_fit=UnivariateSpline(RR,NN,k=3,s=0)
955 
956  RRi=np.linspace(RR[0],RR[-1],N)
957  NNi=spline_fit(RRi)
958  NNi[0]=NN[0]+0.00001
959  NNi[-1]=NN[-1]-0.00001
960  return RRi,NNi
961 
def myInterpol(self, RR, N)
Definition: pyPLUTO.py:926

Here is the caller graph for this function:

def pyPLUTO.pyPLUTO.Tools.RTh2Cyl (   self,
  R,
  Th,
  X1,
  X2 
)
This method does the transformation from spherical coordinates to cylindrical ones.

**Inputs:**

  R - 2D array of spherical radius coordinates.\n
      Th - 2D array of spherical theta-angle coordinates.\n
  X1 - 2D array of radial component of given vector\n
  X2 - 2D array of thetoidal component of given vector\n

**Outputs:**

  This routine outputs two 2D arrays after transformation.
  
**Usage:**

  ``import pyPLUTO as pp``\n
  ``import numpy as np``\n
  ``D = pp.pload(0)``\n
  ``ppt=pp.Tools()``\n
  ``TH,R=np.meshgrid(D.x2,D.x1)``\n
  ``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.

Definition at line 892 of file pyPLUTO.py.

892  def RTh2Cyl(self,R,Th,X1,X2):
893  """ This method does the transformation from spherical coordinates to cylindrical ones.
894 
895  **Inputs:**
896 
897  R - 2D array of spherical radius coordinates.\n
898  Th - 2D array of spherical theta-angle coordinates.\n
899  X1 - 2D array of radial component of given vector\n
900  X2 - 2D array of thetoidal component of given vector\n
901 
902  **Outputs:**
903 
904  This routine outputs two 2D arrays after transformation.
905 
906  **Usage:**
907 
908  ``import pyPLUTO as pp``\n
909  ``import numpy as np``\n
910  ``D = pp.pload(0)``\n
911  ``ppt=pp.Tools()``\n
912  ``TH,R=np.meshgrid(D.x2,D.x1)``\n
913  ``Br,Bz=ppt.RTh2Cyl(R,TH,D.bx1,D.bx2)``
914 
915  D.bx1 and D.bx2 should be vectors in spherical coordinates. After transformation (Br,Bz) corresponds to vector in cilindrical coordinates.
916 
917 
918  """
919  Y1=X1*np.sin(Th)+X2*np.cos(Th)
920  Y2=X1*np.cos(Th)-X2*np.sin(Th)
921  return Y1,Y2
922 
923 
924 
925 
def RTh2Cyl(self, R, Th, X1, X2)
Definition: pyPLUTO.py:892

The documentation for this class was generated from the following file: