ARTS
2.3.1285(git:92a29ea9-dirty)
|
Declarations for the T-Matrix interface. More...
Go to the source code of this file.
Functions | |
void | calcSingleScatteringDataProperties (SingleScatteringData &ssd, ConstMatrixView ref_index_real, ConstMatrixView ref_index_imag, const Numeric equiv_radius=200, const Index np=-1, const Numeric axial_ratio=1.000001, const Numeric precision=0.001, const Index ndgs=2, const Index robust=0, const Index quiet=1) |
Calculate SingleScatteringData properties. More... | |
void | tmatrix_ampld_test (const Verbosity &verbosity) |
T-Matrix validation test. More... | |
void | tmatrix_tmd_test (const Verbosity &verbosity) |
T-Matrix validation test. More... | |
void | calc_ssp_random_test (const Verbosity &verbosity) |
Single scattering properties calculation for randomly oriented particles. More... | |
void | calc_ssp_fixed_test (const Verbosity &verbosity) |
Single scattering properties calculation for particles with fixed orientation. More... | |
void calc_ssp_fixed_test | ( | const Verbosity & | verbosity | ) |
Single scattering properties calculation for particles with fixed orientation.
Two cases are calculated. One with oblate particles which is equivalent to the following PyARTS case:
from PyARTS import arts_types from PyARTS import constants
params = {'ptype': constants.PTYPE_AZIMUTH_RND, 'f_grid': [230e9, 240e9], 'T_grid': [220, 250], 'za_grid': numpy.arange(0, 181, 10), 'aa_grid': numpy.arange(0, 181, 10), 'equiv_radius': 200, # equivalent volume radius 'NP':-1, # -1 for spheroid, -2 for cylinder, positive for chebyshev 'phase':'ice', 'mrr': numpy.array([[1.78031135, 1.78150475], [1.78037238, 1.78147686]]), 'mri': numpy.array([[0.00278706, 0.00507565], [0.00287245, 0.00523012]]), 'axial_ratio': 1.5} s = arts_types.SingleScatteringData(params) s.calc()
And one with prolate particles which is equivalent to the following PyARTS case:
from PyARTS import arts_types from PyARTS import constants
params = {'ptype': constants.PTYPE_AZIMUTH_RND, 'f_grid': [230e9, 240e9], 'T_grid': [220, 250], 'za_grid': numpy.arange(0, 181, 10), 'aa_grid': numpy.arange(0, 181, 10), 'equiv_radius': 200, # equivalent volume radius 'NP':-1, # -1 for spheroid, -2 for cylinder, positive for chebyshev 'phase':'ice', 'mrr': numpy.array([[1.78031135, 1.78150475], [1.78037238, 1.78147686]]), 'mri': numpy.array([[0.00278706, 0.00507565], [0.00287245, 0.00523012]]), 'axial_ratio': 0.7} s = arts_types.SingleScatteringData(params) s.calc()
Definition at line 1507 of file tmatrix.cc.
References SingleScatteringData::aa_grid, calcSingleScatteringDataProperties(), CREATE_OUT0, SingleScatteringData::f_grid, ConstVectorView::nelem(), nlinspace(), SingleScatteringData::ptype, PTYPE_AZIMUTH_RND, SingleScatteringData::T_grid, and SingleScatteringData::za_grid.
Referenced by TMatrixTest().
void calc_ssp_random_test | ( | const Verbosity & | verbosity | ) |
Single scattering properties calculation for randomly oriented particles.
Two cases are calculated. One with oblate particles which is equivalent to the following PyARTS case:
from PyARTS import arts_types
params = {'ptype': constants.PTYPE_AZIMUTH_RND, 'f_grid': [230e9, 240e9], 'T_grid': [220, 250], 'za_grid': numpy.arange(0, 181, 10), 'aa_grid': numpy.arange(0, 181, 10), 'equiv_radius': 200, # equivalent volume radius 'NP':-1, # -1 for spheroid, -2 for cylinder, positive for chebyshev 'phase':'ice', 'mrr': numpy.array([[1.78031135, 1.78150475], [1.78037238, 1.78147686]]), 'mri': numpy.array([[0.00278706, 0.00507565], [0.00287245, 0.00523012]]), 'axial_ratio': 1.5} s = arts_types.SingleScatteringData(params) s.calc()
And one with prolate particles which is equivalent to the following PyARTS case:
from PyARTS import arts_types
params = {'ptype': constants.PTYPE_AZIMUTH_RND, 'f_grid': [230e9, 240e9], 'T_grid': [220, 250], 'za_grid': numpy.arange(0, 181, 10), 'aa_grid': numpy.arange(0, 181, 10), 'equiv_radius': 200, # equivalent volume radius 'NP':-1, # -1 for spheroid, -2 for cylinder, positive for chebyshev 'phase':'ice', 'mrr': numpy.array([[1.78031135, 1.78150475], [1.78037238, 1.78147686]]), 'mri': numpy.array([[0.00278706, 0.00507565], [0.00287245, 0.00523012]]), 'axial_ratio': 0.7} s = arts_types.SingleScatteringData(params) s.calc()
Definition at line 1454 of file tmatrix.cc.
References SingleScatteringData::aa_grid, calcSingleScatteringDataProperties(), CREATE_OUT0, SingleScatteringData::f_grid, ConstVectorView::nelem(), nlinspace(), SingleScatteringData::ptype, PTYPE_TOTAL_RND, SingleScatteringData::T_grid, and SingleScatteringData::za_grid.
Referenced by TMatrixTest().
void calcSingleScatteringDataProperties | ( | SingleScatteringData & | ssd, |
ConstMatrixView | ref_index_real, | ||
ConstMatrixView | ref_index_imag, | ||
const Numeric | equiv_radius = 200 , |
||
const Index | np = -1 , |
||
const Numeric | axial_ratio = 1.000001 , |
||
const Numeric | precision = 0.001 , |
||
const Index | ndgs = 2 , |
||
const Index | robust = 0 , |
||
const Index | quiet = 1 |
||
) |
Calculate SingleScatteringData properties.
Port of calc_SSP function from PyARTS.
[in,out] | ssd | As input all grids must be set (f, T, za and aa). These grids are given by ssd are used to calculate pha_mat_data, ext_mat_data and abs_vec_data. Also ssd.ptype must be set as input. The output ssd has remaining fields set. |
[in] | ref_index_real | Vector with real parts of refractive index Number of rows must match elements in ssd.f_grid Number of cols must match elements in ssd.T_grid |
[in] | ref_index_imag | Vector with imaginary parts of refractive index |
[in] | equiv_radius | equivalent volume radius [micrometer] |
[in] | np | Particle type (-1 for spheroid, -2 for cylinder) |
[in] | axial_ratio | Axial ratio of particles |
[in] | precision | Accuracy of the computations |
[in] | ndgs | The number of division points in computing integrals over the particle surface. |
Definition at line 979 of file tmatrix.cc.
References SingleScatteringData::abs_vec_data, SingleScatteringData::ext_mat_data, SingleScatteringData::f_grid, ConstMatrixView::ncols(), ConstVectorView::nelem(), ConstMatrixView::nrows(), SingleScatteringData::pha_mat_data, PI, SingleScatteringData::ptype, PTYPE_TOTAL_RND, Tensor5::resize(), Tensor7::resize(), SPEED_OF_LIGHT, SingleScatteringData::T_grid, tmatrix_random_orientation(), and SingleScatteringData::za_grid.
Referenced by calc_ssp_fixed_test(), and calc_ssp_random_test().
void tmatrix_ampld_test | ( | const Verbosity & | verbosity | ) |
T-Matrix validation test.
Executes the standard test included with the double precision T-Matrix code for nonspherical particles in a fixed orientation. Should give the same as running the 3rdparty/tmatrix/tmatrix_ampld executable.
Definition at line 1297 of file tmatrix.cc.
References alpha, ampl_(), ampmat_to_phamat(), beta, CREATE_OUT0, and tmatrix_().
Referenced by TMatrixTest().
void tmatrix_tmd_test | ( | const Verbosity & | verbosity | ) |
T-Matrix validation test.
Executes the standard test included with the double precision T-Matrix code for randomly oriented nonspherical particles. Should give the same as running the 3rdparty/tmatrix/tmatrix_tmd executable.
Definition at line 1361 of file tmatrix.cc.
References CREATE_OUT0, VectorView::get_c_array(), and tmd_().
Referenced by TMatrixTest().