CppADCodeGen 2.4.3
A C++ Algorithmic Differentiation Package with Source Code Generation
Loading...
Searching...
No Matches
dummy_deriv.hpp
1#ifndef CPPAD_CG_DUMMY_DERIV_INCLUDED
2#define CPPAD_CG_DUMMY_DERIV_INCLUDED
3/* --------------------------------------------------------------------------
4 * CppADCodeGen: C++ Algorithmic Differentiation with Source Code Generation:
5 * Copyright (C) 2012 Ciengis
6 * Copyright (C) 2018 Joao Leal
7 *
8 * CppADCodeGen is distributed under multiple licenses:
9 *
10 * - Eclipse Public License Version 1.0 (EPL1), and
11 * - GNU General Public License Version 3 (GPL3).
12 *
13 * EPL1 terms and conditions can be found in the file "epl-v10.txt", while
14 * terms and conditions for the GPL3 can be found in the file "gpl3.txt".
15 * ----------------------------------------------------------------------------
16 * Author: Joao Leal
17 */
18#include <Eigen/Dense>
19#include <Eigen/Sparse>
20#include <Eigen/Eigenvalues>
21#include <Eigen/LU>
22#include <Eigen/QR>
23
24#include <cppad/cg/dae_index_reduction/pantelides.hpp>
25#include <cppad/cg/dae_index_reduction/dummy_deriv_util.hpp>
26
27namespace CppAD {
28namespace cg {
29
33template<class Base>
34class DummyDerivatives : public DaeIndexReduction<Base> {
35 using CGBase = CG<Base>;
36 using ADCG = AD<CGBase>;
37 using VectorB = Eigen::Matrix<Base, Eigen::Dynamic, 1>;
38 using VectorCB = Eigen::Matrix<std::complex<Base>, Eigen::Dynamic, 1>;
39 using MatrixB = Eigen::Matrix<Base, Eigen::Dynamic, Eigen::Dynamic>;
40protected:
48 std::vector<Base> x_;
52 std::vector<Base> normVar_;
56 std::vector<Base> normEq_;
60 std::unique_ptr<ADFun<CGBase> > reducedFun_;
65 std::vector<bool> jacSparsity_;
66 // the initial index of time derivatives
67 size_t diffVarStart_;
68 // the initial index of the differentiated equations
69 size_t diffEqStart_;
75 Eigen::SparseMatrix<Base, Eigen::RowMajor> jacobian_;
79 std::vector<Vnode<Base>*> dummyD_;
96 bool avoidConvertAlg2DifVars_;
100 std::set<std::string> avoidAsDummy_;
101public:
102
115 const std::vector<Base>& x,
116 const std::vector<Base>& normVar,
117 const std::vector<Base>& normEq) :
118 DaeIndexReduction<Base>(idxIdentify.getOriginalModel()),
119 idxIdentify_(&idxIdentify),
120 x_(x),
121 normVar_(normVar),
122 normEq_(normEq),
123 diffVarStart_(0),
124 diffEqStart_(idxIdentify.getOriginalModel().Range()),
125 reduceEquations_(true),
127 reorder_(true),
128 avoidConvertAlg2DifVars_(true) {
129
130 for (Vnode<Base>* jj : idxIdentify.getGraph().variables()) {
131 if (jj->antiDerivative() != nullptr) {
132 diffVarStart_ = jj->index();
133 break;
134 }
135 }
136 }
137
138 inline virtual ~DummyDerivatives() {
139 }
140
148 inline bool isAvoidConvertAlg2DifVars() const {
149 return avoidConvertAlg2DifVars_;
150 }
151
159 inline void setAvoidConvertAlg2DifVars(bool avoid) {
160 avoidConvertAlg2DifVars_ = avoid;
161 }
162
167 inline bool isGenerateSemiExplicitDae() const {
169 }
170
178 inline void setGenerateSemiExplicitDae(bool generateSemiExplicitDae) {
179 generateSemiExplicitDae_ = generateSemiExplicitDae;
180 }
181
186 inline bool isReduceEquations() const {
187 return reduceEquations_;
188 }
189
196 }
197
201 inline bool isReorder() const {
202 return reorder_;
203 }
204
212 inline void setReorder(bool reorder) {
213 reorder_ = reorder;
214 }
215
225 inline void setAvoidVarsAsDummies(const std::set<std::string>& avoidAsDummy) {
226 avoidAsDummy_ = avoidAsDummy;
227 }
228
235 inline const std::set<std::string>& getAvoidVarsAsDummies() const {
236 return avoidAsDummy_;
237 }
238
239 inline std::unique_ptr<ADFun<CG<Base>>> reduceIndex(std::vector<DaeVarInfo>& newVarInfo,
240 std::vector<DaeEquationInfo>& newEqInfo) override {
241
245 std::vector<DaeVarInfo> reducedVarInfo;
249 std::vector<DaeEquationInfo> reducedEqInfo;
250
251 reducedFun_ = idxIdentify_->reduceIndex(reducedVarInfo, reducedEqInfo);
252 if (reducedFun_.get() == nullptr)
253 return nullptr; //nothing to do (no index reduction required)
254
255 if (this->verbosity_ >= Verbosity::Low)
256 log() << "######## Dummy derivatives method ########" << std::endl;
257
258 newEqInfo = reducedEqInfo; // copy
259 addDummyDerivatives(reducedVarInfo, reducedEqInfo, newVarInfo);
260
262
263 matchVars2Eqs4Elimination(newVarInfo, newEqInfo);
264
265 if (reduceEquations_) {
266 std::vector<DaeVarInfo> varInfo = newVarInfo; // copy
267 std::vector<DaeEquationInfo> eqInfo = newEqInfo; // copy
268 std::unique_ptr<ADFun<CG<Base> > > funShort = reduceEquations(varInfo, newVarInfo,
269 eqInfo, newEqInfo);
270 reducedFun_.swap(funShort);
271 }
272
274 std::vector<DaeVarInfo> varInfo = newVarInfo; // copy
275 std::vector<DaeEquationInfo> eqInfo = newEqInfo; // copy
276 std::unique_ptr<ADFun<CG<Base> > > semiExplicit = generateSemiExplicitDAE(*reducedFun_,
277 varInfo, newVarInfo,
278 eqInfo, newEqInfo);
279 reducedFun_.swap(semiExplicit);
280 }
281 }
282
283 if (reorder_) {
284 std::vector<DaeVarInfo> varInfo = newVarInfo; // copy
285 std::vector<DaeEquationInfo> eqInfo = newEqInfo; // copy
286 std::unique_ptr<ADFun<CG<Base>>> reorderedFun = reorderModelEqNVars(*reducedFun_,
287 varInfo, newVarInfo,
288 eqInfo, newEqInfo);
289 reducedFun_.swap(reorderedFun);
290 }
291
292 return std::unique_ptr<ADFun<CG<Base>>>(reducedFun_.release());
293 }
294
295protected:
296
297 using DaeIndexReduction<Base>::log;
298
299 virtual inline void addDummyDerivatives(const std::vector<DaeVarInfo>& varInfo,
300 const std::vector<DaeEquationInfo>& eqInfo,
301 std::vector<DaeVarInfo>& newVarInfo) {
302 auto& graph = idxIdentify_->getGraph();
303 auto& vnodes = graph.variables();
304 auto& enodes = graph.equations();
305
307
308 // variables of interest
309 std::vector<Vnode<Base>*> vars;
310 vars.reserve(vnodes.size() - diffVarStart_);
311 typename std::vector<Vnode<Base>*>::const_reverse_iterator rj;
312 for (rj = vnodes.rbegin(); rj != vnodes.rend(); ++rj) {
313 Vnode<Base>* jj = *rj;
314 if (jj->antiDerivative() != nullptr && jj->derivative() == nullptr) {
315 vars.push_back(jj); // highest order time derivatives in the index 1 model
316 }
317 }
318
319 // should be already fairly sorted, but sort anyway
320 std::sort(vars.begin(), vars.end(), sortVnodesByOrder<Base>);
321
322 // equations of interest
323 typename std::vector<Enode<Base>*>::const_reverse_iterator ri;
324 std::vector<Enode<Base>*> eqs;
325 eqs.reserve(enodes.size() - diffEqStart_);
326 for (ri = enodes.rbegin(); ri != enodes.rend(); ++ri) {
327 Enode<Base>* ii = *ri;
328 if (ii->derivativeOf() != nullptr && ii->derivative() == nullptr) {
329 eqs.push_back(ii);
330 }
331 }
332
333
334 MatrixB workJac;
335
336 while (true) {
337
338 if (this->verbosity_ >= Verbosity::High) {
339 log() << "# equation selection: ";
340 for (size_t i = 0; i < eqs.size(); i++)
341 log() << *eqs[i] << "; ";
342 log() << "\n";
343
344 log() << "# variable selection: ";
345 for (size_t j = 0; j < vars.size(); j++)
346 log() << *vars[j] << "; ";
347 log() << "\n";
348 }
349
350 // Exploit the current equations for elimination of candidates
351 selectDummyDerivatives(eqs, vars, workJac);
352
359 std::vector<Enode<Base>*> newEqs;
360 newEqs.reserve(eqs.size());
361
362 for (Enode<Base>* i : eqs) {
363 Enode<Base>* ii = i->derivativeOf();
364 if (ii != nullptr && ii->derivativeOf() != nullptr) {
365 newEqs.push_back(ii);
366 }
367 }
368 eqs.swap(newEqs);
369
370 if (eqs.empty()) {
371 break;
372 }
373
380 std::vector<Vnode<Base>*> varsNew;
381 varsNew.reserve(vars.size());
382 for (Vnode<Base>* j : vars) {
383 Vnode<Base>* v = j->antiDerivative();
384 if (v != nullptr && v->antiDerivative() != nullptr) {
385 varsNew.push_back(v);
386 }
387 }
388 vars.swap(varsNew);
389 }
390
391
395 newVarInfo = varInfo; //copy
396 for (Vnode<Base>* j : dummyD_) {
397 CPPADCG_ASSERT_UNKNOWN(j->tapeIndex() >= 0);
398 CPPADCG_ASSERT_UNKNOWN(j->antiDerivative() != nullptr);
399 CPPADCG_ASSERT_UNKNOWN(j->antiDerivative()->tapeIndex() >= 0);
400
401 newVarInfo[j->tapeIndex()].setAntiDerivative(-1);
402 newVarInfo[j->antiDerivative()->tapeIndex()].setDerivative(-1);
403 }
404
405 if (this->verbosity_ >= Verbosity::Low) {
406 log() << "## dummy derivatives:\n";
407
408 for (Vnode<Base>* j : dummyD_)
409 log() << "# " << *j << " \t" << newVarInfo[j->tapeIndex()].getName() << "\n";
410 log() << "# \n";
411 if (this->verbosity_ >= Verbosity::High) {
412 graph.printModel(log(), *reducedFun_, newVarInfo, eqInfo);
413 }
414 }
415
416 }
417
425 inline std::unique_ptr<ADFun<CGBase> > reduceEquations(const std::vector<DaeVarInfo>& reducedVarInfo,
426 std::vector<DaeVarInfo>& newVarInfo,
427 const std::vector<DaeEquationInfo>& reducedEqInfo,
428 std::vector<DaeEquationInfo>& newEqInfo) {
429 using namespace std;
430 using std::vector;
431 using std::map;
432 using std::map;
433
434 auto& graph = idxIdentify_->getGraph();
435 //auto& vnodes = graph.variables();
436 auto& enodes = graph.equations();
437
438 CPPADCG_ASSERT_UNKNOWN(reducedVarInfo.size() == reducedFun_->Domain());
439 CPPADCG_ASSERT_UNKNOWN(reducedEqInfo.size() == reducedFun_->Range());
440
441 newEqInfo = reducedEqInfo; // copy
442
446 CodeHandler<Base> handler;
447
448 vector<CGBase> indep0(reducedFun_->Domain());
449 handler.makeVariables(indep0);
450
451 vector<CGBase> res0 = graph.forward0(*reducedFun_, indep0);
452
453 map<int, int> assignedVar2Eq;
454 for (size_t i = 0; i < newEqInfo.size(); ++i) {
455 DaeEquationInfo& newEq = newEqInfo[i];
456 if (newEq.getAssignedVarIndex() >= 0)
457 assignedVar2Eq[newEq.getAssignedVarIndex()] = i;
458 }
459
465 vector<int> eqIndexReduced2Short(enodes.size());
466 for (size_t i = 0; i < eqIndexReduced2Short.size(); i++) {
467 eqIndexReduced2Short[i] = i;
468 }
469
475 vector<int> tapeIndexReduced2Short(reducedVarInfo.size());
476 for (size_t j = 0; j < tapeIndexReduced2Short.size(); j++) {
477 tapeIndexReduced2Short[j] = j;
478 }
479
484 set<size_t> erasedVariables;
485 set<size_t> erasedEquations;
486
487 if (this->verbosity_ >= Verbosity::High)
488 log() << "Reducing total number of equations by symbolic manipulation:" << std::endl;
489
490 for (Vnode<Base>* dummy : dummyD_) {
491
495 map<int, int>::const_iterator ita = assignedVar2Eq.find(dummy->tapeIndex());
496 if (ita == assignedVar2Eq.end()) {
497 if (this->verbosity_ >= Verbosity::High)
498 log() << "unable to solve for variable " << dummy->name() << "." << std::endl;
499
500 continue; // unable to solve for a dummy variable: keep the equation and variable
501 }
502 int bestEquation = ita->second;
503
504 try {
505 // eliminate all references to the dummy variable by substitution
506 handler.substituteIndependent(indep0[dummy->tapeIndex()], res0[bestEquation]);
507 tapeIndexReduced2Short[dummy->tapeIndex()] = -1;
508 eqIndexReduced2Short[bestEquation] = -1;
509
510 if (this->verbosity_ >= Verbosity::High) {
511 log() << "######### use equation " << *enodes[newEqInfo[bestEquation].getId()] << " to solve for variable " << dummy->name() << std::endl;
512 erasedVariables.insert(dummy->tapeIndex());
513 erasedEquations.insert(bestEquation);
514 printModel(log(), handler, res0, reducedVarInfo, erasedVariables, erasedEquations);
515 }
516
517 } catch (const CGException& ex) {
518 // unable to solve for a dummy variable: keep the equation and variable
519 if (this->verbosity_ >= Verbosity::High)
520 log() << "unable to use equation " << *enodes[newEqInfo[bestEquation].getId()] << " to solve for variable " << dummy->name() << ": " << ex.what() << std::endl;
521 }
522 }
523
524 // determine the new equation indexes
525 for (size_t i = 0; i < eqIndexReduced2Short.size(); i++) {
526 if (eqIndexReduced2Short[i] < 0) { // removed equation
527 for (size_t ii = i + 1; ii < eqIndexReduced2Short.size(); ii++) {
528 if (eqIndexReduced2Short[ii] >= 0) {
529 eqIndexReduced2Short[ii]--;
530 }
531 }
532 }
533 }
534
535 // determine the new indexes in the tape
536 for (size_t p = 0; p < tapeIndexReduced2Short.size(); p++) {
537 if (tapeIndexReduced2Short[p] < 0) {
538 // removed from model
539 for (size_t p2 = p + 1; p2 < tapeIndexReduced2Short.size(); p2++) {
540 if (tapeIndexReduced2Short[p2] >= 0) {
541 tapeIndexReduced2Short[p2]--;
542 }
543 }
544 }
545 }
546
550 CPPADCG_ASSERT_UNKNOWN(tapeIndexReduced2Short.size() == reducedVarInfo.size());
551
552 newVarInfo = reducedVarInfo; // copy
553 for (int p = tapeIndexReduced2Short.size() - 1; p >= 0; p--) {
554 if (tapeIndexReduced2Short[p] < 0) { // removed from model
555 newVarInfo.erase(newVarInfo.begin() + p);
556 for (size_t pp = 0; pp < tapeIndexReduced2Short.size(); pp++) {
557 DaeVarInfo& v = newVarInfo[pp];
558 if (v.getAntiDerivative() > p) {
559 v.setAntiDerivative(v.getAntiDerivative() - 1);
560 } else if (v.getAntiDerivative() == p) {
561 v.setAntiDerivative(-1);
562 }
563 if (v.getDerivative() > p) {
564 v.setDerivative(v.getDerivative() - 1);
565 } else if (v.getDerivative() == p) {
566 v.setDerivative(-1);
567 }
568 }
569 }
570 }
571
572 for (int p = eqIndexReduced2Short.size() - 1; p >= 0; p--) {
573 if (eqIndexReduced2Short[p] < 0) {// removed from model
574 newEqInfo.erase(newEqInfo.begin() + p);
575 } else {
576 DaeEquationInfo& eq = newEqInfo[p];
577 int reducedVIndex = eq.getAssignedVarIndex();
578 if (reducedVIndex >= 0)
579 eq.setAssignedVarIndex(tapeIndexReduced2Short[reducedVIndex]);
580 if (eq.getAntiDerivative() >= 0)
581 eq.setAntiDerivative(eqIndexReduced2Short[eq.getAntiDerivative()]);
582 }
583 }
584
589 std::unique_ptr<ADFun<CGBase> > shortFun(generateReorderedModel(handler, res0,
590 reducedVarInfo, newVarInfo,
591 reducedEqInfo, newEqInfo));
592
593 if (this->verbosity_ >= Verbosity::High) {
594 log() << "DAE with less equations and variables:\n";
595 graph.printModel(log(), *shortFun, newVarInfo, newEqInfo);
596 }
597
598 return shortFun;
599 }
600
611 inline std::unique_ptr<ADFun<CGBase> > generateSemiExplicitDAE(ADFun<CG<Base> >& fun,
612 const std::vector<DaeVarInfo>& varInfo,
613 std::vector<DaeVarInfo>& newVarInfo,
614 const std::vector<DaeEquationInfo>& eqInfo,
615 std::vector<DaeEquationInfo>& newEqInfo) {
616 using namespace std;
617 using std::vector;
618 using std::map;
619
620 auto& graph = idxIdentify_->getGraph();
621
622 newEqInfo = eqInfo; // copy (we will have the same number of equations)
623
627 CodeHandler<Base> handler;
628
629 vector<CGBase> indep0(fun.Domain());
630 handler.makeVariables(indep0);
631
632 vector<CGBase> res0 = graph.forward0(fun, indep0);
633
634 map<int, int> assignedVar2Eq;
635 for (size_t i = 0; i < newEqInfo.size(); ++i) {
636 DaeEquationInfo& newEq = newEqInfo[i];
637 assignedVar2Eq[newEq.getAssignedVarIndex()] = i;
638 }
639
643 for (size_t j = 0; j < varInfo.size(); ++j) {
644 const DaeVarInfo& jj = varInfo[j];
645 if (jj.getAntiDerivative() < 0) {
646 continue; // not a time derivative
647 }
648 CGBase& indep = indep0[j]; // the time derivative
652 map<int, int>::const_iterator ita = assignedVar2Eq.find(j);
653 if (ita == assignedVar2Eq.end()) {
654 throw CGException("Failed to generate semi-explicit DAE: unable to create an explicit equation for ", jj.getName());
655 }
656
657 int bestEquation = ita->second;
658
659 try {
660 CGBase& dep = res0[bestEquation]; // the equation residual
661
662 handler.substituteIndependent(indep, dep); // removes indep from the list of variables
663
664 OperationNode<Base>* alias = indep.getOperationNode();
665 CPPADCG_ASSERT_UNKNOWN(alias != nullptr && alias->getOperationType() == CGOpCode::Alias);
666 dep.getOperationNode()->makeAlias(alias->getArguments()[0]); // not a residual anymore but equal to alias (explicit equation)
667
668 // it is now an explicit differential equation
669 newEqInfo[bestEquation].setExplicit(true);
670 // the derivative variable will disappear, associate the equation with the original variable
671 newEqInfo[bestEquation].setAssignedVarIndex(jj.getAntiDerivative());
672 } catch (const CGException& ex) {
673 // unable to solve for a dummy variable: keep the equation and variable
674 throw CGException("Failed to generate semi-explicit DAE: ", ex.what());
675 }
676 }
677
682 vector<int> varIndexOld2New(varInfo.size(), -1);
683 size_t count = 0;
684 for (size_t j = 0; j != varInfo.size(); ++j) {
685 // exclude derivatives (they will be removed)
686 if (varInfo[j].getAntiDerivative() < 0) {
687 varIndexOld2New[j] = count++;
688 }
689 }
690
691 for (size_t i = 0; i < newEqInfo.size(); ++i) {
692 const DaeEquationInfo& ii = newEqInfo[i];
693 int j = ii.getAssignedVarIndex();
694 if (j >= 0)
695 newEqInfo[i].setAssignedVarIndex(varIndexOld2New[j]);
696 }
697
701 newVarInfo = varInfo;
702 for (int j = newVarInfo.size() - 1; j >= 0; --j) {
703 if (newVarInfo[j].getAntiDerivative() >= 0) {
704 // a derivative
705 newVarInfo.erase(newVarInfo.begin() + j);
706 }
707 }
708 for (size_t j = 0; j < newVarInfo.size(); j++) {
709 newVarInfo[j].setDerivative(-1); // no derivatives in tape
710 }
711
716 std::unique_ptr<ADFun<CGBase> > semiExplicitFun(generateReorderedModel(handler, res0, varInfo, newVarInfo, eqInfo, newEqInfo));
717
718 if (this->verbosity_ >= Verbosity::High) {
719 log() << "Semi-Eplicit DAE:\n";
720 graph.printModel(log(), *semiExplicitFun, newVarInfo, newEqInfo);
721 }
722
723 return semiExplicitFun;
724 }
725
726 inline void matchVars2Eqs4Elimination(std::vector<DaeVarInfo>& varInfo,
727 std::vector<DaeEquationInfo>& eqInfo) {
728 using std::vector;
729 using std::map;
730
731 auto& graph = idxIdentify_->getGraph();
732 auto& vnodes = graph.variables();
733 auto& enodes = graph.equations();
734
735 CPPADCG_ASSERT_UNKNOWN(eqInfo.size() == enodes.size());
736 CPPADCG_ASSERT_UNKNOWN(varInfo.size() == reducedFun_->Domain());
737 CPPADCG_ASSERT_UNKNOWN(eqInfo.size() == reducedFun_->Range());
738
739 CodeHandler<Base> handler;
740
741 vector<CGBase> indep0(reducedFun_->Domain());
742 handler.makeVariables(indep0);
743
744 vector<CGBase> res0 = graph.forward0(*reducedFun_, indep0);
745
746 vector<bool> jacSparsity = jacobianSparsity<vector<bool> >(*reducedFun_);
747
748 vector<Vnode<Base>*> diffVariables;
749 vector<Vnode<Base>*> dummyVariables;
750 vector<Vnode<Base>*> variables;
751 vector<Enode<Base>*> equations(enodes.size(), nullptr);
752 try {
756 // create variable nodes
757 map<Vnode<Base>*, Vnode<Base>*> eliminateOrig2New;
758 for (size_t j = 0; j < vnodes.size(); j++) {
759 Vnode<Base>* v = vnodes[j];
760 if (std::find(dummyD_.begin(), dummyD_.end(), v) != dummyD_.end()) {
761 if (reduceEquations_) {
762 dummyVariables.push_back(new Vnode<Base>(j, v->tapeIndex(), v->name()));
763 eliminateOrig2New[v] = dummyVariables.back();
764 }
765 } else if (v->antiDerivative() != nullptr) {
767 diffVariables.push_back(new Vnode<Base>(j, v->tapeIndex(), v->name()));
768 eliminateOrig2New[v] = diffVariables.back();
769 }
770 }
771 }
772
773 variables.reserve(diffVariables.size() + dummyVariables.size());
774 variables.insert(variables.end(), diffVariables.begin(), diffVariables.end()); // must be added 1st (priority)
775 variables.insert(variables.end(), dummyVariables.begin(), dummyVariables.end());
776
777 // create equation nodes
778 for (size_t i = 0; i < enodes.size(); i++) {
779 equations[i] = new Enode<Base>(i, enodes[i]->name());
780 const vector<Vnode<Base>*>& origVars = enodes[i]->originalVariables();
781 for (size_t p = 0; p < origVars.size(); p++) {
782 Vnode<Base>* jOrig = origVars[p];
783
784 typename map<Vnode<Base>*, Vnode<Base>*>::const_iterator it;
785 it = eliminateOrig2New.find(jOrig);
786 if (it != eliminateOrig2New.end() &&
787 jacSparsity[varInfo.size() * i + jOrig->tapeIndex()]) {
788 Vnode<Base>* j = it->second;
789
790 CGBase& dep = res0[i]; // the equation residual
791 CGBase& indep = indep0[j->tapeIndex()];
792
793 if (handler.isSolvable(*dep.getOperationNode(), *indep.getOperationNode())) {
794 equations[i]->addVariable(j);
795 }
796 }
797 }
798 }
799
800 map<size_t, Vnode<Base>*> tape2FreeVariables;
801 for (Vnode<Base>* j : variables) {
802 tape2FreeVariables[j->tapeIndex()] = j;
803 }
804
805
809 while (true) {
810 size_t assigned;
811 do {
812 do {
816 do {
817 assigned = 0;
818 for (Vnode<Base>* j : diffVariables) {
819 if (!j->isDeleted() && j->equations().size() == 1) {
820 Enode<Base>& i = *j->equations()[0];
821 if (i.assignmentVariable() == nullptr) {
822 if (!assignVar2Equation(i, res0, *j, indep0, handler,
823 jacSparsity, tape2FreeVariables,
824 equations, varInfo)) {
825 throw CGException("Failed to solve equation ", i.name(), " for variable ", j->name());
826 }
827 assigned++;
828 }
829 }
830 }
831 } while (assigned > 0);
832
837 assigned = 0;
838 for (Vnode<Base>* j : dummyVariables) {
839 if (!j->isDeleted() && j->equations().size() == 1) {
840 Enode<Base>& i = *j->equations()[0];
841 if (i.assignmentVariable() == nullptr) {
842 if (assignVar2Equation(i, res0, *j, indep0, handler,
843 jacSparsity, tape2FreeVariables,
844 equations, varInfo))
845 assigned++;
846 }
847 }
848 }
849 } while (assigned > 0);
850
855 assigned = 0;
856 for (Enode<Base>* i : equations) {
857 if (i->assignmentVariable() == nullptr && i->variables().size() == 1) {
858 Vnode<Base>* j = i->variables()[0];
859 if (assignVar2Equation(*i, res0, *j, indep0, handler,
860 jacSparsity, tape2FreeVariables,
861 equations, varInfo))
862 assigned++;
863 }
864 }
865
866 } while (assigned > 0);
867
873 assigned = 0;
874 for (Vnode<Base>* j : variables) {
875 if (!j->isDeleted()) {
876 for (Enode<Base>* i : j->equations()) {
877 if (i->assignmentVariable() == nullptr) {
878 if (assignVar2Equation(*i, res0, *j, indep0, handler,
879 jacSparsity, tape2FreeVariables,
880 equations, varInfo)) {
881 assigned++;
882 break;
883 }
884 }
885 }
886 if (assigned > 0)
887 break;
888 }
889 }
890
891 if (assigned == 0) {
892 break; // done
893 }
894 }
895
896
901 for (Vnode<Base>* j : variables) { // previous assignments must not change!
902 if (j->assignmentEquation() != nullptr) {
903 j->deleteNode();
904 }
905 }
906
908 for(Enode<Base>* i: equations) {
909 if (i->assignmentVariable() == nullptr) {
910 augment.augmentPath(*i);
911 }
912 }
913
917 for (Vnode<Base>* j : variables) {
918 if (j->assignmentEquation() != nullptr) {
919 int i = j->assignmentEquation()->index();
920 DaeEquationInfo& eq = eqInfo[i];
921
922 if (eq.getAssignedVarIndex() != int(j->tapeIndex())) {
923 eq.setAssignedVarIndex(j->tapeIndex());
924 }
925 }
926 }
927
928 // verify results
930 std::string error;
931 for (Vnode<Base>* j : diffVariables) {
932 if (j->assignmentEquation() == nullptr) {
933 // failed!!!
934 if (!error.empty())
935 error += ",";
936 error += " " + j->name();
937 }
938 }
939 if (!error.empty())
940 throw CGException("Failed to generate semi-explicit DAE. Could not solve system for the following variables:", error);
941 }
942
943 } catch (...) {
944 deleteVectorValues(diffVariables);
945 deleteVectorValues(dummyVariables);
946 deleteVectorValues(equations);
947 throw;
948 }
949
950 if (this->verbosity_ >= Verbosity::High) {
951 for (Vnode<Base>* j : variables) {
952 if (j->assignmentEquation() != nullptr)
953 log() << "## Variable " + j->name() << " assigned to equation " << j->assignmentEquation()->name() << "\n";
954 }
955 log() << std::endl;
956 }
957 deleteVectorValues(diffVariables);
958 deleteVectorValues(dummyVariables);
959 deleteVectorValues(equations);
960 }
961
962 inline bool assignVar2Equation(Enode<Base>& i, std::vector<CGBase>& res0,
963 Vnode<Base>& j, std::vector<CGBase>& indep0,
964 CodeHandler<Base>& handler,
965 std::vector<bool>& jacSparsity,
966 const std::map<size_t, Vnode<Base>*>& tape2FreeVariables,
967 std::vector<Enode<Base>*>& equations,
968 std::vector<DaeVarInfo>& varInfo) {
969 using namespace std;
970 using std::vector;
971 using std::map;
972
973 std::vector<bool> localJacSparsity = jacSparsity;
974 const size_t n = varInfo.size();
975
979 CGBase& dep = res0[i.index()];
980 CGBase& indep = indep0[j.tapeIndex()];
981
982 std::string indepName;
983 if (indep.getOperationNode()->getName() != nullptr) {
984 indepName = *indep.getOperationNode()->getName();
985 }
986
987
988 try {
989 handler.substituteIndependent(indep, dep, false); // indep not removed from the list of variables yet
990
991 OperationNode<Base>* alias = indep.getOperationNode();
992 CPPADCG_ASSERT_UNKNOWN(alias != nullptr && alias->getOperationType() == CGOpCode::Alias);
993
994 // it is now an explicit differential equation
995 //newEqInfo[bestEquation].setExplicit(true);
996 // the derivative variable will disappear, associate the equation with the original variable
997 //newEqInfo[bestEquation].setAssignedVarIndex(jj.getAntiDerivative());
998 } catch (const CGException& ex) {
999 // unable to solve for a dummy variable: keep the equation and variable
1000 throw CGException("Failed to solve equation ", i.name(), " for variable ", j.name(), ": ", ex.what());
1001 }
1002
1003
1008 vector<size_t> nnzs;
1009 for (const auto& it : tape2FreeVariables) {
1010 size_t tapeJ = it.first;
1011 if (localJacSparsity[n * i.index() + tapeJ] && tapeJ != j.tapeIndex()) {
1012 nnzs.push_back(tapeJ);
1013 }
1014 }
1015 set<Enode<Base>*> affected;
1016 for (size_t e = 0; e < equations.size(); ++e) {
1017 if (equations[e] != &i && localJacSparsity[n * e + j.tapeIndex()]) {
1018 localJacSparsity[n * e + j.tapeIndex()] = false; // eliminated by substitution
1019 affected.insert(equations[e]);
1020 for (size_t p = 0; p < nnzs.size(); ++p) {
1021 localJacSparsity[n * e + nnzs[p]] = true;
1022 }
1023 }
1024 }
1025
1026 if (affected.size() > 0) {
1030 map<size_t, set<Enode<Base>*> > solvable;
1031 for (size_t e = 0; e < equations.size(); ++e) {
1032 Enode<Base>* eq = equations[e];
1033 if (&i != eq && affected.find(eq) == affected.end()) {
1034 // no change
1035 for (size_t v = 0; v < eq->variables().size(); ++v) {
1036 solvable[eq->variables()[v]->tapeIndex()].insert(eq);
1037 }
1038 }
1039 }
1040
1041 // redetermine solvability
1042 for (Enode<Base>* itAff : affected) {
1043 Enode<Base>& a = *itAff;
1044 for (const auto& it : tape2FreeVariables) {
1045 size_t jj = it.first;
1046 if (localJacSparsity[n * a.index() + jj]) {
1047 if (handler.isSolvable(*res0[a.index()].getOperationNode(), *indep0[jj].getOperationNode())) {
1048 solvable[jj].insert(&a);
1049 }
1050 }
1051 }
1052 }
1053
1054 // check if any variable stops being solvable
1055 bool ok = true;
1056 for (const auto& it : tape2FreeVariables) {
1057 Vnode<Base>* v = it.second;
1058 if (v == &j)
1059 continue;
1060
1061 if (v->assignmentEquation() != nullptr) {
1062 if (affected.count(v->assignmentEquation()) > 0 &&
1063 solvable[v->tapeIndex()].count(v->assignmentEquation()) == 0) {
1064 ok = false;
1065 break;
1066 }
1067 } else if (solvable[v->tapeIndex()].size() == 0) {
1068 ok = false;
1069 break;
1070 }
1071 }
1072
1073 if (!ok) {
1074 handler.undoSubstituteIndependent(*indep.getOperationNode());
1075 if (indepName.size() > 0) {
1076 indep.getOperationNode()->setName(indepName);
1077 }
1078 return false;
1079 }
1080
1084 for (Enode<Base>* itAff : affected) {
1085 Enode<Base>& a = *itAff;
1086 for (const auto& it : tape2FreeVariables) {
1087 size_t v = it.first;
1088 if (localJacSparsity[n * a.index() + v]) {
1089 if (solvable[v].count(&a) > 0) {
1090 a.addVariable(it.second);
1091 } else {
1092 // not solvable anymore
1093 a.deleteNode(it.second);
1094 }
1095 }
1096 }
1097 }
1098
1099 }
1100
1104 handler.removeIndependent(*indep.getOperationNode());
1105
1109 j.setAssignmentEquation(i, log(), this->verbosity_);
1110 j.deleteNode(log(), this->verbosity_);
1111
1112 jacSparsity = localJacSparsity;
1113
1114 return true;
1115 }
1116
1117 inline std::unique_ptr<ADFun<CGBase> > reorderModelEqNVars(ADFun<CG<Base> >& fun,
1118 const std::vector<DaeVarInfo>& varInfo,
1119 std::vector<DaeVarInfo>& newVarInfo,
1120 const std::vector<DaeEquationInfo>& eqInfo,
1121 std::vector<DaeEquationInfo>& newEqInfo) {
1122
1123 using namespace std;
1124 using std::vector;
1125
1126 auto& graph = idxIdentify_->getGraph();
1127
1131 std::set<size_t> oldVarWithDerivatives; // indexes of old variables (before reordering) with derivatives
1132 for (size_t i = 0; i < eqInfo.size(); i++) {
1133 if (eqInfo[i].isExplicit() && eqInfo[i].getAssignedVarIndex() >= 0) {
1134 oldVarWithDerivatives.insert(eqInfo[i].getAssignedVarIndex());
1135 }
1136 }
1137
1138 if (oldVarWithDerivatives.empty()) {
1139 // no semi-explicit model generated
1140 for (size_t j = 0; j < varInfo.size(); j++) {
1141 int index = j;
1142 bool differential = false;
1143 while (varInfo[index].getAntiDerivative() >= 0) {
1144 index = varInfo[index].getAntiDerivative();
1145 differential = true;
1146 }
1147
1148 if (differential) {
1149 oldVarWithDerivatives.insert(index);
1150 }
1151 }
1152 }
1153
1157 std::vector<DaeVarOrderInfo> varOrder(varInfo.size());
1158 for (size_t j = 0; j < varInfo.size(); j++) {
1159 size_t j0;
1160 int derivOrder = graph.determineVariableDiffOrder(varInfo, j, j0);
1161 if (varInfo[j].isIntegratedVariable()) {
1162 derivOrder = -2; // so that it goes last
1163 }
1164 bool hasDerivatives = oldVarWithDerivatives.find(j) != oldVarWithDerivatives.end();
1165 varOrder[j] = DaeVarOrderInfo(j, j0, hasDerivatives, derivOrder);
1166 }
1167
1168 std::sort(varOrder.begin(), varOrder.end(), sortVariablesByOrder);
1169
1173 std::vector<size_t> varIndexOld2New(varInfo.size(), -1);
1174 for (size_t j = 0; j < varOrder.size(); ++j) {
1175 varIndexOld2New[varOrder[j].originalIndex] = j;
1176 }
1177
1178 newVarInfo.resize(varInfo.size());
1179 for (size_t j = 0; j < varOrder.size(); ++j) {
1180 newVarInfo[j] = varInfo[varOrder[j].originalIndex];
1181 int oldDerivOfIndex = newVarInfo[j].getAntiDerivative();
1182 if (oldDerivOfIndex >= 0)
1183 newVarInfo[j].setAntiDerivative(varIndexOld2New[oldDerivOfIndex]);
1184 int oldDerivIndex = newVarInfo[j].getDerivative();
1185 if (oldDerivIndex >= 0)
1186 newVarInfo[j].setDerivative(varIndexOld2New[oldDerivIndex]);
1187 }
1188
1192 newEqInfo = eqInfo; //copy
1193 for (size_t i = 0; i < newEqInfo.size(); i++) {
1194 int oldVIndex = newEqInfo[i].getAssignedVarIndex();
1195 if (oldVIndex >= 0) {
1196 newEqInfo[i].setAssignedVarIndex(varIndexOld2New[oldVIndex]);
1197 }
1198 }
1199
1200 std::vector<DaeEqOrderInfo> eqOrder(newEqInfo.size());
1201 for (size_t i = 0; i < newEqInfo.size(); i++) {
1202 int assignedVar = newEqInfo[i].getAssignedVarIndex();
1203 size_t i0 = i;
1204 while (newEqInfo[i0].getAntiDerivative() >= 0) {
1205 i0 = newEqInfo[i0].getAntiDerivative();
1206 }
1207 bool isDifferential = newEqInfo[i].isExplicit() || (assignedVar >= 0 && newVarInfo[assignedVar].getAntiDerivative() >= 0);
1208 eqOrder[i] = DaeEqOrderInfo(i, i0, isDifferential, assignedVar);
1209 }
1210
1211 std::sort(eqOrder.begin(), eqOrder.end(), sortEquationByAssignedOrder2);
1212
1213 std::vector<DaeEquationInfo> newEqInfo2(newEqInfo.size());
1214 for (size_t i = 0; i < eqOrder.size(); i++) {
1215 newEqInfo2[i] = newEqInfo[eqOrder[i].originalIndex];
1216 }
1217 newEqInfo = newEqInfo2;
1218
1219
1223 CodeHandler<Base> handler;
1224
1225 vector<CGBase> indep0(fun.Domain());
1226 handler.makeVariables(indep0);
1227
1228 const vector<CGBase> res0 = graph.forward0(fun, indep0);
1229
1233 std::unique_ptr<ADFun<CGBase> > reorderedFun(generateReorderedModel(handler, res0, varInfo, newVarInfo, eqInfo, newEqInfo));
1234
1235 if (this->verbosity_ >= Verbosity::High) {
1236 log() << "reordered DAE equations and variables:\n";
1237 graph.printModel(log(), *reorderedFun, newVarInfo, newEqInfo);
1238 }
1239
1240 return reorderedFun;
1241 }
1242
1244 const std::vector<CGBase>& res0,
1245 const std::vector<DaeVarInfo>& varInfo,
1246 const std::vector<DaeVarInfo>& newVarInfo,
1247 const std::vector<DaeEquationInfo>& eqInfo,
1248 const std::vector<DaeEquationInfo>& newEqInfo) const {
1249 using std::vector;
1250
1251 vector<ADCG> indepNewOrder(handler.getIndependentVariableSize());
1252 CPPADCG_ASSERT_UNKNOWN(indepNewOrder.size() == newVarInfo.size());
1253
1254 for (size_t p = 0; p < newVarInfo.size(); p++) {
1255 int origIndex = newVarInfo[p].getOriginalIndex();
1256 if (origIndex >= 0) {
1257 indepNewOrder[p] = x_[origIndex];
1258 }
1259 }
1260
1261 Independent(indepNewOrder);
1262
1269 std::set<size_t> newIds;
1270 for (size_t j = 0; j < newVarInfo.size(); j++) {
1271 newIds.insert(newVarInfo[j].getId());
1272 }
1273
1274 std::map<size_t, size_t> varId2HandlerIndex;
1275 size_t handlerIndex = 0; // start the variable count again since some variable might have been removed
1276 for (size_t j = 0; j < varInfo.size(); j++) {
1277 int id = varInfo[j].getId();
1278 if (newIds.find(id) != newIds.end()) {
1279 varId2HandlerIndex[id] = handlerIndex++; // not removed from model
1280 }
1281 }
1282
1283 vector<ADCG> indepHandlerOrder(handler.getIndependentVariableSize());
1284 for (size_t p = 0; p < newVarInfo.size(); p++) {
1285 size_t id = newVarInfo[p].getId();
1286 indepHandlerOrder[varId2HandlerIndex[id]] = indepNewOrder[p];
1287 }
1288
1289 // reorder equations
1290 std::map<size_t, size_t> eqId2OldIndex;
1291 for (size_t i = 0; i < eqInfo.size(); i++) {
1292 eqId2OldIndex[eqInfo[i].getId()] = i;
1293 }
1294
1295 vector<CGBase> resNewOrder(newEqInfo.size());
1296 for (size_t i = 0; i < newEqInfo.size(); i++) {
1297 size_t oldIndex = eqId2OldIndex[newEqInfo[i].getId()];
1298 resNewOrder[i] = res0[oldIndex];
1299 }
1300
1301 // evaluate the model
1302 Evaluator<Base, CGBase> evaluator0(handler);
1303 evaluator0.setPrintFor(idxIdentify_->getGraph().isPreserveNames()); // variable names saved with CppAD::PrintFor
1304 vector<ADCG> depNewOrder = evaluator0.evaluate(indepHandlerOrder, resNewOrder);
1305
1306 return new ADFun<CGBase>(indepNewOrder, depNewOrder);
1307 }
1308
1313 inline void determineJacobian() {
1314 using namespace std;
1315 using std::vector;
1316
1317 const size_t n = reducedFun_->Domain();
1318 const size_t m = reducedFun_->Range();
1319
1320 auto& graph = idxIdentify_->getGraph();
1321 auto& vnodes = graph.variables();
1322 auto& enodes = graph.equations();
1323
1324 jacSparsity_ = jacobianReverseSparsity<vector<bool>, CGBase>(*reducedFun_); // in the original variable order
1325
1326 vector<size_t> row, col;
1327 row.reserve((vnodes.size() - diffVarStart_) * (m - diffEqStart_));
1328 col.reserve(row.capacity());
1329
1330 for (size_t i = diffEqStart_; i < m; i++) {
1331 for (size_t j = diffVarStart_; j < vnodes.size(); j++) {
1332 CPPADCG_ASSERT_UNKNOWN(vnodes[j]->antiDerivative() != nullptr);
1333 size_t t = vnodes[j]->tapeIndex();
1334 if (jacSparsity_[i * n + t]) {
1335 row.push_back(i);
1336 col.push_back(t);
1337 }
1338 }
1339 }
1340
1341 vector<CG<Base> > jac(row.size());
1342
1343 vector<CG<Base> > indep(n);
1344 std::copy(x_.begin(), x_.end(), indep.begin());
1345 std::fill(indep.begin() + x_.size(), indep.end(), 0);
1346
1347 CppAD::sparse_jacobian_work work; // temporary structure for CPPAD
1348 reducedFun_->SparseJacobianReverse(indep, jacSparsity_,
1349 row, col, jac, work);
1350
1351 // resize and zero matrix
1352 jacobian_.resize(m - diffEqStart_, vnodes.size() - diffVarStart_);
1353
1354 map<size_t, Vnode<Base>*> origIndex2var;
1355 for (size_t j = diffVarStart_; j< vnodes.size(); j++) {
1356 Vnode<Base>* jj = vnodes[j];
1357 origIndex2var[jj->tapeIndex()] = jj;
1358 }
1359
1360 // normalize values
1361 for (size_t e = 0; e < jac.size(); e++) {
1362 Enode<Base>* eqOrig = enodes[row[e]]->originalEquation();
1363 Vnode<Base>* vOrig = origIndex2var[col[e]]->originalVariable(graph.getOrigTimeDependentCount());
1364
1365 // normalized jacobian value
1366 Base normVal = jac[e].getValue() * normVar_[vOrig->tapeIndex()]
1367 / normEq_[eqOrig->index()];
1368
1369 size_t i = row[e]; // same order
1370 size_t j = origIndex2var[col[e]]->index(); // different order than in model/tape
1371
1372 jacobian_.coeffRef(i - diffEqStart_, j - diffVarStart_) = normVal;
1373 }
1374
1375 jacobian_.makeCompressed();
1376
1377 if (this->verbosity_ >= Verbosity::High) {
1378 log() << "\npartial jacobian:\n" << jacobian_ << "\n\n";
1379 //cout << jacobian_.triangularView<Eigen::Lower > () << "\n\n";
1380 }
1381 }
1382
1383 inline void selectDummyDerivatives(const std::vector<Enode<Base>* >& eqs,
1384 const std::vector<Vnode<Base>* >& vars,
1385 MatrixB& work) {
1386
1387 if (eqs.size() == vars.size()) {
1388 dummyD_.insert(dummyD_.end(), vars.begin(), vars.end());
1389 if (this->verbosity_ >= Verbosity::High) {
1390 log() << "# new dummy derivatives: ";
1391 for (size_t j = 0; j < vars.size(); j++)
1392 log() << *vars[j] << "; ";
1393 log() << " \n";
1394 }
1395#ifndef NDEBUG
1396 for (Vnode<Base>* it : vars) {
1397 CPPADCG_ASSERT_UNKNOWN(std::find(dummyD_.begin(), dummyD_.end(), it) == dummyD_.end());
1398 }
1399#endif
1400 return;
1401 }
1402
1406 std::set<size_t> excludeCols;
1407 std::set<size_t> avoidCols;
1408 for (size_t j = 0; j < vars.size(); j++) {
1409 Vnode<Base>* jj = vars[j];
1410 bool notZero = false;
1411 for (size_t i = 0; i < eqs.size(); i++) {
1412 Enode<Base>* ii = eqs[i];
1413 Base val = jacobian_.coeff(ii->index() - diffEqStart_, jj->index() - diffVarStart_);
1414 if (val != Base(0.0)) {
1415 notZero = true;
1416 break;
1417 }
1418 }
1419 if (!notZero) {
1420 // all zeros: must not choose this column/variable
1421 excludeCols.insert(j);
1422 } else if (avoidAsDummy_.find(vars[j]->name()) != avoidAsDummy_.end()) {
1423 // avoid using this column/variable if there are other alternatives
1424 avoidCols.insert(j);
1425 }
1426 }
1427
1428 if (eqs.size() <= vars.size() - (excludeCols.size() + avoidCols.size())) {
1429 // there might be sufficient alternatives to avoid using the columns/variables disabled by the user
1430 excludeCols.insert(avoidCols.begin(), avoidCols.end());
1431 } else {
1432 if (this->verbosity_ >= Verbosity::Low)
1433 log() << "Must use variables defined to be avoided by the user!\n";
1434 }
1435
1436 std::vector<Vnode<Base>* > varsLocal;
1437
1438 Eigen::ColPivHouseholderQR<MatrixB> qr;
1439
1440 auto orderColumns = [&]() {
1441 varsLocal.reserve(vars.size() - excludeCols.size());
1442 for (size_t j = 0; j < vars.size(); j++) {
1443 if (excludeCols.find(j) == excludeCols.end()) {
1444 varsLocal.push_back(vars[j]);
1445 }
1446 }
1447
1448
1449 work.setZero(eqs.size(), varsLocal.size());
1450
1451 // determine the rows that only contain a single nonzero (a single column)
1452 for (size_t i = 0; i < eqs.size(); i++) {
1453 Enode<Base>* ii = eqs[i];
1454 for (size_t j = 0; j < varsLocal.size(); j++) {
1455 const Vnode<Base>* jj = varsLocal[j];
1456 Base val = jacobian_.coeff(ii->index() - diffEqStart_, jj->index() - diffVarStart_);
1457 if (val != Base(0.0)) {
1458 work(i, j) = val;
1459 }
1460 }
1461 }
1462
1463 if (this->verbosity_ >= Verbosity::High)
1464 log() << "subset Jac:\n" << work << "\n";
1465
1466 qr.compute(work);
1467
1468 if (qr.info() != Eigen::Success) {
1469 throw CGException("Failed to select dummy derivatives! "
1470 "QR decomposition of a submatrix of the Jacobian failed!");
1471 } else if (qr.rank() < work.rows()) {
1472 throw CGException("Failed to select dummy derivatives! "
1473 "The resulting system is probably singular for the provided data.");
1474 }
1475
1476 using PermutationMatrix = typename Eigen::ColPivHouseholderQR<MatrixB>::PermutationType;
1477 using Indices = typename PermutationMatrix::IndicesType;
1478
1479 const PermutationMatrix& p = qr.colsPermutation();
1480 const Indices& indices = p.indices();
1481
1482 if (this->verbosity_ >= Verbosity::High) {
1483 log() << "## matrix Q:\n";
1484 MatrixB q = qr.matrixQ();
1485 log() << q << "\n";
1486 log() << "## matrix R:\n";
1487 MatrixB r = qr.matrixR().template triangularView<Eigen::Upper>();
1488 log() << r << "\n";
1489 log() << "## matrix P: " << indices.transpose() << "\n";
1490 }
1491
1492 if (indices.size() < work.rows()) {
1493 throw CGException("Failed to select dummy derivatives! "
1494 "The resulting system is probably singular for the provided data.");
1495 }
1496 };
1497
1498 try {
1499 orderColumns();
1500 } catch (const CGException& ex) {
1501 if (avoidCols.empty()) {
1502 throw;
1503 }
1504
1505 if (this->verbosity_ >= Verbosity::Low)
1506 log() << "Failed to determine dummy derivatives without the variables defined to be avoided by the user\n";
1507
1508 // try again with the variables selected to be avoided by the user
1509 for (size_t j : avoidCols)
1510 excludeCols.erase(j);
1511
1512 varsLocal.clear();
1513
1514 orderColumns();
1515 }
1516
1517 const auto& p = qr.colsPermutation();
1518 const auto& indices = p.indices();
1519
1520 std::vector<Vnode<Base>* > newDummies;
1521 if (avoidConvertAlg2DifVars_) {
1522 auto& graph = idxIdentify_->getGraph();
1523 const auto& varInfo = graph.getOriginalVariableInfo();
1524
1525 // add algebraic first
1526 for (int i = 0; newDummies.size() < size_t(work.rows()) && i < qr.rank(); i++) {
1527 Vnode<Base>* v = varsLocal[indices(i)];
1528 CPPADCG_ASSERT_UNKNOWN(v->originalVariable() != nullptr);
1529 size_t tape = v->originalVariable()->tapeIndex();
1530 CPPADCG_ASSERT_UNKNOWN(tape < varInfo.size());
1531 if (varInfo[tape].getDerivative() < 0) {
1532 // derivative of a variable which was originally algebraic only
1533 newDummies.push_back(v);
1534 }
1535 }
1536 // add remaining
1537 for (int i = 0; newDummies.size() < size_t(work.rows()); i++) {
1538 Vnode<Base>* v = varsLocal[indices(i)];
1539 CPPADCG_ASSERT_UNKNOWN(v->originalVariable() != nullptr);
1540 size_t tape = v->originalVariable()->tapeIndex();
1541 CPPADCG_ASSERT_UNKNOWN(tape < varInfo.size());
1542 if (varInfo[tape].getDerivative() >= 0) {
1543 // derivative of a variable which was already differential
1544 newDummies.push_back(v);
1545 }
1546 }
1547
1548 } else {
1549 // use order provided by the householder column pivoting
1550 for (int i = 0; i < work.rows(); i++) {
1551 newDummies.push_back(varsLocal[indices(i)]);
1552 }
1553 }
1554
1555 if (this->verbosity_ >= Verbosity::High) {
1556 log() << "## new dummy derivatives: "; //"(condition = " << bestCond << "): ";
1557 for (Vnode<Base>* it : newDummies)
1558 log() << *it << "; ";
1559 log() << " \n\n";
1560 }
1561#ifndef NDEBUG
1562 for (Vnode<Base>* it : newDummies) {
1563 CPPADCG_ASSERT_UNKNOWN(std::find(dummyD_.begin(), dummyD_.end(), it) == dummyD_.end());
1564 }
1565#endif
1566
1567 dummyD_.insert(dummyD_.end(), newDummies.begin(), newDummies.end());
1568 }
1569
1570 inline static void printModel(std::ostream& out,
1571 CodeHandler<Base>& handler,
1572 const std::vector<CGBase>& res,
1573 const std::vector<DaeVarInfo>& varInfo,
1574 const std::set<size_t>& erasedVariables,
1575 const std::set<size_t>& erasedEquations) {
1576 using std::vector;
1577
1578 std::vector<std::string> indepNames;
1579 for (size_t p = 0; p < varInfo.size(); p++) {
1580 if (erasedVariables.find(p) == erasedVariables.end()) {
1581 // not erased from model
1582 indepNames.push_back(varInfo[p].getName());
1583 }
1584 }
1585 CPPADCG_ASSERT_UNKNOWN(handler.getIndependentVariableSize() == indepNames.size());
1586
1587 LanguageC<Base> lang("double");
1588 vector<CGBase> resAux;
1589 for (size_t p = 0; p < res.size(); ++p) {
1590 if (erasedEquations.find(p) == erasedEquations.end()) {
1591 resAux.push_back(res[p]);
1592 }
1593 }
1594 std::vector<std::string> depNames;
1595 LangCCustomVariableNameGenerator<Base> nameGen(depNames, indepNames);
1596 handler.generateCode(out, lang, resAux, nameGen);
1597 }
1598
1599 inline static void printGraphSparsity(std::ostream& out,
1600 const std::vector<bool>& jacSparsity,
1601 const std::map<size_t, Vnode<Base>*>& tape2FreeVariables,
1602 const std::vector<Enode<Base>*>& equations,
1603 const size_t n) {
1604 for (size_t e = 0; e < equations.size(); ++e) {
1605 Enode<Base>* eq = equations[e];
1606 size_t count = 0;
1607 for (const auto& it : tape2FreeVariables) {
1608 if (jacSparsity[n * eq->index() + it.first]) {
1609 if (count == 0)
1610 out << "# Equation " << e << ": \t";
1611 out << " " << it.second->name();
1612 count++;
1613 }
1614 }
1615 if (count > 0)
1616 out << "\n";
1617 }
1618
1619 out << std::endl;
1620 }
1621
1622 template<class T>
1623 inline static void deleteVectorValues(std::vector<T*>& v) {
1624 for (size_t i = 0; i < v.size(); i++) {
1625 delete v[i];
1626 }
1627 v.clear();
1628 }
1629
1630};
1631
1632} // END cg namespace
1633} // END CppAD namespace
1634
1635#endif
bool augmentPath(Enode< Base > &i) override final
void removeIndependent(Node &indep)
Definition graph_mod.hpp:76
size_t getIndependentVariableSize() const
void makeVariables(VectorCG &variables)
void substituteIndependent(const CGB &indep, const CGB &dep, bool removeFromIndeps=true)
Definition graph_mod.hpp:22
virtual void generateCode(std::ostream &out, Language< Base > &lang, CppAD::vector< CGB > &dependent, VariableNameGenerator< Base > &nameGen, const std::string &jobName="source")
void undoSubstituteIndependent(Node &indep)
Definition graph_mod.hpp:65
ADFun< CG< Base > > & getOriginalModel() const
void setDerivative(int derivative)
const std::string & getName() const
int getAntiDerivative() const
std::vector< bool > jacSparsity_
ADFun< CGBase > * generateReorderedModel(CodeHandler< Base > &handler, const std::vector< CGBase > &res0, const std::vector< DaeVarInfo > &varInfo, const std::vector< DaeVarInfo > &newVarInfo, const std::vector< DaeEquationInfo > &eqInfo, const std::vector< DaeEquationInfo > &newEqInfo) const
std::unique_ptr< ADFun< CGBase > > reorderModelEqNVars(ADFun< CG< Base > > &fun, const std::vector< DaeVarInfo > &varInfo, std::vector< DaeVarInfo > &newVarInfo, const std::vector< DaeEquationInfo > &eqInfo, std::vector< DaeEquationInfo > &newEqInfo)
bool isGenerateSemiExplicitDae() const
void setReorder(bool reorder)
std::vector< Base > normEq_
std::unique_ptr< ADFun< CGBase > > reducedFun_
std::vector< Vnode< Base > * > dummyD_
void selectDummyDerivatives(const std::vector< Enode< Base > * > &eqs, const std::vector< Vnode< Base > * > &vars, MatrixB &work)
std::vector< Base > x_
void setGenerateSemiExplicitDae(bool generateSemiExplicitDae)
void matchVars2Eqs4Elimination(std::vector< DaeVarInfo > &varInfo, std::vector< DaeEquationInfo > &eqInfo)
std::vector< Base > normVar_
void setAvoidVarsAsDummies(const std::set< std::string > &avoidAsDummy)
void setAvoidConvertAlg2DifVars(bool avoid)
std::unique_ptr< ADFun< CGBase > > reduceEquations(const std::vector< DaeVarInfo > &reducedVarInfo, std::vector< DaeVarInfo > &newVarInfo, const std::vector< DaeEquationInfo > &reducedEqInfo, std::vector< DaeEquationInfo > &newEqInfo)
bool isAvoidConvertAlg2DifVars() const
std::unique_ptr< ADFun< CG< Base > > > reduceIndex(std::vector< DaeVarInfo > &newVarInfo, std::vector< DaeEquationInfo > &newEqInfo) override
DummyDerivatives(DaeStructuralIndexReduction< Base > &idxIdentify, const std::vector< Base > &x, const std::vector< Base > &normVar, const std::vector< Base > &normEq)
std::set< std::string > avoidAsDummy_
const std::set< std::string > & getAvoidVarsAsDummies() const
std::unique_ptr< ADFun< CGBase > > generateSemiExplicitDAE(ADFun< CG< Base > > &fun, const std::vector< DaeVarInfo > &varInfo, std::vector< DaeVarInfo > &newVarInfo, const std::vector< DaeEquationInfo > &eqInfo, std::vector< DaeEquationInfo > &newEqInfo)
bool assignVar2Equation(Enode< Base > &i, std::vector< CGBase > &res0, Vnode< Base > &j, std::vector< CGBase > &indep0, CodeHandler< Base > &handler, std::vector< bool > &jacSparsity, const std::map< size_t, Vnode< Base > * > &tape2FreeVariables, std::vector< Enode< Base > * > &equations, std::vector< DaeVarInfo > &varInfo)
DaeStructuralIndexReduction< Base > *const idxIdentify_
void setReduceEquations(bool reduceEquations)
Eigen::SparseMatrix< Base, Eigen::RowMajor > jacobian_
virtual void addDummyDerivatives(const std::vector< DaeVarInfo > &varInfo, const std::vector< DaeEquationInfo > &eqInfo, std::vector< DaeVarInfo > &newVarInfo)
Enode< Base > * derivative() const
std::vector< ActiveOut > evaluate(ArrayView< const ActiveOut > indepNew, ArrayView< const CG< ScalarIn > > depOld)
Definition evaluator.hpp:93
const std::vector< Argument< Base > > & getArguments() const
CGOpCode getOperationType() const
Vnode< Base > * derivative() const
Vnode< Base > * antiDerivative() const