PISM, A Parallel Ice Sheet Model 2.3.0-79cae578d committed by Constantine Khrulev on 2026-03-22
Loading...
Searching...
No Matches
IP_SSATaucTaoTikhonovProblemLCL.cc
Go to the documentation of this file.
1// Copyright (C) 2012, 2014, 2015, 2016, 2017, 2020, 2021, 2022, 2023, 2025 David Maxwell and Constantine Khroulev
2//
3// This file is part of PISM.
4//
5// PISM is free software; you can redistribute it and/or modify it under the
6// terms of the GNU General Public License as published by the Free Software
7// Foundation; either version 3 of the License, or (at your option) any later
8// version.
9//
10// PISM is distributed in the hope that it will be useful, but WITHOUT ANY
11// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12// FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
13// details.
14//
15// You should have received a copy of the GNU General Public License
16// along with PISM; if not, write to the Free Software
17// Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18
19#include "pism/inverse/IP_SSATaucTaoTikhonovProblemLCL.hh"
20#include "pism/util/Grid.hh"
21#include "pism/util/Config.hh"
22#include "pism/util/Context.hh"
23#include "pism/util/petscwrappers/DM.hh"
24
25#include <memory>
26
27namespace pism {
28namespace inverse {
29
32
33// typedef TikhonovProblemListener<InverseProblem> Listener;
34// typedef typename Listener::Ptr ListenerPtr;
35
39 double eta,
40 IPFunctional<DesignVec> &designFunctional,
41 IPFunctional<StateVec> &stateFunctional)
42: m_ssaforward(ssaforward),
43 m_dGlobal(d0.grid(), "design variable (global)"),
44 m_d0(d0),
45 m_dzeta(d0.grid(), "dzeta"),
46 m_u(d0.grid(), "state variable"),
47 m_du(d0.grid(), "du"),
48 m_u_obs(u_obs),
49 m_eta(eta),
50 m_d_Jdesign(d0.grid(), "Jdesign design variable"),
51 m_u_Jdesign(d0.grid(), "Jdesign state variable"),
52 m_designFunctional(designFunctional),
53 m_stateFunctional(stateFunctional)
54{
55
56 PetscErrorCode ierr;
57 std::shared_ptr<const Grid> grid = m_d0.grid();
58
59 double stressScale = grid->ctx()->config()->get_number("inverse.design.param_tauc_scale");
60 m_constraintsScale = grid->Lx()*grid->Ly()*4*stressScale;
61
62 m_velocityScale = grid->ctx()->config()->get_number("inverse.ssa.velocity_scale", "m second-1");
63
64 m_d.reset(new DesignVecGhosted(grid, "design variable"));
65
67
68 m_uGlobal.reset(new StateVec(grid, "state variable (global)"));
69
70 m_u_diff.reset(new StateVec1(grid, "state residual"));
71
72 m_d_diff.reset(new DesignVecGhosted(grid, "design residual"));
73
74 m_grad_state.reset(new StateVec(grid, "state gradient"));
75
76 m_grad_design.reset(new DesignVec(grid, "design gradient"));
77
78 m_constraints.reset(new StateVec(grid, "PDE constraints"));
79
81
82 ierr = DMSetMatType(da, MATBAIJ);
83 PISM_CHK(ierr, "DMSetMatType");
84
85 ierr = DMCreateMatrix(da, m_Jstate.rawptr());
86 PISM_CHK(ierr, "DMCreateMatrix");
87
88 int nLocalNodes = grid->xm()*grid->ym();
89 int nGlobalNodes = grid->Mx()*grid->My();
90 ierr = MatCreateShell(grid->com, 2*nLocalNodes, nLocalNodes, 2*nGlobalNodes, nGlobalNodes,
91 this, m_Jdesign.rawptr());
92 PISM_CHK(ierr, "MatCreateShell");
93
94 ierr = MatShellSetOperation(m_Jdesign, MATOP_MULT,
95 (void(*)(void))jacobian_design_callback);
96 PISM_CHK(ierr, "MatShellSetOperation");
97
98 ierr = MatShellSetOperation(m_Jdesign, MATOP_MULT_TRANSPOSE,
100 PISM_CHK(ierr, "MatShellSetOperation");
101
102 m_x.reset(new IPTwoBlockVec(m_dGlobal.vec(),m_uGlobal->vec()));
103}
104
108
109std::shared_ptr<IP_SSATaucTaoTikhonovProblemLCL::StateVec>
111
112 m_x->scatterToB(m_uGlobal->vec());
114
115 return m_uGlobal;
116}
117
118std::shared_ptr<IP_SSATaucTaoTikhonovProblemLCL::DesignVec>
120 m_x->scatterToA(m_d->vec());
121 return m_d;
122}
123
125 PetscErrorCode ierr;
126 ierr = TaoSetStateDesignIS(tao,
127 m_x->blockBIndexSet() /*state*/,
128 m_x->blockAIndexSet() /*design*/);
129 PISM_CHK(ierr, "TaoSetStateDesignIS");
130
133
135 m_constraints->vec(),
137
139}
140
142 PetscErrorCode ierr;
143
144 // Has to be a PetscInt because of the TaoGetSolutionStatus call.
145 PetscInt its;
146 ierr = TaoGetSolutionStatus(tao, &its, NULL, NULL, NULL, NULL, NULL);
147 PISM_CHK(ierr, "TaoGetSolutionStatus");
148
149 int nListeners = m_listeners.size();
150 for (int k = 0; k < nListeners; k++) {
151 m_listeners[k]->iteration(*this, m_eta,
155 m_u_diff,
158 }
159}
160
162 double *value, Vec gradient) {
163
164 m_x->scatter(x,m_dGlobal.vec(),m_uGlobal->vec());
166
167 // Variable 'm_dGlobal' has no ghosts. We need ghosts for computation with the design variable.
168 m_d->copy_from(m_dGlobal);
169
170 m_d_diff->copy_from(*m_d);
171 m_d_diff->add(-1,m_d0);
173 m_grad_design->scale(1/m_eta);
174
175 m_u_diff->copy_from(*m_uGlobal);
176 m_u_diff->add(-1, m_u_obs);
179
180 m_x->gather(m_grad_design->vec(), m_grad_state->vec(), gradient);
181
184
185 *value = m_val_design / m_eta + m_val_state;
186}
187
188std::shared_ptr<TerminationReason> IP_SSATaucTaoTikhonovProblemLCL::formInitialGuess(Vec *x) {
189 m_d->copy_from(m_dGlobal);
190 auto reason = m_ssaforward.linearize_at(*m_d);
191 if (reason->failed()) {
192 return reason;
193 }
194
195 m_uGlobal->copy_from(*m_ssaforward.solution());
196 m_uGlobal->scale(1.0 / m_velocityScale);
197
198 m_x->gather(m_dGlobal.vec(), m_uGlobal->vec());
199
200 // This is probably irrelevant.
202
203 *x = *m_x;
205}
206
208 PetscErrorCode ierr;
209
210 m_x->scatter(x,m_dGlobal.vec(),m_uGlobal->vec());
212
213 m_d->copy_from(m_dGlobal);
215
217
219
220 ierr = VecScale(r,1./m_constraintsScale);
221 PISM_CHK(ierr, "VecScale");
222}
223
225 Mat Jstate,
226 Mat /*Jpc*/,
227 Mat /*Jinv*/,
228 MatStructure *s) {
229 PetscErrorCode ierr;
230
231 m_x->scatter(x, m_dGlobal.vec(), m_uGlobal->vec());
233
234 m_d->copy_from(m_dGlobal);
236
238
240
241 *s = SAME_NONZERO_PATTERN;
242
243 ierr = MatScale(Jstate, m_velocityScale / m_constraintsScale);
244 PISM_CHK(ierr, "MatScale");
245}
246
248 // I'm not sure if the following are necessary (i.e. will the copies that happen
249 // in evaluateObjectiveAndGradient be sufficient) but we'll do them here
250 // just in case.
251 m_x->scatter(x,m_dGlobal.vec(),m_uGlobal->vec());
255}
256
258
259 PetscErrorCode ierr;
260 {
261 ierr = DMGlobalToLocalBegin(*m_dzeta.dm(), x, INSERT_VALUES, m_dzeta.vec());
262 PISM_CHK(ierr, "DMGlobalToLocalBegin");
263
264 ierr = DMGlobalToLocalEnd(*m_dzeta.dm(), x, INSERT_VALUES, m_dzeta.vec());
265 PISM_CHK(ierr, "DMGlobalToLocalEnd");
266 }
267
269
271
272 ierr = VecScale(y, 1.0 / m_constraintsScale); PISM_CHK(ierr, "VecScale");
273}
274
276
277 PetscErrorCode ierr;
278 {
279 ierr = DMGlobalToLocalBegin(*m_du.dm(), x, INSERT_VALUES, m_du.vec());
280 PISM_CHK(ierr, "DMGlobalToLocalBegin");
281
282 ierr = DMGlobalToLocalEnd(*m_du.dm(), x, INSERT_VALUES, m_du.vec());
283 PISM_CHK(ierr, "DMGlobalToLocalEnd");
284 }
285
287
289
290 ierr = VecScale(y, 1.0 / m_constraintsScale); PISM_CHK(ierr, "VecScale");
291}
292
294 try {
296 PetscErrorCode ierr = MatShellGetContext(A,&ctx);
297 PISM_CHK(ierr, "MatShellGetContext");
298
300 } catch (...) {
301 MPI_Comm com = MPI_COMM_SELF;
302 PetscErrorCode ierr = PetscObjectGetComm((PetscObject)A, &com); CHKERRQ(ierr);
304 SETERRQ(com, 1, "A PISM callback failed");
305 }
306 return 0;
307}
308
310 Vec y) {
311 try {
313 PetscErrorCode ierr = MatShellGetContext(A,&ctx);
314 PISM_CHK(ierr, "MatShellGetContext");
315
317 } catch (...) {
318 MPI_Comm com = MPI_COMM_SELF;
319 PetscErrorCode ierr = PetscObjectGetComm((PetscObject)A, &com); CHKERRQ(ierr);
321 SETERRQ(com, 1, "A PISM callback failed");
322 }
323 return 0;
324}
325
326} // end of namespace inverse
327} // end of namespace pism
static std::shared_ptr< TerminationReason > success()
T * rawptr()
Definition Wrapper.hh:34
void copy_from(const Array2D< T > &source)
Definition Array2D.hh:101
petsc::Vec & vec() const
Definition Array.cc:313
std::shared_ptr< const Grid > grid() const
Definition Array.cc:134
std::shared_ptr< petsc::DM > dm() const
Definition Array.cc:327
virtual void valueAt(IMVecType &x, double *OUTPUT)=0
Computes the value of the functional at the vector x.
virtual void gradientAt(IMVecType &x, IMVecType &gradient)=0
Computes the gradient of the functional at the vector x.
Abstract base class for functions from ice model vectors to .
virtual void apply_jacobian_design(array::Vector &u, array::Scalar &dzeta, array::Vector &du)
Applies the design Jacobian matrix to a perturbation of the design variable.
virtual void apply_jacobian_design_transpose(array::Vector &u, array::Vector &du, array::Scalar &dzeta)
Applies the transpose of the design Jacobian matrix to a perturbation of the state variable.
virtual std::shared_ptr< array::Vector > solution()
Returns the last solution of the SSA as computed by linearize_at.
virtual void assemble_residual(array::Vector &u, array::Vector &R)
Computes the residual function as defined in the class-level documentation.
virtual std::shared_ptr< TerminationReason > linearize_at(array::Scalar &zeta)
Sets the current value of the design variable and solves the SSA to find the associated .
petsc::DM & get_da() const
Exposes the DMDA of the underlying grid for the benefit of TAO.
virtual void set_design(array::Scalar &zeta)
Sets the current value of of the design paramter .
virtual void assemble_jacobian_state(array::Vector &u, Mat J)
Assembles the state Jacobian matrix.
Implements the forward problem of the map taking to the corresponding solution of the SSA.
virtual void evaluateObjectiveAndGradient(Tao tao, Vec x, double *value, Vec gradient)
static PetscErrorCode jacobian_design_transpose_callback(Mat A, Vec x, Vec y)
virtual void evaluateConstraintsJacobianState(Tao tao, Vec x, Mat Jstate, Mat Jpc, Mat Jinv, MatStructure *s)
virtual std::shared_ptr< TerminationReason > formInitialGuess(Vec *x)
IP_SSATaucTaoTikhonovProblemLCL(IP_SSATaucForwardProblem &ssaforward, DesignVec &d0, StateVec &u_obs, double eta, IPFunctional< DesignVec > &designFunctional, IPFunctional< StateVec > &stateFunctional)
virtual void evaluateConstraintsJacobianDesign(Tao tao, Vec x, Mat Jdesign)
static PetscErrorCode jacobian_design_callback(Mat A, Vec x, Vec y)
Defines a Tikhonov minimization problem of determining from SSA velocities to be solved with a TaoBa...
static void connect(Tao tao, Problem &p, Vec c, Mat Jc, Mat Jd, Mat Jcpc=NULL, Mat Jcinv=NULL)
Definition TaoUtil.hh:463
static void connect(Tao tao, Problem &p)
Definition TaoUtil.hh:232
Adaptor to connect a TAO objective and gradient function callback to a C++ object method.
Definition TaoUtil.hh:411
#define PISM_CHK(errcode, name)
static const double k
Definition exactTestP.cc:42
void handle_fatal_errors(MPI_Comm com)