Flow123d
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
flow123d
src
mesh
accessors.hh
Go to the documentation of this file.
1
/*
2
* Accessors.hh
3
*
4
* Created on: Dec 4, 2012
5
* Author: jb
6
*/
7
8
#ifndef ACCESSORS_HH_
9
#define ACCESSORS_HH_
10
11
#include "
mesh/bounding_box.hh
"
12
#include "
mesh/mesh_types.hh
"
13
#include "
mesh/region.hh
"
14
#include "
mesh/elements.h
"
15
#include "
mesh/mesh.h
"
16
#include <armadillo>
17
18
/**
19
* Element accessor templated just by dimension of the embedding space, used by Fields.
20
* This should allow algorithms over elements where dimension of particular element is runtime parameter.
21
*
22
* This class suites as interface of Fields to the mesh elements, in particular this accessor knows directly
23
* the region, and also can be used as an accessor that works on the whole region if used by Fields that do not depend on
24
* particular elements as FieldConstant, FiledFormula, and FieldPython.
25
*
26
* TODO:
27
* - make this kind of accessor subclass of FieldCommonBase or at least move it into src/fields
28
* since it has functionality particular for Fields
29
*
30
* Ideas:
31
* need function to calculate intersection (object) of two ElementAccessors, but this definitely should be templated by
32
* dimension of the ref. element (or rather shape of ref. element), here we can have case dispatch
33
*
34
*/
35
template
<
int
spacedim>
36
class
ElementAccessor
{
37
public
:
38
/**
39
* Default invalid accessor.
40
*/
41
ElementAccessor
()
42
:
mesh_
(NULL)
43
{}
44
45
/**
46
* Regional accessor.
47
*/
48
ElementAccessor
(
const
Mesh
*mesh,
RegionIdx
r_idx)
49
:
mesh_
(mesh),
dim_
(
undefined_dim_
),
r_idx_
(r_idx)
50
{}
51
52
/**
53
* Element accessor.
54
*/
55
ElementAccessor
(
const
Mesh
*mesh,
unsigned
int
idx
,
bool
boundary)
56
:
mesh_
(mesh),
boundary_
(boundary),
element_idx_
(idx),
r_idx_
(
element
()->
region_idx
())
57
{
58
dim_
=
element
()->
dim
();
59
}
60
61
inline
bool
is_regional
()
const
{
62
return
dim_
==
undefined_dim_
;
63
}
64
65
inline
bool
is_elemental
()
const
{
66
return
(
is_valid
() && !
is_regional
() );
67
}
68
69
inline
bool
is_valid
()
const
{
70
return
mesh_
!= NULL;
71
}
72
73
inline
unsigned
int
dim
()
const
74
{
return
dim_
; }
75
76
inline
const
Element
*
element
()
const
{
77
if
(
boundary_
)
return
(
Element
*)(
mesh_
->
bc_elements
(
element_idx_
)) ;
78
else
return
(
Element
*)(
mesh_
->
element
(
element_idx_
)) ;
79
}
80
81
inline
arma::vec::fixed<spacedim>
centre
()
const
{
82
ASSERT
(
is_valid
(),
"Invalid element accessor."
);
83
if
(
is_regional
() )
return
arma::vec::fixed<spacedim>();
84
else
return
element
()->
centre
();
85
}
86
87
88
inline
Region
region
()
const
89
{
return
Region
(
r_idx_
,
mesh_
->
region_db
()); }
90
91
inline
RegionIdx
region_idx
()
const
92
{
return
r_idx_
; }
93
94
/// We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular mesh
95
//inline unsigned int region_id() const {
96
// return region().id();
97
//}
98
99
//const BoundingBox &bounding_box();
100
101
inline
bool
is_boundary
()
const
{
102
return
boundary_
;
103
}
104
105
inline
unsigned
int
idx
()
const
{
106
return
element_idx_
;
107
}
108
private
:
109
/**
110
* When dim_ == undefined_dim_ ; the value of element_idx_ is invalid.
111
* Is used for ElementAccessors for whole region
112
*/
113
static
const
unsigned
int
undefined_dim_
= 100;
114
115
/// Dimension of reference element.
116
unsigned
int
dim_
;
117
118
/// Pointer to the mesh owning the element.
119
const
Mesh
*
mesh_
;
120
/// True if the element is boundary, i.e. stored in Mesh::bc_elements, bulk elements are stored in Mesh::element
121
bool
boundary_
;
122
123
/// Index into Mesh::bc_elements or Mesh::element array.
124
unsigned
int
element_idx_
;
125
126
/// Region index.
127
RegionIdx
r_idx_
;
128
};
129
130
131
132
133
/******************************************************************* implementations
134
*
135
*
136
*/
137
/*
138
template<int spacedim>
139
const BoundingBox &ElementAccessor<spacedim>::bounding_box() {
140
return box_;
141
}
142
*/
143
144
#endif
/* ACCESSORS_HH_ */
Generated on Thu May 29 2014 23:14:48 for Flow123d by
1.8.4