Go to the documentation of this file.
73 out1 <<
" Executing for loop body, index: " << i <<
"\n";
93 Index first_ybatch_index = 0;
96 bool do_abort =
false;
108 Index job_counter = 0;
128 #pragma omp parallel for schedule(dynamic) if (!arts_omp_in_parallel() && \
130 firstprivate(l_ws, l_ybatch_calc_agenda)
135 if (do_abort)
continue;
136 #pragma omp critical(ybatchCalc_job_counter)
137 { l_job_counter = ++job_counter; }
141 os <<
" Job " << l_job_counter <<
" of " <<
ybatch_n <<
", Index "
157 l_ybatch_calc_agenda);
160 #pragma omp critical(ybatchCalc_assign_y)
162 #pragma omp critical(ybatchCalc_assign_y_aux)
169 if (Knr != 0 || Knc != 0) {
170 if (Knr !=
y.nelem()) {
172 os <<
"First dimension of Jacobian must have same length as the measurement *y*.\n"
173 <<
"Length of *y*: " <<
y.nelem() <<
"\n"
174 <<
"Dimensions of *jacobian*: (" << Knr <<
", " << Knc
180 #pragma omp critical(ybatchCalc_setabort)
183 throw runtime_error(os.str());
192 }
catch (
const std::exception& e) {
193 if (robust && !do_abort) {
197 <<
"y Vector in output variable ybatch will be empty for this job.\n"
198 <<
"The runtime error produced was:\n"
204 #pragma omp critical(ybatchCalc_setabort)
209 <<
" failed. Aborting...\n";
216 #pragma omp critical(ybatchCalc_push_fail_msg)
217 fail_msg.push_back(os.str());
221 if (fail_msg.
nelem()) {
224 if (!do_abort) os <<
"\nError messages from failed batch cases:\n";
225 for (ArrayOfString::const_iterator it = fail_msg.begin();
226 it != fail_msg.end();
231 throw runtime_error(os.str());
253 const Index& nelem_p_grid,
254 const String& met_profile_path,
255 const String& met_profile_pnd_path,
299 ybatch.resize(no_profiles);
302 for (
Index i = 0; i < no_profiles; ++i) {
303 ostringstream lat_os, lon_os;
306 if (
lat[i] < 0) lat_prec--;
309 if (
abs(
lat[i]) >= 100) lat_prec--;
312 lat_os.setf(ios::showpoint | ios::fixed);
313 lat_os << setprecision((
int)lat_prec) <<
lat[i];
316 if (
lon[i] < 0) lon_prec--;
319 if (
abs(
lon[i]) >= 100) lon_prec--;
321 lon_os.setf(ios::showpoint | ios::fixed);
322 lon_os << setprecision((
int)lon_prec) <<
lon[i];
331 ".lon_" + lon_os.str() +
".t.xml",
337 ".lon_" + lon_os.str() +
".z.xml",
343 ".lon_" + lon_os.str() +
".H2O.xml",
353 lat_os.str() +
".lon_" + lon_os.str() +
".pnd15.xml",
386 const Vector& tfr_p_grid =
399 Numeric cl_grid_min, cl_grid_max;
403 cl_grid_min = tfr_p_grid[0];
406 Index level_counter = 0;
409 for (
Index ip = 0; ip < N_p; ++ip) {
419 cl_grid_max = tfr_p_grid[ip + 1];
429 if (level_counter == 0) {
487 const Index& nelem_p_grid,
488 const String& met_profile_path,
509 ybatch.resize(no_profiles);
511 Vector sat_za_from_profile;
525 for (
Index i = 0; i < no_profiles; ++i) {
526 ostringstream lat_os, lon_os;
529 if (
lat[i] < 0) lat_prec--;
532 if (
abs(
lat[i]) >= 100) lat_prec--;
535 lat_os.setf(ios::showpoint | ios::fixed);
536 lat_os << setprecision((
int)lat_prec) <<
lat[i];
539 if (
lon[i] < 0) lon_prec--;
542 if (
abs(
lon[i]) >= 100) lon_prec--;
544 lon_os.setf(ios::showpoint | ios::fixed);
545 lon_os << setprecision((
int)lon_prec) <<
lon[i];
546 cout << lat_os.str() << endl;
547 cout << lon_os.str() << endl;
549 sat_za = sat_za_from_profile[i];
555 cout <<
"sensor_los" << sat_za_from_profile[i] << endl;
556 cout <<
"sensor_los" << sat_za << endl;
561 ".lon_" + lon_os.str() +
".t.xml",
566 ".lon_" + lon_os.str() +
".z.xml",
574 ".lon_" + lon_os.str() +
".H2O.xml",
581 <<
"--------------------------------------------------------------------------"
584 << met_profile_path +
"profile.lat_" + lat_os.str() +
".lon_" +
586 <<
"is executed now" << endl;
588 <<
"--------------------------------------------------------------------------"
596 cout <<
"z_surface" <<
z_surface << endl;
597 const Vector& tfr_p_grid =
627 cout <<
"t_field_raw[0](0,0,0)" << tfr_p_grid[0] << endl;
628 cout <<
"t_field_raw[0](N_p -1,0,0)" << tfr_p_grid[N_p - 1] << endl;
664 Index first_ybatch_index = 0;
667 bool do_abort =
false;
679 Index job_counter = 0;
695 #pragma omp parallel for schedule(dynamic) if (!arts_omp_in_parallel() && \
697 firstprivate(l_ws, l_dobatch_calc_agenda)
702 if (do_abort)
continue;
703 #pragma omp critical(dobatchCalc_job_counter)
704 { l_job_counter = ++job_counter; }
708 os <<
" Job " << l_job_counter <<
" of " <<
ybatch_n <<
", Index "
726 l_dobatch_calc_agenda);
728 #pragma omp critical(dobatchCalc_assign_cloudbox_field)
730 #pragma omp critical(dobatchCalc_assign_radiance_field)
732 #pragma omp critical(dobatchCalc_assign_irradiance_field)
734 #pragma omp critical(dobatchCalc_assign_spectral_irradiance_field)
738 }
catch (
const std::exception& e) {
739 if (robust && !do_abort) {
743 <<
"element in output variables will be empty for this job.\n"
744 <<
"The runtime error produced was:\n"
750 #pragma omp critical(dobatchCalc_setabort)
755 <<
" failed. Aborting...\n";
762 #pragma omp critical(dobatchCalc_push_fail_msg)
763 fail_msg.push_back(os.str());
767 if (fail_msg.
nelem()) {
770 if (!do_abort) os <<
"\nError messages from failed batch cases:\n";
771 for (ArrayOfString::const_iterator it = fail_msg.begin();
772 it != fail_msg.end();
777 throw runtime_error(os.str());
Index atmosphere_dim(Workspace &ws) noexcept
ArrayOfVector y_aux(Workspace &ws) noexcept
int arts_omp_get_thread_num()
Wrapper for omp_get_thread_num.
Tensor5 radiance_field(Workspace &ws) noexcept
Agenda forloop_agenda(Workspace &ws) noexcept
void ybatchMetProfilesClear(Workspace &ws, ArrayOfVector &ybatch, const ArrayOfArrayOfSpeciesTag &abs_species, const Agenda &met_profile_calc_agenda, const Vector &f_grid, const Matrix &met_amsu_data, const Matrix &sensor_pos, const Vector &refellipsoid, const Index &nelem_p_grid, const String &met_profile_path, const Verbosity &verbosity)
WORKSPACE METHOD: ybatchMetProfilesClear.
Numeric lat(Workspace &ws) noexcept
Index ybatch_n(Workspace &ws) noexcept
ArrayOfTensor4 dobatch_irradiance_field(Workspace &ws) noexcept
void cloudboxSetManually(Index &cloudbox_on, ArrayOfIndex &cloudbox_limits, const Index &atmosphere_dim, const Vector &p_grid, const Vector &lat_grid, const Vector &lon_grid, const Numeric &p1, const Numeric &p2, const Numeric &lat1, const Numeric &lat2, const Numeric &lon1, const Numeric &lon2, const Verbosity &)
WORKSPACE METHOD: cloudboxSetManually.
Vector lat_grid(Workspace &ws) noexcept
ArrayOfGriddedField3 pnd_field_raw(Workspace &ws) noexcept
Verbosity verbosity(Workspace &ws) noexcept
Tensor5 spectral_irradiance_field(Workspace &ws) noexcept
void DOBatchCalc(Workspace &ws, ArrayOfTensor7 &dobatch_cloudbox_field, ArrayOfTensor5 &dobatch_radiance_field, ArrayOfTensor4 &dobatch_irradiance_field, ArrayOfTensor5 &dobatch_spectral_irradiance_field, const Index &ybatch_start, const Index &ybatch_n, const Agenda &dobatch_calc_agenda, const Index &robust, const Verbosity &verbosity)
WORKSPACE METHOD: DOBatchCalc.
void std(VectorView std, const Vector &y, const ArrayOfVector &ys, const Index start=0, const Index end=-1)
Compute the standard deviation of the ranged ys.
ArrayOfIndex cloudbox_limits(Workspace &ws) noexcept
Vector y(Workspace &ws) noexcept
G0 G2 FVC Y DV Numeric Numeric Numeric Zeeman LowerQuantumNumbers void * data
Tensor4 irradiance_field(Workspace &ws) noexcept
Index TotalNumberOfElements(const Array< Array< base > > &aa)
Determine total number of elements in an ArrayOfArray.
Numeric lon(Workspace &ws) noexcept
Agenda dobatch_calc_agenda(Workspace &ws) noexcept
ArrayOfTensor5 dobatch_spectral_irradiance_field(Workspace &ws) noexcept
void ybatch_calc_agendaExecute(Workspace &ws, Vector &y, ArrayOfVector &y_aux, Matrix &jacobian, const Index ybatch_index, const Agenda &input_agenda)
void forloop_agendaExecute(Workspace &ws, const Index forloop_index, const Agenda &input_agenda)
void xml_read_from_file(const String &filename, T &type, const Verbosity &verbosity)
Reads data from XML file.
The implementation for String, the ARTS string class.
ArrayOfTensor7 dobatch_cloudbox_field(Workspace &ws) noexcept
GriddedField3 t_field_raw(Workspace &ws) noexcept
Index nelem() const
Returns the number of elements.
This file contains declerations of functions of physical character.
ArrayOfArrayOfSpeciesTag abs_species(Workspace &ws) noexcept
Vector p_grid(Workspace &ws) noexcept
NUMERIC Numeric
The type to use for all floating point numbers.
Vector refellipsoid(Workspace &ws) noexcept
Vector f_grid(Workspace &ws) noexcept
Matrix sensor_pos(Workspace &ws) noexcept
ArrayOfArrayOfSingleScatteringData scat_data(Workspace &ws) noexcept
Vector lon_grid(Workspace &ws) noexcept
Agenda met_profile_calc_agenda(Workspace &ws) noexcept
ArrayOfTensor5 dobatch_radiance_field(Workspace &ws) noexcept
void xml_write_to_file(const String &filename, const T &type, const FileType ftype, const Index no_clobber, const Verbosity &verbosity)
Write data to XML file.
Index ybatch_index(Workspace &ws) noexcept
constexpr Rational start(Rational Ju, Rational Jl, Polarization type) noexcept
Gives the lowest M for a polarization type of this transition.
void ybatchMetProfiles(Workspace &ws, ArrayOfVector &ybatch, const ArrayOfArrayOfSpeciesTag &abs_species, const Agenda &met_profile_calc_agenda, const Vector &f_grid, const Matrix &met_amsu_data, const Matrix &sensor_pos, const Vector &refellipsoid, const Vector &lat_grid, const Vector &lon_grid, const Index &atmosphere_dim, const ArrayOfArrayOfSingleScatteringData &scat_data, const Index &nelem_p_grid, const String &met_profile_path, const String &met_profile_pnd_path, const Verbosity &verbosity)
WORKSPACE METHOD: ybatchMetProfiles.
void ForLoop(Workspace &ws, const Agenda &forloop_agenda, const Index &start, const Index &stop, const Index &step, const Verbosity &verbosity)
WORKSPACE METHOD: ForLoop.
Declaration of functions in rte.cc.
ArrayOfMatrix ybatch_jacobians(Workspace &ws) noexcept
ArrayOfArrayOfVector ybatch_aux(Workspace &ws) noexcept
void VectorNLogSpace(Vector &x, const Index &n, const Numeric &start, const Numeric &stop, const Verbosity &verbosity)
WORKSPACE METHOD: VectorNLogSpace.
Matrix met_amsu_data(Workspace &ws) noexcept
Tensor7 cloudbox_field(Workspace &ws) noexcept
Matrix jacobian(Workspace &ws) noexcept
const Index GFIELD3_P_GRID
INDEX Index
The type to use for all integer numbers and indices.
ArrayOfGriddedField3 vmr_field_raw(Workspace &ws) noexcept
Header file for helper functions for OpenMP.
GriddedField3 z_field_raw(Workspace &ws) noexcept
void dobatch_calc_agendaExecute(Workspace &ws, Tensor7 &spectral_radiance_field, Tensor5 &radiance_field, Tensor4 &irradiance_field, Tensor5 &spectral_irradiance_field, const Index ybatch_index, const Agenda &input_agenda)
Agenda ybatch_calc_agenda(Workspace &ws) noexcept
Index ybatch_start(Workspace &ws) noexcept
ArrayOfVector ybatch(Workspace &ws) noexcept
A constant view of a Vector.
Index nelem() const
Number of elements.
Matrix sensor_los(Workspace &ws) noexcept
void ybatchCalc(Workspace &ws, ArrayOfVector &ybatch, ArrayOfArrayOfVector &ybatch_aux, ArrayOfMatrix &ybatch_jacobians, const Index &ybatch_start, const Index &ybatch_n, const Agenda &ybatch_calc_agenda, const Index &robust, const Verbosity &verbosity)
WORKSPACE METHOD: ybatchCalc.
Matrix z_surface(Workspace &ws) noexcept
The global header file for ARTS.
Index cloudbox_on(Workspace &ws) noexcept
This file contains basic functions to handle XML data files.
void met_profile_calc_agendaExecute(Workspace &ws, Vector &y, const GriddedField3 &t_field_raw, const ArrayOfGriddedField3 &vmr_field_raw, const GriddedField3 &z_field_raw, const ArrayOfGriddedField3 &pnd_field_raw, const Vector &p_grid, const Matrix &sensor_los, const Index cloudbox_on, const ArrayOfIndex &cloudbox_limits, const Matrix &z_surface, const Agenda &input_agenda)