Loading [MathJax]/extensions/tex2jax.js
MAIA bb96820c
Multiphysics at AIA
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fvzonalstg.cpp
Go to the documentation of this file.
1// Copyright (C) 2024 The m-AIA AUTHORS
2//
3// This file is part of m-AIA (https://git.rwth-aachen.de/aia/m-AIA/m-AIA)
4//
5// SPDX-License-Identifier: LGPL-3.0-only
6
7#include "fvzonalstg.h"
8
9#include <algorithm>
10#include <iostream>
11#include <stack>
12#include <vector>
15#include "MEMORY/alloc.h"
16#include "UTIL/functions.h"
17#include "UTIL/kdtree.h"
18#include "globals.h"
19#include "globalvariables.h"
20
21using namespace std;
22
23
24template <MInt nDim, class SysEqn>
25FvZonalSTG<nDim, SysEqn>::FvZonalSTG(const MInt couplingId, RANS* R, LES* L)
26 : Coupling(couplingId), FvZonal<nDim, SysEqn>(couplingId, R, L) {
27 // Coupling(couplingId), m_RANSSolver(R), m_LESSolver(L) {
28
29 m_log << "Create Zonal coupler for RANS Solver (" << RANSSolver().m_solverId << ") and LES Solver ("
30 << LESSolver().m_solverId << ")" << endl;
31
39 m_uvRANSFactor = F1;
40 if(Context::propertyExists("uvRANSFactor")) {
41 m_uvRANSFactor = Context::getBasicProperty<MFloat>("uvRANSFactor", AT_, &m_uvRANSFactor);
42 }
43}
44
45
46template <MInt nDim, class SysEqn>
48 TRACE();
49
50 m_RANSSolverId = RANSSolver().m_solverId;
51 m_LESSolverId = LESSolver().m_solverId;
52
53 LESSolver().m_LESNoVarAverage = LESVarAverageData::noAvgVars;
54
55 // necessary for STG Method
56 initRANSValues();
57 initLESValues();
58
59 determineZonalPositions();
60
61 if(m_cylindricCommunication) {
62 initCylinderExchange();
63 // calcLESSectorAverage();
64 calcRANSSectorValues();
65 }
66
67 // necessary for boundary condition which is called in finalizeInitSolver
68 transferSolverData();
69}
70
71
72template <MInt nDim, class SysEqn>
74 TRACE();
75
76 if(m_cylindricCommunication) {
77 initCylinderExchange();
78 }
79
80 if(m_cylindricCommunication) {
81 // calcLESSectorAverage();
82 calcRANSSectorValues();
83 }
84
85 transferSolverData();
86}
87
88
89template <MInt nDim, class SysEqn>
91 TRACE();
92
93 if(solverId != m_RANSSolverId || solverId != m_LESSolverId) return;
94
95 transferSolverData();
96}
97
98
99template <MInt nDim, class SysEqn>
101 TRACE();
102
103 transferSolverData();
104}
105
106
107template <MInt nDim, class SysEqn>
109 TRACE();
110
111 if(globalTimeStep % m_zonalTransferInterval == 0 || globalTimeStep == LESSolver().m_restartTimeStep + 1) {
112 if(m_cylindricCommunication) {
113 // initCylinderExchange();
114 // calcLESSectorAverage();
115 calcRANSSectorValues();
116 }
117
118 transferSolverData();
119 }
120}
121
122
123template <MInt nDim, class SysEqn>
125 TRACE();
126
127 LESSolver().m_averagePos.clear();
128 LESSolver().m_averageDir.clear();
129
130 // if(m_preliminarySTGSponge) {
131 // m_7901Position = m_preliminarySTGSpongePosition;
132 // m_averagePos = m_7901Position;
133 // m_averageDir = m_preliminarySTGSpongeDirection;
134 // } else {
135
136 m_7901faceNormalDir = -1;
137 m_7909faceNormalDir = -1;
138 m_7901Position = std::numeric_limits<MFloat>::infinity();
139 m_7909Position = std::numeric_limits<MFloat>::infinity();
140
141 // read cutoff positions to determine average and sponge cells
142 if(Context::propertyExists("cutOffBndryIds")) {
143 MInt propertyLength = Context::propertyLength("cutOffBndryIds", m_RANSSolverId);
144 for(MInt i = 0; i < propertyLength; i++) {
145 MInt bcId = Context::getSolverProperty<MFloat>("cutOffBndryIds", m_RANSSolverId, AT_, i);
146 if(bcId == 7901) {
147 m_7901faceNormalDir = Context::getSolverProperty<MFloat>("cutOffDirections", m_RANSSolverId, AT_, i);
148 }
149 }
150 }
151 if(m_7901faceNormalDir != -1) {
152 m_7901faceNormalDir = (m_7901faceNormalDir % 2 == 0) ? m_7901faceNormalDir + 1 : m_7901faceNormalDir - 1;
153 }
154
155 if(Context::propertyExists("cutOffBndryIds")) {
156 MInt propertyLength = Context::propertyLength("cutOffBndryIds", m_LESSolverId);
157 for(MInt i = 0; i < propertyLength; i++) {
158 MInt bcId = Context::getSolverProperty<MFloat>("cutOffBndryIds", m_LESSolverId, AT_, i);
159 if(bcId == 7909) {
160 m_7909faceNormalDir = Context::getSolverProperty<MFloat>("cutOffDirections", m_LESSolverId, AT_, i);
161 }
162 }
163 }
164 if(m_7909faceNormalDir != -1) {
165 m_7909faceNormalDir = (m_7909faceNormalDir % 2 == 0) ? m_7909faceNormalDir + 1 : m_7909faceNormalDir - 1;
166 }
167
168 if(m_STGSponge && m_7901faceNormalDir == -1) {
169 TERMM(1, "m_STGSponge true, bc7901 has to be set to determine spongeCells!");
170 }
171
172 m_7901wallDir = 1;
173 m_7901periodicDir = 2;
174
175 if(Context::propertyExists("bc7901wallDir")) {
176 m_7901wallDir = Context::getSolverProperty<MInt>("bc7901wallDir", m_RANSSolverId, AT_);
177 }
178 if(Context::propertyExists("bc7901periodicDir")) {
179 m_7901periodicDir = Context::getSolverProperty<MInt>("bc7901periodicDir", m_RANSSolverId, AT_);
180 }
181
182 // determine real m_7909Position based on cell centers
183 MFloat tmpPos7901 = ((m_7901faceNormalDir + 1) % 2 > 0) ? -std::numeric_limits<MFloat>::infinity()
184 : std::numeric_limits<MFloat>::infinity();
185 MFloat tmpPos7909 = ((m_7909faceNormalDir + 1) % 2 > 0) ? -std::numeric_limits<MFloat>::infinity()
186 : std::numeric_limits<MFloat>::infinity();
187 MFloat compFactor7901 = ((m_7901faceNormalDir + 1) % 2 > 0) ? F1 : -F1;
188 MFloat compFactor7909 = ((m_7909faceNormalDir + 1) % 2 > 0) ? F1 : -F1;
189
190 if(LESSolver().grid().isActive()) {
191 for(MInt cellId = 0; cellId < a_noFvGridCellsLES(); cellId++) {
192 MInt RANSId = convertIdParent(LESSolver(), RANSSolver(), cellId);
193 if(RANSId != -1) {
194 MFloat pos = LESSolver().a_coordinate(cellId, 0);
195 if(compFactor7901 * pos > compFactor7901 * tmpPos7901) {
196 tmpPos7901 = pos;
197 }
198 if(compFactor7909 * pos > compFactor7909 * tmpPos7909) {
199 tmpPos7909 = pos;
200 }
201 }
202 }
203
204 if(compFactor7901 > 0) {
205 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7901, 1, MPI_DOUBLE, MPI_MAX, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE",
206 "tmpPos7901");
207 } else {
208 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7901, 1, MPI_DOUBLE, MPI_MIN, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE",
209 "tmpPos7901");
210 }
211 if(compFactor7909 > 0) {
212 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7909, 1, MPI_DOUBLE, MPI_MAX, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE",
213 "tmpPos7909");
214 } else {
215 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7909, 1, MPI_DOUBLE, MPI_MIN, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE",
216 "tmpPos7909");
217 }
218
219 if(m_7901faceNormalDir != -1) {
220 m_7901Position = tmpPos7901;
221 }
222 if(m_7909faceNormalDir != -1) {
223 m_7909Position = tmpPos7909;
224 }
225 }
226
227 if(RANSSolver().grid().isActive()) {
228 if(compFactor7901 > 0) {
229 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7901, 1, MPI_DOUBLE, MPI_MAX, RANSSolver().mpiComm(), AT_, "MPI_IN_PLACE",
230 "tmpPos7901");
231 } else {
232 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7901, 1, MPI_DOUBLE, MPI_MIN, RANSSolver().mpiComm(), AT_, "MPI_IN_PLACE",
233 "tmpPos7901");
234 }
235 if(compFactor7909 > 0) {
236 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7909, 1, MPI_DOUBLE, MPI_MAX, RANSSolver().mpiComm(), AT_, "MPI_IN_PLACE",
237 "tmpPos7909");
238 } else {
239 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7909, 1, MPI_DOUBLE, MPI_MIN, RANSSolver().mpiComm(), AT_, "MPI_IN_PLACE",
240 "tmpPos7909");
241 }
242 MPI_Allreduce(MPI_IN_PLACE, &tmpPos7909, 1, MPI_DOUBLE, MPI_MAX, RANSSolver().mpiComm(), AT_, "MPI_IN_PLACE",
243 "tmpPos7909");
244
245 if(m_7901faceNormalDir != -1) {
246 m_7901Position = tmpPos7901;
247 }
248 if(m_7909faceNormalDir != -1) {
249 m_7909Position = tmpPos7909;
250 }
251 }
252
253 m_averagePos = m_7901Position;
254 m_averageDir = abs(m_7901wallDir + m_7901periodicDir - nDim);
255
256 if(!std::isinf(m_7901Position)) {
257 // positions to find cells for LES average
258 LESSolver().m_7901Position = m_7901Position;
259 LESSolver().m_averagePos.push_back(m_averagePos);
260 LESSolver().m_averageDir.push_back(m_averageDir);
261 LESSolver().m_averageReconstructNut.push_back(false);
262 }
263}
264
265
266template <MInt nDim, class SysEqn>
268 TRACE();
269
270 // write LES data to RANS Solver
271 if(RANSSolver().grid().isActive()) {
272 if(!m_cylindricCommunication) {
273 for(MInt cellId = 0; cellId < a_noFvGridCellsRANS(); cellId++) {
274 MInt LESId = convertIdParent(RANSSolver(), LESSolver(), cellId);
275 if(LESId != -1) {
276 vector<MInt>::iterator findAvgId =
277 find(LESSolver().m_LESAverageCells.begin(), LESSolver().m_LESAverageCells.end(), LESId);
278 if(findAvgId != LESSolver().m_LESAverageCells.end()) {
279 MInt avgId = distance(LESSolver().m_LESAverageCells.begin(), findAvgId);
280 for(MInt var = 0; var < noRANSVariables(); var++) {
281 ASSERT(cellId < (MInt)RANSSolver().m_LESValues[var].size(),
282 "Trying to access data [" + to_string(var) + "][" + to_string(cellId)
283 + "] in LESSolver().m_LESVarAverage(RANS) with length "
284 + to_string(RANSSolver().m_LESValues[var].size())
285 + ", domainId: " + to_string(RANSSolver().domainId()));
286
287 RANSSolver().m_LESValues[var][cellId] = LESSolver().m_LESVarAverage[var][avgId];
288 }
289 }
290 }
291 }
292 //}
293 } else {
294 if(m_cylindricCommunication && m_cylinderExchangeIdsOffset > -1) {
295 for(MInt exchangeIndex = 0; exchangeIndex < m_noRANSExchangeCells; exchangeIndex++) {
296 MInt exchangeIndexOffset = exchangeIndex + m_cylinderExchangeIdsOffset;
297 MInt RANSId = m_globalCylinderExchangeIds[exchangeIndexOffset];
298 if(RANSSolver().a_isBndryGhostCell(RANSId)) continue;
299 for(MInt var = 0; var < noRANSVariables(); var++) {
300 ASSERT(RANSId < (MInt)RANSSolver().m_LESValues[var].size(),
301 "Trying to access data [" + to_string(var) + "][" + to_string(RANSId)
302 + "] in LESSolver().m_LESVarAverage(RANS) with length "
303 + to_string(RANSSolver().m_LESValues[var].size())
304 + ", domainId: " + to_string(RANSSolver().domainId()) + " | "
305 << RANSSolver().c_noCells());
306
307 RANSSolver().m_LESValues[var][RANSId] =
308 m_globalCylinderLESExchangeValues[exchangeIndexOffset * LESSolver().m_LESNoVarAverage + var];
309 }
310 }
311 }
312 }
313 }
314
315 if(!m_cylindricCommunication) {
316 // write RANS data to LES Solver
317 if(LESSolver().grid().isActive()) {
318 for(MInt LESId = 0; LESId < a_noFvGridCellsLES(); LESId++) {
319 MInt RANSId = convertIdParent(LESSolver(), RANSSolver(), LESId);
320
321 if(RANSId != -1) {
322 for(MInt var = 0; var < noRANSVariables(); var++) {
323 ASSERT(LESId < (MInt)LESSolver().m_RANSValues[var].size(),
324 "Trying to access data [" + to_string(var) + "][" + to_string(LESId)
325 + "] in m_RANSValues with length " + to_string(LESSolver().m_RANSValues[var].size())
326 + ", domainId: " + to_string(LESSolver().domainId()) + " | "
327 << LESSolver().c_noCells());
328
329 LESSolver().m_RANSValues[var][LESId] = RANSSolver().a_pvariable(RANSId, var);
330 }
331
332 if(m_STGSponge) {
333 // calculate m_uvRans
334 const MFloat fre = 1.0 / RANSSolver().sysEqn().m_Re0;
335 const MFloat nu_t = RANSSolver().a_pvariable(RANSId, RANSSolver().m_sysEqn.PV->N);
336
337 MFloat du[3][3]{{F0, F0, F0}, {F0, F0, F0}, {F0, F0, F0}};
338 const MInt recData = RANSSolver().a_reconstructionData(RANSId);
339 const MFloat u[3] = {RANSSolver().a_pvariable(RANSId, RANSSolver().m_sysEqn.PV->U),
340 RANSSolver().a_pvariable(RANSId, RANSSolver().m_sysEqn.PV->V),
341 RANSSolver().a_pvariable(RANSId, RANSSolver().m_sysEqn.PV->W)};
342
343 for(MInt nghbr = 0; nghbr < RANSSolver().a_noReconstructionNeighbors(RANSId); nghbr++) {
344 const MInt recNghbrId = RANSSolver().a_reconstructionNeighborId(RANSId, nghbr);
345 if(recNghbrId > -1) {
346 const MFloat recConst_x = RANSSolver().m_reconstructionConstants[nDim * (recData + nghbr) + 0];
347 const MFloat recConst_y = RANSSolver().m_reconstructionConstants[nDim * (recData + nghbr) + 1];
348 const MFloat recConst_z = RANSSolver().m_reconstructionConstants[nDim * (recData + nghbr) + 2];
349 for(MInt dim = 0; dim < nDim; ++dim) {
350 const MFloat delta_u =
351 RANSSolver().a_pvariable(recNghbrId, RANSSolver().m_sysEqn.PV->VV[dim]) - u[dim];
352 du[dim][0] += recConst_x * delta_u;
353 du[dim][1] += recConst_y * delta_u;
354 du[dim][2] += recConst_z * delta_u;
355 }
356 }
357 }
358 MFloat sij[3][3]{{F0, F0, F0}, {F0, F0, F0}, {F0, F0, F0}};
359 MFloat SijSij = F0;
360 for(MInt d1 = 0; d1 < nDim; d1++) {
361 for(MInt d2 = 0; d2 < nDim; d2++) {
362 sij[d1][d2] = 0.5 * (du[d1][d2] + du[d2][d1]);
363 SijSij += sij[d1][d2] * sij[d1][d2];
364 }
365 }
366 const MFloat sr1 = (sij[0][1] + sij[1][0]) * (sij[0][1] + sij[1][0]);
367 const MFloat sr2 = (sij[1][2] + sij[2][1]) * (sij[1][2] + sij[2][1]);
368 const MFloat sr3 = (sij[0][2] + sij[2][0]) * (sij[0][2] + sij[2][0]);
369 const MFloat srt = std::max(sqrt(sr1 + sr2 + sr3), epss);
370 const MFloat rr1 = sqrt(sr1) / srt;
371
372 MFloat uvRans = -sqrt(2.0 * SijSij) * rr1 * nu_t * fre;
373 LESSolver().m_RANSValues[noRANSVariables()][LESId] = uvRans;
374 }
375 }
376 }
377 }
378 } else {
379 cylinderExchange();
380 }
381}
382
383template <MInt nDim, class SysEqn>
385 TRACE();
386
387 // write STG Sponge data to LES block
388 if(LESSolver().grid().isActive()) {
389 if(m_cylindricCommunication) {
390 for(MInt cellId = 0; cellId < a_noFvGridCellsLES(); cellId++) {
391 // write value to cell
392 if(LESSolver().a_isHalo(cellId)) continue;
393 MFloat y = LESSolver().a_coordinate(cellId, 1);
394 MFloat z = LESSolver().a_coordinate(cellId, 2);
395 MFloat sector = maia::math::getSector(y, z, m_azimuthalAngle);
396
397 for(MInt var = 0; var < nDim + 3; var++) {
398 ASSERT(cellId < (MInt)LESSolver().m_STGSpongeFactor[var].size(),
399 "Trying to access data [" + to_string(var) + "][" + to_string(cellId)
400 + "] in m_STGSpongeFactor with length " + to_string(LESSolver().m_STGSpongeFactor[var].size())
401 + ", domainId: " + to_string(LESSolver().domainId()));
402
403 if(m_periodicSpongeInterpolationIndex[cellId].size() == 1) {
404 MInt interpolationIndex = m_periodicSpongeInterpolationIndex[cellId][0];
405 MInt exchangeIndex = m_periodicSpongeCylinderExchangeIndex[interpolationIndex];
406 MFloat data = F0;
407 if(var < nDim) {
408 data = LESSolver().a_pvariable(cellId, var)
409 - m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + var];
410 }
411 if(var == nDim) {
412 data = m_uvRans[interpolationIndex];
413 }
414 if(var == nDim + 1) {
415 data = alpha * m_uvErr[interpolationIndex] + beta * m_uvInt[interpolationIndex];
416 }
417 if(var == nDim + 2) {
418 MFloat uv_LES =
419 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + noLESVariables()
420 + 1]
421 - m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + 0]
422 * m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + 1];
423 data = uv_LES;
424 }
425 LESSolver().m_STGSpongeFactor[var][cellId] = data;
426 }
427
428 if(m_periodicSpongeInterpolationIndex[cellId].size() > 1) {
429 // Interpolation - Least square method
430 MFloat n = F0;
431 MFloat Eyi = F0;
432 MFloat Eyi2 = F0;
433 MFloat Ezi = F0;
434 MFloat Ezi2 = F0;
435 MFloat Eyizi = F0;
436 MFloat Eui = F0;
437 MFloat Euiyi = F0;
438 MFloat Euizi = F0;
439
440 for(MInt i = 0; i < (MInt)m_periodicSpongeInterpolationIndex[cellId].size(); i++) {
441 // MBool tmp = false;
442 MInt interpolationIndex = m_periodicSpongeInterpolationIndex[cellId][i];
443 MInt exchangeIndex = m_periodicSpongeCylinderExchangeIndex[interpolationIndex];
444
445 MFloat yExchange = m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 1];
446 MFloat zExchange = m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 2];
447 MFloat rotationAngle =
448 -((sector - m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 3]) * m_azimuthalAngle)
449 / 180.0 * PI;
450 MFloat yRotated = cos(rotationAngle) * yExchange + sin(rotationAngle) * zExchange;
451 MFloat zRotated = -sin(rotationAngle) * yExchange + cos(rotationAngle) * zExchange;
452
453 MFloat delta_y = LESSolver().a_coordinate(cellId, 1) - yRotated;
454 MFloat delta_z = LESSolver().a_coordinate(cellId, 2) - zRotated;
455
456
457 n += 1;
458 Eyi += delta_y;
459 Eyi2 += POW2(delta_y);
460 Ezi += delta_z;
461 Ezi2 += POW2(delta_z);
462 Eyizi += delta_y * delta_z;
463
464 MFloat data = F0;
465 if(var < nDim) {
466 data = LESSolver().a_pvariable(cellId, var)
467 - m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + var];
468 }
469 if(var == nDim) {
470 data = m_uvRans[interpolationIndex];
471 }
472 if(var == nDim + 1) {
473 data = alpha * m_uvErr[interpolationIndex] + beta * m_uvInt[interpolationIndex];
474 }
475 if(var == nDim + 2) {
476 MFloat uv_LES =
477 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + noLESVariables()
478 + 1]
479 - m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + 0]
480 * m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + 1];
481 data = uv_LES;
482 }
483 Eui += data;
484 Euiyi += data * delta_y;
485 Euizi += data * delta_z;
486 }
487
488 MFloat detA =
489 n * Eyi2 * Ezi2 + 2 * Eyi * Eyizi * Ezi - (Ezi * Eyi2 * Ezi + n * Eyizi * Eyizi + Ezi2 * Eyi * Eyi);
490
491 if(detA < 1e-12) {
492 MInt interpolationIndex = m_periodicSpongeInterpolationIndex[cellId][0];
493 LESSolver().m_STGSpongeFactor[var][cellId] =
494 m_globalCylinderRANSExchangeValues[interpolationIndex * m_noRANSCylinderExchangeVariables + var];
495 } else {
496 LESSolver().m_STGSpongeFactor[var][cellId] =
497 1 / detA
498 * (Eui * (Eyi2 * Ezi2 - Eyizi * Eyizi) + Euiyi * (Ezi * Eyizi - Eyi * Ezi2)
499 + Euizi * (Eyi * Eyizi - Eyi2 * Ezi));
500 }
501 }
502 }
503 }
504 } else {
505 for(MInt LESId = 0; LESId < a_noFvGridCellsLES(); LESId++) {
506 MInt RANSId = convertIdParent(LESSolver(), RANSSolver(), LESId);
507 if(RANSId != -1) {
508 for(MInt i = 0; i < LESSolver().m_globalNoSpongeLocations; i++) {
509 if(abs(LESSolver().m_globalSpongeLocations[i].first - RANSSolver().a_coordinate(RANSId, m_7901wallDir))
510 < eps) {
511 LESSolver().m_uvRans[i] = m_uvRans[i];
512 }
513 }
514 }
515 }
516 }
517 }
518}
519
520
524template <MInt nDim, class SysEqn>
526 TRACE();
527
528 if(LESSolver().grid().isActive()) {
529 std::vector<MFloat> RANSSectors;
530 RANSSectors.push_back(std::numeric_limits<MFloat>::max());
531 RANSSectors.push_back(-std::numeric_limits<MFloat>::max());
532 std::vector<MFloat> RANSSectorLimits;
533 RANSSectorLimits.push_back(std::numeric_limits<MFloat>::max());
534 RANSSectorLimits.push_back(-std::numeric_limits<MFloat>::max());
535 std::vector<MInt> cylinderExchangeIds;
536 std::vector<MFloat> cylinderExchangeLocations;
537
538 mAlloc(m_RANSSectors, 2, "m_RANSSectors", F0, FUN_);
539 mAlloc(m_RANSSectorLimits, 2, "m_RANSSectorLimits", F0, FUN_);
540
541 MInt noRANSExchangeCells = 0;
542 for(MInt RANSId = 0; RANSId < a_noFvGridCellsRANS(); RANSId++) {
543 MFloat halfCellLength = RANSSolver().grid().halfCellLength(RANSId);
544 MFloat pos = RANSSolver().a_coordinate(RANSId, 0);
545 if(approx(m_7901Position, pos, halfCellLength) || approx(m_7909Position, pos, halfCellLength)) {
546 // MInt RANSId = convertIdParent(RANSSolver(), RANSSolver(), RANSId);
547 // if(RANSId > -1){
548 if(RANSSolver().c_isLeafCell(RANSId) && !RANSSolver().a_isHalo(RANSId)) {
549 MFloat x = RANSSolver().a_coordinate(RANSId, 0);
550 MFloat y = RANSSolver().a_coordinate(RANSId, 1);
551 MFloat z = RANSSolver().a_coordinate(RANSId, 2);
552 MFloat sector = maia::math::getSector(y, z, m_azimuthalAngle);
553 MFloat angle = maia::math::getAngle(y, z);
554
555 if(sector < RANSSectors[0]) {
556 RANSSectors[0] = sector;
557 }
558 if(sector > RANSSectors[1]) {
559 RANSSectors[1] = sector;
560 }
561 if(angle < RANSSectorLimits[0]) {
562 RANSSectorLimits[0] = angle;
563 }
564 if(angle > RANSSectorLimits[1]) {
565 RANSSectorLimits[1] = angle;
566 }
567
568 noRANSExchangeCells++;
569 cylinderExchangeIds.push_back(RANSId);
570 cylinderExchangeLocations.push_back(x);
571 cylinderExchangeLocations.push_back(y);
572 cylinderExchangeLocations.push_back(z);
573 cylinderExchangeLocations.push_back(maia::math::getSector(y, z, m_azimuthalAngle));
574 // }
575 }
576 }
577 }
578
579 if(noRANSExchangeCells > 0) {
580 MInt noRANSExchangeCells_ = noRANSExchangeCells;
581 // add ghostCells for better boundary interpolation
582 for(MInt i = 0; i < noRANSExchangeCells_; i++) {
583 MInt RANSId = cylinderExchangeIds[i];
584 MInt noNghbrIds = RANSSolver().a_noReconstructionNeighbors(RANSId);
585 for(MInt n = 0; n < noNghbrIds; n++) {
586 MInt nghbrId = RANSSolver().a_reconstructionNeighborId(RANSId, n);
587 if(RANSSolver().a_isBndryGhostCell(nghbrId)) {
588 MFloat x = RANSSolver().a_coordinate(nghbrId, 0);
589 MFloat y = RANSSolver().a_coordinate(nghbrId, 1);
590 MFloat z = RANSSolver().a_coordinate(nghbrId, 2);
591 noRANSExchangeCells++;
592 cylinderExchangeIds.push_back(nghbrId);
593 cylinderExchangeLocations.push_back(x);
594 cylinderExchangeLocations.push_back(y);
595 cylinderExchangeLocations.push_back(z);
596 cylinderExchangeLocations.push_back(maia::math::getSector(y, z, m_azimuthalAngle));
597 }
598 }
599 }
600 }
601
602 //--------------------------------------------------------------------------
603#ifdef _MB_DEBUG_
604 if(noRANSExchangeCells > 0) {
605 stringstream f;
606 f.clear();
607 f << "globalCylinderExchange_" << to_string(RANSSolver().domainId()) << ".txt";
608 MString fname = f.str();
609 MString header = "#id isHalo x1 x2 x3 sector";
610 ofstream data;
611 data.precision(16);
612 data.open(fname);
613 data << header << endl;
614 for(MInt i = 0; i < noRANSExchangeCells; i++) {
615 MString line = "";
616 line.append(to_string(cylinderExchangeIds[i]) + " "
617 + to_string(RANSSolver().a_isBndryGhostCell(cylinderExchangeIds[i])) + " "
618 + to_string(cylinderExchangeLocations[i * (nDim + 1) + 0]) + " "
619 + to_string(cylinderExchangeLocations[i * (nDim + 1) + 1]) + " "
620 + to_string(cylinderExchangeLocations[i * (nDim + 1) + 2]) + " "
621 + to_string(cylinderExchangeLocations[i * (nDim + 1) + 3]));
622 data << line << endl;
623 }
624 data.close();
625 }
626#endif
627 //--------------------------------------------------------------------------
628
629 // determine global limiting angles of RANS sector, necessary for finding equivalent solver cells for interpolation
630 MPI_Allreduce(MPI_IN_PLACE, &RANSSectors[0], 1, MPI_DOUBLE, MPI_MIN, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE",
631 "RANSSectors");
632 MPI_Allreduce(MPI_IN_PLACE, &RANSSectors[1], 1, MPI_DOUBLE, MPI_MAX, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE",
633 "RANSSectors");
634 MPI_Allreduce(MPI_IN_PLACE, &RANSSectorLimits[0], 1, MPI_DOUBLE, MPI_MIN, LESSolver().mpiComm(), AT_,
635 "MPI_IN_PLACE", "RANSSectorLimits");
636 MPI_Allreduce(MPI_IN_PLACE, &RANSSectorLimits[1], 1, MPI_DOUBLE, MPI_MAX, LESSolver().mpiComm(), AT_,
637 "MPI_IN_PLACE", "RANSSectorLimits");
638
639 m_RANSSectors[0] = RANSSectors[0];
640 m_RANSSectors[1] = RANSSectors[1];
641 m_RANSSectorLimits[0] = RANSSectorLimits[0];
642 m_RANSSectorLimits[1] = RANSSectorLimits[1];
643
644 m_noCylindricalGlobalExchangeLocations = 0;
645 m_noCylindricalGlobalRANSExchangeValues = 0;
646 m_noCylindricalGlobalLESExchangeValues = 0;
647 m_noCylindricalGlobalExchangeIds = 0;
648 m_noGlobalRANSExchangeCells = 0;
649
650 // set up global arrays for exchange
651 MPI_Allreduce(&noRANSExchangeCells, &m_noGlobalRANSExchangeCells, 1, MPI_INT, MPI_SUM, LESSolver().mpiComm(), AT_,
652 "noRANSExchangeCells", "m_noGlobalRANSExchangeCells");
653
654 m_noRANSCylinderExchangeVariables = noRANSVariables() + (MInt)m_STGSponge;
655 m_noCylindricalGlobalExchangeLocations = m_noGlobalRANSExchangeCells * (nDim + 1);
656 m_noCylindricalGlobalRANSExchangeValues = m_noGlobalRANSExchangeCells * m_noRANSCylinderExchangeVariables;
657 m_noCylindricalGlobalLESExchangeValues = m_noGlobalRANSExchangeCells * LESSolver().m_LESNoVarAverage;
658 m_noCylindricalGlobalExchangeIds = m_noGlobalRANSExchangeCells;
659
660 m_cylinderExchangeIdsOffset = -1;
661 m_noRANSExchangeCells = noRANSExchangeCells;
662 mAlloc(m_globalCylinderExchangeLocations, m_noCylindricalGlobalExchangeLocations,
663 "m_globalCylinderExchangeLocations", F0, FUN_);
664 mAlloc(m_globalCylinderRANSExchangeValues, m_noCylindricalGlobalRANSExchangeValues,
665 "m_globalCylinderRANSExchangeValues", F0, FUN_);
666 mAlloc(m_globalCylinderLESExchangeValues, m_noCylindricalGlobalLESExchangeValues,
667 "m_globalCylinderLESExchangeValues", F0, FUN_);
668 mAlloc(m_globalCylinderExchangeIds, m_noCylindricalGlobalExchangeIds, "m_globalCylinderExchangeIds", 0, FUN_);
669
670 // Create communicator for cylindic exchange communication
671 //--------------------------------------------------------------------------
672 MPI_Comm_size(LESSolver().mpiComm(), &m_commSizeCylExchange);
673 std::vector<MInt> cellsPerDomain(m_commSizeCylExchange);
674 std::vector<MInt> cylRanks(m_commSizeCylExchange);
675 MInt myCylRank = LESSolver().domainId();
676
677 // if(m_commSizeCylExchange > 1) {
678 MPI_Allgather(&noRANSExchangeCells, 1, MPI_INT, &cellsPerDomain[0], 1, MPI_INT, LESSolver().mpiComm(), AT_,
679 "noBc7909Cells ", "bc7909CellsperDomain");
680 MPI_Allgather(&myCylRank, 1, MPI_INT, &cylRanks[0], 1, MPI_INT, LESSolver().mpiComm(), AT_, "myCylRank",
681 "cylRanks");
682
683 MInt noInvolvedRanks = 0;
684 std::vector<MInt> involvedRanks(m_commSizeCylExchange);
685 for(MInt i = 0; i < m_commSizeCylExchange; i++) {
686 if(cellsPerDomain[i]) {
687 involvedRanks[noInvolvedRanks] = i;
688 ++noInvolvedRanks;
689 }
690 }
691
692 MPI_Comm commCyl;
693 MPI_Group group;
694 MPI_Group groupCyl;
695 MPI_Comm_group(LESSolver().mpiComm(), &group, AT_, "group");
696 MPI_Group_incl(group, noInvolvedRanks, &involvedRanks[0], &groupCyl, AT_);
697
698 MPI_Comm_create(LESSolver().mpiComm(), groupCyl, &commCyl, AT_, "commCyl");
699 m_commCyl = commCyl;
700
701 MInt myCommRank = -1;
702 for(MInt r = 0; r < noInvolvedRanks; r++) {
703 if(myCylRank == involvedRanks[r]) {
704 myCommRank = r;
705 break;
706 }
707 }
708
709 m_cylRoot = involvedRanks[0];
710
711 // defined here for global exchange for to determine m_cylinderExchangeIdsOffset
712 ScratchSpace<MInt> displsIds(noInvolvedRanks, "displsIds", FUN_);
713 //--------------------------------------------------------------------------
714
715 if(noRANSExchangeCells > 0) {
716 m_cylCommActive = true;
717
718 // exchange arrays
719 MInt noLocations = noRANSExchangeCells * (nDim + 1);
720 MInt noIds = noRANSExchangeCells;
721 if(noInvolvedRanks > 0) {
722 ScratchSpace<MInt> recvbufLocations(noInvolvedRanks, "recvbuf", FUN_);
723 ScratchSpace<MInt> recvbufIds(noInvolvedRanks, "recvbuf", FUN_);
724 recvbufLocations.fill(0);
725 // recvbufValues.fill(0);
726 recvbufIds.fill(0);
727
728 MPI_Gather(&noLocations, 1, MPI_INT, &recvbufLocations[0], 1, MPI_INT, 0, commCyl, AT_, "noLocations",
729 "recvbufLocations");
730 // MPI_Gather(&noValues, 1, MPI_INT, &recvbufValues[0], 1, MPI_INT, 0, commCyl,
731 // AT_, "noRANSExchangeCells", "recvbufValues");
732 MPI_Gather(&noIds, 1, MPI_INT, &recvbufIds[0], 1, MPI_INT, 0, commCyl, AT_, "noIds", "recvbufIds");
733
734 ScratchSpace<MInt> displsLocations(noInvolvedRanks, "displsLocations", FUN_);
735 if(myCommRank == 0) {
736 MInt offsetLocations = 0;
737 // MInt offsetValues = 0;
738 MInt offsetIds = 0;
739 for(MInt dom = 0; dom < noInvolvedRanks; dom++) {
740 displsLocations[dom] = offsetLocations;
741 displsIds[dom] = offsetIds;
742 offsetLocations += recvbufLocations[dom];
743 offsetIds += recvbufIds[dom];
744 }
745 }
746
747 MPI_Gatherv(&cylinderExchangeLocations[0], noLocations, MPI_DOUBLE, &m_globalCylinderExchangeLocations[0],
748 &recvbufLocations[myCommRank], &displsLocations[myCommRank], MPI_DOUBLE, 0, commCyl, AT_,
749 "cylinderExchangeLocations", "m_globalCylinderExchangeLocations");
750 MPI_Gatherv(&cylinderExchangeIds[0], noIds, MPI_INT, &m_globalCylinderExchangeIds[0], &recvbufIds[myCommRank],
751 &displsIds[myCommRank], MPI_INT, 0, commCyl, AT_, "cylinderExchangeIds",
752 "m_globalCylinderExchangeIds");
753
754 } else {
755 // JANNIK:is this correct, check this
756 // MInt index = 0;
757 // for(MInt RANSId = 0; RANSId < a_noFvGridCellsRANS(); RANSId++) {
758 // if(!RANSSolver().a_isHalo(RANSId) && RANSSolver().c_isLeafCell(RANSId)) {
759 // MInt LESId = convertIdParent(RANSSolver(), LESSolver(), RANSId);
760 // if(LESId != -1) {
761 // for(MInt pos = 0; pos < (nDim + 1); pos++) {
762 // m_globalCylinderExchangeLocations[index * (nDim + 1) + pos] =
763 // cylinderExchangeLocations[index * (nDim + 1) + pos];
764 // }
765 // m_globalCylinderExchangeIds[index] = cylinderExchangeIds[index];
766 // index++;
767 // }
768 // }
769 // }
770 }
771 }
772
773 MPI_Bcast(&displsIds[0], noInvolvedRanks, MPI_INT, m_cylRoot, LESSolver().mpiComm(), AT_, "displsIds");
774 for(MInt dom = 0; dom < noInvolvedRanks; dom++) {
775 if(dom == myCommRank) {
776 m_cylinderExchangeIdsOffset = displsIds[dom];
777 }
778 }
779
780 MPI_Allreduce(MPI_IN_PLACE, &m_cylRoot, 1, MPI_INT, MPI_MAX, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE",
781 "m_cylRoot");
782 MPI_Bcast(&m_globalCylinderExchangeLocations[0], m_noCylindricalGlobalExchangeLocations, MPI_DOUBLE, m_cylRoot,
783 LESSolver().mpiComm(), AT_, "m_globalCylinderExchangeLocations");
784 MPI_Bcast(&m_globalCylinderExchangeIds[0], m_noCylindricalGlobalExchangeIds, MPI_INT, m_cylRoot,
785 LESSolver().mpiComm(), AT_, "m_globalCylinderExchangeIds");
786
787 // debug output
788 //--------------------------------------------------------------------------
789#ifdef _MB_DEBUG_
790 if(LESSolver().domainId() == 0) {
791 stringstream fn2;
792 fn2.clear();
793 fn2 << "globalCylinderExchange.txt";
794 MString fname2 = fn2.str();
795 MString header2 = "#exchangeIndex id x1 x2 x3 sector";
796 ofstream data2;
797 data2.precision(16);
798 data2.open(fname2);
799 data2 << header2 << endl;
800 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
801 MString line = "";
802 line.append(to_string(exchangeIndex) + " " + to_string(m_globalCylinderExchangeIds[exchangeIndex]) + " "
803 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 0]) + " "
804 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 1]) + " "
805 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 2]) + " "
806 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 3]));
807 data2 << line << endl;
808 }
809 data2.close();
810 }
811#endif
812 //--------------------------------------------------------------------------
813
814 //
815 // Find equivalent LES cells to RANSExchange cells
816 //
817 mAlloc(m_globalCylinderInterpolationIndex, a_noFvGridCellsLES(), "m_globalCylinderInterpolationIndex", FUN_);
818 mAlloc(m_cylinderInterpolationAngle, a_noFvGridCellsLES(), "m_cylinderInterpolationAngle", FUN_);
819 mAlloc(m_globalCylinderInterpolationCell, m_noGlobalRANSExchangeCells, "m_globalCylinderInterpolationCell", FUN_);
820
821 for(MInt exchangeIndex = 0; exchangeIndex < m_noCylindricalGlobalLESExchangeValues; exchangeIndex++) {
822 m_globalCylinderInterpolationNumber.push_back(0);
823 }
824
825 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
826 m_globalCylinderInterpolationCell[exchangeIndex].clear();
827 }
828
829 MInt maxInterpCells = 4;
830 vector<vector<pair<MFloat, MInt>>> minDistIndex(
831 a_noFvGridCellsLES(),
832 vector<pair<MFloat, MInt>>(maxInterpCells, make_pair(std::numeric_limits<MFloat>::infinity(), -1)));
833 vector<vector<pair<MFloat, MInt>>> minDistCell(
834 m_noGlobalRANSExchangeCells,
835 vector<pair<MFloat, MInt>>(maxInterpCells, make_pair(std::numeric_limits<MFloat>::infinity(), -1)));
836
837 // find corresponding cell in exchangeValues
838 for(MInt LESId = 0; LESId < a_noFvGridCellsLES(); LESId++) {
839 m_globalCylinderInterpolationIndex[LESId].clear();
840 MFloat x = LESSolver().a_coordinate(LESId, 0);
841 MFloat y = LESSolver().a_coordinate(LESId, 1);
842 MFloat z = LESSolver().a_coordinate(LESId, 2);
843 MFloat cellHalfLength = F1B2 * LESSolver().c_cellLengthAtLevel(LESSolver().a_level(LESId));
844 MFloat sector = maia::math::getSector(y, z, m_azimuthalAngle);
845 MFloat angle = maia::math::getAngle(y, z);
846
847 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
848 MFloat xExchange = m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 0];
849 MFloat yExchange = m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 1];
850 MFloat zExchange = m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 2];
851
852 // fill data structure for calculation of LES sector average
853 for(MInt r = 0; r < 2; r++) { // JANNIK: generalize this
854 MFloat rotationAngle_ = -((sector - m_RANSSectors[r]) * m_azimuthalAngle) / 180.0 * PI;
855 MFloat yRotated_ = cos(rotationAngle_) * y - sin(rotationAngle_) * z;
856 MFloat zRotated_ = sin(rotationAngle_) * y + cos(rotationAngle_) * z;
857 if(approx(x, xExchange, F3 * cellHalfLength) && approx(yRotated_, yExchange, F3 * cellHalfLength)
858 && approx(zRotated_, zExchange, F3 * cellHalfLength)) {
859 // rotated cell is in vicinity of current cell
860 // save as interpolation cell
861 MInt RANSId_ = m_globalCylinderExchangeIds[exchangeIndex];
862 if(!RANSSolver().a_isBndryGhostCell(RANSId_) && LESSolver().c_isLeafCell(LESId)) {
863 m_globalCylinderInterpolationCell[exchangeIndex].push_back(
864 make_pair<MInt, MFloat>((MInt)LESId, (MFloat)rotationAngle_));
865 }
866 }
867 }
868
869 // fill data structure for exchange of sector RANS data to LES domain
870 MFloat rotationAngle = F0;
871 MFloat rotationAngle1 = -((sector - m_RANSSectors[0]) * m_azimuthalAngle) / 180.0 * PI;
872 MFloat rotationAngle2 = -((sector - m_RANSSectors[1]) * m_azimuthalAngle) / 180.0 * PI;
873 MFloat minAngleDist1 = F0;
874 MFloat minAngleDist2 = F0;
875 if(angle + rotationAngle1 > m_RANSSectorLimits[0] || angle + rotationAngle1 < m_RANSSectorLimits[1]) {
876 if(abs(angle + rotationAngle1 - m_RANSSectorLimits[0])
877 < abs(angle + rotationAngle1 - m_RANSSectorLimits[1])) {
878 minAngleDist1 = abs(angle + rotationAngle1 - m_RANSSectorLimits[0]);
879 } else {
880 minAngleDist1 = abs(angle + rotationAngle1 - m_RANSSectorLimits[1]);
881 }
882 }
883 if(angle + rotationAngle2 > m_RANSSectorLimits[0] || angle + rotationAngle2 < m_RANSSectorLimits[1]) {
884 if(abs(angle + rotationAngle2 - m_RANSSectorLimits[0])
885 < abs(angle + rotationAngle2 - m_RANSSectorLimits[1])) {
886 minAngleDist2 = abs(angle + rotationAngle2 - m_RANSSectorLimits[0]);
887 } else {
888 minAngleDist2 = abs(angle + rotationAngle2 - m_RANSSectorLimits[1]);
889 }
890 }
891 if(minAngleDist2 > minAngleDist1) {
892 rotationAngle = rotationAngle2;
893 } else {
894 rotationAngle = rotationAngle1;
895 }
896 MInt RANSId = convertId(LESSolver(), RANSSolver(), LESId);
897 // equivalent RANS cell exists, no rotation necessary
898 if(RANSId != -1) rotationAngle = 0;
899 m_cylinderInterpolationAngle[LESId] = rotationAngle;
900 MFloat yRotated_ = cos(rotationAngle) * y - sin(rotationAngle) * z;
901 MFloat zRotated_ = sin(rotationAngle) * y + cos(rotationAngle) * z;
902 if(approx(x, xExchange, F3 * cellHalfLength) && approx(yRotated_, yExchange, F4 * cellHalfLength)
903 && approx(zRotated_, zExchange, F4 * cellHalfLength)) {
904 MFloat delta_y = yRotated_ - yExchange;
905 MFloat delta_z = zRotated_ - zExchange;
906 MFloat dist = sqrt(POW2(delta_y) + POW2(delta_z));
907
908 MBool checkWindowHalo = false;
909 if(!checkWindowHalo) {
910 if(dist < minDistIndex[LESId][maxInterpCells - 1].first) {
911 minDistIndex[LESId][maxInterpCells - 1] = make_pair(dist, exchangeIndex);
912 // sort by dist
913 sort(minDistIndex[LESId].begin(), minDistIndex[LESId].end());
914 }
915 }
916 }
917 }
918 }
919
920 for(MInt LESId = 0; LESId < a_noFvGridCellsLES(); LESId++) {
921 MFloat minY = std::numeric_limits<MFloat>::infinity();
922 MFloat maxY = -std::numeric_limits<MFloat>::infinity();
923 MFloat minZ = std::numeric_limits<MFloat>::infinity();
924 MFloat maxZ = -std::numeric_limits<MFloat>::infinity();
925 MFloat rotationAngle = m_cylinderInterpolationAngle[LESId];
926 MFloat* coordLES = &LESSolver().a_coordinate(LESId, 0);
927 MFloat yRotated = cos(rotationAngle) * coordLES[1] - sin(rotationAngle) * coordLES[2];
928 MFloat zRotated = sin(rotationAngle) * coordLES[1] + cos(rotationAngle) * coordLES[2];
929 // if(minDistIndex[LESId][0].second != -1 && minDistIndex[LESId][0].first < 1e-16){
930 // MInt exchangeIndex = minDistIndex[LESId][0].second;
931 // m_globalCylinderInterpolationIndex[LESId].push_back(exchangeIndex);
932 // } else {
933 for(MInt i = 0; i < maxInterpCells; i++) {
934 if(minDistIndex[LESId][i].second != -1) {
935 MInt exchangeIndex = minDistIndex[LESId][i].second;
936 MBool usePoint = true;
937 MFloat* coordRANS = &m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1)];
938 if(i == maxInterpCells - 1) {
939 if(yRotated > minY && yRotated < maxY && zRotated > minZ && zRotated < maxZ) usePoint = false;
940 }
941 if(coordRANS[1] > yRotated && coordRANS[1] > maxY) maxY = coordRANS[1];
942 if(coordRANS[1] < yRotated && coordRANS[1] < minY) minY = coordRANS[1];
943 if(coordRANS[2] > zRotated && coordRANS[2] > maxZ) maxZ = coordRANS[2];
944 if(coordRANS[2] < zRotated && coordRANS[2] < minZ) minZ = coordRANS[2];
945 if(usePoint) m_globalCylinderInterpolationIndex[LESId].push_back(exchangeIndex);
946 }
947 }
948
949 // no proper 2D interpolation can be performed, use nearest neighbor extrapolation
950 if(yRotated < minY || yRotated > maxY || zRotated < minZ || zRotated > maxZ) {
951 if(m_globalCylinderInterpolationIndex[LESId].size() > 1) {
952 m_globalCylinderInterpolationIndex[LESId].erase(m_globalCylinderInterpolationIndex[LESId].begin() + 1,
953 m_globalCylinderInterpolationIndex[LESId].end());
954 }
955 }
956 // }
957 }
958
959 // delete interpolation cells if none is not a halo cell
960 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
961 MInt testForOnlyHalo = 0;
962 for(MInt i = 0; i < (MInt)m_globalCylinderInterpolationCell[exchangeIndex].size(); i++) {
963 MInt cellId = m_globalCylinderInterpolationCell[exchangeIndex][i].first;
964 testForOnlyHalo += (MInt)LESSolver().a_isHalo(cellId);
965 }
966
967 if(testForOnlyHalo == (MInt)m_globalCylinderInterpolationCell[exchangeIndex].size()) {
968 m_globalCylinderInterpolationCell[exchangeIndex].clear();
969 }
970 }
971
972 // debug
973 //--------------------------------------------------------------------------
974#ifdef _MB_DEBUG_
975 stringstream fn3;
976 fn3.clear();
977 fn3 << "globalCylinderLESCell_" << to_string(LESSolver().domainId()) << ".txt";
978 MString fname3 = fn3.str();
979 MString header3 = "#exchangeIndex ...";
980 ofstream data3;
981 data3.precision(16);
982 data3.open(fname3);
983 data3 << header3 << endl;
984 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
985 MString line = to_string(exchangeIndex) + " | ";
986 for(MInt i = 0; i < (MInt)m_globalCylinderInterpolationCell[exchangeIndex].size(); i++) {
987 MInt cellId = m_globalCylinderInterpolationCell[exchangeIndex][i].first;
988 MInt rotationAngle = m_globalCylinderInterpolationCell[exchangeIndex][i].second;
989 line.append("(" + to_string(cellId) + " " + to_string(rotationAngle) + " "
990 + to_string(LESSolver().a_isHalo(cellId)) + " ," + to_string(LESSolver().a_coordinate(cellId, 0))
991 + " " + to_string(LESSolver().a_coordinate(cellId, 1)) + " "
992 + to_string(LESSolver().a_coordinate(cellId, 2)) + ") ");
993 }
994 data3 << line << endl;
995 }
996 data3.close();
997
998 stringstream fn4;
999 fn4.clear();
1000 fn4 << "globalCylinderLESIndex_" << to_string(LESSolver().domainId()) << ".txt";
1001 MString fname4 = fn4.str();
1002 MString header4 = "#LESId exchangeIndex ";
1003 ofstream data4;
1004 data4.precision(16);
1005 data4.open(fname4);
1006 data4 << header4 << endl;
1007 for(MInt LESId = 0; LESId < a_noFvGridCellsLES(); LESId++) {
1008 if(m_globalCylinderInterpolationIndex[LESId].empty()) continue;
1009 MString line =
1010 to_string(LESId) + " " + to_string(LESSolver().a_isHalo(LESId)) + " "
1011 + to_string(LESSolver().c_isLeafCell(LESId)) + " | " + to_string(LESSolver().a_coordinate(LESId, 0)) + " "
1012 + to_string(LESSolver().a_coordinate(LESId, 1)) + " " + to_string(LESSolver().a_coordinate(LESId, 2)) + " ";
1013 for(MInt i = 0; i < (MInt)m_globalCylinderInterpolationIndex[LESId].size(); i++) {
1014 MInt exchangeIndex = m_globalCylinderInterpolationIndex[LESId][i];
1015 line.append("(" + to_string(exchangeIndex) + ", "
1016 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 0]) + " "
1017 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 1]) + " "
1018 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 2]) + " "
1019 + to_string(m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 3]) + ") ");
1020 }
1021 data4 << line << endl;
1022 }
1023 data4.close();
1024#endif
1025 //--------------------------------------------------------------------------
1026 }
1027}
1028
1029
1033template <MInt nDim, class SysEqn>
1035 TRACE();
1036
1037 if(LESSolver().grid().isActive()) {
1038 // interpolation
1039 for(MInt LESId = 0; LESId < a_noFvGridCellsLES(); LESId++) {
1040 // write value to cell
1041 if(LESSolver().a_isHalo(LESId)) continue;
1042 MFloat y = LESSolver().a_coordinate(LESId, 1);
1043 MFloat z = LESSolver().a_coordinate(LESId, 2);
1044 MFloat rotationAngle = m_cylinderInterpolationAngle[LESId];
1045 MFloat yRotated_ = cos(rotationAngle) * y - sin(rotationAngle) * z;
1046 MFloat zRotated_ = sin(rotationAngle) * y + cos(rotationAngle) * z;
1047
1048 if(m_globalCylinderInterpolationIndex[LESId].empty()) continue;
1049
1050 for(MInt var = 0; var < noRANSVariables(); var++) {
1051 ASSERT(LESId < (MInt)LESSolver().m_RANSValues[var].size(),
1052 "Trying to access data [" + to_string(var) + "][" + to_string(LESId) + "] in m_RANSValues with length "
1053 + to_string(LESSolver().m_RANSValues[var].size())
1054 + ", domainId: " + to_string(LESSolver().domainId()));
1055
1056 if(m_globalCylinderInterpolationIndex[LESId].size() == 1) {
1057 MInt interpolationIndex = m_globalCylinderInterpolationIndex[LESId][0];
1058 LESSolver().m_RANSValues[var][LESId] =
1059 m_globalCylinderRANSExchangeValues[interpolationIndex * m_noRANSCylinderExchangeVariables + var];
1060 }
1061
1062 if(m_globalCylinderInterpolationIndex[LESId].size() > 1) {
1063 MInt interpolationIndex = m_globalCylinderInterpolationIndex[LESId][0];
1064 MFloat yExchange = m_globalCylinderExchangeLocations[interpolationIndex * (nDim + 1) + 1];
1065 MFloat zExchange = m_globalCylinderExchangeLocations[interpolationIndex * (nDim + 1) + 2];
1066 MFloat delta_y = yRotated_ - yExchange;
1067 MFloat delta_z = zRotated_ - zExchange;
1068
1069 // no interpolation necessary, since equivalent RANS cell exists
1070 if(sqrt(pow(delta_y, 2) + pow(delta_z, 2)) < 1e-16 /*halfCellLength/20.0*/) {
1071 LESSolver().m_RANSValues[var][LESId] =
1072 m_globalCylinderRANSExchangeValues[interpolationIndex * m_noRANSCylinderExchangeVariables + var];
1073
1074 // Interpolation - Least square method
1075 } else {
1076 MFloat n = F0;
1077 MFloat Eyi = F0;
1078 MFloat Eyi2 = F0;
1079 MFloat Ezi = F0;
1080 MFloat Ezi2 = F0;
1081 MFloat Eyizi = F0;
1082 MFloat Eui = F0;
1083 MFloat Euiyi = F0;
1084 MFloat Euizi = F0;
1085
1086 for(MInt i = 0; i < (MInt)m_globalCylinderInterpolationIndex[LESId].size(); i++) {
1087 interpolationIndex = m_globalCylinderInterpolationIndex[LESId][i];
1088 yExchange = m_globalCylinderExchangeLocations[interpolationIndex * (nDim + 1) + 1];
1089 zExchange = m_globalCylinderExchangeLocations[interpolationIndex * (nDim + 1) + 2];
1090 delta_y = yRotated_ - yExchange;
1091 delta_z = zRotated_ - zExchange;
1092 n += 1;
1093 Eyi += delta_y;
1094 Eyi2 += POW2(delta_y);
1095 Ezi += delta_z;
1096 Ezi2 += POW2(delta_z);
1097 Eyizi += delta_y * delta_z;
1098 Eui += m_globalCylinderRANSExchangeValues[interpolationIndex * m_noRANSCylinderExchangeVariables + var];
1099 Euiyi += m_globalCylinderRANSExchangeValues[interpolationIndex * m_noRANSCylinderExchangeVariables + var]
1100 * delta_y;
1101 Euizi += m_globalCylinderRANSExchangeValues[interpolationIndex * m_noRANSCylinderExchangeVariables + var]
1102 * delta_z;
1103 }
1104
1105 MFloat detA =
1106 n * Eyi2 * Ezi2 + 2 * Eyi * Eyizi * Ezi - (Ezi * Eyi2 * Ezi + n * Eyizi * Eyizi + Ezi2 * Eyi * Eyi);
1107
1108 if(detA > 1e-16) {
1109 LESSolver().m_RANSValues[var][LESId] =
1110 1 / detA
1111 * (Eui * (Eyi2 * Ezi2 - Eyizi * Eyizi) + Euiyi * (Ezi * Eyizi - Eyi * Ezi2)
1112 + Euizi * (Eyi * Eyizi - Eyi2 * Ezi));
1113 }
1114 }
1115 }
1116 }
1117
1118 // rotate veloctiy components
1119 MInt vId = LESSolver().sysEqn().PV->V;
1120 MInt wId = LESSolver().sysEqn().PV->W;
1121 MFloat v_ = LESSolver().m_RANSValues[vId][LESId];
1122 MFloat w_ = LESSolver().m_RANSValues[wId][LESId];
1123 LESSolver().m_RANSValues[vId][LESId] = v_ * cos(-rotationAngle) - w_ * sin(-rotationAngle);
1124 LESSolver().m_RANSValues[wId][LESId] = v_ * sin(-rotationAngle) + w_ * cos(-rotationAngle);
1125 }
1126 }
1127}
1128
1129
1133template <MInt nDim, class SysEqn>
1135 TRACE();
1136
1137 if(LESSolver().grid().isActive()) {
1138 ScratchSpace<MInt> numberCylinderExchangeValues(m_noCylindricalGlobalRANSExchangeValues,
1139 "numberCylinderExchangeValues", FUN_);
1140 numberCylinderExchangeValues.fill(0);
1141
1142 if(m_cylCommActive) {
1143 // reset vector
1144 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
1145 for(MInt var = 0; var < m_noRANSCylinderExchangeVariables; var++) {
1146 m_globalCylinderRANSExchangeValues[exchangeIndex * m_noRANSCylinderExchangeVariables + var] = F0;
1147 }
1148 }
1149
1150 // fill vector
1151 for(MInt exchangeIndex_ = 0; exchangeIndex_ < m_noRANSExchangeCells; exchangeIndex_++) {
1152 MInt exchangeIndex = exchangeIndex_ + m_cylinderExchangeIdsOffset;
1153 MInt RANSId = m_globalCylinderExchangeIds[exchangeIndex];
1154 if(RANSId > -1) {
1155 for(MInt var = 0; var < noRANSVariables(); var++) {
1156 m_globalCylinderRANSExchangeValues[exchangeIndex * m_noRANSCylinderExchangeVariables + var] =
1157 RANSSolver().a_pvariable(RANSId, var);
1158 numberCylinderExchangeValues[exchangeIndex * m_noRANSCylinderExchangeVariables + var] = 1;
1159 }
1160 }
1161 }
1162
1163 MPI_Allreduce(MPI_IN_PLACE, &m_globalCylinderRANSExchangeValues[0], m_noCylindricalGlobalRANSExchangeValues,
1164 MPI_DOUBLE, MPI_SUM, m_commCyl, AT_, "MPI_IN_PLACE", "m_globalCylinderRANSExchangeValues");
1165 MPI_Allreduce(MPI_IN_PLACE, &numberCylinderExchangeValues[0], m_noCylindricalGlobalRANSExchangeValues, MPI_INT,
1166 MPI_SUM, m_commCyl, AT_, "MPI_IN_PLACE", "numberCylinderExchangeValues");
1167 }
1168
1169 MPI_Bcast(&m_globalCylinderRANSExchangeValues[0], m_noCylindricalGlobalRANSExchangeValues, MPI_DOUBLE, m_cylRoot,
1170 LESSolver().mpiComm(), AT_, "m_globalCylinderRANSExchangeValues");
1171 MPI_Bcast(&numberCylinderExchangeValues[0], m_noCylindricalGlobalRANSExchangeValues, MPI_INT, m_cylRoot,
1172 LESSolver().mpiComm(), AT_, "numberCylinderExchangeValues");
1173
1174
1175 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
1176 for(MInt var = 0; var < noRANSVariables(); var++) {
1177 if(numberCylinderExchangeValues[exchangeIndex * m_noRANSCylinderExchangeVariables + var] > 0) {
1178 m_globalCylinderRANSExchangeValues[exchangeIndex * m_noRANSCylinderExchangeVariables + var] /=
1179 numberCylinderExchangeValues[exchangeIndex * m_noRANSCylinderExchangeVariables + var];
1180 }
1181 }
1182 }
1183 }
1184}
1185
1186
1190template <MInt nDim, class SysEqn>
1192 TRACE();
1193
1194 if(LESSolver().grid().isActive()) {
1195 MInt vId = LESSolver().sysEqn().PV->V;
1196 MInt wId = LESSolver().sysEqn().PV->W;
1197
1198 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
1199 if((MInt)m_globalCylinderInterpolationCell[exchangeIndex].size() == 0) continue;
1200
1201 // No interpolation necessary
1202 if((MInt)m_globalCylinderInterpolationCell[exchangeIndex].size() == 1) {
1203 MInt cellId = m_globalCylinderInterpolationCell[exchangeIndex][0].first;
1204
1205 MFloat rotationAngle = m_globalCylinderInterpolationCell[exchangeIndex][0].second;
1206
1207 vector<MInt>::iterator findAvgId =
1208 find(LESSolver().m_LESAverageCells.begin(), LESSolver().m_LESAverageCells.end(), cellId);
1209 if(findAvgId != LESSolver().m_LESAverageCells.end()) {
1210 MInt avgId = distance(LESSolver().m_LESAverageCells.begin(), findAvgId);
1211 for(MInt var = 0; var < LESSolver().m_LESNoVarAverage; var++) {
1212 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + var] =
1213 LESSolver().m_LESVarAverage[var][avgId];
1214 m_globalCylinderInterpolationNumber[exchangeIndex * LESSolver().m_LESNoVarAverage + var] = 1;
1215 }
1216
1217 // rotate veloctiy components
1218 MFloat v_ = LESSolver().m_LESVarAverage[vId][avgId];
1219 MFloat w_ = LESSolver().m_LESVarAverage[wId][avgId];
1220 MFloat v_rot = v_ * cos(rotationAngle) - w_ * sin(rotationAngle);
1221 MFloat w_rot = v_ * sin(rotationAngle) + w_ * cos(rotationAngle);
1222 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + vId] = v_rot;
1223 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + wId] = w_rot;
1224 }
1225 // continue;
1226 }
1227
1228 // Interpolation - Least square method
1229 if(m_globalCylinderInterpolationCell[exchangeIndex].size() > 1) {
1230 for(MInt var = 0; var < LESSolver().m_LESNoVarAverage; var++) {
1231 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + var] = F0;
1232
1233 MFloat n = F0;
1234 MFloat Eyi = F0;
1235 MFloat Eyi2 = F0;
1236 MFloat Ezi = F0;
1237 MFloat Ezi2 = F0;
1238 MFloat Eyizi = F0;
1239 MFloat Eui = F0;
1240 MFloat Euiyi = F0;
1241 MFloat Euizi = F0;
1242
1243 for(MInt i = 0; i < (MInt)m_globalCylinderInterpolationCell[exchangeIndex].size(); i++) {
1244 MInt cellId = m_globalCylinderInterpolationCell[exchangeIndex][i].first;
1245 MFloat y = LESSolver().a_coordinate(cellId, 1);
1246 MFloat z = LESSolver().a_coordinate(cellId, 2);
1247 MFloat rotationAngle = m_globalCylinderInterpolationCell[exchangeIndex][i].second;
1248 MFloat yRotated_ = cos(rotationAngle) * y - sin(rotationAngle) * z;
1249 MFloat zRotated_ = sin(rotationAngle) * y + cos(rotationAngle) * z;
1250 MFloat yExchange = m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 1];
1251 MFloat zExchange = m_globalCylinderExchangeLocations[exchangeIndex * (nDim + 1) + 2];
1252 MFloat delta_y = yExchange - yRotated_;
1253 MFloat delta_z = zExchange - zRotated_;
1254
1255 vector<MInt>::iterator findAvgId =
1256 find(LESSolver().m_LESAverageCells.begin(), LESSolver().m_LESAverageCells.end(), cellId);
1257 if(findAvgId != LESSolver().m_LESAverageCells.end()) {
1258 MInt avgId = distance(LESSolver().m_LESAverageCells.begin(), findAvgId);
1259
1260 // rotate veloctiy components
1261 MFloat v_ = LESSolver().m_LESVarAverage[vId][avgId];
1262 MFloat w_ = LESSolver().m_LESVarAverage[wId][avgId];
1263 MFloat v_rot = v_ * cos(rotationAngle) - w_ * sin(rotationAngle);
1264 MFloat w_rot = v_ * sin(rotationAngle) + w_ * cos(rotationAngle);
1265
1266 n += 1;
1267 Eyi += delta_y;
1268 Eyi2 += POW2(delta_y);
1269 Ezi += delta_z;
1270 Ezi2 += POW2(delta_z);
1271 Eyizi += delta_y * delta_z;
1272
1273 if(var == vId) {
1274 Eui += v_rot;
1275 Euiyi += v_rot * delta_y;
1276 Euizi += v_rot * delta_z;
1277 } else if(var == wId) {
1278 Eui += w_rot;
1279 Euiyi += w_rot * delta_y;
1280 Euizi += w_rot * delta_z;
1281 } else {
1282 Eui += LESSolver().m_LESVarAverage[var][avgId];
1283 Euiyi += LESSolver().m_LESVarAverage[var][avgId] * delta_y;
1284 Euizi += LESSolver().m_LESVarAverage[var][avgId] * delta_z;
1285 }
1286 }
1287 }
1288
1289 MFloat detA =
1290 n * Eyi2 * Ezi2 + 2 * Eyi * Eyizi * Ezi - (Ezi * Eyi2 * Ezi + n * Eyizi * Eyizi + Ezi2 * Eyi * Eyi);
1291
1292 if(detA >= 1e-12) {
1293 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + var] =
1294 F1 / detA
1295 * (Eui * (Eyi2 * Ezi2 - Eyizi * Eyizi) + Euiyi * (Ezi * Eyizi - Eyi * Ezi2)
1296 + Euizi * (Eyi * Eyizi - Eyi2 * Ezi));
1297 }
1298
1299 if(!approx(m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + var], F0,
1300 1e-12)) {
1301 m_globalCylinderInterpolationNumber[exchangeIndex * LESSolver().m_LESNoVarAverage + var] = 1;
1302 }
1303 }
1304 }
1305 }
1306
1307 MPI_Allreduce(MPI_IN_PLACE, &m_globalCylinderInterpolationNumber[0], m_noCylindricalGlobalLESExchangeValues,
1308 MPI_INT, MPI_SUM, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE", "m_globalCylinderInterpolationNumber");
1309
1310 MPI_Allreduce(MPI_IN_PLACE, &m_globalCylinderLESExchangeValues[0], m_noCylindricalGlobalLESExchangeValues,
1311 MPI_DOUBLE, MPI_SUM, LESSolver().mpiComm(), AT_, "MPI_IN_PLACE", "m_globalCylinderLESExchangeValues");
1312
1313
1314 for(MInt exchangeIndex = 0; exchangeIndex < m_noGlobalRANSExchangeCells; exchangeIndex++) {
1315 for(MInt var = 0; var < LESSolver().m_LESNoVarAverage; var++) {
1316 if(m_globalCylinderInterpolationNumber[exchangeIndex * LESSolver().m_LESNoVarAverage + var] > 0) {
1317 m_globalCylinderLESExchangeValues[exchangeIndex * LESSolver().m_LESNoVarAverage + var] /=
1318 m_globalCylinderInterpolationNumber[exchangeIndex * LESSolver().m_LESNoVarAverage + var];
1319 }
1320 }
1321 }
1322 }
1323}
1324
1325
void mAlloc(T *&a, const MLong N, const MString &objectName, MString function)
allocates memory for one-dimensional array 'a' of size N
Definition: alloc.h:173
static MInt propertyLength(const MString &name, MInt solverId=m_noSolvers)
Returns the number of elements of a property.
Definition: context.cpp:538
static MBool propertyExists(const MString &name, MInt solver=m_noSolvers)
This function checks if a property exists in general.
Definition: context.cpp:494
void initCylinderExchange()
Create mapping for LES cells to equivalent RANS cells.
Definition: fvzonalstg.cpp:525
void balancePost() override
Definition: fvzonalstg.cpp:100
void determineZonalPositions()
Definition: fvzonalstg.cpp:124
FvZonalSTG(const MInt couplingId, RANS *R, LES *L)
void calcRANSSectorValues()
fill m_globalCylinderRANSExchangeValues with RANS values for exchange with LES Solver
void preCouple(MInt)
Definition: fvzonalstg.cpp:108
void finalizeCouplerInit()
Definition: fvzonalstg.cpp:73
void transferSolverData()
Definition: fvzonalstg.cpp:267
void finalizeAdaptation(const MInt solverId) override
Definition: fvzonalstg.cpp:90
void init()
Definition: fvzonalstg.cpp:47
void cylinderExchange()
Interpolate RANS sector data to LES domain using mapping created in initCylinderExchange.
void calcLESSectorAverage()
fill m_globalCylinderLESExchangeValues with sector averaged LES values for exchange with RANS Solver
void transferSpongeData()
Definition: fvzonalstg.cpp:384
This class is a ScratchSpace.
Definition: scratch.h:758
void fill(T val)
fill the scratch with a given value
Definition: scratch.h:311
MInt convertId(SolverA &solverA, SolverB &solverB, const MInt solverAId)
Conversion from solverA id to the solverB id on the same-level only!
Definition: couplingutils.h:21
MInt convertIdParent(SolverA &solverA, SolverB &solverB, const MInt solverAId)
Conversion from solverA id to the solverB id If no cell on the same level is found,...
Definition: couplingutils.h:46
constexpr Real POW2(const Real x)
Definition: functions.h:119
MBool approx(const T &, const U &, const T)
Definition: functions.h:272
MInt globalTimeStep
InfoOutFile m_log
int32_t MInt
Definition: maiatypes.h:62
std::basic_string< char > MString
Definition: maiatypes.h:55
double MFloat
Definition: maiatypes.h:52
bool MBool
Definition: maiatypes.h:58
int MPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm *newcomm, const MString &name, const MString &varname)
same as MPI_Comm_create, but updates the number of MPI communicators
int MPI_Group_incl(MPI_Group group, int n, const int ranks[], MPI_Group *newgroup, const MString &name)
same as MPI_Group_incl
int MPI_Comm_group(MPI_Comm comm, MPI_Group *group, const MString &name, const MString &varname)
same as MPI_Comm_group
int MPI_Gatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int recvcounts[], const int displs[], MPI_Datatype recvtype, int root, MPI_Comm comm, const MString &name, const MString &sndvarname, const MString &rcvvarname)
same as MPI_Gatherv
int MPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, const MString &name, const MString &sndvarname, const MString &rcvvarname)
same as MPI_Allreduce
int MPI_Bcast(void *buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm, const MString &name, const MString &varname)
same as MPI_Bcast
int MPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm, const MString &name, const MString &sndvarname, const MString &rcvvarname)
same as MPI_Gather
int MPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm, const MString &name, const MString &sndvarname, const MString &rcvvarname)
same as MPI_Allgather
MFloat getSector(MFloat y, MFloat z, MFloat azimuthalAngle)
Definition: maiamath.h:912
MFloat getAngle(MFloat y, MFloat z)
Definition: maiamath.h:925
MFloat dist(const Point< DIM > &p, const Point< DIM > &q)
Definition: pointbox.h:54
define array structures