LatticeYangMills
parallel.h
Go to the documentation of this file.
1 /******************************************************************************
2 *
3 * MIT License
4 *
5 * Copyright (c) 2018 Giovanni Pederiva
6 *
7 * Permission is hereby granted, free of charge, to any person obtaining a copy
8 * of this software and associated documentation files (the "Software"), to deal
9 * in the Software without restriction, including without limitation the rights
10 * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
11 * copies of the Software, and to permit persons to whom the Software is
12 * furnished to do so, subject to the following conditions:
13 *
14 * The above copyright notice and this permission notice shall be included in
15 * all copies or substantial portions of the Software.
16 *
17 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
18 * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
19 * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
20 * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
21 * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
22 * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
23 * SOFTWARE.
24 ******************************************************************************/
25 
34 #pragma once
35 #include "Utils/clusterspecifier.h"
36 #ifndef LACONIA
37 #include <mpi/mpi.h>
38 #else
39 #include <mpi.h>
40 #endif
41 #include <array>
42 
56 class Parallel{
57 // Public interface of the Parallel class
58 public:
59  // Init-Final methods
60  static void initialize();
61  static void createGeometry(std::array<int, 4> latticeSize,
62  std::array<int, 4> subLatticeSize);
63  static void finalize();
64 
65  // Data members
66  static int rank() { return m_rank; }
67  static int numProcs() { return m_numProcs; }
68  static int activeProcs() { return m_activeProcs; }
69  static bool isActive() { return m_isActive; }
70  static MPI_Comm cartCoordComm() { return m_cartCoordComm; }
71  static std::array<int, 4>& subBlocks() { return m_subBlocks; }
72  static std::array<int, 4>& rankCoord() { return m_rankCoord; }
73  static std::array<int, 4>& latticeSubSize() { return m_latticeSubSize; }
74  static std::array<int, 4>& latticeFullSize() { return m_latticeFullSize; }
75  static std::array<int, 4>& parity(){ return m_parity; }
76 
77  // Neighbor getters
78  static int getNeighbor(int direction, int sign);
79  static int getSecondNeighbor(int direction1, int sign1,
80  int direction2, int sign2);
81 
82  // File IO
83  static void openFile (MPI_File& file, const char *fileName);
84  static void closeFile(MPI_File& file);
85 
86 private:
87  // Members made private for safety
88  static int m_rank;
89  static int m_numProcs;
90  static int m_activeProcs;
91  static bool m_isActive;
92  static MPI_Comm m_cartCoordComm;
93  static std::array<int, 4> m_subBlocks;
94  static std::array<int, 4> m_rankCoord;
95  static std::array<int, 4> m_latticeSubSize;
96  static std::array<int, 4> m_latticeFullSize;
97  static std::array<int, 4> m_parity;
98  static MPI_Comm m_cartCoords;
99  static std::array< std::array<int, 2>, 4> m_neighbor;
100  static std::array< std::array<std::array<std::array<int, 2>, 4>, 2>, 4> m_secondNeighbor;
101 
102  // Init-Only methods
103  static void createNeighborLists();
104  static void assignNeighbor(int direction);
105  static void assignSecondNeighbor(int dir1, int dir2);
106  static void MyMPI_Cart_shift2(MPI_Comm comm, int dir1, int shift1,
107  int dir2, int shift2,
108  int &source, int &dest);
109 };
static void createNeighborLists()
Create castesian coordinates and creates neighbor lists for every processor in every direction...
Definition: parallel.cpp:100
static void assignSecondNeighbor(int dir1, int dir2)
use MPI utilities to find the neighbors in two direction
Definition: parallel.cpp:172
static void initialize()
initializes MPI, gets the rank number and the total numer of processors
Definition: parallel.cpp:61
static void createGeometry(std::array< int, 4 > latticeSize, std::array< int, 4 > subLatticeSize)
creates the parallel gemometry of the lattice
Definition: parallel.cpp:77
Class to manage the parallelization scheme.
Definition: parallel.h:56
static void assignNeighbor(int direction)
use MPI utilities to find the neighbors in one direction
Definition: parallel.cpp:163
static void MyMPI_Cart_shift2(MPI_Comm comm, int dir1, int shift1, int dir2, int shift2, int &source, int &dest)
use MPI utilities to find the neighbors in two direction
Definition: parallel.cpp:186
static void closeFile(MPI_File &file)
closes a file with MPI
Definition: parallel.cpp:156
static int getSecondNeighbor(int direction1, int sign1, int direction2, int sign2)
return the rank of the neighbor along the two given directions and signs
Definition: parallel.cpp:210
static int getNeighbor(int direction, int sign)
return the rank of the neighbor along the given direction and sign
Definition: parallel.cpp:203
static void finalize()
finalizes MPI
Definition: parallel.cpp:92
static void openFile(MPI_File &file, const char *fileName)
opens a file with MPI
Definition: parallel.cpp:148