Free cookie consent management tool by TermsFeed Policy Generator

source: branches/2925_AutoDiffForDynamicalModels/HeuristicLab.Problems.DynamicalSystemsModelling/3.3/sundials/include/sunlinsol/sunlinsol_sptfqmr.h @ 16222

Last change on this file since 16222 was 16222, checked in by gkronber, 5 years ago

#2925:

  • added comments about parameter identification for differential equation models
  • added source code of cvodes library (part of sundials) which provides functionality to calculate gradients for the parameters of partial differential equation models efficiently using the 'adjoint state method'.
  • added compiled version of cvodes
File size: 9.3 KB
RevLine 
[16222]1/*
2 * -----------------------------------------------------------------
3 * Programmer(s): Daniel Reynolds @ SMU
4 * Based on code sundials_sptfqmr.h by: Aaron Collier @ LLNL
5 * -----------------------------------------------------------------
6 * LLNS/SMU Copyright Start
7 * Copyright (c) 2017, Southern Methodist University and
8 * Lawrence Livermore National Security
9 *
10 * This work was performed under the auspices of the U.S. Department
11 * of Energy by Southern Methodist University and Lawrence Livermore
12 * National Laboratory under Contract DE-AC52-07NA27344.
13 * Produced at Southern Methodist University and the Lawrence
14 * Livermore National Laboratory.
15 *
16 * All rights reserved.
17 * For details, see the LICENSE file.
18 * LLNS/SMU Copyright End
19 * -----------------------------------------------------------------
20 * This is the header file for the SPTFQMR implementation of the
21 * SUNLINSOL module.  The SPTFQMR algorithm is based on the
22 * Scaled Preconditioned Transpose-free Quasi-Minimum Residual method.
23 *
24 * The SPTFQMR algorithm solves a linear system A x = b.
25 * Preconditioning is allowed on the left, right, or both.
26 * Scaling is allowed on both sides.
27 * We denote the preconditioner and scaling matrices as follows:
28 *   P1 = left preconditioner
29 *   P2 = right preconditioner
30 *   S1 = diagonal matrix of scale factors for P1-inverse b
31 *   S2 = diagonal matrix of scale factors for P2 x
32 * The matrices A, P1, and P2 are not required explicitly; only
33 * routines that provide A, P1-inverse, and P2-inverse as
34 * operators are required.
35 *
36 * In this notation, SPTFQMR applies the underlying TFQMR method to
37 * the equivalent transformed system
38 *   Abar xbar = bbar , where
39 *   Abar = S1 (P1-inverse) A (P2-inverse) (S2-inverse) ,
40 *   bbar = S1 (P1-inverse) b , and   xbar = S2 P2 x .
41 *
42 * The scaling matrices must be chosen so that vectors S1
43 * P1-inverse b and S2 P2 x have dimensionless components.
44 * If preconditioning is done on the left only (P2 = I), by a
45 * matrix P, then S2 must be a scaling for x, while S1 is a
46 * scaling for P-inverse b, and so may also be taken as a scaling
47 * for x.  Similarly, if preconditioning is done on the right only
48 * (P1 = I, P2 = P), then S1 must be a scaling for b, while S2 is
49 * a scaling for P x, and may also be taken as a scaling for b.
50 *
51 * The stopping test for the SPTFQMR iterations is on the L2 norm of
52 * the scaled preconditioned residual:
53 *      || bbar - Abar xbar ||_2  <  delta
54 * with an input test constant delta.
55 *
56 * The usage of this SPTFQMR solver involves supplying up to three
57 * routines and making a variety of calls.  The user-supplied routines are
58 *    atimes (A_data, x, y) to compute y = A x, given x,
59 *    psolve (P_data, y, x, lr) to solve P1 x = y or P2 x = y for
60 *           x, given y,
61 *    psetup (P_data) to perform any 'setup' operations in
62 *           preparation for calling psolve.
63 * The three user calls are:
64 *    SUNLinearSolver LS = SUNSPTFQMR(y, pretype, maxl);
65 *           to create the linear solver structure,
66 *    flag = SUNLinSolSetATimes(LS, A_data, atimes);
67 *           to set the matrix-vector product setup/apply routines,
68 *    flag = SUNLinSolSetPreconditioner(LS, P_data, psetup, psolve);
69 *           to *optionally* set the preconditioner setup/apply routines,
70 *    flag = SUNLinSolInitialize(LS);
71 *           to perform internal solver memory allocations,
72 *    flag = SUNLinSolSetup(LS, NULL);
73 *           to call the psetup routine (if non-NULL);
74 *    flag = SUNLinSolSolve(LS, NULL, x, b, w, tol);
75 *           to solve the linear system to the tolerance 'tol'
76 *    long int nli = SUNLinSolNumIters(LS);
77 *           to *optionally* retrieve the number of linear iterations
78 *           performed by the solver,
79 *    long int lastflag = SUNLinSolLastFlag(LS);
80 *           to *optionally* retrieve the last internal solver error flag,
81 *    flag = SUNLinSolFree(LS);
82 *           to free the solver memory.
83 * Complete details for specifying atimes, psetup and psolve
84 * and for the usage calls are given below.
85 *
86 * -----------------------------------------------------------------
87 *
88 * Part I contains declarations specific to the SPTFQMR implementation
89 * of the supplied SUNLINSOL module.
90 *
91 * Part II contains the prototype for the constructor
92 * SUNSPTFQMR as well as implementation-specific prototypes
93 * for various useful solver operations.
94 *
95 * Notes:
96 *
97 *   - The definition of the generic SUNLinearSolver structure can
98 *     be found in the header file sundials_linearsolver.h.
99 *
100 * -----------------------------------------------------------------
101 */
102
103#ifndef _SUNLINSOL_SPTFQMR_H
104#define _SUNLINSOL_SPTFQMR_H
105
106#include <sundials/sundials_linearsolver.h>
107#include <sundials/sundials_matrix.h>
108#include <sundials/sundials_nvector.h>
109#include <sundials/sundials_sptfqmr.h>
110
111#ifdef __cplusplus  /* wrapper to enable C++ usage */
112extern "C" {
113#endif
114
115/* Default SPTFQMR solver parameters */
116#define SUNSPTFQMR_MAXL_DEFAULT    5
117
118/*
119 * -----------------------------------------------------------------
120 * PART I: SPTFQMR implementation of SUNLinearSolver
121 *
122 * The SPTFQMR implementation of the SUNLinearSolver 'content'
123 * structure contains:
124 *     maxl -- number of BiCGStab iterations to allow
125 *     pretype -- flag for type of preconditioning to employ
126 *     numiters -- number of iterations from most-recent solve
127 *     resnorm -- final linear residual norm from most-recent solve
128 *     last_flag -- last error return flag from internal setup/solve
129 *     ATimes -- function pointer to ATimes routine
130 *     ATData -- pointer to structure for ATimes
131 *     Psetup -- function pointer to preconditioner setup routine
132 *     Psolve -- function pointer to preconditioner solve routine
133 *     PData -- pointer to structure for Psetup/Psolve
134 *     s1, s2 -- vector pointers for supplied scaling matrices
135 *     r_star -- a vector (type N_Vector) which holds the initial
136 *         scaled, preconditioned linear system residual
137 *     q, d, v, p and u -- vectors (type N_Vector) used for
138 *         workspace by the SPTFQMR algorithm
139 *     r -- array of vectors (type N_Vector) used for workspace
140 *         within the SPTFQMR algorithm
141 *     vtemp1/vtemp2/vtemp3 -- scratch vectors (type N_Vector) used
142 *         as temporary vector storage during calculations.
143 * -----------------------------------------------------------------
144 */
145 
146struct _SUNLinearSolverContent_SPTFQMR {
147  int maxl;
148  int pretype;
149  int numiters;
150  realtype resnorm;
151  long int last_flag;
152
153  ATimesFn ATimes;
154  void* ATData;
155  PSetupFn Psetup;
156  PSolveFn Psolve;
157  void* PData;
158
159  N_Vector s1;
160  N_Vector s2;
161  N_Vector r_star;
162  N_Vector q;
163  N_Vector d;
164  N_Vector v;
165  N_Vector p;
166  N_Vector *r;
167  N_Vector u;
168  N_Vector vtemp1;
169  N_Vector vtemp2;
170  N_Vector vtemp3;
171};
172
173typedef struct _SUNLinearSolverContent_SPTFQMR *SUNLinearSolverContent_SPTFQMR;
174
175/*
176 * -----------------------------------------------------------------
177 * PART II: functions exported by sunlinsol_sptfqmr
178 *
179 * CONSTRUCTOR:
180 *    SUNSPTFQMR creates and allocates memory for a SPTFQMR solver
181 *
182 * "SET" ROUTINES:
183 *    SUNSPTFQMRSSetPrecType updates the type of preconditioning to
184 *       use.  Supported values are PREC_NONE, PREC_LEFT, PREC_RIGHT
185 *       and PREC_BOTH.
186 *    SUNSPTFQMRSetMaxl updates the maximum number of iterations to
187 *       allow in the solver.
188 * -----------------------------------------------------------------
189 */
190
191SUNDIALS_EXPORT SUNLinearSolver SUNSPTFQMR(N_Vector y, int pretype, int maxl);
192SUNDIALS_EXPORT int SUNSPTFQMRSetPrecType(SUNLinearSolver S, int pretype);
193SUNDIALS_EXPORT int SUNSPTFQMRSetMaxl(SUNLinearSolver S, int maxl);
194
195/*
196 * -----------------------------------------------------------------
197 * SPTFQMR implementations of various useful linear solver operations
198 * -----------------------------------------------------------------
199 */
200
201SUNDIALS_EXPORT SUNLinearSolver_Type SUNLinSolGetType_SPTFQMR(SUNLinearSolver S);
202SUNDIALS_EXPORT int SUNLinSolInitialize_SPTFQMR(SUNLinearSolver S);
203SUNDIALS_EXPORT int SUNLinSolSetATimes_SPTFQMR(SUNLinearSolver S, void* A_data,
204                                               ATimesFn ATimes);
205SUNDIALS_EXPORT int SUNLinSolSetPreconditioner_SPTFQMR(SUNLinearSolver S,
206                                                       void* P_data,
207                                                       PSetupFn Pset,
208                                                       PSolveFn Psol);
209SUNDIALS_EXPORT int SUNLinSolSetScalingVectors_SPTFQMR(SUNLinearSolver S,
210                                                       N_Vector s1,
211                                                       N_Vector s2);
212SUNDIALS_EXPORT int SUNLinSolSetup_SPTFQMR(SUNLinearSolver S, SUNMatrix A);
213SUNDIALS_EXPORT int SUNLinSolSolve_SPTFQMR(SUNLinearSolver S, SUNMatrix A,
214                                           N_Vector x, N_Vector b, realtype tol);
215SUNDIALS_EXPORT int SUNLinSolNumIters_SPTFQMR(SUNLinearSolver S);
216SUNDIALS_EXPORT realtype SUNLinSolResNorm_SPTFQMR(SUNLinearSolver S);
217SUNDIALS_EXPORT N_Vector SUNLinSolResid_SPTFQMR(SUNLinearSolver S);
218SUNDIALS_EXPORT long int SUNLinSolLastFlag_SPTFQMR(SUNLinearSolver S);
219SUNDIALS_EXPORT int SUNLinSolSpace_SPTFQMR(SUNLinearSolver S,
220                                           long int *lenrwLS,
221                                           long int *leniwLS);
222SUNDIALS_EXPORT int SUNLinSolFree_SPTFQMR(SUNLinearSolver S);
223
224
225#ifdef __cplusplus
226}
227#endif
228
229#endif
Note: See TracBrowser for help on using the repository browser.