-
Notifications
You must be signed in to change notification settings - Fork 296
/
ida_solver.c
1940 lines (1665 loc) · 60.4 KB
/
ida_solver.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
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
* This file is part of OpenModelica.
*
* Copyright (c) 1998-2016, 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 ida_solver.c
*/
#include <string.h>
#include <setjmp.h>
#include "omc_config.h"
#include "openmodelica.h"
#include "openmodelica_func.h"
#include "simulation_data.h"
#include "util/omc_error.h"
#include "gc/omc_gc.h"
#include "simulation/options.h"
#include "simulation/simulation_runtime.h"
#include "simulation/results/simulation_result.h"
#include "simulation/solver/solver_main.h"
#include "simulation/solver/model_help.h"
#include "simulation/solver/external_input.h"
#include "simulation/solver/epsilon.h"
#include "simulation/solver/omc_math.h"
#include "simulation/solver/ida_solver.h"
#include "simulation/solver/dassl.h"
#ifdef WITH_SUNDIALS
/* readability */
#define SCALE_MODE 0
#define RESCALE_MODE 1
#define MINIMAL_SCALE_FACTOR 1e-8
#include <sundials/sundials_nvector.h>
#include <nvector/nvector_serial.h>
#include <idas/idas.h>
#include <idas/idas_dense.h>
#include <idas/idas_klu.h>
#include <idas/idas_spgmr.h>
#include <idas/idas_spbcgs.h>
#include <idas/idas_sptfqmr.h>
static int callDenseJacobian(long int Neq, double tt, double cj,
N_Vector yy, N_Vector yp, N_Vector rr,
DlsMat Jac, void *user_data,
N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
static int callSparseJacobian(realtype tt, realtype cj,
N_Vector yy, N_Vector yp, N_Vector rr, SlsMat Jac, void *user_data,
N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
static int residualFunctionIDA(double time, N_Vector yy, N_Vector yp, N_Vector res, void* userData);
int rootsFunctionIDA(double time, N_Vector yy, N_Vector yp, double *gout, void* userData);
static int getScalingFactors(DATA* data, IDA_SOLVER *idaData, SlsMat scaleMatrix);
static int idaScaleData(IDA_SOLVER *idaData);
static int idaReScaleData(IDA_SOLVER *idaData);
static int idaScaleVector(N_Vector vec, double* factors, unsigned int size);
static int idaReScaleVector(N_Vector vec, double* factors, unsigned int size);
static IDA_SOLVER *idaDataGlobal;
static int initializedSolver = 0;
int ida_event_update(DATA* data, threadData_t *threadData);
int checkIDAflag(int flag)
{
TRACE_PUSH
int retVal;
switch(flag)
{
case IDA_SUCCESS:
case IDA_TSTOP_RETURN:
retVal = 0;
break;
default:
retVal = 1;
break;
}
TRACE_POP
return retVal;
}
void errOutputIDA(int error_code, const char *module, const char *function,
char *msg, void *userData)
{
TRACE_PUSH
DATA* data = (DATA*)(((IDA_USERDATA*)((IDA_SOLVER*)userData)->simData)->data);
infoStreamPrint(LOG_SOLVER, 1, "#### IDA error message #####");
infoStreamPrint(LOG_SOLVER, 0, " -> error code %d\n -> module %s\n -> function %s", error_code, module, function);
infoStreamPrint(LOG_SOLVER, 0, " Message: %s", msg);
messageClose(LOG_SOLVER);
TRACE_POP
}
/* initial main ida data */
int
ida_solver_initial(DATA* data, threadData_t *threadData, SOLVER_INFO* solverInfo, IDA_SOLVER *idaData){
TRACE_PUSH
int flag;
long int i;
double* tmp;
/* default max order */
int maxOrder = 5;
/* sim data */
idaData->simData = (IDA_USERDATA*)malloc(sizeof(IDA_USERDATA));
idaData->simData->data = data;
idaData->simData->threadData = threadData;
/* initialize constants */
idaData->setInitialSolution = 0;
/* start initialization routines of sundials */
idaData->ida_mem = IDACreate();
if (idaData->ida_mem == NULL){
throwStreamPrint(threadData, "##IDA## Initialization of IDA solver failed!");
}
idaData->daeMode = 0;
idaData->residualFunction = residualFunctionIDA;
idaData->N = (long int)data->modelData->nStates;
/* change parameter for DAE mode */
if (omc_flag[FLAG_DAE_MODE])
{
if (compiledInDAEMode)
{
idaData->daeMode = 1;
idaData->N = data->modelData->nStates + data->simulationInfo->daeModeData->nAlgebraicDAEVars;
}
else
{
warningStreamPrint(LOG_STDOUT, 0, "-daeMode flag is used, the model is not compiled in DAE mode. See compiler flag: +daeMode. Continue as usual.");
}
}
/* initialize states and der(states) */
if (idaData->daeMode)
{
idaData->states = (double*) malloc(idaData->N*sizeof(double));
idaData->statesDer = (double*) calloc(idaData->N,sizeof(double));
memcpy(idaData->states, data->localData[0]->realVars, sizeof(double)*data->modelData->nStates);
// and also algebraic vars
data->simulationInfo->daeModeData->getAlgebraicDAEVars(data, threadData, idaData->states + data->modelData->nStates);
memcpy(idaData->statesDer, data->localData[0]->realVars + data->modelData->nStates, sizeof(double)*data->modelData->nStates);
idaData->y = N_VMake_Serial(idaData->N, idaData->states);
idaData->yp = N_VMake_Serial(idaData->N, idaData->statesDer);
}
else
{
idaData->y = N_VMake_Serial(idaData->N, data->localData[0]->realVars);
idaData->yp = N_VMake_Serial(idaData->N, data->localData[0]->realVars + data->modelData->nStates);
}
flag = IDAInit(idaData->ida_mem,
idaData->residualFunction,
data->simulationInfo->startTime,
idaData->y,
idaData->yp);
/* allocate memory for jacobians calculation */
idaData->ysave = (double*) malloc(idaData->N*sizeof(double));
idaData->ypsave = (double*) malloc(idaData->N*sizeof(double));
idaData->delta_hh = (double*) malloc(idaData->N*sizeof(double));
idaData->errwgt = N_VNew_Serial(idaData->N);
idaData->newdelta = N_VNew_Serial(idaData->N);
/* allocate memory for initialization process */
tmp = (double*) malloc(idaData->N*sizeof(double));
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Something goes wrong while initialize IDA solver!");
}
flag = IDASetUserData(idaData->ida_mem, idaData);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Something goes wrong while initialize IDA solver!");
}
flag = IDASetErrHandlerFn(idaData->ida_mem, errOutputIDA, idaData);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Something goes wrong while set error handler!");
}
/* set nominal values of the states for absolute tolerances */
infoStreamPrint(LOG_SOLVER, 1, "The relative tolerance is %g. Following absolute tolerances are used for the states: ", data->simulationInfo->tolerance);
for(i=0; i < data->modelData->nStates; ++i)
{
tmp[i] = fmax(fabs(data->modelData->realVarsData[i].attribute.nominal), 1e-32);
infoStreamPrint(LOG_SOLVER_V, 0, "%ld. %s -> %g", i+1, data->modelData->realVarsData[i].info.name, tmp[i]);
}
/* daeMode: set nominal values for algebraic variables */
if (idaData->daeMode)
{
data->simulationInfo->daeModeData->getAlgebraicDAEVarNominals(data, threadData, tmp + data->modelData->nStates);
}
/* multiply by tolerance to obtain a relative tolerace */
for(i=0; i < idaData->N; ++i)
{
tmp[i] *= data->simulationInfo->tolerance;
}
messageClose(LOG_SOLVER);
flag = IDASVtolerances(idaData->ida_mem,
data->simulationInfo->tolerance,
N_VMake_Serial(idaData->N,tmp));
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Setting tolerances fails while initialize IDA solver!");
}
if (omc_flag[FLAG_IDA_SCALING]) // idaNoScaling
{
/* allocate memory for scaling */
idaData->yScale = (double*) malloc(idaData->N*sizeof(double));
idaData->ypScale = (double*) malloc(idaData->N*sizeof(double));
idaData->resScale = (double*) malloc(idaData->N*sizeof(double));
/* set yScale from nominal values */
for(i=0; i < data->modelData->nStates; ++i)
{
idaData->yScale[i] = fabs(data->modelData->realVarsData[i].attribute.nominal);
idaData->ypScale[i] = 1.0;
}
/* daeMode: set nominal values for algebraic variables */
if (idaData->daeMode)
{
data->simulationInfo->daeModeData->getAlgebraicDAEVarNominals(data, threadData, idaData->yScale + data->modelData->nStates);
for(i=data->modelData->nStates; i < idaData->N; ++i)
{
idaData->ypScale[i] = 1.0;
}
}
infoStreamPrint(LOG_SOLVER_V, 1, "The scale factors for all ida states: ");
for(i=0; i < idaData->N; ++i)
{
infoStreamPrint(LOG_SOLVER_V, 0, "%ld. scaleFactor: %g", i+1, idaData->yScale[i]);
}
messageClose(LOG_SOLVER_V);
}
/* initialize */
idaData->disableScaling = 0;
/* set root function */
/* if FLAG_NO_ROOTFINDING is set, solve perform without internal root finding*/
if(!omc_flag[FLAG_NO_ROOTFINDING])
{
solverInfo->solverRootFinding = 1;
flag = IDARootInit(idaData->ida_mem, data->modelData->nZeroCrossings, rootsFunctionIDA);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Setting root function fails while initialize IDA solver!");
}
}
infoStreamPrint(LOG_SOLVER, 0, "ida uses internal root finding method %s", solverInfo->solverRootFinding?"YES":"NO");
/* define maximum integration order of dassl */
if (omc_flag[FLAG_MAX_ORDER])
{
maxOrder = atoi(omc_flagValue[FLAG_MAX_ORDER]);
flag = IDASetMaxOrd(idaData->ida_mem, maxOrder);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Failed to set max integration order!");
}
}
infoStreamPrint(LOG_SOLVER, 0, "maximum integration order %d", maxOrder);
/* if FLAG_NOEQUIDISTANT_GRID is set, choose ida step method */
if (omc_flag[FLAG_NOEQUIDISTANT_GRID])
{
idaData->internalSteps = 1; /* TRUE */
solverInfo->solverNoEquidistantGrid = 1;
}
else
{
idaData->internalSteps = 0; /* FALSE */
}
infoStreamPrint(LOG_SOLVER, 0, "use equidistant time grid %s", idaData->internalSteps?"NO":"YES");
/* check if Flags FLAG_NOEQUIDISTANT_OUT_FREQ or FLAG_NOEQUIDISTANT_OUT_TIME are set */
if (idaData->internalSteps){
if (omc_flag[FLAG_NOEQUIDISTANT_OUT_FREQ])
{
idaData->stepsFreq = atoi(omc_flagValue[FLAG_NOEQUIDISTANT_OUT_FREQ]);
}
else if (omc_flag[FLAG_NOEQUIDISTANT_OUT_TIME])
{
idaData->stepsTime = atof(omc_flagValue[FLAG_NOEQUIDISTANT_OUT_TIME]);
flag = IDASetMaxStep(idaData->ida_mem, idaData->stepsTime);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Setting max steps of the IDA solver!");
}
infoStreamPrint(LOG_SOLVER, 0, "maximum step size %g", idaData->stepsTime);
} else {
idaData->stepsFreq = 1;
idaData->stepsTime = 0.0;
}
if (omc_flag[FLAG_NOEQUIDISTANT_OUT_FREQ] && omc_flag[FLAG_NOEQUIDISTANT_OUT_TIME]){
warningStreamPrint(LOG_STDOUT, 0, "The flags are \"noEquidistantOutputFrequency\" "
"and \"noEquidistantOutputTime\" are in opposition "
"to each other. The flag \"noEquidistantOutputFrequency\" superiors.");
}
infoStreamPrint(LOG_SOLVER, 0, "as the output frequency control is used: %d", idaData->stepsFreq);
infoStreamPrint(LOG_SOLVER, 0, "as the output frequency time step control is used: %f", idaData->stepsTime);
}
/* if FLAG_IDA_LS is set, choose ida linear solver method */
if (omc_flag[FLAG_IDA_LS])
{
for(i=1; i< IDA_LS_MAX;i++)
{
if(!strcmp((const char*)omc_flagValue[FLAG_IDA_LS], IDA_LS_METHOD[i])){
idaData->linearSolverMethod = (int)i;
break;
}
}
if(idaData->linearSolverMethod == IDA_LS_UNKNOWN)
{
if (ACTIVE_WARNING_STREAM(LOG_SOLVER))
{
warningStreamPrint(LOG_SOLVER, 1, "unrecognized ida linear solver method %s, current options are:", (const char*)omc_flagValue[FLAG_IDA_LS]);
for(i=1; i < IDA_LS_MAX; ++i)
{
warningStreamPrint(LOG_SOLVER, 0, "%-15s [%s]", IDA_LS_METHOD[i], IDA_LS_METHOD_DESC[i]);
}
messageClose(LOG_SOLVER);
}
throwStreamPrint(threadData,"unrecognized ida linear solver method %s", (const char*)omc_flagValue[FLAG_IDA_LS]);
}
}
else
{
idaData->linearSolverMethod = IDA_LS_KLU;
}
switch (idaData->linearSolverMethod){
case IDA_LS_SPGMR:
flag = IDASpgmr(idaData->ida_mem, 0);
break;
case IDA_LS_SPBCG:
flag = IDASpbcg(idaData->ida_mem, 0);
break;
case IDA_LS_SPTFQMR:
flag = IDASptfqmr(idaData->ida_mem, 0);
break;
case IDA_LS_DENSE:
flag = IDADense(idaData->ida_mem, idaData->N);
break;
case IDA_LS_KLU:
/* Set KLU after initialized sparse pattern of the jacobian for nnz */
break;
default:
throwStreamPrint(threadData,"unrecognized linear solver method %s", (const char*)omc_flagValue[FLAG_IDA_LS]);
break;
}
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Setting linear solver method fails while initialize IDA solver!");
} else {
infoStreamPrint(LOG_SOLVER, 0, "IDA linear solver method selected %s", IDA_LS_METHOD_DESC[idaData->linearSolverMethod]);
}
/* if FLAG_JACOBIAN is set, choose jacobian calculation method */
if (omc_flag[FLAG_JACOBIAN])
{
for(i=1; i< JAC_MAX;i++)
{
if(!strcmp((const char*)omc_flagValue[FLAG_JACOBIAN], JACOBIAN_METHOD[i])){
idaData->jacobianMethod = (int)i;
break;
}
}
if(idaData->daeMode && ( idaData->jacobianMethod == COLOREDSYMJAC || idaData->jacobianMethod == SYMJAC ))
{
warningStreamPrint(LOG_STDOUT, 1, "Symbolic Jacobians are currently not supported by DAE Mode. Switch back to numerical Jacobian!");
idaData->jacobianMethod = COLOREDNUMJAC;
}
if(idaData->jacobianMethod == JAC_UNKNOWN)
{
if (ACTIVE_WARNING_STREAM(LOG_SOLVER))
{
warningStreamPrint(LOG_SOLVER, 1, "unrecognized jacobian calculation method %s, current options are:", (const char*)omc_flagValue[FLAG_JACOBIAN]);
for(i=1; i < JAC_MAX; ++i)
{
warningStreamPrint(LOG_SOLVER, 0, "%-15s [%s]", JACOBIAN_METHOD[i], JACOBIAN_METHOD_DESC[i]);
}
messageClose(LOG_SOLVER);
}
throwStreamPrint(threadData,"unrecognized jacobian calculation method %s", (const char*)omc_flagValue[FLAG_JACOBIAN]);
}
/* default case colored numerical jacobian */
}
else
{
idaData->jacobianMethod = COLOREDNUMJAC;
}
/* selects the calculation method of the jacobian */
/* in daeMode sparse pattern is already initialized in DAEMODE_DATA */
if(!idaData->daeMode && (idaData->jacobianMethod == COLOREDNUMJAC ||
idaData->jacobianMethod == COLOREDSYMJAC ||
idaData->jacobianMethod == SYMJAC))
{
if (data->callback->initialAnalyticJacobianA(data, threadData))
{
infoStreamPrint(LOG_STDOUT, 0, "Jacobian or SparsePattern is not generated or failed to initialize! Switch back to normal.");
idaData->jacobianMethod = INTERNALNUMJAC;
if (idaData->linearSolverMethod == IDA_LS_KLU)
{
idaData->linearSolverMethod = IDA_LS_DENSE;
flag = IDADense(idaData->ida_mem, idaData->N);
warningStreamPrint(LOG_STDOUT, 0, "IDA linear solver method also switched back to %s", IDA_LS_METHOD_DESC[idaData->linearSolverMethod]);
}
}
}
/* set up the appropriate function pointer */
if (idaData->linearSolverMethod == IDA_LS_KLU)
{
if (idaData->daeMode)
{
idaData->NNZ = data->simulationInfo->daeModeData->sparsePattern->numberOfNoneZeros;
flag = IDAKLU(idaData->ida_mem, idaData->N, idaData->NNZ);
idaData->jacobianMethod = COLOREDNUMJAC; /* In DAE mode only numerical matrix is supported */
}
else
{
idaData->NNZ = data->simulationInfo->analyticJacobians[data->callback->INDEX_JAC_A].sparsePattern.numberOfNoneZeros;
flag = IDAKLU(idaData->ida_mem, idaData->N, idaData->NNZ);
/* to add a cj identety matrix */
idaData->tmpJac = NewSparseMat(idaData->N, idaData->N, idaData->N);
}
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Setup of linear solver KLU failed!");
}
switch (idaData->jacobianMethod){
case SYMJAC:
case NUMJAC:
case COLOREDSYMJAC:
case COLOREDNUMJAC:
flag = IDASlsSetSparseJacFn(idaData->ida_mem, callSparseJacobian);
break;
default:
throwStreamPrint(threadData,"For the klu solver jacobian calculation method has to be %s or %s", JACOBIAN_METHOD[COLOREDSYMJAC], JACOBIAN_METHOD[COLOREDNUMJAC]);
break;
}
}
else
{
switch (idaData->jacobianMethod){
case SYMJAC:
case COLOREDSYMJAC:
case COLOREDNUMJAC:
case NUMJAC:
/* set jacobian function */
flag = IDADlsSetDenseJacFn(idaData->ida_mem, callDenseJacobian);
break;
case INTERNALNUMJAC:
break;
default:
throwStreamPrint(threadData,"unrecognized jacobian calculation method %s", (const char*)omc_flagValue[FLAG_JACOBIAN]);
break;
}
}
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Setting jacobian function fails while initialize IDA solver!");
} else {
infoStreamPrint(LOG_SOLVER, 0, "jacobian is calculated by %s", JACOBIAN_METHOD_DESC[idaData->jacobianMethod]);
}
/* set max error test fails */
if (omc_flag[FLAG_IDA_MAXERRORTESTFAIL])
{
int maxErrorTestFails = atoi(omc_flagValue[FLAG_IDA_MAXERRORTESTFAIL]);
flag = IDASetMaxErrTestFails(idaData->ida_mem, maxErrorTestFails);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASetMaxErrTestFails failed!");
}
}
/* set maximum number of nonlinear solver iterations at one step */
if (omc_flag[FLAG_IDA_MAXNONLINITERS])
{
int maxNonlinIters = atoi(omc_flagValue[FLAG_IDA_MAXNONLINITERS]);
flag = IDASetMaxNonlinIters(idaData->ida_mem, maxNonlinIters);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASetMaxNonlinIters failed!");
}
}
/* maximum number of nonlinear solver convergence failures at one step */
if (omc_flag[FLAG_IDA_MAXCONVFAILS])
{
int maxConvFails = atoi(omc_flagValue[FLAG_IDA_MAXCONVFAILS]);
flag = IDASetMaxConvFails(idaData->ida_mem, maxConvFails);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASetMaxConvFails failed!");
}
}
/* safety factor in the nonlinear convergence test */
if (omc_flag[FLAG_IDA_NONLINCONVCOEF])
{
double nonlinConvCoef = atof(omc_flagValue[FLAG_IDA_NONLINCONVCOEF]);
flag = IDASetNonlinConvCoef(idaData->ida_mem, nonlinConvCoef);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASetNonlinConvCoef failed!");
}
}
/* configure algebraic variables as such */
if (idaData->daeMode)
{
if (omc_flag[FLAG_NO_SUPPRESS_ALG])
{
flag = IDASetSuppressAlg(idaData->ida_mem, TRUE);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Suppress algebraic variables in the local error test failed");
}
}
for(i=0; i<idaData->N; ++i)
{
tmp[i] = (i<data->modelData->nStates)? 1.0: 0.0;
}
flag = IDASetId(idaData->ida_mem, N_VMake_Serial(idaData->N,tmp));
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Mark algebraic variables as such failed!");
}
}
/* define initial step size */
if (omc_flag[FLAG_INITIAL_STEP_SIZE])
{
double initialStepSize = atof(omc_flagValue[FLAG_INITIAL_STEP_SIZE]);
assertStreamPrint(threadData, initialStepSize >= DASSL_STEP_EPS, "Selected initial step size %e is too small.", initialStepSize);
flag = IDASetInitStep(idaData->ida_mem, initialStepSize);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Set initial step size failed!");
}
infoStreamPrint(LOG_SOLVER, 0, "initial step size: %g", initialStepSize);
}
else
{
infoStreamPrint(LOG_SOLVER, 0, "initial step size is set automatically.");
}
/* configure the sensitivities part */
idaData->idaSmode = omc_flag[FLAG_IDAS] ? 1 : 0;
if (idaData->idaSmode)
{
idaData->Np = data->modelData->nSensitivityParamVars;
idaData->yS = N_VCloneVectorArray_Serial(idaData->Np, idaData->y);
idaData->ySp = N_VCloneVectorArray_Serial(idaData->Np, idaData->yp);
for(i=0; i<idaData->Np; ++i)
{
int j;
for(j=0; j<idaData->N; ++j)
{
NV_Ith_S(idaData->yS[i],j) = 0;
NV_Ith_S(idaData->ySp[i],j) = 0;
}
}
flag = IDASensInit(idaData->ida_mem, idaData->Np, IDA_SIMULTANEOUS, NULL, idaData->yS, idaData->ySp);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASensInit failed!");
}
flag = IDASetSensParams(idaData->ida_mem, data->simulationInfo->realParameter, NULL, data->simulationInfo->sensitivityParList);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASetSensParams failed!");
}
flag = IDASetSensDQMethod(idaData->ida_mem, IDA_FORWARD, 0);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASetSensDQMethod failed!");
}
flag = IDASensEEtolerances(idaData->ida_mem);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASensEEtolerances failed!");
}
/*
flag = IDASetSensErrCon(idaData->ida_mem, TRUE);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASensEEtolerances failed!");
}
*/
/* allocate result workspace */
idaData->ySResult = N_VCloneVectorArrayEmpty_Serial(idaData->Np, idaData->y);
for(i = 0; i < idaData->Np; ++i)
{
N_VSetArrayPointer_Serial((data->simulationInfo->sensitivityMatrix + i*idaData->N), idaData->ySResult[i]);
}
}
if (compiledInDAEMode == 3){
idaDataGlobal = idaData;
initializedSolver = 1;
data->callback->functionDAE = ida_event_update;
}
free(tmp);
TRACE_POP
return 0;
}
/* deinitial ida data */
int
ida_solver_deinitial(IDA_SOLVER *idaData){
TRACE_PUSH
free(idaData->simData);
free(idaData->ysave);
free(idaData->ypsave);
free(idaData->delta_hh);
if (!idaData->daeMode && idaData->linearSolverMethod == IDA_LS_KLU){
DestroySparseMat(idaData->tmpJac);
}
if (idaData->daeMode)
{
free(idaData->states);
free(idaData->statesDer);
}
if (idaData->idaSmode)
{
N_VDestroyVectorArray_Serial(idaData->yS, idaData->Np);
N_VDestroyVectorArray_Serial(idaData->ySp, idaData->Np);
N_VDestroyVectorArray_Serial(idaData->ySResult, idaData->Np);
}
N_VDestroy_Serial(idaData->errwgt);
N_VDestroy_Serial(idaData->newdelta);
IDAFree(&idaData->ida_mem);
TRACE_POP
return 0;
}
int
ida_event_update(DATA* data, threadData_t *threadData)
{
IDA_SOLVER *idaData = idaDataGlobal;
int flag;
long nonLinIters;
double init_h;
if (compiledInDAEMode == 3){
if (initializedSolver){
data->simulationInfo->needToIterate = 0;
/* get new values from data -> TODO: update all discrete */
data->simulationInfo->discreteCall = 1;
memcpy(idaData->states, data->localData[0]->realVars, sizeof(double)*data->modelData->nStates);
data->simulationInfo->daeModeData->getAlgebraicDAEVars(data, threadData, idaData->states + data->modelData->nStates);
memcpy(idaData->statesDer, data->localData[0]->realVars + data->modelData->nStates, sizeof(double)*data->modelData->nStates);
/* update inner algebraic get new values from data */
idaData->residualFunction(data->localData[0]->timeValue, idaData->y, idaData->yp, idaData->newdelta, idaData);
data->simulationInfo->daeModeData->getAlgebraicDAEVars(data, threadData, idaData->states + data->modelData->nStates);
data->simulationInfo->discreteCall = 0;
infoStreamPrint(LOG_SOLVER, 0, "##IDA## do event update at %.15g", data->localData[0]->timeValue);
flag = IDAReInit(idaData->ida_mem,
data->localData[0]->timeValue,
idaData->y,
idaData->yp);
/* get initial step to provide a direction of the solution */
IDAGetActualInitStep(idaData->ida_mem, &init_h);
/* provide a feasible step-size if it's too small */
if (init_h < DBL_EPSILON){
init_h = DBL_EPSILON;
IDASetInitStep(idaData->ida_mem, init_h);
infoStreamPrint(LOG_SOLVER, 0, "##IDA## corrected step-size at %.15g", init_h);
}
/* increase limits of the non-linear solver */
IDASetMaxNumStepsIC(idaData->ida_mem, 2*idaData->N*10);
IDASetMaxNumJacsIC(idaData->ida_mem, 2*idaData->N*10);
IDASetMaxNumItersIC(idaData->ida_mem, 2*idaData->N*10);
/* Calc Consistent y_algebraic and y_prime with current y */
flag = IDACalcIC(idaData->ida_mem, IDA_YA_YDP_INIT, data->localData[0]->timeValue+init_h);
/* debug */
IDAGetNumNonlinSolvIters(idaData->ida_mem, &nonLinIters);
infoStreamPrint(LOG_SOLVER, 0, "##IDA## IDACalcIC run status %d.\nIterations : %ld\n", flag, nonLinIters);
/* try again without line search if first try fails */
if (checkIDAflag(flag)){
infoStreamPrint(LOG_SOLVER, 0, "##IDA## first event iteration failed. Start next try without line search!");
IDASetLineSearchOffIC(idaData->ida_mem, TRUE);
flag = IDACalcIC(idaData->ida_mem, IDA_YA_YDP_INIT, data->localData[0]->timeValue+data->simulationInfo->tolerance);
IDAGetNumNonlinSolvIters(idaData->ida_mem, &nonLinIters);
infoStreamPrint(LOG_SOLVER, 0, "##IDA## IDACalcIC run status %d.\nIterations : %ld\n", flag, nonLinIters);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## discrete update failed flag %d!", flag);
}
}
/* obtain consistent values of y_algebraic and y_prime */
IDAGetConsistentIC(idaData->ida_mem, idaData->y, idaData->yp);
/* update inner algebraic variables */
data->simulationInfo->discreteCall = 1;
idaData->residualFunction(data->localData[0]->timeValue, idaData->y, idaData->yp, idaData->newdelta, idaData);
data->simulationInfo->discreteCall = 0;
memcpy(data->localData[0]->realVars, idaData->states, sizeof(double)*data->modelData->nStates);
// and also algebraic vars
data->simulationInfo->daeModeData->setAlgebraicDAEVars(data, threadData, idaData->states + data->modelData->nStates);
memcpy(data->localData[0]->realVars + data->modelData->nStates, idaData->statesDer, sizeof(double)*data->modelData->nStates);
/* reset initial step size again to default */
IDASetInitStep(idaData->ida_mem, 0.0);
}
}
else{
data->callback->functionDAE(data, threadData);
}
}
/* main ida function to make a step */
int
ida_solver_step(DATA* data, threadData_t *threadData, SOLVER_INFO* solverInfo)
{
TRACE_PUSH
double tout = 0;
int i = 0, flag;
int retVal = 0, finished = FALSE;
int saveJumpState;
long int tmp;
static unsigned int stepsOutputCounter = 1;
int stepsMode;
int restartAfterLSFail = 0;
IDA_SOLVER *idaData = (IDA_SOLVER*) solverInfo->solverData;
SIMULATION_DATA *sData = data->localData[0];
SIMULATION_DATA *sDataOld = data->localData[1];
MODEL_DATA *mData = (MODEL_DATA*) data->modelData;
/* alloc all work arrays */
if (!idaData->daeMode)
{
N_VSetArrayPointer(data->localData[0]->realVars, idaData->y);
N_VSetArrayPointer(data->localData[1]->realVars + data->modelData->nStates, idaData->yp);
}
if (solverInfo->didEventStep)
{
idaData->setInitialSolution = 0;
}
/* reinit solver */
if (!idaData->setInitialSolution)
{
debugStreamPrint(LOG_SOLVER, 0, "Re-initialized IDA Solver");
/* initialize states and der(states) */
if (idaData->daeMode)
{
memcpy(idaData->states, data->localData[0]->realVars, sizeof(double)*data->modelData->nStates);
/* and also algebraic vars */
data->simulationInfo->daeModeData->getAlgebraicDAEVars(data, threadData, idaData->states + data->modelData->nStates);
memcpy(idaData->statesDer, data->localData[0]->realVars + data->modelData->nStates, sizeof(double)*data->modelData->nStates);
}
/* calculate matrix for residual scaling */
if (omc_flag[FLAG_IDA_SCALING])
{
getScalingFactors(data, idaData, NULL);
/* scale idaData->y and idaData->yp */
infoStreamPrint(LOG_SOLVER_V, 1, "Scale y and yp");
idaScaleData(idaData);
messageClose(LOG_SOLVER_V);
}
flag = IDAReInit(idaData->ida_mem,
solverInfo->currentTime,
idaData->y,
idaData->yp);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Something goes wrong while reinit IDA solver after event!");
}
/* calculate matrix for residual scaling */
if (omc_flag[FLAG_IDA_SCALING])
{
/* scale idaData->y and idaData->yp */
idaReScaleData(idaData);
}
if (idaData->idaSmode)
{
for(i=0; i<idaData->Np; ++i)
{
int j;
for(j=0; j<idaData->N; ++j)
{
NV_Ith_S(idaData->yS[i],j) = 0;
NV_Ith_S(idaData->ySp[i],j) = 0;
}
}
flag = IDASensReInit(idaData->ida_mem, IDA_SIMULTANEOUS, idaData->yS, idaData->ySp);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## set IDASensInit failed!");
}
}
idaData->setInitialSolution = 1;
}
saveJumpState = threadData->currentErrorStage;
threadData->currentErrorStage = ERROR_INTEGRATOR;
/* try */
#if !defined(OMC_EMCC)
MMC_TRY_INTERNAL(simulationJumpBuffer)
#endif
/* Check that tout is not less than timeValue otherwise the solver
* will come in trouble.
* If that is the case we skip the current step. */
if (solverInfo->currentStepSize < DASSL_STEP_EPS)
{
infoStreamPrint(LOG_SOLVER, 0, "Desired step to small try next one");
infoStreamPrint(LOG_SOLVER, 0, "Interpolate linear");
/* linear extrapolation */
for(i = 0; i < idaData->N; i++)
{
NV_Ith_S(idaData->y, i) = NV_Ith_S(idaData->y, i) + NV_Ith_S(idaData->yp, i) * solverInfo->currentStepSize;
}
sData->timeValue = solverInfo->currentTime + solverInfo->currentStepSize;
data->callback->functionODE(data, threadData);
solverInfo->currentTime = sData->timeValue;
TRACE_POP
return 0;
}
/* Calculate steps until TOUT is reached */
if (idaData->internalSteps)
{
/* If internalSteps are selected, let IDA run to stopTime or next sample event */
if (data->simulationInfo->nextSampleEvent < data->simulationInfo->stopTime)
{
tout = data->simulationInfo->nextSampleEvent;
}
else
{
tout = data->simulationInfo->stopTime;
}
stepsMode = IDA_ONE_STEP;
flag = IDASetStopTime(idaData->ida_mem, tout);
if (checkIDAflag(flag)){
throwStreamPrint(threadData, "##IDA## Something goes wrong while set stopTime!");
}
}
else
{
tout = solverInfo->currentTime + solverInfo->currentStepSize;
stepsMode = IDA_NORMAL;
}
do
{
infoStreamPrint(LOG_SOLVER, 1, "##IDA## new step from %.15g to %.15g", solverInfo->currentTime, tout);
/* read input vars */
externalInputUpdate(data);
data->callback->input_function(data, threadData);
if (omc_flag[FLAG_IDA_SCALING])
{
/* scale idaData->y and idaData->yp */
idaScaleData(idaData);
}
flag = IDASolve(idaData->ida_mem, tout, &solverInfo->currentTime, idaData->y, idaData->yp, stepsMode);
if (omc_flag[FLAG_IDA_SCALING])
{
/* rescale idaData->y and idaData->yp */
idaReScaleData(idaData);
}
/* set time to current time */
sData->timeValue = solverInfo->currentTime;
/* error handling */
if ( !checkIDAflag(flag) && solverInfo->currentTime >= tout)
{
infoStreamPrint(LOG_SOLVER, 0, "##IDA## step done to time = %.15g", solverInfo->currentTime);
finished = TRUE;
}
else
{
if (!checkIDAflag(flag))
{
infoStreamPrint(LOG_SOLVER, 0, "##IDA## continue integration time = %.15g", solverInfo->currentTime);
}
else if (flag == IDA_ROOT_RETURN)
{
infoStreamPrint(LOG_SOLVER, 0, "##IDA## root found at time = %.15g", solverInfo->currentTime);
finished = TRUE;
}
else if (flag == IDA_TOO_MUCH_WORK)
{
warningStreamPrint(LOG_SOLVER, 0, "##IDA## has done too much work with small steps at time = %.15g", solverInfo->currentTime);
}
else if (flag == IDA_LSETUP_FAIL && !restartAfterLSFail )
{
flag = IDAReInit(idaData->ida_mem,
solverInfo->currentTime,
idaData->y,
idaData->yp);
restartAfterLSFail = 1;
warningStreamPrint(LOG_SOLVER, 0, "##IDA## linear solver failed try once again = %.15g", solverInfo->currentTime);
}
else
{
infoStreamPrint(LOG_STDOUT, 0, "##IDA## %d error occurred at time = %.15g", flag, solverInfo->currentTime);
finished = TRUE;
retVal = flag;
}
}
/* closing new step message */
messageClose(LOG_SOLVER);
/* emit step, if step mode is selected */
if (idaData->internalSteps)
{
infoStreamPrint(LOG_SOLVER, 0, "##IDA## noEquadistant stepsOutputCounter %d by freq %d at time = %.15g", stepsOutputCounter, idaData->stepsFreq, solverInfo->currentTime);
if (omc_flag[FLAG_NOEQUIDISTANT_OUT_FREQ]){
/* output every n-th time step */
if (stepsOutputCounter >= idaData->stepsFreq){
stepsOutputCounter = 1; /* next line set it to one */
infoStreamPrint(LOG_SOLVER, 0, "##IDA## noEquadistant output %d by freq at time = %.15g", stepsOutputCounter, solverInfo->currentTime);