This file is indexed.

/usr/include/sundials/sundials_iterative.h is in libsundials-serial-dev 2.5.0-3+b3.

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
/*
 * -----------------------------------------------------------------
 * $Revision: 1.2 $
 * $Date: 2006/11/29 00:05:07 $
 * ----------------------------------------------------------------- 
 * Programmer(s): Scott D. Cohen and Alan C. Hindmarsh @ LLNL
 * -----------------------------------------------------------------
 * Copyright (c) 2002, The Regents of the University of California.
 * Produced at the Lawrence Livermore National Laboratory.
 * All rights reserved.
 * For details, see the LICENSE file.
 * -----------------------------------------------------------------
 * This header file contains declarations intended for use by
 * generic iterative solvers of Ax = b. The enumeration gives
 * symbolic names for the type  of preconditioning to be used.
 * The function type declarations give the prototypes for the
 * functions to be called within an iterative linear solver, that
 * are responsible for
 *    multiplying A by a given vector v (ATimesFn), and
 *    solving the preconditioner equation Pz = r (PSolveFn).
 * -----------------------------------------------------------------
 */

#ifndef _ITERATIVE_H
#define _ITERATIVE_H

#ifdef __cplusplus  /* wrapper to enable C++ usage */
extern "C" {
#endif

#include <sundials/sundials_nvector.h>


/*
 * -----------------------------------------------------------------
 * enum : types of preconditioning                                
 * -----------------------------------------------------------------
 * PREC_NONE  : The iterative linear solver should not use             
 *              preconditioning.                                       
 *                                                                
 * PREC_LEFT  : The iterative linear solver uses preconditioning on    
 *              the left only.                                         
 *                                                                
 * PREC_RIGHT : The iterative linear solver uses preconditioning on    
 *              the right only.                                        
 *                                                                
 * PREC_BOTH  : The iterative linear solver uses preconditioning on    
 *              both the left and the right.                           
 * -----------------------------------------------------------------
 */

enum { PREC_NONE, PREC_LEFT, PREC_RIGHT, PREC_BOTH };

/*
 * -----------------------------------------------------------------
 * enum : types of Gram-Schmidt routines                          
 * -----------------------------------------------------------------
 * MODIFIED_GS  : The iterative solver uses the modified          
 *                Gram-Schmidt routine ModifiedGS listed in this  
 *                file.                                           
 *                                                                
 * CLASSICAL_GS : The iterative solver uses the classical         
 *                Gram-Schmidt routine ClassicalGS listed in this 
 *                file.                                           
 * -----------------------------------------------------------------
 */

enum { MODIFIED_GS = 1, CLASSICAL_GS = 2 };

/*
 * -----------------------------------------------------------------
 * Type: ATimesFn                                                 
 * -----------------------------------------------------------------
 * An ATimesFn multiplies Av and stores the result in z. The      
 * caller is responsible for allocating memory for the z vector.  
 * The parameter A_data is a pointer to any information about A   
 * which the function needs in order to do its job. The vector v  
 * is unchanged. An ATimesFn returns 0 if successful and a        
 * non-zero value if unsuccessful.                                
 * -----------------------------------------------------------------
 */

typedef int (*ATimesFn)(void *A_data, N_Vector v, N_Vector z);

/*
 * -----------------------------------------------------------------
 * Type: PSolveFn                                                 
 * -----------------------------------------------------------------
 * A PSolveFn solves the preconditioner equation Pz = r for the   
 * vector z. The caller is responsible for allocating memory for  
 * the z vector. The parameter P_data is a pointer to any         
 * information about P which the function needs in order to do    
 * its job. The parameter lr is input, and indicates whether P    
 * is to be taken as the left preconditioner or the right         
 * preconditioner: lr = 1 for left and lr = 2 for right.          
 * If preconditioning is on one side only, lr can be ignored.     
 * The vector r is unchanged.                                     
 * A PSolveFn returns 0 if successful and a non-zero value if     
 * unsuccessful.  On a failure, a negative return value indicates 
 * an unrecoverable condition, while a positive value indicates   
 * a recoverable one, in which the calling routine may reattempt  
 * the solution after updating preconditioner data.               
 * -----------------------------------------------------------------
 */

typedef int (*PSolveFn)(void *P_data, N_Vector r, N_Vector z, int lr);

/*
 * -----------------------------------------------------------------
 * Function: ModifiedGS                                           
 * -----------------------------------------------------------------
 * ModifiedGS performs a modified Gram-Schmidt orthogonalization  
 * of the N_Vector v[k] against the p unit N_Vectors at           
 * v[k-1], v[k-2], ..., v[k-p].                                   
 *                                                                
 * v is an array of (k+1) N_Vectors v[i], i=0, 1, ..., k.         
 * v[k-1], v[k-2], ..., v[k-p] are assumed to have L2-norm        
 * equal to 1.                                                    
 *                                                                
 * h is the output k by k Hessenberg matrix of inner products.    
 * This matrix must be allocated row-wise so that the (i,j)th     
 * entry is h[i][j]. The inner products (v[i],v[k]),              
 * i=i0, i0+1, ..., k-1, are stored at h[i][k-1]. Here            
 * i0=MAX(0,k-p).                                                 
 *                                                                
 * k is the index of the vector in the v array that needs to be   
 * orthogonalized against previous vectors in the v array.        
 *                                                                
 * p is the number of previous vectors in the v array against     
 * which v[k] is to be orthogonalized.                            
 *                                                                
 * new_vk_norm is a pointer to memory allocated by the caller to  
 * hold the Euclidean norm of the orthogonalized vector v[k].     
 *                                                                
 * If (k-p) < 0, then ModifiedGS uses p=k. The orthogonalized     
 * v[k] is NOT normalized and is stored over the old v[k]. Once   
 * the orthogonalization has been performed, the Euclidean norm   
 * of v[k] is stored in (*new_vk_norm).                           
 *                                                                
 * ModifiedGS returns 0 to indicate success. It cannot fail.      
 * -----------------------------------------------------------------
 */                                                                

SUNDIALS_EXPORT int ModifiedGS(N_Vector *v, realtype **h, int k, int p, 
			       realtype *new_vk_norm);

/*
 * -----------------------------------------------------------------
 * Function: ClassicalGS                                          
 * -----------------------------------------------------------------
 * ClassicalGS performs a classical Gram-Schmidt                  
 * orthogonalization of the N_Vector v[k] against the p unit      
 * N_Vectors at v[k-1], v[k-2], ..., v[k-p]. The parameters v, h, 
 * k, p, and new_vk_norm are as described in the documentation    
 * for ModifiedGS.                                                
 *                                                                
 * temp is an N_Vector which can be used as workspace by the      
 * ClassicalGS routine.                                           
 *                                                                
 * s is a length k array of realtype which can be used as         
 * workspace by the ClassicalGS routine.                          
 *
 * ClassicalGS returns 0 to indicate success. It cannot fail.     
 * -----------------------------------------------------------------
 */

SUNDIALS_EXPORT int ClassicalGS(N_Vector *v, realtype **h, int k, int p, 
				realtype *new_vk_norm, N_Vector temp, realtype *s);

/*
 * -----------------------------------------------------------------
 * Function: QRfact                                               
 * -----------------------------------------------------------------
 * QRfact performs a QR factorization of the Hessenberg matrix H. 
 *                                                                
 * n is the problem size; the matrix H is (n+1) by n.             
 *                                                                
 * h is the (n+1) by n Hessenberg matrix H to be factored. It is  
 * stored row-wise.                                               
 *                                                                
 * q is an array of length 2*n containing the Givens rotations    
 * computed by this function. A Givens rotation has the form:     
 * | c  -s |                                                      
 * | s   c |.                                                     
 * The components of the Givens rotations are stored in q as      
 * (c, s, c, s, ..., c, s).                                       
 *                                                                
 * job is a control flag. If job==0, then a new QR factorization  
 * is performed. If job!=0, then it is assumed that the first     
 * n-1 columns of h have already been factored and only the last  
 * column needs to be updated.                                    
 *                                                                
 * QRfact returns 0 if successful. If a zero is encountered on    
 * the diagonal of the triangular factor R, then QRfact returns   
 * the equation number of the zero entry, where the equations are 
 * numbered from 1, not 0. If QRsol is subsequently called in     
 * this situation, it will return an error because it could not   
 * divide by the zero diagonal entry.                             
 * -----------------------------------------------------------------
 */                                                                

SUNDIALS_EXPORT int QRfact(int n, realtype **h, realtype *q, int job);

/*                                                                
 * -----------------------------------------------------------------
 * Function: QRsol                                                
 * -----------------------------------------------------------------
 * QRsol solves the linear least squares problem                  
 *                                                                
 * min (b - H*x, b - H*x), x in R^n,                              
 *                                                                
 * where H is a Hessenberg matrix, and b is in R^(n+1).           
 * It uses the QR factors of H computed by QRfact.                
 *                                                                
 * n is the problem size; the matrix H is (n+1) by n.             
 *                                                                
 * h is a matrix (computed by QRfact) containing the upper        
 * triangular factor R of the original Hessenberg matrix H.       
 *                                                                
 * q is an array of length 2*n (computed by QRfact) containing    
 * the Givens rotations used to factor H.                         
 *                                                                
 * b is the (n+1)-vector appearing in the least squares problem   
 * above.                                                         
 *                                                                
 * On return, b contains the solution x of the least squares      
 * problem, if QRsol was successful.                              
 *                                                                
 * QRsol returns a 0 if successful.  Otherwise, a zero was        
 * encountered on the diagonal of the triangular factor R.        
 * In this case, QRsol returns the equation number (numbered      
 * from 1, not 0) of the zero entry.                              
 * -----------------------------------------------------------------
 */                                                                

SUNDIALS_EXPORT int QRsol(int n, realtype **h, realtype *q, realtype *b);

#ifdef __cplusplus
}
#endif

#endif