29 #include "libmesh/numeric_vector.h" 30 #include "libmesh/dof_map.h" 31 #include "libmesh/sparse_matrix.h" 32 #include "libmesh/linear_solver.h" 41 _first_sensitivity_step (true),
43 _n_iters_to_store (n),
44 _assembly_ops (nullptr),
45 _if_highest_derivative_solution (false) {
62 libmesh_assert(_assembly_ops);
64 _assembly_ops->set_assembly(assembly);
73 _assembly_ops->clear_assembly();
80 libmesh_assert(_assembly_ops);
82 return *_assembly_ops;
99 _assembly_ops = &assembly_ops;
103 for (
unsigned int i=0; i<_n_iters_to_store; i++) {
104 std::ostringstream iter;
111 nm =
"transient_solution_";
113 sys.add_vector(nm,
true, libMesh::GHOSTED);
116 nm =
"transient_solution_sensitivity_";
118 sys.add_vector(nm,
true, libMesh::GHOSTED);
121 nm =
"transient_velocity_";
123 sys.add_vector(nm,
true, libMesh::GHOSTED);
124 nm =
"transient_velocity_sensitivity_";
126 sys.add_vector(nm,
true, libMesh::GHOSTED);
128 if (_ode_order > 1) {
130 nm =
"transient_acceleration_";
132 sys.add_vector(nm,
true, libMesh::GHOSTED);
133 nm =
"transient_acceleration_sensitivity_";
135 sys.add_vector(nm,
true, libMesh::GHOSTED);
157 for (
unsigned int i=0; i<_n_iters_to_store; i++) {
158 std::ostringstream iter;
162 nm =
"transient_solution_";
164 if (sys.have_vector(nm)) sys.remove_vector(nm);
165 nm =
"transient_solution_sensitivity_";
167 if (sys.have_vector(nm)) sys.remove_vector(nm);
170 nm =
"transient_velocity_";
172 if (sys.have_vector(nm)) sys.remove_vector(nm);
173 nm =
"transient_velocity_sensitivity_";
175 if (sys.have_vector(nm)) sys.remove_vector(nm);
177 if (_ode_order > 1) {
179 nm =
"transient_acceleration_";
181 if (sys.have_vector(nm)) sys.remove_vector(nm);
182 nm =
"transient_acceleration_sensitivity_";
184 if (sys.have_vector(nm)) sys.remove_vector(nm);
189 _assembly_ops =
nullptr;
196 libMesh::NumericVector<Real>&
200 libmesh_assert_less(prev_iter, _n_iters_to_store);
203 std::ostringstream oss;
209 nm =
"transient_solution_" + oss.str();
219 libMesh::NumericVector<Real>&
223 libmesh_assert_less(prev_iter, _n_iters_to_store);
226 std::ostringstream oss;
231 nm =
"transient_solution_sensitivity_" + oss.str();
238 libMesh::NumericVector<Real>&
242 libmesh_assert_less(prev_iter, _n_iters_to_store);
245 std::ostringstream oss;
250 nm =
"transient_velocity_" + oss.str();
257 libMesh::NumericVector<Real>&
261 libmesh_assert_less(prev_iter, _n_iters_to_store);
264 std::ostringstream oss;
269 nm =
"transient_velocity_sensitivity_" + oss.str();
277 libMesh::NumericVector<Real>&
281 libmesh_assert_less(prev_iter, _n_iters_to_store);
284 std::ostringstream oss;
289 nm =
"transient_acceleration_" + oss.str();
296 libMesh::NumericVector<Real>&
300 libmesh_assert_less(prev_iter, _n_iters_to_store);
303 std::ostringstream oss;
308 nm =
"transient_acceleration_sensitivity_" + oss.str();
319 libmesh_assert_msg(
_system,
"System pointer is nullptr.");
333 libmesh_assert_msg(
_system,
"System pointer is nullptr.");
346 bool increment_time) {
348 libmesh_assert(_first_step);
355 _if_highest_derivative_solution =
true;
366 _if_highest_derivative_solution =
false;
368 std::pair<unsigned int, Real>
372 libMesh::SparseMatrix<Real> *
373 pc = sys.request_matrix(
"Preconditioner");
375 std::unique_ptr<libMesh::NumericVector<Real> >
376 dvec(
velocity().zero_clone().release());
381 solver_params.second,
382 solver_params.first);
384 libMesh::NumericVector<Real> *vec =
nullptr;
386 switch (_ode_order) {
401 vec->add(-1., *dvec);
404 #ifdef LIBMESH_ENABLE_CONSTRAINTS 405 sys.get_dof_map().enforce_constraints_exactly(sys, vec,
true);
412 for (
unsigned int i=_n_iters_to_store-1; i>0; i--) {
421 if (increment_time) sys.time +=
dt;
434 libmesh_assert(_first_sensitivity_step);
441 _if_highest_derivative_solution =
true;
454 _if_highest_derivative_solution =
false;
456 std::pair<unsigned int, Real>
460 libMesh::SparseMatrix<Real> *
461 pc = sys.request_matrix(
"Preconditioner");
463 std::unique_ptr<libMesh::NumericVector<Real> >
464 dvec(
velocity().zero_clone().release());
469 solver_params.second,
470 solver_params.first);
472 libMesh::NumericVector<Real> *vec =
nullptr;
474 switch (_ode_order) {
489 vec->add(-1., *dvec);
492 #ifdef LIBMESH_ENABLE_CONSTRAINTS 493 sys.get_dof_map().enforce_constraints_exactly(sys, vec,
true);
498 for (
unsigned int i=_n_iters_to_store-1; i>0; i--) {
509 _first_sensitivity_step =
false;
517 std::vector<libMesh::NumericVector<Real>*>& sol) {
526 libmesh_assert(!sol.size());
529 sol.resize(_ode_order+1);
531 const std::vector<libMesh::dof_id_type>&
532 send_list = sys.get_dof_map().get_send_list();
535 for (
unsigned int i=0; i<=_ode_order; i++) {
537 sol[i] = libMesh::NumericVector<Real>::build(sys.comm()).release();
538 sol[i]->init(sys.n_dofs(),
548 if (!_if_highest_derivative_solution)
550 current_sol.localize(*sol[i], send_list);
552 solution().localize(*sol[i], send_list);
559 libMesh::NumericVector<Real>& vel = this->
velocity();
561 if (!_if_highest_derivative_solution)
567 vel.localize(*sol[i], send_list);
574 libMesh::NumericVector<Real>& acc = this->
acceleration();
576 if (!_if_highest_derivative_solution)
582 acc.localize(*sol[i], send_list);
599 std::vector<libMesh::NumericVector<Real>*>& sol) {
603 libmesh_assert_less_equal(prev_iter, _n_iters_to_store);
609 libmesh_assert(!sol.size());
612 sol.resize(_ode_order+1);
614 const std::vector<libMesh::dof_id_type>&
615 send_list = sys.get_dof_map().get_send_list();
618 for (
unsigned int i=0; i<=_ode_order; i++) {
620 sol[i] = libMesh::NumericVector<Real>::build(sys.comm()).release();
621 sol[i]->init(sys.n_dofs(),
654 MAST::TransientSolverBase::
655 build_perturbed_local_quantities(
const libMesh::NumericVector<Real>& current_dsol,
656 std::vector<libMesh::NumericVector<Real>*>& sol) {
665 libmesh_assert(!sol.size());
668 sol.resize(_ode_order+1);
670 const std::vector<libMesh::dof_id_type>&
671 send_list = sys.get_dof_map().get_send_list();
674 std::unique_ptr<libMesh::NumericVector<Real> >
675 tmp(this->
solution().zero_clone().release());
677 for (
unsigned int i=0; i<=_ode_order; i++) {
679 sol[i] = libMesh::NumericVector<Real>::build(sys.comm()).release();
680 sol[i]->init(sys.n_dofs(),
690 current_dsol.localize(*sol[i], send_list);
691 if (_if_highest_derivative_solution)
700 if (_ode_order == 1 && _if_highest_derivative_solution)
702 current_dsol.localize(*sol[i], send_list);
703 else if (_if_highest_derivative_solution)
707 tmp->localize(*sol[i], send_list);
714 if (_ode_order == 2 && _if_highest_derivative_solution)
716 current_dsol.localize(*sol[i], send_list);
717 else if (_if_highest_derivative_solution)
721 tmp->localize(*sol[i], send_list);
755 for (
unsigned int i=_n_iters_to_store-1; i>0; i--) {
764 if (increment_time) sys.time +=
dt;
784 for (
unsigned int i=_n_iters_to_store-1; i>0; i--) {
794 _first_sensitivity_step =
false;
801 const libMesh::Elem& ref_elem,
804 libmesh_assert(_assembly_ops);
805 _assembly_ops->set_elem_data(dim, ref_elem, elem);
812 libmesh_assert(_assembly_ops);
813 _assembly_ops->init(elem);
822 libmesh_assert(_assembly_ops);
823 _assembly_ops->clear_elem();
const MAST::NonlinearSystem & system() const
virtual void build_local_quantities(const libMesh::NumericVector< Real > ¤t_sol, std::vector< libMesh::NumericVector< Real > *> &qtys)
localizes the relevant solutions for system assembly.
MAST::NonlinearSystem & system()
virtual void advance_time_step_with_sensitivity()
advances the time step and copies the current sensitivity solution to old sensitivity solution...
virtual void advance_time_step(bool increment_time=true)
advances the time step and copies the current solution to old solution, and so on.
virtual void sensitivity_solve(MAST::AssemblyBase &assembly, const MAST::FunctionBase &f)
solvers the current time step for sensitivity wrt f
virtual void update_sensitivity_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0
update the transient sensitivity acceleration based on the current sensitivity solution ...
This class implements a system for solution of nonlinear systems.
virtual void set_assembly(MAST::AssemblyBase &assembly)
sets the assembly object
virtual MAST::TransientAssemblyElemOperations & get_elem_operation_object()
virtual void clear_elem_operation_object()
clears the association of this object with the assembly element operation object. ...
virtual void update_delta_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0
update the perturbation in transient acceleration based on the current perturbed solution ...
virtual void set_assembly(MAST::AssemblyBase &assembly)
sets the assembly object
virtual void update_delta_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0
update the perturbation in transient velocity based on the current perturbed solution ...
virtual ~TransientSolverBase()
libMesh::NumericVector< Real > & solution(unsigned int prev_iter=0) const
virtual void set_elem_operation_object(MAST::AssemblyElemOperations &elem_ops)
attaches a element operation to this object, and associated this with the element operation object...
MAST::PhysicsDisciplineBase * _discipline
virtual void set_elem_operation_object(MAST::TransientAssemblyElemOperations &elem_ops)
Attaches the assembly elem operations object that provides the x_dot, M and J quantities for the elem...
libMesh::NumericVector< Real > & solution_sensitivity(unsigned int prev_iter=0) const
libMesh::NumericVector< Real > & acceleration_sensitivity(unsigned int prev_iter=0) const
virtual void set_elem_data(unsigned int dim, const libMesh::Elem &ref_elem, MAST::GeomElem &elem) const =0
some analyses may want to set additional element data before initialization of the GeomElem...
virtual void update_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0
update the transient velocity based on the current solution
libMesh::NumericVector< Real > & velocity(unsigned int prev_iter=0) const
virtual void clear_assembly()
clears the assembly object
libMesh::NumericVector< Real > & acceleration(unsigned int prev_iter=0) const
virtual std::pair< unsigned int, Real > get_linear_solve_parameters()
calls NonlinearImplicitSystem::set_solver_parameters() before accessing the values.
virtual void init(const MAST::GeomElem &elem)=0
initializes the object for calculation of element quantities for the specified elem.
virtual void update_acceleration(libMesh::NumericVector< Real > &acc, const libMesh::NumericVector< Real > &sol)=0
update the transient acceleration based on the current solution
void solve_highest_derivative_and_advance_time_step(MAST::AssemblyBase &assembly, bool increment_time=true)
To be used only for initial conditions.
virtual void clear_assembly()
clears the assembly object
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
virtual void solve(MAST::AssemblyBase &assembly)
solves the current time step for solution and velocity
MAST::AssemblyBase * _assembly
std::unique_ptr< libMesh::LinearSolver< Real > > linear_solver
The LinearSolver for solution of the linear equations.
virtual void clear_elem()
clears the element initialization
virtual void residual_and_jacobian(const libMesh::NumericVector< Real > &X, libMesh::NumericVector< Real > *R, libMesh::SparseMatrix< Real > *J, libMesh::NonlinearImplicitSystem &S)
function that assembles the matrices and vectors quantities for nonlinear solution ...
virtual void build_sensitivity_local_quantities(unsigned int prev_iter, std::vector< libMesh::NumericVector< Real > *> &qtys)
localizes the relevant solutions for system assembly.
libMesh::NumericVector< Real > & velocity_sensitivity(unsigned int prev_iter=0) const
virtual void sensitivity_solve(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly, const MAST::FunctionBase &p, bool if_assemble_jacobian=true)
Solves the sensitivity problem for the provided parameter.
TransientSolverBase(unsigned int o, unsigned int n)
constructor requires the number of iterations to store for the derived solver.
virtual void update_sensitivity_velocity(libMesh::NumericVector< Real > &vel, const libMesh::NumericVector< Real > &sol)=0
update the transient sensitivity velocity based on the current sensitivity solution ...
void solve_highest_derivative_and_advance_time_step_with_sensitivity(MAST::AssemblyBase &assembly, const MAST::FunctionBase &f)
solves for the sensitivity of highest derivative and advances the time-step.
virtual bool sensitivity_assemble(const libMesh::NumericVector< Real > &X, bool if_localize_sol, const MAST::FunctionBase &f, libMesh::NumericVector< Real > &sensitivity_rhs, bool close_vector=true)
Assembly function.
virtual void clear_elem_operation_object()
Clears the assembly elem operations object.
virtual void solve(MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly)
solves the nonlinear problem with the specified assembly operation object
MAST::SystemInitialization * _system