PISM, A Parallel Ice Sheet Model  stable v2.1.1 committed by Constantine Khrulev on 2024-12-04 13:36:58 -0900
ssa_test_const.cc
Go to the documentation of this file.
1 // Copyright (C) 2010--2018, 2021, 2022 Ed Bueler, Constantine Khroulev, and David Maxwell
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 /* This file implements a test case for the ssa: constant flow. The rheology is
20  nonlinear (i.e. n=3 in the Glen flow law) and the basal shear stress is a
21  nonlinear function of velocity (peseudo-plastic flow with parameter q
22  specified at runtime).
23 
24  The geometry consists of a constant surface slope in the positive
25  x-direction, and a constant velocity is specified as a Dirichlet condition
26  on the boundary that should lead to a constant solution in the interior.
27  Because the solution is constant, the nonzero terms in the SSA are only the
28  basal shear stress and the driving stress.
29  */
30 
31 static char help[] =
32  "\nSSA_TEST_CONST\n"
33  " Testing program for the finite element implementation of the SSA.\n"
34  " Does a time-independent calculation. Does not run IceModel or a derived\n"
35  " class thereof.Also may be used in a PISM\n"
36  " software (regression) test.\n\n";
37 
38 #include <cmath>
39 
40 #include "pism/basalstrength/basal_resistance.hh" // IceBasalResistancePlasticLaw
41 #include "pism/stressbalance/ssa/SSAFD.hh"
42 #include "pism/stressbalance/ssa/SSAFEM.hh"
43 #include "pism/stressbalance/ssa/SSATestCase.hh"
44 #include "pism/util/Mask.hh"
45 #include "pism/util/Context.hh"
46 #include "pism/util/VariableMetadata.hh"
47 #include "pism/util/error_handling.hh"
48 #include "pism/util/io/File.hh"
49 #include "pism/util/petscwrappers/PetscInitializer.hh"
50 #include "pism/util/pism_utilities.hh"
51 #include "pism/util/pism_options.hh"
52 #include "pism/verification/tests/exactTestsIJ.h"
53 
54 namespace pism {
55 namespace stressbalance {
56 
58 {
59 public:
60  SSATestCaseConst(std::shared_ptr<Context> ctx, int Mx, int My, double q,
61  SSAFactory ssafactory):
62  SSATestCase(ctx, Mx, My, 50e3, 50e3, grid::CELL_CORNER, grid::NOT_PERIODIC),
63  basal_q(q)
64  {
65  L = units::convert(ctx->unit_system(), 50.0, "km", "m"); // 50km half-width
66  H0 = 500; // m
67  dhdx = 0.005; // pure number
68  nu0 = units::convert(ctx->unit_system(), 30.0, "MPa year", "Pa s");
69  tauc0 = 1.e4; // Pa
70 
71  m_config->set_flag("basal_resistance.pseudo_plastic.enabled", true);
72  m_config->set_number("basal_resistance.pseudo_plastic.q", basal_q);
73 
74  // Use a pseudo-plastic law with a constant q determined at run time
75  m_config->set_flag("basal_resistance.pseudo_plastic.enabled", true);
76 
77  // The following is irrelevant because we will force linear rheology later.
79 
80  m_ssa = ssafactory(m_grid);
81  };
82 
83 protected:
84  virtual void initializeSSACoefficients();
85 
86  virtual void exactSolution(int i, int j,
87  double x, double y, double *u, double *v);
88 
89  double basal_q,
90  L, H0, dhdx, nu0, tauc0;
91 };
92 
94 
95  // Force linear rheology
98 
99  // The finite difference code uses the following flag to treat the non-periodic grid correctly.
100  m_config->set_flag("stress_balance.ssa.compute_surface_gradient_inward", true);
101 
102  // Set constant thickness, tauc
103  m_bc_mask.set(0);
105  m_tauc.set(tauc0);
106 
109 
110  for (auto p = m_grid->points(); p; p.next()) {
111  const int i = p.i(), j = p.j();
112 
113  double u, v;
114  const double x = m_grid->x(i), y=m_grid->y(j);
115 
116  m_geometry.bed_elevation(i, j) = -x*(dhdx);
118 
119  bool edge = ((j == 0) || (j == (int)m_grid->My() - 1) ||
120  (i == 0) || (i == (int)m_grid->Mx() - 1));
121  if (edge) {
122  m_bc_mask(i,j) = 1;
123  exactSolution(i, j, x, y, &u, &v);
124  m_bc_values(i, j).u = u;
125  m_bc_values(i, j).v = v;
126  }
127  }
128 
133 }
134 
135 
136 void SSATestCaseConst::exactSolution(int /*i*/, int /*j*/,
137  double /*x*/, double /*y*/,
138  double *u, double *v) {
139  double earth_grav = m_config->get_number("constants.standard_gravity"),
140  tauc_threshold_velocity = m_config->get_number("basal_resistance.pseudo_plastic.u_threshold",
141  "m second-1"),
142  ice_rho = m_config->get_number("constants.ice.density");
143 
144  *u = pow(ice_rho * earth_grav * H0 * dhdx / tauc0, 1./basal_q)*tauc_threshold_velocity;
145  *v = 0;
146 }
147 
148 } // end of namespace stressbalance
149 } // end of namespace pism
150 
151 int main(int argc, char *argv[]) {
152 
153  using namespace pism;
154  using namespace pism::stressbalance;
155 
156  MPI_Comm com = MPI_COMM_WORLD; // won't be used except for rank,size
157  petsc::Initializer petsc(argc, argv, help);
158 
159  com = PETSC_COMM_WORLD;
160 
161  /* This explicit scoping forces destructors to be called before PetscFinalize() */
162  try {
163  std::shared_ptr<Context> ctx = context_from_options(com, "ssa_test_const");
164  Config::Ptr config = ctx->config();
165 
166  std::string usage = "\n"
167  "usage of SSA_TEST_CONST:\n"
168  " run ssa_test_const -Mx <number> -My <number> -ssa_method <fd|fem>\n"
169  "\n";
170 
171  bool stop = show_usage_check_req_opts(*ctx->log(), "ssa_test_const", {}, usage);
172 
173  if (stop) {
174  return 0;
175  }
176 
177  // Parameters that can be overridden by command line options
178  unsigned int Mx = config->get_number("grid.Mx");
179  unsigned int My = config->get_number("grid.My");
180 
181  config->set_number("basal_resistance.pseudo_plastic.q", 1.0);
182  double basal_q = config->get_number("basal_resistance.pseudo_plastic.q");
183 
184  auto method = config->get_string("stress_balance.ssa.method");
185  auto output_file = config->get_string("output.file");
186 
187  bool write_output = config->get_string("output.size") != "none";
188 
189  // Determine the kind of solver to use.
190  SSAFactory ssafactory = NULL;
191  if (method == "fem") {
192  ssafactory = SSAFEMFactory;
193  } else if (method == "fd") {
194  ssafactory = SSAFDFactory;
195  } else {
196  /* can't happen */
197  }
198 
199  SSATestCaseConst testcase(ctx, Mx, My, basal_q, ssafactory);
200  testcase.init();
201  testcase.run();
202  testcase.report("const");
203  if (write_output) {
204  testcase.write(output_file);
205  }
206  }
207  catch (...) {
208  handle_fatal_errors(com);
209  return 1;
210  }
211 
212  return 0;
213 }
std::shared_ptr< Config > Ptr
std::shared_ptr< EnthalpyConverter > Ptr
Converts between specific enthalpy and temperature or liquid content.
array::Scalar2 ice_surface_elevation
Definition: Geometry.hh:57
array::Scalar2 ice_thickness
Definition: Geometry.hh:51
array::Scalar2 bed_elevation
Definition: Geometry.hh:47
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition: Array.hh:65
void set(double c)
Result: v[j] <- c for all j.
Definition: Array.cc:707
void update_ghosts()
Updates ghost points.
Definition: Array.cc:693
void set_min_thickness(double my_min_thickness)
Set minimum thickness to trigger use of extension.
Definition: SSA.cc:50
void set_notional_strength(double my_nuH)
Set strength = (viscosity times thickness).
Definition: SSA.cc:41
virtual void initializeSSACoefficients()
Set up the coefficient variables as appropriate for the test case.
virtual void exactSolution(int i, int j, double x, double y, double *u, double *v)
SSATestCaseConst(std::shared_ptr< Context > ctx, int Mx, int My, double q, SSAFactory ssafactory)
virtual void init()
Initialize the test case at the start of a run.
Definition: SSATestCase.cc:94
EnthalpyConverter::Ptr m_enthalpyconverter
Definition: SSATestCase.hh:102
std::shared_ptr< Grid > m_grid
Definition: SSATestCase.hh:94
virtual void report(const std::string &testname)
Report on the generated solution.
Definition: SSATestCase.cc:121
virtual void run()
Solve the SSA.
Definition: SSATestCase.cc:103
virtual void write(const std::string &filename)
Save the computation and data to a file.
Definition: SSATestCase.cc:296
SSAStrengthExtension * strength_extension
Definition: SSA.hh:115
@ NOT_PERIODIC
Definition: Grid.hh:51
@ CELL_CORNER
Definition: Grid.hh:53
SSA * SSAFDFactory(std::shared_ptr< const Grid > g)
Constructs a new SSAFD.
Definition: SSAFD.cc:49
SSA *(* SSAFactory)(std::shared_ptr< const Grid >)
Definition: SSA.hh:77
SSA * SSAFEMFactory(std::shared_ptr< const Grid > g)
Factory function for constructing a new SSAFEM.
Definition: SSAFEM.cc:123
Stress balance models and related diagnostics.
Definition: Isochrones.hh:31
double convert(System::Ptr system, double input, const std::string &spec1, const std::string &spec2)
Convert a quantity from unit1 to unit2.
Definition: Units.cc:70
std::shared_ptr< Context > context_from_options(MPI_Comm com, const std::string &prefix, bool print)
Create a default context using options.
Definition: Context.cc:170
void handle_fatal_errors(MPI_Comm com)
bool show_usage_check_req_opts(const Logger &log, const std::string &execname, const std::vector< std::string > &required_options, const std::string &usage)
In a single call a driver program can provide a usage string to the user and check if required option...
Definition: pism_options.cc:51
int main(int argc, char *argv[])
static char help[]