/* This file is part of Lemma, a geophysical modelling and inversion API */ /* This Source Code Form is subject to the terms of the Mozilla Public * License, v. 2.0. If a copy of the MPL was not distributed with this * file, You can obtain one at http://mozilla.org/MPL/2.0/. */ /** @file @author Trevor Irons @date 10/11/2010 @version $Id: logbarriercg.h 265 2015-03-27 16:05:21Z tirons $ **/ #ifndef LOGBARRIERCG_INC #define LOGBARRIERCG_INC #include #include #include #include #include #include "cg.h" #include "bicgstab.h" #include "lemma.h" #ifdef LEMMAUSEVTK #include "matplot.h" #endif namespace Lemma { template < typename cgScalar > cgScalar PhiB (const cgScalar& mux, const cgScalar& muy, const cgScalar& minVal, const cgScalar& maxVal, const VectorXr x) { cgScalar phib = std::abs((x.array() - minVal).log().sum()*mux); // phib += std::abs((maxVal - x.array()).log().sum()*muy); return phib; } // PhiB for block log barrier template < typename cgScalar > cgScalar PhiB2 (const cgScalar& mux, const cgScalar& muy, const cgScalar& minVal, const cgScalar& maxVal, const VectorXr x, const int& block, const int &nblocks) { cgScalar phib = std::abs((x.array() - minVal).log().sum()*mux); //phib += std::abs((maxVal - x.array()).log().sum()*muy); for (int ib=0; ib cgScalar PhiB2 (const cgScalar& minVal, const cgScalar& maxVal, const VectorXr x, const int& block, const int &nblocks) { cgScalar phib = std::abs((x.array() - minVal).log().sum()); //phib += std::abs((maxVal - x.array()).log().sum()*muy); for (int ib=0; ib cgScalar PhiB2_NN (const cgScalar& mux, const cgScalar& minVal, const VectorXr x) { cgScalar phib = std::abs((x.array() - minVal).log().sum()*mux); return phib; } /** Impliments a logarithmic barrier CG solution of a Real linear system of * the form \f[ \mathbf{A} \mathbf{x} = \mathbf{b} \f] s.t. \f$ x \in * (minVal, maxVal) \f$. Note that this method optimized the complete * solution, using the large matrix ATA. If you have a system with a huge * number of columns, see the implicit version of this routine. Solution of * the dual problem (interior-point) follows "Tikhonov Regularization with * Nonnegativity constraint, (Calavetti et. al. 2004)" * @param[in] A is a Real Matrix. * @param[in] b is a Real vector * @param[in] x0 is current value. * @param[in] xr is reference model * @param[in] Wd is a Sparse Real matrix, that specifies data objective * function. * @param[in] Wm is a Sparse Real matrix, that sepcifies the model * objective function. * @param[in] minVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ * @param[in] maxVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ */ template < typename cgScalar > VectorXr LogBarrierCG_4(const MatrixXr &A, const VectorXr &b, const VectorXr &xr, const VectorXr &x0, const Eigen::SparseMatrix& WdTWd, const Eigen::SparseMatrix& WmTWm, const cgScalar &minVal, const cgScalar &maxVal, const Real& lambda) { // TODO pass correct phid and phim into this so that we are working with that. // Check that everything is aligned OK if (A.rows() != b.size() ) { std::cerr << "Matrix A is not aligned with Vector b" << "\n"; exit(1); } // std::cout << "A\n" << A << std::endl; // std::cout << "WdTWd\n" << WdTWd << std::endl; // std::cout << "WmTWm\n" << WmTWm << std::endl; MatrixXr ATWdTWdA = A.transpose()*WdTWd*A; ATWdTWdA += lambda*WmTWm; int N = A.cols(); // number of model int M = A.rows(); // number of data cgScalar SIGMA = .25; //1e-2; // 1e-1; cgScalar epsilon = 1e-55; // Ratio between phib and phim+phid Real delta = minVal + 1e-8; // initial guess, start with reference model VectorXr x = VectorXr::Zero(N).array() + epsilon; VectorXr b_pre = VectorXr::Zero(M); cgScalar phid = (b_pre - b).norm(); cgScalar phim = x.norm(); cgScalar mux = (phid + lambda*phim) * 1e4; cgScalar muy = 1e-30; cgScalar phib = PhiB(mux, muy, minVal, maxVal, x); int iloop(0); do { //////////////////////////////////////////////////////////// // Log barrier terms VectorXr X1 = VectorXr::Ones(N).array() / ((x0+x).array()-minVal) ; VectorXr Y1 = VectorXr::Ones(N).array() / (maxVal-(x0+x).array()) ; VectorXr X2 = VectorXr::Ones(N).array() / (((x0+x).array()-minVal)*((x0+x).array()-minVal)); VectorXr Y2 = VectorXr::Ones(N).array() / ((maxVal-(x0+x).array())*(maxVal-(x0+x).array())); ///////////////////////////////////////////////////////////// // Solve system ///////////////////////////////////////////////////////////// // CG solution of complete system VectorXr b2 = (A.transpose()*WdTWd*(b)).array() // b_pre-b_obs //- (lambda*WmTWm*(x-xr)).array() + 2.*mux*X1.array() + 2.*muy*Y1.array(); // Eigen requires these temporaries :( MatrixXr A2 = ATWdTWdA; A2 += lambda*WmTWm; A2.diagonal().array() += mux*X2.array() + muy*Y2.array(); //std::cout << "all set up" << std::endl; Eigen::ConjugateGradient CG; //CG.setTolerance(1e-30); CG.compute(A2); VectorXr ztilde = CG.solve(b2); // full soln // update x VectorXr h = ztilde - x; cgScalar d(.25); for (int ix=0; ix maxVal) { std::cout << "overstepp BIG\t" << x(i) << "\t" << x0(i) << "\t" << maxVal<< std::endl; x(i) = maxVal - x0(i) - delta; exit(1); } } //std::cout << "|h|=" << h.norm() << " |x|=" << x.norm() << std::endl; // Determine mu steps to take VectorXr s1 = mux * (X2.asDiagonal()*h - 2.*X1); VectorXr s2 = muy * (Y2.asDiagonal()*h - 2.*Y1); //std::cout << "s1 = " << s1.transpose() << "\ns2 = " << s2.transpose() << std::endl; //std::cout << "ztilde = " << ztilde.transpose() << std::endl; //std::cout << "x = " << x.transpose() << std::endl; //std::cout << "h = " << h.transpose() << std::endl; // determine mu for next step mux = SIGMA/((cgScalar)(N)) * std::abs( s1.dot(x) ) ; muy = SIGMA/((cgScalar)(N)) * std::abs( s2.dot(x) ) ; b_pre = A*x; phid = (b_pre - b).norm(); //std::sqrt( (WdTWd*(b_pre-b)).norm() ) ; //(b - (b_pre)).norm(); //Real phid_n = (b_pre-b).norm() ; //(b - (b_pre)).norm(); phim = (WmTWm*(x-xr)).norm(); phib = PhiB(mux, muy, minVal, maxVal, x0+x); ++iloop; std::cout.precision(14); std::cout << iloop << "\t" << phid << "\t" << phim << "\t" << mux << "\t" << muy << "\t" << phib<< std::endl; //} while (std::abs(phib) > epsilon && phid > 1e2); } while (std::abs(phib / (phid+lambda*phim)) > epsilon); std::cout << "final x\t" << x.array().exp().transpose() << std::endl; return x; } /** Impliments a logarithmic barrier CG solution of a Real linear system of * the form \f[ \mathbf{A} \mathbf{x} = \mathbf{b} \f] s.t. \f$ x \in * (minVal, maxVal) \f$. Note that this method optimized the complete * solution, using the large matrix ATA. If you have a system with a huge * number of columns, see the implicit version of this routine. Solution of * the dual problem (interior-point) follows "Tikhonov Regularization with * Nonnegativity constraint, (Calavetti et. al. 2004)" * @param[in] A is a Real Matrix. * @param[in] b is a Real vector * @param[in] minVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ * @param[in] maxVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ */ template < typename cgScalar > VectorXr LogBarrierCG(const MatrixXr &A, const VectorXr &b, const cgScalar &minVal, const cgScalar &maxVal) { // Check that everything is aligned OK if (A.rows() != b.size() ) { std::cerr << "Matrix A is not aligned with Vector b" << "\n"; exit(1); } // TODO make ATA implicit MatrixXr ATA = A.transpose()*A; int N = A.cols(); // number of model int M = A.rows(); // number of data int MAXITER = 100; // M*N; cgScalar SIGMA = 1e-2; //1e-1; cgScalar delta = 1e-4; // Cholesky preconditioner //Eigen::FullPivLU Pre; //Eigen::ColPivHouseholderQR Pre; //Pre.compute(ATA); // Determine starting lambda_0, find lim_sup of the norm of impulse responses and scale cgScalar limsup = 1e10; for (int i=0; i MM = // Eigen::SparseMatrix(A2.rows(), A2.cols()); // for (int i=0; i maxVal) { // x.segment(ib*block, block).array() *= (.9); // } // } // for (int i=0; i epsilon); // report phireport.precision(12); std::cout << ii << "\t" << x.norm() << "\t" << (b-(A*x)).norm() << "\t" << lambda << "\t" << mux << "\t" << muy << "\t" << iter << "\t" << iloop << std::endl; phireport << ii << "\t" << x.norm() << "\t" << (b-(A*x)).norm() << "\t" << lambda << "\t" << mux << "\t" << muy << "\t" << iter << "\t" << iloop << std::endl; // write model file std::fstream modfile; std::string fname = "iteration" + to_string(ii) + ".dat"; modfile.open( fname.c_str(), std::ios::out); modfile << ii << "\t" << x.norm() << "\t" << (b-(A*x)).norm() << "\t" << lambda << "\t" << mux << "\t" << muy << "\t" << iter << "\n"; modfile << x << "\n"; modfile.close(); // write predicted data file std::fstream predata; fname = "iteration" + to_string(ii) + "pre.dat"; predata.open(fname.c_str(), std::ios::out); predata << b_pre << std::endl; predata.close(); // update lambda // @todo smarter lambda change lambda *= .92; } phireport.close(); // TODO, determine optimal solution return x; } /** Impliments a logarithmic barrier CG solution of a Real linear system of * the form \f[ \mathbf{A} \mathbf{x} = \mathbf{b} \f] s.t. \f$ x \in * (minVal, maxVal) \f$. Note that this method optimized the complete * solution, using the large matrix ATA. If you have a system with a huge * number of columns, see the implicit version of this routine. Solution of * the dual problem (interior-point) follows "Tikhonov Regularization with * Nonnegativity constraint, (Calavetti et. al. 2004)" * @param[in] A is a Real Matrix. * @param[in] b is a Real vector * @param[in] minVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ * @param[in] maxVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ * @param[in] block is the number of parameters to sum together for the * upper barrier term. So block block number parameters are kept under maxVal. * as such A.rows() / block must be evenly divisible. */ template VectorXr LogBarrierCG(const MatrixXr &A, const VectorXr &b, const cgScalar &minVal, const cgScalar &maxVal, const int& block) { // Check that everything is aligned OK if (A.rows() != b.size() ) { std::cerr << "Matrix A is not aligned with Vector b" << "\n"; exit(1); } // write predicted data file std::fstream obsdata; std::string fname = "obsdata.dat"; obsdata.open(fname.c_str(), std::ios::out); obsdata << b << std::endl; obsdata.close(); // #ifdef LEMMAUSEVTK // double blue[3] = {0.0,0.0,1.0}; // double red[3] = {1.0,0.0,0.0}; // VectorXr ind = VectorXr::LinSpaced(b.size(), 1, b.size()); // #endif // TODO make ATA implicit MatrixXr ATA = A.transpose()*A; int N = A.cols(); // number of model int M = A.rows(); // number of data int MAXITER = 100; // M*N; cgScalar SIGMA = 1e-2; //1e-1; cgScalar delta = 1e-4; // Cholesky preconditioner //Eigen::FullPivLU Pre; //Eigen::ColPivHouseholderQR Pre; //Pre.compute(ATA); // Determine starting lambda_0, find lim_sup of the norm of impulse responses and scale cgScalar limsup = 1e10; for (int i=0; i MM = // Eigen::SparseMatrix(A2.rows(), A2.cols()); // for (int i=0; i maxVal) { // x.segment(ib*block, block).array() *= (.9); // } // } // for (int i=0; i epsilon); // report phireport.precision(12); std::cout << ii << "\t" << x.norm() << "\t" << (b-(A*x)).norm() << "\t" << lambda << "\t" << mux << "\t" << muy << "\t" << iter << "\t" << iloop << std::endl; phireport << ii << "\t" << x.norm() << "\t" << (b-(A*x)).norm() << "\t" << lambda << "\t" << mux << "\t" << muy << "\t" << iter << "\t" << iloop << std::endl; std::fstream modfile; std::string fname = "iteration" + to_string(ii) + ".dat"; modfile.open( fname.c_str(), std::ios::out); modfile << ii << "\t" << x.norm() << "\t" << (b-(A*x)).norm() << "\t" << lambda << "\t" << mux << "\t" << muy << "\t" << iter << "\n"; modfile << x << "\n"; modfile.close(); // write predicted data file std::fstream predata; fname = "iteration" + to_string(ii) + "pre.dat"; predata.open(fname.c_str(), std::ios::out); predata << b_pre << std::endl; predata.close(); // update lambda // @todo smarter lambda change lambda *= .85; } phireport.close(); // TODO, determine optimal solution return x; } /** Impliments a logarithmic barrier CG solution of a Real linear system of * the form \f[ \mathbf{A} \mathbf{x} = \mathbf{b} \f] s.t. \f$ x \in * (minVal, maxVal) \f$. Note that this method optimized the complete * solution, using the large matrix ATA. If you have a system with a huge * number of columns, see the implicit version of this routine. Solution of * the dual problem (interior-point) follows "Tikhonov Regularization with * Nonnegativity constraint, (Calavetti et. al. 2004)" * @param[in] A is a Real Matrix. * @param[in] xref is a reference model * @param[in] b is a Real vector * @param[in] minVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ * @param[in] maxVal is the minimum allowed value \f$ x \in (minVal, maxVal) \f$ * @param[in] block is the number of parameters to sum together for the * upper barrier term. So block block number parameters are kept under maxVal. * as such A.rows() / block must be evenly divisible. * @param[in] WdTWd is the data objective function * @param[in] WmTWm is the model objective function */ template VectorXr LogBarrierCG(const MatrixXr &A, const VectorXr &xr, const VectorXr &b, const cgScalar &minVal, const cgScalar &maxVal, const int& block, const Eigen::SparseMatrix& WdTWd, const Eigen::SparseMatrix& WmTWm, Real lambda0=1e1) { // Check that everything is aligned OK if (A.rows() != b.size() ) { std::cerr << "Matrix A is not aligned with Vector b" << "\n"; std::cerr << "A.rows() " << A.rows() << "\n"; std::cerr << "A.cols() " << A.cols() << "\n"; std::cerr << "b.size() " << b.size() << "\n"; exit(1); } // write predicted data file std::fstream obsdata; std::string fname = "obsdata.dat"; obsdata.open(fname.c_str(), std::ios::out); obsdata << b << std::endl; obsdata.close(); // TODO make ATA implicit MatrixXr ATWdTWdA = A.transpose()*WdTWd*A; int N = A.cols(); // number of model int M = A.rows(); // number of data int MAXITER = 175; // M*N; cgScalar SIGMA = .25;//5.85; //1e-2; // .25; //1e-2; // 1e-1; cgScalar delta = 1e-4; // // Determine starting lambda_0, find lim_sup of the norm of impulse responses and scale // cgScalar limsup = 1e10; // for (int i=0; i cg; /// @todo add stopping criteria. for (int ii=0; ii MM = // Eigen::SparseMatrix(A2.rows(), A2.cols()); // for (int i=0; i= maxVal) { x.segment(ib*block, block).array() *= .99; } } for (int i=0; i VectorXr LogBarrierCG_NN(const MatrixXr &A, const VectorXr &xr, const VectorXr &b, const cgScalar &minVal, const Eigen::SparseMatrix& WdTWd, const Eigen::SparseMatrix& WmTWm, Real lambda0=1e1, int MAXITER=175) { // Check that everything is aligned OK if (A.rows() != b.size() ) { std::cerr << "Matrix A is not aligned with Vector b" << "\n"; std::cerr << "A.rows() " << A.rows() << "\n"; std::cerr << "A.cols() " << A.cols() << "\n"; std::cerr << "b.size() " << b.size() << "\n"; exit(1); } // write predicted data file std::fstream obsdata; std::string fname = "obsdata.dat"; obsdata.open(fname.c_str(), std::ios::out); obsdata << b << std::endl; obsdata.close(); // TODO make ATA implicit, or at least only compute half MatrixXr ATWdTWdA = A.transpose()*WdTWd*A; int N = A.cols(); // number of model int M = A.rows(); // number of data //int MAXITER = 175; // M*N; cgScalar SIGMA = .25;//5.85; //1e-2; // .25; //1e-2; // 1e-1; cgScalar delta = 1e-4; // // Determine starting lambda_0, find lim_sup of the norm of impulse responses and scale // cgScalar limsup = 1e10; // for (int i=0; i cg; // Use ILUT preconditioner Eigen::ConjugateGradient< MatrixXr, Eigen::Upper, Eigen::DiagonalPreconditioner > cg; //Eigen::ConjugateGradient< MatrixXr, Eigen::Upper, Eigen::IncompleteLUT > cg; cgScalar phid = (WdTWd*(b - b_pre)).norm(); cgScalar phim = (WmTWm*(x - xr)).norm(); cgScalar phib = PhiB2_NN(1., minVal, x); cgScalar mux = (phid + lambda*phim) / phib; //cgScalar tol; // = 1e-5*phid; // 1e-14; std::fstream phireport("phimphid.dat", std::ios::out); std::cout << "Starting CG iterations 4" << std::endl; /// @todo add stopping criteria. for (int ii=0; ii