ARTS
2.2.66
|
Implementation of the T-Matrix interface. More...
#include "tmatrix.h"
#include <stdexcept>
#include <cstring>
#include <cmath>
#include "complex.h"
#include "messages.h"
#include "matpackI.h"
#include "make_vector.h"
#include "math_funcs.h"
#include "optproperties.h"
Go to the source code of this file.
Functions | |
void | calc_phamat (Matrix &z, const Index &nmax, const Numeric &lam, const Numeric &thet0, const Numeric &thet, const Numeric &phi0, const Numeric &phi, const Numeric &beta, const Numeric &alpha) |
void | ampmat_to_phamat (Matrix &z, const Complex &s11, const Complex &s12, const Complex &s21, const Complex &s22) |
Calculate phase matrix from amplitude matrix. More... | |
void | integrate_phamat_alpha10 (Matrix &phamat, const Index &nmax, const Numeric &lam, const Numeric &thet0, const Numeric &thet, const Numeric &phi0, const Numeric &phi, const Numeric &beta, const Numeric &alpha_1, const Numeric &alpha_2) |
Integrate phase matrix over particle orientation angle. More... | |
void | integrate_phamat_alpha6 (Matrix &phamat, const Index &nmax, const Numeric &lam, const Numeric &thet0, const Numeric &thet, const Numeric &phi0, const Numeric &phi, const Numeric &beta, const Numeric &alpha_1, const Numeric &alpha_2) |
Integrate phase matrix over particle orientation angle. More... | |
void | integrate_phamat_theta0_phi10 (Matrix &phamat, const Index &nmax, const Numeric &lam, const Numeric &thet0_1, const Numeric &thet0_2, const Numeric &thet, const Numeric &phi0, const Numeric &phi_1, const Numeric &phi_2, const Numeric &beta, const Numeric &alpha) |
Integrate phase matrix over angles thet0 and phi. More... | |
void | integrate_phamat_theta0_phi_alpha6 (Matrix &phamat, const Index &nmax, const Numeric &lam, const Numeric &thet0_1, const Numeric &thet0_2, const Numeric &thet, const Numeric &phi0, const Numeric &phi_1, const Numeric &phi_2, const Numeric &beta, const Numeric &alpha_1, const Numeric &alpha_2) |
Integrate phase matrix over angles thet0, phi and alpha. More... | |
void | tmd_ (const Numeric &rat, const Index &ndistr, const Numeric &axmax, const Index &npnax, const Numeric &b, const Numeric &gam, const Index &nkmax, const Numeric &eps, const Index &np, const Numeric &lam, const Numeric &mrr, const Numeric &mri, const Numeric &ddelt, const Index &npna, const Index &ndgs, const Numeric &r1rat, const Numeric &r2rat, const Index &quiet, Numeric &reff, Numeric &veff, Numeric &cext, Numeric &csca, Numeric &walb, Numeric &asymm, Numeric *f11, Numeric *f22, Numeric *f33, Numeric *f44, Numeric *f12, Numeric *f34, char *errmsg) |
T-matrix code for randomly oriented nonspherical particles. More... | |
void | tmatrix_ (const Numeric &rat, const Numeric &axi, const Index &np, const Numeric &lam, const Numeric &eps, const Numeric &mrr, const Numeric &mri, const Numeric &ddelt, const Index &quiet, Index &nmax, Numeric &csca, Numeric &cext, char *errmsg) |
T-matrix code for nonspherical particles in a fixed orientation. More... | |
void | ampl_ (const Index &nmax, const Numeric &lam, const Numeric &thet0, const Numeric &thet, const Numeric &phi0, const Numeric &phi, const Numeric &alpha, const Numeric &beta, Complex &s11, Complex &s12, Complex &s21, Complex &s22) |
T-matrix code for nonspherical particles in a fixed orientation. More... | |
void | avgtmatrix_ (const Index &nmax) |
Perform orientation averaging. More... | |
void | tmatrix_random_orientation (Numeric &cext, Numeric &csca, Vector &f11, Vector &f22, Vector &f33, Vector &f44, Vector &f12, Vector &f34, const Numeric equiv_radius, const Numeric aspect_ratio, const Index np, const Numeric lam, const Numeric ref_index_real, const Numeric ref_index_imag, const Numeric precision, const Index nza, const Index quiet=1) |
Calculate properties for randomly oriented particles. More... | |
void | tmatrix_fixed_orientation (Numeric &cext, Numeric &csca, Index &nmax, const Numeric equiv_radius, const Numeric aspect_ratio, const Index np, const Numeric lam, const Numeric ref_index_real, const Numeric ref_index_imag, const Numeric precision, const Index quiet=1) |
Calculate properties for particles in a fixed orientation. More... | |
void | calcSingleScatteringDataProperties (SingleScatteringData &ssd, ConstMatrixView ref_index_real, ConstMatrixView ref_index_imag, const Numeric equiv_radius, const Index np, const Numeric aspect_ratio, const Numeric precision) |
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 ampl_ | ( | const Index & | nmax, |
const Numeric & | lam, | ||
const Numeric & | thet0, | ||
const Numeric & | thet, | ||
const Numeric & | phi0, | ||
const Numeric & | phi, | ||
const Numeric & | alpha, | ||
const Numeric & | beta, | ||
Complex & | s11, | ||
Complex & | s12, | ||
Complex & | s21, | ||
Complex & | s22 | ||
) |
T-matrix code for nonspherical particles in a fixed orientation.
This is the interface to the T-Matrix ampl Fortran subroutine. It calculates the amplitude matrix. The T-Matrix is passed from tmatrix_() internally via common blocks.
See 3rdparty/tmatrix/tmd.lp.f for the complete documentation of the T-Matrix codes.
[in] | nmax | Iteration count. Calculated by tmatrix_() |
[in] | lam | Wavelength of incident light [microns] |
[in] | thet0 | Zenith angle of the incident beam [degrees] |
[in] | thet | zenith angle of the scattered beam in degrees |
[in] | phi0 | azimuth angle of the incident beam in degrees |
[in] | phi | azimuth angle of the scattered beam in degrees |
[in] | alpha | Euler angle (in degrees) specifying the orientation of the scattering particle relative to the laboratory reference frame |
[in] | beta | Euler angle (in degrees) specifying the orientation of the scattering particle relative to the laboratory reference frame |
[out] | s11 |
Definition at line 402 of file tmatrix.cc.
Referenced by calc_phamat(), calcSingleScatteringDataProperties(), and tmatrix_ampld_test().
void ampmat_to_phamat | ( | Matrix & | z, |
const Complex & | s11, | ||
const Complex & | s12, | ||
const Complex & | s21, | ||
const Complex & | s22 | ||
) |
Calculate phase matrix from amplitude matrix.
Ported from the T-Matrix Fortran code in 3rdparty/tmatrix/ampld.lp.f
[out] | z | Phase Matrix |
[in] | s11 | Calculated by ampl_() |
[in] | s12 | Calculated by ampl_() |
[in] | s21 | Calculated by ampl_() |
[in] | s22 | Calculated by ampl_() |
Definition at line 459 of file tmatrix.cc.
References Matrix::resize().
Referenced by calc_phamat(), and tmatrix_ampld_test().
void avgtmatrix_ | ( | const Index & | nmax | ) |
Perform orientation averaging.
This should be called after tmatrix_() for prolate particles. Data is passed from the tmatrix_() subroutine internally in the Fortran code via common blocks.
[in] | nmax | Iteration count. Calculated by tmatrix_() |
Definition at line 418 of file tmatrix.cc.
Referenced by calcSingleScatteringDataProperties().
void calc_phamat | ( | Matrix & | z, |
const Index & | nmax, | ||
const Numeric & | lam, | ||
const Numeric & | thet0, | ||
const Numeric & | thet, | ||
const Numeric & | phi0, | ||
const Numeric & | phi, | ||
const Numeric & | beta, | ||
const Numeric & | alpha | ||
) |
Definition at line 426 of file tmatrix.cc.
References ampl_(), and ampmat_to_phamat().
Referenced by calcSingleScatteringDataProperties().
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.PARTICLE_TYPE_HORIZ_AL, '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]]), 'aspect_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.PARTICLE_TYPE_HORIZ_AL, '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]]), 'aspect_ratio': 0.7} s = arts_types.SingleScatteringData(params) s.calc()
Definition at line 1409 of file tmatrix.cc.
References SingleScatteringData::aa_grid, SingleScatteringData::abs_vec_data, calcSingleScatteringDataProperties(), CREATE_OUT0, SingleScatteringData::ext_mat_data, SingleScatteringData::f_grid, joker, ConstVectorView::nelem(), nlinspace(), SingleScatteringData::particle_type, PARTICLE_TYPE_HORIZ_AL, SingleScatteringData::pha_mat_data, 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.PARTICLE_TYPE_HORIZ_AL, '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]]), 'aspect_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.PARTICLE_TYPE_HORIZ_AL, '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]]), 'aspect_ratio': 0.7} s = arts_types.SingleScatteringData(params) s.calc()
Definition at line 1351 of file tmatrix.cc.
References SingleScatteringData::aa_grid, SingleScatteringData::abs_vec_data, calcSingleScatteringDataProperties(), CREATE_OUT0, SingleScatteringData::ext_mat_data, SingleScatteringData::f_grid, joker, ConstVectorView::nelem(), nlinspace(), SingleScatteringData::particle_type, PARTICLE_TYPE_MACROS_ISO, SingleScatteringData::pha_mat_data, 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 | aspect_ratio = 1.000001 , |
||
const Numeric | precision = 0.001 |
||
) |
Calculate SingleScatteringData properties.
Port of calc_SSP function from PyARTS.
[in,out] | ssd | Grids given by ssd are used to calculate pha_mat_data, ext_mat_data and abs_vec_data |
[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] | phase | Particle phase ("ice"), currently unused |
[in] | aspect_ratio | Aspect ratio of particles |
[in] | precision | Accuracy of the computations |
Definition at line 927 of file tmatrix.cc.
References SingleScatteringData::aa_grid, SingleScatteringData::abs_vec_data, ampl_(), avgtmatrix_(), beta, calc_phamat(), SingleScatteringData::ext_mat_data, SingleScatteringData::f_grid, integrate_phamat_alpha10(), integrate_phamat_theta0_phi10(), integrate_phamat_theta0_phi_alpha6(), joker, ConstMatrixView::ncols(), ConstVectorView::nelem(), ConstMatrixView::nrows(), SingleScatteringData::particle_type, PARTICLE_TYPE_HORIZ_AL, PARTICLE_TYPE_MACROS_ISO, SingleScatteringData::pha_mat_data, PI, Tensor5::resize(), Tensor7::resize(), SPEED_OF_LIGHT, SingleScatteringData::T_grid, tmatrix_fixed_orientation(), tmatrix_random_orientation(), and SingleScatteringData::za_grid.
Referenced by calc_ssp_fixed_test(), calc_ssp_random_test(), and scat_data_arrayFromMeta().
void integrate_phamat_alpha10 | ( | Matrix & | phamat, |
const Index & | nmax, | ||
const Numeric & | lam, | ||
const Numeric & | thet0, | ||
const Numeric & | thet, | ||
const Numeric & | phi0, | ||
const Numeric & | phi, | ||
const Numeric & | beta, | ||
const Numeric & | alpha_1, | ||
const Numeric & | alpha_2 | ||
) |
Integrate phase matrix over particle orientation angle.
Performs a ten point Gauss–Legendre integration over the orientation angle (first Euler angle) of the particles from alpha1 to alpha2.
[out] | phamat | Integrated phase matrix |
[in] | nmax | FIXME OLE |
[in] | lam | See ampl_() |
[in] | thet0 | See ampl_() |
[in] | thet | See ampl_() |
[in] | phi0 | See ampl_() |
[in] | phi | See ampl_() |
[in] | beta | See ampl_() |
[in] | alpha_1 | See alpha in ampl_() |
[in] | alpha_2 | See alpha in ampl_() |
Definition at line 523 of file tmatrix.cc.
References Matrix::resize().
Referenced by calcSingleScatteringDataProperties().
void integrate_phamat_alpha6 | ( | Matrix & | phamat, |
const Index & | nmax, | ||
const Numeric & | lam, | ||
const Numeric & | thet0, | ||
const Numeric & | thet, | ||
const Numeric & | phi0, | ||
const Numeric & | phi, | ||
const Numeric & | beta, | ||
const Numeric & | alpha_1, | ||
const Numeric & | alpha_2 | ||
) |
Integrate phase matrix over particle orientation angle.
Performs a six point Gauss–Legendre integration over the orientation angle (first Euler angle) of the particles from alpha1 to alpha2.
[out] | phamat | Integrated phase matrix |
[in] | nmax | FIXME OLE |
[in] | lam | See ampl_() |
[in] | thet0 | See ampl_() |
[in] | thet | See ampl_() |
[in] | phi0 | See ampl_() |
[in] | phi | See ampl_() |
[in] | beta | See ampl_() |
[in] | alpha_1 | See alpha in ampl_() |
[in] | alpha_2 | See alpha in ampl_() |
Definition at line 576 of file tmatrix.cc.
References Matrix::resize().
void integrate_phamat_theta0_phi10 | ( | Matrix & | phamat, |
const Index & | nmax, | ||
const Numeric & | lam, | ||
const Numeric & | thet0_1, | ||
const Numeric & | thet0_2, | ||
const Numeric & | thet, | ||
const Numeric & | phi0, | ||
const Numeric & | phi_1, | ||
const Numeric & | phi_2, | ||
const Numeric & | beta, | ||
const Numeric & | alpha | ||
) |
Integrate phase matrix over angles thet0 and phi.
Performs a ten point Gauss–Legendre integration over the zenith angle of the incident beam (thet0) and the azimuth angle of the scattered beam (phi).
[out] | phamat | Integrated phase matrix |
[in] | nmax | FIXME OLE |
[in] | lam | See ampl_() |
[in] | thet0_1 | See thet0 in ampl_() |
[in] | thet0_2 | See thet0 in ampl_() |
[in] | thet | See ampl_() |
[in] | phi0 | See ampl_() |
[in] | phi_1 | See phi in ampl_() |
[in] | phi_2 | See phi in ampl_() |
[in] | beta | See ampl_() |
[in] | alpha | See in ampl_() |
Definition at line 631 of file tmatrix.cc.
References PI, and Matrix::resize().
Referenced by calcSingleScatteringDataProperties().
void integrate_phamat_theta0_phi_alpha6 | ( | Matrix & | phamat, |
const Index & | nmax, | ||
const Numeric & | lam, | ||
const Numeric & | thet0_1, | ||
const Numeric & | thet0_2, | ||
const Numeric & | thet, | ||
const Numeric & | phi0, | ||
const Numeric & | phi_1, | ||
const Numeric & | phi_2, | ||
const Numeric & | beta, | ||
const Numeric & | alpha_1, | ||
const Numeric & | alpha_2 | ||
) |
Integrate phase matrix over angles thet0, phi and alpha.
Performs a six point Gauss–Legendre integration over the zenith angle of the incident beam (thet0), the azimuth angle of the scattered beam (phi) and the orientation angle (first Euler angle) of the particles.
[out] | phamat | Integrated phase matrix |
[in] | nmax | FIXME OLE |
[in] | lam | See ampl_() |
[in] | thet0_1 | See thet0 in ampl_() |
[in] | thet0_2 | See thet0 in ampl_() |
[in] | thet | See ampl_() |
[in] | phi0 | See ampl_() |
[in] | phi_1 | See phi in ampl_() |
[in] | phi_2 | See phi in ampl_() |
[in] | beta | See ampl_() |
[in] | alpha_1 | See alpha in ampl_() |
[in] | alpha_2 | See alpha in ampl_() |
Definition at line 714 of file tmatrix.cc.
References PI, and Matrix::resize().
Referenced by calcSingleScatteringDataProperties().
void tmatrix_ | ( | const Numeric & | rat, |
const Numeric & | axi, | ||
const Index & | np, | ||
const Numeric & | lam, | ||
const Numeric & | eps, | ||
const Numeric & | mrr, | ||
const Numeric & | mri, | ||
const Numeric & | ddelt, | ||
const Index & | quiet, | ||
Index & | nmax, | ||
Numeric & | csca, | ||
Numeric & | cext, | ||
char * | errmsg | ||
) |
T-matrix code for nonspherical particles in a fixed orientation.
This is the interface to the T-Matrix tmatrix Fortran subroutine. It calculates extinction and scattering cross section per particle. The T-Matrix is calculated internally and used accessed later by ampl_() via common blocks.
See 3rdparty/tmatrix/ampld.lp.f for the complete documentation of the T-Matrix codes.
[in] | rat | 1 - particle size is specified in terms of the equal-volume-sphere radius != 1 - particle size is specified in terms of the equal-surface-area-sphere radius |
[in] | axi | Equivalent-sphere radius [microns] |
[in] | np | Shape of the particles For spheroids NP=-1 and EPS is the ratio of the horizontal to rotational axes. EPS is larger than 1 for oblate spheroids and smaller than 1 for prolate spheroids. For cylinders NP=-2 and EPS is the ratio of the diameter to the length. |
[in] | lam | Wavelength of light [microns] |
[in] | eps | Shape of the particles (see np above) |
[in] | mrr | Vector with real parts of refractive index |
[in] | mri | Vector with imaginary parts of refractive index |
[in] | ddelt | Accuracy of the computations |
[in] | quiet | 0 = Verbose output from Fortran code, 1 = silent |
[out] | nmax | Iteration count |
[out] | cext | Extinction cross section per particle |
[out] | csca | Scattering cross section per particle |
[out] | errmsg | Error message string from Fortran code |
Definition at line 383 of file tmatrix.cc.
Referenced by tmatrix_ampld_test(), and tmatrix_fixed_orientation().
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 1189 of file tmatrix.cc.
References ampl_(), ampmat_to_phamat(), beta, CREATE_OUT0, and tmatrix_().
Referenced by TMatrixTest().
void tmatrix_fixed_orientation | ( | Numeric & | cext, |
Numeric & | csca, | ||
Index & | nmax, | ||
const Numeric | equiv_radius, | ||
const Numeric | aspect_ratio, | ||
const Index | np, | ||
const Numeric | lam, | ||
const Numeric | ref_index_real, | ||
const Numeric | ref_index_imag, | ||
const Numeric | precision, | ||
const Index | quiet = 1 |
||
) |
Calculate properties for particles in a fixed orientation.
This is a simplified interface to the tmatrix_() function for randomly oriented particles based on the PyARTS function tmat_fxd
[out] | cext | Extinction cross section per particle |
[out] | csca | Scattering cross section per particle |
[out] | nmax | FIXME OLE |
[in] | equiv_radius | See parameter axmax in tmd_() |
[in] | aspect_ratio | See parameter eps in tmd_() |
[in] | np | See tmd_() |
[in] | lam | See tmd_() |
[in] | ref_index_real | See parameter mrr in tmd_() |
[in] | ref_index_imag | See parameter mri in tmd_() |
[in] | precision | See parameter ddelt in tmd_() |
[in] | quiet | See tmd_() |
Definition at line 896 of file tmatrix.cc.
References tmatrix_().
Referenced by calcSingleScatteringDataProperties().
void tmatrix_random_orientation | ( | Numeric & | cext, |
Numeric & | csca, | ||
Vector & | f11, | ||
Vector & | f22, | ||
Vector & | f33, | ||
Vector & | f44, | ||
Vector & | f12, | ||
Vector & | f34, | ||
const Numeric | equiv_radius, | ||
const Numeric | aspect_ratio, | ||
const Index | np, | ||
const Numeric | lam, | ||
const Numeric | ref_index_real, | ||
const Numeric | ref_index_imag, | ||
const Numeric | precision, | ||
const Index | nza, | ||
const Index | quiet = 1 |
||
) |
Calculate properties for randomly oriented particles.
This is a simplified interface to the tmd_() function for randomly oriented particles based on the PyARTS function tmat_rnd
[out] | cext | Extinction cross section per particle |
[out] | csca | Scattering cross section per particle |
[out] | f11 | See tmd_() |
[out] | f22 | See tmd_() |
[out] | f33 | See tmd_() |
[out] | f44 | See tmd_() |
[out] | f12 | See tmd_() |
[out] | f34 | See tmd_() |
[in] | equiv_radius | See parameter axmax in tmd_() |
[in] | aspect_ratio | See parameter eps in tmd_() |
[in] | np | See tmd_() |
[in] | lam | See tmd_() |
[in] | ref_index_real | See parameter mrr in tmd_() |
[in] | ref_index_imag | See parameter mri in tmd_() |
[in] | precision | See parameter ddelt in tmd_() |
[in] | nza | See parameter npna in tmd_() |
[in] | quiet | See tmd_() |
Definition at line 802 of file tmatrix.cc.
References VectorView::get_c_array(), Vector::resize(), and tmd_().
Referenced by calcSingleScatteringDataProperties().
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 1256 of file tmatrix.cc.
References CREATE_OUT0, VectorView::get_c_array(), and tmd_().
Referenced by TMatrixTest().
void tmd_ | ( | const Numeric & | rat, |
const Index & | ndistr, | ||
const Numeric & | axmax, | ||
const Index & | npnax, | ||
const Numeric & | b, | ||
const Numeric & | gam, | ||
const Index & | nkmax, | ||
const Numeric & | eps, | ||
const Index & | np, | ||
const Numeric & | lam, | ||
const Numeric & | mrr, | ||
const Numeric & | mri, | ||
const Numeric & | ddelt, | ||
const Index & | npna, | ||
const Index & | ndgs, | ||
const Numeric & | r1rat, | ||
const Numeric & | r2rat, | ||
const Index & | quiet, | ||
Numeric & | reff, | ||
Numeric & | veff, | ||
Numeric & | cext, | ||
Numeric & | csca, | ||
Numeric & | walb, | ||
Numeric & | asymm, | ||
Numeric * | f11, | ||
Numeric * | f22, | ||
Numeric * | f33, | ||
Numeric * | f44, | ||
Numeric * | f12, | ||
Numeric * | f34, | ||
char * | errmsg | ||
) |
T-matrix code for randomly oriented nonspherical particles.
This is the interface to the T-Matrix tmd Fortran subroutine.
See 3rdparty/tmatrix/tmd.lp.f for the complete documentation of the T-Matrix codes.
[in] | rat | 1 - particle size is specified in terms of the equal-volume-sphere radius != 1 - particle size is specified in terms of the equal-surface-area-sphere radius |
[in] | ndistr | Specifies the distribution of equivalent-sphere radii NDISTR = 1 - modified gamma distribution [Eq. (40) of Ref. 7] AXI=alpha B=r_c GAM=gamma NDISTR = 2 - log-normal distribution [Eq. 41) of Ref. 7] AXI=r_g B=[ln(sigma_g)]**2 NDISTR = 3 - power law distribution [Eq. (42) of Ref. 7] AXI=r_eff (effective radius) B=v_eff (effective variance) Parameters R1 and R2 (see below) are calculated automatically for given AXI and B NDISTR = 4 - gamma distribution [Eq. (39) of Ref. 7] AXI=a B=b NDISTR = 5 - modified power law distribution [Eq. (24) in M. I. Mishchenko et al., Bidirectional reflectance of flat, optically thick particulate laters: an efficient radiative transfer solution and applications to snow and soil surfaces, J. Quant. Spectrosc. Radiat. Transfer, Vol. 63, 409-432 (1999)]. B=alpha |
[in] | axmax | The code computes NPNAX size distributions of the same type and with the same values of B and GAM in one run. The parameter AXI varies from AXMAX to AXMAX/NPNAX in steps of AXMAX/NPNAX. To compute a single size distribution, use NPNAX=1 and AXMAX equal to AXI of this size distribution. |
[in] | npnax | See axmax above |
[in] | b | See ndistr above |
[in] | gam | See ndistr above |
[in] | nkmax | NKMAX<=988 is such that NKMAX+2 is the number of Gaussian quadrature points used in integrating over the size distribution for particles with AXI=AXMAX. For particles with AXI=AXMAX-AXMAX/NPNAX, AXMAX-2*AXMAX/NPNAX, etc. the number of Gaussian points linearly decreases. For the modified power law distribution, the number of integration points on the interval [0,R1] is also equal to NKMAX. |
[in] | eps | Shape of the particles For spheroids NP=-1 and EPS is the ratio of the horizontal to rotational axes. EPS is larger than 1 for oblate spheroids and smaller than 1 for prolate spheroids. For cylinders NP=-2 and EPS is the ratio of the diameter to the length. For Chebyshev particles NP must be positive and is the degree of the Chebyshev polynomial, while EPS is the deformation parameter. |
[in] | np | Shape of the particles (see eps above) |
[in] | lam | Wavelength of light [microns] |
[in] | mrr | Vector with real parts of refractive index |
[in] | mri | Vector with imaginary parts of refractive index |
[in] | ddelt | Accuracy of the computations |
[in] | npna | Number of equidistant scattering angles (from 0 to 180 deg) for which the scattering matrix is calculated |
[in] | ndgs | Parameter controlling the number of division points in computing integrals over the particle surface. For compact particles, the recommended value is 2. For highly aspherical particles larger values (3, 4,...) may be necessary to obtain convergence. The code does not check convergence over this parameter. Therefore, control comparisons of results obtained with different NDGS-values are recommended. |
[in] | r1rat | |
[in] | r2rat | |
[in] | quiet | 0 = Verbose output from Fortran code, 1 = silent |
[out] | reff | Effective radius of the size distribution |
[out] | veff | Effective variance of the size distribution |
[out] | cext | Extinction cross section per particle |
[out] | csca | Scattering cross section per particle |
[out] | walb | Single scattering albedo |
[out] | asymm | Asymmetry parameter of the phase function |
[out] | f11 | Elements of the normalized scattering matrix versus scattering angle |
[out] | f22 | Elements of the normalized scattering matrix versus scattering angle |
[out] | f33 | Elements of the normalized scattering matrix versus scattering angle |
[out] | f44 | Elements of the normalized scattering matrix versus scattering angle |
[out] | f12 | Elements of the normalized scattering matrix versus scattering angle |
[out] | f34 | Elements of the normalized scattering matrix versus scattering angle |
[out] | errmsg | Error message string from Fortran code |
Definition at line 346 of file tmatrix.cc.
Referenced by tmatrix_random_orientation(), and tmatrix_tmd_test().