void EvaluationMX::evaluateMX(const MXPtrV& arg, MXPtrV& res, const MXPtrVV& fseed, MXPtrVV& fsens, const MXPtrVV& aseed, MXPtrVV& asens, bool output_given) { // Number of sensitivity directions int nfdir = fsens.size(); casadi_assert(nfdir==0 || fcn_.spCanEvaluate(true)); int nadir = aseed.size(); casadi_assert(nadir==0 || fcn_.spCanEvaluate(false)); // Get/generate the derivative function FX d = fcn_.derivative(nfdir, nadir); // Temporary vector<MX> tmp; // Assemble inputs vector<MX> d_arg; d_arg.reserve(d.getNumInputs()); // Nondifferentiated inputs tmp = getVector(arg); d_arg.insert(d_arg.end(), tmp.begin(), tmp.end()); for (MXPtrVV::const_iterator i = fseed.begin(); i != fseed.end(); ++i) { tmp = getVector(*i); d_arg.insert(d_arg.end(), tmp.begin(), tmp.end()); } for (MXPtrVV::const_iterator i = aseed.begin(); i != aseed.end(); ++i) { tmp = getVector(*i); d_arg.insert(d_arg.end(), tmp.begin(), tmp.end()); } // Evaluate symbolically vector<MX> d_res = d.call(d_arg); vector<MX>::const_iterator d_res_it = d_res.begin(); // Collect the nondifferentiated results for (MXPtrV::iterator i = res.begin(); i != res.end(); ++i, ++d_res_it) { if (!output_given && *i) **i = *d_res_it; } // Collect the forward sensitivities for (MXPtrVV::iterator j = fsens.begin(); j != fsens.end(); ++j) { for (MXPtrV::iterator i = j->begin(); i != j->end(); ++i, ++d_res_it) { if (*i) **i = *d_res_it; } } // Collect the adjoint sensitivities for (MXPtrVV::iterator j = asens.begin(); j != asens.end(); ++j) { for (MXPtrV::iterator i = j->begin(); i != j->end(); ++i, ++d_res_it) { if(*i && !d_res_it->isNull()){ **i += *d_res_it; } } } // Make sure that we've got to the end of the outputs casadi_assert(d_res_it==d_res.end()); }
EvaluationMX::EvaluationMX(const FX& fcn, const std::vector<MX> &arg) : fcn_(fcn) { // Number inputs and outputs int num_in = fcn.getNumInputs(); // All dependencies of the function vector<MX> d = arg; d.resize(num_in); setDependencies(d); setSparsity(CRSSparsity(1, 1, true)); }
MXFunction vec (const FX &a_) { FX a = a_; // Pass null if input is null if (a.isNull()) return MXFunction(); // Get the MX inputs, only used for shape const std::vector<MX> &symbolicInputMX = a.symbolicInput(); // Have a vector with MX that have the shape of vec(symbolicInputMX ) std::vector<MX> symbolicInputMX_vec(a.getNumInputs()); // Make vector valued MX's out of them std::vector<MX> symbolicInputMX_vec_reshape(a.getNumInputs()); // Apply the vec-transformation to the inputs for (int i=0;i<symbolicInputMX.size();++i) { std::stringstream s; s << "X_flat_" << i; symbolicInputMX_vec[i] = MX(s.str(),vec(symbolicInputMX[i].sparsity())); symbolicInputMX_vec_reshape[i] = trans(reshape(symbolicInputMX_vec[i],trans(symbolicInputMX[i].sparsity()))); } // Call the original function with the vecced inputs std::vector<MX> symbolicOutputMX = a.call(symbolicInputMX_vec_reshape); // Apply the vec-transformation to the outputs for (int i=0;i<symbolicOutputMX.size();++i) symbolicOutputMX[i] = vec(symbolicOutputMX[i]); // Make a new function with the vecced input/outputs MXFunction ret(symbolicInputMX_vec,symbolicOutputMX); // Initialize it if a was if (a.isInit()) ret.init(); return ret; }
void EvaluationMX::create(const FX& fcn, const std::vector<MX> &arg, std::vector<MX> &res, const std::vector<std::vector<MX> > &fseed, std::vector<std::vector<MX> > &fsens, const std::vector<std::vector<MX> > &aseed, std::vector<std::vector<MX> > &asens, bool output_given) { // Number inputs and outputs int num_in = fcn.getNumInputs(); int num_out = fcn.getNumOutputs(); // Number of directional derivatives int nfdir = fseed.size(); int nadir = aseed.size(); // Create the evaluation node MX ev; if(nfdir>0 || nadir>0){ // Create derivative function Derivative dfcn(fcn,nfdir,nadir); stringstream ss; ss << "der_" << fcn.getOption("name") << "_" << nfdir << "_" << nadir; dfcn.setOption("verbose",fcn.getOption("verbose")); dfcn.setOption("name",ss.str()); dfcn.init(); // All inputs vector<MX> darg; darg.reserve(num_in*(1+nfdir) + num_out*nadir); darg.insert(darg.end(),arg.begin(),arg.end()); // Forward seeds for(int dir=0; dir<nfdir; ++dir){ darg.insert(darg.end(),fseed[dir].begin(),fseed[dir].end()); } // Adjoint seeds for(int dir=0; dir<nadir; ++dir){ darg.insert(darg.end(),aseed[dir].begin(),aseed[dir].end()); } ev.assignNode(new EvaluationMX(dfcn, darg)); } else { ev.assignNode(new EvaluationMX(fcn, arg)); } // Output index int ind = 0; // Create the output nodes corresponding to the nondifferented function res.resize(num_out); for (int i = 0; i < num_out; ++i, ++ind) { if(!output_given){ if(!fcn.output(i).empty()){ res[i].assignNode(new OutputNode(ev, ind)); } else { res[i] = MX(); } } } // Forward sensitivities fsens.resize(nfdir); for(int dir = 0; dir < nfdir; ++dir){ fsens[dir].resize(num_out); for (int i = 0; i < num_out; ++i, ++ind) { if (!fcn.output(i).empty()){ fsens[dir][i].assignNode(new OutputNode(ev, ind)); } else { fsens[dir][i] = MX(); } } } // Adjoint sensitivities asens.resize(nadir); for (int dir = 0; dir < nadir; ++dir) { asens[dir].resize(num_in); for (int i = 0; i < num_in; ++i, ++ind) { if (!fcn.input(i).empty()) { asens[dir][i].assignNode(new OutputNode(ev, ind)); } else { asens[dir][i] = MX(); } } } }