diff --git a/python/py_model.cc b/python/py_model.cc index 025d114e6..69aaeb72b 100644 --- a/python/py_model.cc +++ b/python/py_model.cc @@ -1,74 +1,75 @@ /* -------------------------------------------------------------------------- */ #include "py_aka_array.hh" /* -------------------------------------------------------------------------- */ #include #include #include /* -------------------------------------------------------------------------- */ #include #include #include /* -------------------------------------------------------------------------- */ namespace py = pybind11; /* -------------------------------------------------------------------------- */ namespace akantu { /* -------------------------------------------------------------------------- */ void register_model(py::module & mod) { py::class_(mod, "SparseMatrix") .def("getMatrixType", &SparseMatrix::getMatrixType) .def("size", &SparseMatrix::size); py::class_(mod, "SparseMatrixAIJ") .def("getIRN", &SparseMatrixAIJ::getIRN) .def("getJCN", &SparseMatrixAIJ::getJCN) .def("getA", &SparseMatrixAIJ::getA); py::class_(mod, "DOFManager") .def("getMatrix", [](DOFManager & self, const std::string & name) { return dynamic_cast( self.getMatrix(name)); }, py::return_value_policy::reference); py::class_(mod, "NonLinearSolver") .def( "set", [](NonLinearSolver & self, const std::string & id, const Real & val) { if (id == "max_iterations") self.set(id, int(val)); else self.set(id, val); }) .def("set", [](NonLinearSolver & self, const std::string & id, const SolveConvergenceCriteria & val) { self.set(id, val); }); py::class_(mod, "ModelSolver", py::multiple_inheritance()) .def("getNonLinearSolver", (NonLinearSolver & (ModelSolver::*)(const ID &)) & ModelSolver::getNonLinearSolver, py::arg("solver_id") = "", py::return_value_policy::reference) - .def("solveStep", &ModelSolver::solveStep, py::arg("solver_id") = ""); + .def("solveStep", [](ModelSolver & self, const ID & solver_id = "") { + self.solveStep(solver_id); + }); py::class_(mod, "Model", py::multiple_inheritance()) .def("setBaseName", &Model::setBaseName) .def("getFEEngine", &Model::getFEEngine, py::arg("name") = "", py::return_value_policy::reference) .def("addDumpFieldVector", &Model::addDumpFieldVector) .def("addDumpField", &Model::addDumpField) .def("setBaseNameToDumper", &Model::setBaseNameToDumper) .def("addDumpFieldVectorToDumper", &Model::addDumpFieldVectorToDumper) .def("addDumpFieldToDumper", &Model::addDumpFieldToDumper) .def("dump", &Model::dump) .def("initNewSolver", &Model::initNewSolver) .def("getDOFManager", &Model::getDOFManager, py::return_value_policy::reference); - } } // namespace akantu diff --git a/src/model/common/model_solver.cc b/src/model/common/model_solver.cc index c7d50b993..93994a487 100644 --- a/src/model/common/model_solver.cc +++ b/src/model/common/model_solver.cc @@ -1,365 +1,370 @@ /** * @file model_solver.cc * * @author Nicolas Richart * * @date creation: Tue Aug 18 2015 * @date last modification: Wed Feb 21 2018 * * @brief Implementation of ModelSolver * * @section LICENSE * * Copyright (©) 2015-2018 EPFL (Ecole Polytechnique Fédérale de Lausanne) * Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides) * * Akantu is free software: you can redistribute it and/or modify it under the * terms of the GNU Lesser General Public License as published by the Free * Software Foundation, either version 3 of the License, or (at your option) any * later version. * * Akantu is distributed in the hope that it will be useful, but WITHOUT ANY * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR * A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more * details. * * You should have received a copy of the GNU Lesser General Public License * along with Akantu. If not, see . * */ /* -------------------------------------------------------------------------- */ #include "model_solver.hh" #include "dof_manager.hh" #include "dof_manager_default.hh" #include "mesh.hh" #include "non_linear_solver.hh" #include "time_step_solver.hh" #if defined(AKANTU_USE_PETSC) #include "dof_manager_petsc.hh" #endif /* -------------------------------------------------------------------------- */ namespace akantu { /* -------------------------------------------------------------------------- */ template static T getOptionToType(const std::string & opt_str) { std::stringstream sstr(opt_str); T opt; sstr >> opt; return opt; } /* -------------------------------------------------------------------------- */ ModelSolver::ModelSolver(Mesh & mesh, const ModelType & type, const ID & id, UInt memory_id) : Parsable(ParserType::_model, id), SolverCallback(), model_type(type), parent_id(id), parent_memory_id(memory_id), mesh(mesh), dof_manager(nullptr), default_solver_id("") {} /* -------------------------------------------------------------------------- */ ModelSolver::~ModelSolver() = default; /* -------------------------------------------------------------------------- */ std::tuple ModelSolver::getParserSection() { auto sub_sections = getStaticParser().getSubSections(ParserType::_model); auto it = std::find_if( sub_sections.begin(), sub_sections.end(), [&](auto && section) { ModelType type = getOptionToType(section.getName()); // default id should be the model type if not defined std::string name = section.getParameter("name", this->parent_id); return type == model_type and name == this->parent_id; }); if (it == sub_sections.end()) return std::make_tuple(ParserSection(), true); return std::make_tuple(*it, false); } /* -------------------------------------------------------------------------- */ void ModelSolver::initDOFManager() { // default without external solver activated at compilation same as mumps that // is the historical solver but with only the lumped solver ID solver_type = "default"; #if defined(AKANTU_USE_MUMPS) solver_type = "default"; #elif defined(AKANTU_USE_PETSC) solver_type = "petsc"; #endif ParserSection section; bool is_empty; std::tie(section, is_empty) = this->getParserSection(); if (not is_empty) { solver_type = section.getOption(solver_type); this->initDOFManager(section, solver_type); } else { this->initDOFManager(solver_type); } } /* -------------------------------------------------------------------------- */ void ModelSolver::initDOFManager(const ID & solver_type) { try { this->dof_manager = DOFManagerFactory::getInstance().allocate( solver_type, mesh, this->parent_id + ":dof_manager_" + solver_type, this->parent_memory_id); } catch (...) { AKANTU_EXCEPTION( "To use the solver " << solver_type << " you will have to code it. This is an unknown solver type."); } this->setDOFManager(*this->dof_manager); } /* -------------------------------------------------------------------------- */ void ModelSolver::initDOFManager(const ParserSection & section, const ID & solver_type) { this->initDOFManager(solver_type); auto sub_sections = section.getSubSections(ParserType::_time_step_solver); // parsing the time step solvers for (auto && section : sub_sections) { ID type = section.getName(); ID solver_id = section.getParameter("name", type); auto tss_type = getOptionToType(type); auto tss_options = this->getDefaultSolverOptions(tss_type); auto sub_solvers_sect = section.getSubSections(ParserType::_non_linear_solver); auto nb_non_linear_solver_section = section.getNbSubSections(ParserType::_non_linear_solver); auto nls_type = tss_options.non_linear_solver_type; if (nb_non_linear_solver_section == 1) { auto && nls_section = *(sub_solvers_sect.first); nls_type = getOptionToType(nls_section.getName()); } else if (nb_non_linear_solver_section > 0) { AKANTU_EXCEPTION("More than one non linear solver are provided for the " "time step solver " << solver_id); } this->getNewSolver(solver_id, tss_type, nls_type); if (nb_non_linear_solver_section == 1) { const auto & nls_section = *(sub_solvers_sect.first); this->dof_manager->getNonLinearSolver(solver_id).parseSection( nls_section); } auto sub_integrator_sections = section.getSubSections(ParserType::_integration_scheme); for (auto && is_section : sub_integrator_sections) { const auto & dof_type_str = is_section.getName(); ID dof_id; try { ID tmp = is_section.getParameter("name"); dof_id = tmp; } catch (...) { AKANTU_EXCEPTION("No degree of freedom name specified for the " "integration scheme of type " << dof_type_str); } auto it_type = getOptionToType(dof_type_str); IntegrationScheme::SolutionType s_type = is_section.getParameter( "solution_type", tss_options.solution_type[dof_id]); this->setIntegrationScheme(solver_id, dof_id, it_type, s_type); } for (auto & is_type : tss_options.integration_scheme_type) { if (!this->hasIntegrationScheme(solver_id, is_type.first)) { this->setIntegrationScheme(solver_id, is_type.first, is_type.second, tss_options.solution_type[is_type.first]); } } } if (section.hasParameter("default_solver")) { ID default_solver = section.getParameter("default_solver"); if (this->hasSolver(default_solver)) { this->setDefaultSolver(default_solver); } else AKANTU_EXCEPTION( "The solver \"" << default_solver << "\" was not created, it cannot be set as default solver"); } } /* -------------------------------------------------------------------------- */ TimeStepSolver & ModelSolver::getSolver(const ID & solver_id) { ID tmp_solver_id = solver_id; if (tmp_solver_id == "") tmp_solver_id = this->default_solver_id; TimeStepSolver & tss = this->dof_manager->getTimeStepSolver(tmp_solver_id); return tss; } /* -------------------------------------------------------------------------- */ const TimeStepSolver & ModelSolver::getSolver(const ID & solver_id) const { ID tmp_solver_id = solver_id; if (solver_id == "") tmp_solver_id = this->default_solver_id; const TimeStepSolver & tss = this->dof_manager->getTimeStepSolver(tmp_solver_id); return tss; } /* -------------------------------------------------------------------------- */ TimeStepSolver & ModelSolver::getTimeStepSolver(const ID & solver_id) { return this->getSolver(solver_id); } /* -------------------------------------------------------------------------- */ const TimeStepSolver & ModelSolver::getTimeStepSolver(const ID & solver_id) const { return this->getSolver(solver_id); } /* -------------------------------------------------------------------------- */ NonLinearSolver & ModelSolver::getNonLinearSolver(const ID & solver_id) { return this->getSolver(solver_id).getNonLinearSolver(); } /* -------------------------------------------------------------------------- */ const NonLinearSolver & ModelSolver::getNonLinearSolver(const ID & solver_id) const { return this->getSolver(solver_id).getNonLinearSolver(); } /* -------------------------------------------------------------------------- */ bool ModelSolver::hasSolver(const ID & solver_id) const { ID tmp_solver_id = solver_id; if (solver_id == "") tmp_solver_id = this->default_solver_id; if (not this->dof_manager) { AKANTU_EXCEPTION("No DOF manager was initialized"); } return this->dof_manager->hasTimeStepSolver(tmp_solver_id); } /* -------------------------------------------------------------------------- */ void ModelSolver::setDefaultSolver(const ID & solver_id) { AKANTU_DEBUG_ASSERT( this->hasSolver(solver_id), "Cannot set the default solver to a solver that does not exists"); this->default_solver_id = solver_id; } /* -------------------------------------------------------------------------- */ -void ModelSolver::solveStep(const ID & solver_id) { +void ModelSolver::solveStep(SolverCallback & callback, const ID & solver_id) { AKANTU_DEBUG_IN(); TimeStepSolver & tss = this->getSolver(solver_id); // make one non linear solve - tss.solveStep(*this); + tss.solveStep(callback); AKANTU_DEBUG_OUT(); } +/* -------------------------------------------------------------------------- */ +void ModelSolver::solveStep(const ID & solver_id) { + solveStep(*this, solver_id); +} + /* -------------------------------------------------------------------------- */ void ModelSolver::getNewSolver(const ID & solver_id, TimeStepSolverType time_step_solver_type, NonLinearSolverType non_linear_solver_type) { if (this->default_solver_id == "") { this->default_solver_id = solver_id; } if (non_linear_solver_type == NonLinearSolverType::_auto) { switch (time_step_solver_type) { case TimeStepSolverType::_dynamic: case TimeStepSolverType::_static: non_linear_solver_type = NonLinearSolverType::_newton_raphson; break; case TimeStepSolverType::_dynamic_lumped: non_linear_solver_type = NonLinearSolverType::_lumped; break; case TimeStepSolverType::_not_defined: AKANTU_EXCEPTION(time_step_solver_type << " is not a valid time step solver type"); break; } } this->initSolver(time_step_solver_type, non_linear_solver_type); NonLinearSolver & nls = this->dof_manager->getNewNonLinearSolver( solver_id, non_linear_solver_type); this->dof_manager->getNewTimeStepSolver(solver_id, time_step_solver_type, nls, *this); } /* -------------------------------------------------------------------------- */ Real ModelSolver::getTimeStep(const ID & solver_id) const { const TimeStepSolver & tss = this->getSolver(solver_id); return tss.getTimeStep(); } /* -------------------------------------------------------------------------- */ void ModelSolver::setTimeStep(Real time_step, const ID & solver_id) { TimeStepSolver & tss = this->getSolver(solver_id); return tss.setTimeStep(time_step); } /* -------------------------------------------------------------------------- */ void ModelSolver::setIntegrationScheme( const ID & solver_id, const ID & dof_id, const IntegrationSchemeType & integration_scheme_type, IntegrationScheme::SolutionType solution_type) { TimeStepSolver & tss = this->dof_manager->getTimeStepSolver(solver_id); tss.setIntegrationScheme(dof_id, integration_scheme_type, solution_type); } /* -------------------------------------------------------------------------- */ bool ModelSolver::hasDefaultSolver() const { return (this->default_solver_id != ""); } /* -------------------------------------------------------------------------- */ bool ModelSolver::hasIntegrationScheme(const ID & solver_id, const ID & dof_id) const { TimeStepSolver & tss = this->dof_manager->getTimeStepSolver(solver_id); return tss.hasIntegrationScheme(dof_id); } /* -------------------------------------------------------------------------- */ void ModelSolver::predictor() {} /* -------------------------------------------------------------------------- */ void ModelSolver::corrector() {} /* -------------------------------------------------------------------------- */ TimeStepSolverType ModelSolver::getDefaultSolverType() const { return TimeStepSolverType::_dynamic_lumped; } /* -------------------------------------------------------------------------- */ ModelSolverOptions ModelSolver::getDefaultSolverOptions(__attribute__((unused)) const TimeStepSolverType & type) const { ModelSolverOptions options; options.non_linear_solver_type = NonLinearSolverType::_auto; return options; } } // namespace akantu diff --git a/src/model/common/model_solver.hh b/src/model/common/model_solver.hh index e1116cafb..13d45410d 100644 --- a/src/model/common/model_solver.hh +++ b/src/model/common/model_solver.hh @@ -1,194 +1,199 @@ /** * @file model_solver.hh * * @author Nicolas Richart * * @date creation: Fri Jun 18 2010 * @date last modification: Wed Jan 31 2018 * * @brief Class regrouping the common solve interface to the different models * * @section LICENSE * * Copyright (©) 2010-2018 EPFL (Ecole Polytechnique Fédérale de Lausanne) * Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides) * * Akantu is free software: you can redistribute it and/or modify it under the * terms of the GNU Lesser General Public License as published by the Free * Software Foundation, either version 3 of the License, or (at your option) any * later version. * * Akantu is distributed in the hope that it will be useful, but WITHOUT ANY * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR * A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more * details. * * You should have received a copy of the GNU Lesser General Public License * along with Akantu. If not, see . * */ /* -------------------------------------------------------------------------- */ #include "aka_common.hh" #include "integration_scheme.hh" #include "parsable.hh" #include "solver_callback.hh" #include "synchronizer_registry.hh" /* -------------------------------------------------------------------------- */ #include /* -------------------------------------------------------------------------- */ #ifndef __AKANTU_MODEL_SOLVER_HH__ #define __AKANTU_MODEL_SOLVER_HH__ namespace akantu { class Mesh; class DOFManager; class TimeStepSolver; class NonLinearSolver; struct ModelSolverOptions; } namespace akantu { class ModelSolver : public Parsable, public SolverCallback, public SynchronizerRegistry { /* ------------------------------------------------------------------------ */ /* Constructors/Destructors */ /* ------------------------------------------------------------------------ */ public: ModelSolver(Mesh & mesh, const ModelType & type, const ID & id, UInt memory_id); ~ModelSolver() override; /// initialize the dof manager based on solver type passed in the input file void initDOFManager(); /// initialize the dof manager based on the used chosen solver type void initDOFManager(const ID & solver_type); protected: /// initialize the dof manager based on the used chosen solver type void initDOFManager(const ParserSection & section, const ID & solver_type); /// Callback for the model to instantiate the matricees when needed virtual void initSolver(TimeStepSolverType /*time_step_solver_type*/, NonLinearSolverType /*non_linear_solver_type*/) {} /// get the section in the input file (if it exsits) corresponding to this /// model std::tuple getParserSection(); /* ------------------------------------------------------------------------ */ /* Methods */ /* ------------------------------------------------------------------------ */ public: /// solve a step using a given pre instantiated time step solver and /// nondynamic linear solver virtual void solveStep(const ID & solver_id = ""); + /// solve a step using a given pre instantiated time step solver and + /// non linear solver with a user defined callback instead of the + /// model itself /!\ This can mess up everything + virtual void solveStep(SolverCallback & callback, const ID & solver_id = ""); + /// Initialize a time solver that can be used afterwards with its id void getNewSolver( const ID & solver_id, TimeStepSolverType time_step_solver_type, NonLinearSolverType non_linear_solver_type = NonLinearSolverType::_auto); /// set an integration scheme for a given dof and a given solver void setIntegrationScheme(const ID & solver_id, const ID & dof_id, const IntegrationSchemeType & integration_scheme_type, IntegrationScheme::SolutionType solution_type = IntegrationScheme::_not_defined); /// set an externally instantiated integration scheme void setIntegrationScheme(const ID & solver_id, const ID & dof_id, IntegrationScheme & integration_scheme, IntegrationScheme::SolutionType solution_type = IntegrationScheme::_not_defined); /* ------------------------------------------------------------------------ */ /* SolverCallback interface */ /* ------------------------------------------------------------------------ */ public: /// Predictor interface for the callback void predictor() override; /// Corrector interface for the callback void corrector() override; /* ------------------------------------------------------------------------ */ /* Accessors */ /* ------------------------------------------------------------------------ */ public: /// Default time step solver to instantiate for this model virtual TimeStepSolverType getDefaultSolverType() const; /// Default configurations for a given time step solver virtual ModelSolverOptions getDefaultSolverOptions(const TimeStepSolverType & type) const; /// get access to the internal dof manager DOFManager & getDOFManager() { return *this->dof_manager; } /// get the time step of a given solver Real getTimeStep(const ID & solver_id = "") const; /// set the time step of a given solver virtual void setTimeStep(Real time_step, const ID & solver_id = ""); /// set the parameter 'param' of the solver 'solver_id' // template // void set(const ID & param, const T & value, const ID & solver_id = ""); /// get the parameter 'param' of the solver 'solver_id' // const Parameter & get(const ID & param, const ID & solver_id = "") const; /// answer to the question "does the solver exists ?" bool hasSolver(const ID & solver_id) const; /// changes the current default solver void setDefaultSolver(const ID & solver_id); /// is a default solver defined bool hasDefaultSolver() const; /// is an integration scheme set for a given solver and a given dof bool hasIntegrationScheme(const ID & solver_id, const ID & dof_id) const; TimeStepSolver & getTimeStepSolver(const ID & solver_id = ""); NonLinearSolver & getNonLinearSolver(const ID & solver_id = ""); const TimeStepSolver & getTimeStepSolver(const ID & solver_id = "") const; const NonLinearSolver & getNonLinearSolver(const ID & solver_id = "") const; private: TimeStepSolver & getSolver(const ID & solver_id); const TimeStepSolver & getSolver(const ID & solver_id) const; /* ------------------------------------------------------------------------ */ /* Class Members */ /* ------------------------------------------------------------------------ */ protected: ModelType model_type; private: ID parent_id; UInt parent_memory_id; /// Underlying mesh Mesh & mesh; /// Underlying dof_manager (the brain...) std::unique_ptr dof_manager; /// Default time step solver to use ID default_solver_id; }; struct ModelSolverOptions { NonLinearSolverType non_linear_solver_type; std::map integration_scheme_type; std::map solution_type; }; } // akantu #endif /* __AKANTU_MODEL_SOLVER_HH__ */