tgreens.hh - pism - [fork] customized build of PISM, the parallel ice sheet model (tillflux branch)
HTML git clone git://src.adamsgaard.dk/pism
DIR Log
DIR Files
DIR Refs
DIR LICENSE
---
tgreens.hh (2143B)
---
1 // Copyright (C) 2007--2009, 2014, 2015, 2017, 2019 Ed Bueler
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 #ifndef __greens_hh
20 #define __greens_hh
21
22 #include <gsl/gsl_spline.h>
23
24 namespace pism {
25 namespace bed {
26
27 //! @brief The integrand in defining the elastic Green's function from
28 //! the [@ref Farrell] earth model.
29 /*!
30 * For G^E(r), the Green's function of spherical layered elastic earth
31 * model. From data in \ref LingleClark. See also \ref BLKfastearth.
32 */
33 double ge_integrand(unsigned ndim, const double* args, void* data);
34
35 class greens_elastic {
36 public:
37 greens_elastic();
38 ~greens_elastic();
39 double operator()(double r);
40 private:
41 static const int N = 42;
42 static const double rmkm[N];
43 static const double GE[N];
44
45 gsl_interp_accel* acc;
46 gsl_spline* spline;
47
48 // disable copy constructor and the assignment operator:
49 greens_elastic(const greens_elastic &other);
50 greens_elastic& operator=(const greens_elastic&);
51 };
52
53 //! @brief Parameters used to access elastic Green's function from the
54 //! [@ref Farrell] earth model.
55 struct ge_data {
56 double dx, dy;
57 int p, q;
58 greens_elastic *G;
59 };
60
61 //! @brief Actually compute the response of the viscous half-space
62 //! model in \ref LingleClark, to a disc load.
63 double viscDisc(double t, double H0, double R0, double r,
64 double rho, double rho_ice, double grav, double D, double eta);
65
66 } // end of namespace bed
67 } // end of namespace pism
68
69 #endif /* __greens_hh */