4#ifndef _COUPLING_SCENARIO_COUETTESCENARIO_H_
5#define _COUPLING_SCENARIO_COUETTESCENARIO_H_
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"
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"
31#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
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"
48#define SYNTHETICMD_SEQUENCE "SYNTHETIC-MD"
84 _timeIntegrationService->run(
_cfg.couplingCycles);
91#if defined(LS1_MARDYN)
92 Log::global_log = std::make_unique<Log::Logger>(Log::Error);
93#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
94 global_log->set_mpi_output_root(0);
112 if (
_cfg.totalNumberMDSimulations < 0)
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);
139 const int EQUI_CYCLES = 3;
140 for (
int i = 0; i < EQUI_CYCLES; i++) {
156 for (
int cycle = 0; cycle <
_cfg.couplingCycles; cycle++)
165 std::string filename(
"couette.xml");
171 std::cout <<
"ERROR CouetteScenario: Invalid SimpleMD config!" << std::endl;
178 std::cout <<
"ERROR CouetteScenario: Invalid MaMiCo config!" << std::endl;
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!"
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]);
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);
209 _timeIntegrationService = std::make_unique<coupling::services::ParallelTimeIntegrationService<3>>(
_mamicoConfig,
this);
210 _rank = _timeIntegrationService->getRank();
213 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[0],
_simpleMDConfig.getDomainConfiguration().getGlobalDomainOffset(),
215 _mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap());
219 simplemd::configurations::DomainDecompConfiguration::DecompositionType::STATIC_IRREG_RECT_GRID) {
220 coupling::indexing::IndexingService<3>::getInstance().initWithMDSize(
223 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize(),
_mamicoConfig.getParallelTopologyConfiguration().getParallelTopologyType(),
224 _mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap(), (
unsigned int)
_rank
225#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
227 _timeIntegrationService->getPintComm()
231 coupling::indexing::IndexingService<3>::getInstance().initWithMDSize(
234 _mamicoConfig.getParallelTopologyConfiguration().getParallelTopologyType(),
_mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap(),
236#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
238 _timeIntegrationService->getPintComm()
251 unsigned int totNumMD;
252 if (
_cfg.totalNumberMDSimulations > 0)
253 totNumMD =
_cfg.totalNumberMDSimulations;
255 totNumMD =
_cfg.lowerBoundNumberMDSimulations;
258#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
260 _timeIntegrationService->getPintComm()
265 if (_instanceHandling ==
nullptr) {
267 std::cout <<
"ERROR CouetteScenario::initSolvers() : _instanceHandling == NULL!" << std::endl;
269 std::exit(EXIT_FAILURE);
274 gettimeofday(&
_tv.start, NULL);
275 gettimeofday(&
_tv.output, NULL);
280 _instanceHandling->switchOffCoupling();
282 _instanceHandling->switchOnCoupling();
287 _instanceHandling->setMDSolverInterface();
312 _MDBoundarySetupDone =
false;
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!");
365 ->getSequence(SYNTHETICMD_SEQUENCE)
367 new std::function<std::vector<double>(std::vector<double>, std::vector<std::array<unsigned int, 3>>)>{
369 [
this](std::vector<double> inputScalars,
373 std::vector<std::array<unsigned int, 3>> cellIndices
376 gettimeofday(&
_tv.start, NULL);
381 const unsigned int size = inputScalars.size();
382 auto dx = coupling::indexing::IndexingService<3>::getInstance().getCouplingCellSize();
384 const double mass = (
_cfg.density) * couplingCellSize[0] * couplingCellSize[1] * couplingCellSize[2];
386 std::vector<double> syntheticMasses;
387 for (
unsigned int i = 0; i < size; i++) {
388 syntheticMasses.push_back(mass);
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);
397 return syntheticMasses;
399 new std::function<std::vector<std::array<double, 3>>(std::vector<std::array<double, 3>>, std::vector<std::array<unsigned int, 3>>)>{
401 [
this](std::vector<std::array<double, 3>> inputVectors,
402 std::vector<std::array<unsigned int, 3>> cellIndices
405 gettimeofday(&
_tv.start, NULL);
410 const unsigned int size = inputVectors.size();
411 auto dx = coupling::indexing::IndexingService<3>::getInstance().getCouplingCellSize();
413 const double mass = (
_cfg.density) * couplingCellSize[0] * couplingCellSize[1] * couplingCellSize[2];
415 using coupling::indexing::IndexTrait;
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],
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++) {
428 CellIndex<3, IndexTrait::vector> globalIndex(
429 CellIndex<3, IndexTrait::local, IndexTrait::md2macro, IndexTrait::noGhost>{i});
432 for (
unsigned int d = 0; d < 3; d++) {
433 cellMidPoint[d] = cellMidPoint[d] + ((double)globalIndex.get()[d]) * couplingCellSize[d];
441 syntheticMomenta.push_back({momentum[0], momentum[1], momentum[2]});
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);
450 return syntheticMomenta;
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 "
459 << SYNTHETICMD_SEQUENCE <<
"' in config." << std::endl;
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;
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"
491 gettimeofday(&
_tv.start_total, NULL);
494 std::cout <<
"Finish CouetteScenario::initSolvers() " << std::endl;
504 gettimeofday(&
_tv.start, NULL);
508 if (
_cfg.wallInitCycles > 0 && cycle ==
_cfg.wallInitCycles) {
512 if (
_cfg.wallOscillations != 0) {
514 vel = vel * cos(2 * M_PI *
_cfg.wallOscillations * cycle /
_cfg.couplingCycles);
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);
529#ifdef MAMICO_USE_COLLECTIVE_MPI
539 void varyMD(
int cycle) {
548 int addMDInstances = 0;
556 std::tie(couplingCell, idx) = pair;
564 (1 / std::sqrt(3)) * (
_mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[0] /
566 double cellVolume =
_mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[0] *
567 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[1] *
568 _mamicoConfig.getCouplingCellConfiguration().getCouplingCellSize()[2];
570 coupling::error::ErrorEstimation errorControl(vel[0],
_cfg.temp, mass,
_simpleMDConfig.getMoleculeConfiguration().getMass(), soundSpeed,
573 double simtime = (double)cycle /
_cfg.couplingCycles;
574 double targetError =
_cfg.absVelErrEnd * simtime +
_cfg.absVelErrStart * (1 - simtime);
575 errorControl.setAbsVelocityError(targetError);
578 addMDInstances = NoMD -
_multiMDService->getTotalNumberOfMDSimulations();
581#if (COUPLING_MD_PARALLEL == COUPLING_MD_YES)
582 MPI_Bcast(&addMDInstances, 1, MPI_INT, 0, _timeIntegrationService->getPintComm());
584 if (addMDInstances < 0) {
585 if ((
int)_multiMDMediator->getNumberOfActiveMDSimulations() + addMDInstances <
_cfg.lowerBoundNumberMDSimulations) {
586 addMDInstances = (
_cfg.lowerBoundNumberMDSimulations - (int)_multiMDMediator->getNumberOfActiveMDSimulations());
589 if (addMDInstances < 0) {
591 std::cout <<
"Removing " << -addMDInstances <<
" of " <<
_multiMDService->getTotalNumberOfMDSimulations() <<
" MD simulations" << std::endl;
592 _multiMDMediator->rmNMDSimulations(-addMDInstances);
593 }
else if (addMDInstances > 0) {
595 std::cout <<
"Adding " << addMDInstances <<
" to " <<
_multiMDService->getTotalNumberOfMDSimulations() <<
" MD simulations" << std::endl;
596 _multiMDMediator->addNMDSimulations(addMDInstances);
605 gettimeofday(&
_tv.start, NULL);
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);
622#ifdef MAMICO_USE_COLLECTIVE_MPI
637#ifdef MAMICO_USE_COLLECTIVE_MPI
648 if (
_cfg.computeSNR && cycle >=
_cfg.filterInitCycles) {
649 std::cout << cycle -
_cfg.filterInitCycles <<
", ";
650 using namespace coupling::indexing;
652 const double mass =
_cfg.density * dx[0] * dx[1] * dx[2];
659 std::tie(couplingCell, idx) = pair;
660 auto midPoint = idx.getCellMidPoint();
663 _sum_noise += (vx_macro - vx_filter) * (vx_macro - vx_filter);
665 std::cout << vx_macro <<
", " << vx_filter << std::endl;
677 if (
_cfg.twoWayCoupling) {
678#if (BUILD_WITH_OPENFOAM)
679 if ((
_cfg.maSolverType == CouetteConfig::COUETTE_FOAM) && cycle ==
_cfg.filterInitCycles &&
_couetteSolver != NULL) {
683 if ((
_cfg.maSolverType == CouetteConfig::COUETTE_LB ||
_cfg.maSolverType == CouetteConfig::COUETTE_FD) && cycle >=
_cfg.filterInitCycles) {
684 if (!_MDBoundarySetupDone) {
687 _mamicoConfig.getMomentumInsertionConfiguration().getInnerOverlap());
688 _MDBoundarySetupDone =
true;
692#if (BUILD_WITH_OPENFOAM)
693 else if (
_cfg.maSolverType == CouetteConfig::COUETTE_FOAM && cycle >=
_cfg.filterInitCycles &&
_couetteSolver != NULL) {
705 if (
_cfg.computeSNR) {
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;
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;
721 if (_instanceHandling !=
nullptr) {
722 delete _instanceHandling;
732 if (couetteSolverInterface != NULL) {
733 delete couetteSolverInterface;
734 couetteSolverInterface = NULL;
744 if (_multiMDMediator !=
nullptr) {
745 delete _multiMDMediator;
746 _multiMDMediator =
nullptr;
750 std::cout <<
"Finish CouetteScenario::shutdown() " << std::endl;
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;
764 return globalNumberCouplingCells;
771 for (
auto idx : I08()) {
776 if (couplingCell ==
nullptr)
777 throw std::runtime_error(std::string(
"ERROR CouetteScenario::allocateMacro2mdBuffer: couplingCell==NULL!"));
786 for (
auto idx : I12()) {
790 if (couplingCell ==
nullptr)
791 throw std::runtime_error(std::string(
"ERROR CouetteScenario::allocateMacro2mdBuffer: couplingCell==NULL!"));
800 template <
class Container_T>
void write2CSV(Container_T& md2macroBuffer,
int couplingCycle)
const {
801 if (md2macroBuffer.size() == 0)
803 if (
_cfg.csvEveryTimestep < 1 || couplingCycle %
_cfg.csvEveryTimestep > 0)
806 std::stringstream ss;
807 ss <<
"CouetteAvgMultiMDCells_d" << _timeIntegrationService->getPintDomain() <<
"_r" <<
_rank <<
"_c" << couplingCycle;
808 if (_timeIntegrationService->isPintEnabled())
809 ss <<
"_i" << _timeIntegrationService->getIteration();
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;
819 file <<
"I01_x,I01_y,I01_z,vel_x,vel_y,vel_z,mass" << std::endl;
820 for (
auto pair : md2macroBuffer) {
823 std::tie(couplingCell, idx) = pair;
828 file << idx <<
"," << vel[0] <<
"," << vel[1] <<
"," << vel[2] <<
"," << couplingCell->
getMacroscopicMass() << std::endl;
844 using namespace coupling::indexing;
846 double mass = density * dx[0] * dx[1] * dx[2];
847 for (
auto pair : macro2MDBuffer) {
850 std::tie(couplingCell, idx) = pair;
851 auto midPoint = idx.getCellMidPoint();
852 if (
_cfg.maSolverType == CouetteConfig::COUETTE_LB ||
_cfg.maSolverType == CouetteConfig::COUETTE_FD)
869 if (
_cfg.maSolverType == CouetteConfig::COUETTE_ANALYTICAL) {
870 if (
_rank == 0 ||
_cfg.miSolverType == CouetteConfig::SYNTHETIC) {
872 if (solver == NULL) {
873 std::cout <<
"ERROR CouetteScenario::getCouetteSolver(): Analytic solver==NULL!" << std::endl;
878#if (BUILD_WITH_OPENFOAM)
879 else if (
_cfg.maSolverType == CouetteConfig::COUETTE_FOAM) {
881 _cfg.foam.boundariesWithMD,
_cfg.wallVelocity);
882 if (solver == NULL) {
883 std::cout <<
"ERROR CouetteScenario::getCouetteSolver(): IcoFoamInterface solver==NULL!" << std::endl;
889 else if (
_cfg.maSolverType == CouetteConfig::COUETTE_LB) {
891 "LBCouette",
_cfg.lbNumberProcesses, 1,
this);
892 if (solver == NULL) {
893 std::cout <<
"ERROR CouetteScenario::getCouetteSolver(): LB solver==NULL!" << std::endl;
896 }
else if (
_cfg.maSolverType == CouetteConfig::COUETTE_FD) {
898 _cfg.lbNumberProcesses);
899 if (solver == NULL) {
900 std::cout <<
"ERROR CouetteScenario::getCouetteSolver(): FD solver==NULL!" << std::endl;
904 std::cout <<
"ERROR CouetteScenario::getCouetteSolver(): Unknown solver type!" << std::endl;
921 unsigned int outerRegion) {
924 if (
_cfg.maSolverType == CouetteConfig::COUETTE_ANALYTICAL) {
926 } else if (_cfg.maSolverType == CouetteConfig::COUETTE_LB) {
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;
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;
947#if (BUILD_WITH_OPENFOAM)
948 else if (_cfg.maSolverType == CouetteConfig::COUETTE_FOAM) {
952 else if (_cfg.maSolverType == CouetteConfig::COUETTE_FD) {
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;
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;
974 if (interface == NULL) {
975 std::cout <<
"ERROR CouetteScenario::getCouetteSolverInterface(...), rank=" <<
_rank <<
": interface==NULL!" << std::endl;
1041 bool _MDBoundarySetupDone;
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 > ¯o2MDBuffer) 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
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 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