DOLFIN
DOLFIN C++ interface
SubDomain.h
1 // Copyright (C) 2007-2013 Anders Logg
2 //
3 // This file is part of DOLFIN.
4 //
5 // DOLFIN is free software: you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published by
7 // the Free Software Foundation, either version 3 of the License, or
8 // (at your option) any later version.
9 //
10 // DOLFIN is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17 //
18 // First added: 2007-04-10
19 // Last changed: 2013-04-12
20 
21 #ifndef __SUB_DOMAIN_H
22 #define __SUB_DOMAIN_H
23 
24 #include <cstddef>
25 #include <map>
26 #include <dolfin/common/constants.h>
27 #include <Eigen/Dense>
28 
29 namespace dolfin
30 {
31 
32  // Forward declarations
33  class Mesh;
34  template <typename T> class MeshFunction;
35  template <typename T> class MeshValueCollection;
36  template <typename T> class Array;
37 
41 
42  class SubDomain
43  {
44  public:
45 
51  SubDomain(const double map_tol=1.0e-10);
52 
54  virtual ~SubDomain();
55 
65  virtual bool inside(const Array<double>& x, bool on_boundary) const;
66 
76  virtual bool inside(Eigen::Ref<const Eigen::VectorXd> x, bool on_boundary) const;
77 
85  virtual void map(const Array<double>& x, Array<double>& y) const;
86 
87 
95  virtual void map(Eigen::Ref<const Eigen::VectorXd> x, Eigen::Ref<Eigen::VectorXd> y) const;
96 
97 
102  virtual void snap(Array<double>& x) const {}
103 
104  //--- Marking of Mesh ---
105 
114  void mark_cells(Mesh& mesh,
115  std::size_t sub_domain,
116  bool check_midpoint=true) const;
117 
126  void mark_facets(Mesh& mesh,
127  std::size_t sub_domain,
128  bool check_midpoint=true) const;
129 
141  void mark(Mesh& mesh,
142  std::size_t dim,
143  std::size_t sub_domain,
144  bool check_midpoint=true) const;
145 
146  //--- Marking of MeshFunction ---
147 
156  void mark(MeshFunction<std::size_t>& sub_domains,
157  std::size_t sub_domain,
158  bool check_midpoint=true) const;
159 
168  void mark(MeshFunction<int>& sub_domains,
169  int sub_domain,
170  bool check_midpoint=true) const;
171 
180  void mark(MeshFunction<double>& sub_domains,
181  double sub_domain,
182  bool check_midpoint=true) const;
183 
192  void mark(MeshFunction<bool>& sub_domains,
193  bool sub_domain,
194  bool check_midpoint=true) const;
195 
196  //--- Marking of MeshValueCollection ---
197 
208  void mark(MeshValueCollection<std::size_t>& sub_domains,
209  std::size_t sub_domain,
210  const Mesh& mesh,
211  bool check_midpoint=true) const;
212 
223  void mark(MeshValueCollection<int>& sub_domains,
224  int sub_domain,
225  const Mesh& mesh,
226  bool check_midpoint=true) const;
227 
238  void mark(MeshValueCollection<double>& sub_domains,
239  double sub_domain,
240  const Mesh& mesh,
241  bool check_midpoint=true) const;
242 
253  void mark(MeshValueCollection<bool>& sub_domains,
254  bool sub_domain,
255  const Mesh& mesh,
256  bool check_midpoint=true) const;
257 
262  std::size_t geometric_dimension() const;
263 
268  virtual void set_property(std::string name, double value);
269 
274  virtual double get_property(std::string name) const;
275 
280  const double map_tolerance;
281 
282  private:
283 
286  template<typename S, typename T>
287  void apply_markers(S& sub_domains,
288  T sub_domain,
289  const Mesh& mesh,
290  bool check_midpoint) const;
291 
292  template<typename T>
293  void apply_markers(std::map<std::size_t, std::size_t>& sub_domains,
294  std::size_t dim,
295  T sub_domain,
296  const Mesh& mesh,
297  bool check_midpoint) const;
298 
299  // Friends
300  friend class DirichletBC;
301  friend class PeriodicBC;
302 
303  // Geometric dimension, needed for SWIG interface, will be set before
304  // calls to inside() and map()
305  mutable std::size_t _geometric_dimension;
306 
307  };
308 
309 }
310 
311 #endif
virtual ~SubDomain()
Destructor.
Definition: SubDomain.cpp:46
Definition: SubDomain.h:42
void mark(Mesh &mesh, std::size_t dim, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:94
Definition: adapt.h:29
virtual bool inside(const Array< double > &x, bool on_boundary) const
Definition: SubDomain.cpp:51
Definition: Array.h:41
void mark_facets(Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:87
const double map_tolerance
Definition: SubDomain.h:280
virtual void map(const Array< double > &x, Array< double > &y) const
Definition: SubDomain.cpp:65
SubDomain(const double map_tol=1.0e-10)
Definition: SubDomain.cpp:40
virtual void snap(Array< double > &x) const
Definition: SubDomain.h:102
Interface for setting (strong) Dirichlet boundary conditions.
Definition: DirichletBC.h:124
Definition: GenericFile.h:38
virtual void set_property(std::string name, double value)
Definition: SubDomain.cpp:381
std::size_t geometric_dimension() const
Definition: SubDomain.cpp:165
void mark_cells(Mesh &mesh, std::size_t sub_domain, bool check_midpoint=true) const
Definition: SubDomain.cpp:80
virtual double get_property(std::string name) const
Definition: SubDomain.cpp:388
Definition: Mesh.h:82