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 : #include "AdjointTransientSolve.h" 11 : 12 : #include "FEProblemBase.h" 13 : #include "NonlinearSystemBase.h" 14 : #include "MooseUtils.h" 15 : 16 : #include "libmesh/sparse_matrix.h" 17 : #include "libmesh/numeric_vector.h" 18 : 19 : using namespace libMesh; 20 : 21 : InputParameters 22 126 : AdjointTransientSolve::validParams() 23 : { 24 126 : InputParameters params = AdjointSolve::validParams(); 25 126 : return params; 26 : } 27 : 28 63 : AdjointTransientSolve::AdjointTransientSolve(Executioner & ex) 29 : : AdjointSolve(ex), 30 : Restartable(this, "Executioners"), 31 126 : _old_time_residual(_nl_adjoint.getResidualTimeVector()), 32 126 : _forward_solutions(declareRecoverableData<std::vector<std::string>>("forward_solutions")) 33 : { 34 63 : } 35 : 36 : bool 37 3384 : AdjointTransientSolve::solve() 38 : { 39 3384 : bool converged = AdjointSolve::solve(); 40 : 41 : // Gather the contribution of this timestep to add to the next solve's source 42 3384 : if (converged) 43 3384 : evaluateTimeResidual(_nl_adjoint.solution(), _old_time_residual); 44 : 45 3384 : return converged; 46 : } 47 : 48 : void 49 3741 : AdjointTransientSolve::insertForwardSolution(int tstep) 50 : { 51 : // Time step should not be negative 52 3741 : if (tstep < 0) 53 0 : mooseError("Negative time step occurred."); 54 3741 : auto t_step = cast_int<std::size_t>(tstep); // Avoid compiler warnings 55 : // Should not be inserting a time greater the last one inserted 56 3741 : if (t_step > _forward_solutions.size()) 57 0 : mooseError("Trying to insert a solution at a time-step greater than one past the previously " 58 : "inserted time step. Previous time step = ", 59 0 : (int)_forward_solutions.size() - 1, 60 : ", inserted time step = ", 61 : t_step, 62 : "."); 63 : 64 : // Add vector name to member variable if this time has not occurred 65 3741 : if (t_step == _forward_solutions.size()) 66 7482 : _forward_solutions.push_back(getForwardSolutionName(tstep)); 67 : 68 : // Get/add vector from/to adjoint system 69 3741 : auto & solution = _nl_adjoint.addVector(_forward_solutions[t_step], false, PARALLEL); 70 : 71 : // Set the vector to the inserted solution 72 3741 : solution = _nl_forward.solution(); 73 3741 : } 74 : 75 : void 76 3384 : AdjointTransientSolve::setForwardSolution(int tstep) 77 : { 78 : // Make sure the time step was saved 79 3384 : if (tstep < 0 || tstep >= cast_int<int>(_forward_solutions.size())) 80 0 : mooseError("Could not find forward solution at time step ", tstep, "."); 81 : auto t_step = cast_int<std::size_t>(tstep); // Avoid compiler warnings 82 : 83 : // Copy the solutions to states that exist in the system 84 : unsigned int state = 0; 85 10152 : while (_nl_forward.hasSolutionState(state)) 86 : { 87 6768 : const auto & name = t_step > state ? _forward_solutions[t_step - state] : _forward_solutions[0]; 88 6768 : _nl_forward.solutionState(state) = _nl_adjoint.getVector(name); 89 6768 : state++; 90 : } 91 : 92 3384 : _nl_forward.update(); 93 3384 : } 94 : 95 : void 96 3384 : AdjointTransientSolve::assembleAdjointSystem(SparseMatrix<Number> & matrix, 97 : const NumericVector<Number> & solution, 98 : NumericVector<Number> & rhs) 99 : { 100 : // Assemble the steady-state version of the adjoint problem 101 3384 : AdjointSolve::assembleAdjointSystem(matrix, solution, rhs); 102 : // Add the contribution from old solutions 103 3384 : rhs += _old_time_residual; 104 3384 : } 105 : 106 : void 107 3384 : AdjointTransientSolve::evaluateTimeResidual(const NumericVector<Number> & solution, 108 : NumericVector<Number> & residual) 109 : { 110 : // This tag should exist, but the matrix might not necessarily be added 111 3384 : auto time_matrix_tag = _problem.getMatrixTagID("TIME"); 112 : // Use the adjoint system matrix to hold the time Jacobian 113 3384 : auto & time_matrix = static_cast<ImplicitSystem &>(_nl_adjoint.system()).get_system_matrix(); 114 : 115 : // Make sure we tell the problem which system we are evaluating 116 3384 : _problem.setCurrentNonlinearSystem(_forward_sys_num); 117 : // Accumulate the time part of the Jacobian into the adjoint matrix 118 3384 : _problem.computeJacobianTag(*_nl_forward.currentSolution(), time_matrix, time_matrix_tag); 119 : 120 : // Perform multiplication of the adjoint solution on the transpose of the time Jacobian 121 3384 : residual.zero(); 122 3384 : residual.add_vector_transpose(solution, time_matrix); 123 3384 : }