These functions are not methods of mglData
class. However it provide additional functionality to handle data. So I put it in this chapter.
DAT 'type' real imag
¶mglData
mglTransform (const mglDataA &
real, const mglDataA &
imag, const char *
type)
¶HMDT
mgl_transform (HCDT
real, HCDT
imag, const char *
type)
¶Does integral transformation of complex data real, imag on specified direction. The order of transformations is specified in string type: first character for x-dimension, second one for y-dimension, third one for z-dimension. The possible character are: ‘f’ is forward Fourier transformation, ‘i’ is inverse Fourier transformation, ‘s’ is Sine transform, ‘c’ is Cosine transform, ‘h’ is Hankel transform, ‘n’ or ‘ ’ is no transformation.
DAT 'type' ampl phase
¶mglData
mglTransformA const mglDataA &
ampl, const mglDataA &
phase, const char *
type)
¶HMDT
mgl_transform_a HCDT
ampl, HCDT
phase, const char *
type)
¶The same as previous but with specified amplitude ampl and phase phase of complex numbers.
reDat imDat 'dir'
¶complexDat 'dir'
¶void
mglFourier const mglDataA &
re, const mglDataA &
im, const char *
dir)
¶mglDataC
: void
FFT (const char *
dir)
¶void
mgl_data_fourier HCDT
re, HCDT
im, const char *
dir)
¶void
mgl_datac_fft (HADT
dat, const char *
dir)
¶Does Fourier transform of complex data re+i*im in directions dir. Result is placed back into re and im data arrays. If dir contain ‘i’ then inverse Fourier is used.
RES real imag dn
['dir'='x']
¶mglData
mglSTFA (const mglDataA &
real, const mglDataA &
imag, int
dn, char
dir='x'
)
¶HMDT
mgl_data_stfa (HCDT
real, HCDT
imag, int
dn, char
dir)
¶Short time Fourier transformation for real and imaginary parts. Output is amplitude of partial Fourier of length dn. For example if dir=‘x’, result will have size {int(nx/dn), dn, ny} and it will contain res[i,j,k]=|\sum_d^dn exp(I*j*d)*(real[i*dn+d,k]+I*imag[i*dn+d,k])|/dn.
dat xdat ydat
¶mglData
mglTriangulation (const mglDataA &
x, const mglDataA &
y)
¶void
mgl_triangulation_2d (HCDT
x, HCDT
y)
¶Do Delone triangulation for 2d points and return result suitable for triplot and tricont. See Making regular data, for sample code and picture.
RES ADAT BDAT CDAT DDAT 'how'
¶mglData
mglTridMat (const mglDataA &
A, const mglDataA &
B, const mglDataA &
C, const mglDataA &
D, const char *
how)
¶mglDataC
mglTridMatC (const mglDataA &
A, const mglDataA &
B, const mglDataA &
C, const mglDataA &
D, const char *
how)
¶HMDT
mgl_data_tridmat (HCDT
A, HCDT
B, HCDT
C, HCDT
D, const char*
how)
¶HADT
mgl_datac_tridmat (HCDT
A, HCDT
B, HCDT
C, HCDT
D, const char*
how)
¶Get array as solution of tridiagonal system of equations A[i]*x[i-1]+B[i]*x[i]+C[i]*x[i+1]=D[i]. String how may contain:
Data dimensions of arrays A, B, C should be equal. Also their dimensions need to be equal to all or to minor dimension(s) of array D. See PDE solving hints, for sample code and picture.
RES 'ham' ini_re ini_im [dz=0.1 k0=100
]
¶mglData
mglPDE (HMGL
gr, const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, mreal
dz=0.1
, mreal
k0=100
, const char *
opt=""
)
¶mglDataC
mglPDEc (HMGL
gr, const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, mreal
dz=0.1
, mreal
k0=100
, const char *
opt=""
)
¶HMDT
mgl_pde_solve (HMGL
gr, const char *
ham, HCDT
ini_re, HCDT
ini_im, mreal
dz, mreal
k0, const char *
opt)
¶HADT
mgl_pde_solve_c (HMGL
gr, const char *
ham, HCDT
ini_re, HCDT
ini_im, mreal
dz, mreal
k0, const char *
opt)
¶Solves equation du/dz = i*k0*ham(p,q,x,y,z,|u|)[u], where p=-i/k0*d/dx, q=-i/k0*d/dy are pseudo-differential operators. Parameters ini_re, ini_im specify real and imaginary part of initial field distribution. Parameters Min, Max set the bounding box for the solution. Note, that really this ranges are increased by factor 3/2 for purpose of reducing reflection from boundaries. Parameter dz set the step along evolutionary coordinate z. At this moment, simplified form of function ham is supported – all “mixed” terms (like ‘x*p’->x*d/dx) are excluded. For example, in 2D case this function is effectively ham = f(p,z) + g(x,z,u). However commutable combinations (like ‘x*q’->x*d/dy) are allowed. Here variable ‘u’ is used for field amplitude |u|. This allow one solve nonlinear problems – for example, for nonlinear Shrodinger equation you may set ham="p^2 + q^2 - u^2"
. You may specify imaginary part for wave absorption, like ham = "p^2 + i*x*(x>0)"
. See also apde, qo2d, qo3d. See PDE solving hints, for sample code and picture.
RES 'ham' ini_re ini_im [dz=0.1 k0=100
]
¶mglData
mglAPDE (HMGL
gr, const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, mreal
dz=0.1
, mreal
k0=100
, const char *
opt=""
)
¶mglDataC
mglAPDEc (HMGL
gr, const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, mreal
dz=0.1
, mreal
k0=100
, const char *
opt=""
)
¶HMDT
mgl_pde_solve_adv (HMGL
gr, const char *
ham, HCDT
ini_re, HCDT
ini_im, mreal
dz, mreal
k0, const char *
opt)
¶HADT
mgl_pde_solve_adv_c (HMGL
gr, const char *
ham, HCDT
ini_re, HCDT
ini_im, mreal
dz, mreal
k0, const char *
opt)
¶Solves equation du/dz = i*k0*ham(p,q,x,y,z,|u|)[u], where p=-i/k0*d/dx, q=-i/k0*d/dy are pseudo-differential operators. Parameters ini_re, ini_im specify real and imaginary part of initial field distribution. Parameters Min, Max set the bounding box for the solution. Note, that really this ranges are increased by factor 3/2 for purpose of reducing reflection from boundaries. Parameter dz set the step along evolutionary coordinate z. The advanced and rather slow algorithm is used for taking into account both spatial dispersion and inhomogeneities of media [see A.A. Balakin, E.D. Gospodchikov, A.G. Shalashov, JETP letters v.104, p.690-695 (2016)]. Variable ‘u’ is used for field amplitude |u|. This allow one solve nonlinear problems – for example, for nonlinear Shrodinger equation you may set ham="p^2 + q^2 - u^2"
. You may specify imaginary part for wave absorption, like ham = "p^2 + i*x*(x>0)"
. See also pde. See PDE solving hints, for sample code and picture.
RES 'ham' x0 y0 z0 p0 q0 v0 [dt=0.1 tmax=10]
¶mglData
mglRay (const char *
ham, mglPoint
r0, mglPoint
p0, mreal
dt=0.1
, mreal
tmax=10
)
¶HMDT
mgl_ray_trace (const char *
ham, mreal
x0, mreal
y0, mreal
z0, mreal
px, mreal
py, mreal
pz, mreal
dt, mreal
tmax)
¶Solves GO ray equation like dr/dt = d ham/dp, dp/dt = -d ham/dr. This is Hamiltonian equations for particle trajectory in 3D case. Here ham is Hamiltonian which may depend on coordinates ‘x’, ‘y’, ‘z’, momentums ‘p’=px, ‘q’=py, ‘v’=pz and time ‘t’: ham = H(x,y,z,p,q,v,t). The starting point (at t=0
) is defined by variables r0, p0. Parameters dt and tmax specify the integration step and maximal time for ray tracing. Result is array of {x,y,z,p,q,v,t} with dimensions {7 * int(tmax/dt+1) }.
RES 'df' 'var' ini [dt=0.1 tmax=10
'jump'='']
¶mglData
mglODE (const char *
df, const char *
var, const mglDataA &
ini, mreal
dt=0.1
, mreal
tmax=10
, const char *
jump=""
)
¶mglDataC
mglODEc (const char *
df, const char *
var, const mglDataA &
ini, mreal
dt=0.1
, mreal
tmax=10
, const char *
jump=""
)
¶HMDT
mgl_ode_solve_str (const char *
df, const char *
var, HCDT
ini, mreal
dt, mreal
tmax)
¶HMDT
mgl_ode_solve_str_b (const char *
df, const char *
var, HCDT
ini, mreal
dt, mreal
tmax, const char *
jump)
¶HADT
mgl_ode_solve_str_c (const char *
df, const char *
var, HCDT
ini, mreal
dt, mreal
tmax)
¶HADT
mgl_ode_solve_str_cb (const char *
df, const char *
var, HCDT
ini, mreal
dt, mreal
tmax, const char *
jump)
¶HMDT
mgl_ode_solve (void (*
df)(const mreal *x, mreal *dx, void *par)
, int
n, const mreal *
ini, mreal
dt, mreal
tmax)
¶HMDT
mgl_ode_solve_ex (void (*
df)(const mreal *x, mreal *dx, void *par)
, int
n, const mreal *
ini, mreal
dt, mreal
tmax, void (*
jump)(mreal *x, const mreal *xprev, void *par)
)
¶Solves ODE equations dx/dt = df(x). The functions df can be specified as string of ’;’-separated textual formulas (argument var set the character ids of variables x[i]) or as callback function, which fill dx
array for give x
’s. Parameters ini, dt, tmax set initial values, time step and maximal time of the calculation. Function stop execution if NAN
or INF
values appears. Result is data array with dimensions {n * Nt}, where Nt <= int(tmax/dt+1).
If dt*tmax<0 then regularization is switched on, which change equations to dx/ds = df(x)/max(|df(x)|) to allow accurately passes region of strong df variation or quickly bypass region of small df. Here s is the new "time". At this, real time is determined as dt/ds=max(|df(x)|). If you need real time, then add it into equations manually, like ‘ode res 'y;-sin(x);1' 'xyt' [3,0] 0.3 -100’. This also preserve accuracy at stationary points (i.e. at small df in periodic case).
Functions jump are called each steps to handle special conditions (like reflection at boundary or jumps).
RES 'df' 'var' 'brd' ini [dt=0.1 tmax=10
]
¶mglData
mglODEs (const char *
df, const char *
var, char
brd, const mglDataA &
ini, mreal
dt=0.1
, mreal
tmax=10
)
¶mglDataC
mglODEcs (const char *
df, const char *
var, char
brd, const mglDataA &
ini, mreal
dt=0.1
, mreal
tmax=10
)
¶HMDT
mgl_ode_solve_set (const char *
df, const char *
var, char
brd, HCDT
ini, mreal
dt, mreal
tmax)
¶HADT
mgl_ode_solve_set_c (const char *
df, const char *
var, char
brd, HCDT
ini, mreal
dt, mreal
tmax)
¶Solves difference approximation of PDE as a set of ODE dx/dt = df(x,j). Functions df can be specified as string of ’;’-separated textual formulas, which can depend on index j and current time ‘t’. Argument var set the character ids of variables x[i]. Parameter brd sets the kind of boundary conditions on j: ‘0’ or ‘z’ – zero at border, ‘1’ or ‘c’ – constant at border, ‘2’ or ‘l’ – linear at border (laplacian is zero), ‘3’ or ‘s’ – square at border, ‘4’ or ‘e’ – exponential at border, ‘5’ or ‘g’ – gaussian at border. The cases ‘e’ and ‘g’ are applicable for the complex variant only. Parameters ini, dt, tmax set initial values, time step and maximal time of the calculation. Function stop execution if NAN
or INF
values appears. Result is data array with dimensions {n * Nt}, where Nt <= int(tmax/dt+1). For example, difference aprroximation of diffusion equation with zero boundary conditions can be solved by call: ‘ode res 'u(j+1)-2*u(j)+u(j-1)' 'u' '0' u0’, where ‘u0’ is an initial data array.
If dt*tmax<0 then regularization is switched on, which change equations to dx/ds = df(x)/max(|df(x)|) to allow accurately passes region of strong df variation or quickly bypass region of small df. Here s is the new "time". At this, real time is determined as dt/ds=max(|df(x)|). If you need real time, then add it into equations manually, like ‘ode res 'y;-sin(x);1' 'xyt' [3,0] 0.3 -100’. This also preserve accuracy at stationary points (i.e. at small df in periodic case).
RES 'ham' ini_re ini_im ray [r=1 k0=100
xx yy]
¶mglData
mglQO2d (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mreal
r=1
, mreal
k0=100
)
¶mglData
mglQO2d (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mglData &
xx, mglData &
yy, mreal
r=1
, mreal
k0=100
)
¶mglDataC
mglQO2dc (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mreal
r=1
, mreal
k0=100
)
¶mglDataC
mglQO2dc (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mglData &
xx, mglData &
yy, mreal
r=1
, mreal
k0=100
)
¶HMDT
mgl_qo2d_solve (const char *
ham, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy)
¶HADT
mgl_qo2d_solve_c (const char *
ham, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy)
¶HMDT
mgl_qo2d_func (dual (*
ham)(mreal u, mreal x, mreal y, mreal px, mreal py, void *par)
, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy)
¶HADT
mgl_qo2d_func_c (dual (*
ham)(mreal u, mreal x, mreal y, mreal px, mreal py, void *par)
, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy)
¶Solves equation du/dt = i*k0*ham(p,q,x,y,|u|)[u], where p=-i/k0*d/dx, q=-i/k0*d/dy are pseudo-differential operators (see mglPDE()
for details). Parameters ini_re, ini_im specify real and imaginary part of initial field distribution. Parameters ray set the reference ray, i.e. the ray around which the accompanied coordinate system will be maked. You may use, for example, the array created by ray function. Note, that the reference ray must be smooth enough to make accompanied coodrinates unambiguity. Otherwise errors in the solution may appear. If xx and yy are non-zero then Cartesian coordinates for each point will be written into them. See also pde, qo3d. See PDE solving hints, for sample code and picture.
RES 'ham' ini_re ini_im ray [r=1 k0=100
xx yy zz]
¶mglData
mglQO3d (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mreal
r=1
, mreal
k0=100
)
¶mglData
mglQO3d (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mglData &
xx, mglData &
yy, mglData &
zz, mreal
r=1
, mreal
k0=100
)
¶mglDataC
mglQO3dc (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mreal
r=1
, mreal
k0=100
)
¶mglDataC
mglQO3dc (const char *
ham, const mglDataA &
ini_re, const mglDataA &
ini_im, const mglDataA &
ray, mglData &
xx, mglData &
yy, mglData &
zz, mreal
r=1
, mreal
k0=100
)
¶HMDT
mgl_qo3d_solve (const char *
ham, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy, HMDT
zz)
¶HADT
mgl_qo3d_solve_c (const char *
ham, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy, HMDT
zz)
¶HMDT
mgl_qo3d_func (dual (*
ham)(mreal u, mreal x, mreal y, mreal z, mreal px, mreal py, mreal pz, void *par)
, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy, HMDT
zz)
¶HADT
mgl_qo3d_func_c (dual (*
ham)(mreal u, mreal x, mreal y, mreal z, mreal px, mreal py, mreal pz, void *par)
, HCDT
ini_re, HCDT
ini_im, HCDT
ray, mreal
r, mreal
k0, HMDT
xx, HMDT
yy, HMDT
zz)
¶Solves equation du/dt = i*k0*ham(p,q,v,x,y,z,|u|)[u], where p=-i/k0*d/dx, q=-i/k0*d/dy, v=-i/k0*d/dz are pseudo-differential operators (see mglPDE()
for details). Parameters ini_re, ini_im specify real and imaginary part of initial field distribution. Parameters ray set the reference ray, i.e. the ray around which the accompanied coordinate system will be maked. You may use, for example, the array created by ray function. Note, that the reference ray must be smooth enough to make accompanied coodrinates unambiguity. Otherwise errors in the solution may appear. If xx and yy and zz are non-zero then Cartesian coordinates for each point will be written into them. See also pde, qo2d. See PDE solving hints, for sample code and picture.
RES xdat ydat [zdat]
¶mglData
mglJacobian (const mglDataA &
x, const mglDataA &
y)
¶mglData
mglJacobian (const mglDataA &
x, const mglDataA &
y, const mglDataA &
z)
¶HMDT
mgl_jacobian_2d (HCDT
x, HCDT
y)
¶HMDT
mgl_jacobian_3d (HCDT
x, HCDT
y, HCDT
z)
¶Computes the Jacobian for transformation {i,j,k} to {x,y,z} where initial coordinates {i,j,k} are data indexes normalized in range [0,1]. The Jacobian is determined by formula det||dr_\alpha/d\xi_\beta|| where r={x,y,z} and \xi={i,j,k}. All dimensions must be the same for all data arrays. Data must be 3D if all 3 arrays {x,y,z} are specified or 2D if only 2 arrays {x,y} are specified.
RES xdat ydat
¶mglData
mglTriangulation (const mglDataA &
x, const mglDataA &
y)
¶HMDT
mgl_triangulation_2d (HCDT
x, HCDT
y)
¶Computes triangulation for arbitrary placed points with coordinates {x,y} (i.e. finds triangles which connect points). MathGL use s-hull code for triangulation. The sizes of 1st dimension must be equal for all arrays x.nx=y.nx
. Resulting array can be used in triplot or tricont functions for visualization of reconstructed surface. See Making regular data, for sample code and picture.
mglData
mglGSplineInit (const mglDataA &
x, const mglDataA &
y)
¶mglDataC
mglGSplineCInit (const mglDataA &
x, const mglDataA &
y)
¶HMDT
mgl_gspline_init (HCDT
x, HCDT
y)
¶HADT
mgl_gsplinec_init (HCDT
x, HCDT
y)
¶Prepare coefficients for global cubic spline interpolation.
mreal
mglGSpline (const mglDataA &
coef, mreal
dx, mreal *
d1=0
, mreal *
d2=0
)
¶dual
mglGSplineC (const mglDataA &
coef, mreal
dx, dual *
d1=0
, dual *
d2=0
)
¶mreal
mgl_gspline (HCDT
coef, mreal
dx, mreal *
d1, mreal *
d2)
¶dual
mgl_gsplinec (HCDT
coef, mreal
dx, dual *
d1, dual *
d2)
¶Evaluate global cubic spline (and its 1st and 2nd derivatives d1, d2 if they are not NULL
) using prepared coefficients coef at point dx+x0 (where x0 is 1st element of data x provided to mglGSpline*Init()
function).
RES dat num
[skip=20
]
¶mglData
mglIFS2d (const mglDataA &
dat, long
num, long
skip=20
)
¶HMDT
mgl_data_ifs_2d (HCDT
dat, long
num, long
skip)
¶Computes num points {x[i]=res[0,i], y[i]=res[1,i]} for fractal using iterated function system. Matrix dat is used for generation according the formulas
x[i+1] = dat[0,i]*x[i] + dat[1,i]*y[i] + dat[4,i]; y[i+1] = dat[2,i]*x[i] + dat[3,i]*y[i] + dat[5,i];
Value dat[6,i]
is used as weight factor for i-th row of matrix dat. At this first skip iterations will be omitted. Data array dat must have x-size greater or equal to 7. See also ifs3d, flame2d. See Sample ‘ifs2d’, for sample code and picture.
RES dat num
[skip=20
]
¶mglData
mglIFS3d (const mglDataA &
dat, long
num, long
skip=20
)
¶HMDT
mgl_data_ifs_3d (HCDT
dat, long
num, long
skip)
¶Computes num points {x[i]=res[0,i], y[i]=res[1,i], z[i]=res[2,i]} for fractal using iterated function system. Matrix dat is used for generation according the formulas
x[i+1] = dat[0,i]*x[i] + dat[1,i]*y[i] + dat[2,i]*z[i] + dat[9,i]; y[i+1] = dat[3,i]*x[i] + dat[4,i]*y[i] + dat[5,i]*z[i] + dat[10,i]; z[i+1] = dat[6,i]*x[i] + dat[7,i]*y[i] + dat[8,i]*z[i] + dat[11,i];
Value dat[12,i]
is used as weight factor for i-th row of matrix dat. At this first skip iterations will be omitted. Data array dat must have x-size greater or equal to 13. See also ifs2d. See Sample ‘ifs3d’, for sample code and picture.
RES 'fname' 'name' num
[skip=20
]
¶mglData
mglIFSfile (const char *
fname, const char *
name, long
num, long
skip=20
)
¶HMDT
mgl_data_ifs_file (const char *
fname, const char *
name, long
num, long
skip)
¶Reads parameters of IFS fractal named name from file fname and computes num points for this fractal. At this first skip iterations will be omitted. See also ifs2d, ifs3d.
IFS file may contain several records. Each record contain the name of fractal (‘binary’ in the example below) and the body of fractal, which is enclosed in curly braces {}. Symbol ‘;’ start the comment. If the name of fractal contain ‘(3D)’ or ‘(3d)’ then the 3d IFS fractal is specified. The sample below contain two fractals: ‘binary’ – usual 2d fractal, and ‘3dfern (3D)’ – 3d fractal. See also ifs2d, ifs3d.
binary { ; comment allowed here ; and here .5 .0 .0 .5 -2.563477 -0.000003 .333333 ; also comment allowed here .5 .0 .0 .5 2.436544 -0.000003 .333333 .0 -.5 .5 .0 4.873085 7.563492 .333333 } 3dfern (3D) { .00 .00 0 .0 .18 .0 0 0.0 0.00 0 0.0 0 .01 .85 .00 0 .0 .85 .1 0 -0.1 0.85 0 1.6 0 .85 .20 -.20 0 .2 .20 .0 0 0.0 0.30 0 0.8 0 .07 -.20 .20 0 .2 .20 .0 0 0.0 0.30 0 0.8 0 .07 }
RES dat func num
[skip=20
]
¶mglData
mglFlame2d (const mglDataA &
dat, const mglDataA &
func, long
num, long
skip=20
)
¶HMDT
mgl_data_flame_2d (HCDT
dat, HCDT
func, long
num, long
skip)
¶Computes num points {x[i]=res[0,i], y[i]=res[1,i]} for "flame" fractal using iterated function system. Array func define "flame" function identificator (func[0,i,j]), its weight (func[0,i,j]) and arguments (func[2 ... 5,i,j]). Matrix dat set linear transformation of coordinates before applying the function. The resulting coordinates are
xx = dat[0,i]*x[j] + dat[1,j]*y[i] + dat[4,j]; yy = dat[2,i]*x[j] + dat[3,j]*y[i] + dat[5,j]; x[j+1] = sum_i @var{func}[1,i,j]*@var{func}[0,i,j]_x(xx, yy; @var{func}[2,i,j],...,@var{func}[5,i,j]); y[j+1] = sum_i @var{func}[1,i,j]*@var{func}[0,i,j]_y(xx, yy; @var{func}[2,i,j],...,@var{func}[5,i,j]);
The possible function ids are: mglFlame2d_linear=0, mglFlame2d_sinusoidal, mglFlame2d_spherical, mglFlame2d_swirl, mglFlame2d_horseshoe,
mglFlame2d_polar, mglFlame2d_handkerchief,mglFlame2d_heart, mglFlame2d_disc, mglFlame2d_spiral,
mglFlame2d_hyperbolic, mglFlame2d_diamond, mglFlame2d_ex, mglFlame2d_julia, mglFlame2d_bent,
mglFlame2d_waves, mglFlame2d_fisheye, mglFlame2d_popcorn, mglFlame2d_exponential, mglFlame2d_power,
mglFlame2d_cosine, mglFlame2d_rings, mglFlame2d_fan, mglFlame2d_blob, mglFlame2d_pdj,
mglFlame2d_fan2, mglFlame2d_rings2, mglFlame2d_eyefish, mglFlame2d_bubble, mglFlame2d_cylinder,
mglFlame2d_perspective, mglFlame2d_noise, mglFlame2d_juliaN, mglFlame2d_juliaScope, mglFlame2d_blur,
mglFlame2d_gaussian, mglFlame2d_radialBlur, mglFlame2d_pie, mglFlame2d_ngon, mglFlame2d_curl,
mglFlame2d_rectangles, mglFlame2d_arch, mglFlame2d_tangent, mglFlame2d_square, mglFlame2d_blade,
mglFlame2d_secant, mglFlame2d_rays, mglFlame2d_twintrian, mglFlame2d_cross, mglFlame2d_disc2,
mglFlame2d_supershape, mglFlame2d_flower, mglFlame2d_conic, mglFlame2d_parabola, mglFlame2d_bent2,
mglFlame2d_bipolar, mglFlame2d_boarders, mglFlame2d_butterfly, mglFlame2d_cell, mglFlame2d_cpow,
mglFlame2d_curve, mglFlame2d_edisc, mglFlame2d_elliptic, mglFlame2d_escher, mglFlame2d_foci,
mglFlame2d_lazySusan, mglFlame2d_loonie, mglFlame2d_preBlur, mglFlame2d_modulus, mglFlame2d_oscope,
mglFlame2d_polar2, mglFlame2d_popcorn2, mglFlame2d_scry, mglFlame2d_separation, mglFlame2d_split,
mglFlame2d_splits, mglFlame2d_stripes, mglFlame2d_wedge, mglFlame2d_wedgeJulia, mglFlame2d_wedgeSph,
mglFlame2d_whorl, mglFlame2d_waves2, mglFlame2d_exp, mglFlame2d_log, mglFlame2d_sin,
mglFlame2d_cos, mglFlame2d_tan, mglFlame2d_sec, mglFlame2d_csc, mglFlame2d_cot,
mglFlame2d_sinh, mglFlame2d_cosh, mglFlame2d_tanh, mglFlame2d_sech, mglFlame2d_csch,
mglFlame2d_coth, mglFlame2d_auger, mglFlame2d_flux.
Value dat[6,i]
is used as weight factor for i-th row of matrix dat. At this first skip iterations will be omitted. Sizes of data arrays must be: dat.nx>=7, func.nx>=2 and func.nz=dat.ny. See also ifs2d, ifs3d. See Sample ‘flame2d’, for sample code and picture.