This file is indexed.

/usr/include/trilinos/mrtr_utils.H is in libtrilinos-moertel-dev 12.4.2-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
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
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
/*
#@HEADER
# ************************************************************************
#
#                          Moertel FE Package
#                 Copyright (2006) Sandia Corporation
#
# Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
# license for use of this work by or on behalf of the U.S. Government.
#
# Redistribution and use in source and binary forms, with or without
# modification, are permitted provided that the following conditions are
# met:
#
# 1. Redistributions of source code must retain the above copyright
# notice, this list of conditions and the following disclaimer.
#
# 2. Redistributions in binary form must reproduce the above copyright
# notice, this list of conditions and the following disclaimer in the
# documentation and/or other materials provided with the distribution.
#
# 3. Neither the name of the Corporation nor the names of the
# contributors may be used to endorse or promote products derived from
# this software without specific prior written permission.
#
# THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
# EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
# IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
# PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
# CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
# EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
# PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
# PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
# LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
# NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
# SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#
# Questions? Contact Glen Hansen (gahanse@sandia.gov)
#
# ************************************************************************
#@HEADER
*/
/* ******************************************************************** */
/* See the file COPYRIGHT for a complete copyright notice, contact      */
/* person and disclaimer.                                               */
/* ******************************************************************** */
/*!
 * \file mrtr_utils.H
 *
 * \brief A couple utility methods for the Moertel package
 *
 * \date Last update do Doxygen: 16-Dec-05
 *
 */
#ifndef MOERTEL_UTILS_H
#define MOERTEL_UTILS_H

#include <ctime>
#include <iostream>

// external package headers
#include "Teuchos_RCP.hpp"
#include <Epetra_Comm.h>   // Epetra_Comm.h defines std:: as a global namespace
#ifdef HAVE_MPI
#include <Epetra_MpiComm.h>
#else
#include <Epetra_SerialComm.h>
#endif
#include "Epetra_SerialDenseMatrix.h"
#include "Epetra_SerialDenseSolver.h"
#include "Epetra_CrsMatrix.h"
#include "Epetra_Vector.h"
#include "Epetra_Export.h"
#include "ml_common.h"
#include "ml_include.h"
#include "ml_epetra_utils.h"
#include "ml_epetra.h"
#include "ml_epetra_operator.h"

// Moertel package headers
#include "mrtr_segment.H"
#include "mrtr_functions.H"
#include "mrtr_node.H"
#include "mrtr_point.H"

/*!
\brief MOERTEL: namespace of the Moertel package

The Moertel package depends on \ref Epetra, \ref EpetraExt, \ref Teuchos,
\ref Amesos, \ref ML and \ref AztecOO:<br>
Use at least the following lines in the configure of Trilinos:<br>
\code
--enable-moertel 
--enable-epetra 
--enable-epetraext
--enable-teuchos 
--enable-ml
--enable-aztecoo --enable-aztecoo-teuchos 
--enable-amesos
\endcode

*/
namespace MOERTEL
{

// forward declarations
class Segment;
class Node;

/*!
\brief Allocate a function of the correct type

 For communication reasons, every single derived function class needs
 to have a unique typ-id. This type Id can be communicated easily.
 So when introducing a new derived Function class, one needs to add
 it's type to the enum FunctionType in the virtual base class in
 mrtr_function.H and one needs to add a case to this method
 MOERTEL::AllocateFunction in mrtr_utils.cpp 

\param type : Type of Function to allocate and return pointer to
*/
MOERTEL::Function* AllocateFunction(MOERTEL::Function::FunctionType type, int out);

/*!
\brief Allocate a Segment of the correct type

 For communication reasons, every single derived segment class needs
 to have a unique typ-id. This type Id can be communicated easily.
 So when introducing a new segment class, one needs to add
 it's type to the enum SegmentType in the virtual base class in
 mrtr_segment.H and one needs to add a case to this method
 MOERTEL::AllocateSegment in mrtr_utils.cpp 

\param type : Type of segment to allocate and return pointer to
\param out : Level of output to be generated to stdout ( 0 - 10 )
*/
MOERTEL::Segment* AllocateSegment(int type, int out);


/*!
\brief Cross product

Perform cross product out = g1 x g2 for vectors of dimension 3

*/
bool cross(double* out, const double* g1, const double* g2);

/*!
\brief Dot product

Perform dot product g1 dot g2 for vectors of dimension dim and return result

*/
double dot(const double* g1, const double* g2, const int dim);

/*!
\brief Length of a vector

Return L2 norm of a vector of dimension dim

*/
double length(const double* g, const int dim);

/*!
\brief Solve dense 2x2 system of equations

Ax=b

*/
bool solve22(const double A[][2], double* x, const double* b);

/*!
\brief Solve dense 3x3 system of equations

Ax=b

*/
bool solve33(const double A[][3], double* x, const double* b);


/*!
\brief Return the '10' digit from an integer number
*/
int digit_ten(int i);

/*!
\brief Sort dlist of length N in ascending, sort list2 according to dlist 

This piece of code was lend from the Trilinos package ML
*/
void sort(double* dlist, int N, int* list2);

/*!
\brief Template to swap 2 <kind> instances

<kind> has to implement the assignment operator = 
*/
template<typename kind>
void swap(kind& a, kind& b)
{
  kind tmp = a;
  a = b;
  b = tmp;
  return;
}

/*!
\brief Add matrices A+B

Perform B = scalarB * B + scalarA * A ^ transposeA
If scalarB is 0.0, then B = scalarA * A ^ transposeA isperformed.

This is a modified version of EpetraExt's MatrixMatrixAdd.
FillComplete() must not be called on B upon entry, 
FillComplete() will not be called on B upon exit by this method.

\param A : Matrix A to add to B
\param transposeA : flag indicating whether A*T shall be added
\param scalarA : scalar factor for A
\param B : Matrix B to be added to
\param scalarB : scalar factor for B
\return Zero upon success
*/
int MatrixMatrixAdd(const Epetra_CrsMatrix& A, bool transposeA,double scalarA,
                    Epetra_CrsMatrix& B,double scalarB);

#if 0
/*!
\brief Multiply matrices A*B

matrices A and B are mutliplied and the result is allocated and returned.
The method makes use of the ML matrix-matrix mutliply functions.
The user is responsible for freeing the returned result.

\param A : Matrix A to multiply
\param transA : flag indicating whether A*T shall be used
\param B : Matrix B to multiply
\param transB : flag indicating whether B*T shall be used
\return Result upon success and NULL upon failure
\warning The method fails if any of A or B have empty columns. 
See ML bug 1913 for details.
*/
Epetra_CrsMatrix* MatMatMult(Epetra_CrsMatrix& A, bool transA, 
                             Epetra_CrsMatrix& B, bool transB);
#endif

/*!
\brief Multiply matrices A*B

matrices A and B are mutliplied and the result is allocated and returned.
The method makes uses EpetraExt for multiplication
The user is responsible for freeing the returned result.

\param A : Matrix A to multiply
\param transA : flag indicating whether A*T shall be used
\param B : Matrix B to multiply
\param transB : flag indicating whether B*T shall be used
\return Result upon success and NULL upon failure
*/
Epetra_CrsMatrix* MatMatMult(const Epetra_CrsMatrix& A, bool transA, 
                             const Epetra_CrsMatrix& B, bool transB,
                             int outlevel);


/*!
\brief Allocate and return a matrix padded with val on the diagonal. 
       FillComplete() is NOT called on exit.
*/
Epetra_CrsMatrix* PaddedMatrix(const Epetra_Map rowmap, double val,const int numentriesperrow);


/*!
\brief Strip out values from a matrix below a certain tolerance

Allocates and returns a new matrix and copies A to it where entries
with an absoute value smaller then eps are negelected.
The method calls FillComplete(A.OperatorDomainMap(),A.OperatorRangeMap())
on the result.

\param A : Matrix A to strip
\param eps : tolerance
\return The new matrix upon success, NULL otherwise
*/
Epetra_CrsMatrix* StripZeros(Epetra_CrsMatrix& A, double eps);

/*!
\brief split a matrix into a 2x2 block system where the rowmap of one of the blocks is given

Splits a given matrix into a 2x2 block system where the rowmap of one of the blocks is given
on input. Blocks A11 and A22 are assumed to be square.
All values on entry have to be Teuchos::null except the given rowmap and matrix A.
Note that either A11rowmap or A22rowmap or both have to be nonzero. In case
both rowmaps are supplied they have to be an exact and nonoverlapping split of A->RowMap().
Matrix blocks are FillComplete() on exit.

\param A         : Matrix A on input
\param A11rowmap : rowmap of A11 or null 
\param A22rowmap : rowmap of A22 or null 
\param A11       : on exit matrix block A11 
\param A12       : on exit matrix block A12 
\param A21       : on exit matrix block A21 
\param A22       : on exit matrix block A22 
*/
bool SplitMatrix2x2(Teuchos::RCP<Epetra_CrsMatrix> A,
                    Teuchos::RCP<Epetra_Map>& A11rowmap,
                    Teuchos::RCP<Epetra_Map>& A22rowmap,
                    Teuchos::RCP<Epetra_CrsMatrix>& A11,
                    Teuchos::RCP<Epetra_CrsMatrix>& A12,
                    Teuchos::RCP<Epetra_CrsMatrix>& A21,
                    Teuchos::RCP<Epetra_CrsMatrix>& A22);

/*!
\brief split a rowmap of matrix A

splits A->RowMap() into 2 maps and returns them, where one of the rowmaps has
to be given on input

\param Amap      : Map to split on input
\param Agiven    : on entry submap that is given and part of Amap 
\return the remainder map of Amap that is not overlapping with Agiven 
*/
Epetra_Map* SplitMap(const Epetra_Map& Amap,
                     const Epetra_Map& Agiven);

/*!
\brief split a vector into 2 non-overlapping pieces

*/
bool SplitVector(const Epetra_Vector& x,
                 const Epetra_Map& x1map,
                 Epetra_Vector*&   x1,
                 const Epetra_Map& x2map,
                 Epetra_Vector*&   x2);

/*!
\brief merge results from 2 vectors into one (assumes matching submaps)

*/
bool MergeVector(const Epetra_Vector& x1,
                 const Epetra_Vector& x2,
                 Epetra_Vector& xresult);

/*!
\brief Print matrix to file

Prints an Epetra_CrsMatrix to file in serial and parallel.
Will create several files with process id appended to the name in parallel.
Index base can either be 0 or 1.
The first row of the file gives the global size of the range and domain map, 
the sond row gives the local size of the row- and column map. 

\param name : Name of file without appendix, appendix will be .mtx
\param A : Matrix to print
\param ibase : Index base, should be either 1 or 0 
*/
bool Print_Matrix(std::string name, Epetra_CrsMatrix& A, int ibase);

/*!
\brief Print graph to file

Prints an Epetra_CrsGraph to file in serial and parallel.
Will create several files with process id appended to the name in parallel.
Index base can either be 0 or 1.
The first row of the file gives the global size of the range and domain map, 
the second row gives the local size of the row- and column map. 

\param name : Name of file without appendix, appendix will be .mtx
\param A : Graph to print
\param ibase : Index base, should be either 1 or 0 
*/
bool Print_Graph(std::string name, Epetra_CrsGraph& A, int ibase);

/*!
\brief Print vector to file

Prints an Epetra_Vector to file in serial and parallel.
Will create several files with process id appended to the name in parallel.
Index base can either be 0 or 1.

\param name : Name of file without appendix, appendix will be .vec
\param v : Vector to print
\param ibase : Index base, should be either 1 or 0 
*/
bool Print_Vector(std::string name, Epetra_Vector& v, int ibase);

//! Error reporting method
int ReportError(const std::stringstream &Message);

} // namespace MOERTEL
#endif // MOERTEL_UTILS_H