-
Notifications
You must be signed in to change notification settings - Fork 28
/
Copy pathTractorConvectionDispersionOperator.cpp
1526 lines (1286 loc) · 59.6 KB
/
TractorConvectionDispersionOperator.cpp
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
// =============================================================================
// CADET
//
// Copyright © 2008-2021: The CADET Authors
// Please see the AUTHORS and CONTRIBUTORS file.
//
// All rights reserved. This program and the accompanying materials
// are made available under the terms of the GNU Public License v3.0 (or, at
// your option, any later version) which accompanies this distribution, and
// is available at http://www.gnu.org/licenses/gpl.html
// =============================================================================
#include "model/parts/TractorConvectionDispersionOperator.hpp"
#include "cadet/Exceptions.hpp"
#include "Stencil.hpp"
#include "ParamReaderHelper.hpp"
#include "AdUtils.hpp"
#include "SensParamUtil.hpp"
#include "SimulationTypes.hpp"
#include "linalg/CompressedSparseMatrix.hpp"
#include "model/parts/ConvectionDispersionKernel.hpp"
#ifdef SUPERLU_FOUND
#include "linalg/SuperLUSparseMatrix.hpp"
#endif
#ifdef UMFPACK_FOUND
#include "linalg/UMFPackSparseMatrix.hpp"
#endif
#include "linalg/BandMatrix.hpp"
#include "linalg/Gmres.hpp"
#include "LoggingUtils.hpp"
#include "Logging.hpp"
#include <algorithm>
namespace
{
cadet::model::MultiplexMode readAndRegisterMultiplexParam(cadet::IParameterProvider& paramProvider, std::unordered_map<cadet::ParameterId, cadet::active*>& parameters, std::vector<cadet::active>& values, const std::string& name, unsigned int nComp, unsigned int nRad, cadet::UnitOpIdx uoi)
{
cadet::model::MultiplexMode mode = cadet::model::MultiplexMode::Independent;
readParameterMatrix(values, paramProvider, name, nComp * nRad, 1);
unsigned int nSec = 1;
if (paramProvider.exists(name + "_MULTIPLEX"))
{
const int modeConfig = paramProvider.getInt(name + "_MULTIPLEX");
if (modeConfig == 0)
{
mode = cadet::model::MultiplexMode::Independent;
if (values.size() > 1)
throw cadet::InvalidParameterException("Number of elements in field " + name + " inconsistent with " + name + "_MULTIPLEX (should be 1)");
}
else if (modeConfig == 1)
{
mode = cadet::model::MultiplexMode::Radial;
if (values.size() != nRad)
throw cadet::InvalidParameterException("Number of elements in field " + name + " inconsistent with " + name + "_MULTIPLEX (should be " + std::to_string(nRad) + ")");
}
else if (modeConfig == 2)
{
mode = cadet::model::MultiplexMode::Component;
if (values.size() != nComp)
throw cadet::InvalidParameterException("Number of elements in field " + name + " inconsistent with " + name + "_MULTIPLEX (should be " + std::to_string(nComp) + ")");
}
else if (modeConfig == 3)
{
mode = cadet::model::MultiplexMode::ComponentRadial;
if (values.size() != nComp * nRad)
throw cadet::InvalidParameterException("Number of elements in field " + name + " inconsistent with " + name + "_MULTIPLEX (should be " + std::to_string(nComp * nRad) + ")");
}
else if (modeConfig == 4)
{
mode = cadet::model::MultiplexMode::Section;
nSec = values.size();
}
else if (modeConfig == 5)
{
mode = cadet::model::MultiplexMode::RadialSection;
if (values.size() % nRad != 0)
throw cadet::InvalidParameterException("Number of elements in field " + name + " is not a positive multiple of NRAD (" + std::to_string(nRad) + ")");
nSec = values.size() / nRad;
}
else if (modeConfig == 6)
{
mode = cadet::model::MultiplexMode::ComponentSection;
if (values.size() % nComp != 0)
throw cadet::InvalidParameterException("Number of elements in field " + name + " is not a positive multiple of NCOMP (" + std::to_string(nComp) + ")");
nSec = values.size() / nComp;
}
else if (modeConfig == 7)
{
mode = cadet::model::MultiplexMode::ComponentRadialSection;
if (values.size() % (nComp * nRad) != 0)
throw cadet::InvalidParameterException("Number of elements in field " + name + " is not a positive multiple of NCOMP * NRAD (" + std::to_string(nComp * nRad) + ")");
nSec = values.size() / (nComp * nRad);
}
}
else
{
if (values.size() == 1)
mode = cadet::model::MultiplexMode::Independent;
else if (values.size() == nComp)
mode = cadet::model::MultiplexMode::Component;
else if (values.size() == nRad)
mode = cadet::model::MultiplexMode::Radial;
else if (values.size() == nRad * nComp)
mode = cadet::model::MultiplexMode::ComponentRadial;
else if (values.size() % nComp == 0)
{
mode = cadet::model::MultiplexMode::ComponentSection;
nSec = values.size() / nComp;
}
else if (values.size() % nRad == 0)
{
mode = cadet::model::MultiplexMode::RadialSection;
nSec = values.size() / nRad;
}
else if (values.size() % (nRad * nComp) == 0)
{
mode = cadet::model::MultiplexMode::ComponentRadialSection;
nSec = values.size() / (nComp * nRad);
}
else
throw cadet::InvalidParameterException("Could not infer multiplex mode of field " + name + ", set " + name + "_MULTIPLEX or change number of elements");
// Do not infer cadet::model::MultiplexMode::Section in case of no matches (might hide specification errors)
}
const cadet::StringHash nameHash = cadet::hashStringRuntime(name);
switch (mode)
{
case cadet::model::MultiplexMode::Independent:
case cadet::model::MultiplexMode::Section:
{
std::vector<cadet::active> p(nComp * nRad * nSec);
for (unsigned int s = 0; s < nSec; ++s)
std::fill(p.begin() + s * nRad * nComp, p.begin() + (s+1) * nRad * nComp, values[s]);
values = std::move(p);
for (unsigned int s = 0; s < nSec; ++s)
parameters[cadet::makeParamId(nameHash, uoi, cadet::CompIndep, cadet::ParTypeIndep, cadet::BoundStateIndep, cadet::ReactionIndep, (mode == cadet::model::MultiplexMode::Independent) ? cadet::SectionIndep : s)] = &values[s * nRad * nComp];
}
break;
case cadet::model::MultiplexMode::Component:
case cadet::model::MultiplexMode::ComponentSection:
{
std::vector<cadet::active> p(nComp * nRad * nSec);
for (unsigned int s = 0; s < nSec; ++s)
{
for (unsigned int i = 0; i < nComp; ++i)
std::copy(values.begin() + s * nComp, values.begin() + (s+1) * nComp, p.begin() + i * nComp + s * nComp * nRad);
}
values = std::move(p);
for (unsigned int s = 0; s < nSec; ++s)
{
for (unsigned int i = 0; i < nComp; ++i)
parameters[cadet::makeParamId(nameHash, uoi, i, cadet::ParTypeIndep, cadet::BoundStateIndep, cadet::ReactionIndep, (mode == cadet::model::MultiplexMode::Component) ? cadet::SectionIndep : s)] = &values[s * nRad * nComp + i];
}
}
break;
case cadet::model::MultiplexMode::Radial:
case cadet::model::MultiplexMode::RadialSection:
{
std::vector<cadet::active> p(nComp * nRad * nSec);
for (unsigned int i = 0; i < nRad * nSec; ++i)
std::fill(p.begin() + i * nComp, p.begin() + (i+1) * nComp, values[i]);
values = std::move(p);
for (unsigned int s = 0; s < nSec; ++s)
{
for (unsigned int i = 0; i < nRad; ++i)
parameters[cadet::makeParamId(nameHash, uoi, cadet::CompIndep, i, cadet::BoundStateIndep, cadet::ReactionIndep, (mode == cadet::model::MultiplexMode::Radial) ? cadet::SectionIndep : s)] = &values[s * nRad * nComp + i * nComp];
}
}
break;
case cadet::model::MultiplexMode::ComponentRadial:
case cadet::model::MultiplexMode::ComponentRadialSection:
cadet::registerParam3DArray(parameters, values, [=](bool multi, unsigned int sec, unsigned int compartment, unsigned int comp) { return cadet::makeParamId(nameHash, uoi, comp, compartment, cadet::BoundStateIndep, cadet::ReactionIndep, multi ? sec : cadet::SectionIndep); }, nComp, nRad);
break;
case cadet::model::MultiplexMode::Axial:
case cadet::model::MultiplexMode::AxialRadial:
case cadet::model::MultiplexMode::Type:
case cadet::model::MultiplexMode::ComponentType:
case cadet::model::MultiplexMode::ComponentSectionType:
cadet_assert(false);
break;
}
return mode;
}
bool multiplexParameterValue(const cadet::ParameterId& pId, cadet::StringHash nameHash, cadet::model::MultiplexMode mode, std::vector<cadet::active>& data, unsigned int nComp, unsigned int nRad, double value, std::unordered_set<cadet::active*> const* sensParams)
{
if (pId.name != nameHash)
return false;
switch (mode)
{
case cadet::model::MultiplexMode::Independent:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[0]))
return false;
for (std::size_t i = 0; i < data.size(); ++i)
data[i].setValue(value);
return true;
}
case cadet::model::MultiplexMode::Section:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[pId.section * nComp * nRad]))
return false;
for (unsigned int i = 0; i < nComp * nRad; ++i)
data[i + pId.section * nComp * nRad].setValue(value);
return true;
}
case cadet::model::MultiplexMode::Component:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[pId.component]))
return false;
for (unsigned int i = 0; i < nRad; ++i)
data[i * nComp + pId.component].setValue(value);
return true;
}
case cadet::model::MultiplexMode::ComponentSection:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[pId.component + pId.section * nComp * nRad]))
return false;
for (unsigned int i = 0; i < nRad; ++i)
data[i * nComp + pId.component + pId.section * nComp * nRad].setValue(value);
return true;
}
case cadet::model::MultiplexMode::Radial:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[pId.particleType * nComp]))
return false;
for (unsigned int i = 0; i < nComp; ++i)
data[i + pId.particleType * nComp].setValue(value);
return true;
}
case cadet::model::MultiplexMode::RadialSection:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[pId.particleType * nComp + pId.section * nComp * nRad]))
return false;
for (unsigned int i = 0; i < nComp; ++i)
data[i + pId.particleType * nComp + pId.section * nComp * nRad].setValue(value);
return true;
}
case cadet::model::MultiplexMode::ComponentRadial:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[pId.component + pId.particleType * nComp]))
return false;
data[pId.component + pId.particleType * nComp].setValue(value);
return true;
}
case cadet::model::MultiplexMode::ComponentRadialSection:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
if (sensParams && !cadet::contains(*sensParams, &data[pId.component + pId.particleType * nComp + pId.section * nComp * nRad]))
return false;
data[pId.component + pId.particleType * nComp + pId.section * nComp * nRad].setValue(value);
return true;
}
case cadet::model::MultiplexMode::Axial:
case cadet::model::MultiplexMode::AxialRadial:
case cadet::model::MultiplexMode::Type:
case cadet::model::MultiplexMode::ComponentType:
case cadet::model::MultiplexMode::ComponentSectionType:
cadet_assert(false);
break;
}
return false;
}
bool multiplexParameterAD(const cadet::ParameterId& pId, cadet::StringHash nameHash, cadet::model::MultiplexMode mode, std::vector<cadet::active>& data, unsigned int nComp, unsigned int nRad, unsigned int adDirection, double adValue, std::unordered_set<cadet::active*>& sensParams)
{
if (pId.name != nameHash)
return false;
switch (mode)
{
case cadet::model::MultiplexMode::Independent:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
sensParams.insert(&data[0]);
for (std::size_t i = 0; i < data.size(); ++i)
data[i].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::Section:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
sensParams.insert(&data[pId.section * nComp * nRad]);
for (unsigned int i = 0; i < nComp * nRad; ++i)
data[i + pId.section * nComp * nRad].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::Component:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
sensParams.insert(&data[pId.component]);
for (unsigned int i = 0; i < nRad; ++i)
data[i * nComp + pId.component].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::ComponentSection:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType != cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
sensParams.insert(&data[pId.component + pId.section * nComp * nRad]);
for (unsigned int i = 0; i < nRad; ++i)
data[i * nComp + pId.component + pId.section * nComp * nRad].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::Radial:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
sensParams.insert(&data[pId.particleType * nComp]);
for (unsigned int i = 0; i < nComp; ++i)
data[i + pId.particleType * nComp].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::RadialSection:
{
if ((pId.component != cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
sensParams.insert(&data[pId.particleType * nComp + pId.section * nComp * nRad]);
for (unsigned int i = 0; i < nComp; ++i)
data[i + pId.particleType * nComp + pId.section * nComp * nRad].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::ComponentRadial:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section != cadet::SectionIndep))
return false;
sensParams.insert(&data[pId.component + pId.particleType * nComp]);
data[pId.component + pId.particleType * nComp].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::ComponentRadialSection:
{
if ((pId.component == cadet::CompIndep) || (pId.particleType == cadet::ParTypeIndep) || (pId.boundState != cadet::BoundStateIndep)
|| (pId.reaction != cadet::ReactionIndep) || (pId.section == cadet::SectionIndep))
return false;
sensParams.insert(&data[pId.component + pId.particleType * nComp + pId.section * nComp * nRad]);
data[pId.component + pId.particleType * nComp + pId.section * nComp * nRad].setADValue(adDirection, adValue);
return true;
}
case cadet::model::MultiplexMode::Axial:
case cadet::model::MultiplexMode::AxialRadial:
case cadet::model::MultiplexMode::Type:
case cadet::model::MultiplexMode::ComponentType:
case cadet::model::MultiplexMode::ComponentSectionType:
cadet_assert(false);
break;
}
return false;
}
} // namespace
namespace cadet
{
namespace model
{
namespace parts
{
class TractorConvectionDispersionOperator::LinearSolver
{
public:
virtual ~LinearSolver() CADET_NOEXCEPT { }
virtual bool initialize(IParameterProvider& paramProvider, unsigned int nComp, unsigned int nCol, unsigned int nRad, const Weno& weno) = 0;
virtual void setSparsityPattern(const cadet::linalg::SparsityPattern& pattern) = 0;
virtual void assembleDiscretizedJacobian(double alpha) = 0;
virtual bool factorize() = 0;
virtual bool solveDiscretizedJacobian(double* rhs, double const* weight, double const* init, double outerTol) const = 0;
};
int schurComplementMultiplierTractorCDO(void* userData, double const* x, double* z);
class TractorConvectionDispersionOperator::GmresSolver : public TractorConvectionDispersionOperator::LinearSolver
{
public:
GmresSolver(linalg::CompressedSparseMatrix const* jacC) : _jacC(jacC) { }
virtual ~GmresSolver() CADET_NOEXCEPT { }
virtual bool initialize(IParameterProvider& paramProvider, unsigned int nComp, unsigned int nCol, unsigned int nRad, const Weno& weno)
{
_gmres.initialize(nCol * nComp * nRad, 0, linalg::toOrthogonalization(1), 0);
_gmres.matrixVectorMultiplier(&schurComplementMultiplierTractorCDO, this);
_cache.resize(nCol * nComp * nRad, 0.0);
return true;
}
virtual void setSparsityPattern(const linalg::SparsityPattern& pattern) { }
virtual void assembleDiscretizedJacobian(double alpha)
{
_alpha = alpha;
}
virtual bool factorize() { return true; }
virtual bool solveDiscretizedJacobian(double* rhs, double const* weight, double const* init, double outerTol) const
{
if (init)
std::copy(init, init + _cache.size(), _cache.begin());
const int gmresResult = _gmres.solve(0.05 * outerTol * std::sqrt(_jacC->rows()), weight, rhs, _cache.data());
std::copy(_cache.begin(), _cache.end(), rhs);
return gmresResult == 0;
}
protected:
linalg::CompressedSparseMatrix const* const _jacC;
double _alpha;
mutable linalg::Gmres _gmres; //!< GMRES algorithm for the Schur-complement in linearSolve()
mutable std::vector<double> _cache; //!< GMRES cache for result
int schurComplementMatrixVector(double const* x, double* z) const
{
std::fill(z, z + _jacC->rows(), _alpha);
_jacC->multiplyVector(x, 1.0, 1.0, z);
return 0;
}
// Wrapper for calling the corresponding function in this class
friend int schurComplementMultiplierTractorCDO(void* userData, double const* x, double* z);
};
int schurComplementMultiplierTractorCDO(void* userData, double const* x, double* z)
{
TractorConvectionDispersionOperator::GmresSolver* const cdo = static_cast<TractorConvectionDispersionOperator::GmresSolver*>(userData);
return cdo->schurComplementMatrixVector(x, z);
}
#if defined(UMFPACK_FOUND) || defined(SUPERLU_FOUND)
template <typename sparse_t>
class TractorConvectionDispersionOperator::SparseDirectSolver : public TractorConvectionDispersionOperator::LinearSolver
{
public:
SparseDirectSolver(linalg::CompressedSparseMatrix const* jacC) : _jacC(jacC) { }
virtual ~SparseDirectSolver() CADET_NOEXCEPT { }
virtual bool initialize(IParameterProvider& paramProvider, unsigned int nComp, unsigned int nCol, unsigned int nRad, const Weno& weno)
{
return true;
}
virtual void setSparsityPattern(const linalg::SparsityPattern& pattern)
{
_jacCdisc.assignPattern(pattern);
_jacCdisc.prepare();
}
virtual void assembleDiscretizedJacobian(double alpha)
{
// Copy normal matrix over to factorizable matrix
_jacCdisc.copyFromSamePattern(*_jacC);
for (int i = 0; i < _jacC->rows(); ++i)
_jacCdisc.centered(i, 0) += alpha;
}
virtual bool factorize()
{
return _jacCdisc.factorize();
}
virtual bool solveDiscretizedJacobian(double* rhs, double const* weight, double const* init, double outerTol) const
{
return _jacCdisc.solve(rhs);
}
protected:
linalg::CompressedSparseMatrix const* const _jacC;
sparse_t _jacCdisc;
};
#endif
class TractorConvectionDispersionOperator::DenseDirectSolver : public TractorConvectionDispersionOperator::LinearSolver
{
public:
DenseDirectSolver(linalg::CompressedSparseMatrix const* jacC) : _jacC(jacC) { }
virtual ~DenseDirectSolver() CADET_NOEXCEPT { }
virtual bool initialize(IParameterProvider& paramProvider, unsigned int nComp, unsigned int nCol, unsigned int nRad, const Weno& weno)
{
// Note that we have to increase the lower bandwidth by 1 because the WENO stencil is applied to the
// right cell face (lower + 1 + upper) and to the left cell face (shift the stencil by -1 because influx of cell i
// is outflux of cell i-1)
// We also have to make sure that there's at least one sub and super diagonal for the dispersion term
const unsigned int lb = std::max(weno.lowerBandwidth() + 1u, 1u) * nComp * nRad;
// We have to make sure that there's at least one sub and super diagonal for the dispersion term
const unsigned int ub = std::max(weno.upperBandwidth(), 1u) * nComp * nRad;
// When flow direction is changed, the bandwidths of the Jacobian swap.
// Hence, we have to reserve memory such that the swapped Jacobian can fit into the matrix.
const unsigned int mb = std::max(lb, ub);
// Allocate matrices such that bandwidths can be switched (backwards flow support)
_jacCdisc.resize(nCol * nComp * nRad, mb, mb);
return true;
}
virtual void setSparsityPattern(const linalg::SparsityPattern& pattern) { }
virtual void assembleDiscretizedJacobian(double alpha)
{
// Copy normal matrix over to factorizable matrix
_jacCdisc.setAll(0.0);
linalg::FactorizableBandMatrix::RowIterator jac = _jacCdisc.row(0);
for (int i = 0; i < _jacC->rows(); ++i, ++jac)
{
linalg::sparse_int_t const* const colIdx = _jacC->columnIndicesOfRow(i);
double const* const vals = _jacC->valuesOfRow(i);
const int nnz = _jacC->numNonZerosInRow(i);
// Copy row from sparse matrix to banded matrix
for (int c = 0; c < nnz; ++c)
{
const linalg::sparse_int_t diag = colIdx[c] - i;
jac[diag] = vals[c];
}
// Add time derivative
jac[0] += alpha;
}
}
virtual bool factorize()
{
return _jacCdisc.factorize();
}
virtual bool solveDiscretizedJacobian(double* rhs, double const* weight, double const* init, double outerTol) const
{
return _jacCdisc.solve(rhs);
}
protected:
linalg::CompressedSparseMatrix const* const _jacC;
linalg::FactorizableBandMatrix _jacCdisc;
};
/**
* @brief Creates a TractorConvectionDispersionOperator
*/
TractorConvectionDispersionOperator::TractorConvectionDispersionOperator() : _colPorosities(0), _stencilMemory(sizeof(active) * Weno::maxStencilSize()),
_wenoDerivatives(new double[Weno::maxStencilSize()]), _weno(), _linearSolver(nullptr)
{
}
TractorConvectionDispersionOperator::~TractorConvectionDispersionOperator() CADET_NOEXCEPT
{
delete[] _wenoDerivatives;
delete _linearSolver;
}
/**
* @brief Reads parameters and allocates memory
* @details Has to be called once before the operator is used.
* @param [in] paramProvider Parameter provider for reading parameters
* @param [in] nComp Number of components
* @param [in] nCol Number of axial cells
* @param [in] dynamicReactions Determines whether the sparsity pattern accounts for dynamic reactions
* @return @c true if configuration went fine, @c false otherwise
*/
bool TractorConvectionDispersionOperator::configureModelDiscretization(IParameterProvider& paramProvider, unsigned int nComp, unsigned int nCol, unsigned int nRad, bool dynamicReactions)
{
_nComp = nComp;
_nCol = nCol;
_nRad = nRad;
_hasDynamicReactions = dynamicReactions;
readParameterMatrix(_exchangeMatrix, paramProvider, "EXCHANGE_MATRIX", _nRad * _nRad * _nComp, 1);
paramProvider.pushScope("discretization");
// Read WENO settings and apply them
paramProvider.pushScope("weno");
_weno.order(paramProvider.getInt("WENO_ORDER"));
_weno.boundaryTreatment(paramProvider.getInt("BOUNDARY_MODEL"));
_wenoEpsilon = paramProvider.getDouble("WENO_EPS");
paramProvider.popScope();
// Read solver settings
if (paramProvider.exists("LINEAR_SOLVER_BULK"))
{
const std::string sol = paramProvider.getString("LINEAR_SOLVER_BULK");
if (sol == "DENSE")
_linearSolver = new DenseDirectSolver(&_jacC);
else if (sol == "GMRES")
_linearSolver = new GmresSolver(&_jacC);
#ifdef UMFPACK_FOUND
else if (sol == "UMFPACK")
_linearSolver = new SparseDirectSolver<linalg::UMFPackSparseMatrix>(&_jacC);
#endif
#ifdef SUPERLU_FOUND
else if (sol == "SUPERLU")
_linearSolver = new SparseDirectSolver<linalg::SuperLUSparseMatrix>(&_jacC);
#endif
else
throw InvalidParameterException("Unknown linear solver " + sol + " in field LINEAR_SOLVER_BULK");
}
// Default to sparse solver if available (preferably UMFPACK), fall back to dense
if (!_linearSolver)
{
#if defined(UMFPACK_FOUND)
_linearSolver = new SparseDirectSolver<linalg::UMFPackSparseMatrix>(&_jacC);
#elif defined(SUPERLU_FOUND)
_linearSolver = new SparseDirectSolver<linalg::SuperLUSparseMatrix>(&_jacC);
#else
_linearSolver = new DenseDirectSolver(&_jacC);
#endif
// LOG(Info) << "Default to dense banded linear solver due to invalid or missing LINEAR_SOLVER_BULK setting";
}
paramProvider.popScope();
_radialEdges.resize(nRad + 1);
_radialCenters.resize(nRad);
// _radialCentroids.resize(nRad);
_crossSections.resize(nRad);
_curVelocity.resize(nRad);
setSparsityPattern<active>();
return _linearSolver->initialize(paramProvider, nComp, nCol, nRad, _weno);
}
/**
* @brief Reads model parameters
* @details Only reads parameters that do not affect model structure (e.g., discretization).
* @param [in] unitOpIdx Unit operation id of the owning unit operation model
* @param [in] paramProvider Parameter provider for reading parameters
* @param [out] parameters Map in which local parameters are inserted
* @return @c true if configuration went fine, @c false otherwise
*/
bool TractorConvectionDispersionOperator::configure(UnitOpIdx unitOpIdx, IParameterProvider& paramProvider, std::unordered_map<ParameterId, active*>& parameters)
{
// Read geometry parameters
_colLength = paramProvider.getDouble("COL_LENGTH");
_colRadius = paramProvider.getDouble("COL_RADIUS");
//_exchangeMatrix = paramProvider.getDoubleArray("EXCHANGE_MATRIX");
readScalarParameterOrArray(_colPorosities, paramProvider, "COL_POROSITY", 1);
//readParameterMatrix(_exchangeMatrix, paramProvider, "EXCHANGE_MATRIX", _nRad * _nRad * _nComp, 1);
//_exchangeMatrix = readAndRegisterMultiplexParam(paramProvider, parameters, _radialDispersion, "COL_DISPERSION_RADIAL", _nComp, _nRad, unitOpIdx);
if ((_colPorosities.size() != 1) && (_colPorosities.size() != _nRad))
throw InvalidParameterException("Number of elements in field COL_POROSITY is neither 1 nor NRAD (" + std::to_string(_nRad) + ")");
_singlePorosity = (_colPorosities.size() == 1);
if (_singlePorosity)
_colPorosities = std::vector<active>(_nRad, _colPorosities[0]);
// Read radial discretization mode and default to "EQUIDISTANT"
paramProvider.pushScope("discretization");
const std::string rdt = paramProvider.getString("RADIAL_DISC_TYPE");
if (rdt == "EQUIVOLUME")
_radialDiscretizationMode = RadialDiscretizationMode::Equivolume;
else if (rdt == "USER_DEFINED")
{
_radialDiscretizationMode = RadialDiscretizationMode::UserDefined;
readScalarParameterOrArray(_radialEdges, paramProvider, "RADIAL_COMPARTMENTS", 1);
if (_radialEdges.size() < _nRad + 1)
throw InvalidParameterException("Number of elements in field RADIAL_COMPARTMENTS is less than NRAD + 1 (" + std::to_string(_nRad + 1) + ")");
registerParam1DArray(parameters, _radialEdges, [=](bool multi, unsigned int i) { return makeParamId(hashString("RADIAL_COMPARTMENTS"), unitOpIdx, CompIndep, i, BoundStateIndep, ReactionIndep, SectionIndep); });
}
else
_radialDiscretizationMode = RadialDiscretizationMode::Equidistant;
paramProvider.popScope();
updateRadialDisc();
// Read section dependent parameters (transport)
// Read VELOCITY
_velocity.clear();
if (paramProvider.exists("VELOCITY"))
{
readScalarParameterOrArray(_velocity, paramProvider, "VELOCITY", 1);
if (paramProvider.exists("VELOCITY_MULTIPLEX"))
{
const int mode = paramProvider.getInt("VELOCITY_MULTIPLEX");
if (mode == 0)
// Rad-indep, sec-indep
_singleVelocity = true;
else if (mode == 1)
// Rad-dep, sec-indep
_singleVelocity = false;
else if (mode == 2)
// Rad-indep, sec-dep
_singleVelocity = true;
else if (mode == 3)
// Rad-dep, sec-dep
_singleVelocity = false;
if (!_singleVelocity && (_velocity.size() % _nRad != 0))
throw InvalidParameterException("Number of elements in field VELOCITY is not a positive multiple of NRAD (" + std::to_string(_nRad) + ")");
if ((mode == 0) && (_velocity.size() != 1))
throw InvalidParameterException("Number of elements in field VELOCITY inconsistent with VELOCITY_MULTIPLEX (should be 1)");
if ((mode == 1) && (_velocity.size() != _nRad))
throw InvalidParameterException("Number of elements in field VELOCITY inconsistent with VELOCITY_MULTIPLEX (should be " + std::to_string(_nRad) + ")");
}
else
{
// Infer radial dependence of VELOCITY:
// size not divisible by NRAD -> radial independent
_singleVelocity = ((_velocity.size() % _nRad) != 0);
}
// Expand _velocity to make it component dependent
if (_singleVelocity)
{
std::vector<active> expanded(_velocity.size() * _nRad);
for (std::size_t i = 0; i < _velocity.size(); ++i)
std::fill(expanded.begin() + i * _nRad, expanded.begin() + (i + 1) * _nRad, _velocity[i]);
_velocity = std::move(expanded);
}
}
else
{
_singleVelocity = false;
_velocity.resize(_nRad, 1.0);
}
// Register VELOCITY
if (_singleVelocity)
{
if (_velocity.size() > _nRad)
{
// Register only the first item in each section
for (std::size_t i = 0; i < _velocity.size() / _nRad; ++i)
parameters[makeParamId(hashString("VELOCITY"), unitOpIdx, CompIndep, ParTypeIndep, BoundStateIndep, ReactionIndep, i)] = &_velocity[i * _nRad];
}
else
{
// We have only one parameter
parameters[makeParamId(hashString("VELOCITY"), unitOpIdx, CompIndep, ParTypeIndep, BoundStateIndep, ReactionIndep, SectionIndep)] = &_velocity[0];
}
}
else
registerParam2DArray(parameters, _velocity, [=](bool multi, unsigned int sec, unsigned int compartment) { return makeParamId(hashString("VELOCITY"), unitOpIdx, CompIndep, compartment, BoundStateIndep, ReactionIndep, multi ? sec : SectionIndep); }, _nRad);
_axialDispersionMode = readAndRegisterMultiplexParam(paramProvider, parameters, _axialDispersion, "COL_DISPERSION", _nComp, _nRad, unitOpIdx);
_radialDispersionMode = readAndRegisterMultiplexParam(paramProvider, parameters, _radialDispersion, "COL_DISPERSION_RADIAL", _nComp, _nRad, unitOpIdx);
// Add parameters to map
parameters[makeParamId(hashString("COL_LENGTH"), unitOpIdx, CompIndep, ParTypeIndep, BoundStateIndep, ReactionIndep, SectionIndep)] = &_colLength;
parameters[makeParamId(hashString("COL_RADIUS"), unitOpIdx, CompIndep, ParTypeIndep, BoundStateIndep, ReactionIndep, SectionIndep)] = &_colRadius;
registerParam1DArray(parameters, _colPorosities, [=](bool multi, unsigned int i) { return makeParamId(hashString("COL_POROSITY"), unitOpIdx, CompIndep, multi ? i : ParTypeIndep, BoundStateIndep, ReactionIndep, SectionIndep); });
return true;
}
/**
* @brief Notifies the operator that a discontinuous section transition is in progress
* @details In addition to changing flow direction internally, if necessary, the function returns whether
* the flow direction has changed.
* @param [in] t Current time point
* @param [in] secIdx Index of the new section that is about to be integrated
* @return @c true if flow direction has changed, otherwise @c false
*/
bool TractorConvectionDispersionOperator::notifyDiscontinuousSectionTransition(double t, unsigned int secIdx)
{
bool hasChanged = false;
if (!_velocity.empty())
{
// _curVelocity has already been set to the network flow rate in setFlowRates()
// the direction of the flow (i.e., sign of _curVelocity) is given by _velocity
active const* const dirs = getSectionDependentSlice(_velocity, _nRad, secIdx);
if (secIdx > 0)
{
active const* const dirsPrev = getSectionDependentSlice(_velocity, _nRad, secIdx - 1);
for (unsigned int i = 0; i < _nRad; ++i)
{
if (dirs[i] * dirsPrev[i] < 0.0)
{
hasChanged = true;
break;
}
}
}
for (unsigned int i = 0; i < _nRad; ++i)
{
if (dirs[i] < 0.0)
_curVelocity[i] *= -1.0;
}
}
// Change the sparsity pattern if necessary
if ((secIdx == 0) || hasChanged)
setSparsityPattern<active>();
return hasChanged || (secIdx == 0);
}
/**
* @brief Sets the flow rates for the current time section
* @details The flow rates may change due to valve switches.
* @param [in] compartment Index of the compartment
* @param [in] in Total volumetric inlet flow rate
* @param [in] out Total volumetric outlet flow rate
*/
void TractorConvectionDispersionOperator::setFlowRates(int compartment, const active& in, const active& out) CADET_NOEXCEPT
{
_curVelocity[compartment] = in / (_crossSections[compartment] * _colPorosities[compartment]);
}
void TractorConvectionDispersionOperator::setFlowRates(active const* in, active const* out) CADET_NOEXCEPT
{
for (unsigned int compartment = 0; compartment < _nRad; ++compartment)
_curVelocity[compartment] = in[compartment] / (_crossSections[compartment] * _colPorosities[compartment]);
}
double TractorConvectionDispersionOperator::inletFactor(unsigned int idxSec, int idxRad) const CADET_NOEXCEPT
{
const double h = static_cast<double>(_colLength) / static_cast<double>(_nCol);
return -std::abs(static_cast<double>(_curVelocity[idxRad])) / h;
}
const active& TractorConvectionDispersionOperator::axialDispersion(unsigned int idxSec, int idxRad, int idxComp) const CADET_NOEXCEPT
{
return *(getSectionDependentSlice(_axialDispersion, _nRad * _nComp, idxSec) + idxRad * _nComp + idxComp);
}
const active& TractorConvectionDispersionOperator::radialDispersion(unsigned int idxSec, int idxRad, int idxComp) const CADET_NOEXCEPT
{
return *(getSectionDependentSlice(_radialDispersion, _nRad * _nComp, idxSec) + idxRad * _nComp + idxComp);
}
/**
* @brief Computes the residual of the transport equations
* @param [in] t Current time point
* @param [in] secIdx Index of the current section
* @param [in] y Pointer to unit operation's state vector
* @param [in] yDot Pointer to unit operation's time derivative state vector
* @param [out] res Pointer to unit operation's residual vector
* @param [in] wantJac Determines whether the Jacobian is computed or not
* @return @c 0 on success, @c -1 on non-recoverable error, and @c +1 on recoverable error
*/
int TractorConvectionDispersionOperator::residual(double t, unsigned int secIdx, double const* y, double const* yDot, double* res, bool wantJac, WithoutParamSensitivity)
{
if (wantJac)
return residualImpl<double, double, double, true>(t, secIdx, y, yDot, res);
else
return residualImpl<double, double, double, false>(t, secIdx, y, yDot, res);
}
int TractorConvectionDispersionOperator::residual(double t, unsigned int secIdx, active const* y, double const* yDot, active* res, bool wantJac, WithoutParamSensitivity)
{
if (wantJac)
return residualImpl<active, active, double, true>(t, secIdx, y, yDot, res);
else
return residualImpl<active, active, double, false>(t, secIdx, y, yDot, res);
}
int TractorConvectionDispersionOperator::residual(double t, unsigned int secIdx, double const* y, double const* yDot, active* res, bool wantJac, WithParamSensitivity)
{
if (wantJac)
return residualImpl<double, active, active, true>(t, secIdx, y, yDot, res);
else
return residualImpl<double, active, active, false>(t, secIdx, y, yDot, res);
}
int TractorConvectionDispersionOperator::residual(double t, unsigned int secIdx, active const* y, double const* yDot, active* res, bool wantJac, WithParamSensitivity)
{
if (wantJac)
return residualImpl<active, active, active, true>(t, secIdx, y, yDot, res);
else
return residualImpl<active, active, active, false>(t, secIdx, y, yDot, res);
}
template <typename StateType, typename ResidualType, typename ParamType, bool wantJac>
int TractorConvectionDispersionOperator::residualImpl(double t, unsigned int secIdx, StateType const* y, double const* yDot, ResidualType* res)
{
if (wantJac)
{
// Reset Jacobian
_jacC.setAll(0.0);
}
// Handle convection, axial dispersion (WENO)
const ParamType h = static_cast<ParamType>(_colLength) / static_cast<double>(_nCol);
for (unsigned int i = 0; i < _nRad; ++i)
{