GetFEM  5.5
laplacian.cc
Go to the documentation of this file.
1 /*===========================================================================
2 
3  Copyright (C) 2002-2026 Yves Renard, Julien Pommier.
4 
5  This file is a part of GetFEM
6 
7  GetFEM is free software; you can redistribute it and/or modify it
8  under the terms of the GNU Lesser General Public License as published
9  by the Free Software Foundation; either version 3 of the License, or
10  (at your option) any later version along with the GCC Runtime Library
11  Exception either version 3.1 or (at your option) any later version.
12  This program is distributed in the hope that it will be useful, but
13  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
14  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15  License and GCC Runtime Library Exception for more details.
16  You should have received a copy of the GNU Lesser General Public License
17  along with this program. If not, see https://www.gnu.org/licenses/.
18 
19 ===========================================================================*/
20 
21 /**@file laplacian.cc
22  @brief Laplacian (Poisson) problem.
23 
24  The laplace equation is solved on a regular mesh of the unit
25  square, and is compared to the analytical solution.
26 
27  This program is used to check that getfem++ is working. This is
28  also a good example of use of GetFEM. This program does not use the
29  model bricks intentionally in order to serve as an example of solving
30  a pde directly with the assembly procedures.
31 */
32 
34 #include "getfem/getfem_export.h"
37 #include "gmm/gmm.h"
38 using std::endl; using std::cout; using std::cerr;
39 using std::ends; using std::cin;
40 
41 /* some GetFEM types that we will be using */
42 using bgeot::base_small_vector; /* special class for small (dim<16) vectors */
43 using bgeot::base_node; /* geometrical nodes (derived from base_small_vector)*/
44 using bgeot::scalar_type; /* = double */
45 using bgeot::size_type; /* = unsigned long */
46 
47 /* definition of some matrix/vector types. These ones are built
48  * using the predefined types in Gmm++
49  */
51 typedef gmm::row_matrix<sparse_vector_type> sparse_matrix_type;
52 typedef gmm::col_matrix<sparse_vector_type> col_sparse_matrix_type;
53 typedef std::vector<scalar_type> plain_vector;
54 
55 /* Definitions for the exact solution of the Laplacian problem,
56  * i.e. Delta(sol_u) + sol_f = 0
57  */
58 
59 base_small_vector sol_K; /* a coefficient */
60 /* exact solution */
61 scalar_type sol_u(const base_node &x) { return sin(gmm::vect_sp(sol_K, x)); }
62 /* righ hand side */
63 scalar_type sol_f(const base_node &x)
64 { return gmm::vect_sp(sol_K, sol_K) * sin(gmm::vect_sp(sol_K, x)); }
65 /* gradient of the exact solution */
66 base_small_vector sol_grad(const base_node &x)
67 { return sol_K * cos(gmm::vect_sp(sol_K, x)); }
68 
69 /*
70  structure for the Laplacian problem
71 */
72 struct laplacian_problem {
73 
74  enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
75  getfem::mesh mesh; /* the mesh */
76  getfem::mesh_im mim; /* the integration methods. */
77  getfem::mesh_fem mf_u; /* the main mesh_fem, for the Laplacian solution */
78  getfem::mesh_fem mf_rhs; /* the mesh_fem for the right hand side(f(x),..) */
79  getfem::mesh_fem mf_coef; /* the mesh_fem to represent pde coefficients */
80 
81  scalar_type residual; /* max residual for the iterative solvers */
82  size_type N;
83  bool gen_dirichlet;
84 
85  sparse_matrix_type SM; /* stiffness matrix. */
86  std::vector<scalar_type> U, B; /* main unknown, and right hand side */
87 
88  std::vector<scalar_type> Ud; /* reduced sol. for gen. Dirichlet condition. */
89  col_sparse_matrix_type NS; /* Dirichlet NullSpace
90  * (used if gen_dirichlet is true)
91  */
92  std::string datafilename;
93  bgeot::md_param PARAM;
94 
95  void assembly(void);
96  bool solve(void);
97  void init(void);
98  void compute_error();
99  laplacian_problem(void) : mim(mesh), mf_u(mesh), mf_rhs(mesh),
100  mf_coef(mesh) {}
101 };
102 
103 /* Read parameters from the .param file, build the mesh, set finite element
104  * and integration methods and selects the boundaries.
105  */
106 void laplacian_problem::init(void) {
107 
108  std::string MESH_TYPE = PARAM.string_value("MESH_TYPE","Mesh type ");
109  std::string FEM_TYPE = PARAM.string_value("FEM_TYPE","FEM name");
110  std::string INTEGRATION = PARAM.string_value("INTEGRATION",
111  "Name of integration method");
112 
113  cout << "MESH_TYPE=" << MESH_TYPE << "\n";
114  cout << "FEM_TYPE=" << FEM_TYPE << "\n";
115  cout << "INTEGRATION=" << INTEGRATION << "\n";
116 
117  /* First step : build the mesh */
120  N = pgt->dim();
121  std::vector<size_type> nsubdiv(N);
122  std::fill(nsubdiv.begin(),nsubdiv.end(),
123  PARAM.int_value("NX", "Nomber of space steps "));
124  getfem::regular_unit_mesh(mesh, nsubdiv, pgt,
125  PARAM.int_value("MESH_NOISED") != 0);
126 
127  bgeot::base_matrix M(N,N);
128  for (size_type i=0; i < N; ++i) {
129  static const char *t[] = {"LX","LY","LZ"};
130  M(i,i) = (i<3) ? PARAM.real_value(t[i],t[i]) : 1.0;
131  }
132  if (N>1) { M(0,1) = PARAM.real_value("INCLINE") * PARAM.real_value("LY"); }
133 
134  /* scale the unit mesh to [LX,LY,..] and incline it */
135  mesh.transformation(M);
136 
137  datafilename = PARAM.string_value("ROOTFILENAME","Base name of data files.");
138  scalar_type FT = PARAM.real_value("FT", "parameter for exact solution");
139  residual = PARAM.real_value("RESIDUAL");
140  if (residual == 0.) residual = 1e-10;
141  sol_K.resize(N);
142  for (size_type j = 0; j < N; j++)
143  sol_K[j] = ((j & 1) == 0) ? FT : -FT;
144 
145  /* set the finite element on the mf_u */
146  getfem::pfem pf_u = getfem::fem_descriptor(FEM_TYPE);
147  getfem::pintegration_method ppi = getfem::int_method_descriptor(INTEGRATION);
148 
149  mim.set_integration_method(mesh.convex_index(), ppi);
150  mf_u.set_finite_element(mesh.convex_index(), pf_u);
151 
152  /* set the finite element on mf_rhs (same as mf_u is DATA_FEM_TYPE is
153  not used in the .param file */
154  std::string data_fem_name = PARAM.string_value("DATA_FEM_TYPE");
155  if (data_fem_name.size() == 0) {
156  GMM_ASSERT1(pf_u->is_lagrange(), "You are using a non-lagrange FEM. "
157  << "In that case you need to set "
158  << "DATA_FEM_TYPE in the .param file");
159  mf_rhs.set_finite_element(mesh.convex_index(), pf_u);
160  } else {
161  mf_rhs.set_finite_element(mesh.convex_index(),
162  getfem::fem_descriptor(data_fem_name));
163  }
164 
165  /* set the finite element on mf_coef. Here we use a very simple element
166  * since the only function that need to be interpolated on the mesh_fem
167  * is f(x)=1 ... */
168  mf_coef.set_finite_element(mesh.convex_index(),
169  getfem::classical_fem(pgt,0));
170 
171  /* set boundary conditions
172  * (Neuman on the upper face, Dirichlet elsewhere) */
173  gen_dirichlet = PARAM.int_value("GENERIC_DIRICHLET");
174 
175  if (!pf_u->is_lagrange() && !gen_dirichlet)
176  GMM_WARNING2("With non lagrange fem prefer the generic "
177  "Dirichlet condition option");
178 
179  cout << "Selecting Neumann and Dirichlet boundaries\n";
180  getfem::mesh_region border_faces;
181  getfem::outer_faces_of_mesh(mesh, border_faces);
182  for (getfem::mr_visitor i(border_faces); !i.finished(); ++i) {
183  assert(i.is_face());
184  base_node un = mesh.normal_of_face_of_convex(i.cv(), i.f());
185  un /= gmm::vect_norm2(un);
186  if (gmm::abs(un[N-1] - 1.0) < 1.0E-7) { // new Neumann face
187  mesh.region(NEUMANN_BOUNDARY_NUM).add(i.cv(), i.f());
188  } else {
189  mesh.region(DIRICHLET_BOUNDARY_NUM).add(i.cv(), i.f());
190  }
191  }
192 }
193 
194 void laplacian_problem::assembly(void) {
195  size_type nb_dof = mf_u.nb_dof();
196  size_type nb_dof_rhs = mf_rhs.nb_dof();
197 
198  gmm::resize(B, nb_dof); gmm::clear(B);
199  gmm::resize(U, nb_dof); gmm::clear(U);
200  gmm::resize(SM, nb_dof, nb_dof); gmm::clear(SM);
201 
202  cout << "Number of dof : " << nb_dof << endl;
203  cout << "Assembling stiffness matrix" << endl;
204  getfem::asm_stiffness_matrix_for_laplacian(SM, mim, mf_u, mf_coef,
205  std::vector<scalar_type>(mf_coef.nb_dof(), 1.0));
206 
207  cout << "Assembling source term" << endl;
208  std::vector<scalar_type> F(nb_dof_rhs);
209  getfem::interpolation_function(mf_rhs, F, sol_f);
210  getfem::asm_source_term(B, mim, mf_u, mf_rhs, F);
211 
212  cout << "Assembling Neumann condition" << endl;
213  gmm::resize(F, nb_dof_rhs*N);
214  getfem::interpolation_function(mf_rhs, F, sol_grad);
215  getfem::asm_normal_source_term(B, mim, mf_u, mf_rhs, F,
216  NEUMANN_BOUNDARY_NUM);
217 
218  cout << "take Dirichlet condition into account" << endl;
219  if (!gen_dirichlet) {
220  std::vector<scalar_type> D(nb_dof);
221  getfem::interpolation_function(mf_u, D, sol_u);
222  getfem::assembling_Dirichlet_condition(SM, B, mf_u,
223  DIRICHLET_BOUNDARY_NUM, D);
224  } else {
225  gmm::resize(F, nb_dof_rhs);
226  getfem::interpolation_function(mf_rhs, F, sol_u);
227 
228  gmm::resize(Ud, nb_dof);
229  gmm::resize(NS, nb_dof, nb_dof);
230  col_sparse_matrix_type H(nb_dof_rhs, nb_dof);
231  std::vector<scalar_type> R(nb_dof_rhs);
232  std::vector<scalar_type> RHaux(nb_dof);
233 
234  /* build H and R such that U mush satisfy H*U = R */
235  getfem::asm_dirichlet_constraints(H, R, mim, mf_u, mf_rhs,
236  mf_rhs, F, DIRICHLET_BOUNDARY_NUM);
237 
238  gmm::clean(H, 1e-12);
239 // cout << "H = " << H << endl;
240 // cout << "R = " << R << endl;
241  int nbcols = int(getfem::Dirichlet_nullspace(H, NS, R, Ud));
242  // cout << "Number of irreductible unknowns : " << nbcols << endl;
243  gmm::resize(NS, gmm::mat_ncols(H),nbcols);
244 
245  gmm::mult(SM, Ud, gmm::scaled(B, -1.0), RHaux);
246  gmm::resize(B, nbcols);
247  gmm::resize(U, nbcols);
248  gmm::mult(gmm::transposed(NS), gmm::scaled(RHaux, -1.0), B);
249  sparse_matrix_type SMaux(nbcols, nb_dof);
250  gmm::mult(gmm::transposed(NS), SM, SMaux);
251  gmm::resize(SM, nbcols, nbcols);
252  /* NSaux = NS, but is stored by rows instead of by columns */
253  sparse_matrix_type NSaux(nb_dof, nbcols); gmm::copy(NS, NSaux);
254  gmm::mult(SMaux, NSaux, SM);
255  }
256 }
257 
258 
259 bool laplacian_problem::solve(void) {
260 
261  // see_schmidt(SM, U, B);
262 
263  cout << "Compute preconditionner\n";
264  gmm::iteration iter(residual, 1, 40000);
265  double time = gmm::uclock_sec();
266  if (1) {
267  // gmm::identity_matrix P;
268  // gmm::diagonal_precond<sparse_matrix_type> P(SM);
269  // gmm::mr_approx_inverse_precond<sparse_matrix_type> P(SM, 10, 10E-17);
270  // gmm::ildlt_precond<sparse_matrix_type> P(SM);
271  // gmm::ildltt_precond<sparse_matrix_type> P(SM, 20, 1E-6);
273  // gmm::ilutp_precond<sparse_matrix_type> P(SM, 20, 1E-6);
274  // gmm::ilu_precond<sparse_matrix_type> P(SM);
275  cout << "Time to compute preconditionner : "
276  << gmm::uclock_sec() - time << " seconds\n";
277 
278 
279  //gmm::HarwellBoeing_IO::write("SM", SM);
280 
281  // gmm::cg(SM, U, B, P, iter);
282  gmm::gmres(SM, U, B, P, 50, iter);
283  } else {
284  double rcond;
285 #if defined(GMM_USES_SUPERLU)
286  gmm::SuperLU_solve(SM, U, B, rcond);
287 #endif
288  cout << "cond = " << 1/rcond << "\n";
289  }
290 
291  cout << "Total time to solve : "
292  << gmm::uclock_sec() - time << " seconds\n";
293 
294  if (gen_dirichlet) {
295  std::vector<scalar_type> Uaux(mf_u.nb_dof());
296  gmm::mult(NS, U, Ud, Uaux);
297  gmm::resize(U, mf_u.nb_dof());
298  gmm::copy(Uaux, U);
299  }
300 
301  return (iter.converged());
302 }
303 
304 /* compute the error with respect to the exact solution */
305 void laplacian_problem::compute_error() {
306  std::vector<scalar_type> V(mf_rhs.nb_basic_dof());
307  getfem::interpolation(mf_u, mf_rhs, U, V);
308  for (size_type i = 0; i < mf_rhs.nb_basic_dof(); ++i)
309  V[i] -= sol_u(mf_rhs.point_of_basic_dof(i));
310  cout.precision(16);
311  cout << "L2 error = " << getfem::asm_L2_norm(mim, mf_rhs, V) << endl
312  << "H1 error = " << getfem::asm_H1_norm(mim, mf_rhs, V) << endl
313  << "Linfty error = " << gmm::vect_norminf(V) << endl;
314 }
315 
316 /**************************************************************************/
317 /* main program. */
318 /**************************************************************************/
319 
320 int main(int argc, char *argv[]) {
321 
322  GETFEM_MPI_INIT(argc, argv);
323  GMM_SET_EXCEPTION_DEBUG; // Exceptions make a memory fault, to debug.
324  FE_ENABLE_EXCEPT; // Enable floating point exception for Nan.
325 
326  try {
327  laplacian_problem p;
328  p.PARAM.read_command_line(argc, argv);
329  p.init();
330  p.mesh.write_to_file(p.datafilename + ".mesh");
331  p.assembly();
332  if (!p.solve()) GMM_ASSERT1(false, "Solve procedure has failed");
333  p.compute_error();
334  }
335  GMM_STANDARD_CATCH_ERROR;
336 
337  GETFEM_MPI_FINALIZE;
338 
339  return 0;
340 }
Describe a finite element method linked to a mesh.
Describe an integration method linked to a mesh.
"iterator" class for regions.
structure used to hold a set of convexes and/or convex faces.
Describe a mesh (collection of convexes (elements) and points).
Definition: getfem_mesh.h:98
Incomplete LU with threshold and K fill-in Preconditioner.
The Iteration object calculates whether the solution has reached the desired accuracy,...
Definition: gmm_iter.h:52
sparse vector built upon std::vector.
Definition: gmm_vector.h:995
Miscelleanous assembly routines for common terms. Use the low-level generic assembly....
Compute the gradient of a field on a getfem::mesh_fem.
Export solutions to various formats.
Build regular meshes.
Include common gmm files.
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:976
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
Definition: gmm_blas.h:556
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norminf(const V &v)
Infinity norm of a vector.
Definition: gmm_blas.h:692
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:58
void resize(V &v, size_type n)
*‍/
Definition: gmm_blas.h:209
void clean(L &l, double threshold)
Clean a vector or matrix (replace near-zero entries with zeroes).
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*‍/
Definition: gmm_blas.h:1663
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*‍/
Definition: gmm_blas.h:263
void gmres(const Mat &A, Vec &x, const VecB &b, const Precond &M, int restart, iteration &outer, Basis &KS)
Generalized Minimum Residual.
void asm_stiffness_matrix_for_laplacian(MAT &M, const mesh_im &mim, const mesh_fem &mf, const mesh_fem &mf_data, const VECT &A, const mesh_region &rg=mesh_region::all_convexes())
assembly of , where is scalar.
size_type Dirichlet_nullspace(const MAT1 &H, MAT2 &NS, const VECT1 &R, VECT2 &U0)
Build an orthogonal basis of the kernel of H in NS, gives the solution of minimal norm of H*U = R in ...
void asm_source_term(const VECT1 &B, const mesh_im &mim, const mesh_fem &mf, const mesh_fem &mf_data, const VECT2 &F, const mesh_region &rg=mesh_region::all_convexes())
source term (for both volumic sources and boundary (Neumann) sources).
scalar_type asm_L2_norm(const mesh_im &mim, const mesh_fem &mf, const VEC &U, const mesh_region &rg=mesh_region::all_convexes())
compute , U might be real or complex
scalar_type asm_H1_norm(const mesh_im &mim, const mesh_fem &mf, const VEC &U, const mesh_region &rg=mesh_region::all_convexes())
compute the H1 norm of U.
void asm_dirichlet_constraints(MAT &H, VECT1 &R, const mesh_im &mim, const mesh_fem &mf_u, const mesh_fem &mf_mult, const mesh_fem &mf_r, const VECT2 &r_data, const mesh_region &region, int version=ASMDIR_BUILDALL)
Assembly of Dirichlet constraints in a weak form.
void asm_normal_source_term(VECT1 &B, const mesh_im &mim, const mesh_fem &mf, const mesh_fem &mf_data, const VECT2 &F, const mesh_region &rg)
Normal source term (for boundary (Neumann) condition).
void APIDECL outer_faces_of_mesh(const mesh &m, const dal::bit_vector &cvlst, convex_face_ct &flist)
returns a list of "exterior" faces of a mesh (i.e.
Definition: getfem_mesh.cc:821
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
Definition: getfem_fem.h:243
pfem fem_descriptor(const std::string &name)
get a fem descriptor from its string name.
Definition: getfem_fem.cc:4659
pfem classical_fem(bgeot::pgeometric_trans pgt, short_type k, bool complete=false)
Give a pointer on the structures describing the classical polynomial fem of degree k on a given conve...
Definition: getfem_fem.cc:4566
pgeometric_trans geometric_trans_descriptor(std::string name)
Get the geometric transformation from its string name.
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
void interpolation_function(mesh_fem &mf_target, const VECT &VV, F &f, mesh_region rg=mesh_region::all_convexes())
interpolation of a function f on mf_target.
void regular_unit_mesh(mesh &m, std::vector< size_type > nsubdiv, bgeot::pgeometric_trans pgt, bool noised=false)
Build a regular mesh of the unit square/cube/, etc.
void interpolation(const mesh_fem &mf_source, const mesh_fem &mf_target, const VECTU &U, VECTV &V, int extrapolation=0, double EPS=1E-10, mesh_region rg_source=mesh_region::all_convexes(), mesh_region rg_target=mesh_region::all_convexes())
interpolation/extrapolation of (mf_source, U) on mf_target.
pintegration_method int_method_descriptor(std::string name, bool throw_if_not_found=true)
Get an integration method from its name .