1 #ifndef CPPAD_CG_DUMMY_DERIV_INCLUDED
2 #define CPPAD_CG_DUMMY_DERIV_INCLUDED
18 #include <Eigen/Dense>
19 #include <Eigen/Sparse>
20 #include <Eigen/Eigenvalues>
24 #include <cppad/cg/dae_index_reduction/pantelides.hpp>
25 #include <cppad/cg/dae_index_reduction/dummy_deriv_util.hpp>
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>;
75 Eigen::SparseMatrix<Base, Eigen::RowMajor>
jacobian_;
96 bool avoidConvertAlg2DifVars_;
115 const std::vector<Base>& x,
116 const std::vector<Base>& normVar,
117 const std::vector<Base>& normEq) :
128 avoidConvertAlg2DifVars_(true) {
130 for (
Vnode<Base>* jj : idxIdentify.getGraph().variables()) {
131 if (jj->antiDerivative() !=
nullptr) {
132 diffVarStart_ = jj->index();
149 return avoidConvertAlg2DifVars_;
160 avoidConvertAlg2DifVars_ = avoid;
239 inline std::unique_ptr<ADFun<CG<Base>>>
reduceIndex(std::vector<DaeVarInfo>& newVarInfo,
240 std::vector<DaeEquationInfo>& newEqInfo)
override {
245 std::vector<DaeVarInfo> reducedVarInfo;
249 std::vector<DaeEquationInfo> reducedEqInfo;
255 if (this->verbosity_ >= Verbosity::Low)
256 log() <<
"######## Dummy derivatives method ########" << std::endl;
258 newEqInfo = reducedEqInfo;
266 std::vector<DaeVarInfo> varInfo = newVarInfo;
267 std::vector<DaeEquationInfo> eqInfo = newEqInfo;
268 std::unique_ptr<ADFun<CG<Base> > > funShort =
reduceEquations(varInfo, newVarInfo,
274 std::vector<DaeVarInfo> varInfo = newVarInfo;
275 std::vector<DaeEquationInfo> eqInfo = newEqInfo;
284 std::vector<DaeVarInfo> varInfo = newVarInfo;
285 std::vector<DaeEquationInfo> eqInfo = newEqInfo;
292 return std::unique_ptr<ADFun<CG<Base>>>(
reducedFun_.release());
300 const std::vector<DaeEquationInfo>& eqInfo,
301 std::vector<DaeVarInfo>& newVarInfo) {
303 auto& vnodes = graph.variables();
304 auto& enodes = graph.equations();
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) {
320 std::sort(vars.begin(), vars.end(), sortVnodesByOrder<Base>);
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) {
328 if (ii->derivativeOf() !=
nullptr && ii->
derivative() ==
nullptr) {
338 if (this->verbosity_ >= Verbosity::High) {
339 log() <<
"# equation selection: ";
340 for (
size_t i = 0; i < eqs.size(); i++)
341 log() << *eqs[i] <<
"; ";
344 log() <<
"# variable selection: ";
345 for (
size_t j = 0; j < vars.size(); j++)
346 log() << *vars[j] <<
"; ";
359 std::vector<Enode<Base>*> newEqs;
360 newEqs.reserve(eqs.size());
364 if (ii !=
nullptr && ii->derivativeOf() !=
nullptr) {
365 newEqs.push_back(ii);
380 std::vector<Vnode<Base>*> varsNew;
381 varsNew.reserve(vars.size());
385 varsNew.push_back(v);
395 newVarInfo = varInfo;
397 CPPADCG_ASSERT_UNKNOWN(j->tapeIndex() >= 0);
398 CPPADCG_ASSERT_UNKNOWN(j->antiDerivative() !=
nullptr);
399 CPPADCG_ASSERT_UNKNOWN(j->antiDerivative()->tapeIndex() >= 0);
401 newVarInfo[j->tapeIndex()].setAntiDerivative(-1);
402 newVarInfo[j->antiDerivative()->tapeIndex()].setDerivative(-1);
405 if (this->verbosity_ >= Verbosity::Low) {
406 log() <<
"## dummy derivatives:\n";
409 log() <<
"# " << *j <<
" \t" << newVarInfo[j->tapeIndex()].getName() <<
"\n";
411 if (this->verbosity_ >= Verbosity::High) {
412 graph.printModel(log(), *
reducedFun_, newVarInfo, eqInfo);
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) {
436 auto& enodes = graph.equations();
438 CPPADCG_ASSERT_UNKNOWN(reducedVarInfo.size() ==
reducedFun_->Domain());
439 CPPADCG_ASSERT_UNKNOWN(reducedEqInfo.size() ==
reducedFun_->Range());
441 newEqInfo = reducedEqInfo;
453 map<int, int> assignedVar2Eq;
454 for (
size_t i = 0; i < newEqInfo.size(); ++i) {
456 if (newEq.getAssignedVarIndex() >= 0)
457 assignedVar2Eq[newEq.getAssignedVarIndex()] = i;
466 for (
size_t i = 0; i < eqIndexReduced2Short.size(); i++) {
467 eqIndexReduced2Short[i] = i;
475 vector<int> tapeIndexReduced2Short(reducedVarInfo.size());
476 for (
size_t j = 0; j < tapeIndexReduced2Short.size(); j++) {
477 tapeIndexReduced2Short[j] = j;
484 set<size_t> erasedVariables;
485 set<size_t> erasedEquations;
487 if (this->verbosity_ >= Verbosity::High)
488 log() <<
"Reducing total number of equations by symbolic manipulation:" << std::endl;
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;
502 int bestEquation = ita->second;
507 tapeIndexReduced2Short[dummy->tapeIndex()] = -1;
508 eqIndexReduced2Short[bestEquation] = -1;
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);
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;
525 for (
size_t i = 0; i < eqIndexReduced2Short.size(); i++) {
526 if (eqIndexReduced2Short[i] < 0) {
527 for (
size_t ii = i + 1; ii < eqIndexReduced2Short.size(); ii++) {
528 if (eqIndexReduced2Short[ii] >= 0) {
529 eqIndexReduced2Short[ii]--;
536 for (
size_t p = 0; p < tapeIndexReduced2Short.size(); p++) {
537 if (tapeIndexReduced2Short[p] < 0) {
539 for (
size_t p2 = p + 1; p2 < tapeIndexReduced2Short.size(); p2++) {
540 if (tapeIndexReduced2Short[p2] >= 0) {
541 tapeIndexReduced2Short[p2]--;
550 CPPADCG_ASSERT_UNKNOWN(tapeIndexReduced2Short.size() == reducedVarInfo.size());
552 newVarInfo = reducedVarInfo;
553 for (
int p = tapeIndexReduced2Short.size() - 1; p >= 0; p--) {
554 if (tapeIndexReduced2Short[p] < 0) {
555 newVarInfo.erase(newVarInfo.begin() + p);
556 for (
size_t pp = 0; pp < tapeIndexReduced2Short.size(); pp++) {
561 v.setAntiDerivative(-1);
572 for (
int p = eqIndexReduced2Short.size() - 1; p >= 0; p--) {
573 if (eqIndexReduced2Short[p] < 0) {
574 newEqInfo.erase(newEqInfo.begin() + 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()]);
590 reducedVarInfo, newVarInfo,
591 reducedEqInfo, newEqInfo));
593 if (this->verbosity_ >= Verbosity::High) {
594 log() <<
"DAE with less equations and variables:\n";
595 graph.printModel(log(), *shortFun, newVarInfo, newEqInfo);
612 const std::vector<DaeVarInfo>& varInfo,
613 std::vector<DaeVarInfo>& newVarInfo,
614 const std::vector<DaeEquationInfo>& eqInfo,
615 std::vector<DaeEquationInfo>& newEqInfo) {
634 map<int, int> assignedVar2Eq;
635 for (
size_t i = 0; i < newEqInfo.size(); ++i) {
637 assignedVar2Eq[newEq.getAssignedVarIndex()] = i;
643 for (
size_t j = 0; j < varInfo.size(); ++j) {
648 CGBase& indep = indep0[j];
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());
657 int bestEquation = ita->second;
660 CGBase& dep = res0[bestEquation];
665 CPPADCG_ASSERT_UNKNOWN(alias !=
nullptr && alias->
getOperationType() == CGOpCode::Alias);
666 dep.getOperationNode()->makeAlias(alias->
getArguments()[0]);
669 newEqInfo[bestEquation].setExplicit(
true);
674 throw CGException(
"Failed to generate semi-explicit DAE: ", ex.what());
684 for (
size_t j = 0; j != varInfo.size(); ++j) {
686 if (varInfo[j].getAntiDerivative() < 0) {
687 varIndexOld2New[j] = count++;
691 for (
size_t i = 0; i < newEqInfo.size(); ++i) {
693 int j = ii.getAssignedVarIndex();
695 newEqInfo[i].setAssignedVarIndex(varIndexOld2New[j]);
701 newVarInfo = varInfo;
702 for (
int j = newVarInfo.size() - 1; j >= 0; --j) {
703 if (newVarInfo[j].getAntiDerivative() >= 0) {
705 newVarInfo.erase(newVarInfo.begin() + j);
708 for (
size_t j = 0; j < newVarInfo.size(); j++) {
709 newVarInfo[j].setDerivative(-1);
716 std::unique_ptr<ADFun<CGBase> > semiExplicitFun(
generateReorderedModel(handler, res0, varInfo, newVarInfo, eqInfo, newEqInfo));
718 if (this->verbosity_ >= Verbosity::High) {
719 log() <<
"Semi-Eplicit DAE:\n";
720 graph.printModel(log(), *semiExplicitFun, newVarInfo, newEqInfo);
723 return semiExplicitFun;
727 std::vector<DaeEquationInfo>& eqInfo) {
732 auto& vnodes = graph.variables();
733 auto& enodes = graph.equations();
735 CPPADCG_ASSERT_UNKNOWN(eqInfo.size() == enodes.size());
736 CPPADCG_ASSERT_UNKNOWN(varInfo.size() ==
reducedFun_->Domain());
737 CPPADCG_ASSERT_UNKNOWN(eqInfo.size() ==
reducedFun_->Range());
746 vector<bool> jacSparsity = jacobianSparsity<vector<bool> >(*reducedFun_);
758 for (
size_t j = 0; j < vnodes.size(); j++) {
762 dummyVariables.push_back(
new Vnode<Base>(j, v->tapeIndex(), v->name()));
763 eliminateOrig2New[v] = dummyVariables.back();
767 diffVariables.push_back(
new Vnode<Base>(j, v->tapeIndex(), v->name()));
768 eliminateOrig2New[v] = diffVariables.back();
773 variables.reserve(diffVariables.size() + dummyVariables.size());
774 variables.insert(variables.end(), diffVariables.begin(), diffVariables.end());
775 variables.insert(variables.end(), dummyVariables.begin(), dummyVariables.end());
778 for (
size_t i = 0; i < enodes.size(); i++) {
779 equations[i] =
new Enode<Base>(i, enodes[i]->name());
781 for (
size_t p = 0; p < origVars.size(); p++) {
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()]) {
791 CGBase& indep = indep0[j->tapeIndex()];
793 if (handler.isSolvable(*dep.getOperationNode(), *indep.getOperationNode())) {
794 equations[i]->addVariable(j);
800 map<size_t, Vnode<Base>*> tape2FreeVariables;
802 tape2FreeVariables[j->tapeIndex()] = j;
819 if (!j->isDeleted() && j->equations().size() == 1) {
821 if (i.assignmentVariable() ==
nullptr) {
823 jacSparsity, tape2FreeVariables,
824 equations, varInfo)) {
825 throw CGException(
"Failed to solve equation ", i.name(),
" for variable ", j->name());
831 }
while (assigned > 0);
839 if (!j->isDeleted() && j->equations().size() == 1) {
841 if (i.assignmentVariable() ==
nullptr) {
843 jacSparsity, tape2FreeVariables,
849 }
while (assigned > 0);
857 if (i->assignmentVariable() ==
nullptr && i->variables().size() == 1) {
860 jacSparsity, tape2FreeVariables,
866 }
while (assigned > 0);
875 if (!j->isDeleted()) {
877 if (i->assignmentVariable() ==
nullptr) {
879 jacSparsity, tape2FreeVariables,
880 equations, varInfo)) {
902 if (j->assignmentEquation() !=
nullptr) {
909 if (i->assignmentVariable() ==
nullptr) {
918 if (j->assignmentEquation() !=
nullptr) {
919 int i = j->assignmentEquation()->index();
922 if (eq.getAssignedVarIndex() !=
int(j->tapeIndex())) {
923 eq.setAssignedVarIndex(j->tapeIndex());
932 if (j->assignmentEquation() ==
nullptr) {
936 error +=
" " + j->name();
940 throw CGException(
"Failed to generate semi-explicit DAE. Could not solve system for the following variables:", error);
944 deleteVectorValues(diffVariables);
945 deleteVectorValues(dummyVariables);
946 deleteVectorValues(equations);
950 if (this->verbosity_ >= Verbosity::High) {
952 if (j->assignmentEquation() !=
nullptr)
953 log() <<
"## Variable " + j->name() <<
" assigned to equation " << j->assignmentEquation()->name() <<
"\n";
957 deleteVectorValues(diffVariables);
958 deleteVectorValues(dummyVariables);
959 deleteVectorValues(equations);
965 std::vector<bool>& jacSparsity,
966 const std::map<
size_t,
Vnode<Base>*>& tape2FreeVariables,
968 std::vector<DaeVarInfo>& varInfo) {
973 std::vector<bool> localJacSparsity = jacSparsity;
974 const size_t n = varInfo.size();
979 CGBase& dep = res0[i.index()];
980 CGBase& indep = indep0[j.tapeIndex()];
982 std::string indepName;
983 if (indep.getOperationNode()->getName() !=
nullptr) {
984 indepName = *indep.getOperationNode()->getName();
992 CPPADCG_ASSERT_UNKNOWN(alias !=
nullptr && alias->
getOperationType() == CGOpCode::Alias);
1000 throw CGException(
"Failed to solve equation ", i.name(),
" for variable ", j.name(),
": ", ex.what());
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);
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;
1019 affected.insert(equations[e]);
1020 for (
size_t p = 0; p < nnzs.size(); ++p) {
1021 localJacSparsity[n * e + nnzs[p]] =
true;
1026 if (affected.size() > 0) {
1030 map<size_t, set<Enode<Base>*> > solvable;
1031 for (
size_t e = 0; e < equations.size(); ++e) {
1033 if (&i != eq && affected.find(eq) == affected.end()) {
1035 for (
size_t v = 0; v < eq->variables().size(); ++v) {
1036 solvable[eq->variables()[v]->tapeIndex()].insert(eq);
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);
1056 for (
const auto& it : tape2FreeVariables) {
1061 if (v->assignmentEquation() !=
nullptr) {
1062 if (affected.count(v->assignmentEquation()) > 0 &&
1063 solvable[v->tapeIndex()].count(v->assignmentEquation()) == 0) {
1067 }
else if (solvable[v->tapeIndex()].size() == 0) {
1075 if (indepName.size() > 0) {
1076 indep.getOperationNode()->setName(indepName);
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);
1093 a.deleteNode(it.second);
1109 j.setAssignmentEquation(i, log(), this->verbosity_);
1110 j.deleteNode(log(), this->verbosity_);
1112 jacSparsity = localJacSparsity;
1118 const std::vector<DaeVarInfo>& varInfo,
1119 std::vector<DaeVarInfo>& newVarInfo,
1120 const std::vector<DaeEquationInfo>& eqInfo,
1121 std::vector<DaeEquationInfo>& newEqInfo) {
1123 using namespace std;
1131 std::set<size_t> oldVarWithDerivatives;
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());
1138 if (oldVarWithDerivatives.empty()) {
1140 for (
size_t j = 0; j < varInfo.size(); j++) {
1142 bool differential =
false;
1143 while (varInfo[index].getAntiDerivative() >= 0) {
1144 index = varInfo[index].getAntiDerivative();
1145 differential =
true;
1149 oldVarWithDerivatives.insert(index);
1157 std::vector<DaeVarOrderInfo> varOrder(varInfo.size());
1158 for (
size_t j = 0; j < varInfo.size(); j++) {
1160 int derivOrder = graph.determineVariableDiffOrder(varInfo, j, j0);
1161 if (varInfo[j].isIntegratedVariable()) {
1164 bool hasDerivatives = oldVarWithDerivatives.find(j) != oldVarWithDerivatives.end();
1168 std::sort(varOrder.begin(), varOrder.end(), sortVariablesByOrder);
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;
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]);
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]);
1200 std::vector<DaeEqOrderInfo> eqOrder(newEqInfo.size());
1201 for (
size_t i = 0; i < newEqInfo.size(); i++) {
1202 int assignedVar = newEqInfo[i].getAssignedVarIndex();
1204 while (newEqInfo[i0].getAntiDerivative() >= 0) {
1205 i0 = newEqInfo[i0].getAntiDerivative();
1207 bool isDifferential = newEqInfo[i].isExplicit() || (assignedVar >= 0 && newVarInfo[assignedVar].getAntiDerivative() >= 0);
1211 std::sort(eqOrder.begin(), eqOrder.end(), sortEquationByAssignedOrder2);
1213 std::vector<DaeEquationInfo> newEqInfo2(newEqInfo.size());
1214 for (
size_t i = 0; i < eqOrder.size(); i++) {
1215 newEqInfo2[i] = newEqInfo[eqOrder[i].originalIndex];
1217 newEqInfo = newEqInfo2;
1233 std::unique_ptr<ADFun<CGBase> > reorderedFun(
generateReorderedModel(handler, res0, varInfo, newVarInfo, eqInfo, newEqInfo));
1235 if (this->verbosity_ >= Verbosity::High) {
1236 log() <<
"reordered DAE equations and variables:\n";
1237 graph.printModel(log(), *reorderedFun, newVarInfo, newEqInfo);
1240 return reorderedFun;
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 {
1252 CPPADCG_ASSERT_UNKNOWN(indepNewOrder.size() == newVarInfo.size());
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];
1261 Independent(indepNewOrder);
1269 std::set<size_t> newIds;
1270 for (
size_t j = 0; j < newVarInfo.size(); j++) {
1271 newIds.insert(newVarInfo[j].getId());
1274 std::map<size_t, size_t> varId2HandlerIndex;
1275 size_t handlerIndex = 0;
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++;
1284 for (
size_t p = 0; p < newVarInfo.size(); p++) {
1285 size_t id = newVarInfo[p].getId();
1286 indepHandlerOrder[varId2HandlerIndex[id]] = indepNewOrder[p];
1290 std::map<size_t, size_t> eqId2OldIndex;
1291 for (
size_t i = 0; i < eqInfo.size(); i++) {
1292 eqId2OldIndex[eqInfo[i].getId()] = i;
1296 for (
size_t i = 0; i < newEqInfo.size(); i++) {
1297 size_t oldIndex = eqId2OldIndex[newEqInfo[i].getId()];
1298 resNewOrder[i] = res0[oldIndex];
1303 evaluator0.setPrintFor(
idxIdentify_->getGraph().isPreserveNames());
1304 vector<ADCG> depNewOrder = evaluator0.evaluate(indepHandlerOrder, resNewOrder);
1314 using namespace std;
1321 auto& vnodes = graph.variables();
1322 auto& enodes = graph.equations();
1327 row.reserve((vnodes.size() - diffVarStart_) * (m - diffEqStart_));
1328 col.reserve(row.capacity());
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();
1344 std::copy(
x_.begin(),
x_.end(), indep.begin());
1345 std::fill(indep.begin() +
x_.size(), indep.end(), 0);
1347 CppAD::sparse_jacobian_work work;
1349 row, col, jac, work);
1352 jacobian_.resize(m - diffEqStart_, vnodes.size() - diffVarStart_);
1354 map<size_t, Vnode<Base>*> origIndex2var;
1355 for (
size_t j = diffVarStart_; j< vnodes.size(); j++) {
1357 origIndex2var[jj->tapeIndex()] = jj;
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());
1366 Base normVal = jac[e].getValue() *
normVar_[vOrig->tapeIndex()]
1370 size_t j = origIndex2var[col[e]]->index();
1372 jacobian_.coeffRef(i - diffEqStart_, j - diffVarStart_) = normVal;
1377 if (this->verbosity_ >= Verbosity::High) {
1378 log() <<
"\npartial jacobian:\n" <<
jacobian_ <<
"\n\n";
1387 if (eqs.size() == vars.size()) {
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] <<
"; ";
1406 std::set<size_t> excludeCols;
1407 std::set<size_t> avoidCols;
1408 for (
size_t j = 0; j < vars.size(); j++) {
1410 bool notZero =
false;
1411 for (
size_t i = 0; i < eqs.size(); i++) {
1413 Base val =
jacobian_.coeff(ii->index() - diffEqStart_, jj->index() - diffVarStart_);
1414 if (val != Base(0.0)) {
1421 excludeCols.insert(j);
1424 avoidCols.insert(j);
1428 if (eqs.size() <= vars.size() - (excludeCols.size() + avoidCols.size())) {
1430 excludeCols.insert(avoidCols.begin(), avoidCols.end());
1432 if (this->verbosity_ >= Verbosity::Low)
1433 log() <<
"Must use variables defined to be avoided by the user!\n";
1436 std::vector<Vnode<Base>* > varsLocal;
1438 Eigen::ColPivHouseholderQR<MatrixB> qr;
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]);
1449 work.setZero(eqs.size(), varsLocal.size());
1452 for (
size_t i = 0; i < eqs.size(); i++) {
1454 for (
size_t j = 0; j < varsLocal.size(); j++) {
1456 Base val =
jacobian_.coeff(ii->index() - diffEqStart_, jj->index() - diffVarStart_);
1457 if (val != Base(0.0)) {
1463 if (this->verbosity_ >= Verbosity::High)
1464 log() <<
"subset Jac:\n" << work <<
"\n";
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.");
1476 using PermutationMatrix =
typename Eigen::ColPivHouseholderQR<MatrixB>::PermutationType;
1477 using Indices =
typename PermutationMatrix::IndicesType;
1479 const PermutationMatrix& p = qr.colsPermutation();
1480 const Indices& indices = p.indices();
1482 if (this->verbosity_ >= Verbosity::High) {
1483 log() <<
"## matrix Q:\n";
1484 MatrixB q = qr.matrixQ();
1486 log() <<
"## matrix R:\n";
1487 MatrixB r = qr.matrixR().template triangularView<Eigen::Upper>();
1489 log() <<
"## matrix P: " << indices.transpose() <<
"\n";
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.");
1501 if (avoidCols.empty()) {
1505 if (this->verbosity_ >= Verbosity::Low)
1506 log() <<
"Failed to determine dummy derivatives without the variables defined to be avoided by the user\n";
1509 for (
size_t j : avoidCols)
1510 excludeCols.erase(j);
1517 const auto& p = qr.colsPermutation();
1518 const auto& indices = p.indices();
1520 std::vector<Vnode<Base>* > newDummies;
1521 if (avoidConvertAlg2DifVars_) {
1523 const auto& varInfo = graph.getOriginalVariableInfo();
1526 for (
int i = 0; newDummies.size() <
size_t(work.rows()) && i < qr.rank(); 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) {
1533 newDummies.push_back(v);
1537 for (
int i = 0; newDummies.size() <
size_t(work.rows()); 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) {
1544 newDummies.push_back(v);
1550 for (
int i = 0; i < work.rows(); i++) {
1551 newDummies.push_back(varsLocal[indices(i)]);
1555 if (this->verbosity_ >= Verbosity::High) {
1556 log() <<
"## new dummy derivatives: ";
1558 log() << *it <<
"; ";
1567 dummyD_.insert(
dummyD_.end(), newDummies.begin(), newDummies.end());
1570 inline static void printModel(std::ostream& out,
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) {
1578 std::vector<std::string> indepNames;
1579 for (
size_t p = 0; p < varInfo.size(); p++) {
1580 if (erasedVariables.find(p) == erasedVariables.end()) {
1582 indepNames.push_back(varInfo[p].getName());
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]);
1594 std::vector<std::string> depNames;
1595 LangCCustomVariableNameGenerator<Base> nameGen(depNames, indepNames);
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,
1604 for (
size_t e = 0; e < equations.size(); ++e) {
1605 Enode<Base>* eq = equations[e];
1607 for (
const auto& it : tape2FreeVariables) {
1608 if (jacSparsity[n * eq->index() + it.first]) {
1610 out <<
"# Equation " << e <<
": \t";
1611 out <<
" " << it.second->name();
1623 inline static void deleteVectorValues(std::vector<T*>& v) {
1624 for (
size_t i = 0; i < v.size(); i++) {