ls1-MarDyn
ls1-MarDyn molecular dynamics code
DomainDecompBase.h
1#ifndef DOMAINDECOMPBASE_H_
2#define DOMAINDECOMPBASE_H_
3
4#include "parallel/CollectiveCommBase.h"
5#include <string>
6
7#ifdef ENABLE_MPI
8#include <mpi.h>
9#endif
10#include <particleContainer/RegionParticleIterator.h>
11#include <iostream>
12#include <variant>
13
14#include "HaloRegion.h"
15#include "io/MemoryProfiler.h"
16#include "molecules/MoleculeForwardDeclaration.h"
17#include "utils/Logger.h" // is this used?
18
19class Component;
20class Domain;
22class XMLfileUnits;
23
55public:
58
60 ~DomainDecompBase() override;
61
62 virtual void readXML(XMLfileUnits& xmlconfig);
63
72 void exchangeMolecules(ParticleContainer* moleculeContainer, Domain* domain);
73
79 virtual void exchangeForces(ParticleContainer* moleculeContainer, Domain* domain);
80
87 virtual int getNonBlockingStageCount();
88
97 virtual bool queryBalanceAndExchangeNonBlocking(bool forceRebalancing, ParticleContainer* moleculeContainer, Domain* domain, double etime);
98
110 virtual void balanceAndExchange(double lastTraversalTime, bool forceRebalancing, ParticleContainer* moleculeContainer, Domain* domain);
111
121 virtual bool procOwnsPos(double x, double y, double z, Domain* domain) final;
122
127 void getBoundingBoxMinMax(Domain* domain, double* min, double* max);
128
132 virtual double getBoundingBoxMin(int dimension, Domain* domain);
133
137 virtual double getBoundingBoxMax(int dimension, Domain* domain);
138
143 virtual void printDecomp(const std::string &filename, Domain *domain, ParticleContainer *particleContainer);
144
145
148 virtual int getRank() const;
149
152 virtual int getNumProcs() const;
153
155 virtual void barrier() const;
156
158 virtual double getTime() const;
159
161 virtual unsigned Ndistribution(unsigned localN, float* minrnd, float* maxrnd);
162
164 virtual void assertIntIdentity(int IX);
165 virtual void assertDisjunctivity(ParticleContainer* moleculeContainer) const;
166
175 double getIOCutoffRadius(int dim, Domain* domain, ParticleContainer* moleculeContainer);
176
177
178#ifdef ENABLE_MPI
183 void writeMoleculesToMPIFileBinary(const std::string& filename, ParticleContainer* moleculeContainer) const;
184#endif // ENABLE_MPI
185
193 void writeMoleculesToFile(const std::string& filename, ParticleContainer* moleculeContainer, bool binary = false) const;
194
195
196 void updateSendLeavingWithCopies(bool sendTogether){
197 using Log::global_log;
198 // Count all processes that need to send separately
199 collCommInit(1);
200 collCommAppendInt(!sendTogether);
202 _sendLeavingAndCopiesSeparately = collCommGetInt();
204
205
206 global_log->info() << "Sending leaving particles and halo copies "
207 << (sendLeavingWithCopies() ? "together" : "separately") << std::endl;
208 }
209
210 bool sendLeavingWithCopies() const{
211 // No process needs to send separately => send together
212 return _sendLeavingAndCopiesSeparately == 0;
213 }
214
215
216 //##################################################################
217 // The following methods with prefix "collComm" are all used
218 // in the context of collective communication. Each of the methods
219 // basically has to call the corresponding method from the class
220 // CollectiveCommunication (or CollectiveCommDummy in the sequential
221 // case). To get information about how to use this methods, read
222 // the documentation of the class CollectiveCommunication.
223 //##################################################################
225 virtual void collCommInit(int numValues, int key=0);
227 virtual void collCommFinalize();
229 virtual void collCommAppendInt(int intValue);
231 virtual void collCommAppendUnsLong(unsigned long unsLongValue);
233 virtual void collCommAppendFloat(float floatValue);
235 virtual void collCommAppendDouble(double doubleValue);
237 virtual void collCommAppendLongDouble(long double longDoubleValue);
239 virtual int collCommGetInt();
241 virtual unsigned long collCommGetUnsLong();
243 virtual float collCommGetFloat();
245 virtual double collCommGetDouble();
247 virtual long double collCommGetLongDouble();
249 virtual void collCommAllreduceSum();
253 virtual void collCommAllreduceCustom(ReduceType type);
255 virtual void collCommScanSum();
257 virtual void collCommBroadcast(int root = 0);
258 //returns the ranks of the neighbours
259 virtual std::vector<int> getNeighbourRanks(){
260 return std::vector<int>(0);
261 }
262 //returns the ranks of the neighbours
263 virtual std::vector<int> getNeighbourRanksFullShell(){
264 std::cout << "Not yet implemented";
265 return std::vector<int>(0);
266 }
267
268 //returns the ranks of all ranks
269 virtual std::vector<std::vector<std::vector<int>>> getAllRanks(){
270 std::cout << "Not yet implemented";
271 return std::vector<std::vector<std::vector<int>>>(0);
272 }
273#if defined(ENABLE_MPI)
274 virtual MPI_Comm getCommunicator(){
275 return MPI_COMM_WORLD;
276 }
277#endif
278
279 virtual size_t getTotalSize() override {
280 return _collCommBase.getTotalSize();
281 }
282 virtual void printSubInfo(int offset) override {
283 return;
284 }
285 virtual std::string getName() override {
286 return "DomainDecompBase";
287 }
288
289 virtual void printCommunicationPartners(std::string filename) const {};
290
291protected:
292 void addLeavingMolecules(std::vector<Molecule>&& invalidMolecules, ParticleContainer* moleculeContainer);
293
300 void handleDomainLeavingParticles(unsigned dim, ParticleContainer* moleculeContainer) const;
301
312 void handleDomainLeavingParticlesDirect(const HaloRegion& haloRegion, ParticleContainer* moleculeContainer,
313 std::vector<Molecule>& invalidParticles) const;
314
320 void handleForceExchange(unsigned dim, ParticleContainer* moleculeContainer) const;
321
327 virtual void handleForceExchangeDirect(const HaloRegion& haloRegion, ParticleContainer* moleculeContainer) const;
328
329 void populateHaloLayerWithCopies(unsigned dim, ParticleContainer* moleculeContainer) const;
330
331 void populateHaloLayerWithCopiesDirect(const HaloRegion& haloRegion, ParticleContainer* moleculeContainer,
332 bool positionCheck = true) const;
333
335 int _rank;
336
339
340private:
341 CollectiveCommBase _collCommBase;
342 int _sendLeavingAndCopiesSeparately = 0;
343};
344
345#endif /* DOMAINDECOMPBASE_H_ */
This class is a dummy class which ensures that the collective communication commands also work if the...
Definition: CollectiveCommBase.h:21
virtual size_t getTotalSize() override
Definition: CollectiveCommBase.h:156
Class implementing molecules as rigid rotators consisting out of different interaction sites (LJcente...
Definition: Component.h:14
Definition: NeighbourCommunicationScheme.h:128
handle boundary region and multiple processes
Definition: DomainDecompBase.h:51
virtual void collCommAppendFloat(float floatValue)
has to call appendFloat method of a CollComm class
Definition: DomainDecompBase.cpp:566
void handleDomainLeavingParticlesDirect(const HaloRegion &haloRegion, ParticleContainer *moleculeContainer, std::vector< Molecule > &invalidParticles) const
Definition: DomainDecompBase.cpp:227
virtual double collCommGetDouble()
has to call getDouble method of a CollComm class
Definition: DomainDecompBase.cpp:590
virtual void collCommAppendInt(int intValue)
has to call appendInt method of a CollComm class
Definition: DomainDecompBase.cpp:558
int _numProcs
total number of processes in the simulation
Definition: DomainDecompBase.h:338
void writeMoleculesToMPIFileBinary(const std::string &filename, ParticleContainer *moleculeContainer) const
appends molecule data to the file. The format is the same as that of the input file This version uses...
Definition: DomainDecompBase.cpp:454
int _rank
the id of the current process
Definition: DomainDecompBase.h:335
virtual void balanceAndExchange(double lastTraversalTime, bool forceRebalancing, ParticleContainer *moleculeContainer, Domain *domain)
balance the load (and optimize communication) and exchange boundary particles
Definition: DomainDecompBase.cpp:405
virtual void collCommInit(int numValues, int key=0)
has to call init method of a CollComm class
Definition: DomainDecompBase.cpp:550
virtual bool queryBalanceAndExchangeNonBlocking(bool forceRebalancing, ParticleContainer *moleculeContainer, Domain *domain, double etime)
Checks whether the balance and exchange step can be performed non-blocking.
Definition: DomainDecompBase.cpp:401
double getIOCutoffRadius(int dim, Domain *domain, ParticleContainer *moleculeContainer)
returns an cutoff radius for a dimension for a global linked cells datastructure
Definition: DomainDecompBase.cpp:618
virtual unsigned long collCommGetUnsLong()
has to call getUnsLong method of a CollComm class
Definition: DomainDecompBase.cpp:582
virtual void collCommBroadcast(int root=0)
has to call broadcast method of a CollComm class (none in sequential version)
Definition: DomainDecompBase.cpp:614
virtual int collCommGetInt()
has to call getInt method of a CollComm class
Definition: DomainDecompBase.cpp:578
virtual void collCommAppendDouble(double doubleValue)
has to call appendDouble method of a CollComm class
Definition: DomainDecompBase.cpp:570
void getBoundingBoxMinMax(Domain *domain, double *min, double *max)
get the minimum and maximum coordinate of the bounding box of this process' domain
Definition: DomainDecompBase.cpp:543
virtual double getTime() const
returns the time in seconds since some time in the past
Definition: DomainDecompBase.cpp:423
virtual void collCommAllreduceSum()
has to call allreduceSum method of a CollComm class (none in sequential version)
Definition: DomainDecompBase.cpp:598
virtual long double collCommGetLongDouble()
has to call getLongDouble method of a CollComm class
Definition: DomainDecompBase.cpp:594
virtual double getBoundingBoxMin(int dimension, Domain *domain)
get the minimum of the bounding box of this process' domain in the given dimension (0,...
Definition: DomainDecompBase.cpp:415
virtual void assertIntIdentity(int IX)
checks identity of random number generators
Definition: DomainDecompBase.cpp:433
virtual void barrier() const
synchronizes all processes
Definition: DomainDecompBase.cpp:451
virtual int getNumProcs() const
returns the number of processes
Definition: DomainDecompBase.cpp:447
virtual int getRank() const
returns the own rank
Definition: DomainDecompBase.cpp:443
void handleDomainLeavingParticles(unsigned dim, ParticleContainer *moleculeContainer) const
Definition: DomainDecompBase.cpp:174
void handleForceExchange(unsigned dim, ParticleContainer *moleculeContainer) const
Does the force exchange for each dimension. Will be called for dim=0, 1 and 2.
Definition: DomainDecompBase.cpp:94
void writeMoleculesToFile(const std::string &filename, ParticleContainer *moleculeContainer, bool binary=false) const
appends molecule data to the file. The format is the same as that of the input file If MPI is enabled...
Definition: DomainDecompBase.cpp:505
virtual void printDecomp(const std::string &filename, Domain *domain, ParticleContainer *particleContainer)
writes information about the current decomposition into the given file
Definition: DomainDecompBase.cpp:439
virtual void collCommScanSum()
has to call scanSum method of a CollComm class (none in sequential version)
Definition: DomainDecompBase.cpp:610
void exchangeMolecules(ParticleContainer *moleculeContainer, Domain *domain)
exchange molecules between processes
Definition: DomainDecompBase.cpp:51
virtual float collCommGetFloat()
has to call getFloat method of a CollComm class
Definition: DomainDecompBase.cpp:586
virtual unsigned Ndistribution(unsigned localN, float *minrnd, float *maxrnd)
returns total number of molecules
Definition: DomainDecompBase.cpp:427
virtual void collCommFinalize()
has to call finalize method of a CollComm class
Definition: DomainDecompBase.cpp:554
virtual void collCommAppendLongDouble(long double longDoubleValue)
has to call appendLongDouble method of a CollComm class
Definition: DomainDecompBase.cpp:574
virtual void collCommAllreduceCustom(ReduceType type)
has to call allreduceCustom method of a CollComm class (none in sequential version)
Definition: DomainDecompBase.cpp:606
virtual void collCommAppendUnsLong(unsigned long unsLongValue)
has to call appendUnsLong method of a CollComm class
Definition: DomainDecompBase.cpp:562
virtual void exchangeForces(ParticleContainer *moleculeContainer, Domain *domain)
Exchanges forces at the domain boundaries if it's required by the cell container.
Definition: DomainDecompBase.cpp:88
DomainDecompBase()
The Constructor determines the own rank and the number of the neighbours *‍/.
Definition: DomainDecompBase.cpp:19
virtual int getNonBlockingStageCount()
Definition: DomainDecompBase.cpp:397
virtual bool procOwnsPos(double x, double y, double z, Domain *domain) final
find out whether the given position belongs to the domain of this process
Definition: DomainDecompBase.cpp:409
virtual double getBoundingBoxMax(int dimension, Domain *domain)
get the maximum of the bounding box of this process' domain in the given dimension (0,...
Definition: DomainDecompBase.cpp:419
virtual void handleForceExchangeDirect(const HaloRegion &haloRegion, ParticleContainer *moleculeContainer) const
Definition: DomainDecompBase.cpp:143
~DomainDecompBase() override
The Destructor finalizes MPI.
Definition: DomainDecompBase.cpp:22
virtual void collCommAllreduceSumAllowPrevious()
has to call allreduceSum method of a CollComm class (none in sequential version), allows for values o...
Definition: DomainDecompBase.cpp:602
This class is used to read in the phasespace and to handle macroscopic values.
Definition: Domain.h:47
Definition: NeighbourCommunicationScheme.h:172
Definition: MemoryProfiler.h:13
Definition: NeighbourCommunicationScheme.h:18
This Interface is used to get access to particles and pairs of particles.
Definition: ParticleContainer.h:69
XML file with unit attributes abstraction.
Definition: xmlfileUnits.h:25
Enumeration class corresponding to the type schema type.
Definition: vtk-unstructured.h:1746
::xsd::cxx::tree::string< char, simple_type > string
C++ type corresponding to the string XML Schema built-in type.
Definition: vtk-punstructured.h:270
Definition: HaloRegion.h:10