FlowSieve  3.0.0
Coarse Graining Routines
/home/docs/checkouts/readthedocs.org/user_builds/flowsieve/checkouts/read_the_docs_testing/functions.hpp
Go to the documentation of this file.
1 #ifndef FUNCTIONS_HPP
2 #define FUNCTIONS_HPP 1
3 
4 #include <stdio.h>
5 #include <stdlib.h>
6 #include <vector>
7 #include <string>
8 #include <map>
9 #include <mpi.h>
10 #include "constants.hpp"
11 
24 class dataset {
25 
26  public:
27 
28  // Storage for processor assignments
29  int Nprocs_in_time, Nprocs_in_depth;
30 
31  // Vectors to store the dimension variables
32  std::vector<double> time, depth, latitude, longitude;
33  int Ntime = -1, Ndepth = -1, Nlat = -1, Nlon = -1;
34  int full_Ntime = -1, full_Ndepth = -1;
35 
36  // Store cell areas
37  std::vector<double> areas;
38 
39  // Dictionary for the variables (velocity components, density, etc)
40  std::map< std::string , std::vector<double> > variables;
41 
42  // Vectors and dictionary for storing the regions.
43  // these are what will be used for the post-processing
44  std::vector< std::string > region_names;
45  std::map< std::string, std::vector<bool> > regions;
46  std::vector<double> region_areas, region_areas_water_only;
47 
48  // Store mask data (i.e. land vs water)
49  std::vector<bool> mask, reference_mask;
50 
51  // Store data-chunking info. These keep track of the MPI divisions to ensure
52  // that the output is in the same order as the input.
53  std::vector<int> myCounts, myStarts;
54 
55  // Constructor
56  dataset();
57 
58  // Dimension loaders
59  void load_time( const std::string dim_name, const std::string filename );
60  void load_depth( const std::string dim_name, const std::string filename );
61  void load_latitude( const std::string dim_name, const std::string filename );
62  void load_longitude( const std::string dim_name, const std::string filename );
63 
64  // Compute areas
65  void compute_cell_areas();
66 
67  // Load in variable and store in dictionary
68  void load_variable( const std::string var_name,
69  const std::string var_name_in_file,
70  const std::string filename,
71  const bool read_mask = true,
72  const bool load_counts = true,
73  const bool do_splits = true );
74 
75  // Load in region definitions
76  void load_region_definitions( const std::string filename,
77  const std::string dim_name,
78  const std::string var_name,
79  const MPI_Comm = MPI_COMM_WORLD );
80  void compute_region_areas();
81 
82  // Check the processors divions between dimensions
83  void check_processor_divisions( const int Nprocs_in_time_input, const int Nprocs_in_depth_input, const MPI_Comm = MPI_COMM_WORLD );
84 
85 };
86 
87 void compute_areas(
88  std::vector<double> & areas,
89  const std::vector<double> & longitude,
90  const std::vector<double> & latitude);
91 
92 size_t Index( const int Itime, const int Idepth, const int Ilat, const int Ilon,
93  const int Ntime, const int Ndepth, const int Nlat, const int Nlon );
94 
95 void Index1to4( const size_t index,
96  int & Itime, int & Idepth, int & Ilat, int & Ilon,
97  const int & Ntime, const int & Ndepth, const int & Nlat, const int & Nlon );
98 
99 double distance(const double lon1, const double lat1,
100  const double lon2, const double lat2,
101  const double Llon = 0, const double Llat = 0);
102 
104  std::vector<double> & local_kernel,
105  const double scale,
106  const dataset & source_data,
107  const int Ilat, const int Ilon,
108  const int LAT_lb, const int LAT_ub);
109 
110 void KE_from_vels(
111  std::vector<double> & KE,
112  std::vector<double> * u1,
113  std::vector<double> * u2,
114  std::vector<double> * u3,
115  const std::vector<bool> & mask,
116  const double rho0 = constants::rho0
117  );
118 
119 void vel_Spher_to_Cart(
120  std::vector<double> & u_x,
121  std::vector<double> & u_y,
122  std::vector<double> & u_z,
123  const std::vector<double> & u_r,
124  const std::vector<double> & u_lon,
125  const std::vector<double> & u_lat,
126  const dataset & source_data
127  );
128 
130  double & u_x,
131  double & u_y,
132  double & u_z,
133  const double u_r,
134  const double u_lon,
135  const double u_lat,
136  const double lon,
137  const double lat
138  );
139 
140 
141 void vel_Cart_to_Spher(
142  std::vector<double> & u_r,
143  std::vector<double> & u_lon,
144  std::vector<double> & u_lat,
145  const std::vector<double> & u_x,
146  const std::vector<double> & u_y,
147  const std::vector<double> & u_z,
148  const dataset & source_data
149  );
150 
151 
153  double & u_r,
154  double & u_lon,
155  double & u_lat,
156  const double u_x,
157  const double u_y,
158  const double u_z,
159  const double lon,
160  const double lat
161  );
162 
163 void filtering(const dataset & source_data,
164  const std::vector<double> & scales,
165  const MPI_Comm comm = MPI_COMM_WORLD);
166 
168  const dataset & source_data,
169  const std::vector<double> & scales,
170  const MPI_Comm comm = MPI_COMM_WORLD
171  );
172 
174  std::vector<double*> & coarse_val,
175  const std::vector<const std::vector<double>*> & fields,
176  const dataset & source_data,
177  const int Itime, const int Idepth, const int Ilat, const int Ilon,
178  const int LAT_lb,
179  const int LAT_ub,
180  const double scale,
181  const std::vector<bool> & use_mask,
182  const std::vector<double> & local_kernel,
183  const std::vector<double> * weight = NULL
184  );
185 
186 double kernel(const double distance, const double scale);
187 
188 double kernel_alpha(void);
189 
191  double & vort_r_tmp,
192  double & vort_lon_tmp,
193  double & vort_lat_tmp,
194  double & div_tmp,
195  double & OkuboWeiss_tmp,
196  const std::vector<double> & u_r,
197  const std::vector<double> & u_lon,
198  const std::vector<double> & u_lat,
199  const int Ntime, const int Ndepth, const int Nlat, const int Nlon,
200  const int Itime, const int Idepth, const int Ilat, const int Ilon,
201  const std::vector<double> & longitude,
202  const std::vector<double> & latitude,
203  const std::vector<bool> & mask);
204 
205 void compute_vorticity(
206  std::vector<double> & vort_r,
207  std::vector<double> & vort_lon,
208  std::vector<double> & vort_lat,
209  std::vector<double> & vel_div,
210  std::vector<double> & OkuboWeiss,
211  const std::vector<double> & u_r,
212  const std::vector<double> & u_lon,
213  const std::vector<double> & u_lat,
214  const int Ntime, const int Ndepth, const int Nlat, const int Nlon,
215  const std::vector<double> & longitude,
216  const std::vector<double> & latitude,
217  const std::vector<bool> & mask,
218  const MPI_Comm comm = MPI_COMM_WORLD);
219 
221  double & uxux_tmp, double & uxuy_tmp, double & uxuz_tmp,
222  double & uyuy_tmp, double & uyuz_tmp, double & uzuz_tmp,
223  double & vort_ux_tmp, double & vort_uy_tmp, double & vort_uz_tmp,
224  const std::vector<double> & u_x,
225  const std::vector<double> & u_y,
226  const std::vector<double> & u_z,
227  const std::vector<double> & vort_r,
228  const dataset & source_data,
229  const int Itime, const int Idepth, const int Ilat, const int Ilon,
230  const int LAT_lb, const int LAT_ub,
231  const double scale,
232  const std::vector<double> & local_kernel);
233 
234 void compute_Pi(
235  std::vector<double> & energy_transfer,
236  const dataset & source_data,
237  const std::vector<double> & ux,
238  const std::vector<double> & uy,
239  const std::vector<double> & uz,
240  const std::vector<double> & uxux,
241  const std::vector<double> & uxuy,
242  const std::vector<double> & uxuz,
243  const std::vector<double> & uyuy,
244  const std::vector<double> & uyuz,
245  const std::vector<double> & uzuz,
246  const MPI_Comm comm = MPI_COMM_WORLD);
247 
249  std::vector<double> & energy_transfer,
250  const dataset & source_data,
251  const std::vector<double> & ux,
252  const std::vector<double> & uy,
253  const std::vector<double> & uz,
254  const std::vector<double> & uxux,
255  const std::vector<double> & uxuy,
256  const std::vector<double> & uxuz,
257  const std::vector<double> & uyuy,
258  const std::vector<double> & uyuz,
259  const std::vector<double> & uzuz,
260  const MPI_Comm comm = MPI_COMM_WORLD);
261 
263  std::vector<double> & energy_transfer,
264  const dataset & source_data,
265  const std::vector<double> & ulon,
266  const std::vector<double> & ulat,
267  const std::vector<double> & ulon_ulon,
268  const std::vector<double> & ulon_ulat,
269  const std::vector<double> & ulat_ulat,
270  const MPI_Comm comm = MPI_COMM_WORLD
271  );
272 
273 void compute_Z(
274  std::vector<double> & enstrophy_transfer,
275  const dataset & source_data,
276  const std::vector<double> & ux,
277  const std::vector<double> & uy,
278  const std::vector<double> & uz,
279  const std::vector<double> & coarse_vort_r,
280  const std::vector<double> & vort_ux,
281  const std::vector<double> & vort_uy,
282  const std::vector<double> & vort_uz,
283  const MPI_Comm comm = MPI_COMM_WORLD);
284 
286  std::vector<double> & Lambda_rot,
287  const std::vector<double> & coarse_vort_r,
288  const std::vector<double> & coarse_vort_lon,
289  const std::vector<double> & coarse_vort_lat,
290  const std::vector<double> & coarse_rho,
291  const std::vector<double> & coarse_p,
292  const int Ntime, const int Ndepth, const int Nlat, const int Nlon,
293  const std::vector<double> & longitude,
294  const std::vector<double> & latitude,
295  const std::vector<bool> & mask,
296  const double scale_factor
297  );
298 
299 
301  std::vector<double> & Lambda,
302  const std::vector<double> & coarse_u_r,
303  const std::vector<double> & coarse_u_lon,
304  const std::vector<double> & coarse_u_lat,
305  const std::vector<double> & tilde_u_r,
306  const std::vector<double> & tilde_u_lon,
307  const std::vector<double> & tilde_u_lat,
308  const std::vector<double> & coarse_p,
309  const int Ntime,
310  const int Ndepth,
311  const int Nlat,
312  const int Nlon,
313  const std::vector<double> & longitude,
314  const std::vector<double> & latitude,
315  const std::vector<bool> & mask
316  );
317 
319  std::vector<double> & Lambda_nonlin,
320  const std::vector<double> & coarse_u_r,
321  const std::vector<double> & coarse_u_lon,
322  const std::vector<double> & coarse_u_lat,
323  const std::vector<double> & coarse_rho,
324  const std::vector<double> & coarse_p,
325  const int Ntime, const int Ndepth, const int Nlat, const int Nlon,
326  const std::vector<double> & longitude,
327  const std::vector<double> & latitude,
328  const std::vector<bool> & mask,
329  const double scale_factor
330  );
331 
333  const double p,
334  const double theta);
335 
336 double equation_of_state(
337  const double T,
338  const double S,
339  const double p);
340 
342  std::vector<double> & div_J,
343  const std::vector<double> & u_x,
344  const std::vector<double> & u_y,
345  const std::vector<double> & u_z,
346  const std::vector<double> & uxux,
347  const std::vector<double> & uxuy,
348  const std::vector<double> & uxuz,
349  const std::vector<double> & uyuy,
350  const std::vector<double> & uyuz,
351  const std::vector<double> & uzuz,
352  const std::vector<double> & coarse_p,
353  const std::vector<double> & longitude,
354  const std::vector<double> & latitude,
355  const int Ntime,
356  const int Ndepth,
357  const int Nlat,
358  const int Nlon,
359  const std::vector<bool> & mask);
360 
362  std::vector<double> & means,
363  const std::vector<double> & field,
364  const std::vector<double> & areas,
365  const int Ntime,
366  const int Ndepth,
367  const int Nlat,
368  const int Nlon,
369  const std::vector<bool> & mask);
370 
371 void get_lat_bounds(
372  int & LAT_lb,
373  int & LAT_ub,
374  const std::vector<double> & latitude,
375  const int Ilat,
376  const double scale);
377 
378 void get_lon_bounds(
379  int & LON_lb,
380  int & LON_ub,
381  const std::vector<double> & longitude,
382  const int Ilon,
383  const double centre_lat,
384  const double curr_lat,
385  const double scale);
386 
387 void print_compile_info(
388  const std::vector<double> * scales = NULL);
389 
390 int get_omp_chunksize(const int Nlat, const int Nlon);
391 
392 
394  std::vector<double> & longitude,
395  std::vector<double> & latitude
396  );
397 
398 void mask_out_pole(
399  const std::vector<double> & latitude,
400  std::vector<bool> & mask,
401  const int Ntime,
402  const int Ndepth,
403  const int Nlat,
404  const int Nlon
405  );
406 
407 
408 void roll_field(
409  std::vector<double> & field_to_roll,
410  const std::string dimension,
411  const int roll_count,
412  const int Ntime,
413  const int Ndepth,
414  const int Nlat,
415  const int Nlon
416  );
417 
418 
419 // Extending to the poles
420 
421 void extend_latitude_to_poles(
422  const std::vector<double> & original_latitude,
423  std::vector<double> & extended_latitude,
424  int & orig_Ilat_start,
425  const bool IS_DEGREES = false,
426  const MPI_Comm comm = MPI_COMM_WORLD
427  );
428 
429 void extend_field_to_poles(
430  std::vector<double> & array_to_extend,
431  const dataset & source_data,
432  const std::vector<double> & extended_latitude,
433  const int Ilat_start
434  );
435 
436 void extend_mask_to_poles(
437  std::vector<bool> & mask_to_extend,
438  const dataset & source_data,
439  const std::vector<double> & extended_latitude,
440  const int Ilat_start,
441  const bool extend_val = constants::FILTER_OVER_LAND
442  );
443 
451 
452  public:
454  Timing_Records();
455 
457  void reset();
458 
464  void add_to_record( const double delta, const std::string record_name );
465 
467  void print() const;
468 
469  private:
475  std::map< std::string, double > time_records;
476 };
477 
478 
488 class InputParser {
489 
490  public:
491  InputParser (int &argc, char **argv);
492 
493  const std::string getCmdOption(
494  const std::string &option,
495  const std::string &default_value = ""
496  ) const;
497 
498  bool cmdOptionExists(const std::string &option) const;
499 
500  void getFilterScales( std::vector<double> &filter_scales, const std::string &argname ) const;
501 
502  void getListofStrings( std::vector<std::string> &list_of_strings, const std::string &argname ) const;
503 
504  private:
505  std::vector <std::string> tokens;
506 
507 };
508 
512 bool string_to_bool( std::string str );
513 
514 void print_header_info( );
515 
516 #endif
void compute_Lambda_nonlin_model(std::vector< double > &Lambda_nonlin, const std::vector< double > &coarse_u_r, const std::vector< double > &coarse_u_lon, const std::vector< double > &coarse_u_lat, const std::vector< double > &coarse_rho, const std::vector< double > &coarse_p, const int Ntime, const int Ndepth, const int Nlat, const int Nlon, const std::vector< double > &longitude, const std::vector< double > &latitude, const std::vector< bool > &mask, const double scale_factor)
Compute the non-linear model of the baroclinic transfer term Lambda (see Lees and Aluie 2019) ...
Definition: compute_Lambda_nonlin_model.cpp:26
void KE_from_vels(std::vector< double > &KE, std::vector< double > *u1, std::vector< double > *u2, std::vector< double > *u3, const std::vector< bool > &mask, const double rho0=constants::rho0)
Compute (local) KE density from the provided velocities.
Definition: KE_from_vels.cpp:16
void Index1to4(const size_t index, int &Itime, int &Idepth, int &Ilat, int &Ilon, const int &Ntime, const int &Ndepth, const int &Nlat, const int &Nlon)
Convenience tool to convert logical index to physical index (time, depth, lat, lon).
Definition: Index1to4.cpp:18
const double rho0
Mean fluid density.
Definition: constants.hpp:67
void vel_Cart_to_Spher_at_point(double &u_r, double &u_lon, double &u_lat, const double u_x, const double u_y, const double u_z, const double lon, const double lat)
Convert single Cartesian velocity to spherical velocity.
Definition: vel_Cart_to_Spher_at_point.cpp:28
void get_lon_bounds(int &LON_lb, int &LON_ub, const std::vector< double > &longitude, const int Ilon, const double centre_lat, const double curr_lat, const double scale)
Get longitude integratation bounds for a specific latitude.
Definition: get_lon_bounds.cpp:21
void compute_Pi_shift_deriv(std::vector< double > &energy_transfer, const dataset &source_data, const std::vector< double > &ux, const std::vector< double > &uy, const std::vector< double > &uz, const std::vector< double > &uxux, const std::vector< double > &uxuy, const std::vector< double > &uxuz, const std::vector< double > &uyuy, const std::vector< double > &uyuz, const std::vector< double > &uzuz, const MPI_Comm comm=MPI_COMM_WORLD)
Compute the energy transfer through the current filter scale.
Definition: compute_Pi_shift_deriv.cpp:21
void compute_Z(std::vector< double > &enstrophy_transfer, const dataset &source_data, const std::vector< double > &ux, const std::vector< double > &uy, const std::vector< double > &uz, const std::vector< double > &coarse_vort_r, const std::vector< double > &vort_ux, const std::vector< double > &vort_uy, const std::vector< double > &vort_uz, const MPI_Comm comm=MPI_COMM_WORLD)
Compute the energy transfer through the current filter scale.
Definition: compute_Z.cpp:24
void compute_vorticity(std::vector< double > &vort_r, std::vector< double > &vort_lon, std::vector< double > &vort_lat, std::vector< double > &vel_div, std::vector< double > &OkuboWeiss, const std::vector< double > &u_r, const std::vector< double > &u_lon, const std::vector< double > &u_lat, const int Ntime, const int Ndepth, const int Nlat, const int Nlon, const std::vector< double > &longitude, const std::vector< double > &latitude, const std::vector< bool > &mask, const MPI_Comm comm=MPI_COMM_WORLD)
Wrapper for computing vorticity.
Definition: compute_vorticity.cpp:23
size_t Index(const int Itime, const int Idepth, const int Ilat, const int Ilon, const int Ntime, const int Ndepth, const int Nlat, const int Nlon)
Convenience tool to convert physical index (time, depth, lat, lon) to a logical index.
Definition: Index.cpp:19
double depotential_temperature(const double p, const double theta)
Convert potential temperature to actual temperature following equation 1.154a from Vallis (p...
Definition: depotential_temperature.cpp:13
void get_lat_bounds(int &LAT_lb, int &LAT_ub, const std::vector< double > &latitude, const int Ilat, const double scale)
Get latitude bounds for coarse-graining.
Definition: get_lat_bounds.cpp:19
double kernel(const double distance, const double scale)
Primary kernel function coarse-graining procedure (G in publications)
Definition: kernel.cpp:30
void mask_out_pole(const std::vector< double > &latitude, std::vector< bool > &mask, const int Ntime, const int Ndepth, const int Nlat, const int Nlon)
If the grid includes the poles (lat = 90 degrees), then mask it out.
Definition: mask_out_pole.cpp:21
void filtering(const dataset &source_data, const std::vector< double > &scales, const MPI_Comm comm=MPI_COMM_WORLD)
Main filtering driver.
Definition: filtering.cpp:23
void roll_field(std::vector< double > &field_to_roll, const std::string dimension, const int roll_count, const int Ntime, const int Ndepth, const int Nlat, const int Nlon)
Roll field along dimension.
Definition: roll_field.cpp:19
void vel_Cart_to_Spher(std::vector< double > &u_r, std::vector< double > &u_lon, std::vector< double > &u_lat, const std::vector< double > &u_x, const std::vector< double > &u_y, const std::vector< double > &u_z, const dataset &source_data)
Wrapper that applies vel_Spher_to_Cart_at_point to every point in the domain.
Definition: vel_Cart_to_Spher.cpp:13
Class to store main variables.
Definition: functions.hpp:24
void compute_local_kernel(std::vector< double > &local_kernel, const double scale, const dataset &source_data, const int Ilat, const int Ilon, const int LAT_lb, const int LAT_ub)
Compute an array of the kernel values from a given reference point to every other point in the domain...
Definition: compute_local_kernel.cpp:22
Class to process command-line arguments.
Definition: functions.hpp:488
void compute_vorticity_at_point(double &vort_r_tmp, double &vort_lon_tmp, double &vort_lat_tmp, double &div_tmp, double &OkuboWeiss_tmp, const std::vector< double > &u_r, const std::vector< double > &u_lon, const std::vector< double > &u_lat, const int Ntime, const int Ndepth, const int Nlat, const int Nlon, const int Itime, const int Idepth, const int Ilat, const int Ilon, const std::vector< double > &longitude, const std::vector< double > &latitude, const std::vector< bool > &mask)
Compute the vorticity, divergnce, and OkuboWeiss at a point.
Definition: compute_vorticity_at_point.cpp:26
void convert_coordinates(std::vector< double > &longitude, std::vector< double > &latitude)
Apply degree->radian conversion, if necessary. Applied in-place.
Definition: convert_coordinates.cpp:19
void compute_spatial_average(std::vector< double > &means, const std::vector< double > &field, const std::vector< double > &areas, const int Ntime, const int Ndepth, const int Nlat, const int Nlon, const std::vector< bool > &mask)
Given a field, compute the horizontal average.
Definition: compute_spatial_average.cpp:21
void filtering_helmholtz(const dataset &source_data, const std::vector< double > &scales, const MPI_Comm comm=MPI_COMM_WORLD)
Main filtering driver for Helmholtz decomposed data.
Definition: filtering_helmholtz.cpp:24
void compute_div_transport(std::vector< double > &div_J, const std::vector< double > &u_x, const std::vector< double > &u_y, const std::vector< double > &u_z, const std::vector< double > &uxux, const std::vector< double > &uxuy, const std::vector< double > &uxuz, const std::vector< double > &uyuy, const std::vector< double > &uyuz, const std::vector< double > &uzuz, const std::vector< double > &coarse_p, const std::vector< double > &longitude, const std::vector< double > &latitude, const int Ntime, const int Ndepth, const int Nlat, const int Nlon, const std::vector< bool > &mask)
Compute KE transport caused by div(J)
Definition: compute_div_transport.cpp:94
double kernel_alpha(void)
Compute alpha value for kernel (for baroclinic transfer)
Definition: kernel_alpha.cpp:12
void apply_filter_at_point(std::vector< double *> &coarse_val, const std::vector< const std::vector< double > *> &fields, const dataset &source_data, const int Itime, const int Idepth, const int Ilat, const int Ilon, const int LAT_lb, const int LAT_ub, const double scale, const std::vector< bool > &use_mask, const std::vector< double > &local_kernel, const std::vector< double > *weight=NULL)
Compute filtered field at a single point.
Definition: apply_filter_at_point.cpp:25
void compute_Lambda_full(std::vector< double > &Lambda, const std::vector< double > &coarse_u_r, const std::vector< double > &coarse_u_lon, const std::vector< double > &coarse_u_lat, const std::vector< double > &tilde_u_r, const std::vector< double > &tilde_u_lon, const std::vector< double > &tilde_u_lat, const std::vector< double > &coarse_p, const int Ntime, const int Ndepth, const int Nlat, const int Nlon, const std::vector< double > &longitude, const std::vector< double > &latitude, const std::vector< bool > &mask)
Compute the full baroclinic transfer term Lambda (see Lees and Aluie 2019)
Definition: compute_Lambda_full.cpp:24
Provide namespace for global constants (physical and computational).
void compute_Pi_Helmholtz(std::vector< double > &energy_transfer, const dataset &source_data, const std::vector< double > &ulon, const std::vector< double > &ulat, const std::vector< double > &ulon_ulon, const std::vector< double > &ulon_ulat, const std::vector< double > &ulat_ulat, const MPI_Comm comm=MPI_COMM_WORLD)
Compute the energy transfer through the current filter scale.
Definition: compute_Pi_Helmholtz.cpp:22
Class for storing internal timings.
Definition: functions.hpp:450
void compute_Lambda_rotational(std::vector< double > &Lambda_rot, const std::vector< double > &coarse_vort_r, const std::vector< double > &coarse_vort_lon, const std::vector< double > &coarse_vort_lat, const std::vector< double > &coarse_rho, const std::vector< double > &coarse_p, const int Ntime, const int Ndepth, const int Nlat, const int Nlon, const std::vector< double > &longitude, const std::vector< double > &latitude, const std::vector< bool > &mask, const double scale_factor)
Compute the rotational component of the non-linear model of the baroclinic transfer term Lambda (see ...
Definition: compute_Lambda_rotational.cpp:25
double equation_of_state(const double T, const double S, const double p)
UNESCO 1981 equation of state.
Definition: equation_of_state.cpp:18
int get_omp_chunksize(const int Nlat, const int Nlon)
Compute openmp chunk size for OMP for loops.
Definition: get_omp_chunksize.cpp:18
void vel_Spher_to_Cart(std::vector< double > &u_x, std::vector< double > &u_y, std::vector< double > &u_z, const std::vector< double > &u_r, const std::vector< double > &u_lon, const std::vector< double > &u_lat, const dataset &source_data)
Wrapper that applies vel_Spher_to_Cart_at_point to every point in the domain.
Definition: vel_Spher_to_Cart.cpp:13
void vel_Spher_to_Cart_at_point(double &u_x, double &u_y, double &u_z, const double u_r, const double u_lon, const double u_lat, const double lon, const double lat)
Convert single spherical velocity to Cartesian velocity.
Definition: vel_Spher_to_Cart_at_point.cpp:28
bool string_to_bool(std::string str)
Convert string to boolean.
Definition: string_to_bool.cpp:12
double distance(const double lon1, const double lat1, const double lon2, const double lat2, const double Llon=0, const double Llat=0)
Compute the distance (in metres) between two points in the domain.
Definition: distance.cpp:25
const bool FILTER_OVER_LAND
Boolean to indicate whether or not land values should be filled in with coarse-grained results...
Definition: constants.hpp:118
void apply_filter_at_point_for_quadratics(double &uxux_tmp, double &uxuy_tmp, double &uxuz_tmp, double &uyuy_tmp, double &uyuz_tmp, double &uzuz_tmp, double &vort_ux_tmp, double &vort_uy_tmp, double &vort_uz_tmp, const std::vector< double > &u_x, const std::vector< double > &u_y, const std::vector< double > &u_z, const std::vector< double > &vort_r, const dataset &source_data, const int Itime, const int Idepth, const int Ilat, const int Ilon, const int LAT_lb, const int LAT_ub, const double scale, const std::vector< double > &local_kernel)
Compute filtered quadratic velocities at a single point.
Definition: apply_filter_at_point_for_quadratics.cpp:34
void print_compile_info(const std::vector< double > *scales=NULL)
Print summary of compile-time variables.
Definition: print_compile_info.cpp:16
void compute_Pi(std::vector< double > &energy_transfer, const dataset &source_data, const std::vector< double > &ux, const std::vector< double > &uy, const std::vector< double > &uz, const std::vector< double > &uxux, const std::vector< double > &uxuy, const std::vector< double > &uxuz, const std::vector< double > &uyuy, const std::vector< double > &uyuz, const std::vector< double > &uzuz, const MPI_Comm comm=MPI_COMM_WORLD)
Compute the energy transfer through the current filter scale.
Definition: compute_Pi.cpp:21
void compute_areas(std::vector< double > &areas, const std::vector< double > &longitude, const std::vector< double > &latitude)
Compute the area of each computational cell.
Definition: compute_areas.cpp:19