PolyFEM
Loading...
Searching...
No Matches
CompositeForm.hpp
Go to the documentation of this file.
1#pragma once
2
7
8#include <Eigen/Core>
9
10#include <memory>
11#include <vector>
12
13namespace polyfem::solver
14{
16 {
18
19 public:
20 CompositeForm(const VariableToSimulationGroup &variable_to_simulations, const std::vector<std::shared_ptr<AdjointForm>> &forms) : AdjointForm(variable_to_simulations), forms_(forms) {}
21 CompositeForm(const std::vector<std::shared_ptr<AdjointForm>> &forms) : AdjointForm(forms[0]->get_variable_to_simulations()), forms_(forms) {}
22 virtual ~CompositeForm() {}
23
24 virtual int n_objs() const final { return forms_.size(); }
25
26 virtual Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override final;
27 virtual void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override final;
28
29 virtual double compose(const Eigen::VectorXd &inputs) const = 0;
30 virtual Eigen::VectorXd compose_grad(const Eigen::VectorXd &inputs) const = 0;
31
32 Eigen::VectorXd get_inputs(const Eigen::VectorXd &x) const;
33
34 virtual double value_unweighted(const Eigen::VectorXd &x) const final override;
35 virtual void init(const Eigen::VectorXd &x) final override;
36 virtual bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override;
37 virtual double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override;
38 virtual void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) final override;
39 virtual void line_search_end() final override;
40 virtual void post_step(const polysolve::nonlinear::PostStepData &data) final override;
41 virtual void solution_changed(const Eigen::VectorXd &new_x) final override;
42 virtual bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override;
43
44 private:
45 std::vector<std::shared_ptr<AdjointForm>> forms_;
46 };
47} // namespace polyfem::solver
int x
Storage for additional data required by differntial code.
Definition DiffCache.hpp:21
main class that contains the polyfem solver and all its state
Definition State.hpp:113
const VariableToSimulationGroup & get_variable_to_simulations() const
virtual void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override final
std::vector< std::shared_ptr< AdjointForm > > forms_
CompositeForm(const std::vector< std::shared_ptr< AdjointForm > > &forms)
virtual bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override
Determine if a step from solution x0 to solution x1 is allowed.
virtual double compose(const Eigen::VectorXd &inputs) const =0
virtual Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override final
virtual void init(const Eigen::VectorXd &x) final override
Initialize the form.
virtual double value_unweighted(const Eigen::VectorXd &x) const final override
Compute the value of the form.
virtual double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override
Determine the maximum step size allowable between the current and next solution.
Eigen::VectorXd get_inputs(const Eigen::VectorXd &x) const
CompositeForm(const VariableToSimulationGroup &variable_to_simulations, const std::vector< std::shared_ptr< AdjointForm > > &forms)
virtual int n_objs() const final
virtual void post_step(const polysolve::nonlinear::PostStepData &data) final override
Update fields after a step in the optimization.
virtual Eigen::VectorXd compose_grad(const Eigen::VectorXd &inputs) const =0
virtual bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const final override
Checks if the step is collision free.
virtual void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) final override
Initialize variables used during the line search.
virtual void solution_changed(const Eigen::VectorXd &new_x) final override
Update cached fields upon a change in the solution.
virtual void line_search_end() final override
Clear variables used during the line search.
A collection of VariableToSimulation.