GetFEM  5.5
helmholtz.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 /**
22  @file helmholtz.cc
23  @brief Helmholtz problem (Delta(u) + k^2 u = 0)
24 
25  Diffraction of a plane wave by a circular obstacle.
26 
27  This program is used to check that getfem++ is working. This is also
28  a good example of use of GetFEM.
29 */
30 
31 #include "getfem/getfem_assembling.h" /* import assembly methods (and norms comp.) */
32 #include "getfem/getfem_export.h" /* export functions (save solution in a file) */
35 #include "gmm/gmm.h"
36 using std::endl; using std::cout; using std::cerr;
37 using std::ends; using std::cin;
38 
39 /* some GetFEM types that we will be using */
40 using bgeot::base_small_vector; /* special class for small (dim<16) vectors */
41 using bgeot::base_node; /* geometrical nodes(derived from base_small_vector)*/
42 using bgeot::scalar_type; /* = double */
43 using bgeot::complex_type; /* = std::complex<double> */
44 using bgeot::size_type; /* = unsigned long */
45 using bgeot::short_type;
46 using bgeot::base_matrix; /* small dense matrix. */
47 
48 /* definition of some matrix/vector types. These ones are built
49  using the predefined types in Gmm++ */
51 typedef getfem::modeling_standard_complex_sparse_matrix sparse_matrix;
52 typedef getfem::modeling_standard_complex_plain_vector plain_vector;
53 
54 /*
55  structure for the Helmholtz problem
56 */
57 struct Helmholtz_problem {
58 
59  enum { DIRICHLET_BOUNDARY_NUM = 0, ROBIN_BOUNDARY_NUM = 1};
60  getfem::mesh mesh; /* the mesh */
61  getfem::mesh_im mim; /* the integration methods */
62  getfem::mesh_fem mf_u; /* main mesh_fem, for the elastostatic solution */
63  getfem::mesh_fem mf_rhs; /* mesh_fem for the right hand side (f(x),..) */
64  complex_type wave_number;
65 
66  scalar_type residual; /* max residual for the iterative solvers */
67  int with_mult;
68 
69  std::string datafilename;
70  bgeot::md_param PARAM;
71 
72  bool solve(plain_vector &U);
73  void init(void);
74  void compute_error(plain_vector &U);
75  Helmholtz_problem(void) : mim(mesh), mf_u(mesh), mf_rhs(mesh) {}
76 };
77 
78 complex_type __wave_number;
79 
80 complex_type incoming_field(const base_node &P) {
81  return complex_type(cos(__wave_number.real()*P[1]+.2),
82  sin(__wave_number.real()*P[1]+.2));
83  /*scalar_type s = 0;
84  for (size_type i=1; i < P.size(); ++i) s += P[i]*(1.-P[i]);
85  s = rand()*3. / RAND_MAX;
86  return s;
87  */
88 }
89 
90 /* Read parameters from the .param file, build the mesh, set finite element
91  * and integration methods and selects the boundaries.
92  */
93 void Helmholtz_problem::init(void) {
94  std::string FEM_TYPE = PARAM.string_value("FEM_TYPE","FEM name");
95  std::string INTEGRATION = PARAM.string_value("INTEGRATION",
96  "Name of integration method");
97  cout << "FEM_TYPE=" << FEM_TYPE << "\n";
98  cout << "INTEGRATION=" << INTEGRATION << "\n";
99 
100  /* First step : build the mesh */
101  size_type Nt = PARAM.int_value("NTHETA", "Nomber of space steps "),
102  Nr=PARAM.int_value("NR", "Nomber of space steps ");
103  size_type gt_order = PARAM.int_value("GTDEGREE",
104  "polynomial degree of geometric transformation");
105  scalar_type dtheta=2*M_PI*1./scalar_type(Nt);
106  scalar_type R0 = PARAM.real_value("R0","R0");
107  scalar_type R1 = PARAM.real_value("R1","R1");
108  scalar_type dR = (R1-R0)/scalar_type(Nr-1);
110  pgt = bgeot::parallelepiped_geotrans(2, short_type(gt_order));
111  for (size_type i=0; i < Nt; ++i) {
112  for (size_type j=0; j < Nr-1; ++j) {
113  std::vector<size_type> ipts; ipts.reserve(gmm::sqr(gt_order+1));
114  for (size_type ii=0; ii <= gt_order; ++ii) {
115  for (size_type jj=0; jj <= gt_order; ++jj) {
116  scalar_type r = R0 + scalar_type(j)*dR
117  + scalar_type(jj)*(dR/scalar_type(gt_order));
118  scalar_type t = scalar_type(i)*dtheta
119  + scalar_type(ii)*dtheta/scalar_type(gt_order);
120  ipts.push_back(mesh.add_point(base_node(r*cos(t),r*sin(t))));
121  }
122  }
123  mesh.add_convex(pgt, ipts.begin());
124  }
125  }
126 
127  datafilename = PARAM.string_value("ROOTFILENAME","Base name of data files.");
128  residual = PARAM.real_value("RESIDUAL"); if (residual == 0.) residual = 1e-10;
129 
130  __wave_number = wave_number = complex_type
131  (PARAM.real_value("WAVENUM_R", "Real part of the wave number"),
132  PARAM.real_value("WAVENUM_I", "Imaginary part of the wave number"));
133 
134  with_mult = int(PARAM.int_value("DIRICHLET_VERSION",
135  "Dirichlet condition version"));
136 
137  /* set the finite element on the mf_u */
138  getfem::pfem pf_u =
139  getfem::fem_descriptor(FEM_TYPE);
140  getfem::pintegration_method ppi =
141  getfem::int_method_descriptor(INTEGRATION);
142 
143  mim.set_integration_method(ppi);
144  mf_u.set_finite_element(pf_u);
145 
146  /* set the finite element on mf_rhs (same as mf_u is DATA_FEM_TYPE is
147  not used in the .param file */
148  std::string data_fem_name = PARAM.string_value("DATA_FEM_TYPE");
149  if (data_fem_name.size() == 0) {
150  GMM_ASSERT1(pf_u->is_lagrange(), "You are using a non-lagrange FEM. "
151  << "In that case you need to set "
152  << "DATA_FEM_TYPE in the .param file");
153  mf_rhs.set_finite_element(pf_u);
154  } else {
155  mf_rhs.set_finite_element(getfem::fem_descriptor(data_fem_name));
156  }
157 
158 
159  /* select boundaries */
160  cout << "Selecting Robin and Dirichlet boundaries\n";
161  getfem::mesh_region border_faces;
162  getfem::outer_faces_of_mesh(mesh, border_faces);
163  for (getfem::mr_visitor i(border_faces); !i.finished(); ++i) {
164  assert(i.is_face());
165  if (gmm::vect_norm2(mesh.points_of_face_of_convex(i.cv(),
166  i.f())[0]) > 5.) {
167  mesh.region(ROBIN_BOUNDARY_NUM).add(i.cv(),i.f());
168  } else mesh.region(DIRICHLET_BOUNDARY_NUM).add(i.cv(),i.f());
169  }
170 }
171 
172 /**************************************************************************/
173 /* Model. */
174 /**************************************************************************/
175 
176 
177 bool Helmholtz_problem::solve(plain_vector &U) {
178 
179  // Complex model.
180  getfem::model Helmholtz_model(true);
181 
182  // Main unknown of the problem
183  Helmholtz_model.add_fem_variable("u", mf_u);
184 
185  // Helmholtz term on u.
186  plain_vector K(1); K[0] = wave_number;
187  Helmholtz_model.add_initialized_fixed_size_data("k", K);
188  add_Helmholtz_brick(Helmholtz_model, mim, "u", "k");
189 
190  // Fourier-Robin condition.
191  plain_vector Q(1); Q[0] = wave_number * complex_type(0,1.);
192  Helmholtz_model.add_initialized_fixed_size_data("Q", Q);
193  add_Fourier_Robin_brick(Helmholtz_model, mim, "u", "Q", ROBIN_BOUNDARY_NUM);
194 
195  // Dirichlet condition
196  plain_vector F(mf_rhs.nb_dof());
197  getfem::interpolation_function(mf_rhs, F, incoming_field);
198  Helmholtz_model.add_initialized_fem_data("DirichletData", mf_rhs, F);
200  (Helmholtz_model, mim, "u", mf_u,
201  DIRICHLET_BOUNDARY_NUM, "DirichletData");
202 
203  // Helmholtz_model.listvar(cout);
204 
205  gmm::iteration iter(residual, 1, 40000);
206  getfem::standard_solve(Helmholtz_model, iter);
207 
208  gmm::resize(U, mf_u.nb_dof());
209  gmm::copy(Helmholtz_model.complex_variable("u"), U);
210 
211  // cout << "U = " << U << endl;
212 
213  return (iter.converged());
214 }
215 
216 /**************************************************************************/
217 /* main program. */
218 /**************************************************************************/
219 
220 int main(int argc, char *argv[]) {
221 
222  GETFEM_MPI_INIT(argc, argv);
223  GMM_SET_EXCEPTION_DEBUG; // Exceptions make a memory fault, to debug.
224  FE_ENABLE_EXCEPT; // Enable floating point exception for Nan.
225 
226  Helmholtz_problem p;
227  p.PARAM.read_command_line(argc, argv);
228  p.init();
229  plain_vector U(p.mf_u.nb_dof());
230  if (!p.solve(U)) GMM_ASSERT1(false, "Solve has failed");
231 
232  if (p.PARAM.int_value("VTK_EXPORT")) {
233  cout << "export to " << p.datafilename + ".vtk" << "..\n";
234  getfem::vtk_export vtk_exp(p.datafilename + ".vtk",
235  p.PARAM.int_value("VTK_EXPORT")==1, true);
236  cout << "export to " << p.datafilename + ".vtu" << "..\n";
237  getfem::vtu_export vtu_exp(p.datafilename + ".vtu");
238  getfem::stored_mesh_slice sl(p.mesh, p.mesh.nb_convex() < 2000 ? 8 : 6);
239  vtk_exp.exporting(sl);
240  vtk_exp.write_point_data(p.mf_u, gmm::real_part(U), "helmholtz_rfield");
241  vtk_exp.write_point_data(p.mf_u, gmm::imag_part(U), "helmholtz_ifield");
242  vtu_exp.exporting(sl);
243  vtu_exp.write_point_data(p.mf_u, gmm::real_part(U), "helmholtz_rfield");
244  vtu_exp.write_point_data(p.mf_u, gmm::imag_part(U), "helmholtz_ifield");
245  cout << "export done, you can view the data file with (for example)\n"
246  "mayavi2 -d helmholtz.vtk -f WarpScalar -m Surface -m Outline"
247  "\n";
248  }
249 
250  GETFEM_MPI_FINALIZE;
251 
252  return 0;
253 }
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
`‘Model’' variables store the variables, the data and the description of a model.
The output of a getfem::mesh_slicer which has been recorded.
VTK/VTU export.
Definition: getfem_export.h:67
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....
Export solutions to various formats.
Standard solvers for model bricks.
Build regular meshes.
Include common gmm files.
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:976
void resize(V &v, size_type n)
*‍/
Definition: gmm_blas.h:209
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
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
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
size_type APIDECL add_Dirichlet_condition_with_multipliers(model &md, const mesh_im &mim, const std::string &varname, const std::string &multname, size_type region, const std::string &dataname=std::string())
Add a Dirichlet condition on the variable varname and the mesh region region.
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.
size_type APIDECL add_Helmholtz_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region=size_type(-1))
Add a Helmoltz brick to the model.
pintegration_method int_method_descriptor(std::string name, bool throw_if_not_found=true)
Get an integration method from its name .
size_type APIDECL add_Fourier_Robin_brick(model &md, const mesh_im &mim, const std::string &varname, const std::string &dataexpr, size_type region)
Add a Fourier-Robin brick to the model.
void standard_solve(model &md, gmm::iteration &iter, rmodel_plsolver_type lsolver, abstract_newton_line_search &ls)
A default solver for the model brick system.