|
PISM, A Parallel Ice Sheet Model stable 0.4.1779
|
00001 /* 00002 Copyright (C) 2007--2009 Ed Bueler 00003 00004 This file is part of PISM. 00005 00006 PISM is free software; you can redistribute it and/or modify it under the 00007 terms of the GNU General Public License as published by the Free Software 00008 Foundation; either version 2 of the License, or (at your option) any later 00009 version. 00010 00011 PISM is distributed in the hope that it will be useful, but WITHOUT ANY 00012 WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 00013 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more 00014 details. 00015 00016 You should have received a copy of the GNU General Public License 00017 along with PISM; if not, write to the Free Software 00018 Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA 00019 */ 00020 00021 00022 #include <stdio.h> 00023 #include <stdlib.h> 00024 #include <math.h> 00025 #include <gsl/gsl_errno.h> 00026 #include <gsl/gsl_matrix.h> 00027 #include <gsl/gsl_odeiv.h> 00028 #include "exactTestL.h" 00029 00030 #define pi 3.1415926535897931 00031 #define SperA 31556926.0 /* seconds per year; 365.2422 days */ 00032 00033 #define L 750.0e3 /* m; i.e. 750 km */ 00034 #define b0 500.0 /* m */ 00035 #define z0 1.2 00036 #define g 9.81 00037 #define rho 910.0 00038 #define n 3.0 /* Glen power */ 00039 00040 int funcL(double r, const double u[], double f[], void *params) { 00041 /* 00042 RHS for differential equation: 00043 du 5/8 / a_0 r (L^2 - r^2) \ 1/3 00044 -- = - (8/3) b'(r) u - |---------------------| 00045 dr \ 2 L^2 \tilde\Gamma / 00046 */ 00047 00048 const double Lsqr = L * L; 00049 const double a0 = 0.3 / SperA; /* m/s; i.e. 0.3 m/a */ 00050 const double A = 1.0e-16 / SperA; /* = 3.17e-24 1/(Pa^3 s); EISMINT I flow law */ 00051 const double Gamma = 2 * pow(rho * g,n) * A / (n+2); 00052 const double tilGamma = Gamma * pow(n,n) / (pow(2.0 * n + 2.0, n)); 00053 const double C = a0 / (2.0 * Lsqr * tilGamma); 00054 00055 if (params == NULL) {} /* quash warning "unused parameters" */ 00056 00057 if ((r >= 0.0) && (r <= L)) { 00058 const double freq = z0 * pi / L; 00059 const double bprime = b0 * freq * sin(freq * r); 00060 f[0] = - (8.0/3.0) * bprime * pow(u[0], 5.0/8.0) 00061 - pow(C * r * (Lsqr - r * r), 1.0/3.0); 00062 } else { 00063 f[0] = 0.0; /* no changes outside of defined interval */ 00064 } 00065 return GSL_SUCCESS; 00066 } 00067 00068 00069 /* combination EPS_ABS = 1e-12, EPS_REL=0.0, method = 1 = RK Cash-Karp 00070 is believed to be predictable and accurate */ 00071 int getU(double *r, int N, double *u, 00072 const double EPS_ABS, const double EPS_REL, const int ode_method) { 00073 /* solves ODE for u(r)=H(r)^{8/3}, 0 <= r <= L, for test L 00074 r and u must be allocated vectors of length N; r[] must be decreasing */ 00075 int i, k, count; 00076 const gsl_odeiv_step_type* T; 00077 int status = TESTL_NOT_DONE; 00078 double rr, step; 00079 00080 gsl_odeiv_step* s; 00081 gsl_odeiv_control* c; 00082 gsl_odeiv_evolve* e; 00083 gsl_odeiv_system sys = {funcL, NULL, 1, NULL}; /* Jac-free method and no params */ 00084 00085 /* check first: we have a list, and r is decreasing */ 00086 if (N < 1) return TESTL_NO_LIST; 00087 for (i = 1; i<N; i++) { if (r[i] > r[i-1]) return TESTL_NOT_DECREASING; } 00088 00089 /* setup for GSL ODE solver; following step choices don't need Jacobian, 00090 but should we chose one that does? */ 00091 switch (ode_method) { 00092 case 1: 00093 T = gsl_odeiv_step_rkck; 00094 break; 00095 case 2: 00096 T = gsl_odeiv_step_rk2; 00097 break; 00098 case 3: 00099 T = gsl_odeiv_step_rk4; 00100 break; 00101 case 4: 00102 T = gsl_odeiv_step_rk8pd; 00103 break; 00104 default: 00105 printf("INVALID ode_method in getU(): must be 1,2,3,4\n"); 00106 return TESTL_INVALID_METHOD; 00107 } 00108 00109 s = gsl_odeiv_step_alloc(T, (size_t)1); /* one scalar ode */ 00110 c = gsl_odeiv_control_y_new(EPS_ABS,EPS_REL); 00111 e = gsl_odeiv_evolve_alloc((size_t)1); /* one scalar ode */ 00112 00113 00114 /* outside of ice cap, u = 0 */ 00115 k = 0; 00116 while (r[k] >= L) { 00117 u[k] = 0.0; 00118 k++; 00119 if (k == N) return GSL_SUCCESS; 00120 } 00121 00122 /* initial conditions: (r,u) = (L,0); r decreases from L */ 00123 rr = L; 00124 for (count = k; count < N; count++) { 00125 /* generally use value at end of last interval as initial guess */ 00126 u[count] = (count == 0) ? 0.0 : u[count-1]; 00127 while (rr > r[count]) { 00128 step = r[count] - rr; 00129 status = gsl_odeiv_evolve_apply(e, c, s, &sys, &rr, r[count], &step, &u[count]); 00130 if (status != GSL_SUCCESS) break; 00131 } 00132 } 00133 00134 gsl_odeiv_evolve_free(e); 00135 gsl_odeiv_control_free(c); 00136 gsl_odeiv_step_free(s); 00137 return status; 00138 } 00139 00140 00141 int exactL(double r, double *H, double *b, double *a, 00142 const double EPS_ABS, const double EPS_REL, const int ode_method) { 00143 00144 double u[1] = { 0.0 }; 00145 const double Lsqr = L * L; 00146 const double a0 = 0.3 / SperA; /* m/s; i.e. 0.3 m/a */ 00147 00148 getU(&r,1,u,EPS_ABS,EPS_REL,ode_method); 00149 *H = pow(u[0],3.0/8.0); 00150 *b = - b0 * cos(z0 * pi * r / L); 00151 *a = a0 * (1.0 - (2.0 * r * r / Lsqr)); 00152 return 0; 00153 } 00154 00155 00156 int exactL_list(double *r, int N, double *H, double *b, double *a) { 00157 /* N values r[0] > r[1] > ... > r[N-1] (decreasing) 00158 assumes r, H, b, a are allocated length N arrays */ 00159 00160 const double Lsqr = L * L; 00161 const double a0 = 0.3 / SperA; /* m/s; i.e. 0.3 m/a */ 00162 double *u; 00163 int stat, i; 00164 00165 u = (double *) malloc((size_t)N * sizeof(double)); /* temporary arrays */ 00166 00167 /* combination EPS_ABS = 1e-12, EPS_REL=0.0, method = 1 = RK Cash-Karp 00168 believed to be predictable and accurate */ 00169 stat = getU(r,N,u,1.0e-12,0.0,1); 00170 if (stat != GSL_SUCCESS) { 00171 return stat; 00172 } 00173 00174 for (i = 0; i < N; i++) { 00175 H[i] = pow(u[i],3.0/8.0); 00176 b[i] = - b0 * cos(z0 * pi * r[i] / L); 00177 a[i] = a0 * (1.0 - (2.0 * r[i] * r[i] / Lsqr)); 00178 } 00179 00180 free(u); 00181 return 0; 00182 } 00183
1.7.3