-
Notifications
You must be signed in to change notification settings - Fork 298
/
nonlinearSolverNewton.c
388 lines (331 loc) · 13 KB
/
nonlinearSolverNewton.c
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
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
/*
* This file is part of OpenModelica.
*
* Copyright (c) 1998-2014, Open Source Modelica Consortium (OSMC),
* c/o Linköpings universitet, Department of Computer and Information Science,
* SE-58183 Linköping, Sweden.
*
* All rights reserved.
*
* THIS PROGRAM IS PROVIDED UNDER THE TERMS OF THE BSD NEW LICENSE OR THE
* GPL VERSION 3 LICENSE OR THE OSMC PUBLIC LICENSE (OSMC-PL) VERSION 1.2.
* ANY USE, REPRODUCTION OR DISTRIBUTION OF THIS PROGRAM CONSTITUTES
* RECIPIENT'S ACCEPTANCE OF THE OSMC PUBLIC LICENSE OR THE GPL VERSION 3,
* ACCORDING TO RECIPIENTS CHOICE.
*
* The OpenModelica software and the OSMC (Open Source Modelica Consortium)
* Public License (OSMC-PL) are obtained from OSMC, either from the above
* address, from the URLs: http://www.openmodelica.org or
* http://www.ida.liu.se/projects/OpenModelica, and in the OpenModelica
* distribution. GNU version 3 is obtained from:
* http://www.gnu.org/copyleft/gpl.html. The New BSD License is obtained from:
* http://www.opensource.org/licenses/BSD-3-Clause.
*
* This program is distributed WITHOUT ANY WARRANTY; without even the implied
* warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE, EXCEPT AS
* EXPRESSLY SET FORTH IN THE BY RECIPIENT SELECTED SUBSIDIARY LICENSE
* CONDITIONS OF OSMC-PL.
*
*/
/*! \file nonlinearSolverNewton.c
*/
#ifdef __cplusplus
extern "C" {
#endif
#include <math.h>
#include <stdlib.h>
#include <string.h> /* memcpy */
#include "simulation/simulation_info_json.h"
#include "util/omc_error.h"
#include "util/varinfo.h"
#include "model_help.h"
#include "nonlinearSystem.h"
#include "nonlinearSolverNewton.h"
#include "newtonIteration.h"
#include "external_input.h"
extern double enorm_(int *n, double *x);
int wrapper_fvec_newton(int* n, double* x, double* fvec, void* userdata, int fj);
#ifdef __cplusplus
extern "C" {
#endif
extern int dgesv_(int *n, int *nrhs, doublereal *a, int *lda, int *ipiv, doublereal *b, int *ldb, int *info);
#ifdef __cplusplus
}
#endif
/*! \fn getAnalyticalJacobian
*
* function calculates analytical jacobian
*
* \param [ref] [data]
* \param [out] [jac]
*
* \author wbraun
*
*/
int getAnalyticalJacobianNewton(DATA* data, threadData_t *threadData, double* jac, int sysNumber)
{
int i,j,k,l,ii,currentSys = sysNumber;
NONLINEAR_SYSTEM_DATA* systemData = &(((DATA*)data)->simulationInfo->nonlinearSystemData[currentSys]);
DATA_NEWTON* solverData = (DATA_NEWTON*)(systemData->solverData);
const int index = systemData->jacobianIndex;
ANALYTIC_JACOBIAN* jacobian = &(data->simulationInfo->analyticJacobians[systemData->jacobianIndex]);
memset(jac, 0, (solverData->n)*(solverData->n)*sizeof(double));
for(i=0; i < jacobian->sparsePattern->maxColors; i++)
{
/* activate seed variable for the corresponding color */
for(ii=0; ii < jacobian->sizeCols; ii++)
if(jacobian->sparsePattern->colorCols[ii]-1 == i)
jacobian->seedVars[ii] = 1;
systemData->analyticalJacobianColumn(data, threadData, jacobian, NULL);
for(j = 0; j < jacobian->sizeCols; j++)
{
if(jacobian->seedVars[j] == 1)
{
ii = jacobian->sparsePattern->leadindex[j];
while(ii < jacobian->sparsePattern->leadindex[j+1])
{
l = jacobian->sparsePattern->index[ii];
k = j*jacobian->sizeRows + l;
jac[k] = jacobian->resultVars[l];
ii++;
};
}
/* de-activate seed variable for the corresponding color */
if(jacobian->sparsePattern->colorCols[j]-1 == i)
jacobian->seedVars[j] = 0;
}
}
return 0;
}
/*! \fn wrapper_fvec_newton for the residual Function
* tensolve calls for the subroutine fcn(n, x, fvec, iflag, data)
*
* fj decides whether the function values or the jacobian matrix shall be calculated
* fj = 1 ==> calculate function values
* fj = 0 ==> calculate jacobian matrix
*/
int wrapper_fvec_newton(int* n, double* x, double* fvec, void* userdata, int fj)
{
DATA_USER* uData = (DATA_USER*) userdata;
DATA* data = (DATA*)(uData->data);
void *dataAndThreadData[2] = {data, uData->threadData};
int currentSys = ((DATA_USER*)userdata)->sysNumber;
NONLINEAR_SYSTEM_DATA* systemData = &(data->simulationInfo->nonlinearSystemData[currentSys]);
DATA_NEWTON* solverData = (DATA_NEWTON*)(systemData->solverData);
int flag = 1;
int *iflag=&flag;
if (fj) {
(data->simulationInfo->nonlinearSystemData[currentSys].residualFunc)(dataAndThreadData, x, fvec, iflag);
} else {
/* performance measurement */
rt_ext_tp_tick(&systemData->jacobianTimeClock);
if(systemData->jacobianIndex != -1) {
getAnalyticalJacobianNewton(data, uData->threadData, solverData->fjac, currentSys);
} else {
double delta_h = sqrt(solverData->epsfcn);
double delta_hh;
double xsave;
int i,j,l, linear=0;
for(i = 0; i < *n; i++) {
delta_hh = fmax(delta_h * fmax(fabs(x[i]), fabs(fvec[i])), delta_h);
delta_hh = ((fvec[i] >= 0) ? delta_hh : -delta_hh);
delta_hh = x[i] + delta_hh - x[i];
xsave = x[i];
x[i] += delta_hh;
delta_hh = 1. / delta_hh;
wrapper_fvec_newton(n, x, solverData->rwork, userdata, 1);
solverData->nfev++;
for(j = 0; j < *n; j++) {
l = i * *n + j;
solverData->fjac[l] = (solverData->rwork[j] - fvec[j]) * delta_hh;
}
x[i] = xsave;
}
}
/* performance measurement and statistics */
systemData->jacobianTime += rt_ext_tp_tock(&(systemData->jacobianTimeClock));
systemData->numberOfJEval++;
}
return *iflag;
}
/*! \fn solve non-linear system with newton method
*
* \param [in] [data]
* [sysNumber] index of the corresponding non-linear system
*
* \author wbraun
*/
int solveNewton(DATA *data, threadData_t *threadData, int sysNumber)
{
NONLINEAR_SYSTEM_DATA* systemData = &(data->simulationInfo->nonlinearSystemData[sysNumber]);
DATA_NEWTON* solverData = (DATA_NEWTON*)(systemData->solverData);
int eqSystemNumber = 0;
int i;
double xerror = -1, xerror_scaled = -1;
int success = 0;
int nfunc_evals = 0;
int continuous = 1;
double local_tol = solverData->ftol;
int giveUp = 0;
int retries = 0;
int retries2 = 0;
int nonContinuousCase = 0;
modelica_boolean *relationsPreBackup = NULL;
int casualTearingSet = data->simulationInfo->nonlinearSystemData[sysNumber].strictTearingFunctionCall != NULL;
DATA_USER userdata = {.data = (void*)data, .threadData = threadData, .sysNumber = sysNumber};
/*
* We are given the number of the non-linear system.
* We want to look it up among all equations.
*/
eqSystemNumber = systemData->equationIndex;
relationsPreBackup = (modelica_boolean*) malloc(data->modelData->nRelations*sizeof(modelica_boolean));
solverData->nfev = 0;
/* try to calculate jacobian only once at the beginning of the iteration */
solverData->calculate_jacobian = 0;
// Initialize lambda variable
if (data->simulationInfo->nonlinearSystemData[sysNumber].homotopySupport) {
solverData->x[solverData->n] = 1.0;
solverData->x_new[solverData->n] = 1.0;
}
else {
solverData->x[solverData->n] = 0.0;
solverData->x_new[solverData->n] = 0.0;
}
/* debug output */
if(ACTIVE_STREAM(LOG_NLS_V))
{
int indexes[2] = {1,eqSystemNumber};
infoStreamPrintWithEquationIndexes(LOG_NLS_V, 1, indexes, "Start solving Non-Linear System %d at time %g with Newton Solver",
eqSystemNumber, data->localData[0]->timeValue);
for(i = 0; i < solverData->n; i++) {
infoStreamPrint(LOG_NLS_V, 1, "x[%d] = %.15e", i, data->simulationInfo->discreteCall ? systemData->nlsx[i] : systemData->nlsxExtrapolation[i]);
infoStreamPrint(LOG_NLS_V, 0, "nominal = %g +++ nlsx = %g +++ old = %g +++ extrapolated = %g",
systemData->nominal[i], systemData->nlsx[i], systemData->nlsxOld[i], systemData->nlsxExtrapolation[i]);
messageClose(LOG_NLS_V);
}
messageClose(LOG_NLS_V);
}
/* set x vector */
if(data->simulationInfo->discreteCall) {
memcpy(solverData->x, systemData->nlsx, solverData->n*(sizeof(double)));
} else {
memcpy(solverData->x, systemData->nlsxExtrapolation, solverData->n*(sizeof(double)));
}
/* start solving loop */
while(!giveUp && !success)
{
giveUp = 1;
solverData->newtonStrategy = data->simulationInfo->newtonStrategy;
_omc_newton(wrapper_fvec_newton, solverData, &userdata);
/* check for proper inputs */
if(solverData->info == 0)
printErrorEqSyst(IMPROPER_INPUT, modelInfoGetEquation(&data->modelData->modelDataXml,eqSystemNumber), data->localData[0]->timeValue);
/* reset non-contunuousCase */
if(nonContinuousCase && xerror > local_tol && xerror_scaled > local_tol)
{
memcpy(data->simulationInfo->relationsPre, relationsPreBackup, sizeof(modelica_boolean)*data->modelData->nRelations);
nonContinuousCase = 0;
}
/* check for error */
xerror_scaled = enorm_(&solverData->n, solverData->fvecScaled);
xerror = enorm_(&solverData->n, solverData->fvec);
/* solution found */
if((xerror <= local_tol || xerror_scaled <= local_tol) && solverData->info > 0)
{
success = 1;
nfunc_evals += solverData->nfev;
if(ACTIVE_STREAM(LOG_NLS_V))
{
infoStreamPrint(LOG_NLS_V, 0, "*** System solved ***\n%d restarts", retries);
infoStreamPrint(LOG_NLS_V, 0, "nfunc = %d +++ error = %.15e +++ error_scaled = %.15e", nfunc_evals, xerror, xerror_scaled);
for(i = 0; i < solverData->n; i++)
infoStreamPrint(LOG_NLS_V, 0, "x[%d] = %.15e\n\tresidual = %e", i, solverData->x[i], solverData->fvec[i]);
}
/* take the solution */
memcpy(systemData->nlsx, solverData->x, solverData->n*(sizeof(double)));
/* Then try with old values (instead of extrapolating )*/
}
// If this is the casual tearing set (only exists for dynamic tearing), break after first try
else if(retries < 1 && casualTearingSet)
{
giveUp = 1;
infoStreamPrint(LOG_NLS_V, 0, "### No Solution for the casual tearing set at the first try! ###");
}
else if(retries < 1)
{
memcpy(solverData->x, systemData->nlsxOld, solverData->n*(sizeof(double)));
retries++;
giveUp = 0;
nfunc_evals += solverData->nfev;
infoStreamPrint(LOG_NLS_V, 0, " - iteration making no progress:\t try old values.");
/* try to vary the initial values */
/* evaluate jacobian in every step now */
solverData->calculate_jacobian = 1;
}
else if(retries < 2)
{
for(i = 0; i < solverData->n; i++)
solverData->x[i] += systemData->nominal[i] * 0.01;
retries++;
giveUp = 0;
nfunc_evals += solverData->nfev;
infoStreamPrint(LOG_NLS_V, 0, " - iteration making no progress:\t vary solution point by 1%%.");
/* try to vary the initial values */
}
else if(retries < 3)
{
for(i = 0; i < solverData->n; i++)
solverData->x[i] = systemData->nominal[i];
retries++;
giveUp = 0;
nfunc_evals += solverData->nfev;
infoStreamPrint(LOG_NLS_V, 0, " - iteration making no progress:\t try nominal values as initial solution.");
}
else if(retries < 4 && data->simulationInfo->discreteCall)
{
/* try to solve non-continuous
* work-a-round: since other wise some model does
* stuck in event iteration. e.g.: Modelica.Mechanics.Rotational.Examples.HeatLosses
*/
memcpy(solverData->x, systemData->nlsxOld, solverData->n*(sizeof(double)));
retries++;
/* try to solve a discontinuous system */
continuous = 0;
nonContinuousCase = 1;
memcpy(relationsPreBackup, data->simulationInfo->relationsPre, sizeof(modelica_boolean)*data->modelData->nRelations);
giveUp = 0;
nfunc_evals += solverData->nfev;
infoStreamPrint(LOG_NLS_V, 0, " - iteration making no progress:\t try to solve a discontinuous system.");
}
else if(retries2 < 4)
{
memcpy(solverData->x, systemData->nlsxOld, solverData->n*(sizeof(double)));
/* reduce tolarance */
local_tol = local_tol*10;
retries = 0;
retries2++;
giveUp = 0;
nfunc_evals += solverData->nfev;
infoStreamPrint(LOG_NLS_V, 0, " - iteration making no progress:\t reduce the tolerance slightly to %e.", local_tol);
}
else
{
printErrorEqSyst(ERROR_AT_TIME, modelInfoGetEquation(&data->modelData->modelDataXml,eqSystemNumber), data->localData[0]->timeValue);
if(ACTIVE_STREAM(LOG_NLS_V))
{
infoStreamPrint(LOG_NLS_V, 0, "### No Solution! ###\n after %d restarts", retries);
infoStreamPrint(LOG_NLS_V, 0, "nfunc = %d +++ error = %.15e +++ error_scaled = %.15e", nfunc_evals, xerror, xerror_scaled);
if(ACTIVE_STREAM(LOG_NLS_V))
for(i = 0; i < solverData->n; i++)
infoStreamPrint(LOG_NLS_V, 0, "x[%d] = %.15e\n\tresidual = %e", i, solverData->x[i], solverData->fvec[i]);
}
}
}
if(ACTIVE_STREAM(LOG_NLS_V))
messageClose(LOG_NLS_V);
free(relationsPreBackup);
/* write statistics */
systemData->numberOfFEval = solverData->numberOfFunctionEvaluations;
systemData->numberOfIterations = solverData->numberOfIterations;
return success;
}