ARTS 2.5.11 (git: 6827797f)
m_batch.cc
Go to the documentation of this file.
1/*===========================================================================
2 === File description
3 ===========================================================================*/
4
16/*===========================================================================
17 === External declarations
18 ===========================================================================*/
19
20#include <cmath>
21#include "gridded_fields.h"
22using namespace std;
23
24#include "arts.h"
25#include "arts_omp.h"
26#include "auto_md.h"
27#include "math_funcs.h"
28#include "physics_funcs.h"
29#include "rte.h"
30#include "xml_io.h"
31
32inline constexpr Numeric PI=Constant::pi;
33inline constexpr Numeric DEG2RAD=Conversion::deg2rad(1);
34inline constexpr Numeric RAD2DEG=Conversion::rad2deg(1);
36
37/*===========================================================================
38 === The functions (in alphabetical order)
39 ===========================================================================*/
40
41/* Workspace method: Doxygen documentation will be auto-generated
42
43 2008-07-21 Stefan Buehler */
45 // WS Input:
46 const Agenda& forloop_agenda,
47 // Control Parameters:
48 const Index& start,
49 const Index& stop,
50 const Index& step,
51 const Verbosity& verbosity) {
53
54 for (Index i = start; i <= stop; i += step) {
55 out1 << " Executing for loop body, index: " << i << "\n";
56 forloop_agendaExecute(ws, i, forloop_agenda);
57 }
58}
59
60/* Workspace method: Doxygen documentation will be auto-generated */
62 // WS Output:
63 ArrayOfVector& ybatch,
64 ArrayOfArrayOfVector& ybatch_aux,
65 ArrayOfMatrix& ybatch_jacobians,
66 // WS Input:
67 const Index& ybatch_start,
68 const Index& ybatch_n,
69 const Agenda& ybatch_calc_agenda,
70 // Control Parameters:
71 const Index& robust,
72 const Verbosity& verbosity) {
74
75 Index first_ybatch_index = 0;
76
77 ArrayOfString fail_msg;
78 bool do_abort = false;
79
80 // We allow a start index ybatch_start that is different from 0. We
81 // will calculate ybatch_n jobs starting at the start
82 // index. Internally, we count from zero, which is the right
83 // counting for the output array ybatch. When we call
84 // ybatch_calc_agenda, we add ybatch_start to the internal index
85 // count.
86
87 // We create a counter, so that we can generate nice output about
88 // how many jobs are already done. (All parallel threads later will
89 // increment this, so that we really get an accurate total count!)
90 Index job_counter = 0;
91
92 // Resize the output arrays:
93 ybatch.resize(ybatch_n);
94 ybatch_aux.resize(ybatch_n);
95 ybatch_jacobians.resize(ybatch_n);
96 for (Index i = 0; i < ybatch_n; i++) {
97 ybatch[i].resize(0);
98 ybatch_aux[i].resize(0);
99 ybatch_jacobians[i].resize(0, 0);
100 }
101
102 // Go through the batch:
103
104 if (ybatch_n) {
106
107#pragma omp parallel for schedule(dynamic) if (!arts_omp_in_parallel() && \
108 ybatch_n > 1) firstprivate(wss)
109 for (Index ybatch_index = first_ybatch_index; ybatch_index < ybatch_n;
110 ybatch_index++) {
111 Index l_job_counter; // Thread-local copy of job counter.
112
113 if (do_abort) continue;
114#pragma omp critical(ybatchCalc_job_counter)
115 { l_job_counter = ++job_counter; }
116
117 {
118 ostringstream os;
119 os << " Job " << l_job_counter << " of " << ybatch_n << ", Index "
120 << ybatch_start + ybatch_index << ", Thread-Id "
121 << arts_omp_get_thread_num() << "\n";
122 out2 << os.str();
123 }
124
125 try {
126 Vector y;
127 ArrayOfVector y_aux;
128 Matrix jacobian;
129
131 y,
132 y_aux,
133 jacobian,
134 ybatch_start + ybatch_index,
135 ybatch_calc_agenda);
136
137 if (y.nelem()) {
138#pragma omp critical(ybatchCalc_assign_y)
139 ybatch[ybatch_index] = y;
140#pragma omp critical(ybatchCalc_assign_y_aux)
141 ybatch_aux[ybatch_index] = y_aux;
142
143 // Dimensions of Jacobian:
144 const Index Knr = jacobian.nrows();
145 const Index Knc = jacobian.ncols();
146
147 if (Knr != 0 || Knc != 0) {
148 if (Knr != y.nelem()) {
149 ostringstream os;
150 os << "First dimension of Jacobian must have same length as the measurement *y*.\n"
151 << "Length of *y*: " << y.nelem() << "\n"
152 << "Dimensions of *jacobian*: (" << Knr << ", " << Knc
153 << ")\n";
154 // A mismatch of the Jacobian dimension is a fatal error
155 // and should result in program termination. By setting abort
156 // to true, this will result in a runtime error in the catch
157 // block even if robust == 1
158#pragma omp critical(ybatchCalc_setabort)
159 do_abort = true;
160
161 throw runtime_error(os.str());
162 }
163
164 ybatch_jacobians[ybatch_index] = jacobian;
165
166 // After creation, all individual Jacobi matrices in the array will be
167 // empty (size zero). No need for explicit initialization.
168 }
169 }
170 } catch (const std::exception& e) {
171 if (robust && !do_abort) {
172 ostringstream os;
173 os << "WARNING! Job at ybatch_index " << ybatch_start + ybatch_index
174 << " failed.\n"
175 << "y Vector in output variable ybatch will be empty for this job.\n"
176 << "The runtime error produced was:\n"
177 << e.what() << "\n";
178 out0 << os.str();
179 } else {
180 // The user wants the batch job to fail if one of the
181 // jobs goes wrong.
182#pragma omp critical(ybatchCalc_setabort)
183 do_abort = true;
184
185 ostringstream os;
186 os << " Job at ybatch_index " << ybatch_start + ybatch_index
187 << " failed. Aborting...\n";
188 out1 << os.str();
189 }
190 ostringstream os;
191 os << "Run-time error at ybatch_index " << ybatch_start + ybatch_index
192 << ": \n"
193 << e.what();
194#pragma omp critical(ybatchCalc_push_fail_msg)
195 fail_msg.push_back(os.str());
196 }
197 }
198 }
199
200 if (fail_msg.nelem()) {
201 ostringstream os;
202
203 if (!do_abort) os << "\nError messages from failed batch cases:\n";
204 for (ArrayOfString::const_iterator it = fail_msg.begin();
205 it != fail_msg.end();
206 it++)
207 os << *it << '\n';
208
209 if (do_abort)
210 throw runtime_error(os.str());
211 else
212 out0 << os.str();
213 }
214}
215
216/* Workspace method: Doxygen documentation will be auto-generated */
218 //Output
219 ArrayOfVector& ybatch,
220 //Input
221 const ArrayOfArrayOfSpeciesTag& abs_species,
222 const Agenda& met_profile_calc_agenda,
223 const Vector& f_grid,
224 const Matrix& met_amsu_data,
225 const Matrix& sensor_pos,
226 const Vector& refellipsoid,
227 const Vector& lat_grid,
228 const Vector& lon_grid,
229 const Index& atmosphere_dim,
230 const ArrayOfArrayOfSingleScatteringData& scat_data,
231 //Keyword
232 const Index& nelem_p_grid,
233 const String& met_profile_path,
234 const String& met_profile_pnd_path,
235 const Verbosity& verbosity) {
236 GriddedField3 t_field_raw;
237 GriddedField3 z_field_raw;
238 ArrayOfGriddedField3 vmr_field_raw;
239 ArrayOfGriddedField3 pnd_field_raw;
240 Vector p_grid;
241 Matrix sensor_los;
242 Index cloudbox_on;
243 ArrayOfIndex cloudbox_limits;
244 Matrix z_surface;
245 Vector y;
246 Index no_profiles = met_amsu_data.nrows();
247
248 // *vmr_field_raw* is an ArrayOfArrayOfTensor3 where the first array
249 //holds the gaseous species.
250 //Resize *vmr_field_raw* according to the number of gaseous species
251 //elements
252 vmr_field_raw.resize(abs_species.nelem());
253
254 //pnd_field_raw is an ArrayOfArrayOfTensor3 where the first array
255 //holds the scattering elements.
256 // Number of scattering elements:
257 const Index N_se = TotalNumberOfElements(scat_data);
258
259 pnd_field_raw.resize(N_se);
260
261 // The satellite zenith angle is read in from the amsu data
262 // and converted to arts sensor_los
263 ConstVectorView sat_za_from_data = met_amsu_data(Range(joker), 3);
264
265 sensor_los.resize(1, 1);
266
267 // The lat and lon are extracted to get the proper file names of
268 // profiles
269 ConstVectorView lat = met_amsu_data(Range(joker), 0);
270 ConstVectorView lon = met_amsu_data(Range(joker), 1);
271
272 z_surface.resize(1, 1);
273
274 // The spectra .
275 y.resize(f_grid.nelem());
276
277 // The batch spectra.
278 ybatch.resize(no_profiles);
279
280 // Loop over the number of profiles.
281 for (Index i = 0; i < no_profiles; ++i) {
282 ostringstream lat_os, lon_os;
283
284 Index lat_prec = 3;
285 if (lat[i] < 0) lat_prec--;
286 if (abs(lat[i]) >= 10) {
287 lat_prec--;
288 if (abs(lat[i]) >= 100) lat_prec--;
289 }
290
291 lat_os.setf(ios::showpoint | ios::fixed);
292 lat_os << setprecision((int)lat_prec) << lat[i];
293
294 Index lon_prec = 4;
295 if (lon[i] < 0) lon_prec--;
296 if (abs(lon[i]) >= 10) {
297 lon_prec--;
298 if (abs(lon[i]) >= 100) lon_prec--;
299 }
300 lon_os.setf(ios::showpoint | ios::fixed);
301 lon_os << setprecision((int)lon_prec) << lon[i];
302
303 sensor_los(0, 0) =
304 180.0 - (asin(refellipsoid[0] * sin(sat_za_from_data[i] * DEG2RAD) /
305 sensor_pos(0, 0))) *
306 RAD2DEG;
307
308 //Reads the t_field_raw from file
309 xml_read_from_file(met_profile_path + "profile.lat_" + lat_os.str() +
310 ".lon_" + lon_os.str() + ".t.xml",
311 t_field_raw,
312 verbosity);
313
314 //Reads the z_field_raw from file
315 xml_read_from_file(met_profile_path + "profile.lat_" + lat_os.str() +
316 ".lon_" + lon_os.str() + ".z.xml",
317 z_field_raw,
318 verbosity);
319
320 //Reads the humidity from file - it is only an ArrayofTensor3
321 xml_read_from_file(met_profile_path + "profile.lat_" + lat_os.str() +
322 ".lon_" + lon_os.str() + ".H2O.xml",
323 vmr_field_raw[0],
324 verbosity);
325
326 //Reads the pnd_field_raw for one scattering element
327 //xml_read_from_file("/rinax/storage/users/rekha/uk_data/profiles/new_obs/newest_forecastfields/reff100/profile.lat_"+lat_os.str()+".lon_"+lon_os.str() + ".pnd100.xml", pnd_field_raw[0]);
328
329 //xml_read_from_file(met_profile_pnd_path +"reff100_newformat/profile.lat_"+lat_os.str()+".lon_"+lon_os.str() + ".pnd100.xml", pnd_field_raw[0]);
330
331 xml_read_from_file(met_profile_pnd_path + "lwc_reff15/profile.lat_" +
332 lat_os.str() + ".lon_" + lon_os.str() + ".pnd15.xml",
333 pnd_field_raw[0],
334 verbosity);
335 //Write the profile number into a file.
336 // xml_write_to_file("profile_number.xml", i);
337
338 // Set z_surface from lowest level of z_field
339 z_surface(0, 0) = z_field_raw.data(0, 0, 0);
340
341 /* The vmr_field_raw is an ArrayofArrayofTensor3 where the outer
342 array is for species.
343
344 The oxygen and nitrogen VMRs are set to constant values of 0.209
345 and 0.782, respectively and are used along with humidity field
346 to generate *vmr_field_raw*.*/
347
348 /*The second element of the species. The first 3 Tensors in the
349 array are the same . They are pressure grid, latitude grid and
350 longitude grid. The third tensor which is the vmr is set to a
351 constant value of 0.782, corresponding to N2.*/
352
353 vmr_field_raw[1].resize(vmr_field_raw[0]);
354 vmr_field_raw[1].copy_grids(vmr_field_raw[0]);
355 vmr_field_raw[1] = 0.782; //vmr of N2
356
357 /*the third element of the species. the first 3 Tensors in the
358 array are the same . They are pressure grid, latitude grid and
359 longitude grid. The third tensor which is the vmr is set to a
360 constant value of 0.209, corresponding to O2.*/
361 vmr_field_raw[2].resize(vmr_field_raw[0]);
362 vmr_field_raw[2].copy_grids(vmr_field_raw[0]);
363 vmr_field_raw[2] = 0.209; //vmr of O2
364
365 const Vector& tfr_p_grid =
366 t_field_raw.get_numeric_grid(GFIELD3_P_GRID);
367 // N_p is the number of elements in the pressure grid
368 Index N_p = tfr_p_grid.nelem();
369
370 //Making a p_grid with the first and last element taken from the profile.
372 p_grid, nelem_p_grid, tfr_p_grid[0], tfr_p_grid[N_p - 1], verbosity);
373
374 /*To set the cloudbox limits, the lower and upper cloudbox limits
375 are to be set. The lower cloudbox limit is set to the lowest
376 pressure level. The upper level is the highest level where the
377 ice water content is non-zero.*/
378 Numeric cl_grid_min, cl_grid_max;
379
380 //Lower limit = lowest pressure level of the original grid.
381 //Could it be the interpolated p_grid? FIXME STR
382 cl_grid_min = tfr_p_grid[0];
383
384 // A counter for non-zero ice content
385 Index level_counter = 0;
386
387 // Loop over all pressure levels
388 for (Index ip = 0; ip < N_p; ++ip) {
389 //Checking for non-zero ice content. 0.001 is a threshold for
390 //ice water content.
391 // if((pnd_field_raw[0].data()(ip, 0, 0) > 0.001) || (pnd_field_raw[1](ip, 0, 0) > 0.001))
392 if (pnd_field_raw[0].data(ip, 0, 0) > 0.001) {
393 ++level_counter;
394 //if non-zero ice content is found, it is set to upper
395 //cloudbox limit. Moreover, we take one level higher
396 // than the upper limit because we want the upper limit
397 //to have 0 pnd.
398 cl_grid_max = tfr_p_grid[ip + 1];
399 }
400 }
401
402 //cloudbox limits have dimensions 2*atmosphere_dim
403 cloudbox_limits.resize(atmosphere_dim * 2);
404
405 //if there is no cloud in the considered profile, still we
406 //need to set the upper limit. I here set the first level
407 //for the upper cloudbox limit.
408 if (level_counter == 0) {
409 cl_grid_max = p_grid[1];
410 }
411
412 //Cloudbox is set.
413 cloudboxSetManually(cloudbox_on,
414 cloudbox_limits,
415 atmosphere_dim,
416 p_grid,
417 lat_grid,
418 lon_grid,
419 cl_grid_min,
420 cl_grid_max,
421 0,
422 0,
423 0,
424 0,
425 verbosity);
426
427 /*executing the met_profile_calc_agenda
428 Agenda communication variables are
429 Output of met_profile_calc_agenda : y
430 Input to met_profile_calc_agenda : t_field_raw,
431 z_field_raw, vmr_field_raw, pnd_field_raw, p_grid,
432 sensor_los, cloudbox_on, cloudbox_limits, z_surface, */
433
435 y,
436 t_field_raw,
437 vmr_field_raw,
438 z_field_raw,
439 pnd_field_raw,
440 p_grid,
441 sensor_los,
442 cloudbox_on,
443 cloudbox_limits,
444 z_surface,
445 met_profile_calc_agenda);
446
447 //putting in the spectra *y* for each profile, thus assigning y
448 //to the ith row of ybatch
449 ybatch[i] = y;
450
451 } // closing the loop over profile basenames
452}
453
454/* Workspace method: Doxygen documentation will be auto-generated */
456 //Output
457 ArrayOfVector& ybatch,
458 //Input
459 const ArrayOfArrayOfSpeciesTag& abs_species,
460 const Agenda& met_profile_calc_agenda,
461 const Vector& f_grid,
462 const Matrix& met_amsu_data,
463 const Matrix& sensor_pos,
464 const Vector& refellipsoid,
465 //Keyword
466 const Index& nelem_p_grid,
467 const String& met_profile_path,
468 const Verbosity& verbosity) {
469 GriddedField3 t_field_raw;
470 GriddedField3 z_field_raw;
471 ArrayOfGriddedField3 vmr_field_raw;
472 ArrayOfGriddedField3 pnd_field_raw;
473 Vector p_grid;
474 Matrix sensor_los;
475 Index cloudbox_on = 0;
476 ArrayOfIndex cloudbox_limits;
477 Matrix z_surface;
478 Vector y;
479 Index no_profiles = met_amsu_data.nrows();
480 //Index no_profiles = met_profile_basenames.nelem();
481 // The humidity data is stored as an ArrayOfTensor3 whereas
482 // vmr_field_raw is an ArrayOfArrayOfTensor3
483 GriddedField3 vmr_field_raw_h2o;
484
485 vmr_field_raw.resize(abs_species.nelem());
486
487 y.resize(f_grid.nelem());
488 ybatch.resize(no_profiles);
489
490 Vector sat_za_from_profile;
491 sat_za_from_profile = met_amsu_data(Range(joker), 3);
492 Numeric sat_za;
493
494 sensor_los.resize(1, 1);
495
496 Vector lat, lon;
497 lat = met_amsu_data(Range(joker), 0);
498 lon = met_amsu_data(Range(joker), 1);
499
500// Vector oro_height;
501// oro_height = met_amsu_data(Range(joker), 5);
502
503 z_surface.resize(1, 1);
504 for (Index i = 0; i < no_profiles; ++i) {
505 ostringstream lat_os, lon_os;
506
507 Index lat_prec = 3;
508 if (lat[i] < 0) lat_prec--;
509 if (abs(lat[i]) >= 10) {
510 lat_prec--;
511 if (abs(lat[i]) >= 100) lat_prec--;
512 }
513
514 lat_os.setf(ios::showpoint | ios::fixed);
515 lat_os << setprecision((int)lat_prec) << lat[i];
516
517 Index lon_prec = 4;
518 if (lon[i] < 0) lon_prec--;
519 if (abs(lon[i]) >= 10) {
520 lon_prec--;
521 if (abs(lon[i]) >= 100) lon_prec--;
522 }
523 lon_os.setf(ios::showpoint | ios::fixed);
524 lon_os << setprecision((int)lon_prec) << lon[i];
525 cout << lat_os.str() << endl;
526 cout << lon_os.str() << endl;
527
528 sat_za = sat_za_from_profile[i];
529
530 sensor_los(Range(joker), 0) =
531 180.0 -
532 (asin(refellipsoid[0] * sin(sat_za * PI / 180.) / sensor_pos(0, 0))) *
533 180. / PI;
534 cout << "sensor_los" << sat_za_from_profile[i] << endl;
535 cout << "sensor_los" << sat_za << endl;
536 cout << "sensor_los" << sensor_los << endl;
537 //Reads the t_field_raw from file
538
539 xml_read_from_file(met_profile_path + "profile.lat_" + lat_os.str() +
540 ".lon_" + lon_os.str() + ".t.xml",
541 t_field_raw,
542 verbosity);
543 //Reads the z_field_raw from file
544 xml_read_from_file(met_profile_path + "profile.lat_" + lat_os.str() +
545 ".lon_" + lon_os.str() + ".z.xml",
546 z_field_raw,
547 verbosity);
548
549 //Reads the humidity from file - it is only an ArrayofTensor3
550 // The vmr_field_raw is an ArrayofArrayofTensor3 where the outer
551 // array is for species
552 xml_read_from_file(met_profile_path + "profile.lat_" + lat_os.str() +
553 ".lon_" + lon_os.str() + ".H2O.xml",
554 vmr_field_raw_h2o,
555 verbosity);
556 //xml_read_from_file("/home/home01/rekha/uk/profiles/sat_vmr/profile.lat_"+lat_os.str()//+".lon_"+lon_os.str() + ".H2O_es.xml",
557 // vmr_field_raw_h2o, verbosity);
558
559 cout
560 << "--------------------------------------------------------------------------"
561 << endl;
562 cout << "The file"
563 << met_profile_path + "profile.lat_" + lat_os.str() + ".lon_" +
564 lon_os.str()
565 << "is executed now" << endl;
566 cout
567 << "--------------------------------------------------------------------------"
568 << endl;
569 xml_write_to_file("profile_number.xml", i, FILE_TYPE_ASCII, 0, verbosity);
570 // the first element of the species is water vapour.
571
572 // N_p is the number of elements in the pressure grid
573 //z_surface(0,0) = oro_height[i]+ 0.01;
574 z_surface(0, 0) = z_field_raw.data(0, 0, 0);
575 cout << "z_surface" << z_surface << endl;
576 const Vector& tfr_p_grid =
577 t_field_raw.get_numeric_grid(GFIELD3_P_GRID);
578 Index N_p = tfr_p_grid.nelem();
579
580 vmr_field_raw[0] = vmr_field_raw_h2o;
581
582 // the second element of the species. the first 3 Tensors in the
583 //array are the same . They are pressure grid, latitude grid and
584 // longitude grid. The third tensor which is the vmr is set to a
585 // constant value of 0.782.
586 vmr_field_raw[1].resize(vmr_field_raw[0]);
587 vmr_field_raw[1].copy_grids(vmr_field_raw[0]);
588 vmr_field_raw[1].data(joker, joker, joker) = 0.782;
589
590 // the second element of the species. the first 3 Tensors in the
591 //array are the same . They are pressure grid, latitude grid and
592 // longitude grid. The third tensor which is the vmr is set to a
593 // constant value of 0.209.
594 vmr_field_raw[2].resize(vmr_field_raw[0]);
595 vmr_field_raw[2].copy_grids(vmr_field_raw[0]);
596 vmr_field_raw[2].data(joker, joker, joker) = 0.209;
597
598 //xml_write_to_file(met_profile_basenames[i]+ ".N2.xml", vmr_field_raw[1]);
599 //xml_write_to_file(met_profile_basenames[i]+ ".O2.xml", vmr_field_raw[2]);
600
601 //Making a p_grid with the first and last element taken from the profile.
602 // this is because of the extrapolation problem.
603
605 p_grid, nelem_p_grid, tfr_p_grid[0], tfr_p_grid[N_p - 1], verbosity);
606 cout << "t_field_raw[0](0,0,0)" << tfr_p_grid[0] << endl;
607 cout << "t_field_raw[0](N_p -1,0,0)" << tfr_p_grid[N_p - 1] << endl;
608 xml_write_to_file("p_grid.xml", p_grid, FILE_TYPE_ASCII, 0, verbosity);
609
610 // executing the met_profile_calc_agenda
612 y,
613 t_field_raw,
614 vmr_field_raw,
615 z_field_raw,
616 pnd_field_raw,
617 p_grid,
618 sensor_los,
619 cloudbox_on,
620 cloudbox_limits,
621 z_surface,
622 met_profile_calc_agenda);
623
624 //putting in the spectra *y* for each profile
625 ybatch[i] = y;
626
627 } // closing the loop over profile basenames
628}
629
630/* Workspace method: Doxygen documentation will be auto-generated */
632 ArrayOfTensor7& dobatch_cloudbox_field,
633 ArrayOfTensor5& dobatch_radiance_field,
634 ArrayOfTensor4& dobatch_irradiance_field,
635 ArrayOfTensor5& dobatch_spectral_irradiance_field,
636 const Index& ybatch_start,
637 const Index& ybatch_n,
638 const Agenda& dobatch_calc_agenda,
639 const Index& robust,
640 const Verbosity& verbosity) {
642
643 Index first_ybatch_index = 0;
644
645 ArrayOfString fail_msg;
646 bool do_abort = false;
647
648 // We allow a start index ybatch_start that is different from 0. We
649 // will calculate ybatch_n jobs starting at the start
650 // index. Internally, we count from zero, which is the right
651 // counting for the output array ybatch. When we call
652 // ybatch_calc_agenda, we add ybatch_start to the internal index
653 // count.
654
655 // We create a counter, so that we can generate nice output about
656 // how many jobs are already done. (All parallel threads later will
657 // increment this, so that we really get an accurate total count!)
658 Index job_counter = 0;
659
660 // Resize the output arrays:
661 dobatch_cloudbox_field.resize(ybatch_n);
662 dobatch_radiance_field.resize(ybatch_n);
663 dobatch_irradiance_field.resize(ybatch_n);
664 dobatch_spectral_irradiance_field.resize(ybatch_n);
665
666 // Go through the batch:
667
668 if (ybatch_n) {
670
671#pragma omp parallel for schedule(dynamic) if (!arts_omp_in_parallel() && \
672 ybatch_n > 1) firstprivate(wss)
673 for (Index ybatch_index = first_ybatch_index; ybatch_index < ybatch_n;
674 ybatch_index++) {
675 Index l_job_counter; // Thread-local copy of job counter.
676
677 if (do_abort) continue;
678#pragma omp critical(dobatchCalc_job_counter)
679 { l_job_counter = ++job_counter; }
680
681 {
682 ostringstream os;
683 os << " Job " << l_job_counter << " of " << ybatch_n << ", Index "
684 << ybatch_start + ybatch_index << ", Thread-Id "
685 << arts_omp_get_thread_num() << "\n";
686 out2 << os.str();
687 }
688
689 try {
690 Tensor7 cloudbox_field;
691 Tensor5 radiance_field;
692 Tensor4 irradiance_field;
693 Tensor5 spectral_irradiance_field;
694
696 cloudbox_field,
697 radiance_field,
698 irradiance_field,
699 spectral_irradiance_field,
700 ybatch_start + ybatch_index,
701 dobatch_calc_agenda);
702
703#pragma omp critical(dobatchCalc_assign_cloudbox_field)
704 dobatch_cloudbox_field[ybatch_index] = cloudbox_field;
705#pragma omp critical(dobatchCalc_assign_radiance_field)
706 dobatch_radiance_field[ybatch_index] = radiance_field;
707#pragma omp critical(dobatchCalc_assign_irradiance_field)
708 dobatch_irradiance_field[ybatch_index] = irradiance_field;
709#pragma omp critical(dobatchCalc_assign_spectral_irradiance_field)
710 dobatch_spectral_irradiance_field[ybatch_index] =
711 spectral_irradiance_field;
712
713 } catch (const std::exception& e) {
714 if (robust && !do_abort) {
715 ostringstream os;
716 os << "WARNING! Job at ybatch_index " << ybatch_start + ybatch_index
717 << " failed.\n"
718 << "element in output variables will be empty for this job.\n"
719 << "The runtime error produced was:\n"
720 << e.what() << "\n";
721 out0 << os.str();
722 } else {
723 // The user wants the batch job to fail if one of the
724 // jobs goes wrong.
725#pragma omp critical(dobatchCalc_setabort)
726 do_abort = true;
727
728 ostringstream os;
729 os << " Job at ybatch_index " << ybatch_start + ybatch_index
730 << " failed. Aborting...\n";
731 out1 << os.str();
732 }
733 ostringstream os;
734 os << "Run-time error at ybatch_index " << ybatch_start + ybatch_index
735 << ": \n"
736 << e.what();
737#pragma omp critical(dobatchCalc_push_fail_msg)
738 fail_msg.push_back(os.str());
739 }
740 }
741 }
742
743 if (fail_msg.nelem()) {
744 ostringstream os;
745
746 if (!do_abort) os << "\nError messages from failed batch cases:\n";
747 for (ArrayOfString::const_iterator it = fail_msg.begin();
748 it != fail_msg.end();
749 it++)
750 os << *it << '\n';
751
752 if (do_abort)
753 throw runtime_error(os.str());
754 else
755 out0 << os.str();
756 }
757}
Index TotalNumberOfElements(const Array< Array< base > > &aa)
Determine total number of elements in an ArrayOfArray.
Definition: array.h:224
The global header file for ARTS.
int arts_omp_get_thread_num()
Wrapper for omp_get_thread_num.
Definition: arts_omp.cc:59
Header file for helper functions for OpenMP.
void forloop_agendaExecute(Workspace &ws, const Index forloop_index, const Agenda &input_agenda)
Definition: auto_md.cc:24987
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)
Definition: auto_md.cc:24815
void ybatch_calc_agendaExecute(Workspace &ws, Vector &y, ArrayOfVector &y_aux, Matrix &jacobian, const Index ybatch_index, const Agenda &input_agenda)
Definition: auto_md.cc:26161
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)
Definition: auto_md.cc:25570
The Agenda class.
Definition: agenda_class.h:52
Index nelem() const ARTS_NOEXCEPT
Definition: array.h:75
const Vector & get_numeric_grid(Index i) const
Get a numeric grid.
Workspace class.
Definition: workspace_ng.h:36
Implementation of gridded fields.
void VectorNLogSpace(Vector &x, const Index &n, const Numeric &start, const Numeric &stop, const Verbosity &verbosity)
WORKSPACE METHOD: VectorNLogSpace.
constexpr Numeric DEG2RAD
Definition: m_batch.cc:33
constexpr Numeric RAD2DEG
Definition: m_batch.cc:34
void ForLoop(Workspace &ws, const Agenda &forloop_agenda, const Index &start, const Index &stop, const Index &step, const Verbosity &verbosity)
WORKSPACE METHOD: ForLoop.
Definition: m_batch.cc:44
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.
Definition: m_batch.cc:455
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.
Definition: m_batch.cc:61
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.
Definition: m_batch.cc:217
constexpr Numeric PI
Definition: m_batch.cc:32
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.
Definition: m_batch.cc:631
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.
Definition: m_cloudbox.cc:323
#define CREATE_OUT1
Definition: messages.h:187
#define CREATE_OUTS
Definition: messages.h:191
constexpr Numeric pi
The following mathematical constants are generated in python Decimal package by the code:
constexpr auto deg2rad(auto x) noexcept
Converts degrees to radians.
constexpr auto rad2deg(auto x) noexcept
Converts radians to degrees.
constexpr Index GFIELD3_P_GRID
Global constant, Index of the pressure grid in GriddedField3.
Definition: mystring.h:246
This file contains declerations of functions of physical character.
Declaration of functions in rte.cc.
This file contains basic functions to handle XML data files.
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.
Definition: xml_io.h:155
void xml_read_from_file(const String &filename, T &type, const Verbosity &verbosity)
Reads data from XML file.
Definition: xml_io.h:134
@ FILE_TYPE_ASCII
Definition: xml_io_base.h:26