-
Notifications
You must be signed in to change notification settings - Fork 298
/
kinsolSolver.c
1280 lines (1118 loc) · 43.1 KB
/
kinsolSolver.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-CurrentYear, 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 kinsolSolver.c
*/
#include "kinsolSolver.h"
#include "nonlinearSystem.h"
#include "omc_config.h"
#include "omc_math.h"
#include "simulation/options.h"
#include "simulation/simulation_info_json.h"
#include "sundials_util.h"
#include "util/omc_error.h"
#ifdef WITH_SUNDIALS
#include "events.h"
#include "model_help.h"
#include "openmodelica.h"
#include "openmodelica_func.h"
#include "util/read_matlab4.h"
#include "util/varinfo.h"
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
/* Function prototypes */
static int nlsKinsolResiduals(N_Vector x, N_Vector f, void* userData);
static int nlsSparseJac(N_Vector vecX, N_Vector vecFX, SUNMatrix Jac,
void* userData, N_Vector tmp1, N_Vector tmp2);
int nlsSparseSymJac(N_Vector vecX, N_Vector vecFX, SUNMatrix Jac,
void* userData, N_Vector tmp1, N_Vector tmp2);
static int nlsDenseJac(long int N, N_Vector vecX, N_Vector vecFX,
SUNMatrix Jac, NLS_USERDATA *kinsolUserData,
N_Vector tmp1, N_Vector tmp2);
static void nlsKinsolJacSumSparse(SUNMatrix A);
static void nlsKinsolJacSumDense(SUNMatrix A);
/**
* @brief Set KINSOL configuration.
*
* @param kinsolData Kinsol data with configuration settings.
*/
static void nlsKinsolConfigSetup(NLS_KINSOL_DATA *kinsolData) {
/* Variables */
int flag;
/* configuration */
flag = KINSetFuncNormTol(kinsolData->kinsolMemory,
kinsolData->fnormtol); /* Set function-norm stopping tolerance */
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetFuncNormTol");
kinsolData->resetTol = FALSE;
flag = KINSetScaledStepTol(kinsolData->kinsolMemory,
kinsolData->scsteptol); /* Set scaled-step stopping tolerance */
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetScaledStepTol");
flag = KINSetNumMaxIters(kinsolData->kinsolMemory,
100 * kinsolData->size); /* Set max. number of nonlinear iterations */
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetNumMaxIters");
kinsolData->kinsolStrategy = KIN_LINESEARCH; /* Newton with globalization strategy to solve nonlinear systems */
flag = KINSetNoInitSetup(kinsolData->kinsolMemory, SUNFALSE); /* TODO: This is the default value. Is there a point in calling this function? */
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetNoInitSetup");
kinsolData->retries = 0;
kinsolData->countResCalls = 0;
}
/**
* @brief (Re-) Initialize KINSOL data.
*
* Initialize KINSOL data. If the KINSOL memory block was already initialized
* free it first and then reinitialize.
*
* @param kinsolData KINSOL data.
*/
void resetKinsolMemory(NLS_KINSOL_DATA *kinsolData) {
int flag;
int printLevel;
int size = kinsolData->size;
NONLINEAR_SYSTEM_DATA *nlsData = kinsolData->userData->nlsData;
SPARSE_PATTERN* sparsePattern = nlsData->sparsePattern;
/* Free KINSOL memory block */
if (kinsolData->kinsolMemory) {
KINFree((void *)&kinsolData->kinsolMemory);
}
/* Create KINSOL memory block */
kinsolData->kinsolMemory = KINCreate();
if (kinsolData->kinsolMemory == NULL) {
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: In function KINCreate: An error occurred.");
}
/* Set error handler and print level */
if (!nlsData->logActive) {
printLevel = 0;
} else if (ACTIVE_STREAM(LOG_NLS_V)) {
printLevel = 3;
} else if (ACTIVE_STREAM(LOG_NLS)) {
printLevel = 1;
} else {
printLevel = 0;
}
flag = KINSetPrintLevel(kinsolData->kinsolMemory, printLevel);
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetPrintLevel");
flag = KINSetErrHandlerFn(kinsolData->kinsolMemory, kinsolErrorHandlerFunction, kinsolData);
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetErrHandlerFn");
flag = KINSetInfoHandlerFn(kinsolData->kinsolMemory, kinsolInfoHandlerFunction, NULL);
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetInfoHandlerFn");
flag = KINSetUserData(kinsolData->kinsolMemory, (void*)kinsolData->userData);
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetUserData");
/* Initialize KINSOL object */
flag = KINInit(kinsolData->kinsolMemory, nlsKinsolResiduals,
kinsolData->initialGuess);
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINInit");
/* Create matrix object */
if (kinsolData->linearSolverMethod == NLS_LS_DEFAULT ||
kinsolData->linearSolverMethod == NLS_LS_LAPACK) {
kinsolData->J = SUNDenseMatrix(size, size);
} else if (kinsolData->linearSolverMethod == NLS_LS_KLU) {
if (!sparsePattern) {
kinsolData->nnz = size*size;
} else {
kinsolData->nnz = sparsePattern->numberOfNonZeros;
}
kinsolData->J = SUNSparseMatrix(size, size, kinsolData->nnz, CSC_MAT);
} else {
kinsolData->J = NULL;
}
/* Create linear solver object */
if (kinsolData->linearSolverMethod == NLS_LS_DEFAULT ||
kinsolData->linearSolverMethod == NLS_LS_TOTALPIVOT) {
kinsolData->linSol = SUNLinSol_Dense(kinsolData->y, kinsolData->J);
if (kinsolData->linSol == NULL) {
throwStreamPrint(NULL, "KINSOL: In function SUNLinSol_Dense: Input incompatible.");
}
} else if (kinsolData->linearSolverMethod == NLS_LS_LAPACK) {
kinsolData->linSol = SUNLinSol_LapackDense(kinsolData->y, kinsolData->J);
if (kinsolData->linSol == NULL) {
throwStreamPrint(NULL, "KINSOL: In function SUNLinSol_LapackDense: Input incompatible.");
}
} else if (kinsolData->linearSolverMethod == NLS_LS_KLU) {
kinsolData->linSol = SUNLinSol_KLU(kinsolData->y, kinsolData->J);
if (kinsolData->linSol == NULL) {
throwStreamPrint(NULL, "KINSOL: In function SUNLinSol_KLU: Input incompatible.");
}
} else {
throwStreamPrint(NULL, "KINSOL: Unknown linear solver method.");
}
/* Log used solver */
infoStreamPrint(LOG_NLS, 0, "KINSOL: Using linear solver method %s", NLS_LS_METHOD[kinsolData->linearSolverMethod]);
/* Set linear solver */
flag = KINSetLinearSolver(kinsolData->kinsolMemory, kinsolData->linSol,
kinsolData->J);
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KINLS_FLAG, "KINSetLinearSolver");
/* Set Jacobian for non-linear solver */
if (kinsolData->linearSolverMethod == NLS_LS_KLU) {
if (nlsData->analyticalJacobianColumn != NULL && sparsePattern != NULL) {
flag = KINSetJacFn(kinsolData->kinsolMemory, nlsSparseSymJac); /* Use symbolic Jacobian with sparsity pattern*/
} else if (sparsePattern != NULL) {
flag = KINSetJacFn(kinsolData->kinsolMemory, nlsSparseJac); /* Use numeric Jacobian with sparsity pattern */
} else {
throwStreamPrint(NULL, "KINSOL: In function resetKinsolMemory: Sparse linear solver KLU needs sparse Jacobian, but no sparsity pattern is available. Use a dense non-linear solver instead of KINSOL.");
}
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KINLS_FLAG, "KINSetJacFn");
}
/* Configuration */
nlsKinsolConfigSetup(kinsolData);
}
/**
* @brief Allocate memory for kinsol solver data and initialize KINSOL solver.
*
* @param size Size of non-linear problem.
* @param userData Pointer to set NLS user data.
* @param attemptRetry True if KINSOL should retry with different settings after solution failed.
* @param isPatternAvailable True if sparsity pattern of Jacobian is available. Allocate work vectors for KLU in that case.
* @return NLS_KINSOL_DATA* Pointer to allocated KINSOL data.
*/
NLS_KINSOL_DATA* nlsKinsolAllocate(int size, NLS_USERDATA* userData, modelica_boolean attemptRetry, modelica_boolean isPatternAvailable) {
/* Allocate system data */
NLS_KINSOL_DATA *kinsolData = (NLS_KINSOL_DATA *)calloc(1, sizeof(NLS_KINSOL_DATA));
kinsolData->size = size;
kinsolData->linearSolverMethod = userData->nlsData->nlsLinearSolver;
kinsolData->solved = NLS_FAILED;
kinsolData->fnormtol = newtonFTol; /* function tolerance */
kinsolData->scsteptol = newtonXTol; /* step tolerance */
kinsolData->maxstepfactor = maxStepFactor; /* step tolerance */
kinsolData->nominalJac = 0; /* calculate for scaling the scaled matrix */
kinsolData->attemptRetry = attemptRetry;
kinsolData->initialGuess = N_VNew_Serial(size);
kinsolData->xScale = N_VNew_Serial(size);
kinsolData->fScale = N_VNew_Serial(size);
kinsolData->fRes = N_VNew_Serial(size);
kinsolData->fTmp = N_VNew_Serial(size);
kinsolData->y = N_VNew_Serial(size);
if (isPatternAvailable && kinsolData->linearSolverMethod == NLS_LS_KLU) {
kinsolData->tmp1 = N_VNew_Serial(size);
kinsolData->tmp2 = N_VNew_Serial(size);
} else {
kinsolData->tmp1 = NULL;
kinsolData->tmp2 = NULL;
}
kinsolData->kinsolMemory = NULL;
kinsolData->userData = userData;
resetKinsolMemory(kinsolData);
return kinsolData;
}
/**
* @brief Deallocates memory for KINSOL solver.
*
* Free memory that was allocated with `nlsKinsolAllocate`.
*
* @param kinsolData Pointer to KINSOL data.
*/
void nlsKinsolFree(NLS_KINSOL_DATA* kinsolData) {
KINFree((void *)&kinsolData->kinsolMemory);
N_VDestroy_Serial(kinsolData->initialGuess);
N_VDestroy_Serial(kinsolData->xScale);
N_VDestroy_Serial(kinsolData->fScale);
N_VDestroy_Serial(kinsolData->fRes);
N_VDestroy_Serial(kinsolData->fTmp);
/* Free linear solver data */
SUNLinSolFree(kinsolData->linSol);
SUNMatDestroy(kinsolData->J);
N_VDestroy_Serial(kinsolData->y);
if (kinsolData->tmp1 != NULL) {
N_VDestroy_Serial(kinsolData->tmp1);
N_VDestroy_Serial(kinsolData->tmp2);
}
freeNlsUserData(kinsolData->userData);
free(kinsolData);
return;
}
/**
* @brief Residual function for non-linear problem.
*
* @param x The current value of the variable vector.
* @param f Output vector.
* @param userData Pointer to Kinsol user data.
* @return int Return 0 on success, return 1 on recoverable error.
*/
static int nlsKinsolResiduals(N_Vector x, N_Vector f, void* userData) {
double *xdata = NV_DATA_S(x);
double *fdata = NV_DATA_S(f);
NLS_USERDATA* kinsolUserData = (NLS_USERDATA*)userData;
DATA* data = kinsolUserData->data;
threadData_t* threadData = kinsolUserData->threadData;
NONLINEAR_SYSTEM_DATA* nlsData = kinsolUserData->nlsData;
NLS_KINSOL_DATA* kinsolData = (NLS_KINSOL_DATA*)nlsData->solverData;
RESIDUAL_USERDATA resUserData = {.data=data, .threadData=threadData, .solverData=kinsolUserData->solverData};
int iflag = 1 /* recoverable error */;
/* Update statistics */
kinsolData->countResCalls++;
#ifndef OMC_EMCC
MMC_TRY_INTERNAL(simulationJumpBuffer)
#endif
/* call residual function */
nlsData->residualFunc(&resUserData, xdata, fdata, (const int *)&iflag);
iflag = 0 /* success */;
#ifndef OMC_EMCC
MMC_CATCH_INTERNAL(simulationJumpBuffer)
#endif
return iflag;
}
/**
* @brief Calculate dense Jacobian matrix.
*
* @param N Size of vecX and vecFX.
* @param vecX Vector x.
* @param vecFX Residual vector f(x).
* @param Jac Dense Jacobian matrix J(x).
* @param kinsolUserData Pointer to Kinsol user data.
* @param tmp1 Unused, only to match interface of KINLsJacFn
* @param tmp2 Unused, only to match interface of KINLsJacFn
* @return int Return 0 on success, -1 on failure.
*/
static int nlsDenseJac(long int N,
N_Vector vecX,
N_Vector vecFX,
SUNMatrix Jac,
NLS_USERDATA *kinsolUserData,
N_Vector tmp1,
N_Vector tmp2) {
DATA *data = kinsolUserData->data;
threadData_t *threadData = kinsolUserData->threadData;
NONLINEAR_SYSTEM_DATA *nlsData = kinsolUserData->nlsData;
NLS_KINSOL_DATA *kinsolData = (NLS_KINSOL_DATA *)nlsData->solverData;
if (SUNMatGetID(Jac) != SUNMATRIX_DENSE) {
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: nlsDenseJac illegal input Jac. Matrix is not dense!");
return -1;
}
/* prepare variables */
double *x = N_VGetArrayPointer(vecX);
double *fx = N_VGetArrayPointer(vecFX);
double *xScaling = NV_DATA_S(kinsolData->xScale);
double *fRes = NV_DATA_S(kinsolData->fRes);
double xsave, xscale, sign;
double delta_hh;
const double delta_h = sqrt(DBL_EPSILON * 2e1);
long int i, j;
/* performance measurement */
rt_ext_tp_tick(&nlsData->jacobianTimeClock);
/* Use forward difference quotient to approximate Jacobian */
for (i = 0; i < N; i++) {
xsave = x[i];
delta_hh = delta_h * (fabs(xsave) + 1.0);
if ((xsave + delta_hh >= nlsData->max[i])) {
delta_hh *= -1.0;
}
x[i] += delta_hh;
/* Evaluate Jacobian function */
nlsKinsolResiduals(vecX, kinsolData->fRes, kinsolUserData);
/* Calculate scaled difference quotient */
delta_hh = 1.0 / delta_hh;
for (j = 0; j < N; j++) {
if (kinsolData->nominalJac) {
SM_ELEMENT_D(Jac, j, i) = (fRes[j] - fx[j]) * delta_hh / xScaling[i];
} else {
SM_ELEMENT_D(Jac, j, i) =
(fRes[j] - fx[j]) * delta_hh; /* TODO: Or now Jac(i,j) ??? */
}
}
x[i] = xsave;
}
/* debug */
if (ACTIVE_STREAM(LOG_NLS_JAC)) {
infoStreamPrint(LOG_NLS_JAC, 1, "KINSOL: Dense matrix.");
SUNDenseMatrix_Print(Jac, stdout); /* TODO: Print in LOG_NLS_JAC */
nlsKinsolJacSumDense(Jac);
messageClose(LOG_NLS_JAC);
}
/* performance measurement and statistics */
nlsData->jacobianTime += rt_ext_tp_tock(&(nlsData->jacobianTimeClock));
nlsData->numberOfJEval++;
return 0;
}
/**
* @brief Finish sparse matrix by fixing colprts.
*
* Last value of indexptrs should always be nnz.
* Search for empty rows which would mean the matrix is singular.
*
* @param A CSC matrix
*/
static void finishSparseColPtr(SUNMatrix A, int nnz) {
int i;
/* TODO: Remove this check for performance reasons? */
if (SM_SPARSETYPE_S(A) != CSC_MAT) {
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: In function finishSparseColPtr: Wrong sparse format of SUNMatrix A.");
}
/* Set last value of indexptrs to nnz */
SM_INDEXPTRS_S(A)[SM_COLUMNS_S(A)] = nnz;
/* Check for empty rows */
for (i = 1; i < SM_COLUMNS_S(A) + 1; ++i) {
if (SM_INDEXPTRS_S(A)[i] == SM_INDEXPTRS_S(A)[i - 1]) {
warningStreamPrint(LOG_STDOUT, 0,
"KINSOL: Jacobian row %d singular. See LOG_NLS for "
"more information.",
i);
SM_INDEXPTRS_S(A)[i] = SM_INDEXPTRS_S(A)[i - 1];
}
}
}
/**
* @brief Colored numeric Jacobian evaluation.
*
* Finite differences while using coloring of Jacobian.
* Jacobian matrix format has to be compressed sparse columns (CSC).
*
* @param vecX Input vector x.
* @param vecFX Vector for residual evaluation: f(x)
* @param Jac Jacobian to calculate: J(x)
* @param userData Pointer to user data, tpyecasted to `NLS_USERDATA`.
* @param tmp1 Work vector.
* @param tmp2 Work vector.
* @return int Return 0 on success.
*/
static int nlsSparseJac(N_Vector vecX, N_Vector vecFX, SUNMatrix Jac,
void *userData, N_Vector tmp1, N_Vector tmp2) {
/* Variables */
NLS_USERDATA *kinsolUserData;
DATA *data;
threadData_t *threadData;
NONLINEAR_SYSTEM_DATA *nlsData;
NLS_KINSOL_DATA *kinsolData;
SPARSE_PATTERN *sparsePattern;
if (SUNMatGetID(Jac) != SUNMATRIX_SPARSE || SM_SPARSETYPE_S(Jac) == CSR_MAT) {
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: nlsSparseJac illegal input Jac. Matrix is not sparse!");
return -1;
}
double *x;
double *fx;
double *xsave;
double *delta_hh;
double *xScaling;
double *fRes;
const double delta_h = sqrt(DBL_EPSILON * 2e1);
long int i, j, ii;
int nth;
/* Access userData and nonlinear system data */
kinsolUserData = (NLS_USERDATA *)userData;
data = kinsolUserData->data;
threadData = kinsolUserData->threadData;
nlsData = kinsolUserData->nlsData;
kinsolData = (NLS_KINSOL_DATA *)nlsData->solverData;
sparsePattern = nlsData->sparsePattern;
/* Access N_Vector variables */
x = N_VGetArrayPointer(vecX);
fx = N_VGetArrayPointer(vecFX);
xsave = N_VGetArrayPointer(tmp1);
delta_hh = N_VGetArrayPointer(tmp2);
xScaling = NV_DATA_S(kinsolData->xScale);
fRes = NV_DATA_S(kinsolData->fRes);
nth = 0;
/* performance measurement */
rt_ext_tp_tick(&nlsData->jacobianTimeClock);
/* reset matrix */
SUNMatZero(Jac);
/* Approximate Jacobian */
for (i = 0; i < sparsePattern->maxColors; i++) {
for (ii = 0; ii < kinsolData->size; ii++) {
if (sparsePattern->colorCols[ii] - 1 == i) {
xsave[ii] = x[ii];
delta_hh[ii] = delta_h * (fabs(xsave[ii]) + 1.0);
if ((xsave[ii] + delta_hh[ii] >= nlsData->max[ii])) {
delta_hh[ii] *= -1;
}
x[ii] += delta_hh[ii];
/* Calculate scaled difference quotient */
delta_hh[ii] = 1. / delta_hh[ii];
}
}
/* Evaluate residual function */
nlsKinsolResiduals(vecX, kinsolData->fRes, userData);
/* Save column in Jac and unset seed variables */
for (ii = 0; ii < kinsolData->size; ii++) {
if (sparsePattern->colorCols[ii] - 1 == i) {
nth = sparsePattern->leadindex[ii];
while (nth < sparsePattern->leadindex[ii + 1]) {
j = sparsePattern->index[nth];
if (kinsolData->nominalJac) {
setJacElementSundialsSparse(j, ii, nth, (fRes[j] - fx[j]) * delta_hh[ii] / xScaling[ii], Jac, SM_CONTENT_S(Jac)->M);
} else {
setJacElementSundialsSparse(j, ii, nth, (fRes[j] - fx[j]) * delta_hh[ii], Jac, SM_CONTENT_S(Jac)->M);
}
nth++;
}
x[ii] = xsave[ii];
}
}
}
/* Finish sparse matrix */
finishSparseColPtr(Jac, sparsePattern->numberOfNonZeros);
/* Debug print */
if (ACTIVE_STREAM(LOG_NLS_JAC)) {
infoStreamPrint(LOG_NLS_JAC, 1, "KINSOL: Sparse Matrix.");
SUNSparseMatrix_Print(Jac, stdout);
nlsKinsolJacSumSparse(Jac);
messageClose(LOG_NLS_JAC);
}
if (ACTIVE_STREAM(LOG_DEBUG)) {
sundialsPrintSparseMatrix(Jac, "A", LOG_JAC);
}
/* performance measurement and statistics */
nlsData->jacobianTime += rt_ext_tp_tock(&(nlsData->jacobianTimeClock));
nlsData->numberOfJEval++;
return 0;
}
/**
* @brief Computes symbolic Jacobian matrix Jac(vecX)
*
* @param vecX
* @param vecFX just for interface compatibility, will not be used here
* @param Jac
* @param userData
* @param tmp1 Unused, only to match interface of KINLsJacFn
* @param tmp2 Unused, only to match interface of KINLsJacFn
* @return int
*/
int nlsSparseSymJac(N_Vector vecX, N_Vector vecFX, SUNMatrix Jac,
void *userData, N_Vector tmp1, N_Vector tmp2) {
/* Variables */
NLS_USERDATA *kinsolUserData;
DATA *data;
threadData_t *threadData;
NONLINEAR_SYSTEM_DATA *nlsData;
NLS_KINSOL_DATA *kinsolData;
SPARSE_PATTERN *sparsePattern;
ANALYTIC_JACOBIAN *analyticJacobian;
double *x;
double *fx;
double *xScaling;
long int i, j, ii;
int nth;
if (SUNMatGetID(Jac) != SUNMATRIX_SPARSE || SM_SPARSETYPE_S(Jac) == CSR_MAT) {
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: nlsSparseJac illegal input Jac. Matrix is not sparse!");
return -1;
}
/* Access userData and nonlinear system data */
kinsolUserData = (NLS_USERDATA *)userData;
data = kinsolUserData->data;
threadData = kinsolUserData->threadData;
nlsData = kinsolUserData->nlsData;
analyticJacobian = kinsolUserData->analyticJacobian;
kinsolData = (NLS_KINSOL_DATA *)nlsData->solverData;
sparsePattern = nlsData->sparsePattern;
/* Access N_Vector variables */
x = N_VGetArrayPointer(vecX);
xScaling = NV_DATA_S(kinsolData->xScale);
nth = 0;
/* performance measurement */
rt_ext_tp_tick(&nlsData->jacobianTimeClock);
/* reset matrix */
SUNMatZero(Jac);
/* Evaluate constant equations of Jacobian */
if (analyticJacobian->constantEqns != NULL) {
analyticJacobian->constantEqns(data, threadData, analyticJacobian, NULL);
}
/* Evaluate Jacobian */
for (i = 0; i < sparsePattern->maxColors; i++) {
/* Set seed variables */
for (ii = 0; ii < kinsolData->size; ii++) {
if (sparsePattern->colorCols[ii] - 1 == i) {
analyticJacobian->seedVars[ii] = 1.0;
}
}
/* Evaluate Jacobian column */
((nlsData->analyticalJacobianColumn))(data, threadData, analyticJacobian, NULL);
/* Save column in Jac and unset seed variables */
for (ii = 0; ii < kinsolData->size; ii++) {
if (sparsePattern->colorCols[ii] - 1 == i) {
nth = sparsePattern->leadindex[ii];
while (nth < sparsePattern->leadindex[ii + 1]) {
j = sparsePattern->index[nth];
if (kinsolData->nominalJac) {
setJacElementSundialsSparse(j, ii, nth, analyticJacobian->resultVars[j] / xScaling[ii], Jac, SM_CONTENT_S(Jac)->M);
} else {
setJacElementSundialsSparse(j, ii, nth, analyticJacobian->resultVars[j], Jac, SM_CONTENT_S(Jac)->M);
}
nth++;
}
analyticJacobian->seedVars[ii] = 0;
}
}
}
/* Finish sparse matrix and do a cheap check for singularity */
finishSparseColPtr(Jac, sparsePattern->numberOfNonZeros);
/* Debug print */
if (ACTIVE_STREAM(LOG_NLS_JAC)) {
infoStreamPrint(LOG_NLS_JAC, 1, "KINSOL: Sparse Matrix.");
SUNSparseMatrix_Print(Jac, stdout); /* TODO: Print in LOG_NLS_JAC */
nlsKinsolJacSumSparse(Jac);
messageClose(LOG_NLS_JAC);
}
/* performance measurement and statistics */
nlsData->jacobianTime += rt_ext_tp_tock(&(nlsData->jacobianTimeClock));
nlsData->numberOfJEval++;
return 0;
}
/**
* @brief Check for zero columns of matrix and print absolute sums.
*
* Compute absolute sum for each column and print the result.
* Report a warning if it is zero, since the matrix is singular in that case.
*
* @param A Dense matrix stored columnwise
*/
static void nlsKinsolJacSumDense(SUNMatrix A) {
/* Variables */
int i, j;
double sum;
for (i = 0; i < SM_ROWS_D(A); ++i) {
sum = 0.0;
for (j = 0; j < SM_COLUMNS_D(A); ++j) {
sum += fabs(SM_ELEMENT_D(A, j, i));
}
if (sum == 0.0) { /* TODO: Don't check for equality(!), maybe use DBL_EPSILON */
warningStreamPrint(LOG_NLS_V, 0,
"KINSOL: Column %d of Jacobian is zero. Jacobian is singular.",
i);
} else {
infoStreamPrint(LOG_NLS_JAC, 0, "Column %d of Jacobian absolute sum = %g",
i, sum);
}
}
}
/**
* @brief Check for zero columns of matrix and print absolute sums.
*
* Compute absolute sum for each column and print the result.
* Report a warning if it is zero, since the matrix is singular in that case.
*
* @param A CSC matrix
*/
static void nlsKinsolJacSumSparse(SUNMatrix A) {
/* Variables */
int i, j;
double sum;
/* Check format of A */
if (SM_SPARSETYPE_S(A) != CSC_MAT) {
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: In function nlsKinsolJacSumSparse: Wrong sparse format "
"of SUNMatrix A.");
}
/* Check sums of each column of A */
for (i = 0; i < SM_COLUMNS_S(A); ++i) {
sum = 0.0;
for (j = SM_INDEXPTRS_S(A)[i]; j < SM_INDEXPTRS_S(A)[i + 1]; ++j) {
sum += fabs(SM_DATA_S(A)[j]);
}
if (sum == 0.0) { /* TODO: Don't check for equality(!), maybe use DBL_EPSILON */
warningStreamPrint(LOG_NLS_V, 0,
"KINSOL: Column %d of Jacobian is zero. Jacobian is singular.",
i);
} else {
infoStreamPrint(LOG_NLS_JAC, 0, "Column %d of Jacobian absolute sum = %g",
i, sum);
}
}
}
/**
* @brief Set maximum scaled length of Newton step.
*
* Will be set to the weighted Euclidean l_2 norm of xScale with maxstepfactor
* as weights. maxStep = sqrt(sum_{1=0}^{n-1} (xScale[i]*maxstepfactor)^2)
*
* @param kinsolData
* @param maxstepfactor
*/
static void nlsKinsolSetMaxNewtonStep(NLS_KINSOL_DATA *kinsolData,
double maxstepfactor) {
/* Variables */
int flag;
N_VConst(maxstepfactor, kinsolData->fTmp);
kinsolData->mxnstepin = N_VWL2Norm(kinsolData->xScale, kinsolData->fTmp);
/* Set maximum step size */
flag = KINSetMaxNewtonStep(kinsolData->kinsolMemory, kinsolData->mxnstepin);
checkReturnFlag_SUNDIALS(flag, SUNDIALS_KIN_FLAG, "KINSetMaxNewtonStep");
}
/**
* @brief Set initial guess for KINSOL
*
* Depending on mode extrapolate start value or use old value for
* initialization.
*
* @param data
* @param kinsolData
* @param nlsData
* @param mode Has to be `INITIAL_EXTRAPOLATION` for extrapolation or
* `INITIAL_OLDVALUES` for using old values.
*/
static void nlsKinsolResetInitial(DATA *data, NLS_KINSOL_DATA *kinsolData,
NONLINEAR_SYSTEM_DATA *nlsData,
initialMode mode) {
double *xStart = NV_DATA_S(kinsolData->initialGuess);
/* Set x vector */
switch (mode) {
case INITIAL_EXTRAPOLATION:
if (data->simulationInfo->discreteCall) {
memcpy(xStart, nlsData->nlsx, nlsData->size * (sizeof(double)));
} else {
memcpy(xStart, nlsData->nlsxExtrapolation,
nlsData->size * (sizeof(double)));
}
break;
case INITIAL_OLDVALUES:
memcpy(xStart, nlsData->nlsxOld, nlsData->size * (sizeof(double)));
break;
default:
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: Function nlsKinsolResetInitial: Unknown mode %d.",
(int)mode);
}
}
/**
* @brief Scale x vector.
*
* Scale with 1.0 for mode `SCALING_ONES`.
* Scale with 1/fmax(nominal,|xStart|) for mode `SCALING_NOMINALSTART`.
*
* @param data unused
* @param kinsolData
* @param nlsData
* @param mode Mode for scaling. Use `SCALING_NOMINALSTART` for nominal
* scaling and `SCALING_ONES` for no scaling. Will be
* overwritten by simulation flag `FLAG_NO_SCALING`.
*/
static void nlsKinsolXScaling(DATA *data, NLS_KINSOL_DATA *kinsolData,
NONLINEAR_SYSTEM_DATA *nlsData,
scalingMode mode) {
double *xStart = NV_DATA_S(kinsolData->initialGuess);
double *xScaling = NV_DATA_S(kinsolData->xScale);
int i;
/* if noScaling flag is used overwrite mode */
if (omc_flag[FLAG_NO_SCALING]) {
mode = SCALING_ONES;
}
/* Use nominal value or the actual working point for scaling */
switch (mode) {
case SCALING_NOMINALSTART:
for (i = 0; i < nlsData->size; i++) {
xScaling[i] = 1.0 / fmax(nlsData->nominal[i], fabs(xStart[i]));
}
break;
case SCALING_ONES:
for (i = 0; i < nlsData->size; i++) {
xScaling[i] = 1.0;
}
break;
case SCALING_JACOBIAN:
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: Function nlsKinsolXScaling: Invalid mode SCALING_JACOBIAN.");
default:
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: Function nlsKinsolXScaling: Unknown mode %d.", (int)mode);
}
}
/**
* @brief Scale f(x) vector.
*
* @param data
* @param kinsolData
* @param nlsData
* @param mode
*/
static void nlsKinsolFScaling(DATA *data, NLS_KINSOL_DATA *kinsolData,
NONLINEAR_SYSTEM_DATA *nlsData,
scalingMode mode) {
double *fScaling = NV_DATA_S(kinsolData->fScale);
N_Vector x = kinsolData->initialGuess;
int i, j;
int ret;
/* If noScaling flag is used overwrite mode */
if (omc_flag[FLAG_NO_SCALING]) {
mode = SCALING_ONES;
}
/* Use nominal value or the actual working point for scaling */
switch (mode) {
case SCALING_JACOBIAN:
/* Enable scaled jacobian evaluation */
kinsolData->nominalJac = 1;
/* Calculate the scaled Jacobian */
if (nlsData->isPatternAvailable && kinsolData->linearSolverMethod == NLS_LS_KLU) {
if (kinsolData->solved != NLS_SOLVED) {
kinsolData->nominalJac = 0;
if (nlsData->analyticalJacobianColumn != NULL) {
/* Calculate the sparse Jacobian symbolically */
nlsSparseSymJac(x, kinsolData->fTmp, kinsolData->J, kinsolData->userData, NULL, NULL);
} else {
/* Update f(x) for the numerical jacobian matrix */
nlsKinsolResiduals(x, kinsolData->fTmp, kinsolData->userData);
nlsSparseJac(x, kinsolData->fTmp, kinsolData->J, kinsolData->userData, kinsolData->tmp1, kinsolData->tmp2);
}
}
/* Scale the current Jacobian */
ret = _omc_SUNSparseMatrixVecScaling(kinsolData->J, kinsolData->xScale);
if (ret != 0) {
errorStreamPrint(LOG_STDOUT, 0, "KINSOL: _omc_SUNSparseMatrixVecScaling failed.");
}
} else {
/* Update f(x) for the numerical jacobian matrix */
nlsKinsolResiduals(x, kinsolData->fTmp, kinsolData->userData);
nlsDenseJac(nlsData->size, x, kinsolData->fTmp, kinsolData->J,
kinsolData->userData, NULL, NULL);
}
/* Disable scaled Jacobian evaluation */
kinsolData->nominalJac = 0;
for (i = 0; i < nlsData->size; i++) {
fScaling[i] = 1e-12;
}
switch (SUNMatGetID(kinsolData->J))
{
case SUNMATRIX_SPARSE:
for (i = 0; i < SM_NNZ_S(kinsolData->J); ++i) {
if (fScaling[SM_INDEXVALS_S(kinsolData->J)[i]] < fabs(SM_DATA_S(kinsolData->J)[i])) {
fScaling[SM_INDEXVALS_S(kinsolData->J)[i]] = fabs(SM_DATA_S(kinsolData->J)[i]);
}
}
break;
case SUNMATRIX_DENSE:
for (i = 0; i < nlsData->size; i++) {
for (j = 0; j < nlsData->size; j++) {
if (fScaling[i] < fabs(SM_ELEMENT_D(kinsolData->J, j, i))) {
fScaling[i] = fabs(SM_ELEMENT_D(kinsolData->J, j, i));
}
}
}
break;
default:
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: Function nlsKinsolFScaling: Unknown matrix type.");
}
/* inverse fScale */
N_VInv(kinsolData->fScale, kinsolData->fScale);
break;
case SCALING_ONES:
for (i = 0; i < nlsData->size; i++) {
fScaling[i] = 1.0;
}
break;
case SCALING_NOMINALSTART:
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: Function nlsKinsolFScaling: Invalid mode SCALING_NOMINALSTART.");
default:
errorStreamPrint(LOG_STDOUT, 0,
"KINSOL: Function nlsKinsolFScaling: Unknown mode %d.", (int)mode);
}
}
/**
* @brief Print KINSOL configuration.
*
* Only prints if stream `LOG_NLS_V` is active.
*
* @param kinsolData
* @param nlsData
*/
static void nlsKinsolConfigPrint(NLS_KINSOL_DATA *kinsolData,
NONLINEAR_SYSTEM_DATA *nlsData) {
int retValue;
double fNorm;
DATA *data = kinsolData->userData->data;
int eqSystemNumber = nlsData->equationIndex;
_omc_vector vecStart, vecXScaling, vecFScaling;
if (!useStream[LOG_NLS_V]) {
return;
}
_omc_initVector(&vecStart, kinsolData->size,
NV_DATA_S(kinsolData->initialGuess));
_omc_initVector(&vecXScaling, kinsolData->size,
NV_DATA_S(kinsolData->xScale));
_omc_initVector(&vecFScaling, kinsolData->size,
NV_DATA_S(kinsolData->fScale));
infoStreamPrint(LOG_NLS_V, 1, "Kinsol Configuration");
_omc_printVectorWithEquationInfo(
&vecStart, "Initial guess values", LOG_NLS_V,
modelInfoGetEquation(&data->modelData->modelDataXml, eqSystemNumber));
_omc_printVectorWithEquationInfo(
&vecXScaling, "xScaling", LOG_NLS_V,
modelInfoGetEquation(&data->modelData->modelDataXml, eqSystemNumber));
_omc_printVector(&vecFScaling, "fScaling", LOG_NLS_V);
infoStreamPrint(LOG_NLS_V, 0, "KINSOL F tolerance: %g", kinsolData->fnormtol);
infoStreamPrint(LOG_NLS_V, 0, "KINSOL minimal step size %g",
kinsolData->scsteptol);
infoStreamPrint(LOG_NLS_V, 0, "KINSOL max iterations %d",
20 * kinsolData->size);
infoStreamPrint(LOG_NLS_V, 0, "KINSOL strategy %d",
kinsolData->kinsolStrategy);
infoStreamPrint(LOG_NLS_V, 0, "KINSOL current retry %d", kinsolData->retries);
infoStreamPrint(LOG_NLS_V, 0, "KINSOL max step %g", kinsolData->mxnstepin);
infoStreamPrint(LOG_NLS_V, 0, "KINSOL linear solver %d",
kinsolData->linearSolverMethod);