This file is indexed.

/usr/include/deal.II/lac/sparse_mic.templates.h is in libdeal.ii-dev 6.3.1-1.1.

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
//---------------------------------------------------------------------------
//    Copyright (C) 2002, 2003, 2004, 2005, 2006 by the deal.II authors
//    by the deal.II authors and Stephen "Cheffo" Kolaroff
//
//    This file is subject to QPL and may not be  distributed
//    without copyright and license information. Please refer
//    to the file deal.II/doc/license.html for the  text  and
//    further information on this license.
//
//---------------------------------------------------------------------------
#ifndef __deal2__sparse_mic_templates_h
#define __deal2__sparse_mic_templates_h


#include <base/memory_consumption.h>
#include <lac/sparse_mic.h>
#include <lac/vector.h>

DEAL_II_NAMESPACE_OPEN

template <typename number>
SparseMIC<number>::SparseMIC ()
		:
		diag(0),
		inv_diag(0),
		inner_sums(0)
{}



template <typename number>
SparseMIC<number>::SparseMIC (const SparsityPattern &sparsity)
		:
		SparseLUDecomposition<number> (sparsity),
                diag(0),
                inv_diag(0),
                inner_sums(0)
{}


template <typename number>
SparseMIC<number>::~SparseMIC()
{
  clear();
}


template <typename number>
void SparseMIC<number>::clear()
{
  if (true)
    {
      std::vector<number> tmp;
      tmp.swap (diag);
    };
  if (true)
    {
      std::vector<number> tmp;
      tmp.swap (inv_diag);
    };
  if (true)
    {
      std::vector<number> tmp;
      tmp.swap (inner_sums);
    };

  SparseLUDecomposition<number>::clear();
}


template <typename number>
template <typename somenumber>
inline
void SparseMIC<number>::initialize (const SparseMatrix<somenumber> &matrix,
				    const AdditionalData data)
{
  SparseLUDecomposition<number>::initialize(matrix, data);

  decompose(matrix, data.strengthen_diagonal);
}



template <typename number>
void SparseMIC<number>::reinit (const SparsityPattern& sparsity)
{
  if (true)
    {
      std::vector<number> tmp;
      tmp.swap (diag);
    };
  if (true)
    {
      std::vector<number> tmp;
      tmp.swap (inv_diag);
    };
  if (true)
    {
      std::vector<number> tmp;
      tmp.swap (inner_sums);
    };
  SparseLUDecomposition<number>::reinit (sparsity);
}



template <typename number>
template <typename somenumber>
void SparseMIC<number>::decompose (const SparseMatrix<somenumber> &matrix,
				   const double                    strengthen_diagonal)
{

  SparseLUDecomposition<number>::decompose(matrix, strengthen_diagonal);

  Assert (matrix.m()==matrix.n(), ExcNotQuadratic ());
  Assert (this->m()==this->n(),   ExcNotQuadratic ());
  Assert (matrix.m()==this->m(),  ExcDimensionMismatch(matrix.m(), this->m()));

  Assert (strengthen_diagonal>=0, ExcInvalidStrengthening (strengthen_diagonal));

  if (strengthen_diagonal > 0)
    this->strengthen_diagonal_impl ();

				   // MIC implementation: (S. Margenov lectures)
				   // x[i] = a[i][i] - sum(k=1, i-1,
				   //              a[i][k]/x[k]*sum(j=k+1, N, a[k][j]))

				   // TODO: for sake of simplicity,
				   // those are placed here. A better
				   // implementation would store this
				   // values in the underlying sparse
				   // matrix itself.
  diag.resize (this->m());
  inv_diag.resize (this->m());
  inner_sums.resize (this->m());

				   // precalc sum(j=k+1, N, a[k][j]))
  for(unsigned int row=0; row<this->m(); row++) {
    inner_sums[row] = get_rowsum(row);
  }

  const unsigned int* const col_nums = this->get_sparsity_pattern().get_column_numbers();
  const std::size_t* const rowstarts = this->get_sparsity_pattern().get_rowstart_indices();

  for(unsigned int row=0; row<this->m(); row++)
    {
      number temp = this->diag_element(row);
      number temp1 = 0;
      const unsigned int * const
        first_after_diagonal = this->prebuilt_lower_bound[row];
       
      unsigned int k = 0;
      for (const unsigned int * col=&col_nums[rowstarts[row]+1];
	   col<first_after_diagonal; ++col, ++k)
	temp1 += matrix.global_entry (col-col_nums)/diag[k]*inner_sums[k];
       
      Assert(temp-temp1 > 0, ExcStrengthenDiagonalTooSmall());
      diag[row] = temp - temp1;
       
      inv_diag[row] = 1.0/diag[row];
    }
}



template <typename number>
inline number
SparseMIC<number>::get_rowsum (const unsigned int row) const
{
  Assert(this->m()==this->n(), ExcNotQuadratic());
                                   // get start of this row. skip the
                                   // diagonal element
  const unsigned int * const column_numbers = this->get_sparsity_pattern().get_column_numbers();
  const std::size_t  * const rowstart_indices = this->get_sparsity_pattern().get_rowstart_indices();
  const unsigned int * const rowend = &column_numbers[rowstart_indices[row+1]];

                                   // find the position where the part
                                   // right of the diagonal starts
  const unsigned int * const first_after_diagonal = this->prebuilt_lower_bound[row];
  number rowsum =  0;
  for (const unsigned int * col=first_after_diagonal; col!=rowend; ++col)
    rowsum += this->global_entry (col-column_numbers);

  return rowsum;	
}



template <typename number>
template <typename somenumber>
void
SparseMIC<number>::vmult (Vector<somenumber>       &dst,
                          const Vector<somenumber> &src) const
{
  SparseLUDecomposition<number>::vmult (dst, src);
  Assert (dst.size() == src.size(), ExcDimensionMismatch(dst.size(), src.size()));
  Assert (dst.size() == this->m(), ExcDimensionMismatch(dst.size(), this->m()));

  const unsigned int N=dst.size();
  const std::size_t  * const rowstart_indices = this->get_sparsity_pattern().get_rowstart_indices();
  const unsigned int * const column_numbers   = this->get_sparsity_pattern().get_column_numbers();
                                   // We assume the underlying matrix A is:
                                   // A = X - L - U, where -L and -U are
                                   // strictly lower- and upper- diagonal
                                   // parts of the system.
                                   // 
                                   // Solve (X-L)X{-1}(X-U) x = b
                                   // in 3 steps:
  dst = src;
  for (unsigned int row=0; row<N; ++row)
    {
                                       // Now: (X-L)u = b

                                       // get start of this row. skip
                                       // the diagonal element
      const unsigned int * const rowstart = &column_numbers[rowstart_indices[row]+1];
      const unsigned int * const fad = this->prebuilt_lower_bound[row];
      for (const unsigned int * col=rowstart; col!=fad; ++col)
        dst(row) -= this->global_entry (col-column_numbers) * dst(*col);
      
      dst(row) *= inv_diag[row];
    };

                                   // Now: v = Xu
  for(unsigned int row=0; row<N; row++)
    dst(row) *= diag[row];

                                   // x = (X-U)v
  for (int row=N-1; row>=0; --row)
    {
				       // get end of this row
      const unsigned int * const rowend = &column_numbers[rowstart_indices[row+1]];
      const  unsigned int * const fad = this->prebuilt_lower_bound[row];
      for (const unsigned int * col=fad; col!=rowend; ++col)
        dst(row) -= this->global_entry (col-column_numbers) * dst(*col);

      dst(row) *= inv_diag[row];
    };
}



template <typename number>
unsigned int
SparseMIC<number>::memory_consumption () const
{
  return (SparseLUDecomposition<number>::memory_consumption () +
          MemoryConsumption::memory_consumption(diag) +
          MemoryConsumption::memory_consumption(inv_diag) +
          MemoryConsumption::memory_consumption(inner_sums));
}



DEAL_II_NAMESPACE_CLOSE

#endif // __deal2__sparse_mic_templates_h