source: src/comm/domain_decomposition_mpi.cpp@ 8180d8

Last change on this file since 8180d8 was 8180d8, checked in by Julian Iseringhausen <julian.iseringhausen@…>, 13 years ago

Merge stashed open boundary stuff.

  • Property mode set to 100644
File size: 4.3 KB
Line 
1/*
2 * vmg - a versatile multigrid solver
3 * Copyright (C) 2012 Institute for Numerical Simulation, University of Bonn
4 *
5 * vmg is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * vmg is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19/**
20 * @file domain_decomposition_mpi.cpp
21 * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
22 * @date Mon Jun 27 12:53:50 2011
23 *
24 * @brief Computes a domain decomposition which separates
25 * the finest grid equally for all processes.
26 *
27 */
28
29#ifdef HAVE_CONFIG_H
30#include <config.h>
31#endif
32
33#include "base/interface.hpp"
34#include "comm/comm.hpp"
35#include "comm/domain_decomposition_mpi.hpp"
36#include "grid/grid.hpp"
37#include "grid/multigrid.hpp"
38
39using namespace VMG;
40
41void DomainDecompositionMPI::Compute(Comm* comm, const Interface* interface, std::vector<GlobalIndices>& global)
42{
43 GlobalIndices global_l;
44 Index remainder, procs;
45 Index last_procs = comm->GlobalProcs();
46
47 global.clear();
48
49 for (unsigned int i=0; i<interface->Global().size(); ++i) {
50
51 global_l = interface->Global()[i];
52
53 if (IsActive(comm, global_l.GlobalSizeNew(), procs)) {
54
55 if (i == 0) {
56
57 remainder = global_l.GlobalSizeNew() % procs;
58
59 global_l.LocalSize() = global_l.GlobalSizeNew() / procs;
60 for (int j=0; j<3; ++j)
61 if (comm->GlobalPos()[j] < remainder[j])
62 ++(global_l.LocalSize()[j]);
63
64 global_l.LocalBegin() = global_l.GlobalBegin() + comm->GlobalPos() * global_l.LocalSize();
65 for (int j=0; j<3; ++j)
66 if (comm->GlobalPos()[j] >= remainder[j])
67 global_l.LocalBegin()[j] += remainder[j];
68
69 global_l.LocalEnd() = global_l.LocalBegin() + global_l.LocalSize();
70
71 } else {
72
73 for (int j=0; j<3; ++j) {
74
75 if (procs[j] == last_procs[j]) {
76
77 if (global.back().LocalBegin()[j] == global.back().GlobalBegin()[j])
78 global_l.LocalBegin()[j] = global_l.GlobalBegin()[j];
79 else
80 global_l.LocalBegin()[j] = global.back().LocalBegin()[j] / 2;
81
82 if (global.back().LocalEnd()[j] == global.back().GlobalEnd()[j])
83 global_l.LocalEnd()[j] = global_l.GlobalEnd()[j];
84 else
85 global_l.LocalEnd()[j] = global.back().LocalEnd()[j] / 2;
86
87 global_l.LocalSize()[j] = global_l.LocalEnd()[j] - global_l.LocalBegin()[j];
88
89 } else {
90
91 remainder[j] = global_l.GlobalSizeNew()[j] % procs[j];
92
93 global_l.LocalSize()[j] = global_l.GlobalSizeNew()[j] / procs[j];
94 if (comm->GlobalPos()[j] < remainder[j])
95 ++(global_l.LocalSize()[j]);
96
97 global_l.LocalBegin()[j] = global_l.GlobalBegin()[j] + comm->GlobalPos()[j] * global_l.LocalSize()[j];
98 if (comm->GlobalPos()[j] >= remainder[j])
99 global_l.LocalBegin()[j] += remainder[j];
100
101 global_l.LocalEnd()[j] = global_l.LocalBegin()[j] + global_l.LocalSize()[j];
102
103 }
104 }
105 }
106 }else {
107 global_l.LocalBegin() = 0;
108 global_l.LocalEnd() = 0;
109 global_l.LocalSize() = 0;
110 }
111
112 last_procs = procs;
113
114 global.push_back(global_l);
115
116 }
117}
118
119bool DomainDecompositionMPI::IsActive(Comm* comm, const Index& size_global, Index& procs)
120{
121 bool is_active = true;
122 const int points_min = 5;
123
124 procs = size_global / points_min + 1;
125
126 for (int i=0; i<3; ++i) {
127 procs[i] = std::min(procs[i], comm->GlobalProcs()[i]);
128 is_active &= comm->GlobalPos()[i] < procs[i];
129 }
130
131 return is_active;
132}
133
134void DomainDecompositionMPI::FineToCoarse(Comm* comm, int& begin, int& end, int levels)
135{
136 int last_point = end - 1;
137
138 for (int i=0; i<levels; ++i) {
139
140 if (begin % 2 == 0)
141 begin /= 2;
142 else
143 begin = (begin+1) / 2;
144
145 if (last_point % 2 == 0)
146 last_point /= 2;
147 else
148 last_point = (last_point-1) / 2;
149
150 }
151
152 end = last_point + 1;
153}
Note: See TracBrowser for help on using the repository browser.