/usr/include/dolfin/graph/GraphBuilder.h is in libdolfin-dev 2017.2.0.post0-2.
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 | // Copyright (C) 2010-2013 Garth N. Wells
//
// 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: 2010-02-19
// Last changed: 2013-01-16
#ifndef __GRAPH_BUILDER_H
#define __GRAPH_BUILDER_H
#include <cstdint>
#include <utility>
#include <set>
#include <boost/unordered_map.hpp>
#include <vector>
#include <boost/multi_array.hpp>
#include <dolfin/common/MPI.h>
#include "Graph.h"
namespace dolfin
{
// Forward declarations
class CellType;
class GenericDofMap;
class Mesh;
/// This class builds a Graph corresponding to various objects
class GraphBuilder
{
public:
/// Build local graph from dofmap
static Graph local_graph(const Mesh& mesh, const GenericDofMap& dofmap0,
const GenericDofMap& dofmap1);
/// Build local graph from mesh (general version)
static Graph local_graph(const Mesh& mesh,
const std::vector<std::size_t>& coloring_type);
/// Build local graph (specialized version)
static Graph local_graph(const Mesh& mesh, std::size_t dim0,
std::size_t dim1);
/// Build distributed dual graph (cell-cell connections) from
/// minimal mesh data, and return (num local edges, num
/// non-local edges)
static
std::pair<std::int32_t, std::int32_t>
compute_dual_graph(const MPI_Comm mpi_comm,
const boost::multi_array<std::int64_t, 2>& cell_vertices,
const CellType& cell_type,
const std::int64_t num_global_vertices,
std::vector<std::vector<std::size_t>>& local_graph,
std::set<std::int64_t>& ghost_vertices);
private:
friend class MeshPartitioning;
typedef std::vector<std::pair<std::vector<std::size_t>, std::int32_t>>
FacetCellMap;
// Compute local part of the dual graph, and return number of
// local edges in the graph (undirected)
static std::int32_t
compute_local_dual_graph(const MPI_Comm mpi_comm,
const boost::multi_array<std::int64_t, 2>& cell_vertices,
const CellType& cell_type,
std::vector<std::vector<std::size_t>>& local_graph,
FacetCellMap& facet_cell_map);
// Compute local part of the dual graph, and return number of
// local edges in the graph (undirected)
template<int N>
static std::int32_t
compute_local_dual_graph_keyed(const MPI_Comm mpi_comm,
const boost::multi_array<std::int64_t, 2>& cell_vertices,
const CellType& cell_type,
std::vector<std::vector<std::size_t>>& local_graph,
FacetCellMap& facet_cell_map);
// Build nonlocal part of dual graph for mesh and return number of
// non-local edges. Note: GraphBuilder::compute_local_dual_graph
// should be called before this function is called.
static std::int32_t
compute_nonlocal_dual_graph(const MPI_Comm mpi_comm,
const boost::multi_array<std::int64_t, 2>& cell_vertices,
const CellType& cell_type,
const std::int64_t num_global_vertices,
std::vector<std::vector<std::size_t>>& local_graph,
FacetCellMap& facet_cell_map,
std::set<std::int64_t>& ghost_vertices);
};
}
#endif
|