diff --git a/src/dfpm/solver/parabolic_driver.hpp b/src/dfpm/solver/parabolic_driver.hpp index 0d706fc..ab899fd 100644 --- a/src/dfpm/solver/parabolic_driver.hpp +++ b/src/dfpm/solver/parabolic_driver.hpp @@ -1,187 +1,198 @@ /* ============================================================================= Copyright (c) 2014 - 2016 F. Georget Princeton University All rights reserved. Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met: 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer. 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution. 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission. * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. * ============================================================================= */ #ifndef SPECMICP_DFPMSOLVER_PARABOLICDRIVER_HPP #define SPECMICP_DFPMSOLVER_PARABOLICDRIVER_HPP #include "specmicp_common/types.hpp" #include "specmicp_common/sparse_solvers/sparse_solver.hpp" #include "driver.hpp" #include "parabolic_structs.hpp" namespace specmicp { namespace dfpmsolver { //! \brief The parabolic driver //! //! Parabolic driver for finite element programs template class ParabolicDriver: public Driver { public: using base=Driver; using base::program; using base::get_neq; using base::get_options; using base::get_perfs; using base::scaling; using base::initialize_scaling; ParabolicDriver(Program& the_program): Driver(the_program), m_velocity(Vector::Zero(the_program.get_tot_ndf())), m_solver(nullptr) {} //! \brief Solve a timestep of length dt ParabolicDriverReturnCode solve_timestep(scalar_t dt, Vector& displacement); //! \brief Restart the current timestep ParabolicDriverReturnCode restart_timestep(Vector& displacement); //! \brief Check if the solution has been found ParabolicDriverReturnCode check_convergence(); //! \brief Compute the residuals void compute_residuals(const Vector& displacement, Vector& residual) { compute_residuals(displacement, m_velocity, residual); } //! \brief Compute the residuals void compute_residuals(const Vector& displacement, const Vector& velocity, Vector& residual) { program().compute_residuals(displacement, velocity, residual); get_perfs().nb_call_residuals += 1; } + //! \brief Compute the residuals for use to check convergence + //! + //! For most program, it is the same as calling compute_residuals + void compute_residuals_0(const Vector& displacement, + const Vector& velocity, + Vector& residual) + { + program().compute_residuals_0(displacement, velocity, residual); + get_perfs().nb_call_residuals += 1; + } + //! \brief Compute the jacobian void compute_jacobian(Vector& displacement, Vector& velocity, Eigen::SparseMatrix& jacobian ); //! \brief Return the norm of the current residuals double norm_residuals() {return m_residuals.norm();} //! \brief Read/Write reference to the velocity vector const Vector& get_velocity() const {return m_velocity;} Vector& velocity() {return m_velocity;} void set_velocity(Vector& velocity_vector); //! \brief Reset the Sparse solver //! This function should be called when a new solver need to be used void reset_solver() {return m_solver.reset(nullptr);} //! \brief Call this function if the pattern of the jacobian has changed void jacobian_pattern_has_changed() {reset_solver();} //! \brief Initialize the computation void initialize_timestep(scalar_t dt, Eigen::VectorXd& displacement); //! \brief Strang Linesearch //! //! ref : //! - Matthies et al. (1979) //! - JHP course notes ParabolicLinesearchReturnCode strang_linesearch( Vector& update, Vector& displacements, scalar_t& lambda ); scalar_t residuals_0() {return m_norm_0;} private: void reset_velocity(); //! \brief Backtracking Linesearch //! //! ref : //! - Algo A6.3.1 : Dennis and Schnabel (1983) //! - Nocedal & Wrigth (2006) ParabolicLinesearchReturnCode backtracking_linesearch( Vector& update, Vector& displacements, scalar_t& lambda_out ); //! Update the variables in displacement void update_variable( const Vector& update, scalar_t lambda, Vector& displacement ); //! \brief Set the predictor void set_predictor(Vector& displacement); //! \brief perform the linesearch ParabolicDriverReturnCode linesearch( Vector &update, Vector &displacements ); //! Compute the residuals for the linesearch double compute_residuals_linesearch( Vector& update, scalar_t lambda, Vector& displacement ); double compute_residuals_strang_linesearch( Vector &update, scalar_t lambda, Vector &displacement ); scalar_t update_norm(const Vector& update); double m_norm_0 {0.0}; double m_current_dt {-1.0}; Eigen::VectorXd m_gradient; Eigen::VectorXd m_velocity; Eigen::VectorXd m_residuals; Eigen::VectorXd m_predictor; Eigen::SparseMatrix m_jacobian; sparse_solvers::SparseSolverPtr, Vector, Vector> m_solver; bool m_velocity_is_initialized {false}; }; } // end namespace dfpmsolver } // end namespace specmicp // implementation // ============== #include "parabolic_driver.inl" #endif // SPECMICP_DFPMSOLVER_PARABOLICDRIVER_HPP diff --git a/src/dfpm/solver/parabolic_driver.inl b/src/dfpm/solver/parabolic_driver.inl index 72955db..e343091 100644 --- a/src/dfpm/solver/parabolic_driver.inl +++ b/src/dfpm/solver/parabolic_driver.inl @@ -1,548 +1,548 @@ /* ============================================================================= Copyright (c) 2014 - 2016 F. Georget Princeton University All rights reserved. Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met: 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer. 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution. 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission. * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. * ============================================================================= */ #ifndef SPECMICP_DFPMSOLVER_PARABOLICDRIVER_HPP #include "parabolic_driver.hpp" // for syntaxic coloration only #endif #include "specmicp_common/log.hpp" #include "specmicp_common/eigen/incl_eigen_sparse_core.hpp" #ifdef EIGEN3_UNSUPPORTED_FOUND #include "unsupported/Eigen/SparseExtra" #endif namespace specmicp { namespace dfpmsolver { template void ParabolicDriver::compute_jacobian(Vector& displacement, Vector& velocity, Eigen::SparseMatrix& jacobian ) { program().compute_jacobian(displacement, velocity, jacobian, get_options().alpha*m_current_dt); jacobian = jacobian*scaling().asDiagonal(); jacobian.makeCompressed(); get_perfs().nb_call_jacobian += 1; } template ParabolicDriverReturnCode ParabolicDriver::solve_timestep(scalar_t dt, Eigen::VectorXd& displacement) { initialize_timestep(dt, displacement); return restart_timestep(displacement); } template void ParabolicDriver::initialize_timestep(scalar_t dt, Eigen::VectorXd& displacement) { // 1) scaling initialize_scaling(); // 2) predictor, velocity m_residuals = Eigen::VectorXd::Zero(get_neq()); m_current_dt = dt; set_predictor(displacement); reset_velocity(); program().apply_bc(dt, displacement, m_velocity); // 3) initial residuals - compute_residuals(displacement, m_velocity, m_residuals); + compute_residuals_0(displacement, m_velocity, m_residuals); m_norm_0 = norm_residuals(); // - If the norm is too low, we set it to 1 // This is important in reactive transport when // an equation can be driven out of equilibrium // during a timestep if (m_norm_0 < get_options().absolute_tolerance) m_norm_0 = 1.0; // 4) other initialisations get_perfs().nb_iterations = 0; } template void ParabolicDriver::set_predictor(Vector& displacement) { if (get_options().alpha < 1) m_predictor = displacement + (1.0-get_options().alpha)*m_current_dt*m_velocity; else m_predictor = displacement; } template scalar_t ParabolicDriver::update_norm(const Vector& update) { // l-∞ scaled norm scalar_t norm = 0.0; for (index_t dof=0; dof ParabolicDriverReturnCode ParabolicDriver::restart_timestep(Vector& displacement) { ParabolicDriverReturnCode return_code = ParabolicDriverReturnCode::NotConvergedYet; //m_solver.reset(nullptr); Eigen::VectorXd update(get_neq()); update.setZero(); get_perfs().current_update = 0; bool force_recompute_jacobian = true; while (return_code == ParabolicDriverReturnCode::NotConvergedYet) { compute_residuals(displacement, m_velocity, m_residuals); get_perfs().absolute_residual = m_residuals.norm(); get_perfs().current_residual = m_residuals.norm()/m_norm_0; get_perfs().current_update = update_norm(update); DEBUG << " NB iterations : " << get_perfs().nb_iterations << " - res : " << get_perfs().current_residual << " - update : " << get_perfs().current_update; return_code = check_convergence(); if (return_code != ParabolicDriverReturnCode::NotConvergedYet) break; if (m_solver == nullptr) { m_solver = sparse_solvers::get_sparse_solver< Eigen::SparseMatrix, Vector, Vector>(get_options().sparse_solver); compute_jacobian(displacement, m_velocity, m_jacobian); m_solver->analyse_pattern(m_jacobian); sparse_solvers::SparseSolverReturnCode retcode = m_solver->decompose(m_jacobian); if (retcode != sparse_solvers::SparseSolverReturnCode::Success) { ERROR << "Error when factorizing jacobian : " << (int) retcode; #ifdef EIGEN3_UNSUPPORTED_FOUND WARNING << "Saved matrix in : dfpm_matrix.mkt"; Eigen::saveMarket(m_jacobian, "dfpm_matrix.mkt"); #endif return ParabolicDriverReturnCode::ErrorLinearSystem; } force_recompute_jacobian =false; } else if (force_recompute_jacobian or get_perfs().nb_iterations % get_options().quasi_newton == 0) { compute_jacobian(displacement, m_velocity, m_jacobian); sparse_solvers::SparseSolverReturnCode retcode = m_solver->decompose(m_jacobian); if (retcode != sparse_solvers::SparseSolverReturnCode::Success) { ERROR << "Error when factorizing jacobian : " << (int) retcode; #ifdef EIGEN3_UNSUPPORTED_FOUND WARNING << "Saved matrix in : dfpm_matrix.mkt"; Eigen::saveMarket(m_jacobian, "dfpm_matrix.mkt"); #endif return ParabolicDriverReturnCode::ErrorLinearSystem; } } get_perfs().nb_iterations += 1; m_gradient = m_jacobian.transpose()*m_residuals; sparse_solvers::SparseSolverReturnCode retcode = m_solver->solve_scaling(m_residuals, scaling(), update); if (retcode != sparse_solvers::SparseSolverReturnCode::Success) { ERROR << "Error when solving linear system : " << (int) retcode << std::endl; return ParabolicDriverReturnCode::ErrorLinearSystem; } //if (update.norm() < get_options().step_tolerance) return_code = ParabolicDriverReturnCode::ErrorMinimized; //else return_code = linesearch(update, displacement); } return return_code; } template ParabolicDriverReturnCode ParabolicDriver::check_convergence() { ParabolicDriverReturnCode termcode = ParabolicDriverReturnCode::NotConvergedYet; const int nb_iterations = get_perfs().nb_iterations; const scalar_t norm_residuals = get_perfs().current_residual; const scalar_t norm_update = get_perfs().current_update; //std::cout << "Residuals : " << nb_iterations << " - " << norm_residuals/m_norm_0 << std::endl; DEBUG << "Residuals : " << nb_iterations << " - " << norm_residuals/m_norm_0; if (norm_residuals < get_options().residuals_tolerance) { termcode = ParabolicDriverReturnCode::ResidualMinimized; } else if (get_perfs().absolute_residual < get_options().absolute_tolerance) { termcode = ParabolicDriverReturnCode::ResidualMinimized; } else if (nb_iterations > 0 and norm_update > 0.0 and norm_update < 1.01*get_options().step_tolerance) { if (norm_residuals > get_options().threshold_stationary_point) { ERROR << "Stationary point detected !"; termcode = ParabolicDriverReturnCode::StationaryPoint; } WARNING << "Error is minimized - may indicate a stationnary point"; termcode = ParabolicDriverReturnCode::ErrorMinimized; } else if (nb_iterations > get_options().maximum_iterations) { ERROR << "Maximum number of iterations reached (" << get_options().maximum_iterations << ")"; termcode = ParabolicDriverReturnCode::MaxIterations; } else if (get_perfs().maximum_step_taken) { get_perfs().nb_consecutive_max_step_taken += 1; get_perfs().nb_max_step_taken += 1; if (get_perfs().nb_consecutive_max_step_taken == get_options().max_iterations_at_max_length) { ERROR << "Divergence detected - Maximum step length taken two many times"; termcode = ParabolicDriverReturnCode::MaxStepTakenTooManyTimes; } } else { get_perfs().nb_consecutive_max_step_taken = 0; } get_perfs().return_code = termcode; return termcode; } template double ParabolicDriver::compute_residuals_linesearch( Vector &update, scalar_t lambda, Vector &displacement ) { Eigen::VectorXd velocity(m_velocity); Eigen::VectorXd residual = Eigen::VectorXd::Zero(get_neq()); program().update_solution(update, lambda, get_options().alpha*m_current_dt, m_predictor, displacement, velocity); compute_residuals(displacement, velocity, residual); return 0.5*residual.squaredNorm(); } template double ParabolicDriver::compute_residuals_strang_linesearch( Vector &update, scalar_t lambda, Vector &displacement ) { Eigen::VectorXd velocity(m_velocity); Eigen::VectorXd residual = Eigen::VectorXd::Zero(get_neq()); program().update_solution(update, lambda, get_options().alpha*m_current_dt, m_predictor, displacement, velocity); compute_residuals(displacement, velocity, residual); return update.dot(residual); } template void ParabolicDriver::update_variable( const Vector& update, scalar_t lambda, Vector& displacement ) { program().update_solution(update, lambda, get_options().alpha*m_current_dt, m_predictor, displacement, m_velocity); } template ParabolicDriverReturnCode ParabolicDriver::linesearch( Vector &update, Vector &displacements ) { base::is_step_too_long(update); get_perfs().maximum_step_taken = false; scalar_t lambda; ParabolicLinesearchReturnCode retcode; switch (get_options().linesearch) { case ParabolicLinesearch::Backtracking: retcode = backtracking_linesearch(update, displacements, lambda); break; case ParabolicLinesearch::Strang: retcode = strang_linesearch(update, displacements, lambda); break; default: throw std::runtime_error("Linesearch type for Parabolic driver is not recognized"); break; } if (retcode != ParabolicLinesearchReturnCode::Success) { return ParabolicDriverReturnCode::LinesearchFailed; } update_variable(update, lambda, displacements); update *= lambda; return ParabolicDriverReturnCode::NotConvergedYet; } namespace internal { //! \brief Return true if a and b have the same sign inline bool have_same_sign(double a, double b) { return (std::copysign(1., a) * std::copysign(1., b) > 0); } } // end namespace internal template ParabolicLinesearchReturnCode ParabolicDriver::strang_linesearch( Vector& update, Vector& displacements, scalar_t& lambda ) { DEBUG << "Strang linesearch"; Eigen::VectorXd xp(displacements); const scalar_t s_tol = 0.5; const scalar_t s_max = 2.0; const int lin_max = 10; scalar_t s_b = 0.0; scalar_t g_b = 0.5*m_residuals.squaredNorm(); scalar_t s_a = 1.0; scalar_t g_a = compute_residuals_strang_linesearch(update, 1.0, xp); scalar_t newtlen = (scaling().asDiagonal()*update).norm(); // const scalar_t s_r = s_a; const scalar_t g_r = g_a; if (std::abs(g_a) <= s_tol*std::abs(g_b)) { DEBUG << "Skip linesearch "; lambda = 1.0; if (lambda == 1.0 and (newtlen > 0.99 * get_options().maximum_step_length)) { get_perfs().maximum_step_taken = true; } return ParabolicLinesearchReturnCode::Success; } while (internal::have_same_sign(g_a, g_b) and s_a < s_max) { s_b = s_a; s_a = 2*s_a; g_b = g_a; g_a = compute_residuals_strang_linesearch(update, s_a, xp); } scalar_t g_ = g_a; scalar_t g_0 = g_a; scalar_t s = s_a; int l; for (l=0; l= g_r) { WARNING << "Failed to find better update in Strang linesearch"; //lambda = 0.1*s_r; return backtracking_linesearch(update, displacements, lambda); } lambda = s_a; if (lambda == 1.0 and (newtlen > 0.99 * get_options().maximum_step_length)) { get_perfs().maximum_step_taken = true; } return ParabolicLinesearchReturnCode::Success; } template ParabolicLinesearchReturnCode ParabolicDriver::backtracking_linesearch( Vector& update, Vector& displacements, scalar_t& lambda_out ) { // References // ---------- // - Algo A6.3.1 : Dennis and Schnabel (1983) // - Nocedal & Wrigth (2006) DEBUG << "Linesearch"; Eigen::VectorXd xp(displacements); double fcp; int retcode = 2; // 2 not converged, 1 problem, 0 success const scalar_t alpha = 1e-4; scalar_t newtlen = (scaling().asDiagonal()*update).norm(); scalar_t init_slope = m_gradient.dot(update); const scalar_t rellength = update_norm(update); const scalar_t minlambda = get_options().step_tolerance / rellength; scalar_t lambda = 1.0; scalar_t lambda_prev = lambda; scalar_t merit_value = 0.5*m_residuals.squaredNorm(); // new residual fcp = compute_residuals_linesearch(update, lambda, xp); // Skip linesearch if enough progress is done // ------------------------------------------ if (fcp < get_options().coeff_accept_newton_step *merit_value) { DEBUG << "Skip linesearch "; lambda_out = 1.0; return ParabolicLinesearchReturnCode::Success; } // The linesearch // -------------- scalar_t fc = merit_value; scalar_t fcp_prev; int cnt = 0; do { SPAM << "cnt : " < 0.99 * get_options().maximum_step_length)) { get_perfs().maximum_step_taken = true; } break; } else if (lambda < minlambda) { retcode = 0; lambda = minlambda; //retcode = 1; break; } else { // Select a new step length // - - - - - - - - - - - - double lambdatmp; if (cnt == 0) { // only a quadratic at the first lambdatmp = - init_slope / (2*(fcp - fc -init_slope)); } else { const scalar_t factor = 1.0 /(lambda - lambda_prev); const scalar_t x1 = fcp - fc - lambda*init_slope; const scalar_t x2 = fcp_prev - fc - lambda_prev*init_slope; const scalar_t a = factor * ( x1/(lambda*lambda) - x2/(lambda_prev*lambda_prev)); const scalar_t b = factor * ( -x1*lambda_prev/(lambda*lambda) + x2*lambda/(lambda_prev*lambda_prev)); if (a == 0) { // cubic interpolation is in fact a quadratic interpolation lambdatmp = - init_slope/(2*b); } else { const scalar_t disc = b*b-3*a*init_slope; lambdatmp = (-b+std::sqrt(disc))/(3*a); } if (lambdatmp > 0.5*lambda ) lambdatmp = 0.5*lambda; } lambda_prev = lambda; fcp_prev = fcp; if (lambdatmp < 0.1*lambda) { lambda = 0.1 * lambda; } else { lambda = lambdatmp; } } if (not std::isfinite(lambda)) { ERROR << "Lambda is non finite - we stop \n" << "Proposed update : \n ---- \n" << update << "\n ---- \n" << "Current displacements \n ---- \n" << displacements << "\n ---- \n"; return ParabolicLinesearchReturnCode::Divergence; } fcp = compute_residuals_linesearch(update, lambda, xp); ++cnt; } while(retcode == 2 and cnt < 50); if (cnt == 50) { ERROR << "Too much linesearch iterations ! We stop"; return ParabolicLinesearchReturnCode::MaximumIterations; } lambda_out = lambda; switch (retcode) { case 0: return ParabolicLinesearchReturnCode::Success; break; case 1: return ParabolicLinesearchReturnCode::LambdaTooSmall; default: return ParabolicLinesearchReturnCode::NotSupposedToHappen; break; } } template void ParabolicDriver::set_velocity(Vector& velocity_vector) { m_velocity.resize(program().get_tot_ndf()); m_velocity.setZero(); for (index_t dof=0; dof void ParabolicDriver::reset_velocity() { for (index_t dof=0; dof Princeton University All rights reserved. Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met: 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer. 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution. 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission. * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. * ============================================================================= */ #ifndef SPECMICP_DFPMSOLVER_PARABOLICPROGRAM_HPP #define SPECMICP_DFPMSOLVER_PARABOLICPROGRAM_HPP #include "specmicp_common/types.hpp" #include "specmicp_common/eigen/incl_eigen_sparse_core.hpp" #include #include "dfpm_program.hpp" namespace specmicp { namespace dfpmsolver { //! Concept class for a program template class ParabolicProgram: DFPMProgram> { public: - Derived* derived() {static_cast(this);} + Derived* derived() {return static_cast(this);} //! \brief Compute the residuals void compute_residuals(const Vector& displacement, const Vector& velocity, Vector& residual ) { derived()->compute_residuals(displacement, velocity, residual); } + //! \brief Compute the initial residuals + //! + //! By default, it calls the same residuals as compute_residuals + void compute_residuals_0(const Vector& displacement, + const Vector& velocity, + Vector& residual + ) + { + derived()->compute_residuals(displacement, velocity, residual); + } //! \brief Compute the jacobian void compute_jacobian(Vector& displacement, Vector& velocity, Eigen::SparseMatrix& jacobian, scalar_t alphadt ) { derived()->compute_jacobian(displacement, velocity, jacobian, alphadt); } //! \brief Update the solutions void update_solution(const Vector& update, scalar_t lambda, scalar_t alpha_dt, Vector& predictor, Vector& displacement, Vector& velocity) { derived()->update_solution(update, lambda, alpha_dt, predictor, displacement, velocity); } //! \brief Apply boundary conditions to the velocity vector //! //! by default do nothing. void apply_bc(scalar_t dt, const Vector& displacement, Vector& velocity) {} }; } // end namespace dfpmsolver } // end namespace specmicp #endif // SPECMICP_DFPMSOLVER_PARABOLICPROGRAM_HPP