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