Skip to content
This repository was archived by the owner on May 18, 2019. It is now read-only.

Commit b4aa137

Browse files
sjoelundOpenModelica-Hudson
authored andcommitted
[FMI] Set default dense and sparse linear solvers
Belonging to [master]: - #1907
1 parent c1445a8 commit b4aa137

File tree

5 files changed

+30
-18
lines changed

5 files changed

+30
-18
lines changed

SimulationRuntime/c/simulation/simulation_runtime.cpp

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -268,8 +268,9 @@ static int getlinearSparseSolverMethod()
268268
const char *cflags = omc_flagValue[FLAG_LSS];
269269
const string *method = cflags ? new string(cflags) : NULL;
270270

271-
if(!method)
272-
return LSS_KLU; /* default method */
271+
if (!method) {
272+
return LSS_DEFAULT; /* default method */
273+
}
273274

274275
for(i=1; i<LSS_MAX; ++i)
275276
if(*method == LSS_NAME[i])

SimulationRuntime/c/simulation/solver/linearSystem.c

Lines changed: 22 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -120,36 +120,47 @@ int initializeLinearSystems(DATA *data, threadData_t *threadData)
120120
{
121121
switch(data->simulationInfo->lssMethod)
122122
{
123-
#if !defined(OMC_MINIMAL_RUNTIME)
124-
case LSS_LIS:
125-
linsys[i].setAElement = setAElementLis;
126-
linsys[i].setBElement = setBElementLis;
127-
allocateLisData(size, size, nnz, linsys[i].solverData);
128-
break;
129-
#endif
130123
#ifdef WITH_UMFPACK
131124
case LSS_UMFPACK:
132125
linsys[i].setAElement = setAElementUmfpack;
133126
linsys[i].setBElement = setBElement;
134127
allocateUmfPackData(size, size, nnz, linsys[i].solverData);
135128
break;
129+
case LSS_DEFAULT:
136130
case LSS_KLU:
137131
linsys[i].setAElement = setAElementKlu;
138132
linsys[i].setBElement = setBElement;
139133
allocateKluData(size, size, nnz, linsys[i].solverData);
140134
break;
141135
#else
136+
case LSS_KLU:
142137
case LSS_UMFPACK:
143-
throwStreamPrint(threadData, "OMC is compiled without UMFPACK, if you want use umfpack please compile OMC with UMFPACK.");
138+
throwStreamPrint(threadData, "OMC is compiled without UMFPACK, if you want use klu or umfpack please compile OMC with UMFPACK.");
144139
break;
145140
#endif
146-
141+
#if !defined(OMC_MINIMAL_RUNTIME)
142+
#if !defined(WITH_UMFPACK)
143+
case LSS_DEFAULT:
144+
#endif
145+
case LSS_LIS:
146+
linsys[i].setAElement = setAElementLis;
147+
linsys[i].setBElement = setBElementLis;
148+
allocateLisData(size, size, nnz, linsys[i].solverData);
149+
break;
150+
#endif
151+
#if defined(OMC_MINIMAL_RUNTIME) && !defined(WITH_UMFPACK)
152+
case LSS_DEFAULT:
153+
{
154+
int indexes[2] = {1, linsys[i].index};
155+
infoStreamPrintWithEquationIndexes(LOG_STDOUT, 0, indexes, "The simulation runtime does not have access to sparse solvers. Defaulting to a dense linear system solver instead.");
156+
break;
157+
}
158+
#endif
147159
default:
148160
throwStreamPrint(threadData, "unrecognized linear solver");
149161
}
150162
}
151-
152-
else{
163+
if(linsys[i].useSparseSolver == 0) { /* Not an else-statement because there might not be a sparse linear solver available */
153164
switch(data->simulationInfo->lsMethod)
154165
{
155166
case LS_LAPACK:

SimulationRuntime/c/simulation/solver/model_help.c

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -940,8 +940,8 @@ void initializeDataStruc(DATA *data, threadData_t *threadData)
940940
data->simulationInfo->nlsMethod = NLS_HOMOTOPY;
941941
#endif
942942
data->simulationInfo->nlsLinearSolver = NLS_LS_LAPACK;
943-
data->simulationInfo->lsMethod = LS_LAPACK;
944-
data->simulationInfo->lssMethod = LS_UMFPACK;
943+
data->simulationInfo->lsMethod = LS_DEFAULT;
944+
data->simulationInfo->lssMethod = LSS_DEFAULT;
945945
data->simulationInfo->mixedMethod = MIXED_SEARCH;
946946
data->simulationInfo->newtonStrategy = NEWTON_PURE;
947947
data->simulationInfo->nlsCsvInfomation = 0;

SimulationRuntime/c/util/simulation_options.c

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -679,7 +679,7 @@ const char *LS_DESC[LS_MAX] = {
679679

680680
const char *LSS_NAME[LS_MAX] = {
681681
"LS_UNKNOWN",
682-
682+
"default",
683683
#if !defined(OMC_MINIMAL_RUNTIME)
684684
/* LS_LIS */ "lis",
685685
#else
@@ -691,7 +691,7 @@ const char *LSS_NAME[LS_MAX] = {
691691

692692
const char *LSS_DESC[LS_MAX] = {
693693
"unknown",
694-
694+
"the default sparse linear solver (or a dense solver if there is none available)"
695695
#if !defined(OMC_MINIMAL_RUNTIME)
696696
/* LS_LIS */ "method using iterative solver Lis",
697697
#else

SimulationRuntime/c/util/simulation_options.h

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -214,14 +214,14 @@ enum LINEAR_SPARSE_SOLVER
214214
{
215215
LSS_NONE = 0,
216216

217+
LSS_DEFAULT,
217218
#if !defined(OMC_MINIMAL_RUNTIME)
218219
LSS_LIS,
219220
#else
220221
LSS_LIS_NOT_AVAILABLE,
221222
#endif
222223
LSS_KLU,
223224
LSS_UMFPACK,
224-
225225
LSS_MAX
226226
};
227227

0 commit comments

Comments
 (0)