PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
IP_SSATaucTaoTikhonovProblemLCL.cc
Go to the documentation of this file.
1 // Copyright (C) 2012, 2014, 2015, 2016, 2017, 2020, 2021, 2022, 2023 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/ConfigInterface.hh"
22 #include "pism/util/Context.hh"
23 #include <memory>
24 
25 namespace pism {
26 namespace inverse {
27 
30 
31 // typedef TikhonovProblemListener<InverseProblem> Listener;
32 // typedef typename Listener::Ptr ListenerPtr;
33 
37  double eta,
38  IPFunctional<DesignVec> &designFunctional,
39  IPFunctional<StateVec> &stateFunctional)
40 : m_ssaforward(ssaforward),
41  m_dGlobal(d0.grid(), "design variable (global)"),
42  m_d0(d0),
43  m_dzeta(d0.grid(), "dzeta"),
44  m_u(d0.grid(), "state variable"),
45  m_du(d0.grid(), "du"),
46  m_u_obs(u_obs),
47  m_eta(eta),
48  m_d_Jdesign(d0.grid(), "Jdesign design variable"),
49  m_u_Jdesign(d0.grid(), "Jdesign state variable"),
50  m_designFunctional(designFunctional),
51  m_stateFunctional(stateFunctional)
52 {
53 
54  PetscErrorCode ierr;
55  std::shared_ptr<const Grid> grid = m_d0.grid();
56 
57  double stressScale = grid->ctx()->config()->get_number("inverse.design.param_tauc_scale");
58  m_constraintsScale = grid->Lx()*grid->Ly()*4*stressScale;
59 
60  m_velocityScale = grid->ctx()->config()->get_number("inverse.ssa.velocity_scale", "m second-1");
61 
62  m_d.reset(new DesignVecGhosted(grid, "design variable"));
63 
65 
66  m_uGlobal.reset(new StateVec(grid, "state variable (global)"));
67 
68  m_u_diff.reset(new StateVec1(grid, "state residual"));
69 
70  m_d_diff.reset(new DesignVecGhosted(grid, "design residual"));
71 
72  m_grad_state.reset(new StateVec(grid, "state gradient"));
73 
74  m_grad_design.reset(new DesignVec(grid, "design gradient"));
75 
76  m_constraints.reset(new StateVec(grid, "PDE constraints"));
77 
79 
80  ierr = DMSetMatType(da, MATBAIJ);
81  PISM_CHK(ierr, "DMSetMatType");
82 
83  ierr = DMCreateMatrix(da, m_Jstate.rawptr());
84  PISM_CHK(ierr, "DMCreateMatrix");
85 
86  int nLocalNodes = grid->xm()*grid->ym();
87  int nGlobalNodes = grid->Mx()*grid->My();
88  ierr = MatCreateShell(grid->com, 2*nLocalNodes, nLocalNodes, 2*nGlobalNodes, nGlobalNodes,
89  this, m_Jdesign.rawptr());
90  PISM_CHK(ierr, "MatCreateShell");
91 
92  ierr = MatShellSetOperation(m_Jdesign, MATOP_MULT,
93  (void(*)(void))jacobian_design_callback);
94  PISM_CHK(ierr, "MatShellSetOperation");
95 
96  ierr = MatShellSetOperation(m_Jdesign, MATOP_MULT_TRANSPOSE,
97  (void(*)(void))jacobian_design_transpose_callback);
98  PISM_CHK(ierr, "MatShellSetOperation");
99 
100  m_x.reset(new IPTwoBlockVec(m_dGlobal.vec(),m_uGlobal->vec()));
101 }
102 
104  m_dGlobal.copy_from(d0);
105 }
106 
107 std::shared_ptr<IP_SSATaucTaoTikhonovProblemLCL::StateVec>
109 
110  m_x->scatterToB(m_uGlobal->vec());
111  m_uGlobal->scale(m_velocityScale);
112 
113  return m_uGlobal;
114 }
115 
116 std::shared_ptr<IP_SSATaucTaoTikhonovProblemLCL::DesignVec>
118  m_x->scatterToA(m_d->vec());
119  return m_d;
120 }
121 
123  PetscErrorCode ierr;
124  ierr = TaoSetStateDesignIS(tao,
125  m_x->blockBIndexSet() /*state*/,
126  m_x->blockAIndexSet() /*design*/);
127  PISM_CHK(ierr, "TaoSetStateDesignIS");
128 
131 
133  m_constraints->vec(),
135 
137 }
138 
140  PetscErrorCode ierr;
141 
142  // Has to be a PetscInt because of the TaoGetSolutionStatus call.
143  PetscInt its;
144  ierr = TaoGetSolutionStatus(tao, &its, NULL, NULL, NULL, NULL, NULL);
145  PISM_CHK(ierr, "TaoGetSolutionStatus");
146 
147  int nListeners = m_listeners.size();
148  for (int k = 0; k < nListeners; k++) {
149  m_listeners[k]->iteration(*this, m_eta,
153  m_u_diff,
154  m_grad_state,
155  m_constraints);
156  }
157 }
158 
160  double *value, Vec gradient) {
161 
162  m_x->scatter(x,m_dGlobal.vec(),m_uGlobal->vec());
163  m_uGlobal->scale(m_velocityScale);
164 
165  // Variable 'm_dGlobal' has no ghosts. We need ghosts for computation with the design variable.
166  m_d->copy_from(m_dGlobal);
167 
168  m_d_diff->copy_from(*m_d);
169  m_d_diff->add(-1,m_d0);
171  m_grad_design->scale(1/m_eta);
172 
173  m_u_diff->copy_from(*m_uGlobal);
174  m_u_diff->add(-1, m_u_obs);
177 
178  m_x->gather(m_grad_design->vec(), m_grad_state->vec(), gradient);
179 
182 
183  *value = m_val_design / m_eta + m_val_state;
184 }
185 
186 std::shared_ptr<TerminationReason> IP_SSATaucTaoTikhonovProblemLCL::formInitialGuess(Vec *x) {
187  m_d->copy_from(m_dGlobal);
188  auto reason = m_ssaforward.linearize_at(*m_d);
189  if (reason->failed()) {
190  return reason;
191  }
192 
193  m_uGlobal->copy_from(*m_ssaforward.solution());
194  m_uGlobal->scale(1.0 / m_velocityScale);
195 
196  m_x->gather(m_dGlobal.vec(), m_uGlobal->vec());
197 
198  // This is probably irrelevant.
199  m_uGlobal->scale(m_velocityScale);
200 
201  *x = *m_x;
203 }
204 
206  PetscErrorCode ierr;
207 
208  m_x->scatter(x,m_dGlobal.vec(),m_uGlobal->vec());
209  m_uGlobal->scale(m_velocityScale);
210 
211  m_d->copy_from(m_dGlobal);
213 
215 
217 
218  ierr = VecScale(r,1./m_constraintsScale);
219  PISM_CHK(ierr, "VecScale");
220 }
221 
223  Mat Jstate,
224  Mat /*Jpc*/,
225  Mat /*Jinv*/,
226  MatStructure *s) {
227  PetscErrorCode ierr;
228 
229  m_x->scatter(x, m_dGlobal.vec(), m_uGlobal->vec());
230  m_uGlobal->scale(m_velocityScale);
231 
232  m_d->copy_from(m_dGlobal);
234 
236 
238 
239  *s = SAME_NONZERO_PATTERN;
240 
241  ierr = MatScale(Jstate, m_velocityScale / m_constraintsScale);
242  PISM_CHK(ierr, "MatScale");
243 }
244 
246  // I'm not sure if the following are necessary (i.e. will the copies that happen
247  // in evaluateObjectiveAndGradient be sufficient) but we'll do them here
248  // just in case.
249  m_x->scatter(x,m_dGlobal.vec(),m_uGlobal->vec());
250  m_uGlobal->scale(m_velocityScale);
253 }
254 
256 
257  PetscErrorCode ierr;
258  {
259  ierr = DMGlobalToLocalBegin(*m_dzeta.dm(), x, INSERT_VALUES, m_dzeta.vec());
260  PISM_CHK(ierr, "DMGlobalToLocalBegin");
261 
262  ierr = DMGlobalToLocalEnd(*m_dzeta.dm(), x, INSERT_VALUES, m_dzeta.vec());
263  PISM_CHK(ierr, "DMGlobalToLocalEnd");
264  }
265 
267 
269 
270  ierr = VecScale(y, 1.0 / m_constraintsScale); PISM_CHK(ierr, "VecScale");
271 }
272 
274 
275  PetscErrorCode ierr;
276  {
277  ierr = DMGlobalToLocalBegin(*m_du.dm(), x, INSERT_VALUES, m_du.vec());
278  PISM_CHK(ierr, "DMGlobalToLocalBegin");
279 
280  ierr = DMGlobalToLocalEnd(*m_du.dm(), x, INSERT_VALUES, m_du.vec());
281  PISM_CHK(ierr, "DMGlobalToLocalEnd");
282  }
283 
285 
287 
288  ierr = VecScale(y, 1.0 / m_constraintsScale); PISM_CHK(ierr, "VecScale");
289 }
290 
292  try {
294  PetscErrorCode ierr = MatShellGetContext(A,&ctx);
295  PISM_CHK(ierr, "MatShellGetContext");
296 
298  } catch (...) {
299  MPI_Comm com = MPI_COMM_SELF;
300  PetscErrorCode ierr = PetscObjectGetComm((PetscObject)A, &com); CHKERRQ(ierr);
301  handle_fatal_errors(com);
302  SETERRQ(com, 1, "A PISM callback failed");
303  }
304  return 0;
305 }
306 
308  Vec y) {
309  try {
311  PetscErrorCode ierr = MatShellGetContext(A,&ctx);
312  PISM_CHK(ierr, "MatShellGetContext");
313 
315  } catch (...) {
316  MPI_Comm com = MPI_COMM_SELF;
317  PetscErrorCode ierr = PetscObjectGetComm((PetscObject)A, &com); CHKERRQ(ierr);
318  handle_fatal_errors(com);
319  SETERRQ(com, 1, "A PISM callback failed");
320  }
321  return 0;
322 }
323 
324 } // end of namespace inverse
325 } // end of namespace pism
static std::shared_ptr< TerminationReason > success()
T * rawptr()
Definition: Wrapper.hh:39
void copy_from(const Array2D< T > &source)
Definition: Array2D.hh:73
petsc::Vec & vec() const
Definition: Array.cc:339
std::shared_ptr< const Grid > grid() const
Definition: Array.cc:132
std::shared_ptr< petsc::DM > dm() const
Definition: Array.cc:353
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 .
Definition: IPFunctional.hh:41
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.
petsc::DM & get_da() const
Exposes the DMDA of the underlying grid for the benefit of TAO.
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 .
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:459
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:407
#define PISM_CHK(errcode, name)
static const double k
Definition: exactTestP.cc:42
void handle_fatal_errors(MPI_Comm com)