sparsity_pattern.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 #ifndef LIBMESH_SPARSITY_PATTERN_H
20 #define LIBMESH_SPARSITY_PATTERN_H
21 
22 // Local Includes
23 #include "libmesh/elem_range.h"
26 
27 // C++ includes
28 #include <vector>
29 
30 namespace libMesh
31 {
32 
33  // Forward declaractions
34  class MeshBase;
35  class DofMap;
36  class CouplingMatrix;
37 
38 // ------------------------------------------------------------
39 // Sparsity Pattern
40 
47 namespace SparsityPattern // use a namespace so member classes can be forward-declared.
48 {
49  typedef std::vector<dof_id_type, Threads::scalable_allocator<dof_id_type> > Row;
50  class Graph : public std::vector<Row> {};
51 
52  class NonlocalGraph : public std::map<dof_id_type, Row> {};
53 
63  template<typename BidirectionalIterator>
64  static void sort_row (const BidirectionalIterator begin,
65  BidirectionalIterator middle,
66  const BidirectionalIterator end);
67 
79  class Build : public ParallelObject
80  {
81  private:
82  const MeshBase &mesh;
83  const DofMap &dof_map;
87 
88  public:
89 
92 
93  std::vector<dof_id_type> n_nz;
94  std::vector<dof_id_type> n_oz;
95 
96  Build (const MeshBase &mesh_in,
97  const DofMap &dof_map_in,
98  const CouplingMatrix *dof_coupling_in,
99  const bool implicit_neighbor_dofs_in,
100  const bool need_full_sparsity_pattern_in);
101 
102  Build (Build &other, Threads::split);
103 
104  void operator()(const ConstElemRange &range);
105 
106  void join (const Build &other);
107 
108  void parallel_sync ();
109  };
110 
111 #if defined(__GNUC__) && (__GNUC__ < 4) && !defined(__INTEL_COMPILER)
112 
117  void _dummy_function(void);
118 #endif
119 
120 }
121 
122 
123 
124 // ------------------------------------------------------------
125 // SparsityPattern inline member functions
126 template<typename BidirectionalIterator>
127 inline
128 void SparsityPattern::sort_row (const BidirectionalIterator begin,
129  BidirectionalIterator middle,
130  const BidirectionalIterator end)
131 {
132  if ((begin == middle) || (middle == end)) return;
133 
134  libmesh_assert_greater (std::distance (begin, middle), 0);
135  libmesh_assert_greater (std::distance (middle, end), 0);
136  libmesh_assert (std::unique (begin, middle) == middle);
137  libmesh_assert (std::unique (middle, end) == end);
138 
139  while (middle != end)
140  {
141  BidirectionalIterator
142  b = middle,
143  a = b-1;
144 
145  // Bubble-sort the middle value downward
146  while (!(*a < *b)) // *a & *b are less-than comparable, so use <
147  {
148  std::swap (*a, *b);
149 
150 #if defined(__GNUC__) && (__GNUC__ < 4) && !defined(__INTEL_COMPILER)
151  /* Prohibit optimization at this point since gcc 3.3.5 seems
152  to have a bug. */
154 #endif
155 
156  if (a == begin) break;
157 
158  b=a;
159  --a;
160  }
161 
162  ++middle;
163  }
164 
165  // Assure the algorithm worked if we are in DEBUG mode
166 #ifdef DEBUG
167  {
168  // SGI STL extension!
169  // libmesh_assert (std::is_sorted(begin,end));
170 
171  BidirectionalIterator
172  prev = begin,
173  first = begin;
174 
175  for (++first; first != end; prev=first, ++first)
176  if (*first < *prev)
177  libmesh_assert(false);
178  }
179 #endif
180 
181  // Make sure the two ranges did not contain any common elements
182  libmesh_assert (std::unique (begin, end) == end);
183 }
184 
185 } // namespace libMesh
186 
187 #endif // LIBMESH_SPARSITY_PATTERN_H

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:07 UTC

Hosted By:
SourceForge.net Logo