9 #ifndef Tempus_IntegratorPseudoTransientAdjointSensitivity_decl_hpp
10 #define Tempus_IntegratorPseudoTransientAdjointSensitivity_decl_hpp
13 #include "Tempus_IntegratorBasic.hpp"
14 #include "Tempus_AdjointSensitivityModelEvaluator.hpp"
50 template<
class Scalar>
75 Teuchos::RCP<Teuchos::ParameterList> pList,
76 const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model);
80 const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
81 std::string stepperType);
96 virtual bool advanceTime(
const Scalar timeFinal)
override;
98 virtual Scalar
getTime()
const override;
100 virtual int getIndex()
const override;
104 virtual Teuchos::RCP<Stepper<Scalar> >
getStepper()
const override;
109 virtual Teuchos::RCP<const SolutionHistory<Scalar> >
getSolutionHistory()
const override;
111 virtual Teuchos::RCP<const TimeStepControl<Scalar> >
getTimeStepControl()
const override;
124 Teuchos::RCP<
const Thyra::VectorBase<Scalar> > x0,
125 Teuchos::RCP<
const Thyra::VectorBase<Scalar> > xdot0 = Teuchos::null,
126 Teuchos::RCP<
const Thyra::VectorBase<Scalar> > xdotdot0 = Teuchos::null,
127 Teuchos::RCP<
const Thyra::MultiVectorBase<Scalar> > y0 = Teuchos::null,
128 Teuchos::RCP<
const Thyra::MultiVectorBase<Scalar> > ydot0 = Teuchos::null,
129 Teuchos::RCP<
const Thyra::MultiVectorBase<Scalar> > ydotdot0 = Teuchos::null);
132 virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> >
getX()
const;
134 virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> >
getXdot()
const;
136 virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> >
getXdotdot()
const;
139 virtual Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> >
getDgDp()
const;
155 void describe(Teuchos::FancyOStream & out,
156 const Teuchos::EVerbosityLevel verbLevel)
const override;
160 typedef Thyra::DefaultMultiVectorProductVector<Scalar>
DMVPV;
163 Teuchos::RCP<AdjointSensitivityModelEvaluator<Scalar> >
165 const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
166 const Teuchos::RCP<Teuchos::ParameterList>& inputPL);
170 Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >
model_;
171 Teuchos::RCP<AdjointSensitivityModelEvaluator<Scalar> >
sens_model_;
179 template<
class Scalar>
180 Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> >
182 Teuchos::RCP<Teuchos::ParameterList> pList,
183 const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model);
186 template<
class Scalar>
187 Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> >
189 const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
190 std::string stepperType);
193 template<
class Scalar>
194 Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> >
199 #endif // Tempus_IntegratorPseudoTransientAdjointSensitivity_decl_hpp
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory_
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl) override
virtual Scalar getTime() const override
Get current time.
virtual Teuchos::RCP< Teuchos::Time > getIntegratorTimer() const override
Returns the IntegratorTimer_ for this Integrator.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
virtual void setTempusParameterList(Teuchos::RCP< Teuchos::ParameterList > pl) override
Teuchos::RCP< IntegratorBasic< Scalar > > sens_integrator_
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList() override
std::string description() const override
virtual Teuchos::RCP< Teuchos::Time > getStepperTimer() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return adjoint sensitivity stored in gradient format.
Time integrator suitable for pseudotransient adjoint sensitivity analysis.
Teuchos::RCP< Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar > > integratorPseudoTransientAdjointSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
Non-member constructor.
Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > model_
virtual Teuchos::RCP< TimeStepControl< Scalar > > getNonConstTimeStepControl() override
Status
Status for the Integrator, the Stepper and the SolutionState.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdot() const
Get current the time derivative of the solution, xdot.
virtual void initializeSolutionHistory(Scalar t0, Teuchos::RCP< const Thyra::VectorBase< Scalar > > x0, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xdot0=Teuchos::null, Teuchos::RCP< const Thyra::VectorBase< Scalar > > xdotdot0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > y0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > ydot0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > ydotdot0=Teuchos::null)
Set the initial state from Thyra::VectorBase(s)
virtual ~IntegratorPseudoTransientAdjointSensitivity()
Destructor.
IntegratorPseudoTransientAdjointSensitivity()
Destructor.
Teuchos::RCP< AdjointSensitivityModelEvaluator< Scalar > > createSensitivityModel(const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL)
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdotdot() const
Get current the second time derivative of the solution, xdotdot.
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList() override
virtual Teuchos::RCP< Teuchos::ParameterList > getTempusParameterList() override
Return a copy of the Tempus ParameterList.
virtual int getIndex() const override
Get current index.
Thyra Base interface for time integrators. Time integrators are designed to advance the solution from...
virtual Teuchos::RCP< Stepper< Scalar > > getStepper() const override
Get the Stepper.
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory() const override
Get the SolutionHistory.
virtual Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl() const override
Get the TimeStepControl.
virtual bool advanceTime()
Advance the solution to timeMax, and return true if successful.
void buildSolutionHistory()
Teuchos::RCP< AdjointSensitivityModelEvaluator< Scalar > > sens_model_
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get current the solution, x.
Teuchos::RCP< IntegratorBasic< Scalar > > state_integrator_
Thyra::DefaultMultiVectorProductVector< Scalar > DMVPV
virtual Status getStatus() const override
Get Status.
Teuchos::RCP< DMVPV > dgdp_