-
Notifications
You must be signed in to change notification settings - Fork 9
/
preprocess.hpp
343 lines (315 loc) · 11.3 KB
/
preprocess.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
#ifndef PREPROCESS_HPP
#define PREPROCESS_HPP 1
#include <stdio.h>
#include <stdlib.h>
#include "ALGLIB/linalg.h"
#include <mpi.h>
#include <vector>
/*!
* \file
* \brief Various pre-processing routines related to coarse-graining.
* - Apply ALGLIB interpolation routines to fill in masked areas
* - Apply ALGLIB least-squares solvers to apply toroidal projections
*/
/*!
* \addtogroup InterpolationRoutines
* @{
* \brief Functions directly pertaining to toroidal projection.
*/
/*!
* @ingroup InterpolationRoutines
*/
void interpolate_over_land(
std::vector<double> &interp_field,
const std::vector<double> &field,
const std::vector<double> &time,
const std::vector<double> &depth,
const std::vector<double> &latitude,
const std::vector<double> &longitude,
const std::vector<bool> &mask);
/*!
* @ingroup InterpolationRoutines
*/
void interpolate_over_land_from_coast(
std::vector<double> &interp_field,
const std::vector<double> &field,
const int nlayers,
const std::vector<double> &time,
const std::vector<double> &depth,
const std::vector<double> &latitude,
const std::vector<double> &longitude,
const std::vector<bool> &mask,
const std::vector<int> &myCounts,
const MPI_Comm comm = MPI_COMM_WORLD
);
/*!
* @ingroup InterpolationRoutines
*/
void get_coast(
std::vector<double> &lon_coast,
std::vector<double> &lat_coast,
std::vector<double> &field_coast,
const std::vector<double> &lon_full,
const std::vector<double> &lat_full,
const std::vector<double> &field_full,
const std::vector<bool> &mask,
const int Itime,
const int Idepth,
const int Ntime,
const int Ndepth,
const int Nlat,
const int Nlon);
/*!
* \addtogroup ToroidalProjection
* @{
* \brief Functions directly pertaining to toroidal projection.
*/
/*!
* \brief Applies toroidal projection to the specified velocity field.
* @ingroup ToroidalProjection
*
* Specifically, solves for (scalar field) F that minimizes error in
* \f$ \nabla_H^2 F = \nabla_H \times \vec{u} \cdot \hat{e}_r \f$
*
* This uses sparse differentiation matrices (of the order specified in constants.hpp) and the ALGLIB least-squares solver.
*
* *single_seed*: If true, then only one seed is provided and should be used as the seed for all different times. If false, then the provided seed incorporates a different seed for each time.
*
* @param[in] output_fname Name for the output file
* @param[in,out] u_lon,u_lat velocity field to be projected
* @param[in] time,depth,latitude,longitude dimension vectors (1D)
* @param[in] mask array to distinguish land/water
* @param[in] myCounts Local (to MPI process) dimension sizes
* @param[in] myStarts Vector indicating where the local (to MPI process) region fits in the whole
* @param[in] seed Seed for the least-squares solver
* @param[in] single_seed Indicates if a single seed is used - see notes
* @param[in] comm MPI communicator (default MPI_COMM_WORLD)
*
*/
void Apply_Toroidal_Projection(
const std::string output_fname,
dataset & source_data,
const std::vector<double> & seed,
const bool single_seed,
const double rel_tol,
const int max_iters,
const bool weight_err,
const bool use_mask,
const MPI_Comm comm = MPI_COMM_WORLD
);
void depth_integrate(
std::vector<double> & depth_integral,
const std::vector<double> & field_to_integrate,
const dataset & source_data,
const MPI_Comm comm = MPI_COMM_WORLD
);
void Apply_Potential_Projection(
const std::string output_fname,
dataset & source_data,
const std::vector<double> & seed,
const bool single_seed,
const double rel_tol,
const int max_iters,
const bool weight_err,
const bool use_mask,
const MPI_Comm comm = MPI_COMM_WORLD
);
void Apply_Helmholtz_Projection(
const std::string output_fname,
dataset & source_data,
const std::vector<double> & seed_tor,
const std::vector<double> & seed_pot,
const bool single_seed,
const double rel_tol,
const int max_iters,
const bool weight_err,
const bool use_mask,
const double Tikhov_Laplace,
const MPI_Comm comm = MPI_COMM_WORLD
);
void Apply_Helmholtz_Projection_uiuj(
const std::string output_fname,
dataset & source_data,
const std::vector<double> & seed_v_r,
const std::vector<double> & seed_v_lon,
const std::vector<double> & seed_v_lat,
const bool single_seed,
const double Tikhov_Lambda,
const double Tikhov_Laplace,
const double rel_tol,
const int max_iters,
const bool weight_err,
const bool use_mask,
const MPI_Comm comm = MPI_COMM_WORLD
);
void Apply_Helmholtz_Projection_SymTensor(
const std::string output_fname,
dataset & source_data,
const std::vector<double> & seed_v_r,
const std::vector<double> & seed_v_lon,
const std::vector<double> & seed_v_lat,
const bool single_seed,
const double rel_tol,
const int max_iters,
const bool weight_err,
const bool use_mask,
const MPI_Comm comm = MPI_COMM_WORLD
);
/*!
* \brief Computes the (toroidal) velocity corresponding to field F.
* @ingroup ToroidalProjection
*
* Specifically, \f$ \hat{e}_r \times \frac{1}{r}\nabla_H F \f$
*
* @param[in,out] vel_lon,vel_lat Where to store the toroidal velocities
* @param[in] F Field from which to compute the velocities
* @param[in] longitude,latitude Grid vectors (1D)
* @param[in] Ntime,Ndepth,Nlat,Nlon Dimension sizes
* @param[in] mask Array to distinguish land/water
*
*/
void toroidal_vel_from_F(
std::vector<double> & vel_lon,
std::vector<double> & vel_lat,
const std::vector<double> & F,
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
);
void potential_vel_from_F(
std::vector<double> & vel_lon,
std::vector<double> & vel_lat,
const std::vector<double> & F,
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
);
void uiuj_from_Helmholtz(
std::vector<double> & ulon_ulon,
std::vector<double> & ulon_ulat,
std::vector<double> & ulat_ulat,
const std::vector<double> & v_r,
const std::vector<double> & Phi_v,
const std::vector<double> & Psi_v,
const dataset & source_data
);
/*!
* \brief Computes the curl term (RHS) of the projection operation.
* @ingroup ToroidalProjection
*
* Computes \f$ \nabla_H \times \vec{u} \cdot \hat{e}_r \f$, which
* corresponds to the RHS of the least-squares problem.
*
* Uses the differentation order specified in contants.hpp
*
* *seed*: I couldn't figure out how to provide a seed directly to the solver. Instead,
* seeds are done in the following way. Call the seed \f$ x_0 \f$ and write \f$ x = x' + x_0 \f$.
* If \f$ Ax=b\f$ then \f$Ax' = b - Ax_0 \f$.
* The seed is applied in exactly this way to modify the RHS of the problem.
* The seed is then added back in afterwards.
*
* @param[in,out] Lap Where to store the (sparse) differentiation matrix
* @param[in] longitude,latitude Grid vectors (1D)
* @param[in] Itime,Idepth Indicates current time/depth iteration
* @param[in] Ntime,Ndepth,Nlat,Nlon Dimension sizes
* @param[in] mask Array to distinguish land/water
* @param[in] seed (optional) seed for the solver
*
*/
void toroidal_curl_u_dot_er(
std::vector<double> & out_arr,
const std::vector<double> & u_lon,
const std::vector<double> & u_lat,
const std::vector<double> & longitude,
const std::vector<double> & latitude,
const int Itime,
const int Idepth,
const int Ntime,
const int Ndepth,
const int Nlat,
const int Nlon,
const std::vector<bool> & mask,
const std::vector<double> * seed = NULL
);
void toroidal_sparse_Lap(
alglib::sparsematrix & Lap,
const dataset & source_data,
const int Itime,
const int Idepth,
const std::vector<bool> & mask,
const bool area_weight = false,
const size_t row_skip = 0,
const size_t column_skip = 0
);
void sparse_vel_from_PsiPhi(
alglib::sparsematrix & LHS_matr,
const dataset & source_data,
const int Itime,
const int Idepth,
const std::vector<bool> & mask,
const bool area_weight
);
/*!
* \brief This is just a helper to compute Lap(F). It's provided as an output for diagnostic purposes.
* @ingroup ToroidalProjection
*
* @param[in,out] out_arr Where to store the laplacian
* @param[in] F Field to differentiate
* @param[in] longitude,latitude Grid vectors (1D)
* @param[in] Ntime,Ndepth,Nlat,Nlon Dimension sizes
* @param[in] mask Array to distinguish land/water
*
*/
void toroidal_Lap_F(
std::vector<double> & out_arr,
const std::vector<double> & F,
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
);
/*!
* \brief This is just a helper to compute div(vel). It's provided as an output for diagnostic purposes.
* @ingroup ToroidalProjection
*
* @param[in,out] div Where to store the divergence
* @param[in] vel_lon,vel_lat Velocity fields
* @param[in] longitude,latitude Grid vectors (1D)
* @param[in] Ntime,Ndepth,Nlat,Nlon Dimension sizes
* @param[in] mask Array to distinguish land/water
*
*/
void toroidal_vel_div(
std::vector<double> & div,
const std::vector<double> & vel_lon,
const std::vector<double> & vel_lat,
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
);
void Extract_Beta_Geos_Vel(
std::vector<double> & u_beta,
std::vector<double> & v_beta,
const std::vector<double> & ssh,
const std::vector<bool> & mask,
dataset & source_data,
const double rel_tol,
const int max_iters,
const MPI_Comm comm = MPI_COMM_WORLD
);
#endif