/* 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 88 2013-09-06 17:24:44Z tirons $ **/ #ifndef LOGBARRIERCG_INC #define LOGBARRIERCG_INC #include #include #include #include #include #include "cg.h" #include "bicgstab.h" #include "lemma.h" namespace Lemma { template < typename Scalar > Scalar PhiB (const Scalar& mux, const Scalar& muy, const Scalar& minVal, const Scalar& maxVal, const VectorXr x) { Scalar 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 Scalar > Scalar PhiB2 (const Scalar& mux, const Scalar& muy, const Scalar& minVal, const Scalar& maxVal, const VectorXr x, const int& block, const int &nblocks) { Scalar phib = std::abs((x.array() - minVal).log().sum()*mux); //phib += std::abs((maxVal - x.array()).log().sum()*muy); for (int ib=0; ib Scalar PhiB2 (const Scalar& minVal, const Scalar& maxVal, const VectorXr x, const int& block, const int &nblocks) { Scalar phib = std::abs((x.array() - minVal).log().sum()); //phib += std::abs((maxVal - x.array()).log().sum()*muy); for (int ib=0; ib Scalar PhiB2_NN (const Scalar& mux, const Scalar& minVal, const VectorXr x) { Scalar 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 a 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 Scalar > VectorXr LogBarrierCG(const MatrixXr &A, const VectorXr &b, const VectorXr &x0, const Eigen::SparseMatrix& WdTWd, const Eigen::SparseMatrix& WmTWm, const Scalar &minVal, const Scalar &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); } ///////////////////////////////////////////// // Build all the large static matrices // NOTE, ATA can be large. For some problems an implicit algorithm may // be better suited. //Eigen::SparseMatrix WmTWm = Wm.transpose()*Wm; //Eigen::SparseMatrix WdTWd = Wd.transpose()*Wd; MatrixXr ATWdTWdA = A.transpose()*WdTWd*A; ///////////////////////// // Jacobi Preconditioner Eigen::SparseMatrix MM (ATWdTWdA.rows(), ATWdTWdA.cols()); 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(); // update lambda /// @todo smarter lambda change lambda *= .9; } 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$ */ template < typename Scalar > VectorXr LogBarrierCG(const MatrixXr &A, const VectorXr &b, const Scalar &minVal, const Scalar &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; Scalar SIGMA = 1e-2; //1e-1; Scalar 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 Scalar 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 Scalar &minVal, const Scalar &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; Scalar SIGMA = 1e-2; //1e-1; Scalar 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 Scalar 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 Scalar &minVal, const Scalar &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; Scalar SIGMA = .25;//5.85; //1e-2; // .25; //1e-2; // 1e-1; Scalar delta = 1e-4; // // Determine starting lambda_0, find lim_sup of the norm of impulse responses and scale // Scalar 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 Scalar &minVal, 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, 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; Scalar SIGMA = .25;//5.85; //1e-2; // .25; //1e-2; // 1e-1; Scalar delta = 1e-12; // // Determine starting lambda_0, find lim_sup of the norm of impulse responses and scale // Scalar 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; Scalar phid = (WdTWd*(b - b_pre)).norm(); Scalar phim = (WmTWm*(x - xr)).norm(); Scalar phib = PhiB2_NN(1., minVal, x); Scalar mux = (phid + lambda*phim) / phib; //Scalar tol; // = 1e-5*phid; // 1e-14; std::fstream phireport("phimphid.dat", std::ios::out); mux = (phid + lambda*phim) / phib; /// @todo add stopping criteria. for (int ii=0; ii VectorXr Tikhonov_CG(const MatrixXr &A, const VectorXr &xr, const VectorXr &b, 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.array() << 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 dat int MAXITER = 175; // M*N; //Scalar SIGMA = .25;//5.85; //1e-2; // .25; //1e-2; // 1e-1; Scalar delta = 1e-4; Scalar lambda = lambda0; // * limsup;//e-1;//limsup; Scalar epsilon = 1e-15; // Ratio between phib and phim+phid // initial guess, just near zero for now VectorXr x = VectorXr::Zero(N).array() + delta;// ((minVal + maxVal) / 2.); // predicted b VectorXr b_pre = A*x; Scalar phid = (WdTWd*(b - b_pre)).norm(); Scalar phim = (WmTWm*(x - xr)).norm(); Scalar tol = 1e-5*phid; // 1e-14; std::fstream phireport("phimphid.dat", std::ios::out); //ConjugateGradient > cg; // can't use Eigen, b/c we have mixed types of Re and complex Eigen::ConjugateGradient< MatrixXr > cg; //Eigen::BiCGSTAB< MatrixXcr > solver; /// @todo add stopping criteria. for (int ii=0; ii MM = // Eigen::SparseMatrix(A2.rows(), A2.cols()); // for (int i=0; i())).real(); //iter = cg.iterations(); //VectorXr ztilde = CG(A2, x, b2, iter, tol); //cg.setTolerance(); cg.compute(A2); x = cg.solveWithGuess(b2, x); b_pre = A*x; phid = (WdTWd*(b - b_pre)).norm(); phim = (WmTWm*(x - xr)).norm(); //std::cout << "out cg" << std::endl; ++iloop; itertot += iter; ///////////////////////////////////////////////////////////// // update x, mux, muy //VectorXr h = ztilde; // - x; // update x //h = ztilde - x; // determing steplength //Scalar d = std::min(1., 0.925*(x.array()/h.array().abs()).minCoeff() ); // Scalar d(1.); // for (int ix=0; ix= maxVal) { // x.segment(ib*block, block).array() *= .99; // } // } /* for (int i=0; i VectorXr LogBarrierCGLS(const MatrixXr &A, const VectorXr &b, const Scalar &minVal, const Scalar &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; Scalar SIGMA = .25; // 1e-2; //1e-1; Scalar delta = 1e-8; // 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 Scalar 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; } } #endif // ----- #ifndef LOGBARRIERCG_INC -----