MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
gslib.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "gslib.hpp"
13#include "geom.hpp"
14
15#ifdef MFEM_USE_GSLIB
16
17// Ignore warnings from the gslib header (GCC version)
18#ifdef MFEM_HAVE_GCC_PRAGMA_DIAGNOSTIC
19#pragma GCC diagnostic push
20#pragma GCC diagnostic ignored "-Wunused-function"
21#endif
22
23// External GSLIB header (the MFEM header is gslib.hpp)
24namespace gslib
25{
26#include "gslib.h"
27}
28
29#ifdef MFEM_HAVE_GCC_PRAGMA_DIAGNOSTIC
30#pragma GCC diagnostic pop
31#endif
32
33namespace mfem
34{
35
37 : mesh(NULL),
38 fec_map_lin(NULL),
39 fdata2D(NULL), fdata3D(NULL), cr(NULL), gsl_comm(NULL),
40 dim(-1), points_cnt(0), setupflag(false), default_interp_value(0),
41 avgtype(AvgType::ARITHMETIC), bdr_tol(1e-8)
42{
43 mesh_split.SetSize(4);
44 ir_split.SetSize(4);
45 fes_rst_map.SetSize(4);
46 gf_rst_map.SetSize(4);
47 for (int i = 0; i < mesh_split.Size(); i++)
48 {
49 mesh_split[i] = NULL;
50 ir_split[i] = NULL;
51 fes_rst_map[i] = NULL;
52 gf_rst_map[i] = NULL;
53 }
54
55 gsl_comm = new gslib::comm;
56 cr = new gslib::crystal;
57#ifdef MFEM_USE_MPI
58 int initialized;
59 MPI_Initialized(&initialized);
60 if (!initialized) { MPI_Init(NULL, NULL); }
61 MPI_Comm comm = MPI_COMM_WORLD;
62 comm_init(gsl_comm, comm);
63#else
64 comm_init(gsl_comm, 0);
65#endif
66}
67
69{
70 comm_free(gsl_comm);
71 delete gsl_comm;
72 delete cr;
73 for (int i = 0; i < 4; i++)
74 {
75 if (mesh_split[i]) { delete mesh_split[i]; mesh_split[i] = NULL; }
76 if (ir_split[i]) { delete ir_split[i]; ir_split[i] = NULL; }
77 if (fes_rst_map[i]) { delete fes_rst_map[i]; fes_rst_map[i] = NULL; }
78 if (gf_rst_map[i]) { delete gf_rst_map[i]; gf_rst_map[i] = NULL; }
79 }
80 if (fec_map_lin) { delete fec_map_lin; fec_map_lin = NULL; }
81}
82
83#ifdef MFEM_USE_MPI
85 : mesh(NULL),
86 fec_map_lin(NULL),
87 fdata2D(NULL), fdata3D(NULL), cr(NULL), gsl_comm(NULL),
88 dim(-1), points_cnt(0), setupflag(false), default_interp_value(0),
89 avgtype(AvgType::ARITHMETIC), bdr_tol(1e-8)
90{
91 mesh_split.SetSize(4);
92 ir_split.SetSize(4);
93 fes_rst_map.SetSize(4);
94 gf_rst_map.SetSize(4);
95 for (int i = 0; i < mesh_split.Size(); i++)
96 {
97 mesh_split[i] = NULL;
98 ir_split[i] = NULL;
99 fes_rst_map[i] = NULL;
100 gf_rst_map[i] = NULL;
101 }
102
103 gsl_comm = new gslib::comm;
104 cr = new gslib::crystal;
105 comm_init(gsl_comm, comm_);
106}
107#endif
108
109void FindPointsGSLIB::Setup(Mesh &m, const double bb_t, const double newt_tol,
110 const int npt_max)
111{
112 MFEM_VERIFY(m.GetNodes() != NULL, "Mesh nodes are required.");
113 const int meshOrder = m.GetNodes()->FESpace()->GetMaxElementOrder();
114
115 // call FreeData if FindPointsGSLIB::Setup has been called already
116 if (setupflag) { FreeData(); }
117
118 crystal_init(cr, gsl_comm);
119 mesh = &m;
120 dim = mesh->Dimension();
121 unsigned dof1D = meshOrder + 1;
122
124 if (dim == 2)
125 {
126 if (ir_split[0]) { delete ir_split[0]; ir_split[0] = NULL; }
127 ir_split[0] = new IntegrationRule(3*pow(dof1D, dim));
129
130 if (ir_split[1]) { delete ir_split[1]; ir_split[1] = NULL; }
131 ir_split[1] = new IntegrationRule(pow(dof1D, dim));
133 }
134 else if (dim == 3)
135 {
136 if (ir_split[0]) { delete ir_split[0]; ir_split[0] = NULL; }
137 ir_split[0] = new IntegrationRule(pow(dof1D, dim));
139
140 if (ir_split[1]) { delete ir_split[1]; ir_split[1] = NULL; }
141 ir_split[1] = new IntegrationRule(4*pow(dof1D, dim));
143
144 if (ir_split[2]) { delete ir_split[2]; ir_split[2] = NULL; }
145 ir_split[2] = new IntegrationRule(3*pow(dof1D, dim));
147
148 if (ir_split[3]) { delete ir_split[3]; ir_split[3] = NULL; }
149 ir_split[3] = new IntegrationRule(8*pow(dof1D, dim));
151 }
152
154
155 const int pts_cnt = gsl_mesh.Size()/dim,
156 NEtot = pts_cnt/(int)pow(dof1D, dim);
157
158 if (dim == 2)
159 {
160 unsigned nr[2] = { dof1D, dof1D };
161 unsigned mr[2] = { 2*dof1D, 2*dof1D };
162 double * const elx[2] =
163 {
164 pts_cnt == 0 ? nullptr : &gsl_mesh(0),
165 pts_cnt == 0 ? nullptr : &gsl_mesh(pts_cnt)
166 };
167 fdata2D = findpts_setup_2(gsl_comm, elx, nr, NEtot, mr, bb_t,
168 pts_cnt, pts_cnt, npt_max, newt_tol);
169 }
170 else
171 {
172 unsigned nr[3] = { dof1D, dof1D, dof1D };
173 unsigned mr[3] = { 2*dof1D, 2*dof1D, 2*dof1D };
174 double * const elx[3] =
175 {
176 pts_cnt == 0 ? nullptr : &gsl_mesh(0),
177 pts_cnt == 0 ? nullptr : &gsl_mesh(pts_cnt),
178 pts_cnt == 0 ? nullptr : &gsl_mesh(2*pts_cnt)
179 };
180 fdata3D = findpts_setup_3(gsl_comm, elx, nr, NEtot, mr, bb_t,
181 pts_cnt, pts_cnt, npt_max, newt_tol);
182 }
183 setupflag = true;
184}
185
187 int point_pos_ordering)
188{
189 MFEM_VERIFY(setupflag, "Use FindPointsGSLIB::Setup before finding points.");
190 points_cnt = point_pos.Size() / dim;
196
197 auto xvFill = [&](const double *xv_base[], unsigned xv_stride[])
198 {
199 for (int d = 0; d < dim; d++)
200 {
201 if (point_pos_ordering == Ordering::byNODES)
202 {
203 xv_base[d] = point_pos.GetData() + d*points_cnt;
204 xv_stride[d] = sizeof(double);
205 }
206 else
207 {
208 xv_base[d] = point_pos.GetData() + d;
209 xv_stride[d] = dim*sizeof(double);
210 }
211 }
212 };
213 if (dim == 2)
214 {
215 const double *xv_base[2];
216 unsigned xv_stride[2];
217 xvFill(xv_base, xv_stride);
218 findpts_2(gsl_code.GetData(), sizeof(unsigned int),
219 gsl_proc.GetData(), sizeof(unsigned int),
220 gsl_elem.GetData(), sizeof(unsigned int),
221 gsl_ref.GetData(), sizeof(double) * dim,
222 gsl_dist.GetData(), sizeof(double),
223 xv_base, xv_stride, points_cnt, fdata2D);
224 }
225 else // dim == 3
226 {
227 const double *xv_base[3];
228 unsigned xv_stride[3];
229 xvFill(xv_base, xv_stride);
230 findpts_3(gsl_code.GetData(), sizeof(unsigned int),
231 gsl_proc.GetData(), sizeof(unsigned int),
232 gsl_elem.GetData(), sizeof(unsigned int),
233 gsl_ref.GetData(), sizeof(double) * dim,
234 gsl_dist.GetData(), sizeof(double),
235 xv_base, xv_stride, points_cnt, fdata3D);
236 }
237
238 // Set the element number and reference position to 0 for points not found
239 for (int i = 0; i < points_cnt; i++)
240 {
241 if (gsl_code[i] == 2 ||
242 (gsl_code[i] == 1 && gsl_dist(i) > bdr_tol))
243 {
244 gsl_elem[i] = 0;
245 for (int d = 0; d < dim; d++) { gsl_ref(i*dim + d) = -1.; }
246 gsl_code[i] = 2;
247 }
248 }
249
250 // Map element number for simplices, and ref_pos from [-1,1] to [0,1] for
251 // both simplices and quads. Also sets code to 1 for points found on element
252 // faces/edges.
254}
255
256void FindPointsGSLIB::FindPoints(Mesh &m, const Vector &point_pos,
257 int point_pos_ordering, const double bb_t,
258 const double newt_tol, const int npt_max)
259{
260 if (!setupflag || (mesh != &m) )
261 {
262 Setup(m, bb_t, newt_tol, npt_max);
263 }
264 FindPoints(point_pos, point_pos_ordering);
265}
266
268 const GridFunction &field_in, Vector &field_out,
269 int point_pos_ordering)
270{
271 FindPoints(point_pos, point_pos_ordering);
272 Interpolate(field_in, field_out);
273}
274
275void FindPointsGSLIB::Interpolate(Mesh &m, const Vector &point_pos,
276 const GridFunction &field_in, Vector &field_out,
277 int point_pos_ordering)
278{
279 FindPoints(m, point_pos, point_pos_ordering);
280 Interpolate(field_in, field_out);
281}
282
284{
285 if (!setupflag) { return; }
286 crystal_free(cr);
287 if (dim == 2)
288 {
289 findpts_free_2(fdata2D);
290 }
291 else
292 {
293 findpts_free_3(fdata3D);
294 }
301 for (int i = 0; i < 4; i++)
302 {
303 if (mesh_split[i]) { delete mesh_split[i]; mesh_split[i] = NULL; }
304 if (ir_split[i]) { delete ir_split[i]; ir_split[i] = NULL; }
305 if (fes_rst_map[i]) { delete fes_rst_map[i]; fes_rst_map[i] = NULL; }
306 if (gf_rst_map[i]) { delete gf_rst_map[i]; gf_rst_map[i] = NULL; }
307 }
308 if (fec_map_lin) { delete fec_map_lin; fec_map_lin = NULL; }
309 setupflag = false;
310}
311
313{
315 if (mesh->Dimension() == 2)
316 {
317 int Nvert = 7;
318 int NEsplit = 3;
319 mesh_split[0] = new Mesh(2, Nvert, NEsplit, 0, 2);
320
321 const double quad_v[7][2] =
322 {
323 {0, 0}, {0.5, 0}, {1, 0}, {0, 0.5},
324 {1./3., 1./3.}, {0.5, 0.5}, {0, 1}
325 };
326 const int quad_e[3][4] =
327 {
328 {0, 1, 4, 3}, {1, 2, 5, 4}, {3, 4, 5, 6}
329 };
330
331 for (int j = 0; j < Nvert; j++)
332 {
333 mesh_split[0]->AddVertex(quad_v[j]);
334 }
335 for (int j = 0; j < NEsplit; j++)
336 {
337 int attribute = j + 1;
338 mesh_split[0]->AddQuad(quad_e[j], attribute);
339 }
340 mesh_split[0]->FinalizeQuadMesh(1, 1, true);
341
344 const int npt = gf_rst_map[0]->Size()/dim;
345 for (int k = 0; k < dim; k++)
346 {
347 for (int j = 0; j < npt; j++)
348 {
349 (*gf_rst_map[0])(j+k*npt) = quad_v[j][k];
350 }
351 }
352
355 }
356 else if (mesh->Dimension() == 3)
357 {
358 mesh_split[0] = new Mesh(Mesh::MakeCartesian3D(1, 1, 1,
360 // Tetrahedron
361 {
362 int Nvert = 15;
363 int NEsplit = 4;
364 mesh_split[1] = new Mesh(3, Nvert, NEsplit, 0, 3);
365
366 const double hex_v[15][3] =
367 {
368 {0, 0, 0.}, {1, 0., 0.}, {0., 1., 0.}, {0, 0., 1.},
369 {0.5, 0., 0.}, {0.5, 0.5, 0.}, {0., 0.5, 0.},
370 {0., 0., 0.5}, {0.5, 0., 0.5}, {0., 0.5, 0.5},
371 {1./3., 0., 1./3.}, {1./3., 1./3., 1./3.}, {0, 1./3., 1./3.},
372 {1./3., 1./3., 0}, {0.25, 0.25, 0.25}
373 };
374 const int hex_e[4][8] =
375 {
376 {7, 10, 4, 0, 12, 14, 13, 6},
377 {10, 8, 1, 4, 14, 11, 5, 13},
378 {14, 11, 5, 13, 12, 9, 2, 6},
379 {7, 3, 8, 10, 12, 9, 11, 14}
380 };
381
382 for (int j = 0; j < Nvert; j++)
383 {
384 mesh_split[1]->AddVertex(hex_v[j]);
385 }
386 for (int j = 0; j < NEsplit; j++)
387 {
388 int attribute = j + 1;
389 mesh_split[1]->AddHex(hex_e[j], attribute);
390 }
391 mesh_split[1]->FinalizeHexMesh(1, 1, true);
392
395 const int npt = gf_rst_map[1]->Size()/dim;
396 for (int k = 0; k < dim; k++)
397 {
398 for (int j = 0; j < npt; j++)
399 {
400 (*gf_rst_map[1])(j+k*npt) = hex_v[j][k];
401 }
402 }
403 }
404 // Prism
405 {
406 int Nvert = 14;
407 int NEsplit = 3;
408 mesh_split[2] = new Mesh(3, Nvert, NEsplit, 0, 3);
409
410 const double hex_v[14][3] =
411 {
412 {0, 0, 0}, {0.5, 0, 0}, {1, 0, 0}, {0, 0.5, 0},
413 {1./3., 1./3., 0}, {0.5, 0.5, 0}, {0, 1, 0},
414 {0, 0, 1}, {0.5, 0, 1}, {1, 0, 1}, {0, 0.5, 1},
415 {1./3., 1./3., 1}, {0.5, 0.5, 1}, {0, 1, 1}
416 };
417 const int hex_e[3][8] =
418 {
419 {0, 1, 4, 3, 7, 8, 11, 10},
420 {1, 2, 5, 4, 8, 9, 12, 11},
421 {3, 4, 5, 6, 10, 11, 12, 13}
422 };
423
424 for (int j = 0; j < Nvert; j++)
425 {
426 mesh_split[2]->AddVertex(hex_v[j]);
427 }
428 for (int j = 0; j < NEsplit; j++)
429 {
430 int attribute = j + 1;
431 mesh_split[2]->AddHex(hex_e[j], attribute);
432 }
433 mesh_split[2]->FinalizeHexMesh(1, 1, true);
434
437 const int npt = gf_rst_map[2]->Size()/dim;
438 for (int k = 0; k < dim; k++)
439 {
440 for (int j = 0; j < npt; j++)
441 {
442 (*gf_rst_map[2])(j+k*npt) = hex_v[j][k];
443 }
444 }
445 }
446 // Pyramid
447 {
448 int Nvert = 23;
449 int NEsplit = 8;
450 mesh_split[3] = new Mesh(3, Nvert, NEsplit, 0, 3);
451
452 const double hex_v[23][3] =
453 {
454 {0.0000, 0.0000, 0.0000}, {0.5000, 0.0000, 0.0000},
455 {0.0000, 0.0000, 0.5000}, {0.3333, 0.0000, 0.3333},
456 {0.0000, 0.5000, 0.0000}, {0.3333, 0.3333, 0.0000},
457 {0.0000, 0.3333, 0.3333}, {0.2500, 0.2500, 0.2500},
458 {1.0000, 0.0000, 0.0000}, {0.5000, 0.0000, 0.5000},
459 {0.5000, 0.5000, 0.0000}, {0.3333, 0.3333, 0.3333},
460 {0.0000, 1.0000, 0.0000}, {0.0000, 0.5000, 0.5000},
461 {0.0000, 0.0000, 1.0000}, {1.0000, 0.5000, 0.0000},
462 {0.6667, 0.3333, 0.3333}, {0.6667, 0.6667, 0.0000},
463 {0.5000, 0.5000, 0.2500}, {1.0000, 1.0000, 0.0000},
464 {0.5000, 0.5000, 0.5000}, {0.5000, 1.0000, 0.0000},
465 {0.3333, 0.6667, 0.3333}
466 };
467 const int hex_e[8][8] =
468 {
469 {2, 3, 1, 0, 6, 7, 5, 4}, {3, 9, 8, 1, 7, 11, 10, 5},
470 {7, 11, 10, 5, 6, 13, 12, 4}, {2, 14, 9, 3, 6, 13, 11, 7},
471 {9, 16, 15, 8, 11, 18, 17, 10}, {16, 20, 19, 15, 18, 22, 21, 17},
472 {18, 22, 21, 17, 11, 13, 12, 10}, {9, 14, 20, 16, 11, 13, 22, 18}
473 };
474
475 for (int j = 0; j < Nvert; j++)
476 {
477 mesh_split[3]->AddVertex(hex_v[j]);
478 }
479 for (int j = 0; j < NEsplit; j++)
480 {
481 int attribute = j + 1;
482 mesh_split[3]->AddHex(hex_e[j], attribute);
483 }
484 mesh_split[3]->FinalizeHexMesh(1, 1, true);
485
488 const int npt = gf_rst_map[3]->Size()/dim;
489 for (int k = 0; k < dim; k++)
490 {
491 for (int j = 0; j < npt; j++)
492 {
493 (*gf_rst_map[3])(j+k*npt) = hex_v[j][k];
494 }
495 }
496 }
497 }
498
499 NE_split_total = 0;
502 int NEsplit = 0;
503 for (int e = 0; e < mesh->GetNE(); e++)
504 {
506 if (gt == Geometry::TRIANGLE || gt == Geometry::PRISM)
507 {
508 NEsplit = 3;
509 }
510 else if (gt == Geometry::TETRAHEDRON)
511 {
512 NEsplit = 4;
513 }
514 else if (gt == Geometry::PYRAMID)
515 {
516 NEsplit = 8;
517 }
518 else if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
519 {
520 NEsplit = 1;
521 }
522 else
523 {
524 MFEM_ABORT("Unsupported geometry type.");
525 }
526 NE_split_total += NEsplit;
527 for (int i = 0; i < NEsplit; i++)
528 {
531 }
532 }
533}
534
536 IntegrationRule *irule,
537 int order)
538{
539 H1_FECollection fec(order, dim);
540 FiniteElementSpace nodal_fes(meshin, &fec, dim);
541 meshin->SetNodalFESpace(&nodal_fes);
542 const int NEsplit = meshin->GetNE();
543
544 const int dof_cnt = nodal_fes.GetFE(0)->GetDof(),
545 pts_cnt = NEsplit * dof_cnt;
546 Vector irlist(dim * pts_cnt);
547
548 const TensorBasisElement *tbe =
549 dynamic_cast<const TensorBasisElement *>(nodal_fes.GetFE(0));
550 MFEM_VERIFY(tbe != NULL, "TensorBasis FiniteElement expected.");
551 const Array<int> &dof_map = tbe->GetDofMap();
552
553 DenseMatrix pos(dof_cnt, dim);
554 Vector posV(pos.Data(), dof_cnt * dim);
555 Array<int> xdofs(dof_cnt * dim);
556
557 // Create an IntegrationRule on the nodes of the reference submesh.
558 MFEM_ASSERT(irule->GetNPoints() == pts_cnt, "IntegrationRule does not have"
559 "the correct number of points.");
560 GridFunction *nodesplit = meshin->GetNodes();
561 int pt_id = 0;
562 for (int i = 0; i < NEsplit; i++)
563 {
564 nodal_fes.GetElementVDofs(i, xdofs);
565 nodesplit->GetSubVector(xdofs, posV);
566 for (int j = 0; j < dof_cnt; j++)
567 {
568 for (int d = 0; d < dim; d++)
569 {
570 irlist(pts_cnt * d + pt_id) = pos(dof_map[j], d);
571 }
572 irule->IntPoint(pt_id).x = irlist(pt_id);
573 irule->IntPoint(pt_id).y = irlist(pts_cnt + pt_id);
574 if (dim == 3)
575 {
576 irule->IntPoint(pt_id).z = irlist(2*pts_cnt + pt_id);
577 }
578 pt_id++;
579 }
580 }
581}
582
584 Vector &node_vals)
585{
586 const GridFunction *nodes = gf_in;
587 const FiniteElementSpace *fes = nodes->FESpace();
588 const int NE = mesh->GetNE();
589 const int vdim = fes->GetVDim();
590
591 IntegrationRule *ir_split_temp = NULL;
592
593 const int maxOrder = fes->GetMaxElementOrder();
594 const int dof_1D = maxOrder+1;
595 const int pts_el = std::pow(dof_1D, dim);
596 const int pts_cnt = NE_split_total * pts_el;
597 node_vals.SetSize(vdim * pts_cnt);
598 node_vals = 0.0;
599
600 int gsl_mesh_pt_index = 0;
601
602 for (int e = 0; e < NE; e++)
603 {
604 const FiniteElement *fe = fes->GetFE(e);
605 const Geometry::Type gt = fe->GetGeomType();
606 bool el_to_split = true;
607 if (gt == Geometry::TRIANGLE)
608 {
609 ir_split_temp = ir_split[0];
610 }
611 else if (gt == Geometry::TETRAHEDRON)
612 {
613 ir_split_temp = ir_split[1];
614 }
615 else if (gt == Geometry::PRISM)
616 {
617 ir_split_temp = ir_split[2];
618 }
619 else if (gt == Geometry::PYRAMID)
620 {
621 ir_split_temp = ir_split[3];
622 }
623 else if (gt == Geometry::SQUARE)
624 {
625 ir_split_temp = ir_split[1];
626 el_to_split = gf_in->FESpace()->IsVariableOrder();
627 }
628 else if (gt == Geometry::CUBE)
629 {
630 ir_split_temp = ir_split[0];
631 el_to_split = gf_in->FESpace()->IsVariableOrder();
632 }
633 else
634 {
635 MFEM_ABORT("Unsupported geometry type.");
636 }
637
638 if (el_to_split) // Triangle/Tet/Prism or Quads/Hex but variable order
639 {
640 // Fill gsl_mesh with location of split points.
641 Vector locval(vdim);
642 for (int i = 0; i < ir_split_temp->GetNPoints(); i++)
643 {
644 const IntegrationPoint &ip = ir_split_temp->IntPoint(i);
645 nodes->GetVectorValue(e, ip, locval);
646 for (int d = 0; d < vdim; d++)
647 {
648 node_vals(pts_cnt*d + gsl_mesh_pt_index) = locval(d);
649 }
650 gsl_mesh_pt_index++;
651 }
652 }
653 else // Quad/Hex and constant polynomial order
654 {
655 const int dof_cnt_split = fe->GetDof();
656
657 const TensorBasisElement *tbe =
658 dynamic_cast<const TensorBasisElement *>(fes->GetFE(e));
659 MFEM_VERIFY(tbe != NULL, "TensorBasis FiniteElement expected.");
660 Array<int> dof_map(dof_cnt_split);
661 const Array<int> &dm = tbe->GetDofMap();
662 if (dm.Size() > 0) { dof_map = dm; }
663 else { for (int i = 0; i < dof_cnt_split; i++) { dof_map[i] = i; } }
664
665 DenseMatrix pos(dof_cnt_split, vdim);
666 Vector posV(pos.Data(), dof_cnt_split * vdim);
667 Array<int> xdofs(dof_cnt_split * vdim);
668
669 fes->GetElementVDofs(e, xdofs);
670 nodes->GetSubVector(xdofs, posV);
671 for (int j = 0; j < dof_cnt_split; j++)
672 {
673 for (int d = 0; d < vdim; d++)
674 {
675 node_vals(pts_cnt * d + gsl_mesh_pt_index) = pos(dof_map[j], d);
676 }
677 gsl_mesh_pt_index++;
678 }
679 }
680 }
681}
682
683
685{
688
689 gsl_mfem_ref -= -1.; // map [-1, 1] to
690 gsl_mfem_ref *= 0.5; // [0, 1]
691
692 int nptorig = points_cnt,
693 npt = points_cnt;
694
695 // tolerance for point to be marked as on element edge/face
696 double btol = 1e-12;
697
698 GridFunction *gf_rst_map_temp = NULL;
699 int nptsend = 0;
700
701 for (int index = 0; index < npt; index++)
702 {
703 if (gsl_code[index] != 2 && gsl_proc[index] != gsl_comm->id)
704 {
705 nptsend +=1;
706 }
707 }
708
709 // Pack data to send via crystal router
710 struct gslib::array *outpt = new gslib::array;
711 struct out_pt { double r[3]; uint index, el, proc, code; };
712 struct out_pt *pt;
713 array_init(struct out_pt, outpt, nptsend);
714 outpt->n=nptsend;
715 pt = (struct out_pt *)outpt->ptr;
716 for (int index = 0; index < npt; index++)
717 {
718 if (gsl_code[index] == 2 || gsl_proc[index] == gsl_comm->id)
719 {
720 continue;
721 }
722 for (int d = 0; d < dim; ++d)
723 {
724 pt->r[d]= gsl_mfem_ref(index*dim + d);
725 }
726 pt->index = index;
727 pt->proc = gsl_proc[index];
728 pt->el = gsl_elem[index];
729 pt->code = gsl_code[index];
730 ++pt;
731 }
732
733 // Transfer data to target MPI ranks
734 sarray_transfer(struct out_pt, outpt, proc, 1, cr);
735 // Map received points
736 npt = outpt->n;
737 pt = (struct out_pt *)outpt->ptr;
738 for (int index = 0; index < npt; index++)
739 {
741 ip.Set3(&pt->r[0]);
742 const int elem = pt->el;
743 const int mesh_elem = split_element_map[elem];
744 const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(mesh_elem);
745 const Geometry::Type gt = fe->GetGeomType();
746 pt->el = mesh_elem;
747
748 if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
749 {
750 // check if it is on element boundary
751 pt->code = Geometry::CheckPoint(gt, ip, -btol) ? 0 : 1;
752 ++pt;
753 continue;
754 }
755 else if (gt == Geometry::TRIANGLE)
756 {
757 gf_rst_map_temp = gf_rst_map[0];
758 }
759 else if (gt == Geometry::TETRAHEDRON)
760 {
761 gf_rst_map_temp = gf_rst_map[1];
762 }
763 else if (gt == Geometry::PRISM)
764 {
765 gf_rst_map_temp = gf_rst_map[2];
766 }
767 else if (gt == Geometry::PYRAMID)
768 {
769 gf_rst_map_temp = gf_rst_map[3];
770 }
771
772 int local_elem = split_element_index[elem];
773 Vector mfem_ref(dim);
774 // map to rst of macro element
775 gf_rst_map_temp->GetVectorValue(local_elem, ip, mfem_ref);
776
777 for (int d = 0; d < dim; d++)
778 {
779 pt->r[d] = mfem_ref(d);
780 }
781
782 // check if point is on element boundary
783 ip.Set3(&pt->r[0]);
784 pt->code = Geometry::CheckPoint(gt, ip, -btol) ? 0 : 1;
785 ++pt;
786 }
787
788 // Transfer data back to source MPI rank
789 sarray_transfer(struct out_pt, outpt, proc, 1, cr);
790 npt = outpt->n;
791
792 // First copy mapped information for points on other procs
793 pt = (struct out_pt *)outpt->ptr;
794 for (int index = 0; index < npt; index++)
795 {
796 gsl_mfem_elem[pt->index] = pt->el;
797 for (int d = 0; d < dim; d++)
798 {
799 gsl_mfem_ref(d + pt->index*dim) = pt->r[d];
800 }
801 gsl_code[pt->index] = pt->code;
802 ++pt;
803 }
804 array_free(outpt);
805 delete outpt;
806
807 // Now map information for points on the same proc
808 for (int index = 0; index < nptorig; index++)
809 {
810 if (gsl_code[index] != 2 && gsl_proc[index] == gsl_comm->id)
811 {
812
814 Vector mfem_ref(gsl_mfem_ref.GetData()+index*dim, dim);
815 ip.Set2(mfem_ref.GetData());
816 if (dim == 3) { ip.z = mfem_ref(2); }
817
818 const int elem = gsl_elem[index];
819 const int mesh_elem = split_element_map[elem];
820 const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(mesh_elem);
821 const Geometry::Type gt = fe->GetGeomType();
822 gsl_mfem_elem[index] = mesh_elem;
823 if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
824 {
825 gsl_code[index] = Geometry::CheckPoint(gt, ip, -btol) ? 0 : 1;
826 continue;
827 }
828 else if (gt == Geometry::TRIANGLE)
829 {
830 gf_rst_map_temp = gf_rst_map[0];
831 }
832 else if (gt == Geometry::TETRAHEDRON)
833 {
834 gf_rst_map_temp = gf_rst_map[1];
835 }
836 else if (gt == Geometry::PRISM)
837 {
838 gf_rst_map_temp = gf_rst_map[2];
839 }
840 else if (gt == Geometry::PYRAMID)
841 {
842 gf_rst_map_temp = gf_rst_map[3];
843 }
844
845 int local_elem = split_element_index[elem];
846 gf_rst_map_temp->GetVectorValue(local_elem, ip, mfem_ref);
847
848 // Check if the point is on element boundary
849 ip.Set2(mfem_ref.GetData());
850 if (dim == 3) { ip.z = mfem_ref(2); }
851 gsl_code[index] = Geometry::CheckPoint(gt, ip, -btol) ? 0 : 1;
852 }
853 }
854}
855
857 Vector &field_out)
858{
859 const int gf_order = field_in.FESpace()->GetMaxElementOrder(),
860 mesh_order = mesh->GetNodalFESpace()->GetMaxElementOrder();
861
862 const FiniteElementCollection *fec_in = field_in.FESpace()->FEColl();
863 const H1_FECollection *fec_h1 = dynamic_cast<const H1_FECollection *>(fec_in);
864 const L2_FECollection *fec_l2 = dynamic_cast<const L2_FECollection *>(fec_in);
865
866 if (fec_h1 && gf_order == mesh_order &&
868 field_in.FESpace()->IsVariableOrder() ==
870 {
871 InterpolateH1(field_in, field_out);
872 return;
873 }
874 else
875 {
876 InterpolateGeneral(field_in, field_out);
877 if (!fec_l2 || avgtype == AvgType::NONE) { return; }
878 }
879
880 // For points on element borders, project the L2 GridFunction to H1 and
881 // re-interpolate.
882 if (fec_l2)
883 {
884 Array<int> indl2;
885 for (int i = 0; i < points_cnt; i++)
886 {
887 if (gsl_code[i] == 1) { indl2.Append(i); }
888 }
889 int borderPts = indl2.Size();
890#ifdef MFEM_USE_MPI
891 MPI_Allreduce(MPI_IN_PLACE, &borderPts, 1, MPI_INT, MPI_SUM, gsl_comm->c);
892#endif
893 if (borderPts == 0) { return; } // no points on element borders
894
895 Vector field_out_l2(field_out.Size());
896 VectorGridFunctionCoefficient field_in_dg(&field_in);
897 int gf_order_h1 = std::max(gf_order, 1); // H1 should be at least order 1
898 H1_FECollection fec(gf_order_h1, dim);
899 const int ncomp = field_in.FESpace()->GetVDim();
900 FiniteElementSpace fes(mesh, &fec, ncomp);
901 GridFunction field_in_h1(&fes);
902
904 {
905 field_in_h1.ProjectDiscCoefficient(field_in_dg, GridFunction::ARITHMETIC);
906 }
907 else if (avgtype == AvgType::HARMONIC)
908 {
909 field_in_h1.ProjectDiscCoefficient(field_in_dg, GridFunction::HARMONIC);
910 }
911 else
912 {
913 MFEM_ABORT("Invalid averaging type.");
914 }
915
916 if (gf_order_h1 == mesh_order) // basis is GaussLobatto by default
917 {
918 InterpolateH1(field_in_h1, field_out_l2);
919 }
920 else
921 {
922 InterpolateGeneral(field_in_h1, field_out_l2);
923 }
924
925 // Copy interpolated values for the points on element border
926 for (int j = 0; j < ncomp; j++)
927 {
928 for (int i = 0; i < indl2.Size(); i++)
929 {
930 int idx = field_in.FESpace()->GetOrdering() == Ordering::byNODES ?
931 indl2[i] + j*points_cnt:
932 indl2[i]*ncomp + j;
933 field_out(idx) = field_out_l2(idx);
934 }
935 }
936 }
937}
938
940 Vector &field_out)
941{
942 FiniteElementSpace ind_fes(mesh, field_in.FESpace()->FEColl());
943 if (field_in.FESpace()->IsVariableOrder())
944 {
945 for (int e = 0; e < ind_fes.GetMesh()->GetNE(); e++)
946 {
947 ind_fes.SetElementOrder(e, field_in.FESpace()->GetElementOrder(e));
948 }
949 ind_fes.Update(false);
950 }
951 GridFunction field_in_scalar(&ind_fes);
952 Vector node_vals;
953
954 const int ncomp = field_in.FESpace()->GetVDim(),
955 points_fld = field_in.Size() / ncomp;
956 MFEM_VERIFY(points_cnt == gsl_code.Size(),
957 "FindPointsGSLIB::InterpolateH1: Inconsistent size of gsl_code");
958
959 field_out.SetSize(points_cnt*ncomp);
960 field_out = default_interp_value;
961
962 for (int i = 0; i < ncomp; i++)
963 {
964 const int dataptrin = i*points_fld,
965 dataptrout = i*points_cnt;
966 if (field_in.FESpace()->GetOrdering() == Ordering::byNODES)
967 {
968 field_in_scalar.NewDataAndSize(field_in.GetData()+dataptrin, points_fld);
969 }
970 else
971 {
972 for (int j = 0; j < points_fld; j++)
973 {
974 field_in_scalar(j) = field_in(i + j*ncomp);
975 }
976 }
977 GetNodalValues(&field_in_scalar, node_vals);
978
979 if (dim==2)
980 {
981 findpts_eval_2(field_out.GetData()+dataptrout, sizeof(double),
982 gsl_code.GetData(), sizeof(unsigned int),
983 gsl_proc.GetData(), sizeof(unsigned int),
984 gsl_elem.GetData(), sizeof(unsigned int),
985 gsl_ref.GetData(), sizeof(double) * dim,
986 points_cnt, node_vals.GetData(), fdata2D);
987 }
988 else
989 {
990 findpts_eval_3(field_out.GetData()+dataptrout, sizeof(double),
991 gsl_code.GetData(), sizeof(unsigned int),
992 gsl_proc.GetData(), sizeof(unsigned int),
993 gsl_elem.GetData(), sizeof(unsigned int),
994 gsl_ref.GetData(), sizeof(double) * dim,
995 points_cnt, node_vals.GetData(), fdata3D);
996 }
997 }
998 if (field_in.FESpace()->GetOrdering() == Ordering::byVDIM)
999 {
1000 Vector field_out_temp = field_out;
1001 for (int i = 0; i < ncomp; i++)
1002 {
1003 for (int j = 0; j < points_cnt; j++)
1004 {
1005 field_out(i + j*ncomp) = field_out_temp(j + i*points_cnt);
1006 }
1007 }
1008 }
1009}
1010
1012 Vector &field_out)
1013{
1014 int ncomp = field_in.VectorDim(),
1015 nptorig = points_cnt,
1016 npt = points_cnt;
1017
1018 field_out.SetSize(points_cnt*ncomp);
1019 field_out = default_interp_value;
1020
1021 if (gsl_comm->np == 1) // serial
1022 {
1023 for (int index = 0; index < npt; index++)
1024 {
1025 if (gsl_code[index] == 2) { continue; }
1028 if (dim == 3) { ip.z = gsl_mfem_ref(index*dim + 2); }
1029 Vector localval(ncomp);
1030 field_in.GetVectorValue(gsl_mfem_elem[index], ip, localval);
1031 if (field_in.FESpace()->GetOrdering() == Ordering::byNODES)
1032 {
1033 for (int i = 0; i < ncomp; i++)
1034 {
1035 field_out(index + i*npt) = localval(i);
1036 }
1037 }
1038 else //byVDIM
1039 {
1040 for (int i = 0; i < ncomp; i++)
1041 {
1042 field_out(index*ncomp + i) = localval(i);
1043 }
1044 }
1045 }
1046 }
1047 else // parallel
1048 {
1049 // Determine number of points to be sent
1050 int nptsend = 0;
1051 for (int index = 0; index < npt; index++)
1052 {
1053 if (gsl_code[index] != 2) { nptsend +=1; }
1054 }
1055
1056 // Pack data to send via crystal router
1057 struct gslib::array *outpt = new gslib::array;
1058 struct out_pt { double r[3], ival; uint index, el, proc; };
1059 struct out_pt *pt;
1060 array_init(struct out_pt, outpt, nptsend);
1061 outpt->n=nptsend;
1062 pt = (struct out_pt *)outpt->ptr;
1063 for (int index = 0; index < npt; index++)
1064 {
1065 if (gsl_code[index] == 2) { continue; }
1066 for (int d = 0; d < dim; ++d) { pt->r[d]= gsl_mfem_ref(index*dim + d); }
1067 pt->index = index;
1068 pt->proc = gsl_proc[index];
1069 pt->el = gsl_mfem_elem[index];
1070 ++pt;
1071 }
1072
1073 // Transfer data to target MPI ranks
1074 sarray_transfer(struct out_pt, outpt, proc, 1, cr);
1075
1076 if (ncomp == 1)
1077 {
1078 // Interpolate the grid function
1079 npt = outpt->n;
1080 pt = (struct out_pt *)outpt->ptr;
1081 for (int index = 0; index < npt; index++)
1082 {
1084 ip.Set3(&pt->r[0]);
1085 pt->ival = field_in.GetValue(pt->el, ip, 1);
1086 ++pt;
1087 }
1088
1089 // Transfer data back to source MPI rank
1090 sarray_transfer(struct out_pt, outpt, proc, 1, cr);
1091 npt = outpt->n;
1092 pt = (struct out_pt *)outpt->ptr;
1093 for (int index = 0; index < npt; index++)
1094 {
1095 field_out(pt->index) = pt->ival;
1096 ++pt;
1097 }
1098 array_free(outpt);
1099 delete outpt;
1100 }
1101 else // ncomp > 1
1102 {
1103 // Interpolate data and store in a Vector
1104 npt = outpt->n;
1105 pt = (struct out_pt *)outpt->ptr;
1106 Vector vec_int_vals(npt*ncomp);
1107 for (int index = 0; index < npt; index++)
1108 {
1110 ip.Set3(&pt->r[0]);
1111 Vector localval(vec_int_vals.GetData()+index*ncomp, ncomp);
1112 field_in.GetVectorValue(pt->el, ip, localval);
1113 ++pt;
1114 }
1115
1116 // Save index and proc data in a struct
1117 struct gslib::array *savpt = new gslib::array;
1118 struct sav_pt { uint index, proc; };
1119 struct sav_pt *spt;
1120 array_init(struct sav_pt, savpt, npt);
1121 savpt->n=npt;
1122 spt = (struct sav_pt *)savpt->ptr;
1123 pt = (struct out_pt *)outpt->ptr;
1124 for (int index = 0; index < npt; index++)
1125 {
1126 spt->index = pt->index;
1127 spt->proc = pt->proc;
1128 ++pt; ++spt;
1129 }
1130
1131 array_free(outpt);
1132 delete outpt;
1133
1134 // Copy data from save struct to send struct and send component wise
1135 struct gslib::array *sendpt = new gslib::array;
1136 struct send_pt { double ival; uint index, proc; };
1137 struct send_pt *sdpt;
1138 for (int j = 0; j < ncomp; j++)
1139 {
1140 array_init(struct send_pt, sendpt, npt);
1141 sendpt->n=npt;
1142 spt = (struct sav_pt *)savpt->ptr;
1143 sdpt = (struct send_pt *)sendpt->ptr;
1144 for (int index = 0; index < npt; index++)
1145 {
1146 sdpt->index = spt->index;
1147 sdpt->proc = spt->proc;
1148 sdpt->ival = vec_int_vals(j + index*ncomp);
1149 ++sdpt; ++spt;
1150 }
1151
1152 sarray_transfer(struct send_pt, sendpt, proc, 1, cr);
1153 sdpt = (struct send_pt *)sendpt->ptr;
1154 for (int index = 0; index < static_cast<int>(sendpt->n); index++)
1155 {
1156 int idx = field_in.FESpace()->GetOrdering() == Ordering::byNODES ?
1157 sdpt->index + j*nptorig :
1158 sdpt->index*ncomp + j;
1159 field_out(idx) = sdpt->ival;
1160 ++sdpt;
1161 }
1162 array_free(sendpt);
1163 }
1164 array_free(savpt);
1165 delete sendpt;
1166 delete savpt;
1167 } // ncomp > 1
1168 } // parallel
1169}
1170
1171void OversetFindPointsGSLIB::Setup(Mesh &m, const int meshid,
1172 GridFunction *gfmax,
1173 const double bb_t, const double newt_tol,
1174 const int npt_max)
1175{
1176 MFEM_VERIFY(m.GetNodes() != NULL, "Mesh nodes are required.");
1177 const int meshOrder = m.GetNodes()->FESpace()->GetMaxElementOrder();
1178
1179 // FreeData if OversetFindPointsGSLIB::Setup has been called already
1180 if (setupflag) { FreeData(); }
1181
1182 crystal_init(cr, gsl_comm);
1183 mesh = &m;
1184 dim = mesh->Dimension();
1185 const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
1186 unsigned dof1D = fe->GetOrder() + 1;
1187
1189 if (dim == 2)
1190 {
1191 if (ir_split[0]) { delete ir_split[0]; ir_split[0] = NULL; }
1192 ir_split[0] = new IntegrationRule(3*pow(dof1D, dim));
1194
1195 if (ir_split[1]) { delete ir_split[1]; ir_split[1] = NULL; }
1196 ir_split[1] = new IntegrationRule(pow(dof1D, dim));
1198 }
1199 else if (dim == 3)
1200 {
1201 if (ir_split[0]) { delete ir_split[0]; ir_split[0] = NULL; }
1202 ir_split[0] = new IntegrationRule(pow(dof1D, dim));
1204
1205 if (ir_split[1]) { delete ir_split[1]; ir_split[1] = NULL; }
1206 ir_split[1] = new IntegrationRule(4*pow(dof1D, dim));
1208
1209 if (ir_split[2]) { delete ir_split[2]; ir_split[2] = NULL; }
1210 ir_split[2] = new IntegrationRule(3*pow(dof1D, dim));
1212
1213 if (ir_split[3]) { delete ir_split[3]; ir_split[3] = NULL; }
1214 ir_split[3] = new IntegrationRule(8*pow(dof1D, dim));
1216 }
1217
1219
1220 MFEM_ASSERT(meshid>=0, " The ID should be greater than or equal to 0.");
1221
1222 const int pts_cnt = gsl_mesh.Size()/dim,
1223 NEtot = pts_cnt/(int)pow(dof1D, dim);
1224
1225 distfint.SetSize(pts_cnt);
1226 if (!gfmax)
1227 {
1228 distfint = 0.0;
1229 }
1230 else
1231 {
1232 GetNodalValues(gfmax, distfint);
1233 }
1234 u_meshid = (unsigned int)meshid;
1235
1236 if (dim == 2)
1237 {
1238 unsigned nr[2] = { dof1D, dof1D };
1239 unsigned mr[2] = { 2*dof1D, 2*dof1D };
1240 double * const elx[2] =
1241 {
1242 pts_cnt == 0 ? nullptr : &gsl_mesh(0),
1243 pts_cnt == 0 ? nullptr : &gsl_mesh(pts_cnt)
1244 };
1245 fdata2D = findptsms_setup_2(gsl_comm, elx, nr, NEtot, mr, bb_t,
1246 pts_cnt, pts_cnt, npt_max, newt_tol,
1247 &u_meshid, &distfint(0));
1248 }
1249 else
1250 {
1251 unsigned nr[3] = { dof1D, dof1D, dof1D };
1252 unsigned mr[3] = { 2*dof1D, 2*dof1D, 2*dof1D };
1253 double * const elx[3] =
1254 {
1255 pts_cnt == 0 ? nullptr : &gsl_mesh(0),
1256 pts_cnt == 0 ? nullptr : &gsl_mesh(pts_cnt),
1257 pts_cnt == 0 ? nullptr : &gsl_mesh(2*pts_cnt)
1258 };
1259 fdata3D = findptsms_setup_3(gsl_comm, elx, nr, NEtot, mr, bb_t,
1260 pts_cnt, pts_cnt, npt_max, newt_tol,
1261 &u_meshid, &distfint(0));
1262 }
1263 setupflag = true;
1264 overset = true;
1265}
1266
1268 Array<unsigned int> &point_id,
1269 int point_pos_ordering)
1270{
1271 MFEM_VERIFY(setupflag, "Use OversetFindPointsGSLIB::Setup before "
1272 "finding points.");
1273 MFEM_VERIFY(overset, " Please setup FindPoints for overlapping grids.");
1274 points_cnt = point_pos.Size() / dim;
1275 unsigned int match = 0; // Don't find points in the mesh if point_id = mesh_id
1276
1282
1283 auto xvFill = [&](const double *xv_base[], unsigned xv_stride[])
1284 {
1285 for (int d = 0; d < dim; d++)
1286 {
1287 if (point_pos_ordering == Ordering::byNODES)
1288 {
1289 xv_base[d] = point_pos.GetData() + d*points_cnt;
1290 xv_stride[d] = sizeof(double);
1291 }
1292 else
1293 {
1294 xv_base[d] = point_pos.GetData() + d;
1295 xv_stride[d] = dim*sizeof(double);
1296 }
1297 }
1298 };
1299 if (dim == 2)
1300 {
1301 const double *xv_base[2];
1302 unsigned xv_stride[2];
1303 xvFill(xv_base, xv_stride);
1304 findptsms_2(gsl_code.GetData(), sizeof(unsigned int),
1305 gsl_proc.GetData(), sizeof(unsigned int),
1306 gsl_elem.GetData(), sizeof(unsigned int),
1307 gsl_ref.GetData(), sizeof(double) * dim,
1308 gsl_dist.GetData(), sizeof(double),
1309 xv_base, xv_stride,
1310 point_id.GetData(), sizeof(unsigned int), &match,
1312 }
1313 else // dim == 3
1314 {
1315 const double *xv_base[3];
1316 unsigned xv_stride[3];
1317 xvFill(xv_base, xv_stride);
1318 findptsms_3(gsl_code.GetData(), sizeof(unsigned int),
1319 gsl_proc.GetData(), sizeof(unsigned int),
1320 gsl_elem.GetData(), sizeof(unsigned int),
1321 gsl_ref.GetData(), sizeof(double) * dim,
1322 gsl_dist.GetData(), sizeof(double),
1323 xv_base, xv_stride,
1324 point_id.GetData(), sizeof(unsigned int), &match,
1326 }
1327
1328 // Set the element number and reference position to 0 for points not found
1329 for (int i = 0; i < points_cnt; i++)
1330 {
1331 if (gsl_code[i] == 2 ||
1332 (gsl_code[i] == 1 && gsl_dist(i) > bdr_tol))
1333 {
1334 gsl_elem[i] = 0;
1335 for (int d = 0; d < dim; d++) { gsl_ref(i*dim + d) = -1.; }
1336 gsl_code[i] = 2;
1337 }
1338 }
1339
1340 // Map element number for simplices, and ref_pos from [-1,1] to [0,1] for both
1341 // simplices and quads.
1343}
1344
1346 Array<unsigned int> &point_id,
1347 const GridFunction &field_in,
1348 Vector &field_out,
1349 int point_pos_ordering)
1350{
1351 FindPoints(point_pos, point_id, point_pos_ordering);
1352 Interpolate(field_in, field_out);
1353}
1354
1355
1356} // namespace mfem
1357
1358#endif // MFEM_USE_GSLIB
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
void DeleteAll()
Delete the whole array.
Definition array.hpp:864
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
T * GetData()
Returns the data.
Definition array.hpp:118
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
real_t * Data() const
Returns the matrix data array.
Definition densemat.hpp:111
Geometry::Type GetGeometryType() const
Definition element.hpp:52
virtual ~FindPointsGSLIB()
Definition gslib.cpp:68
struct gslib::findpts_data_2 * fdata2D
Definition gslib.hpp:67
Array< unsigned int > gsl_code
Definition gslib.hpp:72
Array< Mesh * > mesh_split
Definition gslib.hpp:59
void FindPoints(const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
Definition gslib.cpp:186
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition gslib.cpp:109
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition gslib.cpp:856
virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out)
Use GSLIB for communication and interpolation.
Definition gslib.cpp:939
Array< FiniteElementSpace * > fes_rst_map
Definition gslib.hpp:64
virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals)
Get GridFunction value at the points expected by GSLIB.
Definition gslib.cpp:583
Array< int > split_element_index
Definition gslib.hpp:78
FiniteElementCollection * fec_map_lin
Definition gslib.hpp:66
struct gslib::comm * gsl_comm
Definition gslib.hpp:70
Array< unsigned int > gsl_elem
Definition gslib.hpp:72
Array< unsigned int > gsl_proc
Definition gslib.hpp:72
virtual void SetupIntegrationRuleForSplitMesh(Mesh *mesh, IntegrationRule *irule, int order)
Definition gslib.cpp:535
Array< unsigned int > gsl_mfem_elem
Definition gslib.hpp:72
Array< IntegrationRule * > ir_split
Definition gslib.hpp:62
double default_interp_value
Definition gslib.hpp:75
virtual void SetupSplitMeshes()
Definition gslib.cpp:312
virtual void MapRefPosAndElemIndices()
Definition gslib.cpp:684
virtual void FreeData()
Definition gslib.cpp:283
virtual void InterpolateGeneral(const GridFunction &field_in, Vector &field_out)
Definition gslib.cpp:1011
Array< GridFunction * > gf_rst_map
Definition gslib.hpp:65
Array< int > split_element_map
Definition gslib.hpp:77
struct gslib::findpts_data_3 * fdata3D
Definition gslib.hpp:68
struct gslib::crystal * cr
Definition gslib.hpp:69
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:581
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:725
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition fespace.cpp:3439
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
Definition fespace.cpp:176
void SetElementOrder(int i, int p)
Sets the order of the i'th finite element.
Definition fespace.cpp:149
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:727
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
int GetMaxElementOrder() const
Return the maximum polynomial order.
Definition fespace.hpp:577
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:333
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:326
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
Definition geom.cpp:435
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition gridfunc.cpp:449
FiniteElementSpace * FESpace()
Definition gridfunc.hpp:696
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
int VectorDim() const
Definition gridfunc.cpp:323
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition gridfunc.cpp:476
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
int GetBasisType() const
Definition fe_coll.hpp:285
Class for integration point with weight.
Definition intrules.hpp:35
void Set2(const real_t x1, const real_t x2)
Definition intrules.hpp:89
void Set3(const real_t x1, const real_t x2, const real_t x3)
Definition intrules.hpp:79
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Mesh data type.
Definition mesh.hpp:56
const FiniteElementSpace * GetNodalFESpace() const
Definition mesh.cpp:6206
const Element * GetElement(int i) const
Return pointer to the i'th element object.
Definition mesh.hpp:1283
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
static Mesh MakeCartesian3D(int nx, int ny, int nz, Element::Type type, real_t sx=1.0, real_t sy=1.0, real_t sz=1.0, bool sfc_ordering=true)
Creates a mesh for the parallelepiped [0,sx]x[0,sy]x[0,sz], divided into nx*ny*nz hexahedra if type =...
Definition mesh.cpp:4253
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
Definition mesh.cpp:6153
static Mesh MakeCartesian2D(int nx, int ny, Element::Type type, bool generate_edges=false, real_t sx=1.0, real_t sy=1.0, bool sfc_ordering=true)
Creates mesh for the rectangle [0,sx]x[0,sy], divided into nx*ny quadrilaterals if type = QUADRILATER...
Definition mesh.cpp:4243
void Setup(Mesh &m, const int meshid, GridFunction *gfmax=NULL, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition gslib.cpp:1171
void Interpolate(const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out, int point_pos_ordering=Ordering::byNODES)
Definition gslib.cpp:1345
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id, int point_pos_ordering=Ordering::byNODES)
Definition gslib.cpp:1267
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition fe_base.hpp:1222
Vector coefficient defined by a vector GridFunction.
Vector data type.
Definition vector.hpp:80
void Destroy()
Destroy a vector.
Definition vector.hpp:615
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
void NewDataAndSize(real_t *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition vector.hpp:181
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:578
int dim
Definition ex24.cpp:53
int index(int i, int j, int nx, int ny)
Definition life.cpp:235