Skip to content

Commit a0435a7

Browse files
authored
Removing trapezoid solver (#14286)
- Solver trapezoid was deprecated in #9596 and marked for future removal. - Solver is replaced by GBODE with simulation flags '-s=gbode -gbm=trapezoid -gbctrl=const'
1 parent 2291a39 commit a0435a7

File tree

7 files changed

+12
-37
lines changed

7 files changed

+12
-37
lines changed

OMCompiler/SimulationRuntime/c/simulation/solver/gbode_util.c

Lines changed: 1 addition & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -868,15 +868,14 @@ void replacementString(enum GB_METHOD gbMethod, modelica_boolean constant)
868868
/**
869869
* @brief Display deprecation warning for integration methods replaced by GBODE.
870870
*
871-
* Deprecated methods: impeuler, trapezoid, imprungekutta, rungekuttaSsc
871+
* Deprecated methods: impeuler, imprungekutta, rungekuttaSsc
872872
*
873873
* @param solverMethod Integration method.
874874
*/
875875
void deprecationWarningGBODE(enum SOLVER_METHOD method)
876876
{
877877
switch (method) {
878878
case S_IMPEULER:
879-
case S_TRAPEZOID:
880879
case S_IMPRUNGEKUTTA:
881880
case S_ERKSSC:
882881
break;
@@ -889,9 +888,6 @@ void deprecationWarningGBODE(enum SOLVER_METHOD method)
889888
case S_IMPEULER:
890889
replacementString(RK_IMPL_EULER, TRUE);
891890
break;
892-
case S_TRAPEZOID:
893-
replacementString(RK_TRAPEZOID, TRUE);
894-
break;
895891
case S_IMPRUNGEKUTTA:
896892
replacementString(RK_RADAU_IA_2, TRUE);
897893
break;

OMCompiler/SimulationRuntime/c/simulation/solver/solver_main.c

Lines changed: 2 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -144,7 +144,6 @@ int solver_main_step(DATA* data, threadData_t *threadData, SOLVER_INFO* solverIn
144144
#endif
145145
#ifdef WITH_SUNDIALS
146146
case S_IMPEULER:
147-
case S_TRAPEZOID:
148147
case S_IMPRUNGEKUTTA:
149148
retVal = radau_lobatto_step(data, solverInfo);
150149
if(omc_flag[FLAG_SOLVER_STEPS])
@@ -294,14 +293,12 @@ int initializeSolverData(DATA* data, threadData_t *threadData, SOLVER_INFO* solv
294293
#endif
295294
#ifdef WITH_SUNDIALS
296295
case S_IMPEULER:
297-
case S_TRAPEZOID:
298296
case S_IMPRUNGEKUTTA:
299297
{
300298
int usedImpRKOrder = DEFAULT_IMPRK_ORDER;
301-
if (solverInfo->solverMethod == S_IMPEULER)
299+
if (solverInfo->solverMethod == S_IMPEULER) {
302300
usedImpRKOrder = 1;
303-
if (solverInfo->solverMethod == S_TRAPEZOID)
304-
usedImpRKOrder = 2;
301+
}
305302

306303
/* Check the order if set */
307304
if (omc_flag[FLAG_IMPRK_ORDER])
@@ -399,7 +396,6 @@ int freeSolverData(DATA* data, SOLVER_INFO* solverInfo)
399396
#endif
400397
#ifdef WITH_SUNDIALS
401398
case S_IMPEULER:
402-
case S_TRAPEZOID:
403399
case S_IMPRUNGEKUTTA:
404400
/* free work arrays */
405401
freeKinOde((KINODE*)solverInfo->solverData);
@@ -695,7 +691,6 @@ int solver_main(DATA* data, threadData_t *threadData, const char* init_initMetho
695691
{
696692
#ifndef WITH_SUNDIALS
697693
case S_IMPEULER:
698-
case S_TRAPEZOID:
699694
case S_IMPRUNGEKUTTA:
700695
warningStreamPrint(OMC_LOG_STDOUT, 0, "Sundial/kinsol is needed but not available. Please choose other solver.");
701696
TRACE_POP

OMCompiler/SimulationRuntime/c/util/simulation_options.c

Lines changed: 2 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -224,7 +224,7 @@ const char *FLAG_DESC[FLAG_MAX+1] = {
224224
/* FLAG_IIT */ "[double] value specifies a time for the initialization of the model",
225225
/* FLAG_ILS */ "[int (default 3)] number of lambda steps for homotopy methods",
226226
/* FLAG_IMPRK_ORDER */ "[int (default 5)] value specifies the integration order of the implicit Runge-Kutta method. Valid values: 1-6",
227-
/* FLAG_IMPRK_LS */ "selects the linear solver of the integration methods: impeuler, trapezoid and imprungekuta",
227+
/* FLAG_IMPRK_LS */ "selects the linear solver of the integration methods: impeuler and imprungekuta",
228228
/* FLAG_INITIAL_STEP_SIZE */ "value specifies an initial step size for supported solver",
229229
/* FLAG_INPUT_CSV */ "value specifies an csv-file with inputs for the simulation/optimization of the model",
230230
/* FLAG_INPUT_FILE_STATES */ "value specifies an file with states start values for the optimization of the model",
@@ -439,7 +439,7 @@ const char *FLAG_DETAILED_DESC[FLAG_MAX+1] = {
439439
/* FLAG_IMPRK_ORDER */
440440
" Value specifies the integration order of the implicit Runge-Kutta method. Valid values: 1 to 6. Default order is 5.",
441441
/* FLAG_IMPRK_LS */
442-
" Selects the linear solver of the integration methods impeuler, trapezoid and imprungekuta:\n\n"
442+
" Selects the linear solver of the integration methods impeuler and imprungekuta:\n\n"
443443
" * iterativ - default, sparse iterativ linear solver with fallback case to dense solver\n"
444444
" * dense - dense linear solver, SUNDIALS default method",
445445
/* FLAG_INITIAL_STEP_SIZE */
@@ -1121,7 +1121,6 @@ const char *SOLVER_METHOD_NAME[S_MAX] = {
11211121
/* S_EULER */ "euler",
11221122
/* S_RUNGEKUTTA */ "rungekutta",
11231123
/* S_IMPEULER */ "impeuler",
1124-
/* S_TRAPEZOID */ "trapezoid",
11251124
/* S_IMPRUNGEKUTTA */ "imprungekutta",
11261125
/* S_GBODE */ "gbode",
11271126
/* S_DASSL */ "dassl",

OMCompiler/SimulationRuntime/c/util/simulation_options.h

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -325,7 +325,6 @@ enum SOLVER_METHOD
325325
S_EULER,
326326
S_RUNGEKUTTA,
327327
S_IMPEULER,
328-
S_TRAPEZOID,
329328
S_IMPRUNGEKUTTA,
330329
S_GBODE,
331330
S_DASSL,

doc/UsersGuide/source/solving.rst

Lines changed: 0 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -149,9 +149,6 @@ lot more using the following simulation flags:
149149
old: -s=impeuler
150150
new: -s=gbode -gbm=impl_euler -gbctrl=const
151151
152-
old: -s=trapezoid
153-
new: -s=gbode -gbm=trapezoid -gbctrl=const
154-
155152
old: -s=imprungekutta
156153
new -s=gbode -gbm=(one of the lobatto or radau or gauss RK methods) -gbctrl=const
157154
@@ -179,7 +176,6 @@ with the simflag :ref:`-impRKLS <simflag-imprkls>`. The step-size is
179176
determined as for the basic explicit solvers.
180177

181178
- impeuler - order 1
182-
- trapezoid - order 2
183179
- imprungekutta - Based on Radau IIA and Lobatto IIIA defined by its
184180
Butcher tableau where the order can be adjusted by :ref:`-impRKorder <simflag-imprkorder>`.
185181

testsuite/simulation/modelica/solver/problem1-trapezoid.mos

Lines changed: 3 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -4,7 +4,7 @@
44
// cflags: -d=-newInst
55

66
loadFile("testSolverPackage.mo"); getErrorString();
7-
simulate(testSolver.problem1, stopTime=2e-6, numberOfIntervals=1000, method="trapezoid"); getErrorString();
7+
simulate(testSolver.problem1, stopTime=2e-6, numberOfIntervals=1000, method="gbode", simflags="-gbm=trapezoid -gbctrl=const"); getErrorString();
88

99
res := OpenModelica.Scripting.compareSimulationResults("testSolver.problem1_res.mat",
1010
getEnvironmentVar("REFERENCEFILES")+"/solver/testSolver.problem1.mat",
@@ -48,13 +48,8 @@ getErrorString();
4848
// ""
4949
// record SimulationResult
5050
// resultFile = "testSolver.problem1_res.mat",
51-
// simulationOptions = "startTime = 0.0, stopTime = 2e-06, numberOfIntervals = 1000, tolerance = 1e-06, method = 'trapezoid', fileNamePrefix = 'testSolver.problem1', options = '', outputFormat = 'mat', variableFilter = '.*', cflags = '', simflags = ''",
52-
// messages = "LOG_STDOUT | warning | Integration method 'trapezoid' is deprecated and will be removed in a future version of OpenModelica.
53-
// | | info | | Use integration method GBODE with method 'trapezoid' and constant step size instead:
54-
// | | | | | | Choose integration method 'gbode' in Simulation Setup->General and additional simulation flags '-gbm=trapezoid -gbctrl=const' in Simulation Setup->Simulation Flags.
55-
// | | | | | | or
56-
// | | | | | | Simulation flags '-s=gbode -gbm=trapezoid -gbctrl=const'.
57-
// | | | | | See OpenModelica User's Guide section on GBODE for more details: https://www.openmodelica.org/doc/OpenModelicaUsersGuide/latest/solving.html#gbode
51+
// simulationOptions = "startTime = 0.0, stopTime = 2e-06, numberOfIntervals = 1000, tolerance = 1e-06, method = 'gbode', fileNamePrefix = 'testSolver.problem1', options = '', outputFormat = 'mat', variableFilter = '.*', cflags = '', simflags = '-gbm=trapezoid -gbctrl=const'",
52+
// messages = "LOG_STDOUT | warning | Numerical Jacobians without coloring are currently not supported by GBODE. Colored numerical Jacobian will be used.
5853
// LOG_SUCCESS | info | The initialization finished successfully without homotopy method.
5954
// LOG_SUCCESS | info | The simulation finished successfully.
6055
// "

testsuite/simulation/modelica/solver/problem2-trapezoid.mos

Lines changed: 4 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -1,11 +1,11 @@
1-
// name: problem2-lobatto2
1+
// name: problem2-trapezoid
22
// status: correct
33
// teardown_command: rm -f testSolver.problem2* output.log
44
// cflags: -d=-newInst
55

66
stopTime := 321.8122;
77
loadFile("testSolverPackage.mo"); getErrorString();
8-
simulate(testSolver.problem2, stopTime=stopTime, numberOfIntervals=3000, method="trapezoid", tolerance=1e-12); getErrorString();
8+
simulate(testSolver.problem2, stopTime=stopTime, numberOfIntervals=3000, method="gbode", tolerance=1e-12, simflags="-gbm=trapezoid -gbctrl=const"); getErrorString();
99

1010
res := OpenModelica.Scripting.compareSimulationResults("testSolver.problem2_res.mat",
1111
getEnvironmentVar("REFERENCEFILES")+"/solver/testSolver.problem2.mat",
@@ -56,13 +56,8 @@ val(der(y[8]), stopTime);
5656
// ""
5757
// record SimulationResult
5858
// resultFile = "testSolver.problem2_res.mat",
59-
// simulationOptions = "startTime = 0.0, stopTime = 321.8122, numberOfIntervals = 3000, tolerance = 1e-12, method = 'trapezoid', fileNamePrefix = 'testSolver.problem2', options = '', outputFormat = 'mat', variableFilter = '.*', cflags = '', simflags = ''",
60-
// messages = "LOG_STDOUT | warning | Integration method 'trapezoid' is deprecated and will be removed in a future version of OpenModelica.
61-
// | | info | | Use integration method GBODE with method 'trapezoid' and constant step size instead:
62-
// | | | | | | Choose integration method 'gbode' in Simulation Setup->General and additional simulation flags '-gbm=trapezoid -gbctrl=const' in Simulation Setup->Simulation Flags.
63-
// | | | | | | or
64-
// | | | | | | Simulation flags '-s=gbode -gbm=trapezoid -gbctrl=const'.
65-
// | | | | | See OpenModelica User's Guide section on GBODE for more details: https://www.openmodelica.org/doc/OpenModelicaUsersGuide/latest/solving.html#gbode
59+
// simulationOptions = "startTime = 0.0, stopTime = 321.8122, numberOfIntervals = 3000, tolerance = 1e-12, method = 'gbode', fileNamePrefix = 'testSolver.problem2', options = '', outputFormat = 'mat', variableFilter = '.*', cflags = '', simflags = '-gbm=trapezoid -gbctrl=const'",
60+
// messages = "LOG_STDOUT | warning | Numerical Jacobians without coloring are currently not supported by GBODE. Colored numerical Jacobian will be used.
6661
// LOG_SUCCESS | info | The initialization finished successfully without homotopy method.
6762
// LOG_SUCCESS | info | The simulation finished successfully.
6863
// "

0 commit comments

Comments
 (0)