Line data Source code
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 : 15 : /** 16 : * Transient adjoint solve object. @see TransientAndAdjoint for example usage. General algorithm: 17 : * 18 : * 1. Set forward initial conditions 19 : * 2. Save the forward initial solution using insertForwardSolution(0) 20 : * 3. for t_step = 1,...,num_step 21 : * 4. Solve forward time step 22 : * 5. Save forward solution using insertForwardSolution(t_step) 23 : * 6. Set _old_time_residual = 0 24 : * 6. for t_step = num_step,...,1 25 : * 7. Set forward solution using setForwardSolution(t_step) 26 : * 8. Solve adjoint time step: A* u* = b* + _old_time_residual 27 : * 9. Evaluate time residual: _old_time_residual = A_t* u* 28 : * 29 : * where A* is the transpose of the linearized forward operator at the specified timestep, u* is the 30 : * adjoint solution, and A_t* is the transpose of the linearized forward time operator. 31 : */ 32 : class AdjointTransientSolve : public AdjointSolve, public Restartable 33 : { 34 : public: 35 : AdjointTransientSolve(Executioner & ex); 36 : 37 : static InputParameters validParams(); 38 : 39 : /** 40 : * Overriding parent class so the time-derivative residual is stored for the next time step 41 : */ 42 : virtual bool solve() override; 43 : 44 : /** 45 : * This function should be called after every converged forward time step. It adds a new vector 46 : * (if necessary) to the adjoint system representing the forward solution at the given time step. 47 : * 48 : * @param tstep Time step of the forward solution in which to store in _forward_solutions 49 : */ 50 : void insertForwardSolution(int tstep); 51 : 52 : /** 53 : * Takes the previously saved forward solutions residing in the adjoint system and copies them to 54 : * the available solution states in the forward systems. This should be called at each adjoint 55 : * time step. 56 : * 57 : * @param tstep The forward time step index being evaluated during the adjoint solve. 58 : */ 59 : void setForwardSolution(int tstep); 60 : 61 : protected: 62 : /** 63 : * Overriding parent class so the previous time-derivative residual is added to the 64 : * right-hand-side of the linear solve 65 : */ 66 : virtual void assembleAdjointSystem(SparseMatrix<Number> & matrix, 67 : const NumericVector<Number> & solution, 68 : NumericVector<Number> & rhs) override; 69 : 70 : /** 71 : * This evaluates the time part of the adjoint residual. This is used to accumulate the source 72 : * contribution from previous adjoint time steps for the next adjoint time step. It works by 73 : * evaluating the Jacobian of the time-derivative term on the forward system and multiplying the 74 : * transpose by the current adjoint solution. 75 : * 76 : * @param solution Current adjoint solution 77 : * @param residual Vector to save the result into 78 : */ 79 : void evaluateTimeResidual(const NumericVector<Number> & solution, 80 : NumericVector<Number> & residual); 81 : 82 : /** 83 : * Prescribed name of the forward solution at a specified time step 84 : * 85 : * @param tstep The forward time step 86 : * @return std::string The name of the vector representing the forward solution on the adjoint 87 : * system 88 : */ 89 5705 : static std::string getForwardSolutionName(int tstep) 90 : { 91 11410 : return "forward_solution" + std::to_string(tstep); 92 : } 93 : 94 : private: 95 : /// The residual contribution from previous adjoint solutions 96 : NumericVector<Number> & _old_time_residual; 97 : /// Name of the forward solution at each time step residing in the adjoint system 98 : std::vector<std::string> & _forward_solutions; 99 : };