154 #ifndef GETFEM_CONFIG_H__
155 #define GETFEM_CONFIG_H__
165 #if !defined(GETFEM_PARA_LEVEL)
166 # define GETFEM_PARA_LEVEL 0
169 #define GETFEM_MPI_INIT(argc, argv) {GMM_TRACE1("Running sequential Getfem");}
170 #define GETFEM_MPI_FINALIZE {}
172 #if defined(GMM_USES_MPI)
175 # undef GMM_TRACE_MSG_MPI
176 # define GMM_TRACE_MSG_MPI \
177 int mip_rk__; MPI_Comm_rank(MPI_COMM_WORLD, &mip_rk__); \
180 # undef GETFEM_MPI_INIT
181 # define GETFEM_MPI_INIT(argc, argv) { \
182 MPI_Init(&argc, &argv); \
183 GMM_TRACE1("Running parallelized Getfem level " << GETFEM_PARA_LEVEL); \
185 # undef GETFEM_MPI_FINALIZE
186 # define GETFEM_MPI_FINALIZE { MPI_Finalize(); }
198 using std::endl;
using std::cout;
using std::cerr;
199 using std::ends;
using std::cin;
202 #if GETFEM_PARA_LEVEL > 1
207 # define MUMPS_PARA_SOLVER 1
208 # define SCHWARZADD_PARA_SOLVER 2
210 # ifndef GETFEM_PARA_SOLVER
211 # define GETFEM_PARA_SOLVER MUMPS_PARA_SOLVER
214 template <
typename T>
inline T MPI_SUM_SCALAR(T a) {
215 T b; MPI_Allreduce(&a,&b,1,gmm::mpi_type(a), MPI_SUM, MPI_COMM_WORLD);
219 template <
typename VECT>
inline void MPI_SUM_VECTOR(
const VECT &VV) {
220 VECT &V =
const_cast<VECT &
>(VV);
221 typedef typename gmm::linalg_traits<VECT>::value_type T;
222 std::vector<T> W(gmm::vect_size(V));
223 MPI_Allreduce((
void *)(&(V[0])), &(W[0]), gmm::vect_size(V),
224 gmm::mpi_type(T()), MPI_SUM, MPI_COMM_WORLD);
228 template <
typename VECT>
inline void MPI_MAX_VECTOR(
const VECT &VV) {
229 VECT &V =
const_cast<VECT &
>(VV);
230 typedef typename gmm::linalg_traits<VECT>::value_type T;
231 std::vector<T> W(gmm::vect_size(V));
232 MPI_Allreduce((
void *)(&(V[0])), &(W[0]), gmm::vect_size(V),
233 gmm::mpi_type(T()), MPI_MAX, MPI_COMM_WORLD);
237 template <
typename T>
void MPI_BCAST0_SCALAR(T &a) {
238 MPI_Bcast((
void *)(&a), 1, gmm::mpi_type(a), 0, MPI_COMM_WORLD);
241 template <
typename VECT>
inline void MPI_BCAST0_VECTOR(
const VECT &VV) {
242 VECT &V =
const_cast<VECT &
>(VV);
243 typedef typename gmm::linalg_traits<VECT>::value_type T;
244 MPI_Bcast((
void *)(&(V[0])), gmm::vect_size(V), gmm::mpi_type(T()), 0,
248 template <
typename VECT1,
typename VECT2>
249 inline void MPI_SUM_VECTOR(
const VECT1 &VV,
const VECT2 &WW) {
250 VECT1 &V =
const_cast<VECT1 &
>(VV);
251 VECT2 &W =
const_cast<VECT2 &
>(WW);
252 typedef typename gmm::linalg_traits<VECT1>::value_type T;
253 MPI_Allreduce((
void *)(&(V[0])), &(W[0]),
254 gmm::vect_size(V), gmm::mpi_type(T()),
255 MPI_SUM, MPI_COMM_WORLD);
258 inline bool MPI_IS_MASTER(
void)
259 {
int rk; MPI_Comm_rank(MPI_COMM_WORLD, &rk);
return !rk; }
261 template <
typename T>
inline
263 typedef typename gmm::col_matrix< gmm::rsvector<T> > MAT;
264 MAT &MM =
const_cast<MAT &
>(M);
267 MPI_Comm_rank(MPI_COMM_WORLD, &rk);
268 MPI_Comm_size(MPI_COMM_WORLD, &nbp);
270 size_type nr = gmm::mat_nrows(MM), nc = gmm::mat_ncols(MM);
272 gmm::dense_matrix<int> all_nnz(nc, nbp);
273 std::vector<int> my_nnz(nc), owner(nc);
278 my_nnz[i] = gmm::nnz(gmm::mat_col(MM, i));
281 MPI_Allgather((
void *)(&(my_nnz[0])), nc, MPI_INT,
282 (
void *)(&(all_nnz(0,0))), nc, MPI_INT, MPI_COMM_WORLD);
285 std::vector<int> contributors(nc);
286 for (
int i = 0; i < nc; ++i) {
288 int column_master = -1, max_nnz = 0;
289 contributors.resize(0);
292 for (
int j = nbp-1; j >= 0; --j) {
294 if (rk != j) contributors.push_back(j);
295 if (column_master == -1 || all_nnz(i, j) > max_nnz)
296 { column_master = j; max_nnz = all_nnz(i, j); }
300 if (column_master == rk) {
301 for (
int j = 0; j < int(contributors.size()); ++j) {
302 typename gmm::rsvector<T>::base_type_ &VV = V;
303 int si = all_nnz(i, contributors[j]);
305 MPI_Recv((
void *)(&(VV[0])),
306 si*
sizeof(gmm::elt_rsvector_<T>),
307 MPI_BYTE, contributors[j], 0,
308 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
312 typename gmm::rsvector<T>::base_type_ &VV = MM.col(i);
313 MPI_Send((
void *)(&(VV[0])),
314 my_nnz[i]*
sizeof(gmm::elt_rsvector_<T>),
315 MPI_BYTE, column_master, 0,
324 my_nnz[i] =
gmm::nnz(gmm::mat_col(MM, i));
325 owner[i] = (my_nnz[i]) ? rk : 0;
327 MPI_SUM_VECTOR(my_nnz);
328 MPI_SUM_VECTOR(owner);
335 typename gmm::rsvector<T>::base_type_ &VV = MM.col(i);
336 if (owner[i] != rk) VV.resize(my_nnz[i]);
337 MPI_Bcast((
void *)(&(VV[0])), my_nnz[i]*
sizeof(gmm::elt_rsvector_<T>),
338 MPI_BYTE, owner[i], MPI_COMM_WORLD);
343 template <
typename MAT>
inline void MPI_SUM_SPARSE_MATRIX(
const MAT &M) {
344 typedef typename gmm::linalg_traits<MAT>::value_type T;
345 int nbp; MPI_Comm_size(MPI_COMM_WORLD, &nbp);
347 size_type nr = gmm::mat_nrows(M), nc = gmm::mat_ncols(M);
348 gmm::col_matrix< gmm::rsvector<T> > MM(nr, nc);
349 GMM_WARNING2(
"MPI_SUM_SPARSE_MATRIX: A copy of the matrix is done. "
350 "To avoid it, adapt MPI_SUM_SPARSE_MATRIX to "
351 "your matrix type.");
353 MPI_SUM_SPARSE_MATRIX(MM);
357 template <
typename T>
inline T MPI_SUM_SCALAR(T a) {
return a; }
358 template <
typename VECT>
inline void MPI_SUM_VECTOR(
const VECT &) {}
359 template <
typename VECT>
inline void MPI_MAX_VECTOR(
const VECT &) {}
360 template <
typename T>
void MPI_BCAST0_SCALAR(T &) {}
361 template <
typename VECT>
inline void MPI_BCAST0_VECTOR(
const VECT &) {}
362 template <
typename MAT>
inline void MPI_SUM_SPARSE_MATRIX(
const MAT &) {}
363 template <
typename VECT1,
typename VECT2>
364 inline void MPI_SUM_VECTOR(
const VECT1 &V,
const VECT2 &WW)
366 VECT2 &W =
const_cast<VECT2 &
>(WW);
369 inline bool MPI_IS_MASTER(
void) {
return true; }
374 using bgeot::dim_type;
377 using bgeot::scalar_type;
378 using bgeot::complex_type;
379 using bgeot::long_scalar_type;
380 using bgeot::opt_long_scalar_type;
383 using bgeot::base_vector;
384 using bgeot::base_complex_vector;
385 using bgeot::base_matrix;
386 using bgeot::base_complex_matrix;
387 using bgeot::base_tensor;
388 using bgeot::base_complex_tensor;
392 #if defined(__GNUC__)
395 inline bool isnan(scalar_type x) {
return x != x; }
defines and typedefs for namespace bgeot
Multivariate polynomials.
tensor class, used in mat_elem computations.
sparse vector built upon std::vector.
size_type nnz(const L &l)
count the number of non-zero entries of a vector or matrix.
void copy(const L1 &l1, L2 &l2)
*/
void add(const L1 &l1, L2 &l2)
*/
gmm::uint16_type short_type
used as the common short type integer in the library
size_t size_type
used as the common size type in the library
GEneric Tool for Finite Element Methods.