PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
ssa_test_linear.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: linear flow. The rheology is
20  linear (i.e. n=1 in the Glen flow law) and the basal shear stress is also
21  linear viscous flow. The geometry consists of a constant surface slope in
22  the positive x-direction, and dirichlet conditions leading to an exponential
23  solution are imposed along the entire boundary.
24 */
25 
26 
27 static char help[] =
28  "\nSSA_TEST_EXP\n"
29  " Testing program for the finite element implementation of the SSA.\n"
30  " Does a time-independent calculation. Does not run IceModel or a derived\n"
31  " class thereof.Also may be used in a PISM\n"
32  " software (regression) test.\n\n";
33 
34 #include <cmath>
35 
36 #include "pism/basalstrength/basal_resistance.hh" // IceBasalResistancePlasticLaw
37 #include "pism/stressbalance/ssa/SSAFD.hh"
38 #include "pism/stressbalance/ssa/SSAFEM.hh"
39 #include "pism/stressbalance/ssa/SSATestCase.hh"
40 #include "pism/util/Context.hh"
41 #include "pism/util/VariableMetadata.hh"
42 #include "pism/util/error_handling.hh"
43 #include "pism/util/io/File.hh"
44 #include "pism/util/petscwrappers/PetscInitializer.hh"
45 #include "pism/util/pism_utilities.hh"
46 #include "pism/util/pism_options.hh"
47 #include "pism/verification/tests/exactTestsIJ.h"
48 
49 namespace pism {
50 namespace stressbalance {
51 
53 {
54 public:
55  SSATestCaseExp(std::shared_ptr<Context> ctx, int Mx, int My, SSAFactory ssafactory)
56  : SSATestCase(ctx, Mx, My, 50e3, 50e3, grid::CELL_CORNER, grid::NOT_PERIODIC) {
57  L = units::convert(ctx->unit_system(), 50, "km", "m"); // 50km half-width
58  H0 = 500; // meters
59  dhdx = 0.005; // pure number
60  nu0 = units::convert(ctx->unit_system(), 30.0, "MPa year", "Pa s");
61  tauc0 = 1.e4; // 1kPa
62 
63  // Use a pseudo-plastic law with linear till
64  m_config->set_flag("basal_resistance.pseudo_plastic.enabled", true);
65  m_config->set_number("basal_resistance.pseudo_plastic.q", 1.0);
66 
67  // The following is irrelevant because we will force linear rheology later.
69 
70  m_ssa = ssafactory(m_grid);
71  }
72 
73 protected:
74  virtual void initializeSSACoefficients();
75 
76  virtual void exactSolution(int i, int j,
77  double x, double y, double *u, double *v);
78 
79  double L, H0, dhdx, nu0, tauc0;
80 };
81 
83 
84  // Force linear rheology
87 
88  // The finite difference code uses the following flag to treat the non-periodic grid correctly.
89  m_config->set_flag("stress_balance.ssa.compute_surface_gradient_inward", true);
90 
91  // Set constants for most coefficients.
95  m_tauc.set(tauc0);
96 
97 
98  // Set boundary conditions (Dirichlet all the way around).
99  m_bc_mask.set(0.0);
100 
102 
103  for (auto p = m_grid->points(); p; p.next()) {
104  const int i = p.i(), j = p.j();
105 
106  double myu, myv;
107  const double myx = m_grid->x(i), myy=m_grid->y(j);
108 
109  bool edge = ((j == 0) || (j == (int)m_grid->My() - 1) ||
110  (i == 0) || (i == (int)m_grid->Mx() - 1));
111  if (edge) {
112  m_bc_mask(i,j) = 1;
113  exactSolution(i,j,myx,myy,&myu,&myv);
114  m_bc_values(i,j).u = myu;
115  m_bc_values(i,j).v = myv;
116  }
117  }
118 
121 }
122 
123 
124 void SSATestCaseExp::exactSolution(int /*i*/, int /*j*/, double x, double /*y*/,
125  double *u, double *v) {
126  double tauc_threshold_velocity = m_config->get_number("basal_resistance.pseudo_plastic.u_threshold",
127  "m second-1");
128  double v0 = units::convert(m_sys, 100.0, "m year-1", "m second-1");
129  // double alpha=log(2.)/(2*L);
130  double alpha = sqrt((tauc0/tauc_threshold_velocity) / (4*nu0*H0));
131  *u = v0*exp(-alpha*(x-L));
132  *v = 0;
133 }
134 
135 } // end of namespace stressbalance
136 } // end of namespace pism
137 
138 int main(int argc, char *argv[]) {
139 
140  using namespace pism;
141  using namespace pism::stressbalance;
142 
143  MPI_Comm com = MPI_COMM_WORLD; // won't be used except for rank,size
144  petsc::Initializer petsc(argc, argv, help);
145 
146  com = PETSC_COMM_WORLD;
147 
148  /* This explicit scoping forces destructors to be called before PetscFinalize() */
149  try {
150  std::shared_ptr<Context> ctx = context_from_options(com, "ssa_test_linear");
151  Config::Ptr config = ctx->config();
152 
153  std::string usage = "\n"
154  "usage:\n"
155  " run ssa_test_linear -Mx <number> -My <number> -ssa_method <fd|fem>\n"
156  "\n";
157 
158  bool stop = show_usage_check_req_opts(*ctx->log(), "ssa_test_linear", {}, usage);
159 
160  if (stop) {
161  return 0;
162  }
163 
164  // Parameters that can be overridden by command line options
165  unsigned int Mx = config->get_number("grid.Mx");
166  unsigned int My = config->get_number("grid.My");
167 
168  auto method = config->get_string("stress_balance.ssa.method");
169  auto output_file = config->get_string("output.file");
170 
171  bool write_output = config->get_string("output.size") != "none";
172 
173  // Determine the kind of solver to use.
174  SSAFactory ssafactory = NULL;
175  if (method == "fem") {
176  ssafactory = SSAFEMFactory;
177  } else if (method == "fd") {
178  ssafactory = SSAFDFactory;
179  } else {
180  /* can't happen */
181  }
182 
183  SSATestCaseExp testcase(ctx, Mx, My, ssafactory);
184  testcase.init();
185  testcase.run();
186  testcase.report("linear");
187  if (write_output) {
188  testcase.write(output_file);
189  }
190  }
191  catch (...) {
192  handle_fatal_errors(com);
193  return 1;
194  }
195 
196  return 0;
197 }
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
SSATestCaseExp(std::shared_ptr< Context > ctx, int Mx, int My, SSAFactory ssafactory)
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)
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
const units::System::Ptr m_sys
Definition: SSATestCase.hh:96
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:112
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
static const double v0
Definition: exactTestP.cc:50
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[]