MaMiCo 1.2
Loading...
Searching...
No Matches
CouetteScenario.h
1// This file is part of the Mamico project. For conditions of distribution
2// and use, please see the copyright notice in Mamico's main folder, or at
3// www5.in.tum.de/mamico
4#ifndef _COUPLING_SCENARIO_COUETTESCENARIO_H_
5#define _COUPLING_SCENARIO_COUETTESCENARIO_H_
6
7#include "coupling/ErrorEstimation.h"
8#include "coupling/InstanceHandling.h"
9#include "coupling/MultiMDMediator.h"
10#include "coupling/services/ParallelTimeIntegrationService.h"
11#include "coupling/solvers/CouetteSolver.h"
12#include "coupling/solvers/LBCouetteSolver.h"
13#include "simplemd/configurations/MolecularDynamicsConfiguration.h"
14#include "tarch/configuration/ParseConfiguration.h"
15#include "tarch/utils/MultiMDService.h"
16#include "tarch/utils/RandomNumberService.h"
17#include "tarch/utils/Utils.h"
18#if (BUILD_WITH_OPENFOAM)
19#include "coupling/solvers/IcoFoamInterface.h"
20#endif
21#include "coupling/configurations/CouetteConfiguration.h"
22#include "coupling/configurations/MaMiCoConfiguration.h"
23#include "coupling/indexing/IndexingService.h"
24#include "coupling/interface/MDSimulationFactory.h"
25#include "coupling/interface/impl/SimpleMD/SimpleMDSolverInterface.h"
26#include "coupling/services/MultiMDCellService.h"
27#include "coupling/solvers/CouetteSolverInterface.h"
28#include "coupling/solvers/FDCouetteSolver.h"
29#include "coupling/solvers/LBCouetteSolverInterface.h"
30
31#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
32#include <mpi.h>
33#endif
34#include <chrono>
35#include <math.h>
36#include <random>
37#include <sys/time.h>
38
39#if defined(LS1_MARDYN)
40#include "coupling/interface/impl/ls1/LS1MDSolverInterface.h"
41#include "coupling/interface/impl/ls1/LS1StaticCommData.h"
42#include "utils/Logger.h"
43using Log::global_log;
44#endif
45
46// This is ignored if you dont use synthetic MD. For further instructions cf.
47// SYNTHETIC part of initSolvers().
48#define SYNTHETICMD_SEQUENCE "SYNTHETIC-MD"
49
69class CouetteScenario : public Scenario {
70public:
72 CouetteScenario() : Scenario("CouetteScenario"), _generator(std::chrono::system_clock::now().time_since_epoch().count()) {}
74 virtual ~CouetteScenario() {}
75
78 void run() override {
79 init();
80 if (_cfg.twsLoop) {
81 twsLoop();
82 return;
83 }
84 _timeIntegrationService->run(_cfg.couplingCycles);
85 shutdown();
86 }
87
90 void init() override {
91#if defined(LS1_MARDYN)
92 Log::global_log = std::make_unique<Log::Logger>(Log::Error); // Log::Info
93#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
94 global_log->set_mpi_output_root(0);
95#endif
96#endif
99 }
100
106 void runOneCouplingCycle(int cycle) override {
107 advanceMacro(cycle);
108 varyMD(cycle);
109 advanceMicro(cycle);
110 computeSNR(cycle);
111 twoWayCoupling(cycle);
112 if (_cfg.totalNumberMDSimulations < 0) // dynamic MD
113 _multiMDCellService->finishCycle(cycle, *_instanceHandling);
114 if (_isRootRank) {
115 // Output status info only every 10 seconds
116 gettimeofday(&_tv.end, NULL);
117 int runtime_ms = (int)(((_tv.end.tv_sec - _tv.output.tv_sec) * 1000000 + (_tv.end.tv_usec - _tv.output.tv_usec)) / 1000);
118 if (runtime_ms > 10000) {
119 std::cout << "Finish coupling cycle " << cycle << std::endl;
120 gettimeofday(&_tv.output, NULL);
121 }
122 }
123 }
124
127 void equilibrateMicro() override {
129 return;
130
132 for (auto pair : _couplingBuffer.macro2MDBuffer)
133 buffer << pair;
134 for (auto pair : _couplingBuffer.md2macroBuffer)
135 buffer << pair;
136 fillSendBuffer(_cfg.density, *_couetteSolver, buffer);
137
138 _multiMDCellService->setInnerMomentumImposition(true);
139 const int EQUI_CYCLES = 3; // 2 cycles should be sufficient (10 cycles without inner imposition on MD30), 3 cycles is even safer
140 for (int i = 0; i < EQUI_CYCLES; i++) {
141 _multiMDCellService->sendFromMacro2MD(buffer);
142 _instanceHandling->simulateTimesteps(_simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps(), _mdStepCounter, *_multiMDCellService);
143 }
144 _mdStepCounter += EQUI_CYCLES * _simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps();
145
146 _multiMDCellService->setInnerMomentumImposition(false);
147 }
148
149 coupling::solvers::AbstractCouetteSolver<3>* getSolver() override { return _couetteSolver; }
150
151protected:
153 void twsLoop() {
154 for (_tws = _cfg.twsLoopMin; _tws <= _cfg.twsLoopMax; _tws += _cfg.twsLoopStep) {
155 init();
156 for (int cycle = 0; cycle < _cfg.couplingCycles; cycle++)
157 runOneCouplingCycle(cycle);
158 shutdown();
159 }
160 }
161
165 std::string filename("couette.xml");
166
169 if (!_simpleMDConfig.isValid()) {
170 if (_isRootRank) {
171 std::cout << "ERROR CouetteScenario: Invalid SimpleMD config!" << std::endl;
172 }
173 exit(EXIT_FAILURE);
174 }
176 if (!_mamicoConfig.isValid()) {
177 if (_isRootRank) {
178 std::cout << "ERROR CouetteScenario: Invalid MaMiCo config!" << std::endl;
179 }
180 exit(EXIT_FAILURE);
181 }
182
184
186 _simpleMDConfig.getDomainDecompConfiguration().getDecompType() ==
187 simplemd::configurations::DomainDecompConfiguration::DecompositionType::STATIC_IRREG_RECT_GRID) {
188 std::cout << "ERROR Currently, only LS1 supports irregular rectilinear domain decomposition! Please change md type to ls1, or decomposition to default!"
189 << std::endl;
190 exit(EXIT_FAILURE);
191 }
192
193#if defined(LS1_MARDYN)
194 auto offset = _simpleMDConfig.getDomainConfiguration().getGlobalDomainOffset();
195 coupling::interface::LS1StaticCommData::getInstance().setConfigFilename("ls1config.xml");
196 coupling::interface::LS1StaticCommData::getInstance().setBoxOffsetAtDim(0, offset[0]); // temporary till ls1 offset is natively supported
197 coupling::interface::LS1StaticCommData::getInstance().setBoxOffsetAtDim(1, offset[1]);
198 coupling::interface::LS1StaticCommData::getInstance().setBoxOffsetAtDim(2, offset[2]);
199#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
200 auto subdomainWeights = _simpleMDConfig.getDomainDecompConfiguration().getSubdomainWeights();
201 coupling::interface::LS1StaticCommData::getInstance().setSubdomainWeights(subdomainWeights);
202#endif
203#endif
204 }
205
208 void initSolvers() {
209 _timeIntegrationService = std::make_unique<coupling::services::ParallelTimeIntegrationService<3>>(_mamicoConfig, this);
210 _rank = _timeIntegrationService->getRank(); // returns the rank inside local time domain
211
213 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[0], _simpleMDConfig.getDomainConfiguration().getGlobalDomainOffset(),
214 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize(), getGlobalNumberCouplingCells(_simpleMDConfig, _mamicoConfig),
215 _mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap());
216
217 // init indexing
218 if (_simpleMDConfig.getDomainDecompConfiguration().getDecompType() ==
219 simplemd::configurations::DomainDecompConfiguration::DecompositionType::STATIC_IRREG_RECT_GRID) {
220 coupling::indexing::IndexingService<3>::getInstance().initWithMDSize(
221 _simpleMDConfig.getDomainDecompConfiguration().getSubdomainWeights(), _simpleMDConfig.getDomainConfiguration().getGlobalDomainSize(),
222 _simpleMDConfig.getDomainConfiguration().getGlobalDomainOffset(), _simpleMDConfig.getMPIConfiguration().getNumberOfProcesses(),
223 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize(), _mamicoConfig.getParallelTopologyConfiguration().getParallelTopologyType(),
224 _mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap(), (unsigned int)_rank
225#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
226 ,
227 _timeIntegrationService->getPintComm()
228#endif
229 );
230 } else {
231 coupling::indexing::IndexingService<3>::getInstance().initWithMDSize(
232 _simpleMDConfig.getDomainConfiguration().getGlobalDomainSize(), _simpleMDConfig.getDomainConfiguration().getGlobalDomainOffset(),
233 _simpleMDConfig.getMPIConfiguration().getNumberOfProcesses(), _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize(),
234 _mamicoConfig.getParallelTopologyConfiguration().getParallelTopologyType(), _mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap(),
235 (unsigned int)_rank
236#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
237 ,
238 _timeIntegrationService->getPintComm()
239#endif
240 );
241 }
242
243 // for timing measurements
244 _tv.micro = 0;
245 _tv.macro = 0;
246 _tv.filter = 0;
247
248 // even if _cfg.miSolverType == SYNTHETIC then
249 // multiMDService, _mdSimulations, _mdSolverInterface etc need to be initialized
250
251 unsigned int totNumMD;
252 if (_cfg.totalNumberMDSimulations > 0)
253 totNumMD = _cfg.totalNumberMDSimulations;
254 else // dynamic case, start with _cfg.lowerBoundNumberMDSimulations MD
255 totNumMD = _cfg.lowerBoundNumberMDSimulations;
256
257 _multiMDService = new tarch::utils::MultiMDService<3>(_simpleMDConfig.getMPIConfiguration().getNumberOfProcesses(), totNumMD
258#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
259 ,
260 _timeIntegrationService->getPintComm()
261#endif
262 );
263
265 if (_instanceHandling == nullptr) {
266 if (_isRootRank) {
267 std::cout << "ERROR CouetteScenario::initSolvers() : _instanceHandling == NULL!" << std::endl;
268 }
269 std::exit(EXIT_FAILURE);
270 }
271
272 _mdStepCounter = 0;
273 if (_isRootRank) {
274 gettimeofday(&_tv.start, NULL);
275 gettimeofday(&_tv.output, NULL);
276 }
277
279 // equilibrate MD
280 _instanceHandling->switchOffCoupling();
281 _instanceHandling->equilibrate(_cfg.equSteps, _mdStepCounter);
282 _instanceHandling->switchOnCoupling();
283 _mdStepCounter += _cfg.equSteps;
284 }
285
286 // allocate coupling interfaces
287 _instanceHandling->setMDSolverInterface();
288 _mdSolverInterface = _instanceHandling->getMDSolverInterface();
289
290 if (_cfg.twsLoop) {
291 // initialise coupling cell service for multi-MD case and set single
292 // cell services in each MD simulation
294 _mamicoConfig, "couette.xml", *_multiMDService, _tws);
295 } else {
296 // initialise coupling cell service for multi-MD case and set single
297 // cell services in each MD simulation
299 _mamicoConfig, "couette.xml", *_multiMDService);
300 }
301
302 // init filtering for all md instances
303 _multiMDCellService->constructFilterPipelines();
304
305 _multiMDMediator = new coupling::MultiMDMediator<MY_LINKEDCELL, 3>(*_multiMDCellService, *_instanceHandling, *_multiMDService, couetteSolverInterface);
306
307 // allocate solvers
308 _couetteSolver = NULL;
310 getCouetteSolver(_mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[0],
311 _simpleMDConfig.getSimulationConfiguration().getDt() * _simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps());
312 _MDBoundarySetupDone = false;
313
315 // set couette solver interface in MamicoInterfaceProvider
317
318 _instanceHandling->setCouplingCellServices(*_multiMDCellService);
319 // compute and store temperature in coupling cells (temp=1.1
320 // everywhere)
321 _multiMDCellService->computeAndStoreTemperature(_cfg.temp);
322 }
323
324 /*
325 * A synthethic solver is modeled using a dynamically linked filter,
326 * i.e. a lambda function producing artifical data in every filter step.
327 *
328 * This is how to properly instanciate and use a synthethic solver:
329 * - Create a sequence named like the macro SYNTHETICMD_SEQUENCE in xml. Use
330 * whatever input, but make sure the filter system's output is set to this
331 * sequence.
332 * - Set filtered-values = "macro-mass macro-momentum" for that sequence.
333 * - Use that sequence as input for all sequences that want (unfiltered) MD
334 * input.
335 *
336 * Note that this synthethic solver is designed to be used on sequential
337 * mode only and with only one MD instance.
338 *
339 * TODO
340 * - major bug when there is ONLY a FFF in a sequence
341 * - reduce capture: most variables in lambda can be defined beforehand as
342 * they are const (e.g. everything coming from cfg)
343 * - totalNumberMDSimulations > 1 is theoretically possible with this
344 * redesign. test it and remove the restriction for it to be 1
345 */
346
348 // Synthetic MD runs sequentially only, as described above.
349 if (_cfg.macro2Md || _cfg.totalNumberMDSimulations != 1 || _cfg.lbNumberProcesses[0] != 1 || _cfg.lbNumberProcesses[1] != 1 ||
350 _cfg.lbNumberProcesses[2] != 1) {
351 throw std::runtime_error("ERROR: Syntethic MD is only available in sequential mode!");
352 }
353
354 // set couette solver interface in MamicoInterfaceProvider
356
359
360 // Create new FilterFromFunction instance and insert it into Filtering
361 // System.
362 try {
363 _multiMDCellService->getCouplingCellService(0)
364 .getFilterPipeline()
365 ->getSequence(SYNTHETICMD_SEQUENCE)
366 ->addFilter(
367 new std::function<std::vector<double>(std::vector<double>, std::vector<std::array<unsigned int, 3>>)>{
368 // applyScalar
369 [this](std::vector<double> inputScalars, // not actually used as input:
370 // matching MCS's
371 // addFilterToSequence(...)
372 // signature
373 std::vector<std::array<unsigned int, 3>> cellIndices // not used either
374 ) {
375 if (_isRootRank) {
376 gettimeofday(&_tv.start, NULL);
377 }
378
379 // std::cout << "Entering synthetic MD scalar..." << std::endl;
380
381 const unsigned int size = inputScalars.size();
382 auto dx = coupling::indexing::IndexingService<3>::getInstance().getCouplingCellSize();
383 const tarch::la::Vector<3, double> couplingCellSize(dx);
384 const double mass = (_cfg.density) * couplingCellSize[0] * couplingCellSize[1] * couplingCellSize[2];
385
386 std::vector<double> syntheticMasses;
387 for (unsigned int i = 0; i < size; i++) {
388 syntheticMasses.push_back(mass);
389 }
390 // std::cout << "Generated masses!" << std::endl;
391
392 if (_isRootRank) {
393 gettimeofday(&_tv.end, NULL);
394 _tv.micro += (_tv.end.tv_sec - _tv.start.tv_sec) * 1000000 + (_tv.end.tv_usec - _tv.start.tv_usec);
395 }
396
397 return syntheticMasses;
398 }},
399 new std::function<std::vector<std::array<double, 3>>(std::vector<std::array<double, 3>>, std::vector<std::array<unsigned int, 3>>)>{
400 // applyVector
401 [this](std::vector<std::array<double, 3>> inputVectors, // same as above
402 std::vector<std::array<unsigned int, 3>> cellIndices // same as above
403 ) {
404 if (_isRootRank) {
405 gettimeofday(&_tv.start, NULL);
406 }
407
408 // std::cout << "Entering synthetic MD vector." << std::endl;
409
410 const unsigned int size = inputVectors.size();
411 auto dx = coupling::indexing::IndexingService<3>::getInstance().getCouplingCellSize();
412 const tarch::la::Vector<3, double> couplingCellSize(dx);
413 const double mass = (_cfg.density) * couplingCellSize[0] * couplingCellSize[1] * couplingCellSize[2];
414
415 using coupling::indexing::IndexTrait;
417
418 const tarch::la::Vector<3, double> md2MacroDomainOffset = {
419 CellIndex<3, IndexTrait::local, IndexTrait::md2macro, IndexTrait::noGhost>::lowerBoundary.get()[0] * couplingCellSize[0],
420 CellIndex<3, IndexTrait::local, IndexTrait::md2macro, IndexTrait::noGhost>::lowerBoundary.get()[1] * couplingCellSize[1],
421 CellIndex<3, IndexTrait::local, IndexTrait::md2macro, IndexTrait::noGhost>::lowerBoundary.get()[2] * couplingCellSize[2],
422 };
423
424 std::normal_distribution<double> distribution(0.0, _cfg.noiseSigma);
425 std::vector<std::array<double, 3>> syntheticMomenta;
426 for (unsigned int i = 0; i < size; i++) {
427 // determine cell midpoint
428 CellIndex<3, IndexTrait::vector> globalIndex(
429 CellIndex<3, IndexTrait::local, IndexTrait::md2macro, IndexTrait::noGhost>{i}); // construct global CellIndex from buffer
430 // and convert it to vector
431 tarch::la::Vector<3, double> cellMidPoint(md2MacroDomainOffset - 0.5 * couplingCellSize);
432 for (unsigned int d = 0; d < 3; d++) {
433 cellMidPoint[d] = cellMidPoint[d] + ((double)globalIndex.get()[d]) * couplingCellSize[d];
434 }
435
436 // compute momentum
437 const tarch::la::Vector<3, double> noise(distribution(_generator), distribution(_generator), distribution(_generator));
438 const tarch::la::Vector<3, double> momentum(mass * ((*_couetteSolver).getVelocity(cellMidPoint) + noise));
439
440 // conversion from tarch::la::Vector to std::array
441 syntheticMomenta.push_back({momentum[0], momentum[1], momentum[2]});
442 }
443 // std::cout << "Generated momenta!" << std::endl;
444
445 if (_isRootRank) {
446 gettimeofday(&_tv.end, NULL);
447 _tv.micro += (_tv.end.tv_sec - _tv.start.tv_sec) * 1000000 + (_tv.end.tv_usec - _tv.start.tv_usec);
448 }
449
450 return syntheticMomenta;
451 }},
452 0 // filterIndex
453 );
454 } catch (std::runtime_error& e) {
455 auto expectedError = std::string("ERROR: Could not find Filter Sequence named ").append(SYNTHETICMD_SEQUENCE);
456 if (expectedError.compare(e.what()) == 0) {
457 std::cout << "ERROR: Synthetic MD solver selected without providing "
458 "filter sequence '"
459 << SYNTHETICMD_SEQUENCE << "' in config." << std::endl;
460 exit(EXIT_FAILURE);
461 } else
462 throw e;
463 }
464 }
465
466 // allocate buffers for send/recv operations
467 allocateMacro2mdBuffer(*couetteSolverInterface);
468 allocateMd2macroBuffer(*couetteSolverInterface);
469
470 if (_cfg.initAdvanceCycles > 0 && _couetteSolver != NULL)
471 _couetteSolver->advance(_cfg.initAdvanceCycles * _simpleMDConfig.getSimulationConfiguration().getDt() *
472 _simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps());
473
474 // finish time measurement for initialisation
475 if (_isRootRank) {
476 gettimeofday(&_tv.end, NULL);
477 double runtime = (_tv.end.tv_sec - _tv.start.tv_sec) * 1000000 + (_tv.end.tv_usec - _tv.start.tv_usec);
478 std::cout << "Initialization: " << (int)(runtime / 1000) << "ms" << std::endl;
479 }
480
481 if (_cfg.computeSNR) {
482 std::cout << "Output for every coupling cycle, for the cell 87 in recvBuffer:" << std::endl;
483 std::cout << "cycle number (after filter-init-cycles), vel_x macroscopic "
484 "solver, vel_x filter output"
485 << std::endl;
486 _sum_signal = 0;
487 _sum_noise = 0;
488 }
489
490 if (_isRootRank) {
491 gettimeofday(&_tv.start_total, NULL);
492 }
493 if (_isRootRank) {
494 std::cout << "Finish CouetteScenario::initSolvers() " << std::endl;
495 }
496 }
497
501 void advanceMacro(int cycle) {
502 if (_couetteSolver != NULL) {
503 if (_isRootRank) {
504 gettimeofday(&_tv.start, NULL);
505 }
506
507 // run one time step for macroscopic couette solver
508 if (_cfg.wallInitCycles > 0 && cycle == _cfg.wallInitCycles) {
509 _couetteSolver->setWallVelocity(_cfg.wallVelocity);
510 // When using Synthetic MD,
511 }
512 if (_cfg.wallOscillations != 0) {
513 tarch::la::Vector<3, double> vel = cycle < _cfg.wallInitCycles ? _cfg.wallInitVelocity : _cfg.wallVelocity;
514 vel = vel * cos(2 * M_PI * _cfg.wallOscillations * cycle / _cfg.couplingCycles);
515 _couetteSolver->setWallVelocity(vel);
516 }
517 _couetteSolver->advance(_simpleMDConfig.getSimulationConfiguration().getDt() * _simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps());
518 if (_isRootRank) {
519 gettimeofday(&_tv.end, NULL);
520 _tv.macro += (_tv.end.tv_sec - _tv.start.tv_sec) * 1000000 + (_tv.end.tv_usec - _tv.start.tv_usec);
521 // std::cout << "Finish _couetteSolver->advance " << std::endl;
522 }
523
524 // extract data from couette solver and send them to MD (can take any
525 // index-conversion object)
526 fillSendBuffer(_cfg.density, *_couetteSolver, _couplingBuffer.macro2MDBuffer);
527 }
528 if (_cfg.macro2Md) {
529#ifdef MAMICO_USE_COLLECTIVE_MPI
530 _multiMDCellService->bcastFromMacro2MD(_couplingBuffer.macro2MDBuffer);
531#else
532 _multiMDCellService->sendFromMacro2MD(_couplingBuffer.macro2MDBuffer);
533#endif
534 // std::cout << "Finish _multiMDCellService->sendFromMacro2MD " <<
535 // std::endl;
536 }
537 }
538
539 void varyMD(int cycle) {
540 // Non-dynamic case: Nothing to do.
542 return;
543
544 // modify number of instances only once every 10 cycles
545 if (cycle % 10 != 0)
546 return;
547
548 int addMDInstances = 0;
549
550 if (_rank == 0) {
551 tarch::la::Vector<3, double> vel(0, 0, 0);
552 double mass = 0;
553 for (auto pair : _couplingBuffer.md2macroBuffer) {
554 I01 idx;
556 std::tie(couplingCell, idx) = pair;
557 vel += couplingCell->getMacroscopicMomentum();
558 mass += couplingCell->getMacroscopicMass();
559 }
560 vel = vel / (double)_couplingBuffer.md2macroBuffer.size();
561 mass /= _couplingBuffer.md2macroBuffer.size();
562
563 double soundSpeed =
564 (1 / std::sqrt(3)) * (_mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[0] /
565 (_simpleMDConfig.getSimulationConfiguration().getDt() * _simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps()));
566 double cellVolume = _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[0] *
567 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[1] *
568 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[2];
569
570 coupling::error::ErrorEstimation errorControl(vel[0], _cfg.temp, mass, _simpleMDConfig.getMoleculeConfiguration().getMass(), soundSpeed,
571 _multiMDService->getTotalNumberOfMDSimulations(), cellVolume);
572
573 double simtime = (double)cycle / _cfg.couplingCycles;
574 double targetError = _cfg.absVelErrEnd * simtime + _cfg.absVelErrStart * (1 - simtime);
575 errorControl.setAbsVelocityError(targetError);
576 double NoMD = errorControl.getCorrectorNumberOfSamples(coupling::error::ErrorEstimation::Velocity, coupling::error::ErrorEstimation::Absolute);
577
578 addMDInstances = NoMD - _multiMDService->getTotalNumberOfMDSimulations();
579 }
580
581#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
582 MPI_Bcast(&addMDInstances, 1, MPI_INT, 0, _timeIntegrationService->getPintComm());
583#endif
584 if (addMDInstances < 0) {
585 if ((int)_multiMDMediator->getNumberOfActiveMDSimulations() + addMDInstances < _cfg.lowerBoundNumberMDSimulations) {
586 addMDInstances = (_cfg.lowerBoundNumberMDSimulations - (int)_multiMDMediator->getNumberOfActiveMDSimulations());
587 }
588 }
589 if (addMDInstances < 0) {
590 if (_rank == 0)
591 std::cout << "Removing " << -addMDInstances << " of " << _multiMDService->getTotalNumberOfMDSimulations() << " MD simulations" << std::endl;
592 _multiMDMediator->rmNMDSimulations(-addMDInstances);
593 } else if (addMDInstances > 0) {
594 if (_rank == 0)
595 std::cout << "Adding " << addMDInstances << " to " << _multiMDService->getTotalNumberOfMDSimulations() << " MD simulations" << std::endl;
596 _multiMDMediator->addNMDSimulations(addMDInstances);
597 }
598 }
599
603 void advanceMicro(int cycle) {
604 if (_isRootRank) {
605 gettimeofday(&_tv.start, NULL);
606 }
608 // run MD instances
609 _instanceHandling->simulateTimesteps(_simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps(), _mdStepCounter, *_multiMDCellService);
610 // plot macroscopic time step info in multi md service
611 _multiMDCellService->plotEveryMacroscopicTimestep(cycle);
612
613 _mdStepCounter += _simpleMDConfig.getSimulationConfiguration().getNumberOfTimesteps();
614
615 if (_isRootRank) {
616 gettimeofday(&_tv.end, NULL);
617 _tv.micro += (_tv.end.tv_sec - _tv.start.tv_sec) * 1000000 + (_tv.end.tv_usec - _tv.start.tv_usec);
618 }
619
620 // send back data from MD instances and merge it
621 if (_cfg.md2Macro) {
622#ifdef MAMICO_USE_COLLECTIVE_MPI
623 _tv.filter += _multiMDCellService->reduceFromMD2Macro(_couplingBuffer.md2macroBuffer);
624#else
625 _tv.filter += _multiMDCellService->sendFromMD2Macro(_couplingBuffer.md2macroBuffer);
626#endif
627 // std::cout << "Finish _multiMDCellService->sendFromMD2Macro " <<
628 // std::endl;
629 }
630 }
631
633 if (_cfg.md2Macro) {
634 //_couplingBuffer does not get used here: Instead, the synthetic MD in the
635 // SYNTHETICMD_SEQUENCE generates values. To prevent segfaults, it has
636 // to be nonempty, though.
637#ifdef MAMICO_USE_COLLECTIVE_MPI
638 _tv.filter += _multiMDCellService->reduceFromMD2Macro(_couplingBuffer.md2macroBuffer);
639#else
640 _tv.filter += _multiMDCellService->sendFromMD2Macro(_couplingBuffer.md2macroBuffer);
641#endif
642 }
643 }
644 }
645
647 void computeSNR(int cycle) {
648 if (_cfg.computeSNR && cycle >= _cfg.filterInitCycles) {
649 std::cout << cycle - _cfg.filterInitCycles << ", ";
650 using namespace coupling::indexing;
651 const tarch::la::Vector<3, double> dx(IndexingService<3>::getInstance().getCouplingCellSize());
652 const double mass = _cfg.density * dx[0] * dx[1] * dx[2];
653 unsigned int i = 0;
654 for (auto pair : _couplingBuffer.md2macroBuffer) {
656 if (i == 87) {
657 I01 idx;
659 std::tie(couplingCell, idx) = pair;
660 auto midPoint = idx.getCellMidPoint();
661 double vx_macro = _couetteSolver->getVelocity(midPoint)[0];
662 double vx_filter = (1 / mass * couplingCell->getMacroscopicMomentum())[0];
663 _sum_noise += (vx_macro - vx_filter) * (vx_macro - vx_filter);
664 _sum_signal += vx_macro * vx_macro;
665 std::cout << vx_macro << ", " << vx_filter << std::endl;
666 }
667 i++;
668 }
669 }
670 }
671
675 void twoWayCoupling(int cycle) {
677 if (_cfg.twoWayCoupling) {
678#if (BUILD_WITH_OPENFOAM)
679 if ((_cfg.maSolverType == CouetteConfig::COUETTE_FOAM) && cycle == _cfg.filterInitCycles && _couetteSolver != NULL) {
680 static_cast<coupling::solvers::IcoFoamInterface*>(_couetteSolver)->setupMDBoundary();
681 }
682#endif
683 if ((_cfg.maSolverType == CouetteConfig::COUETTE_LB || _cfg.maSolverType == CouetteConfig::COUETTE_FD) && cycle >= _cfg.filterInitCycles) {
684 if (!_MDBoundarySetupDone) {
686 ->setMDBoundary(_simpleMDConfig.getDomainConfiguration().getGlobalDomainOffset(), _simpleMDConfig.getDomainConfiguration().getGlobalDomainSize(),
687 _mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap());
688 _MDBoundarySetupDone = true;
689 }
690 static_cast<coupling::solvers::LBCouetteSolver*>(_couetteSolver)->setMDBoundaryValues(_couplingBuffer.md2macroBuffer);
691 }
692#if (BUILD_WITH_OPENFOAM)
693 else if (_cfg.maSolverType == CouetteConfig::COUETTE_FOAM && cycle >= _cfg.filterInitCycles && _couetteSolver != NULL) {
694 static_cast<coupling::solvers::IcoFoamInterface*>(_couetteSolver)->setMDBoundaryValues(_couplingBuffer.md2macroBuffer);
695 }
696#endif
697 }
698 // write data to csv-compatible file for evaluation
699 write2CSV(_couplingBuffer.md2macroBuffer, cycle + 1);
700 }
701
704 void shutdown() {
705 if (_cfg.computeSNR) {
706 std::cout << "SNR = " << 10 * log10(_sum_signal / _sum_noise) << std::endl;
707 }
708
709 // finish time measurement for coupled simulation
710 if (_isRootRank) {
711 gettimeofday(&_tv.end, NULL);
712 double time_total = (_tv.end.tv_sec - _tv.start_total.tv_sec) * 1000000 + (_tv.end.tv_usec - _tv.start_total.tv_usec);
713 std::cout << "Finished all coupling cycles after " << time_total / 1000000 << " s" << std::endl;
714 if (_cfg.twsLoop)
715 std::cout << "TWS = " << _tws << std::endl;
716 std::cout << "Time percentages Micro, Macro, Filter, Other: " << std::endl;
717 std::cout << _tv.micro / time_total * 100 << ", " << _tv.macro / time_total * 100 << ", " << _tv.filter / time_total * 100 << ", "
718 << (1 - (_tv.micro + _tv.macro + _tv.filter) / time_total) * 100 << std::endl;
719 }
720
721 if (_instanceHandling != nullptr) {
722 delete _instanceHandling;
723 }
724
727
728 if (_multiMDService != NULL) {
729 delete _multiMDService;
730 _multiMDService = NULL;
731 }
732 if (couetteSolverInterface != NULL) {
733 delete couetteSolverInterface;
734 couetteSolverInterface = NULL;
735 }
736 if (_couetteSolver != NULL) {
737 delete _couetteSolver;
738 _couetteSolver = NULL;
739 }
740 if (_multiMDCellService != NULL) {
741 delete _multiMDCellService;
742 _multiMDCellService = NULL;
743 }
744 if (_multiMDMediator != nullptr) {
745 delete _multiMDMediator;
746 _multiMDMediator = nullptr;
747 }
748
749 if (_isRootRank) {
750 std::cout << "Finish CouetteScenario::shutdown() " << std::endl;
751 }
752 }
753
755 tarch::la::Vector<3, unsigned int> getGlobalNumberCouplingCells(const simplemd::configurations::MolecularDynamicsConfiguration& simpleMDConfig,
756 const coupling::configurations::MaMiCoConfiguration<3>& mamicoConfig) const {
757 tarch::la::Vector<3, double> domainSize(simpleMDConfig.getDomainConfiguration().getGlobalDomainSize());
758 tarch::la::Vector<3, double> dx(mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize());
759 tarch::la::Vector<3, unsigned int> globalNumberCouplingCells(0);
760 for (unsigned int d = 0; d < 3; d++) {
761 int buf = floor(domainSize[d] / dx[d] + 0.5);
762 globalNumberCouplingCells[d] = (unsigned int)buf;
763 }
764 return globalNumberCouplingCells;
765 }
766
771 for (auto idx : I08()) {
772 if (!I12::contains(idx)) {
773 if (tarch::utils::contains(msi.getSourceRanks(idx), (unsigned int)_rank)) {
775 _couplingBuffer.macro2MDBuffer << std::make_pair(couplingCell, idx);
776 if (couplingCell == nullptr)
777 throw std::runtime_error(std::string("ERROR CouetteScenario::allocateMacro2mdBuffer: couplingCell==NULL!"));
778 }
779 }
780 }
781 }
782
786 for (auto idx : I12()) {
787 if (tarch::utils::contains(msi.getTargetRanks(idx), (unsigned int)_rank)) {
789 _couplingBuffer.md2macroBuffer << std::make_pair(couplingCell, idx);
790 if (couplingCell == nullptr)
791 throw std::runtime_error(std::string("ERROR CouetteScenario::allocateMacro2mdBuffer: couplingCell==NULL!"));
792 }
793 }
794 }
795
800 template <class Container_T> void write2CSV(Container_T& md2macroBuffer, int couplingCycle) const {
801 if (md2macroBuffer.size() == 0)
802 return;
803 if (_cfg.csvEveryTimestep < 1 || couplingCycle % _cfg.csvEveryTimestep > 0)
804 return;
805 // form file name and open file
806 std::stringstream ss;
807 ss << "CouetteAvgMultiMDCells_d" << _timeIntegrationService->getPintDomain() << "_r" << _rank << "_c" << couplingCycle;
808 if (_timeIntegrationService->isPintEnabled())
809 ss << "_i" << _timeIntegrationService->getIteration();
810 ss << ".csv";
811 std::ofstream file(ss.str().c_str());
812 if (!file.is_open()) {
813 std::cout << "ERROR CouetteScenario::write2CSV(): Could not open file " << ss.str() << "!" << std::endl;
814 exit(EXIT_FAILURE);
815 }
816
817 // loop over md2macro cells; read macroscopic mass+momentum buffers and
818 // write cell index, mass and velocity to one line in the csv-file
819 file << "I01_x,I01_y,I01_z,vel_x,vel_y,vel_z,mass" << std::endl;
820 for (auto pair : md2macroBuffer) {
821 I01 idx;
823 std::tie(couplingCell, idx) = pair;
825 if (couplingCell->getMacroscopicMass() != 0.0) {
826 vel = (1.0 / couplingCell->getMacroscopicMass()) * vel;
827 }
828 file << idx << "," << vel[0] << "," << vel[1] << "," << vel[2] << "," << couplingCell->getMacroscopicMass() << std::endl;
829 }
830
831 // close file
832 file.close();
833 }
834
841 void fillSendBuffer(const double density, const coupling::solvers::AbstractCouetteSolver<3>& couetteSolver,
844 using namespace coupling::indexing;
845 const tarch::la::Vector<3, double> dx(IndexingService<3>::getInstance().getCouplingCellSize());
846 double mass = density * dx[0] * dx[1] * dx[2];
847 for (auto pair : macro2MDBuffer) {
848 I01 idx;
850 std::tie(couplingCell, idx) = pair;
851 auto midPoint = idx.getCellMidPoint();
852 if (_cfg.maSolverType == CouetteConfig::COUETTE_LB || _cfg.maSolverType == CouetteConfig::COUETTE_FD)
853 mass *= static_cast<const coupling::solvers::LBCouetteSolver*>(&couetteSolver)->getDensity(midPoint);
854 // compute momentum
855 tarch::la::Vector<3, double> momentum(mass * couetteSolver.getVelocity(midPoint));
856 couplingCell->setMicroscopicMass(mass);
857 couplingCell->setMicroscopicMomentum(momentum);
858 }
859 }
860
867 tarch::la::Vector<3, double> vel = _cfg.wallInitCycles > 0 ? _cfg.wallInitVelocity : _cfg.wallVelocity;
868 // analytical solver: is only active on rank 0
869 if (_cfg.maSolverType == CouetteConfig::COUETTE_ANALYTICAL) {
870 if (_rank == 0 || _cfg.miSolverType == CouetteConfig::SYNTHETIC) {
871 solver = new coupling::solvers::CouetteSolver<3>(_cfg.channelheight, vel[0], _cfg.kinVisc);
872 if (solver == NULL) {
873 std::cout << "ERROR CouetteScenario::getCouetteSolver(): Analytic solver==NULL!" << std::endl;
874 exit(EXIT_FAILURE);
875 }
876 }
877 }
878#if (BUILD_WITH_OPENFOAM)
879 else if (_cfg.maSolverType == CouetteConfig::COUETTE_FOAM) {
880 solver = new coupling::solvers::IcoFoamInterface(_rank, _cfg.plotEveryTimestep, _cfg.channelheight, _cfg.foam.directory, _cfg.foam.folder,
881 _cfg.foam.boundariesWithMD, _cfg.wallVelocity);
882 if (solver == NULL) {
883 std::cout << "ERROR CouetteScenario::getCouetteSolver(): IcoFoamInterface solver==NULL!" << std::endl;
884 exit(EXIT_FAILURE);
885 }
886 }
887#endif
888 // LB solver: active on lbNumberProcesses
889 else if (_cfg.maSolverType == CouetteConfig::COUETTE_LB) {
890 solver = new coupling::solvers::LBCouetteSolver(_cfg.channelheight, vel, _cfg.kinVisc, dx, dt, _cfg.plotEveryTimestep, _cfg.plotAverageVelocity,
891 "LBCouette", _cfg.lbNumberProcesses, 1, this);
892 if (solver == NULL) {
893 std::cout << "ERROR CouetteScenario::getCouetteSolver(): LB solver==NULL!" << std::endl;
894 exit(EXIT_FAILURE);
895 }
896 } else if (_cfg.maSolverType == CouetteConfig::COUETTE_FD) {
897 solver = new coupling::solvers::FiniteDifferenceSolver(_cfg.channelheight, vel, _cfg.kinVisc, dx, dt, _cfg.plotEveryTimestep, "FDCouette",
898 _cfg.lbNumberProcesses);
899 if (solver == NULL) {
900 std::cout << "ERROR CouetteScenario::getCouetteSolver(): FD solver==NULL!" << std::endl;
901 exit(EXIT_FAILURE);
902 }
903 } else {
904 std::cout << "ERROR CouetteScenario::getCouetteSolver(): Unknown solver type!" << std::endl;
905 exit(EXIT_FAILURE);
906 }
907 return solver;
908 }
909
919 tarch::la::Vector<3, double> mamicoMeshsize,
920 tarch::la::Vector<3, unsigned int> globalNumberCouplingCells,
921 unsigned int outerRegion) {
924 if (_cfg.maSolverType == CouetteConfig::COUETTE_ANALYTICAL) {
925 interface = new coupling::solvers::CouetteSolverInterface<3>(globalNumberCouplingCells, outerRegion);
926 } else if (_cfg.maSolverType == CouetteConfig::COUETTE_LB) {
927 // compute number of cells of MD offset; detect any mismatches!
928 tarch::la::Vector<3, unsigned int> offsetMDDomain(0);
929 for (unsigned int d = 0; d < 3; d++) {
930 if (mdOffset[d] < 0.0) {
931 std::cout << "ERROR CouetteScenario::getCouetteSolverInterface(...): mdOffset[" << d << "]<0.0!" << std::endl;
932 exit(EXIT_FAILURE);
933 }
934 offsetMDDomain[d] = floor(mdOffset[d] / mamicoMeshsize[d] + 0.5);
935 if (fabs((offsetMDDomain[d] * mamicoMeshsize[d] - mdOffset[d]) / mamicoMeshsize[d]) > 1.0e-8) {
936 std::cout << "ERROR CouetteScenario::getCouetteSolverInterface: MD offset and mesh size mismatch!" << std::endl;
937 exit(EXIT_FAILURE);
938 }
939 }
940 tarch::la::Vector<3, unsigned int> cells_per_process{
942 coupling::solvers::NumericalSolver::getAvgDomainSize(_cfg.channelheight, dx, _cfg.lbNumberProcesses, 1),
943 coupling::solvers::NumericalSolver::getAvgDomainSize(_cfg.channelheight, dx, _cfg.lbNumberProcesses, 2)}};
944 interface =
945 new coupling::solvers::LBCouetteSolverInterface(cells_per_process, _cfg.lbNumberProcesses, offsetMDDomain, globalNumberCouplingCells, outerRegion);
946 }
947#if (BUILD_WITH_OPENFOAM)
948 else if (_cfg.maSolverType == CouetteConfig::COUETTE_FOAM) {
949 interface = new coupling::solvers::CouetteSolverInterface<3>(globalNumberCouplingCells, outerRegion);
950 }
951#endif
952 else if (_cfg.maSolverType == CouetteConfig::COUETTE_FD) {
953 // compute number of cells of MD offset; detect any mismatches!
954 tarch::la::Vector<3, unsigned int> offsetMDDomain(0);
955 for (unsigned int d = 0; d < 3; d++) {
956 if (mdOffset[d] < 0.0) {
957 std::cout << "ERROR CouetteScenario::getCouetteSolverInterface(...): mdOffset[" << d << "]<0.0!" << std::endl;
958 exit(EXIT_FAILURE);
959 }
960 offsetMDDomain[d] = floor(mdOffset[d] / mamicoMeshsize[d] + 0.5);
961 if (fabs((offsetMDDomain[d] * mamicoMeshsize[d] - mdOffset[d]) / mamicoMeshsize[d]) > 1.0e-8) {
962 std::cout << "ERROR CouetteScenario::getCouetteSolverInterface: MD offset and mesh size mismatch!" << std::endl;
963 exit(EXIT_FAILURE);
964 }
965 }
966 tarch::la::Vector<3, unsigned int> cells_per_process{
968 coupling::solvers::NumericalSolver::getAvgDomainSize(_cfg.channelheight, dx, _cfg.lbNumberProcesses, 1),
969 coupling::solvers::NumericalSolver::getAvgDomainSize(_cfg.channelheight, dx, _cfg.lbNumberProcesses, 2)}};
970 interface =
971 new coupling::solvers::LBCouetteSolverInterface(cells_per_process, _cfg.lbNumberProcesses, offsetMDDomain, globalNumberCouplingCells, outerRegion);
972 }
973
974 if (interface == NULL) {
975 std::cout << "ERROR CouetteScenario::getCouetteSolverInterface(...), rank=" << _rank << ": interface==NULL!" << std::endl;
976 exit(EXIT_FAILURE);
977 }
978 return interface;
979 }
980
990
995 timeval start_total;
997 timeval start;
999 timeval end;
1000 timeval output;
1002 double micro;
1004 double macro;
1006 double filter;
1007 };
1008
1012 int _tws;
1014 simplemd::configurations::MolecularDynamicsConfiguration _simpleMDConfig;
1020 unsigned int _mdStepCounter;
1030 std::vector<coupling::interface::MDSolverInterface<MY_LINKEDCELL, 3>*> _mdSolverInterface;
1033 std::default_random_engine _generator;
1041 bool _MDBoundarySetupDone;
1042};
1043
1044#endif // _COUPLING_SCENARIO_COUETTESCENARIO_H_
coupling::solvers::AbstractCouetteSolver< 3 > * getCouetteSolver(const double dx, const double dt)
Definition CouetteScenario.h:864
CouetteScenario()
simple constructor
Definition CouetteScenario.h:72
void allocateMd2macroBuffer(coupling::interface::MacroscopicSolverInterface< 3 > &msi)
Definition CouetteScenario.h:785
void computeSNR(int cycle)
used to compute signal-to-noise ratio, stores values in _sum_noise and _sum_signal
Definition CouetteScenario.h:647
void parseConfigurations()
reads the configuration from the xml file and calls parseCouetteScenarioConfiguration()
Definition CouetteScenario.h:164
TimingValues _tv
a instance of the timingValues
Definition CouetteScenario.h:1039
coupling::configurations::MaMiCoConfiguration< 3 > _mamicoConfig
the config data and information for MaMiCo
Definition CouetteScenario.h:1016
void fillSendBuffer(const double density, const coupling::solvers::AbstractCouetteSolver< 3 > &couetteSolver, coupling::datastructures::FlexibleCellContainer< 3 > &macro2MDBuffer) const
fills send buffer with data from macro/continuum solver
Definition CouetteScenario.h:841
void advanceMacro(int cycle)
advances the continuum solver and collects data to send to md (fillSendBuffer())
Definition CouetteScenario.h:501
void twsLoop()
Executes the entire test several times for a range of time-window-size parameters.
Definition CouetteScenario.h:153
void initSolvers()
initialises the macro and micro solver according to the setup from the xml file and pre-proccses them
Definition CouetteScenario.h:208
void allocateMacro2mdBuffer(coupling::interface::MacroscopicSolverInterface< 3 > &msi)
allocates the send buffer (with values for all coupling cells).
Definition CouetteScenario.h:770
CouplingBuffer _couplingBuffer
the current buffer for the coupling
Definition CouetteScenario.h:1028
int _tws
Definition CouetteScenario.h:1012
void equilibrateMicro() override
Definition CouetteScenario.h:127
void advanceMicro(int cycle)
advances the md solver for one coupling time step and collect the data for the coupling
Definition CouetteScenario.h:603
double _sum_signal
Definition CouetteScenario.h:1035
void init() override
initialises everthing necessary for the test
Definition CouetteScenario.h:90
unsigned int _mdStepCounter
the counter for the time steps, which are done within the md
Definition CouetteScenario.h:1020
double _sum_noise
Definition CouetteScenario.h:1037
tarch::utils::MultiMDService< 3 > * _multiMDService
Definition CouetteScenario.h:1024
void shutdown()
finalize the time measurement, and cleans up at the end of the simulation
Definition CouetteScenario.h:704
simplemd::configurations::MolecularDynamicsConfiguration _simpleMDConfig
the config data and information for SimpleMD
Definition CouetteScenario.h:1014
coupling::interface::MacroscopicSolverInterface< 3 > * getCouetteSolverInterface(const double dx, tarch::la::Vector< 3, double > mdOffset, tarch::la::Vector< 3, double > mamicoMeshsize, tarch::la::Vector< 3, unsigned int > globalNumberCouplingCells, unsigned int outerRegion)
Definition CouetteScenario.h:918
virtual ~CouetteScenario()
a dummy destructor
Definition CouetteScenario.h:74
void twoWayCoupling(int cycle)
sets up the boundaries in the macro solver for the coupling and applies the values from the md in the...
Definition CouetteScenario.h:675
void run() override
runs the simulation
Definition CouetteScenario.h:78
void runOneCouplingCycle(int cycle) override
combines the functioniality necessary for a cycle of the coupled simulation
Definition CouetteScenario.h:106
int _rank
the rank of the current MPI process in the local time domain
Definition CouetteScenario.h:1010
std::default_random_engine _generator
Definition CouetteScenario.h:1033
coupling::configurations::CouetteConfig _cfg
the CouetteConfig for the current setup
Definition CouetteScenario.h:1018
tarch::la::Vector< 3, unsigned int > getGlobalNumberCouplingCells(const simplemd::configurations::MolecularDynamicsConfiguration &simpleMDConfig, const coupling::configurations::MaMiCoConfiguration< 3 > &mamicoConfig) const
Definition CouetteScenario.h:755
coupling::services::MultiMDCellService< MY_LINKEDCELL, 3 > * _multiMDCellService
Definition CouetteScenario.h:1026
std::vector< coupling::interface::MDSolverInterface< MY_LINKEDCELL, 3 > * > _mdSolverInterface
the interface to the md solver
Definition CouetteScenario.h:1030
coupling::solvers::AbstractCouetteSolver< 3 > * _couetteSolver
the current macro/continuum solver
Definition CouetteScenario.h:1022
void write2CSV(Container_T &md2macroBuffer, int couplingCycle) const
write coupling cells that have been received from MD to csv file
Definition CouetteScenario.h:800
bool _isRootRank
if this is the world global root process
Definition Scenario.h:50
Simulation slots are managed (i.e., added/removed) via this class. Works and interacts with the class...
Definition InstanceHandling.h:36
Class to handle interaction between MultiMDCellService and InstanceHandling This is currently mainly ...
Definition MultiMDMediator.h:25
parses all sub-tags for MaMiCo configuration.
Definition MaMiCoConfiguration.h:31
const coupling::configurations::CouplingCellConfiguration< dim > & getCouplingCellConfiguration() const
Definition MaMiCoConfiguration.h:68
defines the cell type with cell-averaged quantities only (no linked cells).
Definition CouplingCell.h:29
const tarch::la::Vector< dim, double > & getMacroscopicMomentum() const
Definition CouplingCell.h:64
void setMicroscopicMass(const double &mass)
Definition CouplingCell.h:42
const double & getMacroscopicMass() const
Definition CouplingCell.h:58
void setMicroscopicMomentum(const tarch::la::Vector< dim, double > &momentum)
Definition CouplingCell.h:48
provides access to coupling cells, which may belong to different indexing domains
Definition FlexibleCellContainer.h:30
@ Velocity
Definition ErrorEstimation.h:51
@ Absolute
Definition ErrorEstimation.h:61
Definition CellIndex.h:85
static bool contains(const coupling::indexing::BaseIndex< dim > &index)
Definition CellIndex.h:227
interface for the macroscopic, i.e. continuum solver
Definition MacroscopicSolverInterface.h:23
virtual std::vector< unsigned int > getSourceRanks(I01 idx)
Definition MacroscopicSolverInterface.h:67
virtual std::vector< unsigned int > getTargetRanks(I01 idx)
Definition MacroscopicSolverInterface.h:78
void setMDSolverInterface(coupling::interface::MDSolverInterface< LinkedCell, dim > *mdSolverInterface)
Definition MamicoInterfaceProvider.h:48
void setMacroscopicSolverInterface(coupling::interface::MacroscopicSolverInterface< dim > *macroscopicSolverInterface)
Definition MamicoInterfaceProvider.h:36
static MamicoInterfaceProvider & getInstance()
Definition MamicoInterfaceProvider.h:28
coupling::interface::MacroscopicSolverInterface< dim > * getMacroscopicSolverInterface()
Definition MamicoInterfaceProvider.h:43
void setCouplingCellService(coupling::services::CouplingCellService< dim > *couplingCellService)
Definition MamicoInterfaceProvider.h:58
Definition MultiMDCellService.h:29
interface for continuum/macro fluid solvers for the Couette scenario
Definition CouetteSolver.h:19
virtual tarch::la::Vector< dim, double > getVelocity(tarch::la::Vector< dim, double > pos) const =0
returns the current velocity at the given position
interface to couette solver
Definition CouetteSolverInterface.h:28
implements an analytic Couette flow solver.
Definition CouetteSolver.h:46
implements a simple one-dimensional finite-diffference solver for the Couette flow.
Definition FDCouetteSolver.h:20
Definition IcoFoamInterface.h:34
interface for the LBCouetteSolver
Definition LBCouetteSolverInterface.h:26
implements a three-dimensional Lattice-Boltzmann Couette flow solver.
Definition LBCouetteSolver.h:56
static int getAvgDomainSize(double channelheight, double dx, tarch::la::Vector< 3, unsigned int > processes, int d)
Definition NumericalSolver.h:225
static void parseConfiguration(const std::string filename, const std::string topleveltag, Configuration &config)
Definition ParseConfiguration.h:63
Definition Vector.h:25
Definition MultiMDService.h:30
everything necessary for coupling operations, is defined in here
Definition AdditiveMomentumInsertion.h:15
holds the buffers for the data transfer
Definition CouetteScenario.h:984
coupling::datastructures::FlexibleCellContainer< 3 > macro2MDBuffer
the buffer for data transfer from macro to md
Definition CouetteScenario.h:986
coupling::datastructures::FlexibleCellContainer< 3 > md2macroBuffer
the buffer for data transfer from md to macro
Definition CouetteScenario.h:988
holds all the variables for the time measurement of a simulation
Definition CouetteScenario.h:993
Configuration parameters for Couette flow scenario.
Definition CouetteConfiguration.h:25
int totalNumberMDSimulations
the number of md simulation instances in a multi-instance coupling, -1 = dynamic
Definition CouetteConfiguration.h:326
static CouetteConfig parseCouetteConfiguration(const std::string &filename)
creates CouetteConfig if all elements exist and can be read
Definition CouetteConfiguration.h:51
@ SIMPLEMD
the SimpleMD solver is used
Definition CouetteConfiguration.h:43
@ SYNTHETIC
the synthetic solver is used
Definition CouetteConfiguration.h:44
@ LS1
the LS1 solver is used
Definition CouetteConfiguration.h:45