ARTS  2.3.1285(git:92a29ea9-dirty)
montecarlo.cc File Reference

functions used by MCGeneral More...

#include <cfloat>
#include <sstream>
#include "auto_md.h"
#include "geodetic.h"
#include "mc_interp.h"
#include "montecarlo.h"

Go to the source code of this file.

Functions

Numeric ext_I (const Numeric &I, const Numeric &Q, const Numeric &kI, const Numeric &kQ, const Numeric &s)
 ext_I. More...
 
void brent_zero (Numeric &sb, const Numeric &a, const Numeric &b, const Numeric &t, const Numeric &rn, const Numeric &I, const Numeric &Q, const Numeric &KI, const Numeric &KQ)
 brent_zero. More...
 
void clear_rt_vars_at_gp (Workspace &ws, MatrixView ext_mat_mono, VectorView abs_vec_mono, Numeric &temperature, const Agenda &propmat_clearsky_agenda, const Numeric &f_mono, const GridPos &gp_p, const GridPos &gp_lat, const GridPos &gp_lon, ConstVectorView p_grid, ConstTensor3View t_field, ConstTensor4View vmr_field)
 clear_rt_vars_at_gp. More...
 
void cloudy_rt_vars_at_gp (Workspace &ws, MatrixView ext_mat_mono, VectorView abs_vec_mono, VectorView pnd_vec, Numeric &temperature, const Agenda &propmat_clearsky_agenda, const Index stokes_dim, const Index f_index, const Vector &f_grid, const GridPos &gp_p, const GridPos &gp_lat, const GridPos &gp_lon, ConstVectorView p_grid_cloud, ConstTensor3View t_field_cloud, ConstTensor4View vmr_field_cloud, const Tensor4 &pnd_field, const ArrayOfArrayOfSingleScatteringData &scat_data, const ArrayOfIndex &cloudbox_limits, const Vector &rte_los)
 cloudy_rt_vars_at_gp. More...
 
void cloud_atm_vars_by_gp (VectorView pressure, VectorView temperature, MatrixView vmr, MatrixView pnd, const ArrayOfGridPos &gp_p, const ArrayOfGridPos &gp_lat, const ArrayOfGridPos &gp_lon, const ArrayOfIndex &cloudbox_limits, ConstVectorView p_grid_cloud, ConstTensor3View t_field_cloud, ConstTensor4View vmr_field_cloud, ConstTensor4View pnd_field)
 cloud_atm_vars_by_gp. More...
 
void get_ppath_transmat (Workspace &ws, MatrixView &trans_mat, const Ppath &ppath, const Agenda &propmat_clearsky_agenda, const Index stokes_dim, const Index f_index, const Vector &f_grid, const Vector &p_grid, const Tensor3 &t_field, const Tensor4 &vmr_field, const ArrayOfIndex &cloudbox_limits, const Tensor4 &pnd_field, const ArrayOfArrayOfSingleScatteringData &scat_data, const Verbosity &verbosity)
 get_ppath_transmat. More...
 
bool is_anyptype_nonTotRan (const ArrayOfArrayOfSingleScatteringData &scat_data)
 is_anyptype_nonTotRan. More...
 
void mcPathTraceGeneral (Workspace &ws, MatrixView evol_op, Vector &abs_vec_mono, Numeric &temperature, MatrixView ext_mat_mono, Rng &rng, Vector &rte_pos, Vector &rte_los, Vector &pnd_vec, Numeric &g, Ppath &ppath_step, Index &termination_flag, bool &inside_cloud, const Agenda &ppath_step_agenda, const Numeric &ppath_lmax, const Numeric &ppath_lraytrace, const Numeric &taustep_limit, const Agenda &propmat_clearsky_agenda, const Index stokes_dim, const Index f_index, const Vector &f_grid, const Vector &p_grid, const Vector &lat_grid, const Vector &lon_grid, const Tensor3 &z_field, const Vector &refellipsoid, const Matrix &z_surface, const Tensor3 &t_field, const Tensor4 &vmr_field, const ArrayOfIndex &cloudbox_limits, const Tensor4 &pnd_field, const ArrayOfArrayOfSingleScatteringData &scat_data, const Verbosity &verbosity)
 mcPathTraceGeneral. More...
 
void mcPathTraceRadar (Workspace &ws, MatrixView evol_op, Vector &abs_vec_mono, Numeric &temperature, MatrixView ext_mat_mono, Rng &rng, Vector &rte_pos, Vector &rte_los, Vector &pnd_vec, Numeric &stot, Numeric &ttot, Ppath &ppath_step, Index &termination_flag, bool &inside_cloud, const Agenda &ppath_step_agenda, const Numeric &ppath_lmax, const Numeric &ppath_lraytrace, const Agenda &propmat_clearsky_agenda, const bool &anyptype_nonTotRan, const Index stokes_dim, const Index f_index, const Vector &f_grid, const Vector &Iprop, const Vector &p_grid, const Vector &lat_grid, const Vector &lon_grid, const Tensor3 &z_field, const Vector &refellipsoid, const Matrix &z_surface, const Tensor3 &t_field, const Tensor4 &vmr_field, const ArrayOfIndex &cloudbox_limits, const Tensor4 &pnd_field, const ArrayOfArrayOfSingleScatteringData &scat_data, const Verbosity &verbosity)
 mcPathTraceRadar. More...
 
void Sample_los (VectorView new_rte_los, Numeric &g_los_csc_theta, MatrixView Z, Rng &rng, ConstVectorView rte_los, const ArrayOfArrayOfSingleScatteringData &scat_data, const Index f_index, const Index stokes_dim, ConstVectorView pnd_vec, ConstVectorView Z11maxvector, const Numeric Csca, const Numeric rtp_temperature, const Index t_interp_order)
 Sample_los. More...
 
void Sample_los_uniform (VectorView new_rte_los, Rng &rng)
 Sample_los_uniform. More...
 

Variables

const Numeric SPEED_OF_LIGHT
 

Detailed Description

functions used by MCGeneral

Author
Cory Davis cory@.nosp@m.met..nosp@m.ed.ac.nosp@m..uk
Date
2003-06-19
Author
Cory Davis cory@.nosp@m.met..nosp@m.ed.ac.nosp@m..uk
Date
2003-06-19

FIXMEDOC ***: mcPathTraceRadar, opt_propExtract, Sample_los_uniform (JUST author) cloudy_rt_vars_at_gp (JUST DATE CHECK), clear_rt_vars_at_gp(JUST DATE CHECK). FIXMEDOC ***: In .cc some root-finding helper functions (for MCRadar) that don't need visibility outside .cc, meaning that are not in .h JUST CLARIFICATION => brent_zero, ext_I. Further docomentation ???

Definition in file montecarlo.cc.

Function Documentation

◆ brent_zero()

void brent_zero ( Numeric sb,
const Numeric a,
const Numeric b,
const Numeric t,
const Numeric rn,
const Numeric I,
const Numeric Q,
const Numeric KI,
const Numeric KQ 
)

brent_zero.

This function seeks the root of a function F(X) in an interval [A,B].

Discussion: The interval [A,B] must be a change of sign interval for F. That is, F(A) and F(B) must be of opposite signs. Then assuming that F is continuous implies the existence of at least one value C between A and B for which F(C) = 0.

The location of the zero is determined to within an accuracy of 6 * MACHEPS * abs ( C ) + 2 * T.

Thanks to Thomas Secretin for pointing out a transcription error in the setting of the value of P, 11 February 2013.

Modifications by Ian S. Adams, U.S. Naval Research Laboratory to conform to ARTS and to hardcode function for root finding while passing in mulitple args for function.

Licensing:

This code is distributed under the GNU LGPL license.

Modified:

11 February 2013, J. Burkardt 15 July 2016 , I. Adams

Author:

Original FORTRAN77 version by Richard Brent. C++ version by John Burkardt.

Reference:

Richard Brent, Algorithms for Minimization Without Derivatives, Dover, 2002, ISBN: 0-486-41998-3, LC: QA402.5.B74.

Parameters:

Input, double A, B, the endpoints of the change of sign interval. Input, double T, a positive error tolerance. Output, double ZERO, the estimated value of a zero of the function F.

Parameters
[out]sbThe estimated value of a zero of the function F.
[in]aThe lower endpoint of the change of sign interval
[in]bThe upper endpoint of the change of sign interval
[in]tA positive error tolerance.
[in]rnA Random number.
[in]I1st Stokes element *** FIXMEDOC ***
[in]Q2nd Stokes element *** FIXMEDOC ***
[in]KIExtinction matrix element 0,0 *** FIXMEDOC ***
[in]KQExtinction matrix element 0,1 *** FIXMEDOC ***
Author
J. Burkhardt
Date
20??-??-??

Definition at line 132 of file montecarlo.cc.

◆ clear_rt_vars_at_gp()

void clear_rt_vars_at_gp ( Workspace ws,
MatrixView  ext_mat_mono,
VectorView  abs_vec_mono,
Numeric temperature,
const Agenda propmat_clearsky_agenda,
const Numeric f_mono,
const GridPos gp_p,
const GridPos gp_lat,
const GridPos gp_lon,
ConstVectorView  p_grid,
ConstTensor3View  t_field,
ConstTensor4View  vmr_field 
)

clear_rt_vars_at_gp.

Calculates a bunch of atmospheric variables at the end of a ppath.

Parameters
[in,out]wsCurrent workspace.
[out]ext_mat_monoTotal monochromatic extinction matrix.
[out]abs_vec_monoTotal monochromatic absorption vector.
[out]temperaturea vector of temperatures
[in]propmat_clearsky_agendaAs the WSA.
[in]f_monoFrequency (single entry vector).
[in]gp_pAn array of pressure gridpoints.
[in]gp_latAn array of latitude gridpoints.
[in]gp_lonAn array of longitude gridpoints.
[in]p_gridThe pressure grid.
[in]t_fieldThe temperature grid.
[in]vmr_fieldVMR field.
Author
Cory Davis
Date
2005-02-19? *** FIXMEDOC ***

Definition at line 244 of file montecarlo.cc.

References interp_atmfield_by_itw(), interp_atmfield_gp2itw(), interpweights(), itw2p(), joker, PropagationMatrix::MatrixAtPosition(), ConstTensor4View::nbooks(), ns, opt_prop_sum_propmat_clearsky(), propmat_clearsky_agendaExecute(), and StokesVector::VectorAtPosition().

Referenced by get_ppath_transmat(), and mcPathTraceRadar().

◆ cloud_atm_vars_by_gp()

void cloud_atm_vars_by_gp ( VectorView  pressure,
VectorView  temperature,
MatrixView  vmr,
MatrixView  pnd,
const ArrayOfGridPos gp_p,
const ArrayOfGridPos gp_lat,
const ArrayOfGridPos gp_lon,
const ArrayOfIndex cloudbox_limits,
ConstVectorView  p_grid_cloud,
ConstTensor3View  t_field_cloud,
ConstTensor4View  vmr_field_cloud,
ConstTensor4View  pnd_field 
)

cloud_atm_vars_by_gp.

Returns pressure, temperature, VMRs and PNDs, at points corresponding to arrays of gridpositions gp_p, gp_lat, and gp_lon. The field and grid input variables all span only the cloudbox

Parameters
[out]pressureA vector of pressures.
[out]temperatureA vector of temperatures.
[out]vmrA n_species by n_p matrix of VMRs.
[out]pndA n_scatelem by n_p matrix of PNDs.
[in]gp_pAn array of pressure gridpoints.
[in]gp_latAn array of latitude gridpoints.
[in]gp_lonAn array of longitude gridpoints.
[in]cloudbox_limitsThe limits of the cloud box.
[in]p_grid_cloudThe subset of p_grid within the cloudbox.
[in]t_field_cloudThe subset of t_field within the cloudbox.
[in]vmr_field_cloudThe subset of vmr_field within the cloudbox.
[in]pnd_fieldParticle number density field.
Author
Cory Davis
Date
2005-06-07

Definition at line 459 of file montecarlo.cc.

References gridpos_upperend_check(), i, interp_atmfield_by_itw(), interp_atmfield_gp2itw(), interpweights(), itw2p(), joker, ConstTensor4View::nbooks(), Array< base >::nelem(), ConstVectorView::nelem(), and ns.

Referenced by cloudy_rt_vars_at_gp(), and iwp_cloud_opt_pathCalc().

◆ cloudy_rt_vars_at_gp()

void cloudy_rt_vars_at_gp ( Workspace ws,
MatrixView  ext_mat_mono,
VectorView  abs_vec_mono,
VectorView  pnd_vec,
Numeric temperature,
const Agenda propmat_clearsky_agenda,
const Index  stokes_dim,
const Index  f_index,
const Vector f_grid,
const GridPos gp_p,
const GridPos gp_lat,
const GridPos gp_lon,
ConstVectorView  p_grid_cloud,
ConstTensor3View  t_field_cloud,
ConstTensor4View  vmr_field_cloud,
const Tensor4 pnd_field,
const ArrayOfArrayOfSingleScatteringData scat_data,
const ArrayOfIndex cloudbox_limits,
const Vector rte_los 
)

cloudy_rt_vars_at_gp.

Calculates a bunch of atmospheric variables at the end of a ppath.

Parameters
[in,out]wsCurrent workspace.
[out]ext_mat_monoTotal monochromatic extinction matrix.
[out]abs_vec_monoTotal monochromatic absorption vector.
[out]pnd_vecVector of particle number densities (one element per scattering element).
[out]temperatureA vector of temperatures
[in]propmat_clearsky_agendaAgenda calculating the absorption coefficient matrices.
[in]stokes_dimThe dimensionality of the Stokes vector (1-4).
[in]f_indexIndex of frequency grid point handeled.
[in]f_gridFrequency grid for monochromatic pencil beam calculations.
[in]gp_pAn array of pressure gridpoints.
[in]gp_latAn array of latitude gridpoints.
[in]gp_lonAn array of longitude gridpoints.
[in]p_grid_cloudThe subset of p_grid within the cloudbox.
[in]t_field_cloudThe subset of t_field within the cloudbox.
[in]vmr_field_cloudThe subset of vmr_field within the cloudbox.
[in]pnd_fieldParticle number density field.
[in]scat_dataArray of single scattering data.
[in]cloudbox_limitsThe limits of the cloud box.
[in]rte_los
Author
Cory Davis
Date
2005-02-19? *** FIXMEDOC ***

Definition at line 326 of file montecarlo.cc.

References cloud_atm_vars_by_gp(), joker, PropagationMatrix::MatrixAtPosition(), mirror_los(), ConstTensor4View::nbooks(), ns, opt_prop_Bulk(), opt_prop_NScatElems(), opt_prop_ScatSpecBulk(), opt_prop_sum_propmat_clearsky(), propmat_clearsky_agendaExecute(), and StokesVector::VectorAtPosition().

Referenced by get_ppath_transmat(), and mcPathTraceRadar().

◆ ext_I()

Numeric ext_I ( const Numeric I,
const Numeric Q,
const Numeric kI,
const Numeric kQ,
const Numeric s 
)

ext_I.

Calculate the extinction of I for a propagating "photon"

Parameters
[in]I1st Stokes element *** FIXMEDOC ***
[in]Q2nd Stokes element *** FIXMEDOC ***
[in]KIExtinction matrix element 0,0 *** FIXMEDOC ***
[in]KQExtinction matrix element 0,1 *** FIXMEDOC ***
[in]sPathlength *** FIXMEDOC ***
Author
Ian Adams
Date
2016-06-15

Definition at line 58 of file montecarlo.cc.

◆ get_ppath_transmat()

void get_ppath_transmat ( Workspace ws,
MatrixView trans_mat,
const Ppath ppath,
const Agenda propmat_clearsky_agenda,
const Index  stokes_dim,
const Index  f_index,
const Vector f_grid,
const Vector p_grid,
const Tensor3 t_field,
const Tensor4 vmr_field,
const ArrayOfIndex cloudbox_limits,
const Tensor4 pnd_field,
const ArrayOfArrayOfSingleScatteringData scat_data_mono,
const Verbosity verbosity 
)

get_ppath_transmat.

Routine to get the transmission matrix along a pre-defined propagation path. This is based on mcPathTraceGeneral using the routines from this source file. Routines from rte.cc require wind and magnetic field data that has not been typically passed to the Monte Carlo routines.

Parameters
[in,out]wsCurrent workspace.
[out]trans_matMatrix defining transmission over the ppath direction multiplied by sin(za)
[in]ppathPropagation path over which transmission matrix is desired
[in]propmat_clearsky_agendaAgenda calculating the absorption coefficient matrices.
[in]stokes_dimThe dimensionality of the Stokes vector (1-4).
[in]f_indexIndex of frequency grid point handeled.
[in]f_gridFrequency grid for monochromatic pencil beam calculations.
[in]p_gridThe pressure grid.
[in]t_fieldThe temperature grid.
[in]vmr_fieldVMR field.
[in]cloudbox_limitsThe limits of the cloud box.
[in]pnd_fieldParticle number density field.
[in]scat_dataArray of single scattering data.
[in]verbosityVerbosity variable to dynamically control the reporting level during runtime.
Author
Ian S. Adams
Date
2015-09-15

Definition at line 544 of file montecarlo.cc.

References clear_rt_vars_at_gp(), cloudy_rt_vars_at_gp(), CREATE_OUT0, ext2trans(), Ppath::gp_lat, Ppath::gp_lon, Ppath::gp_p, id_mat(), is_gp_inside_cloudbox(), joker, Ppath::los, Ppath::lstep, mult(), ConstTensor4View::nbooks(), and Ppath::np.

◆ is_anyptype_nonTotRan()

bool is_anyptype_nonTotRan ( const ArrayOfArrayOfSingleScatteringData scat_data_mono)

is_anyptype_nonTotRan.

Some operations in Monte Carlo simulations are different depending on the ptype of the scattering elements. This function searches scat_data to determine if any of the scattering elements have ptype=30.

Parameters
[in]scat_data_monoMonochromatic single scattering data.
Author
Cory Davis
Date
2004-1-31

Definition at line 694 of file montecarlo.cc.

References is_anyptype_nonTotRan(), Array< base >::nelem(), and PTYPE_TOTAL_RND.

Referenced by is_anyptype_nonTotRan(), and MCRadar().

◆ mcPathTraceGeneral()

void mcPathTraceGeneral ( Workspace ws,
MatrixView  evol_op,
Vector abs_vec_mono,
Numeric temperature,
MatrixView  ext_mat_mono,
Rng rng,
Vector rte_pos,
Vector rte_los,
Vector pnd_vec,
Numeric g,
Ppath ppath_step,
Index termination_flag,
bool &  inside_cloud,
const Agenda ppath_step_agenda,
const Numeric ppath_lmax,
const Numeric ppath_lraytrace,
const Numeric taustep_limit,
const Agenda propmat_clearsky_agenda,
const Index  stokes_dim,
const Index  f_index,
const Vector f_grid,
const Vector p_grid,
const Vector lat_grid,
const Vector lon_grid,
const Tensor3 z_field,
const Vector refellipsoid,
const Matrix z_surface,
const Tensor3 t_field,
const Tensor4 vmr_field,
const ArrayOfIndex cloudbox_limits,
const Tensor4 pnd_field,
const ArrayOfArrayOfSingleScatteringData scat_data,
const Verbosity verbosity 
)

mcPathTraceGeneral.

Performs the tasks of pathlength sampling.

Ray tracing done (but now only as far as determined by pathlength sampling) and calculation of the evolution operator and several atmospheric variables at the new point.

The end point of the ray tracing is returned by ppath_step, where the point of concern has index ppath_step.np-1. However, a somehwat dirty trick is used here to avoid copying of data. Only ppath.np is adjusted, and ppath_step can contain additional points (that should not be used).

Parameters
[in,out]wsCurrent workspace.
[out]evol_opEvolution operator (Stokes attenuation operator; exp(-tau)).
[out]abs_vec_monoTotal monochromatic absorption vector.
[out]temperatureA vector of temperatures.
[out]ext_mat_monoTotal monochromatic extinction matrix.
[out]rngRandom number generator instance.
[out]rte_posPosition for starting radiative transfer simulations.
[out]rte_losIncident line of sight for subsequent ray-tracing.
[out]grandomly chosen extinction path.
[out]ppath_stepA propagation path step.
[out]termination_flagFlag defining whether the path of the photon is terminated
[in]inside_cloudFlag defining inside or not the cloud box
[in]ppath_step_agendaCalculation of a propagation path step.
[in]ppath_lmaxMaximum length between points describing propagation paths.
[in]ppath_lraytraceMaximum length of ray tracing steps when determining propagation paths.
[in]taustep_limitDefines an upper step length in terms of optical thickness for Monte Carlo calculations.
[in]propmat_clearsky_agendaAgenda calculating the absorption coefficient matrices.
[in]stokes_dimThe dimensionality of the Stokes vector (1-4).
[in]f_indexIndex of frequency grid point handeled.
[in]f_gridFrequency grid for monochromatic pencil beam calculations.
[in]p_gridPressure grid.
[in]lat_gridLatitude grid.
[in]lon_gridLongitude grid.
[in]z_fieldThe field of geometrical altitudes.
[in]refellipsoidReference ellipsoid.
[in]z_surfaceThe surface altitude.
[in]t_fieldThe temperature grid.
[in]vmr_fieldVMR field.
[in]cloudbox_limitsThe limits of the cloud box.
[in]pnd_fieldParticle number density field.
[in]scat_dataArray of single scattering data.
[in]verbosityVerbosity variable to dynamically control the reporting level during runtime.
Author
Cory Davis
Date
2005-2-21

Definition at line 711 of file montecarlo.cc.

◆ mcPathTraceRadar()

void mcPathTraceRadar ( Workspace ws,
MatrixView  evol_op,
Vector abs_vec_mono,
Numeric temperature,
MatrixView  ext_mat_mono,
Rng rng,
Vector rte_pos,
Vector rte_los,
Vector pnd_vec,
Numeric stot,
Numeric ttot,
Ppath ppath_step,
Index termination_flag,
bool &  inside_cloud,
const Agenda ppath_step_agenda,
const Numeric ppath_lmax,
const Numeric ppath_lraytrace,
const Agenda propmat_clearsky_agenda,
const bool &  anyptype_nonTotRan,
const Index  stokes_dim,
const Index  f_index,
const Vector f_grid,
const Vector Iprop,
const Vector p_grid,
const Vector lat_grid,
const Vector lon_grid,
const Tensor3 z_field,
const Vector refellipsoid,
const Matrix z_surface,
const Tensor3 t_field,
const Tensor4 vmr_field,
const ArrayOfIndex cloudbox_limits,
const Tensor4 pnd_field,
const ArrayOfArrayOfSingleScatteringData scat_data_mono,
const Verbosity verbosity 
)

mcPathTraceRadar.

Performs the tasks of pathlength sampling.

Ray tracing done (but now only as far as determined by pathlength sampling) and calculation of the evolution operator and several atmospheric variables at the new point.

The end point of the ray tracing is returned by ppath_step, where the point of concern has index ppath_step.np-1. However, a somehwat dirty trick is used here to avoid copying of data. Only ppath.np is adjusted, and ppath_step can contain additional points (that should not be used).

Parameters
[in,out]wsCurrent workspace.
[out]evol_opEvolution operator (Stokes attenuation operator; exp(-tau)).
[out]abs_vec_monoTotal monochromatic absorption vector.
[out]temperatureA vector of temperatures.
[out]ext_mat_monoTotal monochromatic extinction matrix.
[out]rngRandom number generator instance.
[out]rte_posPosition for starting radiative transfer simulations.
[out]rte_losIncident line of sight for subsequent ray-tracing.
[out]stotStarting point of total path lenght
[out]ttotTotal path length *** FIXMEDOC ***
[out]ppath_stepA propagation path step.
[out]termination_flagFlag defining whether the path of the photon is terminated
[in]inside_cloudFlag defining inside or not the cloud box
[in]ppath_step_agendaCalculation of a propagation path step.
[in]ppath_lmaxMaximum length between points describing propagation paths.
[in]ppath_lraytraceMaximum length of ray tracing steps when determining propagation paths.
[in]taustep_limitDefines an upper step length in terms of optical thickness for Monte Carlo calculations.
[in]propmat_clearsky_agendaAgenda calculating the absorption coefficient matrices.
[in]anyptype_nonTotRanFlag definining any particle type, but for totally random oriented.
[in]stokes_dimThe dimensionality of the Stokes vector (1-4).
[in]f_indexIndex of frequency grid point handeled.
[in]f_gridFrequency grid for monochromatic pencil beam calculations.
[in]IpropIncident I component of the Stokes vector *** FIXMEDOC ***
[in]p_gridPressure grid.
[in]lat_gridLatitude grid.
[in]lon_gridLongitude grid.
[in]z_fieldThe field of geometrical altitudes.
[in]refellipsoidReference ellipsoid.
[in]z_surfaceThe surface altitude.
[in]t_fieldThe temperature grid.
[in]vmr_fieldVMR field.
[in]cloudbox_limitsThe limits of the cloud box.
[in]pnd_fieldParticle number density field.
[in]scat_dataArray of single scattering data.
[in]verbosityVerbosity variable to dynamically control the reporting level during runtime.
Author
Cory Davis (mcPathTraceGeneral), Ian S. Adams
Date
2015-09-08

Definition at line 1054 of file montecarlo.cc.

References clear_rt_vars_at_gp(), cloudy_rt_vars_at_gp(), CREATE_OUT0, Rng::draw(), Ppath::end_lstep, ext2trans(), fractional_gp(), Ppath::gp_lat, Ppath::gp_lon, Ppath::gp_p, id_mat(), is_gp_inside_cloudbox(), joker, Ppath::los, Ppath::lstep, mult(), ConstVectorView::nelem(), Ppath::ngroup, Ppath::np, ppath_start_stepping(), ppath_step_agendaExecute(), ppath_what_background(), r, and SPEED_OF_LIGHT.

Referenced by MCRadar().

◆ Sample_los()

void Sample_los ( VectorView  new_rte_los,
Numeric g_los_csc_theta,
MatrixView  Z,
Rng rng,
ConstVectorView  rte_los,
const ArrayOfArrayOfSingleScatteringData scat_data,
const Index  stokes_dim,
const Index  f_index,
ConstVectorView  pnd_vec,
ConstVectorView  Z11maxvector,
const Numeric  Csca,
const Numeric  rtp_temperature,
const Index  t_interp_order = 1 
)

Sample_los.

Sampling the new incident direction according to a rejection method. Calculation of the bulk scattering phase matrix.

Parameters
[out]new_rte_losIncident line of sight for subsequent.
[out]g_los_csc_thetaProbability density for the chosen direction multiplied by sin(za)
[out]ZBulk phase matrix in Stokes notation.
[in,out]rngRng random number generator instance.
[in]rte_losIncident line of sight for subsequent ray-tracing.
[in]scat_dataAs the WSV.
[in]stokes_dimAs the WSV.
[out]pnd_vecVector of particle number densities (one element per scattering element).
[in]Z11maxvectorVector holding the maximum phase function for each scattering element.
[in]CscaScattering cross section
[in]rtp_temperatureAs the WSV.
Author
Cory Davis
Date
2003-06-19

Definition at line 1391 of file montecarlo.cc.

References Rng::draw(), i, joker, mirror_los(), ConstVectorView::nelem(), pha_mat_Bulk(), pha_mat_NScatElems(), pha_mat_ScatSpecBulk(), RAD2DEG, and TotalNumberOfElements().

Referenced by MCGeneral(), and MCIPA().

◆ Sample_los_uniform()

void Sample_los_uniform ( VectorView  new_rte_los,
Rng rng 
)

Sample_los_uniform.

Sampling the new direction uniformly

Parameters
[out]new_rte_losIncident line of sight for subsequent.
[out]rngRng random number generator instance.
Author
*** FIXMEDOC ***
Date
*** FIXMEDOC ***

Definition at line 1471 of file montecarlo.cc.

References Rng::draw(), and RAD2DEG.

Variable Documentation

◆ SPEED_OF_LIGHT

const Numeric SPEED_OF_LIGHT

Referenced by mcPathTraceRadar().