-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathwingwakeinf.cpp
525 lines (479 loc) · 17.3 KB
/
wingwakeinf.cpp
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
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
#include <iostream>
#include <cmath>
#include <vector>
#include <cassert>
#include <atomic>
#include <thread>
#include <mutex>
#include <condition_variable>
#include <chrono>
#ifdef FSI_OPENMP
#include <omp.h>
#endif
#ifdef FSI_MPI
#include <mpi.h>
#endif
#pragma GCC diagnostic ignored "-Wunused-variable"
#pragma GCC diagnostic ignored "-Wunused-but-set-variable"
using namespace std;
using namespace chrono;
static const int X = 0,
Y = 1,
Z = 2;
/*
* wingwakeinf.cpp -
*
* Calculates influence coefficients for 3d VLM method.
*
* This is a MEX-file for MATLAB.
*
*/
#ifdef FSI_MEX_DEBUG
static FILE *fp = nullptr;
#endif
extern condition_variable cv_parent;
extern vector<double **> tandem_output;
extern int world_rank;
extern int min_workinc, max_workinc;
static const double pi=3.14159265358979;
static const double rcut=1.0e-11;
static mutex false_sharing_mtx;
#include <fsi_thread.h>
static inline void cross(double *xcrossy, const double *x, const double *y)
{
xcrossy[X]=x[Y]*y[Z]-x[Z]*y[Y];
xcrossy[Y]=x[Z]*y[X]-x[X]*y[Z];
xcrossy[Z]=x[X]*y[Y]-x[Y]*y[X];
}
static inline void vortex(double *pU,
const double x, const double y, const double z,
const double x1, const double y1, const double z1,
const double x2, const double y2, const double z2,
const double gamma)
{
const double r1[] = {x-x1, y-y1, z-z1};
const double ir1 = sqrt(r1[X] * r1[X] +
r1[Y] * r1[Y] +
r1[Z] * r1[Z]);
if (ir1 < rcut) return;
const double r2[] = {x-x2, y-y2, z-z2};
const double ir2 = sqrt(r2[X] * r2[X] +
r2[Y] * r2[Y] +
r2[Z] * r2[Z]);
if (ir2 < rcut) return;
double r1r2[3];
cross(r1r2, r1, r2);
const double square = r1r2[X] * r1r2[X] +
r1r2[Y] * r1r2[Y] +
r1r2[Z] * r1r2[Z];
if (square < rcut) return;
const double r0[] = {x2-x1, y2-y1, z2-z1};
const double coeff = gamma/4/pi/square*(r0[X] * (r1[X] / ir1 - r2[X] / ir2) +
r0[Y] * (r1[Y] / ir1 - r2[Y] / ir2) +
r0[Z] * (r1[Z] / ir1 - r2[Z] / ir2));
*pU++ += coeff * r1r2[X];
*pU++ += coeff * r1r2[Y];
*pU += coeff * r1r2[Z];
}
// pointers for colocaton points
#define CINIT(C, Cstart, nts3) const double \
*C##ij = Cstart, \
*C##i1j = C##ij + 3, \
*C##ij1 = C##ij + nts3, \
*C##i1j1 = C##i1j + nts3
#define CITER(C) \
C##ij = C##i1j, \
C##i1j += 3, \
C##ij1 = C##i1j1, \
C##i1j1 += 3
#define XYZOFFSET(B, offset) \
B##offset[X], \
B##offset[Y], \
B##offset[Z]
// pointers for boundary points
// assumes that B is never reset
#define BINIT(B, nx2) const double \
*B##i1j = B + 3, \
*B##ij1 = B + (3 * nx2), \
*B##i1j1 = B##i1j + (3 * nx2)
#define BITER(B) \
B = B##i1j, \
B##i1j += 3, \
B##ij1 = B##i1j1, \
B##i1j1 += 3
#define BLASTITER(B) \
B += (3 * 2); \
B##i1j += (3 * 2); \
B##ij1 += (3 * 2); \
B##i1j1 += (3 * 2)
inline void
wingwakeinf_thread::doublet(double * const uw,
const double CXitin, const double CYitin, const double CZitin,
const double *pB,
const double *pMU, const double * const MUE,
const int t) {
int i, j;
uw[X] = 0;
uw[Y] = 0;
uw[Z] = 0;
BINIT(pB, nx2);
for (j = 0; j < ny; j++) {
for (i = 0; i < nx; i++, pMU++) {
const double MUval = *pMU;
vortex(uw,
CXitin, CYitin, CZitin,
pB[X], pB[Y], pB[Z],
XYZOFFSET(pB, ij1), //*pBXij1, *pBYij1, *pBZij1,
MUval);
#if 0
printf("t=%d uvw1[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
t, *uw,
CXitin, CYitin, CZitin,
pB[0], pB[Y], pB[Z],
XYZOFFSET(pB, ij1), //*pBXij1, *pBYij1, *pBZij1,
MUval);
#endif
#ifdef FSI_MEX_DEBUG
// if (it==0 && in == 0)
if (fp != nullptr) fprintf(fp,"t=%d i=%d j=%d uvw1[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
t, i, j, *uw,
CXitin, CYitin, CZitin,
pB[0], pB[Y], pB[Z],
XYZOFFSET(pB, ij1), //*pBXij1, *pBYij1, *pBZij1,
MUval);
#endif
vortex(uw,
CXitin, CYitin, CZitin,
XYZOFFSET(pB, ij1), //*pBXij1, *pBYij1, *pBZij1,
XYZOFFSET(pB, i1j1), //*pBXi1j1, *pBYi1j1, *pBZi1j1,
MUval);
#ifdef FSI_MEX_DEBUG
if (0) if (fp != nullptr) fprintf(fp,"t=%d uvw2[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
t, *uw,
CXitin, CYitin, CZitin,
XYZOFFSET(pB, ij1), //*pBXij1, *pBYij1, *pBZij1,
XYZOFFSET(pB, i1j1), //*pBXi1j1, *pBYi1j1, *pBZi1j1,
MUval);
#endif
vortex(uw,
CXitin, CYitin, CZitin,
XYZOFFSET(pB, i1j1), //*pBXi1j1, *pBYi1j1, *pBZi1j1,
XYZOFFSET(pB, i1j), //*pBXi1j, *pBYi1j, *pBZi1j,
MUval);
#ifdef FSI_MEX_DEBUG
if (0) if (fp != nullptr) fprintf(fp,"t=%d uvw3[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
t, *uw,
CXitin, CYitin, CZitin,
XYZOFFSET(pB, i1j1), //*pBXi1j1, *pBYi1j1, *pBZi1j1,
XYZOFFSET(pB, i1j), //*pBXi1j, *pBYi1j, *pBZi1j,
MUval);
#endif
vortex(uw,
CXitin, CYitin, CZitin,
XYZOFFSET(pB, i1j), //*pBXi1j, *pBYi1j, *pBZi1j,
pB[X], pB[Y], pB[Z],
MUval);
#ifdef FSI_MEX_DEBUG
if (0) if (fp != nullptr) fprintf(fp,"t=%d uvw4[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
t, *uw,
CXitin, CYitin, CZitin,
XYZOFFSET(pB, i1j), //*pBXi1j, *pBYi1j, *pBZi1j,
pB[X], pB[Y], pB[Z],
MUval);
#endif
BITER(pB);
}
// skip over last boundary point and wrapping boundary point
BLASTITER(pB);
}
const double *Cstart = C,
*pMUEstart = MUE;
for (j = 0; j < ny; j++) {
CINIT(C, Cstart, nts3);
//cerr << "C=" << Cstart - C << endl;
Cstart += nts3;
const double *pMUE = pMUEstart;
pMUEstart += nts;
//cerr << "pMUE=" << pMUE - MUE << endl;
for (i = 0; i < t; i++, pMUE++, CITER(C)) {
const double MUEval = *pMUE;
vortex(uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, ij),
XYZOFFSET(C, ij1),
MUEval);
#ifdef FSI_MEX_DEBUG
if (fp != nullptr) fprintf(fp,"uvw1[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
*uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, ij),
XYZOFFSET(C, ij1),
MUEval);
#endif
vortex(uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, ij1),
XYZOFFSET(C, i1j1),
MUEval);
#ifdef FSI_MEX_DEBUG
if (0) if (fp != nullptr) fprintf(fp,"uvw2[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
*uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, ij1),
XYZOFFSET(C, i1j1),
MUEval);
#endif
vortex(uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, i1j1),
XYZOFFSET(C, i1j),
MUEval);
#ifdef FSI_MEX_DEBUG
if (0) if (fp != nullptr) fprintf(fp,"uvw3[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
*uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, i1j1),
XYZOFFSET(C, i1j),
MUEval);
#endif
vortex(uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, i1j), //*pBXi1j, *pBYi1j, *pBZi1j,
XYZOFFSET(C, ij), //*pBXij, *pBYij, *pBZij,
MUEval);
#ifdef FSI_MEX_DEBUG
if (0) if (fp != nullptr) fprintf(fp,"uvw4[0]=%f %f %f %f %f %f %f %f %f %f %f \n",
*uw,
CXitin, CYitin, CZitin,
XYZOFFSET(C, i1j), //*pBXi1j, *pBYi1j, *pBZi1j,
XYZOFFSET(C, ij), //*pBXij, *pBYij, *pBZij,
MUEval);
#endif
}
}
}
#ifdef FSI_MEX_DEBUG
static void openbug(const int t, const char *label) {
if (world_rank != 0) {
#ifdef _WIN32
string buf = "f:/data/wingwakeinfthread_T";
#else
string buf = "/home/wm/ml/WLM_PM_CODE/data/WingWakeInteractionT";
#endif
buf += to_string(t + 1) + ".dat";
fp = fopen(buf.c_str(), world_rank == -1 ? "w" : "a"); // for MPI, this is called repeatedly, so append
if (fp == nullptr) cerr << "Couldn't open %b " << buf << endl;
else {
fprintf(fp, "%s\n", label);
cerr << "Writing " << label << " wingwakeinf debugging file: " << buf << endl;
}
}
}
#endif
void wingwakeinf_thread::mex_thread_init(const int t, const double *sig, const double *sig2, double *uwakeS) {
this->t = t;
this->sig = sig;
this->sig2 = sig2;
this->uwakeS = uwakeS;
// init convenience constants (for current run)
t1 = t + 1;
t1n1 = t1 * n1;
workinc = calc_wake_workinc(t1n1, workfactor);
if (workinc > max_workinc) max_workinc = workinc;
if (workinc < min_workinc) min_workinc = workinc;
// init synchronized variables
order = 0;
threads_finished = 0;
#ifdef FSI_OPENMP
TSTART(mex_time, steady_clock::now()); // start mex timing
omp_set_num_threads(nthreads);
#pragma omp parallel
{
const int thrdid = omp_get_thread_num();
#ifdef FSI_STATS
int ndoublets = 0;
const steady_clock::time_point start = steady_clock::now();
#ifdef FSI_OPENMP
int openmp_denied = 0;
if (thrdid == 0) { // only for one thread
int nthrds = omp_get_num_threads();
if (nthreads != nthrds) openmp_denied += nthreads - nthrds;
}
#endif
#endif
double * const uwakerecv = tandem_output[thrdid][wingwakeinf_index];
// compute, writing to packed private buffer to avoid false sharing caching invalidation
int ord;
while ((ord = order.fetch_add(workinc)) < t1n1) { // while work remains
double *puwake = uwakerecv;
for (int i = 0, iord = ord; i < workinc && iord < t1n1; i++, iord++) {
const double * const pC = &C[coloc(iord)]; // indexes uwakeS(nts, n1)
assert(pC - C <= ((nts * n1) - 1) * 3);
doublet(puwake,
pC[X], pC[Y], pC[Z],
B,
sig, sig2,
t);
puwake += 3;
#ifdef FSI_STATS
ndoublets++;
#endif
}
// distribute computations from packed buffer
const double *puwakerecv = uwakerecv;
{ // start false sharing critical section
#pragma omp critical
for (int i = 0, iord = ord; i < workinc && iord < t1n1; i++, iord++) {
puwake = &uwakeS[coloc(iord)];
assert(puwake - uwakeS < nts * n1 * 3);
*puwake++ += *puwakerecv++; // X
*puwake++ += *puwakerecv++; // Y
*puwake += *puwakerecv++; // Z
}
} // end false sharing critical section
} // while work remains
#ifdef FSI_STATS
// per-thread statistics
mex_stat *pmex_stat = (*thread_stat_map[thrdid])[this];
pmex_stat->call_count++;
pmex_stat->elapsed += duration_cast<mex_thread_time_units_t>(steady_clock::now() - start);
pmex_stat->doublets += ndoublets;
#ifdef FSI_OPENMP
if (openmp_denied != 0) pmex_stat->openmp_denied += openmp_denied;
#endif
#endif
} // end parallel
TSTOP(mex_time, steady_clock::now()); // stop mex timing and add elapsed
#else // low-level threading model
#ifdef FSI_MEX_DEBUG
openbug(t, "threaded");
#endif
mex_parent(); // awaken child threads and wait for them to finish
#endif
}
void wingwakeinf_thread::mex_thread_run (const int partner_rank) {
// WX(nts, ny + 1)
// BX(nx2, n1)
// uwakeS(nts, n1) - row ranges 0:t
// mat MUEA1(2 * m, n); remember that 2*m = nx mat MUEW1doubled(nts, n);
double * const uwakerecv = tandem_output[partner_rank][wingwakeinf_index];
//double * const uwakerecv = new double [outsize];
#ifdef FSI_MPI
bool send_sig = true;
#endif
#ifdef FSI_STATS
int ndoublets = 0;
int nsyncs = 0;
#endif
for (;;) { // until no more data
const int ord = order.fetch_add(workinc);
#ifdef FSI_STATS
nsyncs++;
#endif
if (ord >= t1n1) {
#ifdef FSI_STATS
#ifndef FSI_OPENMP
// statistics
mex_stat *pmex_stat = (*thread_stat_map[this_thread::get_id()])[this];
pmex_stat->doublets += ndoublets;
pmex_stat->nsyncs += nsyncs;
#endif
#endif
//if (itin == ntsn1 + nthreads - 1)...
// handle synchronization with parent in class
if (is_last_thread()) {
#ifdef FSI_MEX_DEBUG
if (fp != nullptr) fclose(fp);
#endif
}
return; // to wait state
}
double *puwake;
#ifdef FSI_MPI
// supply order, indicate INFLUENCE (default)
int tag = wingwakeinf_mask | ord;
//cerr << "wingwakeinf mex_thread_run:sending to partner:" << partner_rank << " ord=" << ord << " workinc="<< workinc<<endl;
int status = MPI_Send(sig, send_sig ? sigsize : 0, MPI_DOUBLE, partner_rank, wingwakeinf_mask | ord, MPI_COMM_WORLD);
assert(status == MPI_SUCCESS);
send_sig = false; // send to remote child only once per time step
try {
status = MPI_Recv(uwakerecv, outsize, MPI_DOUBLE, partner_rank, MPI_ANY_TAG, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
}
catch (exception& e) {
cerr << "exception caught on MPI_Recv in wingwakeinf "
"Exception: " << e.what() << '\n';
}
assert(status == MPI_SUCCESS);
#else
// compute, writing to packed private buffer to avoid false sharing caching invalidation
puwake = uwakerecv;
for (int i = 0, iord = ord; i < workinc && iord < t1n1; i++, iord++, puwake += 3) {
const double * const pC = &C[coloc(iord)]; // indexes uwakeS(nts, n1)
assert(pC - C <= ((nts * n1) - 1) * 3);
//cerr << "it=" << it << " in=" << in << " itin=" << itin << " t=" << t << " n1=" << n1 << endl;
doublet(puwake,
pC[X], pC[Y], pC[Z],
B,
sig, sig2,
t);
}
#endif
// write computations from packed buffer
const double *puwakerecv = uwakerecv;
{ // start false sharing critical section
if (nthreads > 1) unique_lock<mutex> lck(false_sharing_mtx);
for (int i = 0, iord = ord; i < workinc && iord < t1n1; i++, iord++) {
puwake = &uwakeS[coloc(iord)];
assert(puwake - uwakeS < nts * n1 * 3);
//cerr << "puwake=" << puwake - uwakeS << " recv=" << *puwakerecv << endl;
*puwake++ += *puwakerecv++; // X
*puwake++ += *puwakerecv++; // Y
*puwake += *puwakerecv++; // Z
#ifdef FSI_STATS
ndoublets++;
#endif
}
} // end false sharing critical section
}
}
#ifdef FSI_MPI
void wingwakeinf_thread::compute(const int pord, const double * const psig, const int pt) {
#ifdef FSI_MEX_DEBUG
openbug(pt, "computed");
#endif
// wingwakeinf - mat MUEA1(2 * m, n); mat MUEW1doubled(nts, n);
//cube uw(3, nts, n + 1); ZERO(uw);
int ord = pord & ~MEX_MASK;
t = pt;
// init convenience constants (for current run)
t1 = t + 1;
t1n1 = t1 * n1;
workinc = calc_wake_workinc(t1n1, workfactor);
//workinc = t1n1 / nthreads * workfactor + 1; // zero here causes deadlock
sig = psig;
sig2 = &psig[nx * ny];
// reserve space for both A and B
if (uwakeS == nullptr) {
uwakeS = new double [outsize]; //% Influence co-efficient matrix of surface doublet distribution.
}
double *puwake = uwakeS;
for (int nwork = 0; nwork < workinc && ord < t1n1; nwork++, ord++, puwake += 3) {
// decode for row, col
const double * const pC = &C[coloc(ord)]; // indexes uwakeS(nts, n1)
//cerr << "it=" << it << " in=" << in << " itin=" << itin << " t=" << t << " n1=" << n1 << endl;
assert(pC - C <= ((nts * n1) - 1) * 3);
doublet(puwake,
pC[X], pC[Y], pC[Z],
B,
sig, sig2,
t);
}
// operates on only a single XYZ colocation at a time; hence, just multiply work done by 3
assert(puwake - uwakeS <= outsize);
int status = MPI_Send(uwakeS, puwake - uwakeS, MPI_DOUBLE, 0, wingwakeinf_mask | pord, MPI_COMM_WORLD);
assert(status == MPI_SUCCESS);
#ifdef FSI_MEX_DEBUG
if (fp != nullptr) fclose(fp);
#endif
}
#endif