DUNE-FEM (unstable)

partitioning.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3#ifndef DUNE_GRID_YASPGRID_PARTITIONING_HH
4#define DUNE_GRID_YASPGRID_PARTITIONING_HH
5
13#include<array>
14
15#include<dune/common/math.hh>
16
17namespace Dune
18{
19
20 namespace Yasp
21 {
22
36 template<int d>
38 {
39 public:
40 using iTupel = std::array<int, d>;
41 virtual ~Partitioning() = default;
42 virtual void partition(const iTupel&, int, iTupel&, int) const = 0;
43 };
44
45 template<int d>
46 class DefaultPartitioning : public Partitioning<d>
47 {
48 public:
49 using iTupel = std::array<int, d>;
50
56 void partition (const iTupel& size, int P, iTupel& dims, int overlap) const final
57 {
58 double opt=1E100;
59 iTupel trydims;
60
61 trydims.fill(-1);
62 dims.fill(-1);
63
64 optimize_dims(d-1,size,P,dims,trydims,opt,overlap);
65 if (dims[0] == -1)
66 DUNE_THROW(Dune::GridError, "Failed to find a suitable partition");
67 }
68
69 private:
70 void optimize_dims (int i, const iTupel& size, int P, iTupel& dims, iTupel& trydims, double &opt, int overlap ) const
71 {
72 if (i>0) // test all subdivisions recursively
73 {
74 for (int k=1; k<=P; k++)
75 if (
76 P%k==0 // k divides P
77 and (
78 k == 1 // no neighbors
79 or
80 size[i] / k >= 2*overlap // size sufficient for overlap
81 )
82 )
83 {
84 // P divisible by k
85 trydims[i] = k;
86 optimize_dims(i-1,size,P/k,dims,trydims,opt,overlap);
87 }
88 }
89 else
90 {
91 // found a possible combination
92 if (
93 P == 1 // no neighbors
94 or
95 size[0] / P >= 2*overlap // size sufficient for overlap
96 )
97 trydims[0] = P;
98 else
99 return;
100
101 // check for optimality
102 double m = -1.0;
103
104 for (int k=0; k<d; k++)
105 {
106 double mm=((double)size[k])/((double)trydims[k]);
107 if (fmod((double)size[k],(double)trydims[k])>0.0001) mm*=3;
108 if ( mm > m ) m = mm;
109 }
110 //if (_rank==0) std::cout << "optimize_dims: " << size << " | " << trydims << " norm=" << m << std::endl;
111 if (m<opt)
112 {
113 opt = m;
114 dims = trydims;
115 }
116 }
117 }
118 };
119
122 template<int d>
124 {
125 public:
126 typedef std::array<int, d> iTupel;
127 virtual ~PowerDPartitioning() {}
128
129 void partition (const iTupel& size, int P, iTupel& dims, int overlap) const final
130 {
131 for (int i=1; i<=P; ++i)
132 if (Dune::power(i, d) == P) {
133 std::fill(dims.begin(), dims.end(),i);
134 return;
135 }
136
137 DUNE_THROW(GridError, "Power partitioning failed: your number of processes needs to be a " << d << "-th power.");
138 }
139 };
140
145 template<int d>
147 {
148 public:
149 FixedSizePartitioning(const std::array<int,d>& dims) : _dims(dims) {}
150
151 virtual ~FixedSizePartitioning() {}
152
153 void partition(const std::array<int,d>&, int P, std::array<int,d>& dims, int overlap) const final
154 {
155 int prod = 1;
156 for (int i=0; i<d; i++)
157 prod *= _dims[i];
158 if (P != prod)
159 DUNE_THROW(Dune::Exception,"Your processor number doesn't match your partitioning information");
160 dims = _dims;
161 }
162
163 private:
164 std::array<int,d> _dims;
165 };
166
168 }
169}
170
171#endif
Base class for Dune-Exceptions.
Definition: exceptions.hh:96
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:20
Implement partitioner that gets a fixed partitioning from an array If the given partitioning doesn't ...
Definition: partitioning.hh:147
a base class for the yaspgrid partitioning strategy
Definition: partitioning.hh:38
Implement yaspgrid load balance strategy for P=x^{dim} processors.
Definition: partitioning.hh:124
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)