File: sunadaptcontroller_soderlind.h

package info (click to toggle)
sundials 7.1.1%2Bdfsg1-10
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 83,788 kB
  • sloc: ansic: 264,746; f90: 120,987; cpp: 69,481; python: 5,134; sh: 2,483; makefile: 266; perl: 123
file content (132 lines) | stat: -rw-r--r-- 4,599 bytes parent folder | download
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
/* -----------------------------------------------------------------
 * Programmer(s): Daniel R. Reynolds @ SMU
 * -----------------------------------------------------------------
 * SUNDIALS Copyright Start
 * Copyright (c) 2002-2024, Lawrence Livermore National Security
 * and Southern Methodist University.
 * All rights reserved.
 *
 * See the top-level LICENSE and NOTICE files for details.
 *
 * SPDX-License-Identifier: BSD-3-Clause
 * SUNDIALS Copyright End
 * -----------------------------------------------------------------
 * This is the header file for the SUNAdaptController_Soderlind
 * module.
 * -----------------------------------------------------------------*/

#ifndef _SUNADAPTCONTROLLER_SODERLIND_H
#define _SUNADAPTCONTROLLER_SODERLIND_H

#include <stdio.h>
#include <sundials/sundials_adaptcontroller.h>

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

/* ----------------------------------------------------
 * Soderlind implementation of SUNAdaptController
 * ---------------------------------------------------- */

struct _SUNAdaptControllerContent_Soderlind
{
  sunrealtype k1; /* internal controller parameters */
  sunrealtype k2;
  sunrealtype k3;
  sunrealtype k4;
  sunrealtype k5;
  sunrealtype bias; /* error bias factor */
  sunrealtype ep;   /* error from previous step */
  sunrealtype epp;  /* error from 2 steps ago */
  sunrealtype hp;   /* previous step size */
  sunrealtype hpp;  /* step size from 2 steps ago */
  int firststeps;   /* flag to handle first few steps */
};

typedef struct _SUNAdaptControllerContent_Soderlind* SUNAdaptControllerContent_Soderlind;

/* ------------------
 * Exported Functions
 * ------------------ */

SUNDIALS_EXPORT
SUNAdaptController SUNAdaptController_Soderlind(SUNContext sunctx);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetParams_Soderlind(SUNAdaptController C,
                                                  sunrealtype k1, sunrealtype k2,
                                                  sunrealtype k3, sunrealtype k4,
                                                  sunrealtype k5);

SUNDIALS_EXPORT
SUNAdaptController_Type SUNAdaptController_GetType_Soderlind(SUNAdaptController C);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_EstimateStep_Soderlind(SUNAdaptController C,
                                                     sunrealtype h, int p,
                                                     sunrealtype dsm,
                                                     sunrealtype* hnew);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_Reset_Soderlind(SUNAdaptController C);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetDefaults_Soderlind(SUNAdaptController C);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_Write_Soderlind(SUNAdaptController C, FILE* fptr);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetErrorBias_Soderlind(SUNAdaptController C,
                                                     sunrealtype bias);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_UpdateH_Soderlind(SUNAdaptController C,
                                                sunrealtype h, sunrealtype dsm);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_Space_Soderlind(SUNAdaptController C,
                                              long int* lenrw, long int* leniw);

/* Convenience routines to construct subsidiary controllers */

SUNDIALS_EXPORT
SUNAdaptController SUNAdaptController_PID(SUNContext sunctx);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetParams_PID(SUNAdaptController C, sunrealtype k1,
                                            sunrealtype k2, sunrealtype k3);

SUNDIALS_EXPORT
SUNAdaptController SUNAdaptController_PI(SUNContext sunctx);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetParams_PI(SUNAdaptController C, sunrealtype k1,
                                           sunrealtype k2);

SUNDIALS_EXPORT
SUNAdaptController SUNAdaptController_I(SUNContext sunctx);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetParams_I(SUNAdaptController C, sunrealtype k1);

SUNDIALS_EXPORT
SUNAdaptController SUNAdaptController_ExpGus(SUNContext sunctx);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetParams_ExpGus(SUNAdaptController C,
                                               sunrealtype k1, sunrealtype k2);

SUNDIALS_EXPORT
SUNAdaptController SUNAdaptController_ImpGus(SUNContext sunctx);

SUNDIALS_EXPORT
SUNErrCode SUNAdaptController_SetParams_ImpGus(SUNAdaptController C,
                                               sunrealtype k1, sunrealtype k2);

#ifdef __cplusplus
}
#endif

#endif /* _SUNADAPTCONTROLLER_SODERLIND_H */