PISM, A Parallel Ice Sheet Model  stable v0.5
src/earth/greens.cc
Go to the documentation of this file.
00001 // Copyright (C) 2004-2007 Jed Brown and Ed Bueler
00002 //
00003 // This file is part of PISM.
00004 //
00005 // PISM is free software; you can redistribute it and/or modify it under the
00006 // terms of the GNU General Public License as published by the Free Software
00007 // Foundation; either version 2 of the License, or (at your option) any later
00008 // version.
00009 //
00010 // PISM is distributed in the hope that it will be useful, but WITHOUT ANY
00011 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00012 // FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
00013 // details.
00014 //
00015 // You should have received a copy of the GNU General Public License
00016 // along with PISM; if not, write to the Free Software
00017 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00018 
00019 #include <cmath>
00020 #include <gsl/gsl_sf_bessel.h>
00021 #include <gsl/gsl_integration.h>
00022 #include "matlablike.hh"
00023 
00024 #include "greens.hh"
00025 
00026 
00027 double ge_integrand(unsigned ndimMUSTBETWO, const double* xiANDeta, void* paramsIN) {
00028   // Matlab:  function z=integrand(xi,eta,dx,dy,p,q)
00029 
00030   if (ndimMUSTBETWO != 2) { perror("ge_integrand only defined for 2 variables"); }
00031   
00032   // data here is from Lingle & Clark (1985), but claims to give G^E(r) for Farrell (1972)
00033   double rmkm[42] =
00034     {0.0, 0.011,  0.111,  1.112,  2.224,  3.336,  4.448,  6.672,  8.896,  11.12, 17.79,
00035           22.24,  27.80,  33.36,  44.48,  55.60,  66.72,  88.96,  111.2,  133.4, 177.9,
00036           222.4,  278.0,  333.6,  444.8,  556.0,  667.2,  778.4,  889.6, 1001.0, 1112.0,
00037          1334.0, 1779.0, 2224.0, 2780.0, 3336.0, 4448.0, 5560.0, 6672.0, 7784.0, 8896.0,
00038         10008.0};
00039   // rm = rmkm * 1e3 (remember to convert to meters); GE /(10^12 rm) is
00040   //    vertical displacement in meters
00041   // (GE(r=0) has been computed by linear extrapolation:  GE(0) := -33.6488)
00042   double GE[42] =
00043     {-33.6488, -33.64, -33.56, -32.75, -31.86, -30.98, -30.12, -28.44, -26.87, -25.41,
00044                -21.80, -20.02, -18.36, -17.18, -15.71, -14.91, -14.41, -13.69, -13.01,
00045                -12.31, -10.95, -9.757, -8.519, -7.533, -6.131, -5.237, -4.660, -4.272,
00046                -3.999, -3.798, -3.640, -3.392, -2.999, -2.619, -2.103, -1.530, -0.292,
00047                 0.848,  1.676,  2.083,  2.057,  1.643};  
00048 
00049   struct ge_params*  params = (struct ge_params*) paramsIN;
00050   const double  dx = params->dx;
00051   const double  dy = params->dy;
00052   const int     p = params->p;
00053   const int     q = params->q;
00054   const double  xishift = (double) p * dx - xiANDeta[0];
00055   const double  etashift = (double) q * dy - xiANDeta[1];
00056   const double  r = sqrt(xishift * xishift + etashift * etashift);
00057   
00058   double  z;
00059   if (r < 0.01) {
00060     z = GE[0]/(rmkm[1] * 1.0e3 * 1.0e12);
00061   } else if (r > rmkm[41] * 1.0e3) {
00062     z = 0.0;
00063   } else {
00064     z = interp1_linear((double*) rmkm,(double*) GE,42,r / 1.0e3) / (r * 1.0e12);
00065   }
00066   return z;
00067 }
00068 
00069 
00070 double vd_integrand (double kap, void * paramsIN) {
00071 // Matlab:  function y=integrand(kap,rg,D,t,eta,R0,rk)
00072 //            beta=rg + D*kap.^4;
00073 //            expdiff=exp(-beta*t./(2*eta*kap))-ones(size(kap));
00074 //            y=expdiff.*besselj(1.0,kap*R0).*besselj(0.0,kap*rk)./beta;
00075 
00076   struct vd_params*  params = (struct vd_params*) paramsIN;
00077   const double       t = params->t;
00078   const double       R0 = params->R0;
00079   const double       rk = params->rk; 
00080   const double       rho = params->rho; 
00081   const double       grav = params->grav; 
00082   const double       D = params->D;
00083   const double       eta = params->eta;
00084   const double       beta = rho * grav + D * pow(kap,4.0);
00085   const double       expdiff = exp(-beta * t / (2.0 * eta * kap)) - 1.0;
00086   return expdiff * gsl_sf_bessel_J1(kap * R0) * gsl_sf_bessel_J0(kap * rk) / beta;
00087 }
00088 
00089 
00090 double viscDisc(double t, double H0, double R0, double r, 
00091                 double rho, double grav, double D, double eta) {
00092   // t in seconds; H0, R0, r in meters
00093 
00094   const double      ABSTOL = 1.0e-10;
00095   const double      RELTOL = 1.0e-14;
00096   const int         N_gsl_workspace = 1000;
00097   gsl_integration_workspace*
00098                     w = gsl_integration_workspace_alloc(N_gsl_workspace);
00099   double*           pts;
00100   const int         lengthpts = 142;
00101   
00102   // Matlab:  pts=[10.^(-3:-0.05:-10) 1.0e-14];
00103   pts = new double[lengthpts];
00104   for (int j=0; j < lengthpts-1; j++) {
00105     pts[j] = pow(10.0,-3.0 - 0.05 * (double) j);
00106   }
00107   pts[lengthpts-1] = 1.0e-14;
00108 
00109   // result=quadl(@integrand,pts(1),100.0*pts(1),TOL,0,rg,D,t,eta,R0,rk); % kap->infty tail
00110   gsl_function      F;
00111   struct vd_params  params = { t, R0, r, rho, grav, D, eta };
00112   double            result, error;
00113   F.function = &vd_integrand;
00114   F.params = &params;
00115   // regarding tolerance: request is for convergence of all digits and relative tolerance RELTOL
00116   gsl_integration_qag (&F, pts[1], 100.0*pts[1], ABSTOL, RELTOL, N_gsl_workspace, 
00117                        GSL_INTEG_GAUSS21, w, &result, &error);
00118 
00119   double  sum = result; 
00120   // for j=1:length(pts)-1
00121   //   result=result+quadl(@integrand,pts(j+1),pts(j),TOL,0,rg,D,t,eta,R0,rk);
00122   // end
00123   for (int j=0; j < lengthpts-1; j++) {
00124     gsl_integration_qag (&F, pts[j+1], pts[j], ABSTOL, RELTOL, N_gsl_workspace, 
00125                          GSL_INTEG_GAUSS21, w, &result, &error);
00126     sum += result;
00127   }
00128   
00129   delete [] pts;
00130   gsl_integration_workspace_free(w);
00131   // u(k)=rhoi*g*H0*R0*result;
00132   return rho * grav * H0 * R0 * sum;
00133 }
00134 
 All Classes Namespaces Files Functions Variables Typedefs Friends Defines