[ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
v1 ...
Creates new variable with name dat and fills it by numeric values of command arguments v1 ...
. Command can create one-dimensional and two-dimensional arrays with arbitrary values. For creating 2d array the user should use delimiter ‘|’ which means that the following values lie in next row. Array sizes are [maximal of row sizes * number of rows]. For example, command list 1 | 2 3
creates the array [1 0; 2 3]. Note, that the maximal number of arguments is 1000.
Creates new variable with name dat and fills it by data values of arrays of command arguments d1 .... Command can create two-dimensional or three-dimensional (if arrays in arguments are 2d arrays) arrays with arbitrary values. Minor dimensions of all arrays in arguments should be equal to dimensions of first array d1. In the opposite case the argument will be ignored. Note, that the maximal number of arguments is 1000.
mglData
: void
Set (const float *
A, int
NX, int
NY=1
, int
NZ=1
)mglData
: void
Set (const double *
A, int
NX, int
NY=1
, int
NZ=1
)void
mgl_data_set_float (HMDT
dat, const mreal *
A, int
NX, int
NY, int
NZ)void
mgl_data_set_double (HMDT
dat, const double *
A, int
NX, int
NY, int
NZ)mglDataC
: void
Set (const float *
A, int
NX, int
NY=1
, int
NZ=1
)mglDataC
: void
Set (const double *
A, int
NX, int
NY=1
, int
NZ=1
)mglDataC
: void
Set (const dual *
A, int
NX, int
NY=1
, int
NZ=1
)void
mgl_datac_set_float (HADT
dat, const mreal *
A, int
NX, int
NY, int
NZ)void
mgl_datac_set_double (HADT
dat, const double *
A, int
NX, int
NY, int
NZ)void
mgl_datac_set_complex (HADT
dat, const dual *
A, int
NX, int
NY, int
NZ)Allocates memory and copies the data from the flat float*
or double*
array.
mglData
: void
Set (const float **
A, int
N1, int
N2)mglData
: void
Set (const double **
A, int
N1, int
N2)void
mgl_data_set_mreal2 (HMDT
dat, const mreal **
A, int
N1, int
N2)void
mgl_data_set_double2 (HMDT
dat, const double **
A, int
N1, int
N2)Allocates memory and copies the data from the float**
or double**
array with dimensions N1, N2, i.e. from array defined as mreal a[N1][N2];
.
mglData
: void
Set (const float ***
A, int
N1, int
N2)mglData
: void
Set (const double ***
A, int
N1, int
N2)void
mgl_data_set_mreal3 (HMDT
dat, const mreal ***
A, int
N1, int
N2)void
mgl_data_set_double3 (HMDT
dat, const double ***
A, int
N1, int
N2)Allocates memory and copies the data from the float***
or double***
array with dimensions N1, N2, N3, i.e. from array defined as mreal a[N1][N2][N3];
.
mglData
: void
Set (gsl_vector *
v)mglDataC
: void
Set (gsl_vector *
v)void
mgl_data_set_vector (HMDT
dat, gsl_vector *
v)void
mgl_datac_set_vector (HADT
dat, gsl_vector *
v)Allocates memory and copies the data from the gsl_vector *
structure.
mglData
: void
Set (gsl_matrix *
m)mglDataC
: void
Set (gsl_matrix *
m)void
mgl_data_set_matrix (HMDT
dat, gsl_matrix *
m)void
mgl_datac_set_matrix (HADT
dat, gsl_matrix *
m)Allocates memory and copies the data from the gsl_matrix *
structure.
mglData
: void
Set (const mglDataA &
from)mglData
: void
Set (HCDT
from)void
mgl_data_set (HMDT
dat, HCDT
from)mglDataC
: void
Set (const mglDataA &
from)mglDataC
: void
Set (HCDT
from)void
mgl_datac_set (HADT
dat, HCDT
from)Copies the data from mglData
(or mglDataA
) instance from.
mglDataC
: void
Set (const mglDataA &
re, const mglDataA &
im)mglDataC
: void
Set (HCDT
re, HCDT
im)mglDataC
: void
SetAmpl (HCDT
ampl, const mglDataA &
phase)void
mgl_datac_set_ri (HADT
dat, HCDT
re, HCDT
im)void
mgl_datac_set_ap (HADT
dat, HCDT
ampl, HCDT
phase)Copies the data from mglData
instances for real and imaginary parts of complex data arrays.
mglData
: void
Set (const std::vector<int> &
d)mglDataC
: void
Set (const std::vector<int> &
d)mglData
: void
Set (const std::vector<float> &
d)mglDataC
: void
Set (const std::vector<float> &
d)mglData
: void
Set (const std::vector<double> &
d)mglDataC
: void
Set (const std::vector<double> &
d)mglDataC
: void
Set (const std::vector<dual> &
d)Allocates memory and copies the data from the std::vector<T>
array.
mglData
: void
Set (const char *
str, int
NX, int
NY=1
, int
NZ=1
)void
mgl_data_set_values (const char *
str, int
NX, int
NY, int
NZ)mglDataC
: void
Set (const char *
str, int
NX, int
NY=1
, int
NZ=1
)void
mgl_datac_set_values (const char *
str, int
NX, int
NY, int
NZ)Allocates memory and scanf the data from the string.
mglData
: void
Link (mglData &
from)mglData
: void
Link (mreal *
A, int
NX, int
NY=1
, int
NZ=1
)void
mgl_data_link (HMDT
dat, mreal *
A, int
NX, int
NY, int
NZ)mglDataC
: void
Link (mglDataC &
from)mglDataC
: void
Link (dual *
A, int
NX, int
NY=1
, int
NZ=1
)void
mgl_datac_link (HADT
dat, dual *
A, int
NX, int
NY, int
NZ)Links external data array, i.e. don’t delete this array at exit.
num v1 [v2=nan]
Creates new variable with name dat for one-dimensional array of size num. Array elements are equidistantly distributed in range [v1, v2]. If v2=nan
then v2=v1 is used.
mglData
: void
Fill (mreal
v1, mreal
v2, char
dir='x'
)mglDataC
: void
Fill (dual
v1, dual
v2, char
dir='x'
)void
mgl_data_fill (HMDT
dat, mreal
v1, mreal
v2, char
dir)void
mgl_datac_fill (HADT
dat, dual
v1, dual
v2, char
dir)Equidistantly fills the data values to range [v1, v2] in direction dir={‘x’,‘y’,‘z’}.
mglData
: void
Fill (HMGL
gr, const char *
eq, const char *
opt=""
)mglData
: void
Fill (HMGL
gr, const char *
eq, const mglDataA &
vdat, const char *
opt=""
)mglData
: void
Fill (HMGL
gr, const char *
eq, const mglDataA &
vdat, const mglDataA &
wdat, const char *
opt=""
)mglDataC
: void
Fill (HMGL
gr, const char *
eq, const char *
opt=""
)mglDataC
: void
Fill (HMGL
gr, const char *
eq, const mglDataA &
vdat, const char *
opt=""
)mglDataC
: void
Fill (HMGL
gr, const char *
eq, const mglDataA &
vdat, const mglDataA &
wdat, const char *
opt=""
)void
mgl_data_fill_eq (HMGL
gr, HMDT
dat, const char *
eq, HCDT
vdat, HCDT
wdat, const char *
opt)void
mgl_datac_fill_eq (HMGL
gr, HADT
dat, const char *
eq, HCDT
vdat, HCDT
wdat, const char *
opt)Fills the value of array according to the formula in string eq. Formula is an arbitrary expression depending on variables ‘x’, ‘y’, ‘z’, ‘u’, ‘v’, ‘w’. Coordinates ‘x’, ‘y’, ‘z’ are supposed to be normalized in axis range of canvas gr (in difference from Modify
functions). Variable ‘u’ is the original value of the array. Variables ‘v’ and ‘w’ are values of vdat, wdat which can be NULL
(i.e. can be omitted).
dim=0
]mglData
: void
Modify (const char *
eq, int
dim=0
)mglData
: void
Modify (const char *
eq, const mglDataA &
v)mglData
: void
Modify (const char *
eq, const mglDataA &
v, const mglDataA &
w)mglDataC
: void
Modify (const char *
eq, int
dim=0
)mglDataC
: void
Modify (const char *
eq, const mglDataA &
v)mglDataC
: void
Modify (const char *
eq, const mglDataA &
v, const mglDataA &
w)void
mgl_data_modify (HMDT
dat, const char *
eq, int
dim)void
mgl_data_modify_vw (HMDT
dat, const char *
eq, HCDT
v, HCDT
w)void
mgl_datac_modify (HADT
dat, const char *
eq, int
dim)void
mgl_datac_modify_vw (HADT
dat, const char *
eq, HCDT
v, HCDT
w)The same as previous ones but coordinates ‘x’, ‘y’, ‘z’ are supposed to be normalized in range [0,1]. If dim>0 is specified then modification will be fulfilled only for slices >=dim.
mglData
: void
FillSample (const char *
how)void
mgl_data_fill_sample (HMDT
a, const char *
how)Fills data by ’x’ or ’k’ samples for Hankel (’h’) or Fourier (’f’) transform.
mglData
: mglData
Grid (HMGL
gr, const mglDataA &
x, const mglDataA &
y, const mglDataA &
z, const char *
opt=""
)void
mgl_data_grid (HMGL
gr, HMDT
u, HCDT
x, HCDT
y, HCDT
z, const char *
opt)Fills the value of array according to the linear interpolation of triangulated surface, found for arbitrary placed points ‘x’, ‘y’, ‘z’. NAN value is used for grid points placed outside of triangulated surface.
val [i=: j=: k=:]
mglData
: void
Put (mreal
val, int
i=-1
, int
j=-1
, int
k=-1
)void
mgl_data_put_val (HMDT
a, mreal
val, int
i, int
j, int
k)Sets value(s) of array a[i, j, k] = val. Negative indexes i, j, k=-1 set the value val to whole range in corresponding direction(s). For example, Put(val,-1,0,-1);
sets a[i,0,j]=val for i=0...(nx-1), j=0...(nz-1).
i=: j=: k=:
]mglData
: void
Put (const mglDataA &
v, int
i=-1
, int
j=-1
, int
k=-1
)void
mgl_data_put_dat (HMDT
a, HCDT
v, int
i, int
j, int
k)Copies value(s) from array v to the range of original array. Negative indexes i, j, k=-1 set the range in corresponding direction(s). At this minor dimensions of array v should be large than corresponding dimensions of this array. For example, Put(v,-1,0,-1);
sets a[i,0,j]=v.ny>nz ? v[i,j] : v[i], where i=0...(nx-1), j=0...(nz-1) and condition v.nx>=nx is true.
mglData
: void
SetColumnId (const char *
ids)void
mgl_data_set_id (HMDT
a, const char *
ids)mglDataC
: void
SetColumnId (const char *
ids)void
mgl_datac_set_id (HADT
a, const char *
ids)Sets the symbol ids for data columns. The string should contain one symbol ’a’...’z’ per column. These ids are used in column.
[ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] |
This document was generated by Build Daemon user on December 22, 2013 using texi2html 1.82.