CppADCodeGen 2.4.3
A C++ Algorithmic Differentiation Package with Source Code Generation
Loading...
Searching...
No Matches
bipartite_graph.hpp
1#ifndef CPPAD_CG_BIPARTITE_GRAPH_INCLUDED
2#define CPPAD_CG_BIPARTITE_GRAPH_INCLUDED
3/* --------------------------------------------------------------------------
4 * CppADCodeGen: C++ Algorithmic Differentiation with Source Code Generation:
5 * Copyright (C) 2016 Ciengis
6 *
7 * CppADCodeGen is distributed under multiple licenses:
8 *
9 * - Eclipse Public License Version 1.0 (EPL1), and
10 * - GNU General Public License Version 3 (GPL3).
11 *
12 * EPL1 terms and conditions can be found in the file "epl-v10.txt", while
13 * terms and conditions for the GPL3 can be found in the file "gpl3.txt".
14 * ----------------------------------------------------------------------------
15 * Author: Joao Leal
16 */
17
18#include <cppad/cg/dae_index_reduction/bipartite_nodes.hpp>
19#include <cppad/cg/dae_index_reduction/dae_equation_info.hpp>
20#include <cppad/cg/dae_index_reduction/time_diff.hpp>
21
22namespace CppAD {
23namespace cg {
24
29template<class Base>
31protected:
33 using ADCG = CppAD::AD<CGBase>;
34protected:
42 std::vector<DaeVarInfo> varInfo_;
46 std::vector<bool> sparsity_;
47 // Bipartite graph ([equation i][variable j])
48 std::vector<Vnode<Base>*> vnodes_;
49 std::vector<Enode<Base>*> enodes_;
64private:
65 int timeOrigVarIndex_; // time index in the original user model (may not exist)
66 SimpleLogger& logger_;
67public:
68
77 const std::vector<DaeVarInfo>& varInfo,
78 const std::vector<std::string>& eqName,
79 SimpleLogger& logger) :
80 fun_(&fun),
81 varInfo_(varInfo),
84 preserveNames_(false),
85 timeOrigVarIndex_(-1),
86 logger_(logger) {
87
88 using namespace std;
89 using std::vector;
90
91 CPPADCG_ASSERT_UNKNOWN(fun_ != nullptr);
92 const size_t m = fun.Range(); // equation count
93 const size_t n = fun.Domain(); // total variable count
94
95 CPPADCG_ASSERT_UNKNOWN(varInfo_.size() == n);
96 for (size_t j = 0; j < n; ++j) {
97 varInfo_[j].setOriginalIndex(j);
98 varInfo_[j].setId(j);
99 }
100
101 for (size_t j = 0; j < n; ++j) {
102 int deriv = varInfo_[j].getAntiDerivative();
103 CPPADCG_ASSERT_UNKNOWN(deriv < int(varInfo_.size()));
104 if (deriv >= 0) {
105 varInfo_[deriv].setDerivative(j);
106 }
107 }
108
109 for (size_t j = 0; j < n; ++j) {
110 determineVariableOrder(varInfo_[j]);
111 }
112
113 // create equation nodes
114 enodes_.reserve(1.2 * m + 1);
115 enodes_.resize(m);
116 for (size_t i = 0; i < m; i++) {
117 if (i < eqName.size())
118 enodes_[i] = new Enode<Base>(i, eqName[i]);
119 else
120 enodes_[i] = new Enode<Base>(i);
121 }
122
123 // locate the time variable (if present)
124 for (size_t dj = 0; dj < n; dj++) {
125 if (varInfo_[dj].isIntegratedVariable()) {
126 if (timeOrigVarIndex_ >= 0) {
127 throw CGException("More than one time variable (integrated variable) defined");
128 }
129 timeOrigVarIndex_ = dj;
130 }
131 }
132
133 // determine the order of each time derivative
134 vector<int> derivOrder = determineVariableDiffOrder(varInfo_);
135 map<int, vector<size_t> > order2Tape;
136 for (size_t tape = 0; tape < derivOrder.size(); ++tape) {
137 order2Tape[derivOrder[tape]].push_back(tape);
138 }
139 origMaxTimeDivOrder_ = *std::max_element(derivOrder.begin(), derivOrder.end());
140
144 std::string timeVarName;
145 if (timeOrigVarIndex_ < 0) {
146 timeVarName = "t";
147 } else {
148 if (varInfo_[timeOrigVarIndex_].getName().empty()) {
149 varInfo_[timeOrigVarIndex_].setName("t");
150 }
151 timeVarName = varInfo_[timeOrigVarIndex_].getName();
152 }
153
154 stringstream ss;
155 for (int order = -1; order <= origMaxTimeDivOrder_; order++) {
156 //size_t j = 0; j < varInfo_.size(); j++
157 const vector<size_t>& tapeIndexes = order2Tape[order];
158 if (order < 0) {
159 for (size_t p = 0; p < tapeIndexes.size(); ++p) {
160 DaeVarInfo& var = varInfo_[tapeIndexes[p]];
161 if (var.getName().empty()) {
162 ss << "p" << p;
163 var.setName(ss.str());
164 ss.str("");
165 ss.clear();
166 }
167 }
168
169 } else if (order == 0) {
170 for (size_t p = 0; p < tapeIndexes.size(); ++p) {
171 DaeVarInfo& var = varInfo_[tapeIndexes[p]];
172 if (var.getName().empty()) {
173 ss << "x" << p;
174 var.setName(ss.str());
175 ss.str("");
176 ss.clear();
177 }
178 }
179 } else if (order > 0) {
180 for (size_t p = 0; p < tapeIndexes.size(); ++p) {
181 DaeVarInfo& var = varInfo_[tapeIndexes[p]];
182 if (var.getName().empty()) {
183 const DaeVarInfo& deriv = varInfo_[var.getAntiDerivative()];
184 var.setName("d" + deriv.getName() + "d" + timeVarName);
185 }
186 }
187 }
188 }
189
190 // sort the variables according to the time derivative order (constants are kept out)
191 vector<size_t> new2Tape;
192 vector<int> tape2New(n, -1);
193 new2Tape.reserve(n);
194 for (int order = 0; order <= origMaxTimeDivOrder_; order++) {
195 const vector<size_t>& tapeIndexes = order2Tape[order];
196 for (size_t p = 0; p < tapeIndexes.size(); ++p) {
197 size_t tapeIndex = tapeIndexes[p];
198 tape2New[tapeIndex] = new2Tape.size();
199 new2Tape.push_back(tapeIndex);
200 }
201 }
202
203 // create the variable nodes
204 origTimeDependentCount_ = new2Tape.size();
205 vnodes_.resize(origTimeDependentCount_);
206 for (size_t j = 0; j < vnodes_.size(); j++) {
207 size_t tapeIndex = new2Tape[j];
208 int tapeIndex0 = varInfo_[tapeIndex].getAntiDerivative();
209 const std::string& name = varInfo_[tapeIndex].getName();
210
211 CPPADCG_ASSERT_UNKNOWN(varInfo_[tapeIndex].isFunctionOfIntegrated());
212
213 if (tapeIndex0 < 0) {
214 // generate the variable name
215 vnodes_[j] = new Vnode<Base>(j, tapeIndex, name);
216 } else {
217 Vnode<Base>* derivativeOf = vnodes_[tape2New[tapeIndex0]];
218 vnodes_[j] = new Vnode<Base>(j, tapeIndex, derivativeOf, name);
219 }
220 }
221
222 // create the edges
223 sparsity_ = jacobianSparsity<vector<bool>, CGBase>(fun);
224
225 for (size_t i = 0; i < m; i++) {
226 for (size_t p = 0; p < n; p++) {
227 int j = tape2New[p];
228 if (j >= 0 && sparsity_[i * n + p]) {
229 enodes_[i]->addVariable(vnodes_[j]);
230 }
231 }
232 }
233
234 // make sure the system is not under or over determined
235 size_t nvar = 0;
236 for (size_t j = 0; j < vnodes_.size(); j++) {
237 const Vnode<Base>* jj = vnodes_[j];
238 if (!jj->isParameter() && // exclude constants
239 (jj->antiDerivative() != nullptr || // derivatives
240 jj->derivative() == nullptr) // algebraic variables
241 ) {
242 nvar++;
243 }
244 }
245
246 if (nvar != m) {
247 throw CGException("The system is not well determined. "
248 "The of number of equations (", enodes_.size(), ")"
249 " does not match the number of unknown variables "
250 "(", nvar, ").");
251 }
252 }
253
254 BipartiteGraph(const BipartiteGraph& p) = delete;
255
256 BipartiteGraph& operator=(const BipartiteGraph& p) = delete;
257
258 virtual ~BipartiteGraph() {
259 for (size_t i = 0; i < enodes_.size(); i++)
260 delete enodes_[i];
261
262 for (size_t j = 0; j < vnodes_.size(); j++)
263 delete vnodes_[j];
264 }
265
266
267 inline std::vector<Vnode<Base>*>& variables() {
268 return vnodes_;
269 }
270
271 inline const std::vector<Vnode<Base>*>& variables() const {
272 return vnodes_;
273 }
274
275 inline std::vector<Enode<Base>*>& equations() {
276 return enodes_;
277 }
278
279 inline const std::vector<Enode<Base>*>& equations() const {
280 return enodes_;
281 }
282
283 const std::vector<DaeVarInfo>& getOriginalVariableInfo() const {
284 return varInfo_;
285 }
286
287 inline size_t getOrigTimeDependentCount() const {
289 }
290
296 void setPreserveNames(bool p) {
297 preserveNames_ = p;
298 }
299
305 bool isPreserveNames() const {
306 return preserveNames_;
307 }
308
314 inline size_t getStructuralIndex() const {
315 size_t origM = this->fun_->Range();
316 if (origM == enodes_.size()) {
317 // no index reduction performed: it is either an index 1 DAE or an ODE
318 bool isDAE = false;
319 for (size_t j = 0; j < varInfo_.size(); j++) {
320 const DaeVarInfo& jj = varInfo_[j];
321 if (jj.getDerivative() < 0 && !jj.isIntegratedVariable() && jj.isFunctionOfIntegrated()) {
322 isDAE = true; // found algebraic variable
323 break;
324 }
325 }
326 if (!isDAE) {
327 return 0;
328 } else {
329 return 1;
330 }
331 }
332
333 size_t index = 0;
334 for (size_t i = origM; i < enodes_.size(); i++) {
335 Enode<Base>* ii = enodes_[i];
336 size_t eqOrder = 0;
337 if (ii->derivative() == nullptr) {
338 Enode<Base>* eq = ii;
339 while (eq->derivativeOf() != nullptr) {
340 eq = eq->derivativeOf();
341 eqOrder++;
342 }
343 if (eqOrder > index)
344 index = eqOrder;
345 }
346 }
347
348 return index + 1; // one extra differentiation to get an ODE
349 }
350
351 inline void printResultInfo(const std::string& method) {
352 logger_.log() << "\n" << method << " DAE differentiation/structural index reduction:\n\n"
353 " Equations count: " << enodes_.size() << "\n";
354 for (Enode<Base>* ii : enodes_) {
355 logger_.log() << " " << ii->index() << " - " << *ii << "\n";
356 }
357
358 logger_.log() << "\n Variable count: " << vnodes_.size() << "\n";
359
360 for (const Vnode<Base>* jj : vnodes_) {
361 logger_.log() << " " << jj->index() << " - " << *jj;
362 if (jj->assignmentEquation() != nullptr) {
363 logger_.log() << " assigned to " << *jj->assignmentEquation() << "\n";
364 } else if (jj->isParameter()) {
365 logger_.log() << " is a parameter (time independent)\n";
366 } else {
367 logger_.log() << " NOT assigned to any equation\n";
368 }
369 }
370
371 logger_.log() << "\n Degrees of freedom: " << vnodes_.size() - enodes_.size() << std::endl;
372 }
373
374 inline void uncolorAll() {
375 for (Vnode<Base>* j : vnodes_) {
376 j->uncolor();
377 }
378
379 for (Enode<Base>* i : enodes_) {
380 i->uncolor();
381 }
382 }
383
384 inline Vnode<Base>* createDerivate(Vnode<Base>& j) {
385 if (j.derivative() != nullptr)
386 return j.derivative();
387
388 // add new variable derivatives of colored variables
389 size_t newVarCount = vnodes_.size() - origTimeDependentCount_;
390 size_t tapeIndex = varInfo_.size() + newVarCount;
391
392 Vnode<Base>* jDiff = new Vnode<Base> (vnodes_.size(), tapeIndex, &j);
393 vnodes_.push_back(jDiff);
394
395 if (logger_.getVerbosity() >= Verbosity::High)
396 logger_.log() << "Created " << *jDiff << "\n";
397
398 return jDiff;
399 }
400
401 inline Enode<Base>* createDerivate(Enode<Base>& i,
402 bool addOrigVars = true) {
403 if (i.derivative() != nullptr)
404 return i.derivative();
405
406 Enode<Base>* iDiff = new Enode<Base> (enodes_.size(), &i);
407 enodes_.push_back(iDiff);
408
409 // differentiate newI and create edges!!!
410 dirtyDifferentiateEq(i, *iDiff, addOrigVars);
411
412 if (logger_.getVerbosity() >= Verbosity::High)
413 logger_.log() << "Created " << *iDiff << "\n";
414
415 return iDiff;
416 }
417
427 inline void remove(const Enode<Base>& i) {
428 CPPADCG_ASSERT_UNKNOWN(enodes_[i.index()] == &i);
429 CPPADCG_ASSERT_UNKNOWN(i.derivative() == nullptr);
430
431 for (Vnode<Base>* j: i.variables()) {
432 // remove the edges (connections in variables)
433 auto& eqs = j->equations();
434 auto it = std::find(eqs.begin(), eqs.end(), &i);
435 CPPADCG_ASSERT_UNKNOWN(it != eqs.end());
436 eqs.erase(it);
437
441 while(j->equations().empty()) {
442 CPPADCG_ASSERT_UNKNOWN(vnodes_[j->index()] == j);
443
444 if (j->derivative() == nullptr) {
445 vnodes_.erase(vnodes_.cbegin() + j->index());
446
447 // update variable indices
448 for (size_t jj = j->index(); jj < vnodes_.size(); ++jj) {
449 vnodes_[jj]->setTapeIndex(vnodes_[jj]->tapeIndex() - 1);
450 vnodes_[jj]->setIndex(vnodes_[jj]->index() - 1);
451 }
452
453 auto* jOrig = j->antiDerivative();
454 CPPADCG_ASSERT_UNKNOWN(jOrig != nullptr);
455 jOrig->setDerivative(nullptr);
456
457 delete j; // no longer required
458 j = jOrig;
459 }
460 }
461
462 }
463
464 // update equation indices
465 for (size_t ii = i.index() + 1; ii < enodes_.size(); ++ii) {
466 CPPADCG_ASSERT_UNKNOWN(enodes_[ii]->index() > 0);
467 CPPADCG_ASSERT_UNKNOWN(enodes_[ii]->index() == ii);
468 enodes_[ii]->setIndex(enodes_[ii]->index() - 1);
469 }
470
471 if(i.derivativeOf() != nullptr) {
472 i.derivativeOf()->setDerivative(nullptr);
473 }
474
475 auto it = std::find(enodes_.begin(), enodes_.end(), &i);
476 CPPADCG_ASSERT_UNKNOWN(it != enodes_.end());
477 enodes_.erase(it);
478
479 delete &i; // no longer required
480 }
481
499 Enode<Base>& iDiff,
500 bool addOrigVars = true) {
501 for (Vnode<Base>* jj : i.originalVariables()) {
502 if(addOrigVars) {
503 iDiff.addVariable(jj);
504 }
505
506 if (jj->derivative() != nullptr) {
507 iDiff.addVariable(jj->derivative());
508 } else if(!jj->isParameter()) {
509 iDiff.addVariable(createDerivate(*jj));
510 }
511 }
512 }
513
517 inline std::unique_ptr<ADFun<CGBase>> generateNewModel(std::vector<DaeVarInfo>& newVarInfo,
518 std::vector<DaeEquationInfo>& equationInfo,
519 const std::vector<Base>& x) {
520 using std::vector;
521
522 std::unique_ptr<ADFun<CGBase> > reducedFun;
523
524 vector<vector<Enode<Base>*> > newEquations;
525
526 // find new equations that must be generated by differentiation
527 vector<Enode<Base>*> newEqs;
528 size_t origM = this->fun_->Range();
529 for (size_t i = 0; i < origM; i++) {
530 if (enodes_[i]->derivative() != nullptr) {
531 CPPADCG_ASSERT_UNKNOWN(enodes_[i]->derivativeOf() == nullptr);
532 newEqs.push_back(enodes_[i]->derivative());
533 }
534 }
535
536 while (newEqs.size() > 0) {
537 newEquations.push_back(newEqs);
538 newEqs.clear();
539 vector<Enode<Base>*>& eqs = newEquations.back();
540 for (size_t i = 0; i < eqs.size(); i++) {
541 if (eqs[i]->derivative() != nullptr) {
542 newEqs.push_back(eqs[i]->derivative());
543 }
544 }
545 }
546
547 if (newEquations.empty()) {
548 // nothing to do
549 return nullptr;
550 }
551
559 newVarInfo = varInfo_; // copy
560 size_t newVars = vnodes_.size() - origTimeDependentCount_;
561 newVarInfo.reserve(varInfo_.size() + newVars);
562 for (size_t j = origTimeDependentCount_; j < vnodes_.size(); j++) {
563 // new variable derivative added by the Pantelides method
564 Vnode<Base>* jj = vnodes_[j];
565 CPPADCG_ASSERT_UNKNOWN(jj->antiDerivative() != nullptr);
566 size_t antiDeriv = jj->antiDerivative()->tapeIndex();
567 size_t id = newVarInfo.size();
568 newVarInfo.push_back(DaeVarInfo(antiDeriv, jj->name(), id)); // create the new variable
569 DaeVarInfo& newVar = newVarInfo.back();
570 DaeVarInfo& newAntiDeriv = newVarInfo[antiDeriv];
571
572 newAntiDeriv.setDerivative(jj->tapeIndex()); // update the antiderivative
573 newVar.setOrder(newAntiDeriv.getOrder() + 1);
574 newVar.setOriginalAntiDerivative(newVar.getOrder() == 1 ? newAntiDeriv.getOriginalIndex() : newAntiDeriv.getOriginalAntiDerivative());
575 if (jj->derivative() != nullptr) {
576 newVar.setDerivative(jj->derivative()->tapeIndex());
577 }
578 }
579
580 std::map<Enode<Base>*, Vnode<Base>*> assignments;
581 for (size_t j = 0; j < vnodes_.size(); j++) {
582 Vnode<Base>* jj = vnodes_[j];
583 if (jj->assignmentEquation() != nullptr) {
584 assignments[jj->assignmentEquation()] = jj;
585 }
586 }
587
588 equationInfo.resize(enodes_.size());
589 for (size_t i = 0; i < enodes_.size(); i++) {
590 Enode<Base>* ii = enodes_[i];
591 int derivativeOf = ii->derivativeOf() != nullptr ? ii->derivativeOf()->index() : -1;
592 int origIndex = ii->derivativeOf() == nullptr ? i : -1;
593 int assignedVarIndex = assignments.count(ii) > 0 ? assignments[ii]->tapeIndex() : -1;
594
595 equationInfo[i] = DaeEquationInfo(i, origIndex, derivativeOf, assignedVarIndex);
596 }
597
598 size_t timeTapeIndex;
599 {
600 CodeHandler<Base> handler;
601
602 vector<CGBase> indep0(this->fun_->Domain());
603 handler.makeVariables(indep0);
604
605 const vector<CGBase> dep0 = forward0(*this->fun_, indep0);
606
611 vector<ADCG> indepNew;
612 if (timeOrigVarIndex_ >= 0) {
613 indepNew = vector<ADCG>(newVarInfo.size()); // variables + time (vnodes include time)
614 timeTapeIndex = timeOrigVarIndex_;
615 } else {
616 indepNew = vector<ADCG>(newVarInfo.size() + 1); // variables + time (new time variable added)
617 timeTapeIndex = indepNew.size() - 1;
618 }
619
620 // initialize with the user provided values
621 for (size_t j = 0; j < x.size(); j++) {
622 indepNew[j] = x[j];
623 }
624 Independent(indepNew);
625
626 // variables with the relationship between x dxdt and t
627 vector<ADCG> indep2 = prepareTimeDependentVariables(indepNew, newVarInfo, timeTapeIndex);
628 indep2.resize(indep0.size());
629
630 Evaluator<Base, CGBase> evaluator0(handler);
631 evaluator0.setPrintFor(preserveNames_); // variable names saved with CppAD::PrintFor
632 vector<ADCG> depNew = evaluator0.evaluate(indep2, dep0);
633 depNew.resize(enodes_.size());
634
635 try {
636 reducedFun.reset(new ADFun<CGBase>(indepNew, depNew));
637 } catch (const std::exception& ex) {
638 throw CGException("Failed to create ADFun: ", ex.what());
639 }
640
641 if (logger_.getVerbosity() >= Verbosity::High) {
642 logger_.log() << "Original model:\n";
643 printModel(logger_.log(), *reducedFun, newVarInfo, equationInfo);
644 }
645 }
646
647
652 for (size_t d = 0; d < newEquations.size(); d++) {
653 vector<Enode<Base>*>& equations = newEquations[d];
654
655 size_t m = reducedFun->Domain(); // total variable count
656 //size_t n = reducedFun->Range(); // equation count
657
661 CodeHandler<Base> handler0;
662
663 vector<CGBase> indep0(m);
664 handler0.makeVariables(indep0);
665
666 vector<CGBase> dep = forward0(*reducedFun, indep0);
667
671 //forwardTimeDiff(equations, dep, timeTapeIndex);
672 reverseTimeDiff(*reducedFun, equations, dep, timeTapeIndex);
673
677 vector<ADCG> indep2;
678 vector<ADCG> indepNew;
679
680 if (d < newEquations.size() - 1) {
681 indepNew.resize(m);
682 } else if (timeOrigVarIndex_ < 0) {
683 // the very last model creation
684 indepNew.resize(m - 1); // take out time (it was added by this function and not the user)
685 } else {
686 // the very last model creation
687 indepNew.resize(m);
688 }
689
690 for (size_t j = 0; j < x.size(); j++) {
691 indepNew[j] = x[j];
692 }
693 Independent(indepNew);
694
695 if (d < newEquations.size() - 1) {
696 // variables with the relationship between x, dxdt and t
697 indep2 = prepareTimeDependentVariables(indepNew, newVarInfo, timeTapeIndex);
698 } else {
699 indep2 = indepNew;
700 indep2.resize(m);
701 }
702
703 Evaluator<Base, CGBase> evaluator(handler0);
704 evaluator.setPrintFor(preserveNames_); // variable names saved with CppAD::PrintFor
705 vector<ADCG> depNew = evaluator.evaluate(indep2, dep);
706
707 try {
708 reducedFun.reset(new ADFun<CGBase>(indepNew, depNew));
709 } catch (const std::exception& ex) {
710 throw CGException("Failed to create ADFun: ", ex.what());
711 }
712
713 if (logger_.getVerbosity() >= Verbosity::High) {
714 logger_.log() << equations.size() << " new equations:\n";
715 printModel(logger_.log(), *reducedFun, newVarInfo, equationInfo);
716 }
717 }
718
719 return reducedFun;
720 }
721
722 inline static void forwardTimeDiff(ADFun<CGBase>& reducedFun,
723 const std::vector<Enode<Base>*>& equations,
724 std::vector<CG<Base> >& dep,
725 size_t tapeTimeIndex) {
726
727 size_t m = reducedFun.Domain();
728
729 std::vector<CGBase> u(m, CGBase(0));
730 u[tapeTimeIndex] = CGBase(1);
731 std::vector<CGBase> v;
732 try {
733 v = reducedFun.Forward(1, u);
734 } catch (const std::exception& ex) {
735 throw CGException("Failed to determine model Jacobian (forward mode): ", ex.what());
736 }
737
738 for (size_t e = 0; e < equations.size(); e++) {
739 dep[equations[e]->index()] = v[equations[e]->derivativeOf()->index()];
740 }
741 }
742
743 inline static void reverseTimeDiff(ADFun<CGBase>& reducedFun,
744 const std::vector<Enode<Base>*>& equations,
745 std::vector<CG<Base> >& dep,
746 size_t tapeTimeIndex) {
747 size_t m = reducedFun.Domain();
748 size_t n = reducedFun.Range();
749 std::vector<CGBase> u(m);
750 std::vector<CGBase> v(n);
751
752 for (size_t e = 0; e < equations.size(); e++) {
753 size_t i = equations[e]->derivativeOf()->index();
754 if (reducedFun.Parameter(i)) { // return zero for this component of f
755 dep[equations[e]->index()] = 0;
756 } else {
757 // set v to the i-th coordinate direction
758 v[i] = 1;
759
760 // compute the derivative of this component of f
761 try {
762 u = reducedFun.Reverse(1, v);
763 } catch (const std::exception& ex) {
764 throw CGException("Failed to determine model Jacobian (reverse mode): ", ex.what());
765 }
766
767 // reset v to vector of all zeros
768 v[i] = 0;
769
770 // return the result
771 dep[equations[e]->index()] = u[tapeTimeIndex];
772 }
773 }
774 }
775
784 inline std::vector<CppAD::AD<CG<Base> > > prepareTimeDependentVariables(const std::vector<ADCG>& indepOrig,
785 const std::vector<DaeVarInfo>& newVarInfo,
786 size_t timeTapeIndex) const {
787 CPPADCG_ASSERT_UNKNOWN(timeTapeIndex < indepOrig.size());
788
789 using std::vector;
790 using ADCGBase = CppAD::AD<CGBase>;
791
792 vector<ADCGBase> indepOut(indepOrig.size());
793 vector<ADCGBase> ax(3); // function inputs
794 vector<ADCGBase> ay(1); // function output
795
796 ax[2] = indepOrig[timeTapeIndex]; // time
797
798 for (size_t j = 0; j < newVarInfo.size(); j++) {
799 const DaeVarInfo& jj = newVarInfo[j];
800 if (jj.getDerivative() >= 0) {
801 ax[0] = indepOrig[j]; // x
802 ax[1] = indepOrig[jj.getDerivative()]; // dxdt
803 time_var(0, ax, ay);
804 indepOut[j] = ay[0];
805 } else {
806 indepOut[j] = indepOrig[j];
807 }
808 }
809
810 if (newVarInfo.size() < indepOrig.size()) {
811 indepOut[indepOut.size() - 1] = indepOrig[timeTapeIndex];
812 }
813
814 return indepOut;
815 }
816
817 inline void printModel(std::ostream& out,
818 ADFun<CG<Base> >* fun) {
819 printModel(out, fun, varInfo_);
820 }
821
827 inline void printModel(std::ostream& out,
828 ADFun<CG<Base> >& fun,
829 const std::vector<DaeVarInfo>& varInfo,
830 const std::vector<DaeEquationInfo>& eqInfo) const {
831 std::vector<std::string> vnames(varInfo.size());
832 for (size_t i = 0; i < varInfo.size(); ++i) {
833 vnames[i] = varInfo[i].getName();
834 }
835 std::vector<std::string> eqnames(eqInfo.size());
836 for (size_t i = 0; i < eqInfo.size(); ++i) {
837 if(eqInfo[i].isExplicit()) {
838 CPPADCG_ASSERT_UNKNOWN(eqInfo[i].getAssignedVarIndex() >= 0);
839 eqnames[i] = "d" + varInfo[eqInfo[i].getAssignedVarIndex()].getName() + "dt";
840 } else {
841 eqnames[i] = "res[" + std::to_string(i) + "]";
842 }
843 }
844
845 printModel(out, fun, vnames, eqnames);
846 }
847
855 inline void printModel(std::ostream& out,
856 ADFun<CG<Base> >& fun,
857 const std::vector<std::string>& indepNames,
858 const std::vector<std::string>& depNames = std::vector<std::string>()) const {
859 using std::vector;
860
861 CPPADCG_ASSERT_UNKNOWN(fun.Domain() == indepNames.size() || fun.Domain() == indepNames.size() + 1); // with or without time
862
863 CodeHandler<Base> handler;
864
865 vector<CGBase> indep0(fun.Domain());
866 handler.makeVariables(indep0);
867
868 vector<CGBase> dep0 = forward0(fun, indep0);
869
870 LanguageC<double> langC("double");
871
875 LangCCustomVariableNameGenerator<double> nameGen(depNames, indepNames, "res");
876
877 std::ostringstream code;
878 handler.generateCode(code, langC, dep0, nameGen);
879 out << "\n" << code.str() << std::endl;
880 }
881
882 inline void printDot(std::ostream& out) const {
883 out << "digraph {\n";
884 out << " overlap=false\n";
885 out << " rankdir=LR\n";
886 out << " node [style=filled, fillcolor=\"#bdcef5\", color=\"#17128e\"]\n";
887 out << " edge [splines=false, dir=none]\n";
888
889 // variables
890 out << " subgraph variables {\n";
891 out << " rank=min\n";
892 for (const Vnode<Base>* j : vnodes_) {
893 if(!j->isDeleted()) {
894 out << " v" << j->index() << " [label=\"" << j->name() << "\"";
895 if (j->isColored())
896 out << ", color=\"#17c68e\"";
897 out << "]\n";
898 }
899 }
900 out << " }\n";
901
902 // equations
903 out << " subgraph equations {\n";
904 out << " rank=max\n";
905 for (const Enode<Base>* i : enodes_) {
906 out << " e" << i->index() << " [label=\"" << i->name() << "\"";
907 if (i->isColored())
908 out << ", color=\"#17c68e\"";
909 out << "]\n";
910 }
911 out << " }\n";
912
913 // derivatives of equations
914 out << " subgraph eq_derivatives {\n";
915 out << " edge[dir=forward, color=grey]\n";
916 for (const Enode<Base>* i : enodes_) {
917 if (i->derivative() != nullptr && i->derivativeOf() == nullptr) {
918 while (i->derivative() != nullptr) {
919 out << " e" << i->index() << ":e -> e" << i->derivative()->index() << ":e\n";
920 i = i->derivative();
921 }
922 }
923 }
924 out << " }\n";
925
926 // derivatives of variables
927 out << " subgraph var_derivatives {\n";
928 out << " edge[dir=forward, color=grey]\n";
929 for (const Vnode<Base>* j : vnodes_) {
930 if (!j->isDeleted() && j->derivative() != nullptr && (j->antiDerivative() == nullptr || j->antiDerivative()->isDeleted())) {
931 if (!j->derivative()->isDeleted()) {
932 while (j->derivative() != nullptr && !j->derivative()->isDeleted()) {
933 out << " v" << j->index() << ":w -> v" << j->derivative()->index() << ":w\n";
934 j = j->derivative();
935 }
936 }
937 }
938 }
939 out << " }\n";
940
941 // edges
942 for (const Enode<Base>* i : enodes_) {
943 bool added = false;
944 for (const Vnode<Base>* j : i->originalVariables()) {
945 if (!j->isDeleted() && j->assignmentEquation() != i) {
946 if(!added) {
947 out << " ";
948 added = true;
949 }
950 out << "e" << i->index() << " -> v" << j->index() << " ";
951 }
952 }
953 if (added)
954 out << "\n";
955 }
956
957 out << " subgraph assigned {\n";
958 out << " edge[color=blue,penwidth=3.0,style=dashed]\n";
959 for (const Enode<Base>* i : enodes_) {
960 bool added = false;
961
962 for (const Vnode<Base>* j : i->originalVariables()) {
963 if (!j->isDeleted() && j->assignmentEquation() == i) {
964 if(!added) {
965 out << " ";
966 added = true;
967 }
968 out << "e" << i->index() << " -> v" << j->index() << " ";
969 }
970 }
971
972 if (added)
973 out << "\n";
974 }
975
976 out << " }\n";
977 out << "}\n";
978 }
979
980 template<class VectorCGB>
981 inline VectorCGB forward0(ADFun<CGBase>& fun,
982 const VectorCGB& indep0) const {
983
984 if (preserveNames_) {
985 // stream buffer is used to reload names saved with CppAD::PrintFor()
986 OperationNodeNameStreambuf<double> b;
987 std::ostream out(&b);
988
989 return fun.Forward(0, indep0, out);
990 } else {
991 return fun.Forward(0, indep0);
992 }
993 }
994
995 static inline std::vector<int> determineVariableDiffOrder(const std::vector<DaeVarInfo>& varInfo) {
996 size_t n = varInfo.size();
997 // determine the order of each time derivative
998 std::vector<int> derivOrder(n, 0);
999 for (size_t dj = 0; dj < n; dj++) {
1000 size_t j0;
1001 derivOrder[dj] = determineVariableDiffOrder(varInfo, dj, j0);
1002 }
1003
1004 return derivOrder;
1005 }
1006
1007 static inline int determineVariableDiffOrder(const std::vector<DaeVarInfo>& varInfo, size_t index, size_t& j0) {
1008 int derivOrder = -1;
1009 j0 = index;
1010 if (varInfo[index].isFunctionOfIntegrated()) {
1011 derivOrder = 0;
1012 while (varInfo[j0].getAntiDerivative() >= 0) {
1013 CPPADCG_ASSERT_UNKNOWN(j0 < varInfo.size());
1014 CPPADCG_ASSERT_UNKNOWN(varInfo[j0].isFunctionOfIntegrated());
1015 derivOrder++;
1016 j0 = varInfo[j0].getAntiDerivative();
1017 }
1018 }
1019
1020 return derivOrder;
1021 }
1022
1023private:
1024 inline void determineVariableOrder(DaeVarInfo& var) {
1025 if (var.getAntiDerivative() >= 0) {
1026 DaeVarInfo& antiD = varInfo_[var.getAntiDerivative()];
1027 if (antiD.getOriginalAntiDerivative() < 0) {
1028 determineVariableOrder(antiD);
1029 }
1030 var.setOrder(antiD.getOrder() + 1);
1031 var.setOriginalAntiDerivative(var.getOrder() == 1 ? antiD.getOriginalIndex() : antiD.getOriginalAntiDerivative());
1032 }
1033 }
1034};
1035
1036template<class Base>
1037inline std::ostream& operator<<(std::ostream& os,
1038 const BipartiteGraph<Base>& g) {
1039 for (const Enode<Base>* i : g.equations()) {
1040 for (const Vnode<Base>* j : i->originalVariables()) {
1041 os << i->name();
1042 if (j->isDeleted()) {
1043 os << "~~";
1044 } else if (j->assignmentEquation() == i) {
1045 os << "==";
1046 } else {
1047 os << "--";
1048 }
1049 os << j->name() << " ";
1050 }
1051 os << std::endl;
1052 }
1053
1054 return os;
1055}
1056
1057} // END cg namespace
1058} // END CppAD namespace
1059
1060#endif
std::vector< CppAD::AD< CG< Base > > > prepareTimeDependentVariables(const std::vector< ADCG > &indepOrig, const std::vector< DaeVarInfo > &newVarInfo, size_t timeTapeIndex) const
void printModel(std::ostream &out, ADFun< CG< Base > > &fun, const std::vector< DaeVarInfo > &varInfo, const std::vector< DaeEquationInfo > &eqInfo) const
ADFun< CG< Base > > *const fun_
void dirtyDifferentiateEq(Enode< Base > &i, Enode< Base > &iDiff, bool addOrigVars=true)
std::vector< bool > sparsity_
std::unique_ptr< ADFun< CGBase > > generateNewModel(std::vector< DaeVarInfo > &newVarInfo, std::vector< DaeEquationInfo > &equationInfo, const std::vector< Base > &x)
void printModel(std::ostream &out, ADFun< CG< Base > > &fun, const std::vector< std::string > &indepNames, const std::vector< std::string > &depNames=std::vector< std::string >()) const
std::vector< DaeVarInfo > varInfo_
void remove(const Enode< Base > &i)
BipartiteGraph(ADFun< CG< Base > > &fun, const std::vector< DaeVarInfo > &varInfo, const std::vector< std::string > &eqName, SimpleLogger &logger)
void makeVariables(VectorCG &variables)
virtual void generateCode(std::ostream &out, Language< Base > &lang, CppAD::vector< CGB > &dependent, VariableNameGenerator< Base > &nameGen, const std::string &jobName="source")
void setName(const std::string &name)
bool isFunctionOfIntegrated() const
void setDerivative(int derivative)
void setOrder(int order)
void setOriginalAntiDerivative(int originalAntiDerivative)
const std::string & getName() const
bool isIntegratedVariable() const
int getOriginalAntiDerivative() const
int getAntiDerivative() const
Enode< Base > * derivative() const
std::vector< ActiveOut > evaluate(ArrayView< const ActiveOut > indepNew, ArrayView< const CG< ScalarIn > > depOld)
Definition evaluator.hpp:93
Vnode< Base > * derivative() const
Vnode< Base > * antiDerivative() const