Reference documentation for deal.II version 9.0.0
diagonal_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_diagonal_matrix_h
17 #define dealii_diagonal_matrix_h
18 
19 
20 #include <deal.II/lac/vector.h>
21 #include <deal.II/lac/vector_operation.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
46 template <typename VectorType=Vector<double> >
48 {
49 public:
50  typedef typename VectorType::value_type value_type;
51  typedef typename VectorType::size_type size_type;
52 
56  DiagonalMatrix() = default;
57 
62  void reinit (const VectorType &vec);
63 
69  void compress (VectorOperation::values operation);
70 
75  VectorType &get_vector();
76 
80  void clear();
81 
85  const VectorType &get_vector() const;
86 
91  size_type m () const;
92 
97  size_type n () const;
98 
108  value_type operator()(const size_type i,
109  const size_type j) const;
110 
120  value_type &operator()(const size_type i,
121  const size_type j);
122 
134  template <typename number2>
135  void add (const size_type row,
136  const size_type n_cols,
137  const size_type *col_indices,
138  const number2 *values,
139  const bool elide_zero_values = true,
140  const bool col_indices_are_sorted = false);
141 
148  void add (const size_type i,
149  const size_type j,
150  const value_type value);
151 
155  void vmult (VectorType &dst,
156  const VectorType &src) const;
157 
163  void Tvmult (VectorType &dst,
164  const VectorType &src) const;
165 
171  void vmult_add (VectorType &dst,
172  const VectorType &src) const;
173 
179  void Tvmult_add (VectorType &dst,
180  const VectorType &src) const;
181 
189  void initialize_dof_vector(VectorType &dst) const;
190 
194  std::size_t memory_consumption () const;
195 
196 private:
200  VectorType diagonal;
201 };
202 
203 /* ---------------------------------- Inline functions ------------------- */
204 
205 #ifndef DOXYGEN
206 
207 template <typename VectorType>
208 void
210 {
211  diagonal.reinit(0);
212 }
213 
214 
215 
216 template <typename VectorType>
217 std::size_t
219 {
220  return diagonal.memory_consumption();
221 }
222 
223 
224 
225 template <typename VectorType>
226 void
227 DiagonalMatrix<VectorType>::reinit(const VectorType &vec)
228 {
229  diagonal = vec;
230 }
231 
232 
233 
234 template <typename VectorType>
235 void
237 {
238  dst.reinit(diagonal);
239 }
240 
241 
242 
243 template <typename VectorType>
244 void
246 {
247  diagonal.compress(operation);
248 }
249 
250 
251 
252 template <typename VectorType>
253 VectorType &
255 {
256  return diagonal;
257 }
258 
259 
260 
261 template <typename VectorType>
262 const VectorType &
264 {
265  return diagonal;
266 }
267 
268 
269 
270 template <typename VectorType>
271 typename VectorType::size_type
273 {
274  return diagonal.size();
275 }
276 
277 
278 
279 template <typename VectorType>
280 typename VectorType::size_type
282 {
283  return diagonal.size();
284 }
285 
286 
287 
288 template <typename VectorType>
289 typename VectorType::value_type
291  const size_type j) const
292 {
293  Assert (i==j, ExcIndexRange(j,i,i+1));
294  (void) j;
295  return diagonal(i);
296 }
297 
298 
299 
300 template <typename VectorType>
301 typename VectorType::value_type &
303  const size_type j)
304 {
305  Assert (i==j, ExcIndexRange(j,i,i+1));
306  (void)j;
307  return diagonal(i);
308 }
309 
310 
311 
312 template <typename VectorType>
313 template <typename number2>
314 void
315 DiagonalMatrix<VectorType>::add (const size_type row,
316  const size_type n_cols,
317  const size_type *col_indices,
318  const number2 *values,
319  const bool,
320  const bool )
321 {
322  for (size_type i=0; i<n_cols; ++i)
323  if (col_indices[i] == row)
324  diagonal(row) += values[i];
325 }
326 
327 
328 
329 template <typename VectorType>
330 void
331 DiagonalMatrix<VectorType>::add (const size_type i,
332  const size_type j,
333  const value_type value)
334 {
335  if (i == j)
336  diagonal(i) += value;
337 }
338 
339 
340 
341 template <typename VectorType>
342 void
343 DiagonalMatrix<VectorType>::vmult(VectorType &dst,
344  const VectorType &src) const
345 {
346  dst = src;
347  dst.scale(diagonal);
348 }
349 
350 
351 
352 template <typename VectorType>
353 void
354 DiagonalMatrix<VectorType>::Tvmult(VectorType &dst,
355  const VectorType &src) const
356 {
357  vmult(dst, src);
358 }
359 
360 
361 
362 template <typename VectorType>
363 void
365  const VectorType &src) const
366 {
367  VectorType tmp(src);
368  tmp.scale(diagonal);
369  dst += tmp;
370 }
371 
372 
373 
374 template <typename VectorType>
375 void
377  const VectorType &src) const
378 {
379  vmult_add(dst, src);
380 }
381 
382 
383 #endif
384 
385 DEAL_II_NAMESPACE_CLOSE
386 
387 #endif
std::size_t memory_consumption() const
VectorType & get_vector()
void Tvmult(VectorType &dst, const VectorType &src) const
void compress(VectorOperation::values operation)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
value_type operator()(const size_type i, const size_type j) const
void add(const size_type row, const size_type n_cols, const size_type *col_indices, const number2 *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
Matrix is diagonal.
size_type n() const
void vmult(VectorType &dst, const VectorType &src) const
#define Assert(cond, exc)
Definition: exceptions.h:1142
size_type m() const
VectorType diagonal
DiagonalMatrix()=default
void initialize_dof_vector(VectorType &dst) const
void reinit(const VectorType &vec)
void vmult_add(VectorType &dst, const VectorType &src) const
void Tvmult_add(VectorType &dst, const VectorType &src) const