GetFEM  5.5
getfem_config.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2000-2026 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program. If not, see https://www.gnu.org/licenses/.
19 
20  As a special exception, you may use this file as it is a part of a free
21  software library without restriction. Specifically, if other files
22  instantiate templates or use macros or inline functions from this file,
23  or you compile this file and link it with other files to produce an
24  executable, this file does not by itself cause the resulting executable
25  to be covered by the GNU Lesser General Public License. This exception
26  does not however invalidate any other reasons why the executable file
27  might be covered by the GNU Lesser General Public License.
28 
29 ===========================================================================*/
30 
31 /*! @mainpage GetFEM reference documentation.
32 
33  <center><img src="https://getfem.org/_static/logo_getfem_small.png"></center>
34 
35  @section intro Introduction
36 
37  This documentation has been generated from the GetFEM sources, this is not a tutorial.
38 
39  @section Terminology
40 
41  This is just a short summary of the terms employed in GetFEM.
42 
43  The @link mesh mesh@endlink is composed of @link convexes
44  convexes@endlink. What we call convexes can be simple line segments,
45  prisms, tetrahedrons, curved triangles, of even something which is
46  not convex (in the geometrical sense). They all have an associated
47  @link refconvexes reference convex@endlink: for segments, this will
48  be the @f$[0,1]@f$ segment, for triangles this will be the canonical
49  triangle @f$(0,0)-(0,1)-(1,0)@f$ etc... All convexes of the mesh
50  are constructed from the reference convex through a @link
51  bgeot_geometric_trans.h geometric transformation@endlink. In
52  simple cases (when the convexes are simplices for example), this
53  transformation will be linear (hence it is easily inverted, which
54  can be a great advantage). In order to define the geometric
55  transformation, one defines <i>geometrical nodes</i> on the
56  reference convex. The geometrical transformation maps these nodes to
57  the <i>mesh nodes</i>.
58 
59  On the mesh, one defines a set a <i>basis functions</i>: the @link pfem FEM@endlink. A FEM
60  should be associated with each convex. The basis functions are also
61  attached to some geometrical points (which can be arbitrarily
62  chosen). These points are similar to the mesh nodes, but <b>they
63  don't have to be the same</b> (this only happens on very simple cases,
64  such as a classical P1 fem on a triangular mesh). The set of all
65  basis functions on the mesh forms the basis of a vector space, on
66  which the PDE will be solved. These basis functions (and their
67  associated geometrical point) are the <b>degrees of freedom
68  (dof)</b>. The FEM is said to be <i>Lagrangian</i> when each of its basis
69  functions is equal to one at its attached geometrical point, and is
70  zero at the geometrical points of others basis functions. This is an
71  important property as it is very easy to <i>interpolate</i> an arbitrary
72  function on the finite elements space.
73 
74  The finite elements methods involves evaluation of integrals of
75  these basis functions (or product of basis functions etc...) on
76  convexes (and faces of convexes). In simple cases (polynomial basis
77  functions and linear geometrical transformation), one can evaluate
78  analytically these integrals. In other cases, one has to approximate
79  it, using <i>quadrature formulas</i>. Hence, at each convex is attached
80  an @link getfem_integration.h integration method@endlink along with the FEM. If you have to use an
81  approximate integration method, always choose carefully its
82  order(i.e. highest degree of the polynomials who are exactly
83  integrated with the method) : the degree of the FEM, of the
84  polynomial degree of the geometrical transformation, and the nature
85  of the elementary matrix have to be taken into account. If you are
86  unsure about the appropriate degree, always prefer a high order
87  integration method (which will slow down the assembly) to a low
88  order one which will produce a useless linear-system.
89 
90  The process of construction of a global linear system from integrals
91  of basis functions on each convex is the @link asm assembly@endlink.
92 
93  A mesh, with a set of FEM attached to its convexes is called a @link getfem_mesh_fem.h
94  mesh_fem@endlink object in GetFEM.
95 
96  A mesh, with a set of integration methods attached to its convexes
97  is called a @link getfem_mesh_im.h mesh_im@endlink object in GetFEM.
98 
99  A @c mesh_fem can be used to approximate scalar fields (heat, pression,
100  ..), or vector fields (displacement, electric field, ..). A @c mesh_im
101  will be used to perform numerical integrations on these fields.
102  Although GetFEM supports vector FEMs, intrinsic vector FEMs are
103  essentially used in mixed methods (for instance Raviart-Thomas element).
104  Most of the FEM currently available are scalar. Of course,
105  these scalar FEMs can be used to approximate each component of a
106  vector field. This is done by setting the Qdim of the mesh_fem to
107  the dimension of the vector field (i.e. Qdim=1 for scalar field,
108  Qdim=2 for 2D vector field etc...).
109 
110  When solving a PDE, one often has to use more than one FEM. The most
111  important one will be of course the one on which is defined the
112  solution of the PDE. But most PDEs involve various coefficients, for
113  example: @f[ \nabla.(\lambda(x)\nabla u) = f(x). @f] Hence one has
114  to define an FEM for the main unknown @f$u@f$, but also for the data
115  @f$\lambda(x)@f$ and @f$f(x)@f$ if they are not constant. In order
116  to interpolate easily these coefficients in their finite element
117  space, one often choose a Lagrangian FEM.
118 
119  The convexes, mesh nodes, and dof are all numbered. We sometimes
120  refer to the number associated to a convex as its convex id or
121  convex number. Mesh node numbers are also called point id or point
122  number. Faces of convexes do not have a global numbering, but only
123  a local number in each convex.
124 
125  While the @c dof are always numbered consecutively, <b>this is not
126  always the case for point ids and convex ids</b>, especially if you
127  have removed points or convexes from the mesh. To ensure that they
128  form a continuous sequence (starting from 0), you have to use the
129  getfem::getfem_mesh::optimize_structure member function.
130 
131  Historically, GetFEM models where built by combining so called model
132  bricks, which are predefined PDE terms available in GetFEM.
133  Nowadays it is recommended to use the General Weak Form Language (GWFL)
134  instead, which is as easy, more flexible and more transparent.
135 
136  @section Examples
137 
138  @link thermo_elasticity_electrical_coupling.cc tests/thermo_elasticity_electrical_coupling.cc@endlink: solve a multiphysics problem with GWFL.
139  - @link laplacian.cc tests/laplacian.cc@endlink: solve the Laplace equation with the old bricks system.
140  - @link elastostatic.cc tests/elastostatic.cc@endlink: solve a static linear elasticity problem with the old bricks system.
141  - @link helmholtz.cc tests/helmholtz.cc@endlink: solve a scalar wave equation, with robin boundary conditions, with the old bricks system.
142  - @link stokes.cc tests/stokes.cc@endlink: solve the Stokes equation (incompressible viscuous fluid) with the old bricks system.
143  - @link nonlinear_elastostatic.cc tests/nonlinear_elastostatic.cc@endlink: solve a large strain elastostatic problem (torsion of a bar) with the old bricks system.
144  - @link icare.cc contrib/icare/icare.cc@endlink: solve the Navier-Stokes equation (fluid flow around an obstacle) with the old bricks system.
145 */
146 
147 /**@file getfem_config.h
148  @author Yves Renard <[email protected]>
149  @date November 19, 2000.
150  @brief defines and typedefs for namespace getfem
151 */
152 
153 
154 #ifndef GETFEM_CONFIG_H__
155 #define GETFEM_CONFIG_H__
156 
157 #include "bgeot_config.h"
158 
159 // Parallelization options
160 
161 // GETFEM_PARA_LEVEL is the MPI parallelisation level of GetFEM
162 // 0 - Sequential
163 // 1 - Only the assembly procedures are parallelized
164 // 2 - Assembly procedures and linear solves (with MUMPS) are parallelized
165 #if !defined(GETFEM_PARA_LEVEL)
166 # define GETFEM_PARA_LEVEL 0
167 #endif
168 
169 #define GETFEM_MPI_INIT(argc, argv) {GMM_TRACE1("Running sequential Getfem");}
170 #define GETFEM_MPI_FINALIZE {}
171 
172 #if defined(GMM_USES_MPI)
173 # include <mpi.h>
174 
175 # undef GMM_TRACE_MSG_MPI
176 # define GMM_TRACE_MSG_MPI \
177  int mip_rk__; MPI_Comm_rank(MPI_COMM_WORLD, &mip_rk__); \
178  if (mip_rk__ == 0)
179 
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); \
184  }
185 # undef GETFEM_MPI_FINALIZE
186 # define GETFEM_MPI_FINALIZE { MPI_Finalize(); }
187 
188 #endif
189 
190 
191 #include "bgeot_tensor.h"
192 #include "bgeot_poly.h"
193 
194 
195 /// GEneric Tool for Finite Element Methods.
196 namespace getfem {
197 
198  using std::endl; using std::cout; using std::cerr;
199  using std::ends; using std::cin;
200  using gmm::vref;
201 
202 #if GETFEM_PARA_LEVEL > 1
203 
204 // GETFEM_PARA_SOLVER is the parallelisation solver used
205 // MUMPS : use direct parallel solver MUMPS
206 // SCHWARZADD : use a Schwarz additive method
207 # define MUMPS_PARA_SOLVER 1
208 # define SCHWARZADD_PARA_SOLVER 2
209 
210 # ifndef GETFEM_PARA_SOLVER
211 # define GETFEM_PARA_SOLVER MUMPS_PARA_SOLVER
212 # endif
213 
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);
216  return b;
217  }
218 
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);
225  gmm::copy(W, V);
226  }
227 
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);
234  gmm::copy(W, V);
235  }
236 
237  template <typename T> void MPI_BCAST0_SCALAR(T &a) {
238  MPI_Bcast((void *)(&a), 1, gmm::mpi_type(a), 0, MPI_COMM_WORLD);
239  }
240 
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,
245  MPI_COMM_WORLD);
246  }
247 
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);
256  }
257 
258  inline bool MPI_IS_MASTER(void)
259  { int rk; MPI_Comm_rank(MPI_COMM_WORLD, &rk); return !rk; }
260 
261  template <typename T> inline
262  void MPI_SUM_SPARSE_MATRIX(const gmm::col_matrix< gmm::rsvector<T> > &M) {
263  typedef typename gmm::col_matrix< gmm::rsvector<T> > MAT;
264  MAT &MM = const_cast<MAT &>(M);
265  int rk, nbp;
266  MPI_Status status;
267  MPI_Comm_rank(MPI_COMM_WORLD, &rk);
268  MPI_Comm_size(MPI_COMM_WORLD, &nbp);
269  if (nbp < 2) return;
270  size_type nr = gmm::mat_nrows(MM), nc = gmm::mat_ncols(MM);
271 
272  gmm::dense_matrix<int> all_nnz(nc, nbp);
273  std::vector<int> my_nnz(nc), owner(nc);
274  gmm::rsvector<T> V(nr);
275 
276  // Step 1 : each process fill the corresponding nnz line
277  for (size_type i = 0; i < nc; ++i)
278  my_nnz[i] = gmm::nnz(gmm::mat_col(MM, i));
279 
280  // Step 2 : All gather : each proc has all the information
281  MPI_Allgather((void *)(&(my_nnz[0])), nc, MPI_INT,
282  (void *)(&(all_nnz(0,0))), nc, MPI_INT, MPI_COMM_WORLD);
283 
284  // Step 3 : Scan each column and message send/receive
285  std::vector<int> contributors(nc);
286  for (int i = 0; i < nc; ++i) {
287  if (my_nnz[i]) {
288  int column_master = -1, max_nnz = 0;
289  contributors.resize(0);
290 
291  // Determine who is the master for the column
292  for (int j = nbp-1; j >= 0; --j) {
293  if (all_nnz(i, 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); }
297  }
298  }
299 
300  if (column_master == rk) { // receive a column and store
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]);
304  VV.resize(si);
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);
309  gmm::add(V, gmm::mat_col(MM, i));
310  }
311  } else { // send the column to the column master
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,
316  MPI_COMM_WORLD);
317  MM.col(i).clear();
318  }
319  }
320  }
321 
322  // Step 3 : gather the total nnz
323  for (size_type i = 0; i < nc; ++i) {
324  my_nnz[i] = gmm::nnz(gmm::mat_col(MM, i));
325  owner[i] = (my_nnz[i]) ? rk : 0;
326  }
327  MPI_SUM_VECTOR(my_nnz);
328  MPI_SUM_VECTOR(owner);
329 
330  // Step 4 : distribute each column to all the processes
331  // Would it be more efficient to perform a global MPI_SUM on a compressed
332  // storage ?
333  for (size_type i = 0; i < nc; ++i) {
334  if (my_nnz[i]) {
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);
339  }
340  }
341  }
342 
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);
346  if (nbp < 2) return;
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.");
352  gmm::copy(M, MM);
353  MPI_SUM_SPARSE_MATRIX(MM);
354  gmm::copy(MM, const_cast<MAT &>(M));
355  }
356 #else
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)
365  {
366  VECT2 &W = const_cast<VECT2 &>(WW);
367  gmm::copy(V, W);
368  }
369  inline bool MPI_IS_MASTER(void) { return true; }
370 #endif
371 
372  using bgeot::ST_NIL;
373  using bgeot::size_type;
374  using bgeot::dim_type;
375  using bgeot::short_type;
376  using bgeot::short_type;
377  using bgeot::scalar_type;
378  using bgeot::complex_type;
379  using bgeot::long_scalar_type;
380  using bgeot::opt_long_scalar_type;
381 
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;
389  using bgeot::base_poly;
390  using bgeot::base_node;
391 
392 #if defined(__GNUC__)
393  using std::isnan;
394 #else
395  inline bool isnan(scalar_type x) { return x != x; }
396 #endif
397 
398 } /* end of namespace getfem. */
399 
400 
401 #endif /* GETFEM_CONFIG_H__ */
defines and typedefs for namespace bgeot
Multivariate polynomials.
tensor class, used in mat_elem computations.
sparse vector built upon std::vector.
Definition: gmm_vector.h:995
size_type nnz(const L &l)
count the number of non-zero entries of a vector or matrix.
Definition: gmm_blas.h:68
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:976
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
gmm::uint16_type short_type
used as the common short type integer in the library
Definition: bgeot_config.h:72
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
GEneric Tool for Finite Element Methods.