ARTS  2.4.0(git:4fb77825)
m_fluxes.cc
Go to the documentation of this file.
1 /* Copyright (C) 2018
2  Manfred Brath <manfred.brath@uni-hamburg.de>
3 
4  This program is free software; you can redistribute it and/or modify it
5  under the terms of the GNU General Public License as published by the
6  Free Software Foundation; either version 2, or (at your option) any
7  later version.
8 
9  This program is distributed in the hope that it will be useful,
10  but WITHOUT ANY WARRANTY; without even the implied warranty of
11  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  GNU General Public License for more details.
13 
14  You should have received a copy of the GNU General Public License
15  along with this program; if not, write to the Free Software
16  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307,
17  USA. */
18 
19 /*===========================================================================
20  === File description
21  ===========================================================================*/
22 #include <iostream>
23 #include <stdexcept>
24 #include "absorption.h"
25 #include "agenda_class.h"
26 #include "auto_md.h"
27 #include "check_input.h"
28 #include "legendre.h"
29 #include "math_funcs.h"
30 #include "matpackVII.h"
31 #include "messages.h"
32 #include "sorting.h"
33 #include "surface.h"
34 #include "workspace_ng.h"
35 
47 extern const Numeric PI;
48 extern const Numeric DEG2RAD;
49 
50 /*===========================================================================
51  === The functions
52  ===========================================================================*/
53 
54 /* Workspace method: Doxygen documentation will be auto-generated */
56  Vector& aa_grid,
58  // Keywords:
59  const Index& N_za_grid,
60  const Index& N_aa_grid,
61  const String& za_grid_type,
62  const Verbosity&) {
63  // Azimuth angle grid
64  if (N_aa_grid > 1)
65  nlinspace(aa_grid, 0, 360, N_aa_grid);
66  else if (N_aa_grid < 1) {
67  ostringstream os;
68  os << "N_aa_grid must be > 0 (even for 1D).";
69  throw std::runtime_error(os.str());
70  } else {
71  aa_grid.resize(1);
72  aa_grid[0] = 0.;
73  }
74 
75  if (N_za_grid % 2 == 1) {
76  ostringstream os;
77  os << "N_za_grid must be even.";
78  throw runtime_error(os.str());
79  }
80 
81  Index nph = N_za_grid / 2;
82 
83  //calculate zenith angle grid
84  za_grid.resize(N_za_grid);
85  za_grid = 0.;
86  za_grid_weights.resize(N_za_grid);
87  za_grid_weights = 0;
88 
89  if (za_grid_type == "double_gauss") {
90  Vector x;
91  Vector w;
92  Vector xtemp;
93  Vector wtemp;
94  //Numeric theta;
95 
96  //calculate legendre weights and evaluation points
97  gsl_integration_glfixed_table_alloc(xtemp, wtemp, nph);
98 
99  x.resize(nph);
100  w.resize(nph);
101 
102  // reorder and expand weights and abscissa vectors
103  // transform abscissa vector from cos(theta)-space to theta-space
104  // and adjust the domain and weights
105  if (nph % 2 == 1) {
106  x[xtemp.nelem() - 1] = acos((xtemp[0] + 1) / 2) / DEG2RAD;
107  w[wtemp.nelem() - 1] = wtemp[0] / 2;
108 
109  for (Index i = 0; i < xtemp.nelem() - 1; i++) {
110  x[i] = acos((xtemp[xtemp.nelem() - 1 - i] + 1) / 2.) / DEG2RAD;
111  x[xtemp.nelem() + i] = acos(1 - (xtemp[i + 1] + 1) / 2.) / DEG2RAD;
112 
113  w[i] = wtemp[wtemp.nelem() - 1 - i] / 2;
114  w[wtemp.nelem() + i] = wtemp[i + 1] / 2;
115  }
116  } else {
117  for (Index i = 0; i < xtemp.nelem(); i++) {
118  x[i] = acos((xtemp[xtemp.nelem() - 1 - i] + 1) / 2.) / DEG2RAD;
119  x[xtemp.nelem() + i] = acos(1 - (xtemp[i] + 1) / 2.) / DEG2RAD;
120 
121  w[i] = wtemp[wtemp.nelem() - 1 - i] / 2;
122  w[wtemp.nelem() + i] = wtemp[i] / 2;
123  }
124  }
125 
126  for (Index i = 0; i < nph; i++) {
127  //set the angles
128  //theta=x[i];//acos((x[i]+1)/2)/DEG2RAD;
129  za_grid[i] = x[i];
130  za_grid[za_grid.nelem() - 1 - i] = 180 - x[i];
131 
132  // set the weights to the right component
133  za_grid_weights[i] = w[i];
134  za_grid_weights[za_grid_weights.nelem() - 1 - i] = w[i];
135  }
136 
137  } else if (za_grid_type == "linear") {
138  Vector x;
139  Vector w;
141 
142  for (Index i = 0; i < N_za_grid; i++) {
143  za_grid[i] = (x[i] + 1) * 90.;
144 
145  // set the weights to the right component
146  // by adjusting the domain, we also have to adjust the weights
147  za_grid_weights[i] = w[i] * sin(za_grid[i] * DEG2RAD);
148  }
149  } else if (za_grid_type == "linear_mu") {
150  Vector x;
151  Vector w;
152 
153  //calculate weights in cos(theta) space
155 
156  //allocate
157  Vector za_grid_temp;
158  za_grid_temp.resize(x.nelem());
159 
160  for (Index i = 0; i < N_za_grid; i++) {
161  za_grid_temp[i] = acos(x[i]) / DEG2RAD;
162  }
163 
164  //#sort weights and theta in increasing direction of za_grid
165  za_grid = za_grid_temp[Range(x.nelem() - 1, x.nelem(), -1)];
166  za_grid_weights = w[Range(x.nelem() - 1, x.nelem(), -1)];
167 
168  } else {
169  ostringstream os;
170  os << "The selected grid type is not implemented";
171  throw std::runtime_error(os.str());
172  }
173 
174  //be sure that the first and the last angle are within the closed interval
175  //between 0 and 180 deg, because ARTS is picky if the angles are due to numerics
176  // slightly below and above,respectively.
177  if (za_grid[0] < 0) {
178  za_grid[0] = 0.;
179  }
180 
181  if (za_grid[za_grid.nelem() - 1] > 180) {
182  za_grid[za_grid.nelem() - 1] = 180.;
183  }
184 }
185 
186 /* Workspace method: Doxygen documentation will be auto-generated */
188  const Vector& p_grid,
189  const Tensor4& irradiance_field,
191  const Numeric& g0,
192  const Verbosity&) {
193  //allocate
194  heating_rates.resize(irradiance_field.nbooks(),
195  irradiance_field.npages(),
196  irradiance_field.nrows());
197  heating_rates = 0;
198 
199  // allocate some auxiliary variables
200  Numeric net_flux_b; //net_flux bottom
201  Numeric net_flux_c; //net_flux center
202  Numeric net_flux_t; //net_flux top
203  Index idx;
204 
205  // calculate heating rates, we skip the upper and lower boundary here, because to achieve the same
206  // second order accuracy for the derivation of the net flux at the edged, we use
207  // a differentiation based on polynomial interpolation
208  for (Index b = 1; b < irradiance_field.nbooks() - 1; b++) {
209  for (Index p = 0; p < irradiance_field.npages(); p++) {
210  for (Index r = 0; r < irradiance_field.nrows(); r++) {
211  net_flux_b = (irradiance_field(b - 1, p, r, 0) +
212  irradiance_field(b - 1, p, r, 1));
213  net_flux_t = (irradiance_field(b + 1, p, r, 0) +
214  irradiance_field(b + 1, p, r, 1));
215 
216  heating_rates(b, p, r) = (net_flux_t - net_flux_b) /
217  (p_grid[b + 1] - p_grid[b - 1]) * g0 /
218  specific_heat_capacity(b, p, r);
219  }
220  }
221  }
222 
223  idx = irradiance_field.nbooks();
224 
225  // now calculate the heating rates for the upper and lower boundary
226  for (Index p = 0; p < irradiance_field.npages(); p++) {
227  for (Index r = 0; r < irradiance_field.nrows(); r++) {
228  // lower boundary
229  net_flux_b =
230  (irradiance_field(0, p, r, 0) + irradiance_field(0, p, r, 1));
231  net_flux_c =
232  (irradiance_field(1, p, r, 0) + irradiance_field(1, p, r, 1));
233  net_flux_t =
234  (irradiance_field(2, p, r, 0) + irradiance_field(0, p, r, 1));
235 
236  heating_rates(0, p, r) = (-3 * net_flux_b + 4 * net_flux_c - net_flux_t) /
237  (p_grid[2] - p_grid[0]) * g0 /
238  specific_heat_capacity(0, p, r);
239 
240  // lower boundary
241  net_flux_t = (irradiance_field(idx - 1, p, r, 0) +
242  irradiance_field(idx - 1, p, r, 1));
243  net_flux_c = (irradiance_field(idx - 2, p, r, 0) +
244  irradiance_field(idx - 2, p, r, 1));
245  net_flux_b = (irradiance_field(idx - 3, p, r, 0) +
246  irradiance_field(idx - 3, p, r, 1));
247 
248  heating_rates(idx - 1, p, r) =
249  -(-3 * net_flux_t + 4 * net_flux_c - net_flux_b) /
250  (p_grid[2] - p_grid[0]) * g0 / specific_heat_capacity(0, p, r);
251  }
252  }
253 }
254 
255 /* Workspace method: Doxygen documentation will be auto-generated */
257  const Tensor5& radiance_field,
258  const Vector& za_grid,
259  const Vector& aa_grid,
260  const Vector& za_grid_weights,
261  const Verbosity&) {
262  // Number of zenith angles.
263  const Index N_scat_za = za_grid.nelem();
264  const Index N_scat_aa = aa_grid.nelem();
265 
266  Tensor4 radiance_field_aa_integrated;
267 
268  //azimuth integration
269  if (N_scat_aa == 1) //1D case no azimuth dependency
270  {
271  radiance_field_aa_integrated =
273  radiance_field_aa_integrated *= 2 * PI;
274 
275  } else //general case with azimuth dependency
276  {
277  radiance_field_aa_integrated.resize(radiance_field.nshelves(),
278  radiance_field.nbooks(),
279  radiance_field.npages(),
280  radiance_field.nrows());
281  radiance_field_aa_integrated = 0.;
282 
283  for (Index b = 0; b < radiance_field_aa_integrated.nbooks(); b++) {
284  for (Index p = 0; p < radiance_field_aa_integrated.npages(); p++) {
285  for (Index r = 0; r < radiance_field_aa_integrated.nrows(); r++) {
286  for (Index c = 0; c < radiance_field_aa_integrated.ncols(); c++) {
287  for (Index i = 0; i < N_scat_aa - 1; i++) {
288  radiance_field_aa_integrated(b, p, r, c) +=
289  (radiance_field(b, p, r, c, i) +
290  radiance_field(b, p, r, c, i + 1)) /
291  2. * abs(aa_grid[i + 1] - aa_grid[i]) * DEG2RAD;
292  }
293  }
294  }
295  }
296  }
297  }
298 
299  //allocate
300  irradiance_field.resize(radiance_field.nshelves(),
301  radiance_field.nbooks(),
302  radiance_field.npages(),
303  2);
304  irradiance_field = 0;
305 
306  // zenith angle integration
307 
308  for (Index b = 0; b < irradiance_field.nbooks(); b++) {
309  for (Index p = 0; p < irradiance_field.npages(); p++) {
310  for (Index r = 0; r < irradiance_field.nrows(); r++) {
311  for (Index i = 0; i < N_scat_za; i++) {
312  if (za_grid[i] <= 90.) {
313  irradiance_field(b, p, r, 0) +=
314  radiance_field_aa_integrated(b, p, r, i) *
315  cos(za_grid[i] * DEG2RAD) * (-1.) * za_grid_weights[i];
316  } else {
317  irradiance_field(b, p, r, 1) +=
318  radiance_field_aa_integrated(b, p, r, i) *
319  cos(za_grid[i] * DEG2RAD) * (-1.) * za_grid_weights[i];
320  }
321  }
322  }
323  }
324  }
325 }
326 
327 /* Workspace method: Doxygen documentation will be auto-generated */
329  const Vector& f_grid,
330  const Tensor5& spectral_radiation_field,
331  const Verbosity&) {
332  if (f_grid.nelem() != spectral_radiation_field.nshelves()) {
333  throw runtime_error(
334  "The length of f_grid does not match with\n"
335  " the first dimension of the spectral_radiation_field");
336  }
337 
338  //allocate
339  radiation_field.resize(spectral_radiation_field.nbooks(),
340  spectral_radiation_field.npages(),
341  spectral_radiation_field.nrows(),
342  spectral_radiation_field.ncols());
343  radiation_field = 0;
344 
345  // frequency integration
346  for (Index i = 0; i < spectral_radiation_field.nshelves() - 1; i++) {
347  const Numeric df = f_grid[i + 1] - f_grid[i];
348 
349  for (Index b = 0; b < radiation_field.nbooks(); b++) {
350  for (Index p = 0; p < radiation_field.npages(); p++) {
351  for (Index r = 0; r < radiation_field.nrows(); r++) {
352  for (Index c = 0; c < radiation_field.ncols(); c++) {
353  radiation_field(b, p, r, c) +=
354  (spectral_radiation_field(i + 1, b, p, r, c) +
355  spectral_radiation_field(i, b, p, r, c)) /
356  2 * df;
357  }
358  }
359  }
360  }
361  }
362 }
363 
364 /* Workspace method: Doxygen documentation will be auto-generated */
366  const Vector& f_grid,
367  const Tensor7& spectral_radiation_field,
368  const Verbosity&) {
369  if (f_grid.nelem() != spectral_radiation_field.nlibraries()) {
370  throw runtime_error(
371  "The length of f_grid does not match with\n"
372  " the first dimension of the spectral_radiation_field");
373  }
374 
375  //allocate
376  radiation_field.resize(spectral_radiation_field.nvitrines(),
377  spectral_radiation_field.nshelves(),
378  spectral_radiation_field.nbooks(),
379  spectral_radiation_field.npages(),
380  spectral_radiation_field.nrows());
381  radiation_field = 0;
382 
383  // frequency integration
384  for (Index i = 0; i < spectral_radiation_field.nlibraries() - 1; i++) {
385  const Numeric df = f_grid[i + 1] - f_grid[i];
386 
387  for (Index s = 0; s < radiation_field.nshelves(); s++) {
388  for (Index b = 0; b < radiation_field.nbooks(); b++) {
389  for (Index p = 0; p < radiation_field.npages(); p++) {
390  for (Index r = 0; r < radiation_field.nrows(); r++) {
391  for (Index c = 0; c < radiation_field.ncols(); c++) {
392  radiation_field(s, b, p, r, c) +=
393  (spectral_radiation_field(i + 1, s, b, p, r, c, 0) +
394  spectral_radiation_field(i, s, b, p, r, c, 0)) /
395  2 * df;
396  }
397  }
398  }
399  }
400  }
401  }
402 }
403 
404 /* Workspace method: Doxygen documentation will be auto-generated */
408  const Vector& za_grid,
409  const Vector& aa_grid,
410  const Vector& za_grid_weights,
411  const Verbosity&) {
412  // Number of zenith angles.
413  const Index N_scat_za = za_grid.nelem();
414  const Index N_scat_aa = aa_grid.nelem();
415 
416  Tensor5 iy_field_aa_integrated;
417 
418  //azimuth integration
419  if (N_scat_aa == 1) //1D case no azimuth dependency
420  {
421  iy_field_aa_integrated =
423  iy_field_aa_integrated *= 2 * PI;
424 
425  } else //general case with azimuth dependency
426  {
427  iy_field_aa_integrated.resize(spectral_radiance_field.nlibraries(),
428  spectral_radiance_field.nvitrines(),
429  spectral_radiance_field.nshelves(),
430  spectral_radiance_field.nbooks(),
431  spectral_radiance_field.npages());
432  iy_field_aa_integrated = 0.;
433 
434  for (Index s = 0; s < iy_field_aa_integrated.nshelves(); s++) {
435  for (Index b = 0; b < iy_field_aa_integrated.nbooks(); b++) {
436  for (Index p = 0; p < iy_field_aa_integrated.npages(); p++) {
437  for (Index r = 0; r < iy_field_aa_integrated.nrows(); r++) {
438  for (Index c = 0; c < iy_field_aa_integrated.ncols(); c++) {
439  for (Index i = 0; i < N_scat_aa - 1; i++) {
440  iy_field_aa_integrated(s, b, p, r, c) +=
441  (spectral_radiance_field(s, b, p, r, c, i, 0) +
442  spectral_radiance_field(s, b, p, r, c, i + 1, 0)) /
443  2. * abs(aa_grid[i + 1] - aa_grid[i]) * DEG2RAD;
444  }
445  }
446  }
447  }
448  }
449  }
450  }
451 
452  //allocate
454  spectral_radiance_field.nvitrines(),
455  spectral_radiance_field.nshelves(),
456  spectral_radiance_field.nbooks(),
457  2);
459 
460  // zenith angle integration
461  for (Index s = 0; s < spectral_irradiance_field.nshelves(); s++) {
462  for (Index b = 0; b < spectral_irradiance_field.nbooks(); b++) {
463  for (Index p = 0; p < spectral_irradiance_field.npages(); p++) {
464  for (Index r = 0; r < spectral_irradiance_field.nrows(); r++) {
465  for (Index i = 0; i < N_scat_za; i++) {
466  if (za_grid[i] <= 90.) {
467  spectral_irradiance_field(s, b, p, r, 0) +=
468  iy_field_aa_integrated(s, b, p, r, i) *
469  cos(za_grid[i] * DEG2RAD) * (-1.) * za_grid_weights[i];
470  } else {
471  spectral_irradiance_field(s, b, p, r, 1) +=
472  iy_field_aa_integrated(s, b, p, r, i) *
473  cos(za_grid[i] * DEG2RAD) * (-1.) * za_grid_weights[i];
474  }
475  }
476  }
477  }
478  }
479  }
480 }
481 
482 /* Workspace method: Doxygen documentation will be auto-generated */
484  Workspace& ws,
486  Tensor3& trans_field,
488  const Agenda& water_p_eq_agenda,
489  const Agenda& iy_space_agenda,
490  const Agenda& iy_surface_agenda,
491  const Agenda& iy_cloudbox_agenda,
492  const Index& stokes_dim,
493  const Vector& f_grid,
494  const Index& atmosphere_dim,
495  const Vector& p_grid,
496  const Tensor3& z_field,
497  const Tensor3& t_field,
498  const EnergyLevelMap& nlte_field,
499  const Tensor4& vmr_field,
501  const Tensor3& wind_u_field,
502  const Tensor3& wind_v_field,
503  const Tensor3& wind_w_field,
504  const Tensor3& mag_u_field,
505  const Tensor3& mag_v_field,
506  const Tensor3& mag_w_field,
507  const Matrix& z_surface,
508  const Numeric& ppath_lmax,
509  const Numeric& rte_alonglos_v,
511  const Vector& za_grid,
512  const Index& use_parallel_iy,
513  const Verbosity& verbosity) {
514  // Check input
515  if (atmosphere_dim != 1)
516  throw runtime_error("This method only works for atmosphere_dim = 1.");
517 
518  // Sizes
519  const Index nl = p_grid.nelem();
520  const Index nf = f_grid.nelem();
521  const Index nza = za_grid.nelem();
522 
523  // Init spectral_radiance_field and trans_field
524  spectral_radiance_field.resize(nf, nl, 1, 1, nza, 1, stokes_dim);
525  trans_field.resize(nf, nl, nza);
526 
527  // De-activate cloudbox
529  const ArrayOfIndex cloudbox_limits(0);
530 
531  // Various variables
532  const String iy_unit = "1";
533  const ArrayOfString iy_aux_vars(0);
534  const Vector rte_pos2(0);
535  const Index iy_agenda_call1 = 1;
536  const Tensor3 iy_transmission(0, 0, 0);
537  const Index jacobian_do = 0;
539  // Create one altitude just above TOA
540  const Numeric z_space = z_field(nl - 1, 0, 0) + 10;
541 
542  Workspace l_ws(ws);
543  ArrayOfString fail_msg;
544  bool failed = false;
545 
546  // Define iy_main_agenda to be consistent with the assumptions of
547  // this method. This definition of iy_main_agenda will be used to when
548  // calculating the the radiation reflected by the surface
550  iy_main_agenda.append("ppathPlaneParallel", TokVal());
551  iy_main_agenda.append("iyEmissionStandard", TokVal());
552  iy_main_agenda.set_name("iy_main_agenda");
553  iy_main_agenda.check(ws, verbosity);
554 
555  // Index in p_grid where field at surface shall be placed
556  const Index i0 = index_of_zsurface(z_surface(0, 0), z_field(joker, 0, 0));
557 
558  // Loop zenith angles
559  //
560  if (nza)
561 #pragma omp parallel for if (!arts_omp_in_parallel() && nza > 1 && \
562  !use_parallel_iy) firstprivate(l_ws)
563  for (Index i = 0; i < nza; i++) {
564  if (failed) continue;
565  try {
566  // Define output variables
567  Ppath ppath;
575 
576  Index iy_id = i;
577  Vector rte_los(1, za_grid[i]);
578  Vector rte_pos(1, za_grid[i] < 90 ? z_surface(0, 0) : z_space);
579 
582  z_field,
583  z_surface,
584  cloudbox_on,
587  rte_pos,
588  rte_los,
589  ppath_lmax,
590  verbosity);
591  assert(ppath.gp_p[ppath.np - 1].idx == i0 ||
592  ppath.gp_p[ppath.np - 1].idx == nl - 2);
593 
594  if (use_parallel_iy) {
595  iyEmissionStandard(l_ws,
596  iy,
597  iy_aux,
598  diy_dx,
599  ppvar_p,
600  ppvar_t,
601  ppvar_nlte,
602  ppvar_vmr,
603  ppvar_wind,
604  ppvar_mag,
605  ppvar_f,
606  ppvar_iy,
609  iy_id,
610  stokes_dim,
611  f_grid,
613  p_grid,
614  t_field,
615  nlte_field,
616  vmr_field,
617  abs_species,
618  wind_u_field,
619  wind_v_field,
620  wind_w_field,
621  mag_u_field,
622  mag_v_field,
623  mag_w_field,
624  cloudbox_on,
625  iy_unit,
626  iy_aux_vars,
627  jacobian_do,
629  ppath,
630  rte_pos2,
641  verbosity);
642  } else {
644  iy,
645  iy_aux,
646  diy_dx,
647  ppvar_p,
648  ppvar_t,
649  ppvar_nlte,
650  ppvar_vmr,
651  ppvar_wind,
652  ppvar_mag,
653  ppvar_f,
654  ppvar_iy,
657  iy_id,
658  stokes_dim,
659  f_grid,
661  p_grid,
662  t_field,
663  nlte_field,
664  vmr_field,
665  abs_species,
666  wind_u_field,
667  wind_v_field,
668  wind_w_field,
669  mag_u_field,
670  mag_v_field,
671  mag_w_field,
672  cloudbox_on,
673  iy_unit,
674  iy_aux_vars,
675  jacobian_do,
677  ppath,
678  rte_pos2,
689  verbosity);
690  }
691  assert(iy.ncols() == stokes_dim);
692 
693  // First and last points are most easily handled separately
694  if (za_grid[i] < 90) {
695  spectral_radiance_field(joker, i0, 0, 0, i, 0, joker) =
696  ppvar_iy(joker, joker, 0);
697  spectral_radiance_field(joker, nl - 1, 0, 0, i, 0, joker) =
698  ppvar_iy(joker, joker, ppath.np - 1);
699  trans_field(joker, 0, i) = ppvar_trans_partial(0, joker, 0, 0);
700  trans_field(joker, nl - 1, i) =
701  ppvar_trans_partial(ppath.np - 1, joker, 0, 0);
702  } else {
703  spectral_radiance_field(joker, nl - 1, 0, 0, i, 0, joker) =
704  ppvar_iy(joker, joker, 0);
705  spectral_radiance_field(joker, i0, 0, 0, i, 0, joker) =
706  ppvar_iy(joker, joker, ppath.np - 1);
707  trans_field(joker, nl - 1, i) = ppvar_trans_partial(0, joker, 0, 0);
708  trans_field(joker, 0, i) =
709  ppvar_trans_partial(ppath.np - 1, joker, 0, 0);
710  }
711 
712  // Remaining points
713  for (Index p = 1; p < ppath.np - 1; p++) {
714  // We just store values at pressure levels
715  if (ppath.gp_p[p].fd[0] < 1e-2) {
717  joker, ppath.gp_p[p].idx, 0, 0, i, 0, joker) =
718  ppvar_iy(joker, joker, p);
719  trans_field(joker, ppath.gp_p[p].idx, i) =
720  ppvar_trans_partial(p, joker, 0, 0);
721  }
722  }
723 
724  // We don't want undefined values to possibly affect an interpolation,
725  // and to be safe we set the fields for underground levels to equal the
726  // ones at the surface
727  for (Index p = 0; p < i0; p++) {
728  spectral_radiance_field(joker, p, 0, 0, i, 0, joker) =
729  spectral_radiance_field(joker, i0, 0, 0, i, 0, joker);
730  trans_field(joker, p, i) = trans_field(joker, i0, i);
731  }
732 
733  } catch (const std::exception& e) {
734 #pragma omp critical(planep_setabort)
735  failed = true;
736 
737  ostringstream os;
738  os << "Run-time error at nza #" << i << ": \n" << e.what();
739 #pragma omp critical(planep_push_fail_msg)
740  fail_msg.push_back(os.str());
741  }
742  }
743 
744  if (fail_msg.nelem()) {
745  ostringstream os;
746  for (auto& msg : fail_msg) os << msg << '\n';
747  throw runtime_error(os.str());
748  }
749 }
750 
751 /* Workspace method: Doxygen documentation will be auto-generated */
754  const Index& atmosphere_dim,
755  const Vector& p_grid,
756  const Index& cloudbox_on,
758  const Tensor7& cloudbox_field,
759  const Verbosity&) {
760  if (atmosphere_dim > 1)
761  throw runtime_error("This method can only be used for 1D calculations.\n");
762  if (!cloudbox_on)
763  throw runtime_error("Cloudbox is off. This is not handled by this method.");
764  if (cloudbox_limits[0] || cloudbox_limits[1] != p_grid.nelem() - 1)
765  throw runtime_error(
766  "The cloudbox must cover all pressure levels "
767  "to use this method.");
768 
769  // If all OK, it is just to copy
771 }
772 
773 /* Workspace method: Doxygen documentation will be auto-generated */
775  Workspace& ws,
778  const Agenda& water_p_eq_agenda,
779  const Agenda& iy_space_agenda,
780  const Agenda& iy_surface_agenda,
781  const Agenda& iy_cloudbox_agenda,
782  const Index& stokes_dim,
783  const Vector& f_grid,
784  const Index& atmosphere_dim,
785  const Vector& p_grid,
786  const Tensor3& z_field,
787  const Tensor3& t_field,
788  const EnergyLevelMap& nlte_field,
789  const Tensor4& vmr_field,
791  const Tensor3& wind_u_field,
792  const Tensor3& wind_v_field,
793  const Tensor3& wind_w_field,
794  const Tensor3& mag_u_field,
795  const Tensor3& mag_v_field,
796  const Tensor3& mag_w_field,
797  const Matrix& z_surface,
798  const Index& cloudbox_on,
800  const Tensor7& cloudbox_field,
801  const Numeric& ppath_lmax,
802  const Numeric& rte_alonglos_v,
804  const Vector& za_grid,
805  const Index& use_parallel_iy,
806  const Verbosity& verbosity) {
807  // Check input
808  if (atmosphere_dim != 1)
809  throw runtime_error("This method only works for atmosphere_dim = 1.");
810  if (!cloudbox_on)
811  throw runtime_error("No ned to use this method with cloudbox=0.");
812  if (cloudbox_limits[0])
813  throw runtime_error(
814  "The first element of *cloudbox_limits* must be zero "
815  "to use this method.");
816 
817  // Sizes
818  const Index nl = p_grid.nelem();
819  const Index nf = f_grid.nelem();
820  const Index nza = za_grid.nelem();
821 
822  // Init spectral_radiance_field
823  spectral_radiance_field.resize(nf, nl, 1, 1, nza, 1, stokes_dim);
824 
825  // and copy the part taken from cloudbox_field
827  Range(0, cloudbox_limits[1] + 1),
828  joker,
829  joker,
830  joker,
831  joker,
833 
834  // Various variables
836  const String iy_unit = "1";
837  const ArrayOfString iy_aux_vars(0);
838  const Vector rte_pos2(0);
839  const Index iy_agenda_call1 = 1;
840  const Tensor3 iy_transmission(0, 0, 0);
841  const Index jacobian_do = 0;
843  // Create one altitude just above TOA
844  const Numeric z_space = z_field(nl - 1, 0, 0) + 10;
845 
846  Workspace l_ws(ws);
847  ArrayOfString fail_msg;
848  bool failed = false;
849 
850  // Define iy_main_agenda to be consistent with the assumptions of
851  // this method (but the agenda will not be used).
853  iy_main_agenda.append("ppathPlaneParallel", TokVal());
854  iy_main_agenda.append("iyEmissionStandard", TokVal());
855  iy_main_agenda.set_name("iy_main_agenda");
856  iy_main_agenda.check(ws, verbosity);
857 
858  // Variables related to the top of the cloudbox
859  const Index i0 = cloudbox_limits[1];
860  const Numeric z_top = z_field(i0 + 1, 0, 0); // Note i0+1
861 
862  // Loop zenith angles
863  //
864  if (nza)
865 #pragma omp parallel for if (!arts_omp_in_parallel() && nza > 1 && \
866  !use_parallel_iy) firstprivate(l_ws)
867  for (Index i = 0; i < nza; i++) {
868  if (failed) continue;
869  try {
870  // Define output variables
871  Ppath ppath;
879 
880  Index iy_id = i;
881  Vector rte_los(1, za_grid[i]);
882  Vector rte_pos(1, za_grid[i] < 90 ? z_top : z_space);
883 
886  z_field,
887  z_surface,
888  cloudbox_on,
891  rte_pos,
892  rte_los,
893  ppath_lmax,
894  verbosity);
895  assert(ppath.gp_p[ppath.np - 1].idx == i0 ||
896  ppath.gp_p[ppath.np - 1].idx == nl - 2);
897 
898  if (use_parallel_iy) {
899  iyEmissionStandard(l_ws,
900  iy,
901  iy_aux,
902  diy_dx,
903  ppvar_p,
904  ppvar_t,
905  ppvar_nlte,
906  ppvar_vmr,
907  ppvar_wind,
908  ppvar_mag,
909  ppvar_f,
910  ppvar_iy,
913  iy_id,
914  stokes_dim,
915  f_grid,
917  p_grid,
918  t_field,
919  nlte_field,
920  vmr_field,
921  abs_species,
922  wind_u_field,
923  wind_v_field,
924  wind_w_field,
925  mag_u_field,
926  mag_v_field,
927  mag_w_field,
928  cloudbox_on,
929  iy_unit,
930  iy_aux_vars,
931  jacobian_do,
933  ppath,
934  rte_pos2,
945  verbosity);
946  } else {
948  iy,
949  iy_aux,
950  diy_dx,
951  ppvar_p,
952  ppvar_t,
953  ppvar_nlte,
954  ppvar_vmr,
955  ppvar_wind,
956  ppvar_mag,
957  ppvar_f,
958  ppvar_iy,
961  iy_id,
962  stokes_dim,
963  f_grid,
965  p_grid,
966  t_field,
967  nlte_field,
968  vmr_field,
969  abs_species,
970  wind_u_field,
971  wind_v_field,
972  wind_w_field,
973  mag_u_field,
974  mag_v_field,
975  mag_w_field,
976  cloudbox_on,
977  iy_unit,
978  iy_aux_vars,
979  jacobian_do,
981  ppath,
982  rte_pos2,
993  verbosity);
994  }
995  assert(iy.ncols() == stokes_dim);
996 
997  // First and last points are most easily handled separately
998  // But field at top cloudbox already known from copying above
999  if (za_grid[i] < 90) {
1000  spectral_radiance_field(joker, i0 + 1, 0, 0, i, 0, joker) =
1001  ppvar_iy(joker, joker, 0);
1002  spectral_radiance_field(joker, nl - 1, 0, 0, i, 0, joker) =
1003  ppvar_iy(joker, joker, ppath.np - 1);
1004  } else {
1005  spectral_radiance_field(joker, nl - 1, 0, 0, i, 0, joker) =
1006  ppvar_iy(joker, joker, 0);
1007  }
1008 
1009  // Remaining points
1010  for (Index p = 1; p < ppath.np - 1; p++) {
1011  // We just store values at pressure levels
1012  if (ppath.gp_p[p].fd[0] < 1e-2) {
1014  joker, ppath.gp_p[p].idx, 0, 0, i, 0, joker) =
1015  ppvar_iy(joker, joker, p);
1016  }
1017  }
1018 
1019  } catch (const std::exception& e) {
1020 #pragma omp critical(planep_setabort)
1021  failed = true;
1022 
1023  ostringstream os;
1024  os << "Run-time error at nza #" << i << ": \n" << e.what();
1025 #pragma omp critical(planep_push_fail_msg)
1026  fail_msg.push_back(os.str());
1027  }
1028  }
1029 
1030  if (fail_msg.nelem()) {
1031  ostringstream os;
1032  for (auto& msg : fail_msg) os << msg << '\n';
1033  throw runtime_error(os.str());
1034  }
1035 }
ARTS::Var::ppvar_mag
Matrix ppvar_mag(Workspace &ws) noexcept
Definition: autoarts.h:5259
Matrix
The Matrix class.
Definition: matpackI.h:1193
ARTS::Var::ppvar_trans_partial
Tensor4 ppvar_trans_partial(Workspace &ws) noexcept
Definition: autoarts.h:5359
ARTS::Var::atmosphere_dim
Index atmosphere_dim(Workspace &ws) noexcept
Definition: autoarts.h:2510
ARTS::Var::wind_v_field
Tensor3 wind_v_field(Workspace &ws) noexcept
Definition: autoarts.h:7246
ARTS::Var::radiance_field
Tensor5 radiance_field(Workspace &ws) noexcept
Definition: autoarts.h:5476
ARTS::Var::z_field
Tensor3 z_field(Workspace &ws) noexcept
Definition: autoarts.h:7690
Tensor4::resize
void resize(Index b, Index p, Index r, Index c)
Resize function.
Definition: matpackIV.cc:1064
ConstTensor7View::nshelves
Index nshelves() const
Returns the number of shelves.
Definition: matpackVII.cc:48
ConstTensor5View::nbooks
Index nbooks() const
Returns the number of books.
Definition: matpackV.cc:47
spectral_radiance_fieldExpandCloudboxField
void spectral_radiance_fieldExpandCloudboxField(Workspace &ws, Tensor7 &spectral_radiance_field, const Agenda &propmat_clearsky_agenda, const Agenda &water_p_eq_agenda, const Agenda &iy_space_agenda, const Agenda &iy_surface_agenda, const Agenda &iy_cloudbox_agenda, const Index &stokes_dim, const Vector &f_grid, const Index &atmosphere_dim, const Vector &p_grid, const Tensor3 &z_field, const Tensor3 &t_field, const EnergyLevelMap &nlte_field, const Tensor4 &vmr_field, const ArrayOfArrayOfSpeciesTag &abs_species, const Tensor3 &wind_u_field, const Tensor3 &wind_v_field, const Tensor3 &wind_w_field, const Tensor3 &mag_u_field, const Tensor3 &mag_v_field, const Tensor3 &mag_w_field, const Matrix &z_surface, const Index &cloudbox_on, const ArrayOfIndex &cloudbox_limits, const Tensor7 &cloudbox_field, const Numeric &ppath_lmax, const Numeric &rte_alonglos_v, const Tensor3 &surface_props_data, const Vector &za_grid, const Index &use_parallel_iy, const Verbosity &verbosity)
WORKSPACE METHOD: spectral_radiance_fieldExpandCloudboxField.
Definition: m_fluxes.cc:774
ARTS::Var::g0
Numeric g0(Workspace &ws) noexcept
Definition: autoarts.h:3516
ARTS::Var::heating_rates
Tensor3 heating_rates(Workspace &ws) noexcept
Definition: autoarts.h:3568
auto_md.h
ARTS::Var::iy_cloudbox_agenda
Agenda iy_cloudbox_agenda(Workspace &ws) noexcept
Definition: autoarts.h:3742
ARTS::Var::jacobian_quantities
ArrayOfRetrievalQuantity jacobian_quantities(Workspace &ws) noexcept
Definition: autoarts.h:3924
absorption.h
Declarations required for the calculation of absorption coefficients.
w
Complex w(Complex z) noexcept
The Faddeeva function.
Definition: linefunctions.cc:42
TokVal
This stores arbitrary token values and remembers the type.
Definition: token.h:40
Tensor3
The Tensor3 class.
Definition: matpackIII.h:339
surface.h
ARTS::Var::ppvar_vmr
Matrix ppvar_vmr(Workspace &ws) noexcept
Definition: autoarts.h:5372
ConstTensor5View::ncols
Index ncols() const
Returns the number of columns.
Definition: matpackV.cc:56
ARTS::Var::iy_transmission
Tensor3 iy_transmission(Workspace &ws) noexcept
Definition: autoarts.h:3836
joker
const Joker joker
ARTS::Var::iy_id
Index iy_id(Workspace &ws) noexcept
Definition: autoarts.h:3775
ARTS::Var::ppath_inside_cloudbox_do
Index ppath_inside_cloudbox_do(Workspace &ws) noexcept
Definition: autoarts.h:5172
abs
#define abs(x)
Definition: legacy_continua.cc:20626
ARTS::Var::verbosity
Verbosity verbosity(Workspace &ws) noexcept
Definition: autoarts.h:7112
ARTS::Var::spectral_irradiance_field
Tensor5 spectral_irradiance_field(Workspace &ws) noexcept
Definition: autoarts.h:6547
ARTS::Var::ppvar_iy
Tensor3 ppvar_iy(Workspace &ws) noexcept
Definition: autoarts.h:5246
ARTS::Var::cloudbox_limits
ArrayOfIndex cloudbox_limits(Workspace &ws) noexcept
Definition: autoarts.h:2762
ARTS::Var::jacobian_do
Index jacobian_do(Workspace &ws) noexcept
Definition: autoarts.h:3912
ARTS::Var::wind_u_field
Tensor3 wind_u_field(Workspace &ws) noexcept
Definition: autoarts.h:7211
Vector::resize
void resize(Index n)
Resize function.
Definition: matpackI.cc:404
legendre.h
Contains the code to calculate Legendre polynomials.
ARTS::Var::diy_dx
ArrayOfTensor3 diy_dx(Workspace &ws) noexcept
Definition: autoarts.h:3007
ARTS::Var::ppvar_f
Matrix ppvar_f(Workspace &ws) noexcept
Definition: autoarts.h:5233
ARTS::Var::mag_w_field
Tensor3 mag_w_field(Workspace &ws) noexcept
Definition: autoarts.h:4197
gsl_integration_glfixed_table_alloc
bool gsl_integration_glfixed_table_alloc(Vector &x, Vector &w, Index n)
gsl_integration_glfixed_table_alloc
Definition: legendre.cc:2685
ARTS::Var::stokes_dim
Index stokes_dim(Workspace &ws) noexcept
Definition: autoarts.h:6650
spectral_radiance_fieldClearskyPlaneParallel
void spectral_radiance_fieldClearskyPlaneParallel(Workspace &ws, Tensor7 &spectral_radiance_field, Tensor3 &trans_field, const Agenda &propmat_clearsky_agenda, const Agenda &water_p_eq_agenda, const Agenda &iy_space_agenda, const Agenda &iy_surface_agenda, const Agenda &iy_cloudbox_agenda, const Index &stokes_dim, const Vector &f_grid, const Index &atmosphere_dim, const Vector &p_grid, const Tensor3 &z_field, const Tensor3 &t_field, const EnergyLevelMap &nlte_field, const Tensor4 &vmr_field, const ArrayOfArrayOfSpeciesTag &abs_species, const Tensor3 &wind_u_field, const Tensor3 &wind_v_field, const Tensor3 &wind_w_field, const Tensor3 &mag_u_field, const Tensor3 &mag_v_field, const Tensor3 &mag_w_field, const Matrix &z_surface, const Numeric &ppath_lmax, const Numeric &rte_alonglos_v, const Tensor3 &surface_props_data, const Vector &za_grid, const Index &use_parallel_iy, const Verbosity &verbosity)
WORKSPACE METHOD: spectral_radiance_fieldClearskyPlaneParallel.
Definition: m_fluxes.cc:483
ARTS::Var::iy_unit
String iy_unit(Workspace &ws) noexcept
Definition: autoarts.h:3856
index_of_zsurface
Index index_of_zsurface(const Numeric &z_surface, ConstVectorView z_profile)
Lccates the surface with respect to pressure levels.
Definition: surface.cc:54
ARTS::Var::iy
Matrix iy(Workspace &ws) noexcept
Definition: autoarts.h:3690
ARTS::Var::irradiance_field
Tensor4 irradiance_field(Workspace &ws) noexcept
Definition: autoarts.h:3659
ppathPlaneParallel
void ppathPlaneParallel(Ppath &ppath, const Index &atmosphere_dim, const Tensor3 &z_field, const Matrix &z_surface, const Index &cloudbox_on, const ArrayOfIndex &cloudbox_limits, const Index &ppath_inside_cloudbox_do, const Vector &rte_pos, const Vector &rte_los, const Numeric &ppath_lmax, const Verbosity &)
WORKSPACE METHOD: ppathPlaneParallel.
Definition: m_ppath.cc:761
Tensor4
The Tensor4 class.
Definition: matpackIV.h:421
ARTS::Var::iy_space_agenda
Agenda iy_space_agenda(Workspace &ws) noexcept
Definition: autoarts.h:3801
Ppath
The structure to describe a propagation path and releated quantities.
Definition: ppath.h:48
Tensor3::resize
void resize(Index p, Index r, Index c)
Resize function.
Definition: matpackIII.cc:664
Agenda
The Agenda class.
Definition: agenda_class.h:44
ConstTensor5View::npages
Index npages() const
Returns the number of pages.
Definition: matpackV.cc:50
AngularGridsSetFluxCalc
void AngularGridsSetFluxCalc(Vector &za_grid, Vector &aa_grid, Vector &za_grid_weights, const Index &N_za_grid, const Index &N_aa_grid, const String &za_grid_type, const Verbosity &)
WORKSPACE METHOD: AngularGridsSetFluxCalc.
Definition: m_fluxes.cc:55
ARTS::Var::ppvar_trans_cumulat
Tensor4 ppvar_trans_cumulat(Workspace &ws) noexcept
Definition: autoarts.h:5346
ConstTensor7View::nlibraries
Index nlibraries() const
Returns the number of libraries.
Definition: matpackVII.cc:42
ARTS::Var::rte_alonglos_v
Numeric rte_alonglos_v(Workspace &ws) noexcept
Definition: autoarts.h:5628
ARTS::Var::nlte_field
EnergyLevelMap nlte_field(Workspace &ws) noexcept
Definition: autoarts.h:4604
Array
This can be used to make arrays out of anything.
Definition: array.h:108
RadiationFieldSpectralIntegrate
void RadiationFieldSpectralIntegrate(Tensor4 &radiation_field, const Vector &f_grid, const Tensor5 &spectral_radiation_field, const Verbosity &)
WORKSPACE METHOD: RadiationFieldSpectralIntegrate.
Definition: m_fluxes.cc:328
agenda_class.h
Declarations for agendas.
Tensor5::resize
void resize(Index s, Index b, Index p, Index r, Index c)
Resize function.
Definition: matpackV.cc:1743
messages.h
Declarations having to do with the four output streams.
ARTS::Var::rte_pos
Vector rte_pos(Workspace &ws) noexcept
Definition: autoarts.h:5672
PI
const Numeric PI
my_basic_string< char >
heating_ratesFromIrradiance
void heating_ratesFromIrradiance(Tensor3 &heating_rates, const Vector &p_grid, const Tensor4 &irradiance_field, const Tensor3 &specific_heat_capacity, const Numeric &g0, const Verbosity &)
WORKSPACE METHOD: heating_ratesFromIrradiance.
Definition: m_fluxes.cc:187
ConstVectorView::nelem
Index nelem() const
Returns the number of elements.
Definition: matpackI.cc:51
ConstTensor7View::nvitrines
Index nvitrines() const
Returns the number of vitrines.
Definition: matpackVII.cc:45
ARTS::Var::vmr_field
Tensor4 vmr_field(Workspace &ws) noexcept
Definition: autoarts.h:7130
ARTS::Var::rte_los
Vector rte_los(Workspace &ws) noexcept
Definition: autoarts.h:5651
ARTS::Var::ppvar_nlte
EnergyLevelMap ppvar_nlte(Workspace &ws) noexcept
Definition: autoarts.h:5272
ARTS::Var::abs_species
ArrayOfArrayOfSpeciesTag abs_species(Workspace &ws) noexcept
Definition: autoarts.h:2157
ARTS::Var::ppvar_t
Vector ppvar_t(Workspace &ws) noexcept
Definition: autoarts.h:5330
iyEmissionStandard
void iyEmissionStandard(Workspace &ws, Matrix &iy, ArrayOfMatrix &iy_aux, ArrayOfTensor3 &diy_dx, Vector &ppvar_p, Vector &ppvar_t, EnergyLevelMap &ppvar_nlte, Matrix &ppvar_vmr, Matrix &ppvar_wind, Matrix &ppvar_mag, Matrix &ppvar_f, Tensor3 &ppvar_iy, Tensor4 &ppvar_trans_cumulat, Tensor4 &ppvar_trans_partial, const Index &iy_id, const Index &stokes_dim, const Vector &f_grid, const Index &atmosphere_dim, const Vector &p_grid, const Tensor3 &t_field, const EnergyLevelMap &nlte_field, const Tensor4 &vmr_field, const ArrayOfArrayOfSpeciesTag &abs_species, const Tensor3 &wind_u_field, const Tensor3 &wind_v_field, const Tensor3 &wind_w_field, const Tensor3 &mag_u_field, const Tensor3 &mag_v_field, const Tensor3 &mag_w_field, const Index &cloudbox_on, const String &iy_unit, const ArrayOfString &iy_aux_vars, const Index &jacobian_do, const ArrayOfRetrievalQuantity &jacobian_quantities, const Ppath &ppath, const Vector &rte_pos2, const Agenda &propmat_clearsky_agenda, const Agenda &water_p_eq_agenda, const Agenda &iy_main_agenda, const Agenda &iy_space_agenda, const Agenda &iy_surface_agenda, const Agenda &iy_cloudbox_agenda, const Index &iy_agenda_call1, const Tensor3 &iy_transmission, const Numeric &rte_alonglos_v, const Tensor3 &surface_props_data, const Verbosity &verbosity)
WORKSPACE METHOD: iyEmissionStandard.
Definition: m_rte.cc:563
ARTS::Var::p_grid
Vector p_grid(Workspace &ws) noexcept
Definition: autoarts.h:4763
ARTS::Var::ppath
Ppath ppath(Workspace &ws) noexcept
Definition: autoarts.h:5139
ConstTensor4View::ncols
Index ncols() const
Returns the number of columns.
Definition: matpackIV.cc:66
ConstTensor7View::npages
Index npages() const
Returns the number of pages.
Definition: matpackVII.cc:54
ARTS::Var::surface_props_data
Tensor3 surface_props_data(Workspace &ws) noexcept
Definition: autoarts.h:6745
Numeric
NUMERIC Numeric
The type to use for all floating point numbers.
Definition: matpack.h:33
Verbosity
Definition: messages.h:49
iyEmissionStandardSequential
void iyEmissionStandardSequential(Workspace &ws, Matrix &iy, ArrayOfMatrix &iy_aux, ArrayOfTensor3 &diy_dx, Vector &ppvar_p, Vector &ppvar_t, EnergyLevelMap &ppvar_nlte, Matrix &ppvar_vmr, Matrix &ppvar_wind, Matrix &ppvar_mag, Matrix &ppvar_f, Tensor3 &ppvar_iy, Tensor4 &ppvar_trans_cumulat, Tensor4 &ppvar_trans_partial, const Index &iy_id, const Index &stokes_dim, const Vector &f_grid, const Index &atmosphere_dim, const Vector &p_grid, const Tensor3 &t_field, const EnergyLevelMap &nlte_field, const Tensor4 &vmr_field, const ArrayOfArrayOfSpeciesTag &abs_species, const Tensor3 &wind_u_field, const Tensor3 &wind_v_field, const Tensor3 &wind_w_field, const Tensor3 &mag_u_field, const Tensor3 &mag_v_field, const Tensor3 &mag_w_field, const Index &cloudbox_on, const String &iy_unit, const ArrayOfString &iy_aux_vars, const Index &jacobian_do, const ArrayOfRetrievalQuantity &jacobian_quantities, const Ppath &ppath, const Vector &rte_pos2, const Agenda &propmat_clearsky_agenda, const Agenda &water_p_eq_agenda, const Agenda &iy_main_agenda, const Agenda &iy_space_agenda, const Agenda &iy_surface_agenda, const Agenda &iy_cloudbox_agenda, const Index &iy_agenda_call1, const Tensor3 &iy_transmission, const Numeric &rte_alonglos_v, const Tensor3 &surface_props_data, const Verbosity &verbosity)
WORKSPACE METHOD: iyEmissionStandardSequential.
Definition: m_rte.cc:172
ConstTensor4View::npages
Index npages() const
Returns the number of pages.
Definition: matpackIV.cc:60
ARTS::Var::iy_agenda_call1
Index iy_agenda_call1(Workspace &ws) noexcept
Definition: autoarts.h:3704
ARTS::Var::specific_heat_capacity
Tensor3 specific_heat_capacity(Workspace &ws) noexcept
Definition: autoarts.h:6527
ConstTensor4View::nbooks
Index nbooks() const
Returns the number of books.
Definition: matpackIV.cc:57
math_funcs.h
ARTS::Var::iy_aux_vars
ArrayOfString iy_aux_vars(Workspace &ws) noexcept
Definition: autoarts.h:3735
DEG2RAD
const Numeric DEG2RAD
ARTS::Var::f_grid
Vector f_grid(Workspace &ws) noexcept
Definition: autoarts.h:3449
EnergyLevelMap
Definition: energylevelmap.h:60
ARTS::Var::ppvar_wind
Matrix ppvar_wind(Workspace &ws) noexcept
Definition: autoarts.h:5385
irradiance_fieldFromRadiance
void irradiance_fieldFromRadiance(Tensor4 &irradiance_field, const Tensor5 &radiance_field, const Vector &za_grid, const Vector &aa_grid, const Vector &za_grid_weights, const Verbosity &)
WORKSPACE METHOD: irradiance_fieldFromRadiance.
Definition: m_fluxes.cc:256
ARTS::Var::za_grid
Vector za_grid(Workspace &ws) noexcept
Definition: autoarts.h:7771
nlinspace
void nlinspace(Vector &x, const Numeric start, const Numeric stop, const Index n)
nlinspace
Definition: math_funcs.cc:231
ARTS::Var::ppath_lmax
Numeric ppath_lmax(Workspace &ws) noexcept
Definition: autoarts.h:5183
ARTS::Var::iy_surface_agenda
Agenda iy_surface_agenda(Workspace &ws) noexcept
Definition: autoarts.h:3808
Tensor5
The Tensor5 class.
Definition: matpackV.h:506
spectral_radiance_fieldCopyCloudboxField
void spectral_radiance_fieldCopyCloudboxField(Tensor7 &spectral_radiance_field, const Index &atmosphere_dim, const Vector &p_grid, const Index &cloudbox_on, const ArrayOfIndex &cloudbox_limits, const Tensor7 &cloudbox_field, const Verbosity &)
WORKSPACE METHOD: spectral_radiance_fieldCopyCloudboxField.
Definition: m_fluxes.cc:752
ConstTensor4View::nrows
Index nrows() const
Returns the number of rows.
Definition: matpackIV.cc:63
Range
The range class.
Definition: matpackI.h:160
spectral_irradiance_fieldFromSpectralRadianceField
void spectral_irradiance_fieldFromSpectralRadianceField(Tensor5 &spectral_irradiance_field, const Tensor7 &spectral_radiance_field, const Vector &za_grid, const Vector &aa_grid, const Vector &za_grid_weights, const Verbosity &)
WORKSPACE METHOD: spectral_irradiance_fieldFromSpectralRadianceField.
Definition: m_fluxes.cc:405
ConstTensor7View::nbooks
Index nbooks() const
Returns the number of books.
Definition: matpackVII.cc:51
ConstTensor5View::nrows
Index nrows() const
Returns the number of rows.
Definition: matpackV.cc:53
ARTS::Var::iy_aux
ArrayOfMatrix iy_aux(Workspace &ws) noexcept
Definition: autoarts.h:3719
workspace_ng.h
This file contains the Workspace class.
ARTS::Var::iy_main_agenda
Agenda iy_main_agenda(Workspace &ws) noexcept
Definition: autoarts.h:3794
ARTS::Var::rte_pos2
Vector rte_pos2(Workspace &ws) noexcept
Definition: autoarts.h:5696
ConstTensor5View::nshelves
Index nshelves() const
Returns the number of shelves.
Definition: matpackV.cc:44
ARTS::Var::mag_v_field
Tensor3 mag_v_field(Workspace &ws) noexcept
Definition: autoarts.h:4164
ARTS::Var::t_field
Tensor3 t_field(Workspace &ws) noexcept
Definition: autoarts.h:6947
Workspace
Workspace class.
Definition: workspace_ng.h:40
ARTS::Var::za_grid_weights
Vector za_grid_weights(Workspace &ws) noexcept
Definition: autoarts.h:7780
ARTS::Var::wind_w_field
Tensor3 wind_w_field(Workspace &ws) noexcept
Definition: autoarts.h:7279
calculate_weights_linear
void calculate_weights_linear(Vector &x, Vector &w, const Index nph)
calculate_weights_linear
Definition: math_funcs.cc:813
ARTS::Var::cloudbox_field
Tensor7 cloudbox_field(Workspace &ws) noexcept
Definition: autoarts.h:2676
ARTS::Var::propmat_clearsky_agenda
Agenda propmat_clearsky_agenda(Workspace &ws) noexcept
Definition: autoarts.h:5405
ARTS::Var::x
Vector x(Workspace &ws) noexcept
Definition: autoarts.h:7346
Index
INDEX Index
The type to use for all integer numbers and indices.
Definition: matpack.h:39
ConstTensor7View::nrows
Index nrows() const
Returns the number of rows.
Definition: matpackVII.cc:57
check_input.h
ARTS::Var::mag_u_field
Tensor3 mag_u_field(Workspace &ws) noexcept
Definition: autoarts.h:4130
Vector
The Vector class.
Definition: matpackI.h:860
ARTS::Var::spectral_radiance_field
Tensor7 spectral_radiance_field(Workspace &ws) noexcept
Definition: autoarts.h:6573
ARTS::Var::ppvar_p
Vector ppvar_p(Workspace &ws) noexcept
Definition: autoarts.h:5304
ARTS::Var::aa_grid
Vector aa_grid(Workspace &ws) noexcept
Definition: autoarts.h:1717
ARTS::Var::water_p_eq_agenda
Agenda water_p_eq_agenda(Workspace &ws) noexcept
Definition: autoarts.h:7165
sorting.h
Contains sorting routines.
Array::nelem
Index nelem() const
Number of elements.
Definition: array.h:195
ARTS::Var::z_surface
Matrix z_surface(Workspace &ws) noexcept
Definition: autoarts.h:7754
Tensor7
The Tensor7 class.
Definition: matpackVII.h:2382
matpackVII.h
ARTS::Var::cloudbox_on
Index cloudbox_on(Workspace &ws) noexcept
Definition: autoarts.h:2782