MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
kdtree.hpp
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#ifndef MFEM_KDTREE_PROJECTION
13#define MFEM_KDTREE_PROJECTION
14
15#include "../general/kdtree.hpp"
16#include "gridfunc.hpp"
17
18namespace mfem
19{
20
21/// Base class for KDTreeNodalProjection.
23{
24public:
27
28 /// The projection method can be called as many time as necessary with
29 /// different sets of coordinates and corresponding values. For vector
30 /// grid function, users have to specify the data ordering and for all
31 /// cases the user can modify the error tolerance err to smaller or
32 /// bigger value. A node in the target grid function is matching
33 /// a point with coordinates specified in the vector coords if the
34 /// distance between them is smaller than lerr.
35 virtual
36 void Project(const Vector& coords,const Vector& src,
37 int ordering, real_t lerr) = 0;
38
39 /// The project method can be called as many times as necessary with
40 /// different grid functions gf. A node in the target grid function is
41 /// matching a node from the source grid function if the distance
42 /// between them is smaller than lerr.
43 virtual
44 void Project(const GridFunction& gf, real_t lerr) = 0;
45};
46
47/// The class provides methods for projecting function values evaluated on a
48/// set of points to a grid function. The values are directly copied to the
49/// nodal values of the target grid function if any of the points is matching
50/// a node of the grid function. For example, if a parallel grid function is
51/// saved in parallel, every saved chunk can be read on every other process
52/// and mapped to a local grid function that does not have the same structure
53/// as the original one. The functionality is based on a kd-tree search in a
54/// cloud of points.
55template<int kdim=3>
57{
58private:
59 /// Pointer to the KDTree
60 std::unique_ptr<KDTree<int,real_t,kdim>> kdt;
61
62 /// Pointer to the target grid function
63 GridFunction* dest;
64
65 /// Upper corner of the bounding box
66 Vector maxbb;
67
68 /// Lower corner of the bounding box
69 Vector minbb;
70
71public:
72 /// The constructor takes as input an L2 or H1 grid function (it can be
73 /// a vector grid function). The Project method copies a set of values
74 /// to the grid function.
76 {
77 dest=&dest_;
79
80 MFEM_VERIFY(
81 dynamic_cast<const H1_FECollection*>(space->FEColl()) != nullptr ||
82 dynamic_cast<const L2_FECollection*>(space->FEColl()) != nullptr,
83 "Error!");
84
85 Mesh* mesh=space->GetMesh();
86
87 const int dim=mesh->SpaceDimension();
88 MFEM_VERIFY(kdim==dim, "GridFunction dimension does not match!");
89
90 kdt=std::unique_ptr<KDTree<int,real_t,kdim>>(
92
93 std::vector<bool> indt;
94 indt.resize(space->GetVSize()/space->GetVDim(), true);
95
96 minbb.SetSize(dim);
97 maxbb.SetSize(dim);
98
99 //set the loocal coordinates
100 {
102 const IntegrationRule* ir=nullptr;
103 Array<int> vdofs;
104 DenseMatrix elco;
105 int isca=1;
106 if (space->GetOrdering()==Ordering::byVDIM)
107 {
108 isca=space->GetVDim();
109 }
110
111 // intialize the bounding box
112 const FiniteElement* el=space->GetFE(0);
113 trans = space->GetElementTransformation(0);
114 ir=&(el->GetNodes());
115 space->GetElementVDofs(0,vdofs);
116 elco.SetSize(dim,ir->GetNPoints());
117 trans->Transform(*ir,elco);
118 for (int d=0; d<dim; d++)
119 {
120 minbb[d]=elco(d,0);
121 maxbb[d]=elco(d,0);
122 }
123
124 for (int i=0; i<space->GetNE(); i++)
125 {
126 el=space->GetFE(i);
127 // get the element transformation
128 trans = space->GetElementTransformation(i);
129 ir=&(el->GetNodes());
130 space->GetElementVDofs(i,vdofs);
131 elco.SetSize(dim,ir->GetNPoints());
132 trans->Transform(*ir,elco);
133
134 for (int p=0; p<ir->GetNPoints(); p++)
135 {
136 int bind=vdofs[p]/isca;
137 if (indt[bind]==true)
138 {
139 kdt->AddPoint(elco.GetColumn(p),bind);
140 indt[bind]=false;
141
142 for (int d=0; d<kdim; d++)
143 {
144 if (minbb[d]>elco(d,p)) {minbb[d]=elco(d,p);}
145 if (maxbb[d]<elco(d,p)) {maxbb[d]=elco(d,p);}
146 }
147 }
148 }
149 }
150 }
151
152 // build the KDTree
153 kdt->Sort();
154 }
155
156 /// The projection method can be called as many time as necessary with
157 /// different sets of coordinates and corresponding values. For vector
158 /// grid function, users have to specify the data ordering and for all
159 /// cases the user can modify the error tolerance err to smaller or
160 /// bigger value. A node in the target grid function is matching
161 /// a point with coordinates specified in the vector coords if the
162 /// distance between them is smaller than lerr.
163 virtual
164 void Project(const Vector& coords,const Vector& src,
165 int ordering=Ordering::byNODES, real_t lerr=1e-8);
166
167 /// The project method can be called as many times as necessary with
168 /// different grid functions gf. A node in the target grid function is
169 /// matching a node from the source grid function if the distance
170 /// between them is smaller than lerr.
171 virtual
172 void Project(const GridFunction& gf, real_t lerr=1e-8);
173};
174
175} // namespace mfem
176
177#endif // MFEM_KDTREE_PROJECTION
Base class for KDTreeNodalProjection.
Definition kdtree.hpp:23
virtual void Project(const GridFunction &gf, real_t lerr)=0
virtual void Project(const Vector &coords, const Vector &src, int ordering, real_t lerr)=0
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
void GetColumn(int c, Vector &col) const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
Abstract class for all finite elements.
Definition fe_base.hpp:239
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:395
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
FiniteElementSpace * FESpace()
Definition gridfunc.hpp:696
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
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
virtual void Project(const GridFunction &gf, real_t lerr=1e-8)
KDTreeNodalProjection(GridFunction &dest_)
Definition kdtree.hpp:75
virtual void Project(const Vector &coords, const Vector &src, int ordering=Ordering::byNODES, real_t lerr=1e-8)
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Mesh data type.
Definition mesh.hpp:56
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
Vector data type.
Definition vector.hpp:80
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
int dim
Definition ex24.cpp:53
void trans(const Vector &u, Vector &x)
Definition ex27.cpp:412
string space
float real_t
Definition config.hpp:43
real_t p(const Vector &x, real_t t)