File: test_interpolated_fem.cc

package info (click to toggle)
getfem 5.4.4%2Bdfsg1-5
  • links: PTS, VCS
  • area: main
  • in suites: sid, trixie
  • size: 31,640 kB
  • sloc: cpp: 126,151; ansic: 24,798; python: 9,244; sh: 3,648; perl: 1,829; makefile: 1,367
file content (352 lines) | stat: -rw-r--r-- 13,406 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
/*===========================================================================

 Copyright (C) 2002-2020 Yves Renard.

 This file is a part of GetFEM

 GetFEM  is  free software;  you  can  redistribute  it  and/or modify it
 under  the  terms  of the  GNU  Lesser General Public License as published
 by  the  Free Software Foundation;  either version 3 of the License,  or
 (at your option) any later version along with the GCC Runtime Library
 Exception either version 3.1 or (at your option) any later version.
 This program  is  distributed  in  the  hope  that it will be useful,  but
 WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
 or  FITNESS  FOR  A PARTICULAR PURPOSE.  See the GNU Lesser General Public
 License and GCC Runtime Library Exception for more details.
 You  should  have received a copy of the GNU Lesser General Public License
 along  with  this program;  if not, write to the Free Software Foundation,
 Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301, USA.

===========================================================================*/
/***************************************************************************/
/*                                                                         */
/*  Test the interpolated fem method.                                      */
/*                                                                         */
/***************************************************************************/

#include "getfem/getfem_assembling.h"
#include "getfem/getfem_export.h"
#include "getfem/getfem_regular_meshes.h"
#include "gmm/gmm.h"
#include "getfem/getfem_interpolated_fem.h"
#include "getfem/getfem_mesh_fem_sum.h"
using std::endl; using std::cout; using std::cerr;
using std::ends; using std::cin;

using bgeot::base_vector;
using bgeot::base_small_vector;
using bgeot::base_node;
using bgeot::scalar_type;
using bgeot::size_type;
using bgeot::dim_type;

typedef gmm::wsvector<scalar_type> sparse_vector_type;
typedef gmm::col_matrix<sparse_vector_type> sparse_matrix_type;
typedef std::vector<scalar_type> linalg_vector;

/**************************************************************************/
/*  structure representing the problem.                                   */
/**************************************************************************/

struct lap_pb {
  getfem::mesh mesh1, mesh2;
  getfem::mesh_im      mim1, mim2;
  getfem::mesh_fem     mef1,  mef2, mefinterpolated;

  scalar_type LX, LY, LZ;
  int NX1, NX2, N, K, KI, integration;

  bgeot::md_param PARAM;

  void assemble(void);
  void init(void);
  lap_pb(void) : mim1(mesh1), mim2(mesh2), mef1(mesh1), mef2(mesh2),
		 mefinterpolated(mesh1) {}
};

void lap_pb::init(void) {
  dal::bit_vector nn;

  /***********************************************************************/
  /* READING PARAMETER FILE                                              */
  /***********************************************************************/
  
  N = int(PARAM.int_value("N", "Domaine dimension"));
  LX = PARAM.real_value("LX", "Size in X");
  LY = PARAM.real_value("LY", "Size in Y");
  LZ = PARAM.real_value("LZ", "Size in Y");
  NX1 = int(PARAM.int_value("NX1", "Nomber of sace steps "));
  NX2 = int(PARAM.int_value("NX2", "Nomber of sace steps "));
  integration = int(PARAM.int_value("INTEGRATION", "integration method"));
  K = int(PARAM.int_value("K", "Finite element degree"));
  KI = int(PARAM.int_value("KI", "Integration degree"));
  
  /***********************************************************************/
  /*  BUILD MESH.                                                        */
  /***********************************************************************/

  cout << "Mesh generation\n";

  base_node org(N); gmm::clear(org);
  std::vector<base_small_vector> vtab(N);
  std::vector<size_type> ref(N);
  std::fill(ref.begin(), ref.end(), NX1);
  for (dim_type i = 0; i < N; i++)
  { 
    vtab[i] = base_small_vector(N); gmm::clear(vtab[i]);
    (vtab[i])[i] = ((i == 0) ? LX : ((i == 1) ? LY : LZ)) / scalar_type(NX1);
  }
  getfem::parallelepiped_regular_simplex_mesh(mesh1, dim_type(N), org,
					      vtab.begin(), ref.begin());
  mesh1.optimize_structure();

  std::fill(ref.begin(), ref.end(), NX2);
  for (dim_type i = 0; i < N; i++)
  { 
    vtab[i] = base_small_vector(N); gmm::clear(vtab[i]);
    (vtab[i])[i] = ((i == 0) ? LX : ((i == 1) ? LY : LZ)) / scalar_type(NX2);
  }
  getfem::parallelepiped_regular_simplex_mesh(mesh2, dim_type(N), org,
					      vtab.begin(), ref.begin());
  mesh2.optimize_structure();
  cout << "Selecting finite element method.\n";
  char meth[500];
  getfem::pintegration_method ppi;
  switch (integration) {
  case 0  : snprintf(meth, 499, "IM_EXACT_SIMPLEX(%d)", int(N)); break;
  case 1  : snprintf(meth, 499, "IM_NC(%d, %d)", int(N), int(KI)); break;
  case 2  : snprintf(meth, 499, "IM_GAUSS1D(%d)", int(KI)); break;
  case 3  : snprintf(meth, 499, "IM_STRUCTURED_COMPOSITE(IM_NC(%d, %d), %d)",
		    int(N), int(2*K), int(KI)); break;
  case 11 : snprintf(meth, 499, "IM_TRIANGLE(1)"); break;
  case 12 : snprintf(meth, 499, "IM_TRIANGLE(2)"); break;
  case 13 : snprintf(meth, 499, "IM_TRIANGLE(3)"); break;
  case 14 : snprintf(meth, 499, "IM_TRIANGLE(4)"); break;
  case 15 : snprintf(meth, 499, "IM_TRIANGLE(5)"); break;
  case 16 : snprintf(meth, 499, "IM_TRIANGLE(6)"); break;
  case 17 : snprintf(meth, 499, "IM_TRIANGLE(7)"); break;
  case 21 : snprintf(meth, 499, "IM_TETRAHEDRON(1)"); break;
  case 22 : snprintf(meth, 499, "IM_TETRAHEDRON(2)"); break;
  case 23 : snprintf(meth, 499, "IM_TETRAHEDRON(3)"); break;
  case 25 : snprintf(meth, 499, "IM_TETRAHEDRON(5)"); break;
  default : GMM_ASSERT1(false, "Undefined integration method");
  }
  ppi = getfem::int_method_descriptor(meth);
  
  snprintf(meth, 499, "FEM_PK(%d,%d)", int(N), int(K));
  nn = mesh1.convex_index(dim_type(N));
  mim1.set_integration_method(nn, ppi);
  mef1.set_finite_element(nn, getfem::fem_descriptor(meth));
  nn = mesh2.convex_index(dim_type(N));
  mim2.set_integration_method(nn, ppi);
  mef2.set_finite_element(nn, getfem::fem_descriptor(meth));
  nn = mesh1.convex_index(dim_type(N));
  mefinterpolated.set_finite_element
    (nn, getfem::new_interpolated_fem(mef2, mim1));
}

void lap_pb::assemble(void) {
  int nb_dof1 = int(mef1.nb_dof()), nb_dof2 = int(mef2.nb_dof());
  sparse_matrix_type RM1(nb_dof2, nb_dof2);
  double sum, diff;
  
  cout << "Number of dof : " << nb_dof1 << " : " << nb_dof2 << endl;

  cout << "Number of dof of interpolated method: "
       << mefinterpolated.nb_dof() << endl;
 
  cout << "Assembling interpolated mass matrix" << endl;
  getfem::asm_mass_matrix(RM1, mim1, mefinterpolated, mefinterpolated);

  cout << "Mass Matrix\n";
  sum = 0.0;
  for (size_type i = 0; i < RM1.nrows(); i++) { 
    cout << "line " << i << " [ ";
    scalar_type slig = 0;
    for (size_type l = 0; l < RM1.nrows(); l++)
      if (RM1(i, l) != 0.0) {
	cout << "(" << l << "," << RM1(i, l) << ")  ";
	slig += RM1(i, l);
      }
    sum += slig;
    cout << "] -> sum(line)=" << slig << endl;
  }
  cout << endl << " sum: " << sum << endl << endl;

  sparse_matrix_type RM2 = sparse_matrix_type(nb_dof2, nb_dof2);
  cout << "Assembling normal mass matrix" << endl;
  getfem::asm_mass_matrix(RM2, mim2, mef2);
  
  cout << "Mass Matrix\n";
  sum = 0.0; diff = 0.0;
  for (size_type i = 0; i < RM2.nrows(); i++) { 
    cout << "line " << i << " [ ";
    scalar_type slig = 0;
    for (size_type l = 0; l < RM2.nrows(); l++) {
      diff += gmm::abs(RM2(i, l) - RM1(i, l));
      if (RM2(i, l) != 0.0) {
	cout << "(" << l << "," << RM2(i, l) << ")  ";
	slig = slig + RM2(i, l);
      }
    }
    sum += slig;
    cout << "] -> sum(line)=" << slig << endl;
  }
  cout << endl << " sum: " << sum << endl << endl;
  cout << endl << " diff: " << diff << "max_norm="
       << gmm::mat_maxnorm(RM1) << endl << endl;  
}


/* integration of a quarter of circle (approximated with a degree 5 segment) against a planar regular mesh */
void test2() {
  getfem::mesh m1, m2;
  std::vector<size_type> nsubdiv(2); nsubdiv[0] = nsubdiv[1] = 5;

  getfem::regular_unit_mesh
    (m1, nsubdiv, bgeot::geometric_trans_descriptor("GT_QK(2,2)"), false);

  std::vector<base_node> v;
  for (size_type k=0; k < 6; ++k) {
    scalar_type c = scalar_type(k)/5. * M_PI/2;
    v.push_back(base_node(cos(c), sin(c)));
  }
  m2.add_convex_by_points(bgeot::geometric_trans_descriptor("GT_PK(1,5)"),
			  v.begin());
  //m2.add_segment_by_points(bgeot::base_node(.45,.35),
  //                         bgeot::base_node(.75,.65));
  //m2.add_segment_by_points(bgeot::base_node(.8,.7),bgeot::base_node(.23,.3));
  getfem::mesh_fem mf1(m1), mf2(m2);
  getfem::mesh_im mim(m2);
  // getfem::pintegration_method pim1
  //   = getfem::int_method_descriptor("IM_QUAD(17)");
  getfem::pintegration_method pim2
    = getfem::int_method_descriptor("IM_GAUSS1D(10)");
  mim.set_integration_method(m2.convex_index(), pim2);
  mf1.set_finite_element(m1.convex_index(),
			 getfem::fem_descriptor("FEM_QK(2,1)"));
  mf2.set_finite_element(m2.convex_index(),
			 getfem::fem_descriptor("FEM_PK(1,3)"));
  
  getfem::mesh_fem mflnk(m2);
  getfem::pfem ifem = getfem::new_interpolated_fem(mf1, mim);
  mflnk.set_finite_element(m2.convex_index(), ifem);
  sparse_matrix_type MM = sparse_matrix_type(mf2.nb_dof(),
					     mflnk.nb_dof());
  getfem::asm_mass_matrix(MM, mim, mf2, mflnk);
  cout << "MM=" << MM << "\n";
  cout << "mflnk.nb_dof()=" << mflnk.nb_dof() << ", mf2.nb_dof()="
       << mf2.nb_dof() << ", mf1.nb_dof=" << mf1.nb_dof() << "\n";
  cout << "Mass Matrix\n";
  scalar_type sum = 0.0; 
  for (size_type i = 0; i < MM.nrows(); i++) { 
    scalar_type slig = 0;
    for (size_type l = 0; l < MM.ncols(); l++) {
      slig = slig + MM(i, l);
      // cout << "M(" << i << "," << l << ")=" << MM(i,l) << ", slig = "
      //      << slig << "\n";
    }
    sum += slig;
    cout << "sum(line)=" << slig << endl;
  }
  cout << endl << " sum: " << sum << endl << endl;
  GMM_ASSERT1(gmm::abs(sum - M_PI/2.0) < 1e-5,
	      "Test for scalar fem failed : " << gmm::abs(sum - M_PI/2.0));
  sparse_matrix_type MM2 = sparse_matrix_type(mf2.nb_dof(), mf2.nb_dof());
  getfem::asm_mass_matrix(MM2, mim, mf2, mf2);
  cout << "MM2=" << MM2 << "\n";
  sparse_matrix_type MM3 = sparse_matrix_type(mflnk.nb_dof(), mflnk.nb_dof());

  asm_stiffness_matrix_for_homogeneous_laplacian(MM3, mim, mflnk);
}


/* integration of a quarter of circle (approximated with a degree 5 segment) against a planar regular mesh : vector version */
void test3() {
  getfem::mesh m1, m2;
  std::vector<size_type> nsubdiv(2); nsubdiv[0] = nsubdiv[1] = 5;

  getfem::regular_unit_mesh
    (m1, nsubdiv, bgeot::geometric_trans_descriptor("GT_QK(2,2)"), false);

  std::vector<base_node> v;
  for (size_type k=0; k < 6; ++k) {
    scalar_type c = scalar_type(k)/5. * M_PI/2;
    v.push_back(base_node(cos(c), sin(c)));
  }
  m2.add_convex_by_points(bgeot::geometric_trans_descriptor("GT_PK(1,5)"),
			  v.begin());
  //m2.add_segment_by_points(bgeot::base_node(.45,.35),
  //                         bgeot::base_node(.75,.65));
  //m2.add_segment_by_points(bgeot::base_node(.8,.7),bgeot::base_node(.23,.3));
  getfem::mesh_fem mf1(m1), mf2(m2);
  getfem::mesh_im mim(m2);
  // getfem::pintegration_method pim1
  //   = getfem::int_method_descriptor("IM_QUAD(17)");
  getfem::pintegration_method pim2
    = getfem::int_method_descriptor("IM_GAUSS1D(10)");
  mim.set_integration_method(m2.convex_index(), pim2);
  mf1.set_finite_element(m1.convex_index(),
			 getfem::fem_descriptor("FEM_QK(2,1)"));
  mf1.set_qdim(2);
  mf2.set_finite_element(m2.convex_index(),
			 getfem::fem_descriptor("FEM_PK(1,3)"));
  mf2.set_qdim(2);
  
  getfem::mesh_fem mflnk(m2);
  mflnk.set_qdim(2);
  getfem::pfem ifem = getfem::new_interpolated_fem(mf1, mim);
  mflnk.set_finite_element(m2.convex_index(), ifem);
  sparse_matrix_type MM = sparse_matrix_type(mf2.nb_dof(), mflnk.nb_dof());
  getfem::asm_mass_matrix(MM, mim, mf2, mflnk);
  cout << "MM=" << MM << "\n";
  cout << "mflnk.nb_dof()=" << mflnk.nb_dof() << ", mf2.nb_dof()="
       << mf2.nb_dof() << ", mf1.nb_dof=" << mf1.nb_dof() << "\n";
  cout << "Mass matrix\n";
  scalar_type sum = 0.0; 
  for (size_type i = 0; i < MM.nrows(); i++) { 
    scalar_type slig = 0;
    for (size_type l = 0; l < MM.ncols(); l++) {
      slig = slig + MM(i, l);
      // cout << "M(" << i << "," << l << ")=" << MM(i,l) << ", slig = "
      //      << slig << "\n";
    }
    sum += slig;
    cout << "sum(line)=" << slig << endl;
  }
  cout << endl << " sum: " << sum << endl << endl;
  GMM_ASSERT1(gmm::abs(sum - M_PI) < 5e-5,
	      "Test for vectorial fem failed : " << gmm::abs(sum - M_PI));

  // Il faudrait tester les gradients aussi ...
}







/**************************************************************************/
/*  main program.                                                         */
/**************************************************************************/

int main(int argc, char *argv[]) {
  
  GETFEM_MPI_INIT(argc, argv);

  test3(); // test for vectorial fems
  test2(); // test for scalar fems
  lap_pb p;
    
  cout << "initialisation ...\n";
  p.PARAM.read_command_line(argc, argv);
  p.init();
  
  cout << "Assembling \n";
  p.assemble();
  
  GETFEM_MPI_FINALIZE;
  
  return 0; 
}