MaMiCo 1.2
Loading...
Searching...
No Matches
MultiMDCellService.h
1// Copyright (C) 2016 Technische Universitaet Muenchen
2// This file is part of the Mamico project. For conditions of distribution
3// and use, please see the copyright notice in Mamico's main folder, or at
4// www5.in.tum.de/mamico
5#ifndef _MOLECULARDYNAMICS_COUPLING_SERVICES_MULTIMDCELLSERVICE_H_
6#define _MOLECULARDYNAMICS_COUPLING_SERVICES_MULTIMDCELLSERVICE_H_
7
8namespace coupling {
9template <class LinkedCell, unsigned int dim> class InstanceHandling;
10namespace services {
11template <class LinkedCell, unsigned int dim> class MultiMDCellService;
12}
13} // namespace coupling
14
15#include "coupling/InstanceHandling.h"
16#include "coupling/configurations/MaMiCoConfiguration.h"
17#include "coupling/indexing/IndexingService.h"
18#include "coupling/interface/MDSimulationFactory.h"
19#include "coupling/services/CouplingCellService.h"
20#include "coupling/services/CouplingCellServiceDummy.h"
21
29template <class LinkedCell, unsigned int dim> class coupling::services::MultiMDCellService {
30public:
31 MultiMDCellService(std::vector<coupling::interface::MDSolverInterface<LinkedCell, dim>*> mdSolverInterfaces, // MD solver interfaces for each MD simulation
32 // (which uses this rank)
34 simplemd::configurations::MolecularDynamicsConfiguration& mdConfiguration,
35 coupling::configurations::MaMiCoConfiguration<dim>& mamicoConfiguration, const char* filterPipelineConfiguration,
36 tarch::utils::MultiMDService<dim>& multiMDService, int tws = 0)
37 : _multiMDService(multiMDService), _tws(tws), _intNumberProcesses(computeScalarNumberProcesses()), _mdConfiguration(mdConfiguration),
38 _mamicoConfiguration(mamicoConfiguration), _filterPipelineConfiguration(filterPipelineConfiguration),
39 _macroscopicSolverInterface(macroscopicSolverInterface), _postMultiInstanceFilterPipeline(nullptr) {
40 _topologyOffset = computeTopologyOffset();
41 _localNumberMDSimulations = multiMDService.getLocalNumberOfMDSimulations();
42 _totalNumberMDSimulations = multiMDService.getTotalNumberOfMDSimulations();
43
44 _blockOffset = _localNumberMDSimulations * _topologyOffset / _intNumberProcesses;
45
46 _listActiveMDSimulations = std::vector<bool>(_totalNumberMDSimulations, true);
47 _nextFreeBlock = _multiMDService.getNumberLocalComms() - 1;
48 _warmupPhase = std::vector<int>(_totalNumberMDSimulations, 0);
49
50 // If we allow for zero MD instances on ranks, this initializion process
51 // would segfault...
52 if (mdSolverInterfaces.size() == 0) {
53 std::cout << "ERROR: Zero MD instances on rank " << _multiMDService.getGlobalRank() << ". Maybe you forgot to adjust number-md-simulations?" << std::endl;
54 exit(EXIT_FAILURE);
55 }
56
57 // determine globally unique IDs for each coupling cell service.
58 // Assumptions:
59 // - the global parallel topology is subdivided into equally sized blocks of
60 // ranks
61 // - on each block of ranks, the SAME number of MD simulations is executed.
62 // This yields that the variable _topologyOffset together with the local ID
63 // of an MD simulation yields a unique global ID
64 // - an MD simulation executes on one full block of ranks
66 if (_couplingCellServices == NULL) {
67 std::cout << "ERROR "
68 "coupling::services::MultiMDCellService::MultiMDCellService("
69 "...): _couplingCellServices==NULL!"
70 << std::endl;
71 exit(EXIT_FAILURE);
72 }
73
74 // allocate all coupling cell services for macro-only-solver BEFORE
75 // topology offset
76 // -> we have localNumberMDSimulations that run per block of ranks on
77 // intNumberProcesses
78 // -> this yields
79 // localNumberMDSimulations*topologyOffset/intNumberProcesses MD simulations
80 // before the actual block of ranks
81 for (unsigned int i = 0; i < _blockOffset; i++) {
82 _couplingCellServices[i] = createCouplingCellServiceDummy(i, macroscopicSolverInterface, mdConfiguration, _multiMDService, mamicoConfiguration,
83 (i / _localNumberMDSimulations) * _intNumberProcesses);
84 if (_couplingCellServices[i] == NULL) {
85 std::cout << "ERROR "
86 "coupling::services::MultiMDCellService::MultiMDCellServic"
87 "e(...): _couplingCellServices["
88 << i << "]==NULL!" << std::endl;
89 exit(EXIT_FAILURE);
90 }
91 }
92 // allocate all coupling cell services for macro- and micro-interactions
93 for (unsigned int i = _blockOffset; i < _localNumberMDSimulations + _blockOffset; i++) {
95 i, mdSolverInterfaces[i - _blockOffset], macroscopicSolverInterface, mdConfiguration.getMPIConfiguration().getNumberOfProcesses(),
96 _multiMDService.getGlobalRank(), mamicoConfiguration.getParticleInsertionConfiguration(), mamicoConfiguration.getMomentumInsertionConfiguration(),
97 mamicoConfiguration.getBoundaryForceConfiguration(), mamicoConfiguration.getTransferStrategyConfiguration(),
98 mamicoConfiguration.getParallelTopologyConfiguration(), mamicoConfiguration.getThermostatConfiguration(),
99 mdConfiguration.getSimulationConfiguration().getNumberOfTimesteps(), mamicoConfiguration.getCouplingCellConfiguration(),
100 _filterPipelineConfiguration.c_str(), multiMDService, _topologyOffset, _tws);
101 if (_couplingCellServices[i] == NULL) {
102 std::cout << "ERROR "
103 "coupling::services::MultiMDCellService::MultiMDCellServic"
104 "e(...): _couplingCellServices["
105 << i << "]==NULL!" << std::endl;
106 exit(EXIT_FAILURE);
107 }
108 }
109 // allocate all coupling cell services for macro-only-solver AFTER
110 // topology offset
111 for (unsigned int i = _blockOffset + _localNumberMDSimulations; i < _totalNumberMDSimulations; i++) {
112 _couplingCellServices[i] = createCouplingCellServiceDummy(i, macroscopicSolverInterface, mdConfiguration, _multiMDService, mamicoConfiguration,
113 (i / _localNumberMDSimulations) * _intNumberProcesses);
114 if (_couplingCellServices[i] == NULL) {
115 std::cout << "ERROR "
116 "coupling::services::MultiMDCellService::MultiMDCellServic"
117 "e(...): _couplingCellServices["
118 << i << "]==NULL!" << std::endl;
119 exit(EXIT_FAILURE);
120 }
121 }
122
123 _mdConfiguration.getDomainConfigurationNonConst().setInitFromCheckpoint(true);
124 std::stringstream filestem;
125 filestem << "restart_checkpoint_" << (_multiMDService.getGlobalRank()) / computeScalarNumberProcesses() << "_0";
126 _mdConfiguration.getDomainConfigurationNonConst().setCheckpointFilestem(filestem.str());
127 _mdConfiguration.getDomainConfigurationNonConst().setInitFromSequentialCheckpoint(false);
128 }
129
130 ~MultiMDCellService() {
131 for (unsigned int i = 0; i < _totalNumberMDSimulations; i++) {
132 if (_couplingCellServices[i] != nullptr) {
133 delete _couplingCellServices[i];
134 _couplingCellServices[i] = nullptr;
135 }
136 }
137 if (_couplingCellServices != NULL) {
138 delete[] _couplingCellServices;
140 }
141
142 // TODO: fix free bug/possible memory leak here
143 I01 idx;
145 for (auto pair : _couplingCells) {
146 std::tie(cell, idx) = pair;
147 if (cell != nullptr)
148 delete cell;
149 }
152 }
153
158 if (localIndex < _localNumberMDSimulations) {
159 return *(_couplingCellServices[_topologyOffset * _localNumberMDSimulations / _intNumberProcesses + localIndex]);
160 } else {
161 std::cout << "ERROR MultiMDCellService::getCouplingCellService(localIndex): "
162 "localIndex >_localNumberMDSimulations-1!"
163 << std::endl;
164 exit(EXIT_FAILURE);
165 }
166 }
167
169 void plotEveryMacroscopicTimestepforCouplingCellService(unsigned int localIndex, int cycle) {
170 if (_couplingCellServices[_topologyOffset * _localNumberMDSimulations / _intNumberProcesses + localIndex] == nullptr)
171 return;
172 _couplingCellServices[_topologyOffset * _localNumberMDSimulations / _intNumberProcesses + localIndex]->plotEveryMacroscopicTimestep(cycle);
173 }
174
175 void plotEveryMacroscopicTimestep(int cycle) {
176 for (unsigned int i = _blockOffset; i < _blockOffset + _localNumberMDSimulations; ++i) {
177 if (_couplingCellServices[i] == nullptr)
178 continue;
179 _couplingCellServices[i]->plotEveryMacroscopicTimestep(cycle);
180 }
181 }
182
183 void computeAndStoreTemperature(double temp) {
184 for (unsigned int i = _blockOffset; i < _blockOffset + _localNumberMDSimulations; ++i) {
185 if (_couplingCellServices[i] == nullptr)
186 continue;
187 _couplingCellServices[i]->computeAndStoreTemperature(temp);
188 }
189 }
190
193 // just send information to all MD instances. This is currently
194 // sequentialized inside CouplingCellService/SendRecvBuffer
195 for (unsigned int i = 0; i < _totalNumberMDSimulations; i++) {
196 if (_couplingCellServices[i] != nullptr) {
197 _couplingCellServices[i]->sendFromMacro2MD(macro2mdCouplingCellContainer);
198 }
199 }
200 }
201
203 void bcastFromMacro2MD(const std::vector<coupling::datastructures::CouplingCell<dim>*>& couplingCellsFromMacroscopicSolver, const I00* const indices) {
204#if (COUPLING_MD_PARALLEL == COUPLING_MD_NO)
205 // Fall back on sequential operation when MPI is not available (avoids redundant implementation)
206 sendFromMacro2MD(couplingCellsFromMacroscopicSolver, indices);
207 return;
208#else
209
210 std::vector<coupling::sendrecv::DataExchangeFromMacro2MD<dim>*> allDEs(_totalNumberMDSimulations);
211 std::vector<std::vector<coupling::datastructures::CouplingCell<dim>*>> allCouplingCellsFromMamico(_totalNumberMDSimulations);
212 for (unsigned int i = 0; i < _totalNumberMDSimulations; ++i) {
213 if (nullptr == _couplingCellServices[i])
214 continue;
215 allDEs[i] = new coupling::sendrecv::DataExchangeFromMacro2MD<dim>(_macroscopicSolverInterface, _couplingCellServices[i]->getID());
216 if (auto* v = dynamic_cast<CouplingCellServiceImpl<LinkedCell, dim>*>(_couplingCellServices[i])) {
217 allCouplingCellsFromMamico[i] = v->getCouplingCells().getCouplingCells();
218 }
219 }
220
221 for (unsigned int i = 0; i < _totalNumberMDSimulations; ++i) {
222 if (nullptr == _couplingCellServices[i])
223 continue;
224 _couplingCellServices[i]->sendFromMacro2MDPreProcess();
225 }
226
228 fromMacro2MD.bcastFromMacro2MD(allDEs, couplingCellsFromMacroscopicSolver, indices, allCouplingCellsFromMamico);
229
230 for (unsigned int i = 0; i < _totalNumberMDSimulations; ++i) {
231 if (nullptr == _couplingCellServices[i])
232 continue;
233 _couplingCellServices[i]->sendFromMacro2MDPostProcess();
234 }
236 delete de;
237 de = nullptr;
238 }
239#endif
240 }
241
245 /*
246 * On the first coupling step, initialize the reduced couplingCell vector according to `size`
247 */
248 if (!_sumCouplingCells.empty()) {
249 for (unsigned int i = 0; i < _sumCouplingCells.size(); ++i) {
250 _sumCouplingCells[i]->setMacroscopicMass(0);
251 _sumCouplingCells[i]->setMacroscopicMomentum(tarch::la::Vector<dim, double>(0.0));
252 }
253 }
254 for (unsigned int n = 0; n < _totalNumberMDSimulations; ++n) {
255 if (0 != _warmupPhase[n])
256 continue;
257 if (auto* v = dynamic_cast<CouplingCellServiceImpl<LinkedCell, dim>*>(_couplingCellServices[n])) {
258 for (unsigned int i = 0; i < v->getCouplingCells().getCouplingCells().size(); ++i) {
259 if (_sumCouplingCells.size() <= i) {
261 }
262 _sumCouplingCells[i]->addMacroscopicMass(v->getCouplingCells().getCouplingCells()[i]->getMacroscopicMass());
263 _sumCouplingCells[i]->addMacroscopicMomentum(v->getCouplingCells().getCouplingCells()[i]->getMacroscopicMomentum());
264 }
265 }
266 }
267 }
268
275 double reduceFromMD2Macro(const std::vector<coupling::datastructures::CouplingCell<dim>*>& couplingCellsFromMacroscopicSolver, const I00* const indices) {
276#if (COUPLING_MD_PARALLEL == COUPLING_MD_NO)
277 // Fall back on sequential operation when MPI is not available (avoids redundant implementation)
278 return sendFromMD2Macro(couplingCellsFromMacroscopicSolver, indices);
279#else
280 double res = 0;
281 const unsigned int size = (unsigned int)couplingCellsFromMacroscopicSolver.size();
282
283 preprocessingForMD2Macro(couplingCellsFromMacroscopicSolver);
284
285 for (unsigned int i = 0; i < _totalNumberMDSimulations; ++i) {
286 if (nullptr == _couplingCellServices[i] || 0 != _warmupPhase[i])
287 continue;
288 _couplingCellServices[i]->sendFromMD2MacroPreProcess();
289 res += _couplingCellServices[i]->applyFilterPipeline();
290 }
291
296
297 // reset macroscopic data (only those should be used by macroscopic solver anyway) in cells
298 for (coupling::datastructures::CouplingCell<dim>*& couplingCell : _couplingCells) {
299 couplingCell->setMacroscopicMass(0.0);
300 couplingCell->setMacroscopicMomentum(tarch::la::Vector<dim, double>(0.0));
301 }
302
303 std::vector<coupling::sendrecv::DataExchangeFromMD2Macro<dim>*> allDEs(_totalNumberMDSimulations);
304 std::vector<std::vector<coupling::datastructures::CouplingCell<dim>*>> allCouplingCellsFromMamico(_totalNumberMDSimulations);
305 coupling::interface::MacroscopicSolverInterface<dim>* _macroscopicSolverInterface =
307 unsigned int totalNumberEquilibratedMDSimulations = 0;
308 for (unsigned int i = 0; i < _totalNumberMDSimulations; ++i) {
309 if (nullptr == _couplingCellServices[i] || 0 != _warmupPhase[i])
310 continue; // Only reduce using equilibrated MD simulation instances
311 allDEs[i] = new coupling::sendrecv::DataExchangeFromMD2Macro<dim>(_macroscopicSolverInterface, _couplingCellServices[i]->getID());
312 totalNumberEquilibratedMDSimulations += 1;
313 }
314 _fromMD2Macro.reduceFromMD2Macro(allDEs, couplingCellsFromMacroscopicSolver, indices, _sumCouplingCells);
315
316 // average data
317 for (unsigned int i = 0; i < size; i++) {
318 _couplingCells[i]->setMacroscopicMass(couplingCellsFromMacroscopicSolver[i]->getMacroscopicMass() / (double)totalNumberEquilibratedMDSimulations);
319 _couplingCells[i]->setMacroscopicMomentum(1.0 / (double)totalNumberEquilibratedMDSimulations *
320 couplingCellsFromMacroscopicSolver[i]->getMacroscopicMomentum());
321 }
322
323 for (unsigned int i = 0; i < _totalNumberMDSimulations; i++) {
324 if (nullptr == _couplingCellServices[i] || 0 != _warmupPhase[i])
325 continue; // Only equilibrated MD simulation instances
326 _couplingCellServices[i]->sendFromMD2MacroPostProcess();
327 }
328
329 // apply post multi instance FilterPipeline on cell data
330
331#ifdef DEBUG_FILTER_PIPELINE
332 std::cout << "FP: Now applying post-multi-instance filter pipeline" << std::endl;
333#endif
334
335 // FIXME apply filtering to correct region res += (*_postMultiInstanceFilterPipeline)();
336
337 // store data in couplingCellsFromMacroscopicSolver
338 for (unsigned int i = 0; i < size; i++) {
339 couplingCellsFromMacroscopicSolver[i]->setMacroscopicMass(_couplingCells[i]->getMacroscopicMass());
340 couplingCellsFromMacroscopicSolver[i]->setMacroscopicMomentum(_couplingCells[i]->getMacroscopicMomentum());
341 }
342 return res;
343#endif
344 }
345
350 double res = 0;
351
352 preprocessingForMD2Macro(md2macroCouplingCells);
353
354#if (COUPLING_MD_ERROR == COUPLING_MD_YES)
355 if (_couplingCells.size() != md2macroCouplingCells.size())
356 throw std::runtime_error(std::string("Buffers must have the same size!"));
357#endif
358
359 // reset macroscopic data (only those should be used by macroscopic solver
360 // anyway) in duplicate
361 I01 idx;
364 for (auto pair : _couplingCells) {
365 std::tie(cellA, idx) = pair;
366 cellA->setMacroscopicMass(0.0);
368 }
369
370 // receive data from each MD simulation and accumulate information in
371 // duplicate
372 unsigned int totalNumberEquilibratedMDSimulations = 0;
373 for (unsigned int l = 0; l < _totalNumberMDSimulations; l++) {
374 // std::cout << "Rank " <<
375 //_couplingCellServices[l]->getIndexConversion().getThisRank() << ":
376 // Send from MD to Macro for Simulation no. " << l << std::endl;
377 if (_couplingCellServices[l] != nullptr && _warmupPhase[l] == 0) {
378 res += _couplingCellServices[l]->sendFromMD2Macro(md2macroCouplingCells);
379 auto itCouplingCells = _couplingCells.begin();
380 auto itMd2macroCouplingCells = md2macroCouplingCells.begin();
381 while (itCouplingCells != _couplingCells.end()) {
382 std::tie(cellA, idx) = *itCouplingCells;
383 std::tie(cellB, idx) = *itMd2macroCouplingCells;
384 cellA->addMacroscopicMass(cellB->getMacroscopicMass());
386 itCouplingCells++;
387 itMd2macroCouplingCells++;
388 }
389 totalNumberEquilibratedMDSimulations += 1;
390 }
391 }
392
393 // average data
394 for (auto pair : _couplingCells) {
395 std::tie(cellA, idx) = pair;
396 cellA->setMacroscopicMass(cellA->getMacroscopicMass() / (double)totalNumberEquilibratedMDSimulations);
397 cellA->setMacroscopicMomentum(1.0 / (double)totalNumberEquilibratedMDSimulations * cellA->getMacroscopicMomentum());
398 }
399
400 // apply post multi instance FilterPipeline on cell data
401
402#ifdef DEBUG_FILTER_PIPELINE
403 std::cout << "FP: Now applying post-multi-instance filter pipeline" << std::endl;
404#endif
405
406 // FIXME apply filtering to correct region res += (*_postMultiInstanceFilterPipeline)();
407
408 // store data in md2macroCouplingCells
409 auto itCouplingCells = _couplingCells.begin();
410 auto itMd2macroCouplingCells = md2macroCouplingCells.begin();
411 while (itCouplingCells != _couplingCells.end()) {
412 std::tie(cellA, idx) = *itCouplingCells;
413 std::tie(cellB, idx) = *itMd2macroCouplingCells;
414 cellB->setMacroscopicMass(cellA->getMacroscopicMass());
416 itCouplingCells++;
417 itMd2macroCouplingCells++;
418 }
419
420 return res;
421 }
422
427
428 /*
429 * If this is first coupling step, we must allocate space for the
430 * coupling cells we filter and determine averages with.
431 */
432 if (_couplingCells.size() == 0) {
433 // Allocate & init _couplingCells
434 I01 idx;
436 for (auto pair : cells) {
437 std::tie(cell, idx) = pair;
439 _couplingCells << std::make_pair(newCell, idx);
440 }
441 }
442
443 if (_postMultiInstanceFilterPipeline == nullptr) {
444 // Init filter pipeline
445 // FIXME Invalid cell container
446 // _postMultiInstanceFilterPipeline = new coupling::filtering::FilterPipeline<I14, dim>(_couplingCells, coupling::filtering::Scope::postMultiInstance,
447 // _multiMDService, _filterPipelineConfiguration.c_str());
448 }
449 }
450
454 unsigned int rmMDSimulation(coupling::InstanceHandling<LinkedCell, dim>& instanceHandling, const unsigned int& index) {
455#if (COUPLING_MD_DEBUG == COUPLING_MD_YES)
456 if (_couplingCellServices[index] == nullptr) {
457 std::cout << "Rank " << _multiMDService.getGlobalRank() << ": _couplingCellService at " << index << " == NULL" << std::endl;
458 }
459#endif
460
461 delete _couplingCellServices[index];
462 _couplingCellServices[index] = nullptr;
463
464 if (index >= _blockOffset && index < _blockOffset + _localNumberMDSimulations) {
465 unsigned int iSim = index - _blockOffset;
466#if (COUPLING_MD_DEBUG == COUPLING_MD_YES)
467 if (_multiMDService.getGlobalRank() == 0)
468 std::cout << "MultiMDCellService: removing instance " << iSim << "..." << std::endl;
469#endif
470 instanceHandling.rmMDSimulation(iSim);
471 }
472
473 _warmupPhase[index] = -1;
474
475 return index;
476 }
477
482
483 unsigned int newLocalNumberMDSimulations = _localNumberMDSimulations - 1;
484 unsigned int newTotalNumberMDSimulations = _totalNumberMDSimulations - _multiMDService.getNumberLocalComms();
485 unsigned int newBlockOffset = newLocalNumberMDSimulations * _topologyOffset / _intNumberProcesses;
486
487 auto** newCouplingCellServices = new CouplingCellService<dim>*[newTotalNumberMDSimulations];
488
489 for (unsigned int i = 0; i < _multiMDService.getNumberLocalComms(); ++i) {
490 for (unsigned int j = 0; j < newLocalNumberMDSimulations; ++j) {
491 unsigned int index = i * _localNumberMDSimulations + j;
492 unsigned int newIndex = i * newLocalNumberMDSimulations + j;
493
494 newCouplingCellServices[newIndex] = _couplingCellServices[index];
495 }
496
497 auto pos = _warmupPhase.begin() + newLocalNumberMDSimulations * i - i;
498 _warmupPhase.erase(pos);
499 }
500
501 // Update local variables
502 _localNumberMDSimulations = newLocalNumberMDSimulations;
503 _totalNumberMDSimulations = newTotalNumberMDSimulations;
504 _blockOffset = newBlockOffset;
505
506 delete[] _couplingCellServices;
507 _couplingCellServices = newCouplingCellServices;
508
509 //_multiMDService.removeMDSimulationBlock();
510
511 //_listActiveMDSimulations.pop_back();
512 }
513
522
523 unsigned int newLocalNumberMDSimulations = _localNumberMDSimulations + 1;
524 unsigned int newTotalNumberMDSimulations = _totalNumberMDSimulations + _multiMDService.getNumberLocalComms();
525 unsigned int newBlockOffset = newLocalNumberMDSimulations * _topologyOffset / _intNumberProcesses;
526
527 auto** newCouplingCellServices = new CouplingCellService<dim>*[newTotalNumberMDSimulations];
528 for (unsigned int i = 0; i < _multiMDService.getNumberLocalComms(); ++i) {
529 for (unsigned int j = 0; j < _localNumberMDSimulations; ++j) {
530 unsigned int index = i * _localNumberMDSimulations + j;
531 unsigned int newIndex = i * newLocalNumberMDSimulations + j;
532 newCouplingCellServices[newIndex] = _couplingCellServices[index];
533 }
534 newCouplingCellServices[(i + 1) * newLocalNumberMDSimulations - 1] = nullptr;
535
536 auto pos = _warmupPhase.begin() + ((i + 1) * newLocalNumberMDSimulations - 1);
537 _warmupPhase.insert(pos, 0);
538 }
539
540 // Update local variables
541 _localNumberMDSimulations = newLocalNumberMDSimulations;
542 _totalNumberMDSimulations = newTotalNumberMDSimulations;
543 _blockOffset = newBlockOffset;
544 delete[] _couplingCellServices;
545 _couplingCellServices = newCouplingCellServices;
546 }
547
553 coupling::interface::MacroscopicSolverInterface<dim>* macroscopicSolverInterface, const unsigned int& slot) {
554
555 if (slot >= _totalNumberMDSimulations) {
556 std::cout << "ERROR coupling::services::MultiMDCellService::addMDSimulation(): "
557 "Invalid slot "
558 << slot << "!" << std::endl;
559 std::exit(EXIT_FAILURE);
560 }
561 if (_couplingCellServices[slot] != nullptr) {
562 std::cout << "ERROR! "
563 "coupling::services::MultiMDCellService::addMDSimulation(): "
564 "Simulation at "
565 << slot << " already exists!" << std::endl;
566 std::exit(EXIT_FAILURE);
567 }
568
569 std::stringstream filestem;
570 filestem << "restart_checkpoint_" << _multiMDService.getGlobalRank() / computeScalarNumberProcesses();
571 instanceHandling.writeCheckpoint(filestem.str().c_str(), 0);
572
573 if (slot < _blockOffset || slot >= _blockOffset + _localNumberMDSimulations) {
574
576 slot, macroscopicSolverInterface, _mdConfiguration.getMPIConfiguration().getNumberOfProcesses(), _multiMDService.getGlobalRank(),
577 _mdConfiguration.getDomainConfiguration().getGlobalDomainSize(), _mdConfiguration.getDomainConfiguration().getGlobalDomainOffset(),
578 _mamicoConfiguration.getParallelTopologyConfiguration(), _mamicoConfiguration.getCouplingCellConfiguration(),
579 (slot / _localNumberMDSimulations) * _intNumberProcesses);
580 } else {
581
582 unsigned localIndex = slot - _blockOffset;
583 auto* mdSolverInterface = instanceHandling.addMDSimulation(slot, localIndex);
584
586 slot, mdSolverInterface, macroscopicSolverInterface, _mdConfiguration.getMPIConfiguration().getNumberOfProcesses(), _multiMDService.getGlobalRank(),
587 _mamicoConfiguration.getParticleInsertionConfiguration(), _mamicoConfiguration.getMomentumInsertionConfiguration(),
588 _mamicoConfiguration.getBoundaryForceConfiguration(), _mamicoConfiguration.getTransferStrategyConfiguration(),
589 _mamicoConfiguration.getParallelTopologyConfiguration(), _mamicoConfiguration.getThermostatConfiguration(),
590 _mdConfiguration.getSimulationConfiguration().getNumberOfTimesteps(), _mamicoConfiguration.getCouplingCellConfiguration(),
591 _filterPipelineConfiguration.c_str(), _multiMDService, _topologyOffset, _tws);
592 instanceHandling.getSimpleMD()[localIndex]->setCouplingCellService((_couplingCellServices[slot]));
593 _couplingCellServices[slot]->perturbateVelocity();
594 if (_couplingCellServices[slot]->getFilterPipeline() == nullptr) {
595 _couplingCellServices[slot]->initFiltering();
596 }
597 }
598
599 _warmupPhase[slot] = 10;
600
601 return slot;
602 }
603
604 unsigned int getLocalNumberOfMDSimulations() const { return _localNumberMDSimulations; }
605
606 void finishCycle(const unsigned int& cycle, coupling::InstanceHandling<LinkedCell, dim>& instanceHandling) {
607 // TODO this should be move to somewhere else (MDMediator?)
608 for (unsigned int i = 0; i < _totalNumberMDSimulations; ++i) {
609 if (_warmupPhase[i] > 0) {
610 _warmupPhase[i] -= 1;
611 if (_warmupPhase[i] == 0 && i >= _blockOffset && i < _blockOffset + _localNumberMDSimulations) {
612 instanceHandling.switchOnCoupling(i - _blockOffset);
613 }
614 }
615 }
616 // writeCheckpoint(cycle, instanceHandling);
617 // TODO call directly instanceHandling.writeCheckpoint();
618 }
619
620 void writeCheckpoint(const unsigned int& cycle, const coupling::InstanceHandling<LinkedCell, dim>& instanceHandling) const {
621 if (cycle % 10 != 0)
622 return;
623 std::stringstream filestem;
624 filestem << "restart_checkpoint_" << _multiMDService.getGlobalRank() / computeScalarNumberProcesses();
625 instanceHandling.writeCheckpoint(filestem.str().c_str(), 0);
626 }
627
628 // Must be called after construction of MultiMDCellService, but before the
629 // first coupling step.
630 void constructFilterPipelines() {
631 for (unsigned int md = 0; md < _totalNumberMDSimulations; md++) {
632 // get coupling cell service of instance
633 auto& mcs = _couplingCellServices[md];
634
635 // only Impl instances of MCS contain filtering
636 if (dynamic_cast<coupling::services::CouplingCellServiceImpl<LinkedCell, dim>*>(mcs) == nullptr)
637 continue;
638
639 if (mcs->getFilterPipeline() == nullptr) {
640 mcs->initFiltering();
641 }
642 }
643 }
644
645private:
646 unsigned int computeTopologyOffset() const {
647 // determine topology offset of this rank
648 const unsigned int intNumberProcesses = computeScalarNumberProcesses();
649 const unsigned int topologyOffset = (_multiMDService.getGlobalRank() / intNumberProcesses) * intNumberProcesses;
650 return topologyOffset;
651 }
652
653 unsigned int computeScalarNumberProcesses() const {
654 unsigned int np = _multiMDService.getNumberProcessesPerMDSimulation()[0];
655 for (unsigned int d = 1; d < dim; d++) {
656 np = np * _multiMDService.getNumberProcessesPerMDSimulation()[d];
657 }
658 return np;
659 }
660
661 coupling::services::CouplingCellService<dim>*
662 createCouplingCellServiceDummy(unsigned int ID, coupling::interface::MacroscopicSolverInterface<dim>* macroscopicSolverInterface,
663 simplemd::configurations::MolecularDynamicsConfiguration& mdConfiguration, tarch::utils::MultiMDService<dim>& multiMDService,
664 coupling::configurations::MaMiCoConfiguration<dim>& mamicoConfiguration, unsigned int topologyOffset) const {
665 return new coupling::services::CouplingCellServiceDummy<dim>(
666 ID, macroscopicSolverInterface, mdConfiguration.getMPIConfiguration().getNumberOfProcesses(), multiMDService.getGlobalRank(),
667 mdConfiguration.getDomainConfiguration().getGlobalDomainSize(), mdConfiguration.getDomainConfiguration().getGlobalDomainOffset(),
668 mamicoConfiguration.getParallelTopologyConfiguration(), mamicoConfiguration.getCouplingCellConfiguration(), topologyOffset);
669 }
670
671 unsigned int _localNumberMDSimulations;
678 unsigned int _topologyOffset;
679
680 const int _tws;
681 const unsigned int _intNumberProcesses;
682
683 // conversions during filtering TODO after merge with dynamic-md*/
686
687 simplemd::configurations::MolecularDynamicsConfiguration& _mdConfiguration;
689 const std::string _filterPipelineConfiguration;
690 coupling::interface::MacroscopicSolverInterface<dim>* _macroscopicSolverInterface;
691
692 unsigned int _blockOffset;
693 std::vector<bool> _listActiveMDSimulations;
695 unsigned int _nextFreeBlock;
697 std::vector<int> _warmupPhase;
703 /*
704 * Analogon to CouplingCellService's FilterPipeline.
705 * Is applied during this->sendFromMD2Macro.
706 */
708 /* Buffer for copying data from MD to macro */
710};
711#endif // _MOLECULARDYNAMICS_COUPLING_SERVICES_MULTIMDCELLSERVICE_H_
Simulation slots are managed (i.e., added/removed) via this class. Works and interacts with the class...
Definition InstanceHandling.h:35
auto & getSimpleMD() const
Definition InstanceHandling.h:101
void writeCheckpoint(const std::string &filestem, const unsigned int &T) const
Definition InstanceHandling.h:170
void rmMDSimulation(const unsigned int &index)
Definition InstanceHandling.h:253
void switchOnCoupling()
Definition InstanceHandling.h:127
coupling::interface::MDSolverInterface< LinkedCell, dim > * addMDSimulation(unsigned int slot, unsigned int localIndex)
Definition InstanceHandling.h:225
parses all sub-tags for MaMiCo configuration.
Definition MaMiCoConfiguration.h:31
const coupling::configurations::MomentumInsertionConfiguration & getMomentumInsertionConfiguration() const
Definition MaMiCoConfiguration.h:86
const coupling::configurations::ThermostatConfiguration & getThermostatConfiguration() const
Definition MaMiCoConfiguration.h:146
const coupling::configurations::BoundaryForceConfiguration< dim > & getBoundaryForceConfiguration() const
Definition MaMiCoConfiguration.h:99
const coupling::configurations::ParticleInsertionConfiguration & getParticleInsertionConfiguration() const
Definition MaMiCoConfiguration.h:73
const coupling::configurations::ParallelTopologyConfiguration & getParallelTopologyConfiguration() const
Definition MaMiCoConfiguration.h:125
const coupling::configurations::TransferStrategyConfiguration< dim > & getTransferStrategyConfiguration() const
Definition MaMiCoConfiguration.h:112
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
void addMacroscopicMomentum(const tarch::la::Vector< dim, double > &momentum)
Definition CouplingCell.h:85
const tarch::la::Vector< dim, double > & getMacroscopicMomentum() const
Definition CouplingCell.h:64
void addMacroscopicMass(const double &mass)
Definition CouplingCell.h:82
void setMacroscopicMomentum(const tarch::la::Vector< dim, double > &momentum)
Definition CouplingCell.h:61
const double & getMacroscopicMass() const
Definition CouplingCell.h:58
void setMacroscopicMass(const double &mass)
Definition CouplingCell.h:55
provides access to coupling cells, which may belong to different indexing domains
Definition FlexibleCellContainer.h:30
Iterator begin() const
Definition FlexibleCellContainer.h:121
unsigned int size() const
Definition FlexibleCellContainer.h:78
Definition FilterPipeline.h:44
interface to the MD simulation
Definition MDSolverInterface.h:25
interface for the macroscopic, i.e. continuum solver
Definition MacroscopicSolverInterface.h:23
static MamicoInterfaceProvider & getInstance()
Definition MamicoInterfaceProvider.h:28
coupling::interface::MacroscopicSolverInterface< dim > * getMacroscopicSolverInterface()
Definition MamicoInterfaceProvider.h:43
data exchange from the MD solver to the macroscopic solver. Derived from the class coupling::sendrecv...
Definition DataExchangeFromMD2Macro.h:36
data exchange from the macroscopic solver to the MD solver. Derived from the class coupling::sendrecv...
Definition DataExchangeFromMacro2MD.h:36
sends coupling cell information from MaMiCo to the macroscopic solver. Derived from the class couplin...
Definition FromMD2Macro.h:29
SendReceiveBuffer for transfer of quantities from a macroscopic solver to the coupling cells....
Definition FromMacro2MD.h:28
Definition CouplingCellServiceDummy.h:28
Definition CouplingCellService.h:98
Definition CouplingCellService.h:49
Definition MultiMDCellService.h:29
const int _tws
Definition MultiMDCellService.h:680
coupling::filtering::FilterPipeline< I14, dim > * _postMultiInstanceFilterPipeline
Definition MultiMDCellService.h:707
void sendFromMacro2MD(const coupling::datastructures::FlexibleCellContainer< dim > &macro2mdCouplingCellContainer)
Definition MultiMDCellService.h:192
void removeSimulationBlock()
Definition MultiMDCellService.h:481
void bcastFromMacro2MD(const std::vector< coupling::datastructures::CouplingCell< dim > * > &couplingCellsFromMacroscopicSolver, const I00 *const indices)
Definition MultiMDCellService.h:203
void sumUpCouplingCellsFromMamico()
Definition MultiMDCellService.h:244
double reduceFromMD2Macro(const std::vector< coupling::datastructures::CouplingCell< dim > * > &couplingCellsFromMacroscopicSolver, const I00 *const indices)
Definition MultiMDCellService.h:275
std::vector< int > _warmupPhase
Definition MultiMDCellService.h:697
double sendFromMD2Macro(const coupling::datastructures::FlexibleCellContainer< dim > &md2macroCouplingCells)
Definition MultiMDCellService.h:349
simplemd::configurations::MolecularDynamicsConfiguration & _mdConfiguration
Definition MultiMDCellService.h:687
coupling::services::CouplingCellService< dim > & getCouplingCellService(unsigned int localIndex)
Definition MultiMDCellService.h:157
coupling::services::CouplingCellService< dim > ** _couplingCellServices
Definition MultiMDCellService.h:675
unsigned int rmMDSimulation(coupling::InstanceHandling< LinkedCell, dim > &instanceHandling, const unsigned int &index)
Definition MultiMDCellService.h:454
unsigned int addMDSimulation(coupling::InstanceHandling< LinkedCell, dim > &instanceHandling, coupling::interface::MacroscopicSolverInterface< dim > *macroscopicSolverInterface, const unsigned int &slot)
Definition MultiMDCellService.h:552
unsigned int _totalNumberMDSimulations
Definition MultiMDCellService.h:674
void addSimulationBlock()
Definition MultiMDCellService.h:521
void plotEveryMacroscopicTimestepforCouplingCellService(unsigned int localIndex, int cycle)
Definition MultiMDCellService.h:169
unsigned int _nextFreeBlock
Definition MultiMDCellService.h:695
void preprocessingForMD2Macro(const coupling::datastructures::FlexibleCellContainer< dim > &cells)
Definition MultiMDCellService.h:426
tarch::utils::MultiMDService< dim > & _multiMDService
Definition MultiMDCellService.h:677
coupling::datastructures::FlexibleCellContainer< dim > _sumCouplingCells
Definition MultiMDCellService.h:685
Definition Vector.h:24
Definition MultiMDService.h:30
everything necessary for coupling operations, is defined in here
Definition AdditiveMomentumInsertion.h:15