/[MITgcm]/manual/s_phys_pkgs/text/streamice.tex
ViewVC logotype

Annotation of /manual/s_phys_pkgs/text/streamice.tex

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (hide annotations) (download) (as text)
Tue Nov 26 20:35:57 2013 UTC (11 years, 7 months ago) by dgoldberg
Branch: MAIN
File MIME type: application/x-tex
doc for streamice

1 dgoldberg 1.1 \subsection{STREAMICE Package}
2     \label{sec:pkg:streamice}
3     \begin{rawhtml}
4     <!-- CMIREDIR:package_streamice: -->
5     \end{rawhtml}
6    
7     Authors: Daniel Goldberg
8    
9     %----------------------------------------------------------------------
10     \subsubsection{Introduction}
11    
12    
13     Package ``streamice'' provides a dynamic land ice model for MITgcm.
14    
15     %----------------------------------------------------------------------
16    
17     \subsubsection{Equations Solved}
18    
19     As of now, the model tracks only 3 variables: $x$-velocity ($u$), $y$-velocity
20     ($v$), and thickness ($h$). There is also a variable that tracks coverage of
21     fractional cells, discussed...
22    
23     By default the model solves the Shallow Shelf approximation (SSA) for velocity.
24     The SSA is appropriate for floating ice (ice shelf) or ice flowing over a
25     low-friction bed (e.g. MacAyeal, 1989). The SSA consists of the $x$-momentum
26     balance:
27     \begin{equation}
28     \label{eq:xmom}
29     \partial_x(h\nu(4\dot{\varepsilon}_{xx}+2\dot{\varepsilon}_{yy})) +
30     \partial_y(2h\nu\dot{\varepsilon}_{xy}) - \tau_{bx} = \rho g h s_x
31     \end{equation}
32     the $y$-momentum balance:
33     \begin{equation}
34     \label{eq:ymom}
35     \partial_x(2h\nu\dot{\varepsilon}_{xy}) +
36     \partial_y(h\nu(4\dot{\varepsilon}_{yy}+2\dot{\varepsilon}_{xx})) - \tau_{by} =
37     \rho g h s_y.
38     \end{equation}
39     From the velocity field, thickness evolves according to the continuity equation:
40     \begin{equation}
41     \label{eq:cont}
42     h_t + \nabla\cdot(h\vec{u}) = -\dot{b},
43     \end{equation}
44     Where $\dot{b}$ is a basal mass balance (e.g. melting due to contact with the
45     ocean), positive where there is melting. Surface mass balance is not considered,
46     since {\tt pkg/streamice} is intended for localized areas (tens to hundreds of
47     kilometers) over which the integrated effect of surface mass balance is
48     generally small. Where ice is grounded, surface elevation is given by
49     \begin{equation}
50     s = R + h,
51     \end{equation}
52     where $R(x,y)$ is the bathymetry, and the basal elevation $b$ is equal to $R$.
53     If ice is floating, then the assumption of hydrostasy and constant density gives
54     \begin{equation}
55     s = (1-\frac{\rho}{\rho_w} h,
56     \end{equation}
57     where $\rho_w$ is a representative ocean density, and $b=-(\rho/\rho_w)h$. Again
58     by hydrostasy, floation is assumed wherever
59     \begin{equation}
60     h \leq -\frac{\rho_w}{\rho}R
61     \end{equation}
62     is satisfied. Floatation criteria is stored in \texttt{float\_frac\_streamice}, equal to 1 where ice is at floatation.
63    
64     The strain rates $\varepsilon_{ij}$ are generalized to the case of orthogonal
65     curvilinear coordinates, to include the "metric" terms that arise when casting
66     the equations of motion on a sphere or projection on to a sphere (see {\tt
67     pkg/SEAICE}, 6.6.2.4.8 of the MITgcm documentation). Thus
68     \begin{align}
69     \dot{\varepsilon}_{xx} = & u_x + k_1 v, \notag \\
70     \dot{\varepsilon}_{yy} = & v_y + k_1 u, \notag \\
71     \dot{\varepsilon}_{xy} = & \frac{1}{2}(u_y+v_x) + k_1 u + k_2 v. \notag
72     \end{align}
73     $\nu$ has the form arising from Glen's law
74     \begin{equation}
75     \nu =
76     \frac{1}{2}A^{-\frac{1}{n}}\left(\dot{\varepsilon}_{xx}^2+\dot{\varepsilon}_{yy}
77     ^2+\dot{\varepsilon}_{xx}\dot{\varepsilon}_{yy}+\dot{\varepsilon}_{xy}^2+\dot{
78     \varepsilon}_{min}^2\right)^{\frac{1-n}{2n}},
79     \end{equation}
80     though the form is slightly different if a hybrid formulation is used. Whether
81     $\tau_b$ is nonzero depends on whether the floatation condition is satisfied.
82     Currently this is determined simply on an instantaneous cell-by-cell basis
83     (unless subgrid interpolation is used), as is the surface elevation $s$, but
84     possibly this should be rethought if the effects of tides are to be considered.
85     $\vec{\tau}_b$ has the form
86     \begin{equation}
87     \label{eq:sliding_law}
88     \vec{\tau}_b = C (|\vec{u}|^2+u_{min}^2)^{\frac{m-1}{2}}\vec{u}.
89     \end{equation}
90     Again, the form is slightly different if a hybrid formulation is to be used. The
91     scalar term multiplying $\vec{u}$ is referred to as $\beta$ below.
92    
93     The momentum equations are solved together with appropriate boundary conditions,
94     discussed below. In the case of a calving front boundary condition (CFBC), the
95     boundary condition has the following form:
96     \begin{align}
97     \label{eq:cfbcx}
98     (h\nu(4\dot{\varepsilon}_{xx}+2\dot{\varepsilon}_{yy}))n_x +
99     (2h\nu\dot{\varepsilon}_{xy})n_y = & \frac{1}{2}g \left(\rho h^2 - \rho_w
100     b^2\right)n_x \\
101     \label{eq:cfbcy}
102     (2h\nu\dot{\varepsilon}_{xy})n_x +
103     (h\nu(4\dot{\varepsilon}_{yy}+2\dot{\varepsilon}_{xx}))n_y = & \frac{1}{2}g
104     \left(\rho h^2 - \rho_w b^2\right)n_y.
105     \end{align}
106     Here $\vec{n}$ is the normal to the boundary, and $R(x,y)$ is the bathymetry.
107    
108     \paragraph{Hybrid SIA-SSA stress balance}
109    
110     The SSA does not take vertical shear stress or strain rates (e.g.,
111     $\sigma_{xz}$, $\partial u/\partial z$) into account. Although there are other
112     terms in the stress tensor, studies have found that in all but a few cases,
113     vertical shear and longitudinal stresses (represented by the SSA) are sufficient
114     to represent glaciological flow. {\tt streamice} can allow for representation of
115     vertical shear, although the approximation is made that longitudinal stresses
116     are depth-independent. The stress balance is referred to as "hybrid" because it
117     is a joining of the SSA and the Shallow Ice Approximation (SIA), which only
118     accounts only for vertical shear. Such hybrid formulations have been shown to be
119     valid over a larger range of conditions than SSA (\textit{Schoof and Hindmarsh}
120     2010, \textit{Goldberg} 2011).
121    
122     In the hybrid formulation, $\overline{u}$ and $\overline{v}$, the depth-averaged
123     $x-$ and $y-$ velocities, replace $u$ and $v$ in \eqref{eq:xmom},
124     \eqref{eq:ymom}, and \eqref{eq:cont}, and gradients such as $u_x$ are replaced
125     by $(\overline{u})_x$. Viscosity becomes
126     \begin{equation}
127     \nu =
128     \frac{1}{2}A^{-\frac{1}{n}}\left(\dot{\varepsilon}_{xx}^2+\dot{\varepsilon}_{yy}
129     ^2+\dot{\varepsilon}_{xx}\dot{\varepsilon}_{yy}+\dot{\varepsilon}_{xy}^2+\frac{1
130     }{4}u_z^2+\frac{1}{4}v_z^2+\dot{\varepsilon}_{min}^2\right)^{\frac{1-n}{2n}}.
131     \end{equation}
132     In the formulation for $\tau_b$, $u_b$, the horizontal velocity at $u_b$ is used
133     instead. The details are given in \textit{Goldberg} (2011).
134    
135     \subsubsection{Numerical Scheme}
136    
137     %\begin{figure}
138     % \center
139     % \label{fig:setup}
140     % \includegraphics[width=4in]{s_phys_pkgs/figs/halfpipe2.eps}
141     %\end{figure}
142    
143     %\begin{figure}
144     %\center
145     %\includegraphics[width=3in]{s_phys_pkgs/figs/stencil.eps}
146    
147     %\begin{figure}
148     %\center
149     %\label{fig:hmask}
150     %\includegraphics[width=3in]{s_phys_pkgs/figs/mask_cover.eps}
151     %\end{figure}
152    
153     \begin{figure}
154     \begin{center}
155     \label{fig:stencil}
156     \resizebox{5.0in}{3.0in}{\includegraphics{s_phys_pkgs/figs/stencil.eps}}
157     \end{center}
158     \end{figure}
159    
160     \begin{figure}
161     \begin{center}
162     \label{fig:hmask}
163     \resizebox{5.0in}{3.0in}{\includegraphics{s_phys_pkgs/figs/mask_cover.eps}}
164     \end{center}
165     \end{figure}
166    
167    
168     \paragraph{Stress/momentum equations}
169    
170     The stress balance is solved for velocities using a very straightforward,
171     structured-grid finite element method. Generally a finite element method is used
172     to deal with irregularly shaped domains, or to deal with complicated
173     boundary conditions; in this case it is the latter, as explained below. (NOTE:
174     this is not meant as a finite element tutorial, and various mathematical objects
175     are defined nonrigorously. It is simply meant as an overview of the numerical
176     scheme used.)
177     In this method, the numerical velocities
178     $u$ and $v$ (or $\overline{u}$, $\overline{v}$ if a hybrid formulation) have a
179     bilinear shape within each cell. For instance, referring to Fig.
180     \ref{fig:stencil}, at a point ($x,y$) within cell ($i,j$), and assuming a
181     rectangular mesh, the $x-$velocity $u$ would be given by
182     \begin{align}
183     u(x,y) = & (1-\zeta_1)\times(1-\zeta_2)\times u_{i,j}\ + \\
184     & (\zeta_1)\times(1-\zeta_2)\times u_{i+1,j}\ + \\
185     & (1-\zeta_1)\times(\zeta_2)\times u_{i,j+1}\ + \\
186     & (\zeta_1)\times(\zeta_2)\times u_{i+1,j+1},
187     \end{align}
188     where
189     \begin{equation}
190     \zeta_1 = x-x_i, \ \ \zeta_2 = y-y_j
191     \end{equation}
192     and \textit{e.g.} $u_{i,j+1}$ is the nodal value of $u$ at the top right corner
193     of the cell. In finite element terms, the functions multiplying a nodal value of
194     $u$ are the \textit{shape
195     functions} $\phi$ corresponding to that node, e.g. $\phi_{ij}=\zeta_1\zeta_2$ in cell
196     ($i,j$). Note that $h_{ij}$ is defined at the center of each cell. In the code
197     velocities are \texttt{U\_STREAMICE} and \texttt{V\_STREAMICE}, and thickness
198     is \texttt{H\_STREAMICE}.
199    
200     If we take a vector-valued function $\Phi=(\phi,\psi)$ which is in the
201     ``solution space'' of ($u,v$) (i.e. $\phi$ and $\psi$ are bilinear in each cell
202     as
203     defined above, and are zero where velocities are imposed), take the dot product
204     with by (\ref{eq:xmom},\ref{eq:ymom}), and integrate by parts over the domain
205     $D$, we get the weak form of the momentum equation:
206     \begin{align}
207     -&\int_D \left(h\nu\phi_x (4\dot{\varepsilon}_{xx}+2\dot{\varepsilon}_{yy}) +
208     h\nu\phi_y \dot{\varepsilon}_{xy} + \beta (u\phi + v\psi) +
209     h\nu\psi_x\dot{\varepsilon}_{xy} + h\nu\psi_y
210     (4\dot{\varepsilon}_{yy}+2\dot{\varepsilon}_{xx}) \right) dA = \notag \\
211     \label{eq:weak1} & \int_D \Phi\vec{\tau}_d - \int_{\Gamma} \frac{1}{2}g
212     \left(\rho h^2 - \rho_w R^2\right) \Phi\cdot\vec{n} dS,
213     \end{align}
214     where $\Gamma$ is the part of the domain boundary where a calving front stress
215     condition is imposed, and $\vec{\tau}_d$ is the driving stress
216     $\rho g h \nabla s$. Note that the boundary integral does not involve viscosity, or any velocity-dependent terms. This is the advantage of using a finite-element
217     formulation: viscosity need not be represented at boundaries, and so no complicated one-sided derivative expressions need to be used.
218    
219     Let $(U,V)$ be the vector of all nodal values
220     $u_{i,j},v_{i,j}$ that determine $(u,v)$. If we assume that $\nu$ and $\beta$
221     are independent of $(U,V)$,
222     then for any $(\phi,\psi)$, the left hand side of
223     \eqref{eq:weak1} is a linear functional of $(U,V)$, i.e. a linear operator that
224     returns a scalar. We can choose $\phi$ and $\psi$ so that they are nonzero
225     only at a single node; and we can evaluate \eqref{eq:weak1} for each such
226     $\phi_{ij}$ and $\psi_{ij}$, giving a linearly independent set of equations for
227     $U$ and $V$. The left hand side of \eqref{eq:weak1} for each $\phi_{ij}$ and
228     $\psi_{ij}$ is
229     evaluated in the subroutine \texttt{STREAMICE\_CG\_ACTION()} in the source file
230     \texttt{streamice\_cg\_functions.f}.
231     The right hand side is evaluated in \texttt{STREAMICE\_DRIVING\_STRESS()}. The
232     latter is not a straightforward
233     calculation, since the thickness $h$ and the surface $s$ are not expressed by
234     piecewise bilinear functions, and is discussed below.
235    
236     \texttt{STREAMICE\_CG\_ACTION()} represents the \textit{action} of a matrix on
237     the vector represented by $(U,V)$ and \texttt{STREAMICE\_DRIVING\_STRESS()}
238     gives the right hand side of a linear system of equations. This information is
239     enough to solve the linear system, and this is done
240     with the method of conjugate gradients, implemented (with simple Jacobi
241     preconditioning) in \texttt{STREAMICE\_CG\_SOLVE()}.
242    
243     The full nonlinear system is solved in \texttt{STREAMICE\_VEL\_SOLVE()}. This
244     subroutine evaluates driving stress and then calls a loop
245     in which a sequence of linear systems is solved. After each linear solve, the
246     viscosity $\nu$ (the array \texttt{visc\_streamice})
247     and basal traction $\beta$ (\texttt{tau\_beta\_eff\_streamice}) is updated from
248     the new iterate of $u$ and $v$. The iteration
249     continues until a desired tolerance is reached or the maximum number of
250     iterations is reached.
251    
252     Rather than calling \texttt{STREAMICE\_CG\_ACTION()} each iteration of the conjugate gradient,
253     the sparse matrix can be constructed explicitly. This is done if \texttt{STREAMICE\_CONSTRUCT\_MATRIX}
254     is \texttt{\#define}d in \texttt{STREAMICE\_OPTIONS.h}. This speeds up the solution because \texttt{STREAMICE\_CG\_ACTION()}
255     contains a number of nested loops related to quadrature points and interactions between nodes of a cell.
256    
257     The driving stress $\tau_d=\rho g \nabla s$ is not as straightforward to deal
258     with in the weak formulation, as $h$ and $s$ are
259     considered constant in a cell. Furthermore, in the ice shelf the driving stress
260     can be written as
261     \begin{equation}
262     \label{eq:shelf_taud}
263     \vec{\tau}_d = \frac{1}{2}\rho g (1-\frac{\rho}{\rho_w})\nabla(h^2).
264     \end{equation}
265     Among other things, this (along with equations \ref{eq:cfbcx} and
266     \ref{eq:cfbcy}) implies that for an unconfined shelf (one for which the only
267     boundaries are a
268     calving front and a grounding line), the stress state along the grounding line
269     depends
270     only on the location at depth of the grounding line. The numerical
271     representation of velocities
272     must respect this. Thus the driving stress contribution in the $x-$ momentum
273     equation at node $(i,j)$ is given as follows.
274     \begin{equation}
275     \tau_{d,x}(i,j) = \begin{cases}
276     \frac{1}{4}\rho g (\gamma_{ij}
277     h_{ij}+\gamma_{i-1,j}h_{i-1,j})(\gamma_{i-1,j}s_{i-1,j}-\gamma_{ij}s_{ij}) &
278     \mbox{case I} \\
279     \frac{1}{4} \frac{rA_{i-1,j}}{\Delta x_{f,i-1,j}} g (\rho
280     h_{i-1,j}^2 - \rho_w b_{i-1,j}^2) & \mbox{case II} \\
281     -\frac{1}{4} \frac{rA_{i,j}}{\Delta x_{f,ij}} g (\rho
282     h_{ij}^2 - \rho_w b_{ij}^2) & \mbox{case III}.
283     \end{cases}
284     \end{equation}
285     where
286     \begin{equation}
287     \gamma_{ij} = \sqrt{\frac{rA_{ij}}{\Delta x_{f,ij}}},
288     \end{equation}
289     $rA_{ij}$ is the area of cell ($i,j$) and $\Delta x_{f,ij}$ is the width of the
290     cell at its $y-$midpoint. Cases I, II, and III refer to the situations when (I)
291     both cells ($i,j$) and ($i-1,j$) are in the domain, (II) only cell ($i-1,j$) is
292     in the domain, and (III) only cell ($i,j$) is in the domain. (The fact that
293     $\tau_{d,x}(i,j)$ is being calculated means that the implied boundary is a
294     calving stress boundary, as explained below.) The above expression gives the
295     contribution
296     from the cells above node ($i,j$) (in the $y-$direction). A similar expression
297     gives the contribution from cells ($i-1,j-1$) and ($i,j-1$), and
298     corresponding expressions approximate driving stress in the $y-$ momentum
299     equation. In the ice shelf, the case I expressions give approximations for
300     \eqref{eq:shelf_taud} which are discretely integrable; in grounded ice the
301     expression need not be integrable, so any errors associated with varying grid
302     metrics
303     can be subsumed into the discretization error.
304    
305     \paragraph{Orthogonal curvilinear coordinates}
306    
307     If the grid is not rectangular, but given in orthogonal curvilinear coordinates, the evaluation of \eqref{eq:weak1} is not exact, since the gradient of the basis functions are approximated, as we do not have complete knowledge of the coordinate system, only the grid metrics. Still, we consider nodal basis functions that are equal to 1 at a single node, and zero elsewhere in a cell.
308    
309     A cell ($i,j$) has separation $\Delta x_{i,j}$ at the bottom and $\Delta x_{i,j+1}$ at the top (referred to as \texttt{dxG(i,j)} and \texttt{dxG(i,j+1)} in the code). So, for instance, for a cell given by $[\xi_1,\xi_2]\times[\eta_1,\eta_2]$, the $\xi$-derivative of the basis function centered at ($\xi_1,\eta_1$) is approximated by
310     \begin{equation}
311     \left(\frac{1-q}{\Delta x_{i,j}}\right) + \left(\frac{q}{\Delta x_{i,j+1}}\right)
312     \end{equation}
313     where $q$ is the $y$-coordinate of the quadrature point being used. (The code currently uses 2-point gauss quadrature to evaluate the integrals in \eqref{eq:weak1}; with a rectangular grid and cellwise-constant $\beta$ and $\nu$, this is exact.)
314    
315     Basis function derivatives and quadrature weights are stored in \texttt{Dphi} and \texttt{grid\_jacq\_streamice}. Both are initialized in \texttt{STREAMICE\_INIT\_PHI}, called only once.
316    
317     \paragraph{Boundary conditions and masks}
318    
319     The computational domain (which may be smaller than the array/grid as defined by
320     \texttt{SIZE.h} and \texttt{GRID.h}) is determined by a number of mask arrays
321     within the \texttt{STREAMICE} package. They are
322     \begin{itemize}
323     \item $hmask$ (\texttt{STREAMICE\_hmask}): equal to 1 (ice-covered), 0 (open
324     ocean), 2 (partly-covered), or -1 (out of domain)
325     \item $umask$ (\texttt{STREAMICE\_umask}): equal to 1 (an ``active'' velocity
326     node), 3 (a Dirichlet node), or 0 (zero velocity)
327     \item $vmask$ (\texttt{STREAMICE\_vmask}): similar to umask
328     \item $ufacemaskbdry$ (\texttt{STREAMICE\_ufacemaskbdry}): equal to -1
329     (interior face), 0 (no-slip), 1 (no-stress), 2 (calving stress front), 3
330     (Dirichlet), or 4 (flux input boundary); when 4, then \texttt{u\_flux\_bdry\_SI} must be initialized, through binary or parameter file
331     \item $vfacemaskbdry$ (\texttt{STREAMICE\_vfacemaskbdry}): similar to
332     ufacemaskbdry
333     \end{itemize}
334     $hmask$ is defined at cell centers, like $h$. $umask$ and $vmask$ are defined at
335     cell nodes, like velocities. $ufacemaskbdry$ and $vfacemaskbdry$ are defined at
336     cell faces, like velocities in a $C$-grid - but unless
337     \texttt{STREAMICE\_GEOM\_FILE\_SETUP} is \texttt{\#define}d in
338     \texttt{STREAMICE\_OPTIONS.h}, the values are only relevant at the boundaries of
339     the grid.
340    
341     Essentially $hmask$ controls the domain; it is specified at initialization and does not change unless calving front advance is allowed.
342     $ufacemaskbdry$ and $vfacemaskbdry$ are defined at initialization as well. Masks are either initialized through binary files or through the text file \texttt{data.streamice} (see the streamice tutorial)
343    
344     The values of $umask$ and $vmask$ determine which nodal values of $u$ and $v$ are involved in the solve for velocities.
345     Before each call to \texttt{STREAMICE\_VEL\_SOLVE()}, $umask$ and $vmask$ are re-initialized. Values are only relevant if they border a cell where $hmask=1$.
346     Furthermore, if a cell face is a ``no-slip'' face, both $umask$ and $vmask$ are set to zero at the face's nodal endpoints.
347     If the face is a Dirichlet boundary, both nodal endpoints are set to 3. If $ufacemaskbdry=1$ on the $x-$boundary of a cell, i.e. it is a no-stress boundary, then
348     $vmask$ is set to 1 at both endpoints while $umask$ is set to zero (i.e. no normal flow). If the face is a flux boundary, velocities are set to zero (see ???).
349     For a calving stress front, $umask$ and $vmask$ are 1: the nodes represent active degrees of freedom.
350    
351     With $umask$ and $vmask$ appropriately initialized, \texttt{STREAMICE\_VEL\_SOLVE} can proceed rather generally. Contributions to \eqref{eq:weak1} are only
352     evaluated if $hmask=1$ in a given cell, and a given nodal basis function is only considered if $umask=1$ or $vmask=1$ at that node.
353    
354     \paragraph{Thickness evolution}
355    
356     \eqref{eq:cont} is solved in the subroutine \texttt{STREAMICE\_ADVECT\_THICKNESS}, similarly to the advection routines in MITgcm. Mass fluxes are evaluated, first in the $x$-direction. This is done in the generic subroutine \texttt{STREAMICE\_ADV\_FLUX\_FL\_X}. Flux velocity in the $x-$direction at face ($i,j$) are generated by averaging $u_{i,j}$ and $u_{i,j+1}$. Assuming the flux velocity is positive, if $hmask_{i-2,j},\ mask{i-1,j}$ and $hmask_{i,j}$ are equal to 1, then flux thickness, i.e. the interpolation of $h$ at face ($i,j$), is through a minmod limiter as in the \texttt{generic\_advdiff} package. If these values are not available, then a zero-order upwind flux is used. The exception is when \texttt{STREAMICE\_ufacemask(i,j)} is equal to 4; then \texttt{u\_flux\_bdry\_SI(i,j)} is used for the flux. Fluxes are then differenced to update $h$ in cells that are active ($hmask=1$); a similar procedure follows for the $y-$direction.
357    
358     \paragraph{Ice front advance}
359    
360     Flux out of the domain across Dirichlet boundaries is ignored, and flux at no-flow or no-stress boundaries is zero. However, fluxes that cross calving boundaries are stored, and if \texttt{STREAMICE\_move\_front=.true.}, then \texttt{STREAMICE\_ADV\_FRONT} is called after all cells with $hmask=1$ have their thickness updated. In this subroutine, cells with $hmask=0$ or $hmask=2$ with nonzero fluxes entering their boundaries are processed.
361    
362     The algorithm is based on \textit{Albrecht} (2011). In this scheme, if cell ($i,j$) fits the criteria, a \textbf{reference} thickness $h_{ref}(i,j)$ is found, defined as an average over the thickness of all neighboring cells with $hmask=1$ that are contributing mass to ($i,j$). The total mass input over the time step to ($i,j$) is calculated as $V_{in}(ij)$. This is added to the volume of ice already in the cell, which is nonzero only if $hmask_{ij}=2$ at the beginning of the time step, and is equal to
363     \begin{equation}
364     V_{cell}(i,j) = area_{ij} h_{ij}.
365     \end{equation}
366     (If $hmask_{ij}$ is not equal to 1, then $h_{ij}$ has not yet been updated in the current time step.) $area_{ij}$ is stored in the array \texttt{area\_shelf\_streamice}. The total volume, $V_{tot}(i,j)$, is then equal to $V_{cell}(i,j)+V_{in}(i,j)$. Next an effective thickness
367     \begin{equation}
368     h_{eff}(i,j) = \frac{V_{ij}}{rA_{ij}}
369     \end{equation}
370     is found. If $h_{eff}(ij)<h_{ref}(i,j)$ (but nonzero) then
371     \begin{itemize}
372     \item $hmask_{ij}$ is set to 2,
373     \item $h_{ij}$ is set to $h_{ref}(ij)$,
374     \item $area_{ij}$ is set to $\frac{V_{tot}(i,j)}{h_{ref}(i,j)}$.
375     \end{itemize}
376     If, on the other hand, $h_{eff}(ij)\geq h_{ref}(i,j)$,
377     \begin{itemize}
378     \item $hmask_{ij}$ is set to 1,
379     \item $h_{ij}$ is set to $h_{ref}(ij)$,
380     \item $area_{ij}$ is set to $rA_{ij}$,
381     \item Excess flux $q_{ex}(i,j) = V_{tot}(i,j) - h_{eff}(i,j)rA_{ij}$ is found.
382     \end{itemize}
383     If $q_{ex}(i,j)$ is zero, nothing more needs to happen. If it is positive, then
384     \begin{itemize}
385     \item adjacent cells that are not completely covered ($hmask=0$ or 1) are searched for.
386     \item if there are no suitable cells, $q_{ex}(i,j)$ is put back into cell ($i,j$) ($h_{ij}$ is set to $h_{ref}(i,j)+\frac{q_{ex}(i,j)}{rA_{ij}}$).
387     \item if there are suitable cells, $q_{ex}$ is divided \textbf{uniformly} among the corresponding faces of the receiving cells, and the algorithm begins again.
388     \end{itemize}
389     The last step is a recursive one; in theory the recursion could continue indefinitely, but in practice at most one recursive call is done. The decision to spread the excess flux uniformly is not the most physically motivated choice, but with appropriate time steps the excess flux value is expected to be small, and is more a matter of mass/volume conservation. This algorithm, while complicated, is necessary for realistic front advance. If newly-covered cells were given full coverage instead of partial coverage and volume were conserved, the front would quickly diffuse.
390    
391     If \texttt{calve\_to\_mask=.true.}, this sets a limit to how far the front can advance, even if advance is allowed. The front will not advance into cells where the array \texttt{calve\_mask} is not equal to 1. This mask must be set through a binary input file.
392    
393     No calving parameterisation is implemented in \texttt{STREAMICE}. However, front advancement is a precursor for such a development to be added.
394    
395     \paragraph{Surface/basal mass balance}
396    
397     After updating thickness in fully- and partially-covered cells, surface and basal mass balances are applied to nonempty cells. Currently surface mass balance is a uniform value, set through the parameter \texttt{streamice\_adot\_uniform}. Basal melt rates are stored in the array \texttt{BDOT\_streamice}, but are multiplied through by \texttt{float\_frac\_streamice} before being applied.

  ViewVC Help
Powered by ViewVC 1.1.22