/usr/include/dolfin/function/CCFEMFunctionSpace.h is in libdolfin-dev 1.4.0+dfsg-4.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
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 | // Copyright (C) 2013 Anders Logg
//
// This file is part of DOLFIN.
//
// DOLFIN 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.
//
// DOLFIN 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 for more details.
//
// You should have received a copy of the GNU Lesser General Public License
// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
//
// First added: 2013-08-05
// Last changed: 2014-03-03
#ifndef __CCFEM_FUNCTION_SPACE_H
#define __CCFEM_FUNCTION_SPACE_H
#include <vector>
#include <map>
#include <memory>
namespace dolfin
{
// Forward declarations
class FunctionSpace;
class CCFEMDofMap;
class Mesh;
class BoundingBoxTree;
class BoundaryMesh;
// FIXME: Consider moving many of the data structures from this
// class to a new class named CCFEMMesh (or similar), since many of
// them are related to meshes only and not a particular function
// space.
// FIXME: Consider renaming this class and related classes from
// CCFEM to something else, perhaps MultiMeshFunctionSpace?
/// This class represents a cut and composite finite element
/// function space (CCFEM) defined on one or more possibly
/// intersecting meshes.
///
/// A CCFEM function space may be created from a set of standard
/// function spaces by repeatedly calling add(), followed by a call
/// to build(). Note that a CCFEM function space is not useful and
/// its data structures are empty until build() has been called.
class CCFEMFunctionSpace
{
public:
/// Create empty CCFEM function space
CCFEMFunctionSpace();
/// Destructor
~CCFEMFunctionSpace();
/// Return dimension of the CCFEM function space
///
/// *Returns*
/// std::size_t
/// The dimension of the CCFEM function space.
std::size_t dim() const;
/// Return CCFEM dofmap
///
/// *Returns*
/// _CCFEMDofMap_
/// The dofmap.
std::shared_ptr<const CCFEMDofMap> dofmap() const;
/// Return the number function spaces (parts) of the CCFEM function space
///
/// *Returns*
/// std::size_t
/// The number of function spaces (parts) of the CCFEM function space.
std::size_t num_parts() const;
/// Return function space (part) number i
///
/// *Arguments*
/// i (std::size_t)
/// The part number
///
/// *Returns*
/// _FunctionSpace_
/// Function space (part) number i
std::shared_ptr<const FunctionSpace> part(std::size_t i) const;
/// Return the list of uncut cells for given part. The uncut cells
/// are defined as all cells that don't collide with any cells in
/// any other part with higher part number.
///
/// *Arguments*
/// part (std::size_t)
/// The part number
///
/// *Returns*
/// std::vector<unsigned int>
/// List of uncut cell indices for given part
const std::vector<unsigned int>& uncut_cells(std::size_t part) const;
/// Return the list of cut cells for given part. The cut cells are
/// defined as all cells that collide with the boundary of any
/// part with higher part number.
///
/// FIXME: Figure out whether this makes sense; a cell may collide
/// with the boundary of part j but may still be covered
/// completely by the domain of part j + 1. Possible solution is
/// to for each part i check overlapping parts starting from the
/// top and working back down to i + 1.
///
/// *Arguments*
/// part (std::size_t)
/// The part number
///
/// *Returns*
/// std::vector<unsigned int>
/// List of cut cell indices for given part
const std::vector<unsigned int>& cut_cells(std::size_t part) const;
/// Return the list of covered cells for given part. The covered
/// cells are defined as all cells that collide with the domain of
/// any part with higher part number, but not with the boundary of
/// that part; in other words cells that are completely covered by
/// any other part (and which therefore are inactive).
///
/// *Arguments*
/// part (std::size_t)
/// The part number
///
/// *Returns*
/// std::vector<unsigned int>
/// List of covered cell indices for given part
const std::vector<unsigned int>& covered_cells(std::size_t part) const;
/// Return the collision map for cut cells of the given part
///
/// *Arguments*
/// part (std::size_t)
/// The part number
///
/// *Returns*
/// std::map<unsigned int, std::vector<std::pair<std::size_t, unsigned int> > >
/// A map from cell indices of cut cells to a list of
/// cutting cells. Each cutting cell is represented as a
/// pair (part_number, cutting_cell_index).
const std::map<unsigned int,
std::vector<std::pair<std::size_t, unsigned int> > >&
collision_map_cut_cells(std::size_t part) const;
/// Return quadrature rules for cut cells of the given part
///
/// *Arguments*
/// part (std::size_t)
/// The part number
///
/// *Returns*
/// std::map<unsigned int, std::pair<std::vector<double>, std::vector<double> > >
/// A map from cell indices of cut cells to a quadrature
/// rules. Each quadrature rule is represented as a pair
/// of an array of quadrature weights and a corresponding
/// flattened array of quadrature points.
const std::map<unsigned int, std::pair<std::vector<double>, std::vector<double> > > &
quadrature_rule_cut_cells(std::size_t part) const;
/// Add function space (shared pointer version)
///
/// *Arguments*
/// function_space (_FunctionSpace_)
/// The function space.
void add(std::shared_ptr<const FunctionSpace> function_space);
/// Add function space (reference version)
///
/// *Arguments*
/// function_space (_FunctionSpace_)
/// The function space.
void add(const FunctionSpace& function_space);
/// Build CCFEM function space
void build();
private:
// List of function spaces
std::vector<std::shared_ptr<const FunctionSpace> > _function_spaces;
// CCFEM dofmap
std::shared_ptr<CCFEMDofMap> _dofmap;
// List of meshes
std::vector<std::shared_ptr<const Mesh> > _meshes;
// List of boundary meshes
std::vector<std::shared_ptr<BoundaryMesh> > _boundary_meshes;
// List of bounding box trees for meshes
std::vector<std::shared_ptr<BoundingBoxTree> > _trees;
// List of bounding box trees for boundary meshes
std::vector<std::shared_ptr<BoundingBoxTree> > _boundary_trees;
// Cell indices for all uncut cells for all parts. Access data by
//
// c = _uncut_cells[i][j]
//
// where
//
// c = cell index for an uncut cell
// i = the part (mesh) number
// j = the cell number (in the list of uncut cells)
std::vector<std::vector<unsigned int> > _uncut_cells;
// Cell indices for all cut cells for all parts. Access data by
//
// c = _cut_cells[i][j]
//
// where
//
// c = cell index for a cut cell
// i = the part (mesh) number
// j = the cell number (in the list of cut cells)
std::vector<std::vector<unsigned int> > _cut_cells;
// Cell indices for all covered cells for all parts. Access data by
//
// c = _covered_cells[i][j]
//
// where
//
// c = cell index for a covered cell
// i = the part (mesh) number
// j = the cell number (in the list of covered cells)
std::vector<std::vector<unsigned int> > _covered_cells;
// Developer note 1: The data structures _collision_map_cut_cells
// and _quadrature_rules_cut_cells may be changed from maps to
// vectors and indexed by the number of the cut cell (in the list
// of cut cells), instead of indexed by the local cell index as
// here, if we find that this is important for performance.
//
// Developer note 2: Quadrature points are naturally a part of a
// form (or a term in a form) and not a part of a mesh or function
// space. However, for now we use a global (to the function space)
// quadrature rule for all cut cells, for simplicity.
// Collision map for cut cells. Access data by
//
// c = _collision_map_cut_cells[i][j][k]
//
// where
//
// c.first = part number for the cutting mesh
// c.second = cell index for the cutting cell
// i = the part (mesh) number
// j = the cell number (local cell index)
// k = the collision number (in the list of cutting cells)
std::vector<std::map<unsigned int,
std::vector<std::pair<std::size_t, unsigned int> > > >
_collision_maps_cut_cells;
// Quadrature rules for cut cells. Access data by
//
// q = _quadrature_rules_cut_cells[i][j]
//
// where
//
// q.first = quadrature weights, array of length num_points
// q.second = quadrature points, flattened num_points x gdim array
// i = the part (mesh) number
// j = the cell number (local cell index)
std::vector<std::map<unsigned int,
std::pair<std::vector<double>, std::vector<double> > > >
_quadrature_rules_cut_cells;
// Build dofmap
void _build_dofmap();
// Build boundary meshes
void _build_boundary_meshes();
// Build bounding box trees
void _build_bounding_box_trees();
// Build collision maps
void _build_collision_maps();
// Build quadrature rules
void _build_quadrature_rules();
};
}
#endif
|