https://mooseframework.inl.gov
AdjointTransientSolve.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "AdjointSolve.h"
13 #include "Restartable.h"
14 
33 {
34 public:
36 
38 
42  virtual bool solve() override;
43 
50  void insertForwardSolution(int tstep);
51 
59  void setForwardSolution(int tstep);
60 
61 protected:
66  virtual void assembleAdjointSystem(SparseMatrix<Number> & matrix,
67  const NumericVector<Number> & solution,
68  NumericVector<Number> & rhs) override;
69 
79  void evaluateTimeResidual(const NumericVector<Number> & solution,
80  NumericVector<Number> & residual);
81 
89  static std::string getForwardSolutionName(int tstep)
90  {
91  return "forward_solution" + std::to_string(tstep);
92  }
93 
94 private:
98  std::vector<std::string> & _forward_solutions;
99 };
std::vector< std::string > & _forward_solutions
Name of the forward solution at each time step residing in the adjoint system.
Transient adjoint solve object.
AdjointTransientSolve(Executioner &ex)
static InputParameters validParams()
static std::string getForwardSolutionName(int tstep)
Prescribed name of the forward solution at a specified time step.
The solve object is responsible for solving the adjoint version of a forward model.
Definition: AdjointSolve.h:31
virtual void assembleAdjointSystem(SparseMatrix< Number > &matrix, const NumericVector< Number > &solution, NumericVector< Number > &rhs) override
Overriding parent class so the previous time-derivative residual is added to the right-hand-side of t...
void evaluateTimeResidual(const NumericVector< Number > &solution, NumericVector< Number > &residual)
This evaluates the time part of the adjoint residual.
void insertForwardSolution(int tstep)
This function should be called after every converged forward time step.
NumericVector< Number > & _old_time_residual
The residual contribution from previous adjoint solutions.
void setForwardSolution(int tstep)
Takes the previously saved forward solutions residing in the adjoint system and copies them to the av...
virtual bool solve() override
Overriding parent class so the time-derivative residual is stored for the next time step...