Beispiel #1
0
    def _advect(self, dt, E, u, w, h, s, E_inflow, E_surface):
        Q = E.function_space()
        φ, ψ = firedrake.TrialFunction(Q), firedrake.TestFunction(Q)

        U = firedrake.as_vector((u[0], u[1], w))
        flux_cells = -φ * inner(U, grad(ψ)) * h * dx

        mesh = Q.mesh()
        ν = facet_normal_2(mesh)
        outflow = firedrake.max_value(inner(u, ν), 0)
        inflow = firedrake.min_value(inner(u, ν), 0)

        flux_outflow = φ * ψ * outflow * h * ds_v + \
                       φ * ψ * firedrake.max_value(-w, 0) * h * ds_b + \
                       φ * ψ * firedrake.max_value(+w, 0) * h * ds_t
        F = φ * ψ * h * dx + dt * (flux_cells + flux_outflow)

        flux_inflow = -E_inflow * ψ * inflow * h * ds_v \
                      -E_surface * ψ * firedrake.min_value(-w, 0) * h * ds_b \
                      -E_surface * ψ * firedrake.min_value(+w, 0) * h * ds_t
        A = E * ψ * h * dx + dt * flux_inflow

        solver_parameters = {'ksp_type': 'preonly', 'pc_type': 'lu'}
        degree_E = E.ufl_element().degree()
        degree_u = u.ufl_element().degree()
        degree = (3 * degree_E[0] + degree_u[0], 2 * degree_E[1] + degree_u[1])
        form_compiler_parameters = {'quadrature_degree': degree}
        firedrake.solve(F == A,
                        E,
                        solver_parameters=solver_parameters,
                        form_compiler_parameters=form_compiler_parameters)
Beispiel #2
0
    def advective_flux(self, **kwargs):
        keys = (
            "energy",
            "velocity",
            "vertical_velocity",
            "thickness",
            "energy_inflow",
            "energy_surface",
        )
        E, u, w, h, E_inflow, E_surface = itemgetter(*keys)(kwargs)

        Q = E.function_space()
        ψ = firedrake.TestFunction(Q)

        # NOTE: Be careful here going to xz! You might have to separate this into
        # the sum of a horizontal and vertical flux if we're shadowing Firedrake's
        # grad operator with out own specialized one.
        U = firedrake.as_vector((u[0], u[1], w))
        flux_cells = -E * inner(U, grad(ψ)) * h * dx

        ν = FacetNormal(Q.mesh())
        outflow = firedrake.max_value(inner(u, ν), 0)
        inflow = firedrake.min_value(inner(u, ν), 0)

        flux_outflow = (E * outflow * ψ * h * ds_v +
                        E * firedrake.max_value(-w, 0) * ψ * h * ds_b +
                        E * firedrake.max_value(+w, 0) * ψ * h * ds_t)

        flux_inflow = (E_inflow * inflow * ψ * h * ds_v +
                       E_surface * firedrake.min_value(-w, 0) * ψ * h * ds_b +
                       E_surface * firedrake.min_value(+w, 0) * ψ * h * ds_t)

        return flux_cells + flux_outflow + flux_inflow
Beispiel #3
0
def divMelt(h, floating, meltParams, u, Q):
    """ Melt function that is a scaled version of the flux divergence
    h : firedrake function
        ice thickness
    u : firedrake vector function
        surface elevation
    floating : firedrake function
        floating mask
    V : firedrake vector space
        vector space for velocity
    meltParams : dict
        parameters for melt function
    Returns
    -------
    firedrake function
        melt rates
    """

    flux = u * h
    fluxDiv = icepack.interpolate(firedrake.div(flux), Q)
    fluxDivS = firedrakeSmooth(fluxDiv, alpha=8000)
    fluxDivS = firedrake.min_value(
        fluxDivS * floating * meltParams['meltMask'], 0)
    intFluxDiv = firedrake.assemble(fluxDivS * firedrake.dx)
    scale = -1.0 * float(meltParams['intMelt']) / float(intFluxDiv)
    scale = firedrake.Constant(scale)
    melt = icepack.interpolate(
        firedrake.min_value(fluxDivS * scale, meltParams['maxMelt']), Q)

    return melt
Beispiel #4
0
    def advective_flux(self, **kwargs):
        keys = ('energy', 'velocity', 'vertical_velocity', 'thickness',
                'energy_inflow', 'energy_surface')
        keys_alt = ('E', 'u', 'w', 'h', 'E_inflow', 'E_surface')
        E, u, w, h, E_inflow, E_surface = get_kwargs_alt(
            kwargs, keys, keys_alt)

        Q = E.function_space()
        ψ = firedrake.TestFunction(Q)

        U = firedrake.as_vector((u[0], u[1], w))
        flux_cells = -E * inner(U, grad(ψ)) * h * dx

        mesh = Q.mesh()
        ν = facet_normal_2(mesh)
        outflow = firedrake.max_value(inner(u, ν), 0)
        inflow = firedrake.min_value(inner(u, ν), 0)

        flux_outflow = (E * outflow * ψ * h * ds_v +
                        E * firedrake.max_value(-w, 0) * ψ * h * ds_b +
                        E * firedrake.max_value(+w, 0) * ψ * h * ds_t)

        flux_inflow = (E_inflow * inflow * ψ * h * ds_v +
                       E_surface * firedrake.min_value(-w, 0) * ψ * h * ds_b +
                       E_surface * firedrake.min_value(+w, 0) * ψ * h * ds_t)

        return flux_cells + flux_outflow + flux_inflow
Beispiel #5
0
def getModelVelocity(baseName, Q, V, minSigma=5, maxSigma=100):
    """Read in a tiff velocity data set and return
    firedrake interpolate functions.

    Parameters
    ----------
    baseName : str
        baseName should be of the form pattern.*.abc or pattern
        The wildcard (*) will be filled with the suffixes (vx, vy.)
        e.g.,pattern.vx.abc.tif, pattern.vy.abc.tif.
    Q : firedrake function space
        function space
    V : firedrake vector space
        vector space
    Returns
    -------
    uObs firedrake interp function on V
        velocity (m/yr)
    speed firedrake interp function on Q
        speed in (m)
    sigmaX firedrake interp function on Q
        vx error (m)
    sigmaY firedrake interp function on Q
        vy error (m)
    """
    # suffixes for products used
    suffixes = ['vx', 'vy', 'ex', 'ey']
    rasters = {}
    # prep baseName - baseName.*.xyz.tif or baseName.*
    if '*' not in baseName:
        baseName += '.*'
    if '.tif' not in baseName:
        baseName += '.tif'
    # read data
    for suffix in suffixes:
        myBand = baseName.replace('*', suffix)
        if not os.path.exists(myBand):
            u.myerror(f'Velocity/error file - {myBand} - does not exist')
        rasters[suffix] = rasterio.open(myBand, 'r')
    # Firedrake interpolators
    uObs = icepack.interpolate((rasters['vx'], rasters['vy']), V)
    # force error to be at least 1 to avoid 0 or negatives.
    sigmaX = icepack.interpolate(rasters['ex'], Q)
    sigmaX = icepack.interpolate(firedrake.max_value(sigmaX, minSigma), Q)
    sigmaX = icepack.interpolate(firedrake.min_value(sigmaX, maxSigma), Q)
    sigmaY = icepack.interpolate(rasters['ey'], Q)
    sigmaY = icepack.interpolate(firedrake.max_value(sigmaY, minSigma), Q)
    sigmaY = icepack.interpolate(firedrake.min_value(sigmaY, maxSigma), Q)
    speed = icepack.interpolate(firedrake.sqrt(inner(uObs, uObs)), Q)
    # return results
    return uObs, speed, sigmaX, sigmaY
Beispiel #6
0
    def __init__(self, state, accretion=True, accumulation=True):
        super().__init__(state)

        # obtain our fields
        self.water_c = state.fields('water_c')
        self.rain = state.fields('rain')

        # declare function space
        Vt = self.water_c.function_space()

        # define some parameters as attributes
        dt = state.timestepping.dt
        k_1 = Constant(0.001)  # accretion rate in 1/s
        k_2 = Constant(2.2)  # accumulation rate in 1/s
        a = Constant(0.001)  # min cloud conc in kg/kg
        b = Constant(0.875)  # power for rain in accumulation

        # make default rates to be zero
        accr_rate = Constant(0.0)
        accu_rate = Constant(0.0)

        if accretion:
            accr_rate = k_1 * (self.water_c - a)
        if accumulation:
            accu_rate = k_2 * self.water_c * self.rain**b

        # make coalescence rate function, that needs to be the same for all updates in one time step
        coalesce_rate = Function(Vt)

        # adjust coalesce rate using min_value so negative cloud concentration doesn't occur
        self.lim_coalesce_rate = Interpolator(
            conditional(
                self.rain < 0.0,  # if rain is negative do only accretion
                conditional(accr_rate < 0.0, 0.0,
                            min_value(accr_rate, self.water_c / dt)),
                # don't turn rain back into cloud
                conditional(
                    accr_rate + accu_rate < 0.0,
                    0.0,
                    # if accretion rate is negative do only accumulation
                    conditional(
                        accr_rate < 0.0, min_value(accu_rate,
                                                   self.water_c / dt),
                        min_value(accr_rate + accu_rate, self.water_c / dt)))),
            coalesce_rate)

        # tell the prognostic fields what to update to
        self.water_c_new = Interpolator(self.water_c - dt * coalesce_rate, Vt)
        self.rain_new = Interpolator(self.rain + dt * coalesce_rate, Vt)
Beispiel #7
0
def terminus(**kwargs):
    r"""Return the terminal stress part of the ice stream action functional

    The power exerted due to stress at the ice calving terminus :math:`\Gamma`
    is

    .. math::
       E(u) = \frac{1}{2}\int_\Gamma\left(\rho_Igh^2 - \rho_Wgd^2\right)
       u\cdot \nu\, ds

    where :math:`d` is the water depth at the terminus. We assume that sea
    level is at :math:`z = 0` for purposes of calculating the water depth.

    Parameters
    ----------
    velocity : firedrake.Function
    thickness : firedrake.Function
    surface : firedrake.Function
    """
    keys = ('velocity', 'thickness', 'surface')
    keys_alt = ('u', 'h', 's')
    u, h, s = get_kwargs_alt(kwargs, keys, keys_alt)

    d = firedrake.min_value(s - h, 0)
    τ_I = ρ_I * g * h**2 / 2
    τ_W = ρ_W * g * d**2 / 2

    ν = firedrake.FacetNormal(u.ufl_domain())
    return (τ_I - τ_W) * inner(u, ν)
Beispiel #8
0
    def solve(self, dt, **kwargs):
        if not hasattr(self, '_solvers'):
            self._setup(**kwargs)
        else:
            for name, field in kwargs.items():
                if isinstance(field, firedrake.Function):
                    self.fields[name].assign(field)

        δt = self._timestep
        δt.assign(dt)
        D = self.fields.get('damage', self.fields.get('D'))

        solver1, solver2, solver3 = self._solvers
        D1, D2 = self._stages
        dD = self._damage_change

        solver1.solve()
        D1.assign(D + dD)
        solver2.solve()
        D2.assign(3 / 4 * D + 1 / 4 * (D1 + dD))
        solver3.solve()
        D.assign(1 / 3 * D + 2 / 3 * (D2 + dD))

        S = self.model.sources(**self.fields)
        D.project(min_value(max_value(D + δt * S, 0), 1))
        return D.copy(deepcopy=True)
Beispiel #9
0
def reduceNearGLBeta(s, sOrig, zF, grounded, Q, thresh, limit=False):
    """Compute beta reduction in area near grounding line where the height
    above floation is less than thresh.

    Parameters
    ----------
    s : firedrake function
        Current surface.
    sOrig : firedrake function
        Original surface.
    zF : firedrake function
        Flotation height.
    grounded : firedrake function
        grounde mask
    Q : firedrake function space
        scaler function space for model
    thresh : float
        Threshold to determine where to reduce beta (zAbove < thresh)
    """
    # compute original height above flotation for grounded region
    sAboveOrig = (sOrig - zF) * grounded
    # avoid negative/zero values
    sAboveOrig = firedrake.max_value(sAboveOrig, 0.0)
    # Current height above flotation with negative values zeroed out.
    sAbove = firedrake.max_value((s - zF) * grounded, 0)
    # mask so only areas less than thresh but grounded
    sMask = (sAbove < thresh) * grounded
    # print(f'{sAbove.dat.data_ro.min()}, {sAbove.dat.data_ro.max()} {thresh}')
    # Inverse mask
    sMaskInv = sMask < 1
    # scale = fraction of original height above flotation
    # Use 5 to avoid potentially large ratio at small values.
    scaleBeta = sAbove / \
        firedrake.max_value(firedrake.min_value(thresh, sAboveOrig), 3)
    if limit:
        scaleBeta = firedrake.min_value(3, scaleBeta)
    scaleBeta = scaleBeta * sMask + sMaskInv
    # scaleBeta = icepack.interpolate(firedrake.min_value(scaleBeta,1.),Q)
    # sqrt so tau = scale * beta^2
    # scaleBeta = icepack.interpolate(firedrake.sqrt(scaleBeta) * grounded, Q)
    # Removed grounded above because grounded is always applied in friction
    scaleBeta = icepack.interpolate(firedrake.sqrt(scaleBeta), Q)
    # print(f'{scaleBeta.dat.data_ro.min()}, {scaleBeta.dat.data_ro.max()}')
    return scaleBeta
Beispiel #10
0
def readSMB(SMBfile, Q):
    ''' Read SMB file an limit values to +/- 6 to avoid no data values

    Returns water equivalent values.
    '''
    if not os.path.exists:
        myerror(f'readSMB: SMB file  ({SMBfile}) does not exist')
    SMB = mf.getModelVarFromTiff(SMBfile, Q)
    # avoid any unreasonably large value
    SMB = icepack.interpolate(
        firedrake.max_value(firedrake.min_value(SMB, 6), -6), Q)
    return SMB
Beispiel #11
0
    def solve(self, dt, h0, a, u, h_inflow=None):
        r"""Propagate the thickness forward by one timestep

        This function uses the implicit Euler timestepping scheme to avoid
        the stability issues associated to using continuous finite elements
        for advection-type equations. The implicit Euler scheme is stable
        for any timestep; you do not need to ensure that the CFL condition
        is satisfied in order to get an answer. Nonetheless, keeping the
        timestep within the CFL bound is a good idea for accuracy.

        Parameters
        ----------
        dt : float
            Timestep
        h0 : firedrake.Function
            Initial ice thickness
        a : firedrake.Function
            Sum of accumulation and melt rates
        u : firedrake.Function
            Ice velocity
        h_inflow : firedrake.Function
            Thickness of the upstream ice that advects into the domain

        Returns
        -------
        h : firedrake.Function
            Ice thickness at `t + dt`
        """
        grad, ds = self.grad, self.ds

        h_inflow = h_inflow if h_inflow is not None else h0

        Q = h0.function_space()
        h, φ = firedrake.TrialFunction(Q), firedrake.TestFunction(Q)

        n = self.facet_normal(Q.mesh())
        outflow = firedrake.max_value(inner(u, n), 0)
        inflow = firedrake.min_value(inner(u, n), 0)

        flux_cells = -h * inner(u, grad(φ)) * dx
        flux_out = h * φ * outflow * ds
        F = h * φ * dx + dt * (flux_cells + flux_out)

        accumulation = a * φ * dx
        flux_in = -h_inflow * φ * inflow * ds
        A = h0 * φ * dx + dt * (accumulation + flux_in)

        h = h0.copy(deepcopy=True)
        solver_parameters = {'ksp_type': 'preonly', 'pc_type': 'lu'}
        firedrake.solve(F == A, h, solver_parameters=solver_parameters)

        return h
Beispiel #12
0
def reduceNearGLBeta(s, sOrig, zF, grounded, Q, thresh, limit=False):
    """Compute beta reduction in area near grounding line where the height
    above floation is less than thresh.

    Parameters
    ----------
    s : firedrake function
        Current surface.
    sOrig : firedrake function
        Original surface.
    zF : firedrake function
        Flotation height.
    grounded : firedrake function
        grounde mask
    Q : firedrake function space
        scaler function space for model
    thresh : float
        Threshold to determine where to reduce beta (zAbove < thresh)
    """
    # compute original height above flotation for grounded region
    sAboveOrig = (sOrig - zF) * grounded
    # avoid negative/zero values
    sAboveOrig = firedrake.max_value(sAboveOrig, 0.001)
    # Current height above flotation with negative values zeroed out.
    sAbove = firedrake.max_value((s - zF) * grounded, 0)
    # mask so only areas less than thresh but grounded
    sMask = (sAbove <= thresh) * grounded
    # Inverse mask
    sMaskInv = sMask < 1
    # scale = fraction of of original height above flotation
    scaleBeta = sAbove / \
        firedrake.max_value(firedrake.min_value(thresh, sAboveOrig), 5)
    if limit:
        scaleBeta = firedrake.min_value(3, scaleBeta)
    scaleBeta = scaleBeta * sMask + sMaskInv
    # scaleBeta = icepack.interpolate(firedrake.min_value(scaleBeta,1.),Q)
    # sqrt so tau = scale * beta^2
    scaleBeta = icepack.interpolate(firedrake.sqrt(scaleBeta) * grounded, Q)
    return scaleBeta
Beispiel #13
0
    def solve(self, dt, h0, a, u, h_inflow=None):
        r"""Propagate the thickness forward by one timestep

        This function uses an implicit second-order Taylor-Galerkin (also
        known as Lax-Wendroff) scheme to solve the conservative advection
        equation for ice thickness.

        Parameters
        ----------
        dt : float
            Timestep
        h0 : firedrake.Function
            Initial ice thickness
        a : firedrake.Function
            Sum of accumulation and melt rates
        u : firedrake.Function
            Ice velocity
        h_inflow : firedrake.Function
            Thickness of the upstream ice that advects into the domain

        Returns
        -------
        h : firedrake.Function
            Ice thickness at `t + dt`
        """
        grad, div, ds = self.grad, self.div, self.ds

        h_inflow = h_inflow if h_inflow is not None else h0

        Q = h0.function_space()
        h, φ = firedrake.TrialFunction(Q), firedrake.TestFunction(Q)

        n = self.facet_normal(Q.mesh())
        outflow = firedrake.max_value(inner(u, n), 0)
        inflow = firedrake.min_value(inner(u, n), 0)

        flux_cells = -h * inner(u, grad(φ)) * dx
        flux_cells_lax = 0.5 * dt * div(h * u) * inner(u, grad(φ)) * dx
        flux_out = (h - 0.5 * dt * div(h * u)) * φ * outflow * ds
        F = h * φ * dx + dt * (flux_cells + flux_cells_lax + flux_out)

        accumulation = a * φ * dx
        flux_in = -(h_inflow - 0.5 * dt * div(h0 * u)) * φ * inflow * ds
        A = h0 * φ * dx + dt * (accumulation + flux_in)

        h = h0.copy(deepcopy=True)
        solver_parameters = {'ksp_type': 'preonly', 'pc_type': 'lu'}
        firedrake.solve(F == A, h, solver_parameters=solver_parameters)

        return h
Beispiel #14
0
def test_converting_fields():
    δT = 5.0
    T_surface = Tm - δT
    T_expr = firedrake.min_value(Tm, T_surface + 2 * δT * (1 - ζ))
    f_expr = firedrake.max_value(0, 0.0033 * (1 - 2 * ζ))

    model = icepack.models.HeatTransport3D()
    E = firedrake.project(model.energy_density(T_expr, f_expr), Q)
    f = firedrake.project(model.meltwater_fraction(E), Q)
    T = firedrake.project(model.temperature(E), Q)

    avg_meltwater = firedrake.assemble(f * ds_b) / (Lx * Ly)
    assert avg_meltwater > 0

    avg_temp = firedrake.assemble(T * h * dx) / firedrake.assemble(h * dx)
    assert (avg_temp > T_surface) and (avg_temp < Tm)
Beispiel #15
0
    def __call__(self, dt, **kwargs):
        keys = ('thickness', 'velocity', 'accumulation')
        keys_alt = ('h', 'u', 'a')
        h, u, a = utilities.get_kwargs_alt(kwargs, keys, keys_alt)
        h_inflow = kwargs.get('thickness_inflow', kwargs.get('h_inflow', h))

        Q = h.function_space()
        q = firedrake.TestFunction(Q)

        grad, ds, n = self.grad, self.ds, self.facet_normal(Q.mesh())
        u_n = inner(u, n)
        flux_cells = -inner(h * u, grad(q)) * dx
        flux_out = h * firedrake.max_value(u_n, 0) * q * ds
        flux_in = h_inflow * firedrake.min_value(u_n, 0) * q * ds
        accumulation = a * q * dx
        return accumulation - (flux_in + flux_out + flux_cells)
Beispiel #16
0
    def equation(q):
        Q = q.function_space()
        φ = firedrake.TestFunction(Q)

        F = q * u

        mesh = Q.mesh()
        n = firedrake.FacetNormal(mesh)
        c = abs(inner(u, n))

        sources = s * φ * dx
        fluxes = (forms.cell_flux(F, φ) + forms.central_facet_flux(F, φ) +
                  forms.lax_friedrichs_facet_flux(q, c, φ) +
                  q * max_value(0, inner(u, n)) * φ * ds +
                  q_in * min_value(0, inner(u, n)) * φ * ds)

        return sources - fluxes
Beispiel #17
0
    def setup(self, **kwargs):
        r"""Create the internal data structures that help reuse information
        from past prognostic solves"""
        for name, field in kwargs.items():
            if name in self._fields.keys():
                self._fields[name].assign(field)
            else:
                if isinstance(field, firedrake.Constant):
                    self._fields[name] = firedrake.Constant(field)
                elif isinstance(field, firedrake.Function):
                    self._fields[name] = field.copy(deepcopy=True)
                else:
                    raise TypeError(
                        "Input %s field has type %s, must be Constant or Function!"
                        % (name, type(field))
                    )

        dt = firedrake.Constant(1.0)
        h = self._fields["thickness"]
        u = self._fields["velocity"]
        h_0 = h.copy(deepcopy=True)

        Q = h.function_space()
        mesh = Q.mesh()
        n = FacetNormal(mesh)
        outflow = firedrake.max_value(0, inner(u, n))
        inflow = firedrake.min_value(0, inner(u, n))

        # Additional streamlining terms that give 2nd-order accuracy
        q = firedrake.TestFunction(Q)
        ds = firedrake.ds if mesh.layers is None else firedrake.ds_v
        flux_cells = -div(h * u) * inner(u, grad(q)) * dx
        flux_out = div(h * u) * q * outflow * ds
        flux_in = div(h_0 * u) * q * inflow * ds
        d2h_dt2 = flux_cells + flux_out + flux_in

        dh_dt = self._continuity(dt, **self._fields)
        F = (h - h_0) * q * dx - dt * (dh_dt + 0.5 * dt * d2h_dt2)

        problem = firedrake.NonlinearVariationalProblem(F, h)
        self._solver = firedrake.NonlinearVariationalSolver(
            problem, solver_parameters=self._solver_parameters
        )

        self._thickness_old = h_0
        self._timestep = dt
Beispiel #18
0
    def sources(self, **kwargs):
        keys = ("damage", "velocity", "strain_rate", "membrane_stress")
        D, u, ε, M = itemgetter(*keys)(kwargs)

        # Increase/decrease damage depending on stress and strain rates
        ε_1 = eigenvalues(ε)[0]
        σ_e = sqrt(inner(M, M) - det(M))

        ε_h = firedrake.Constant(self.healing_strain_rate)
        σ_d = firedrake.Constant(self.damage_stress)
        γ_h = firedrake.Constant(self.healing_rate)
        γ_d = firedrake.Constant(self.damage_rate)

        healing = γ_h * min_value(ε_1 - ε_h, 0)
        fracture = γ_d * conditional(σ_e - σ_d > 0, ε_1, 0.0) * (1 - D)

        return healing + fracture
Beispiel #19
0
    def __call__(self, dt, **kwargs):
        keys = ("thickness", "velocity", "accumulation")
        h, u, a = itemgetter(*keys)(kwargs)
        h_inflow = kwargs.get("thickness_inflow", h)

        Q = h.function_space()
        q = firedrake.TestFunction(Q)

        mesh = Q.mesh()
        n = FacetNormal(mesh)
        ds = firedrake.ds if mesh.layers is None else firedrake.ds_v

        u_n = inner(u, n)
        flux_cells = -inner(h * u, grad(q)) * dx
        flux_out = h * firedrake.max_value(u_n, 0) * q * ds
        flux_in = h_inflow * firedrake.min_value(u_n, 0) * q * ds
        accumulation = a * q * dx
        return accumulation - (flux_in + flux_out + flux_cells)
Beispiel #20
0
    def flux(self, **kwargs):
        keys = ("damage", "velocity", "damage_inflow")
        D, u, D_inflow = itemgetter(*keys)(kwargs)

        Q = D.function_space()
        φ = firedrake.TestFunction(Q)

        mesh = Q.mesh()
        n = firedrake.FacetNormal(mesh)

        u_n = max_value(0, inner(u, n))
        f = D * u_n
        flux_faces = (f("+") - f("-")) * (φ("+") - φ("-")) * dS
        flux_cells = -D * div(u * φ) * dx
        flux_out = D * max_value(0, inner(u, n)) * φ * ds
        flux_in = D_inflow * min_value(0, inner(u, n)) * φ * ds

        return flux_faces + flux_cells + flux_out + flux_in
Beispiel #21
0
    def sources(self, **kwargs):
        keys = ('damage', 'velocity', 'strain_rate', 'membrane_stress')
        keys_alt = ('D', 'u', 'ε', 'M')
        D, u, ε, M = get_kwargs_alt(kwargs, keys, keys_alt)

        # Increase/decrease damage depending on stress and strain rates
        ε_1 = eigenvalues(ε)[0]
        σ_e = sqrt(inner(M, M) - det(M))

        ε_h = firedrake.Constant(self.healing_strain_rate)
        σ_d = firedrake.Constant(self.damage_stress)
        γ_h = firedrake.Constant(self.healing_rate)
        γ_d = firedrake.Constant(self.damage_rate)

        healing = γ_h * min_value(ε_1 - ε_h, 0)
        fracture = γ_d * conditional(σ_e - σ_d > 0, ε_1, 0.) * (1 - D)

        return healing + fracture
Beispiel #22
0
    def setup(self, **kwargs):
        r"""Create the internal data structures that help reuse information
        from past prognostic solves"""
        for name, field in kwargs.items():
            if name in self._fields.keys():
                self._fields[name].assign(field)
            else:
                if isinstance(field, firedrake.Constant):
                    self._fields[name] = firedrake.Constant(field)
                elif isinstance(field, firedrake.Function):
                    self._fields[name] = field.copy(deepcopy=True)
                else:
                    raise TypeError(
                        'Input fields must be Constant or Function!')

        dt = firedrake.Constant(1.)
        h = self._fields.get('thickness', self._fields.get('h'))
        u = self._fields.get('velocity', self._fields.get('u'))
        h_0 = h.copy(deepcopy=True)

        Q = h.function_space()
        model = self._continuity
        n = model.facet_normal(Q.mesh())
        outflow = firedrake.max_value(0, inner(u, n))
        inflow = firedrake.min_value(0, inner(u, n))

        # Additional streamlining terms that give 2nd-order accuracy
        q = firedrake.TestFunction(Q)
        div, grad, ds = model.div, model.grad, model.ds
        flux_cells = -div(h * u) * inner(u, grad(q)) * dx
        flux_out = div(h * u) * q * outflow * ds
        flux_in = div(h_0 * u) * q * inflow * ds
        d2h_dt2 = flux_cells + flux_out + flux_in

        dh_dt = model(dt, **self._fields)
        F = (h - h_0) * q * dx - dt * (dh_dt + 0.5 * dt * d2h_dt2)

        problem = firedrake.NonlinearVariationalProblem(F, h)
        self._solver = firedrake.NonlinearVariationalSolver(
            problem, solver_parameters=self._solver_parameters)

        self._thickness_old = h_0
        self._timestep = dt
Beispiel #23
0
    def flux(self, **kwargs):
        keys = ('damage', 'velocity', 'damage_inflow')
        keys_alt = ('D', 'u', 'D_inflow')
        D, u, D_inflow = get_kwargs_alt(kwargs, keys, keys_alt)

        Q = D.function_space()
        φ = firedrake.TestFunction(Q)

        mesh = Q.mesh()
        n = firedrake.FacetNormal(mesh)

        u_n = max_value(0, inner(u, n))
        f = D * u_n
        flux_faces = (f('+') - f('-')) * (φ('+') - φ('-')) * dS
        flux_cells = -D * div(u * φ) * dx
        flux_out = D * max_value(0, inner(u, n)) * φ * ds
        flux_in = D_inflow * min_value(0, inner(u, n)) * φ * ds

        return flux_faces + flux_cells + flux_out + flux_in
Beispiel #24
0
    def sources(self, **kwargs):
        keys = ('damage', 'velocity', 'fluidity')
        keys_alt = ('D', 'u', 'A')
        D, u, A = get_kwargs_alt(kwargs, keys, keys_alt)

        # Increase/decrease damage depending on stress and strain rates
        ε = sym(grad(u))
        ε_1 = eigenvalues(ε)[0]

        σ = M(ε, A)
        σ_e = sqrt(inner(σ, σ) - det(σ))

        ε_h = firedrake.Constant(self.healing_strain_rate)
        σ_d = firedrake.Constant(self.damage_stress)
        γ_h = firedrake.Constant(self.healing_rate)
        γ_d = firedrake.Constant(self.damage_rate)

        healing = γ_h * min_value(ε_1 - ε_h, 0)
        fracture = γ_d * conditional(σ_e - σ_d > 0, ε_1, 0.) * (1 - D)

        return healing + fracture
Beispiel #25
0
    def solve(self, dt, h0, a, u, h_inflow=None):
        h_inflow = h_inflow if h_inflow is not None else h0

        Q = h0.function_space()
        h, φ = firedrake.TrialFunction(Q), firedrake.TestFunction(Q)

        ν = facet_normal_2(Q.mesh())
        outflow = firedrake.max_value(inner(u, ν), 0)
        inflow = firedrake.min_value(inner(u, ν), 0)

        flux_cells = -h * inner(u, grad_2(φ)) * dx
        flux_out = h * φ * outflow * ds_v
        F = h * φ * dx + dt * (flux_cells + flux_out)

        accumulation = a * φ * dx
        flux_in = -h_inflow * φ * inflow * ds_v
        A = h0 * φ * dx + dt * (accumulation + flux_in)

        h = h0.copy(deepcopy=True)
        solver_parameters = {'ksp_type': 'preonly', 'pc_type': 'lu'}
        firedrake.solve(F == A, h, solver_parameters=solver_parameters)

        return h
Beispiel #26
0
    def __init__(self, state):
        super().__init__(state)

        # obtain our fields
        self.theta = state.fields('theta')
        self.water_v = state.fields('water_v')
        self.rain = state.fields('rain')
        rho = state.fields('rho')
        try:
            water_c = state.fields('water_c')
            water_l = self.rain + water_c
        except NotImplementedError:
            water_l = self.rain

        # declare function space
        Vt = self.theta.function_space()

        # make rho variables
        # we recover rho into theta space
        if state.vertical_degree == 0 and state.horizontal_degree == 0:
            boundary_method = Boundary_Method.physics
        else:
            boundary_method = None
        Vt_broken = FunctionSpace(state.mesh, BrokenElement(Vt.ufl_element()))
        rho_averaged = Function(Vt)
        self.rho_recoverer = Recoverer(rho,
                                       rho_averaged,
                                       VDG=Vt_broken,
                                       boundary_method=boundary_method)

        # define some parameters as attributes
        dt = state.timestepping.dt
        R_d = state.parameters.R_d
        cp = state.parameters.cp
        cv = state.parameters.cv
        c_pv = state.parameters.c_pv
        c_pl = state.parameters.c_pl
        c_vv = state.parameters.c_vv
        R_v = state.parameters.R_v

        # make useful fields
        Pi = thermodynamics.pi(state.parameters, rho_averaged, self.theta)
        T = thermodynamics.T(state.parameters,
                             self.theta,
                             Pi,
                             r_v=self.water_v)
        p = thermodynamics.p(state.parameters, Pi)
        L_v = thermodynamics.Lv(state.parameters, T)
        R_m = R_d + R_v * self.water_v
        c_pml = cp + c_pv * self.water_v + c_pl * water_l
        c_vml = cv + c_vv * self.water_v + c_pl * water_l

        # use Teten's formula to calculate w_sat
        w_sat = thermodynamics.r_sat(state.parameters, T, p)

        # expression for ventilation factor
        a = Constant(1.6)
        b = Constant(124.9)
        c = Constant(0.2046)
        C = a + b * (rho_averaged * self.rain)**c

        # make appropriate condensation rate
        f = Constant(5.4e5)
        g = Constant(2.55e6)
        h = Constant(0.525)
        dot_r_evap = (((1 - self.water_v / w_sat) * C *
                       (rho_averaged * self.rain)**h) /
                      (rho_averaged * (f + g / (p * w_sat))))

        # make evap_rate function, needs to be the same for all updates in one time step
        evap_rate = Function(Vt)

        # adjust evap rate so negative rain doesn't occur
        self.lim_evap_rate = Interpolator(
            conditional(
                dot_r_evap < 0, 0.0,
                conditional(self.rain < 0.0, 0.0,
                            min_value(dot_r_evap, self.rain / dt))), evap_rate)

        # tell the prognostic fields what to update to
        self.water_v_new = Interpolator(self.water_v + dt * evap_rate, Vt)
        self.rain_new = Interpolator(self.rain - dt * evap_rate, Vt)
        self.theta_new = Interpolator(
            self.theta * (1.0 - dt * evap_rate *
                          (cv * L_v / (c_vml * cp * T) - R_v * cv * c_pml /
                           (R_m * cp * c_vml))), Vt)
Beispiel #27
0
    def __init__(self, state):
        super(Condensation, self).__init__(state)

        # obtain our fields
        self.theta = state.fields('theta')
        self.water_v = state.fields('water_v')
        self.water_c = state.fields('water_c')
        rho = state.fields('rho')

        # declare function space
        Vt = self.theta.function_space()

        param = self.state.parameters

        # define some parameters as attributes
        dt = self.state.timestepping.dt
        R_d = param.R_d
        p_0 = param.p_0
        kappa = param.kappa
        cp = param.cp
        cv = param.cv
        c_pv = param.c_pv
        c_pl = param.c_pl
        c_vv = param.c_vv
        R_v = param.R_v
        L_v0 = param.L_v0
        T_0 = param.T_0
        w_sat1 = param.w_sat1
        w_sat2 = param.w_sat2
        w_sat3 = param.w_sat3
        w_sat4 = param.w_sat4

        # make useful fields
        Pi = ((R_d * rho * self.theta / p_0)**(kappa / (1.0 - kappa)))
        T = Pi * self.theta * R_d / (R_d + self.water_v * R_v)
        p = p_0 * Pi**(1.0 / kappa)
        L_v = L_v0 - (c_pl - c_pv) * (T - T_0)
        R_m = R_d + R_v * self.water_v
        c_pml = cp + c_pv * self.water_v + c_pl * self.water_c
        c_vml = cv + c_vv * self.water_v + c_pl * self.water_c

        # use Teten's formula to calculate w_sat
        w_sat = (w_sat1 / (p * exp(w_sat2 * (T - T_0) /
                                   (T - w_sat3)) - w_sat4))

        # make appropriate condensation rate
        dot_r_cond = ((self.water_v - w_sat) / (dt * (1.0 +
                                                      ((L_v**2.0 * w_sat) /
                                                       (cp * R_v * T**2.0)))))

        # make cond_rate function, that needs to be the same for all updates in one time step
        self.cond_rate = Function(Vt)

        # adjust cond rate so negative concentrations don't occur
        self.lim_cond_rate = Interpolator(
            conditional(dot_r_cond < 0,
                        max_value(dot_r_cond, -self.water_c / dt),
                        min_value(dot_r_cond, self.water_v / dt)),
            self.cond_rate)

        # tell the prognostic fields what to update to
        self.water_v_new = Interpolator(self.water_v - dt * self.cond_rate, Vt)
        self.water_c_new = Interpolator(self.water_c + dt * self.cond_rate, Vt)
        self.theta_new = Interpolator(
            self.theta * (1.0 + dt * self.cond_rate *
                          (cv * L_v / (c_vml * cp * T) - R_v * cv * c_pml /
                           (R_m * cp * c_vml))), Vt)
Beispiel #28
0
def mass_balance(s, max_a, da_ds, ela):
    return min_value((s - ela) * da_ds, max_a)
Beispiel #29
0
    def __init__(self, state, iterations=1):
        super().__init__(state)

        self.iterations = iterations
        # obtain our fields
        self.theta = state.fields('theta')
        self.water_v = state.fields('water_v')
        self.water_c = state.fields('water_c')
        rho = state.fields('rho')
        try:
            rain = state.fields('rain')
            water_l = self.water_c + rain
        except NotImplementedError:
            water_l = self.water_c

        # declare function space
        Vt = self.theta.function_space()

        # make rho variables
        # we recover rho into theta space
        if state.vertical_degree == 0 and state.horizontal_degree == 0:
            boundary_method = Boundary_Method.physics
        else:
            boundary_method = None
        Vt_broken = FunctionSpace(state.mesh, BrokenElement(Vt.ufl_element()))
        rho_averaged = Function(Vt)
        self.rho_recoverer = Recoverer(rho,
                                       rho_averaged,
                                       VDG=Vt_broken,
                                       boundary_method=boundary_method)

        # define some parameters as attributes
        dt = state.timestepping.dt
        R_d = state.parameters.R_d
        cp = state.parameters.cp
        cv = state.parameters.cv
        c_pv = state.parameters.c_pv
        c_pl = state.parameters.c_pl
        c_vv = state.parameters.c_vv
        R_v = state.parameters.R_v

        # make useful fields
        Pi = thermodynamics.pi(state.parameters, rho_averaged, self.theta)
        T = thermodynamics.T(state.parameters,
                             self.theta,
                             Pi,
                             r_v=self.water_v)
        p = thermodynamics.p(state.parameters, Pi)
        L_v = thermodynamics.Lv(state.parameters, T)
        R_m = R_d + R_v * self.water_v
        c_pml = cp + c_pv * self.water_v + c_pl * water_l
        c_vml = cv + c_vv * self.water_v + c_pl * water_l

        # use Teten's formula to calculate w_sat
        w_sat = thermodynamics.r_sat(state.parameters, T, p)

        # make appropriate condensation rate
        dot_r_cond = ((self.water_v - w_sat) / (dt * (1.0 +
                                                      ((L_v**2.0 * w_sat) /
                                                       (cp * R_v * T**2.0)))))

        # make cond_rate function, that needs to be the same for all updates in one time step
        cond_rate = Function(Vt)

        # adjust cond rate so negative concentrations don't occur
        self.lim_cond_rate = Interpolator(
            conditional(dot_r_cond < 0,
                        max_value(dot_r_cond, -self.water_c / dt),
                        min_value(dot_r_cond, self.water_v / dt)), cond_rate)

        # tell the prognostic fields what to update to
        self.water_v_new = Interpolator(self.water_v - dt * cond_rate, Vt)
        self.water_c_new = Interpolator(self.water_c + dt * cond_rate, Vt)
        self.theta_new = Interpolator(
            self.theta * (1.0 + dt * cond_rate *
                          (cv * L_v / (c_vml * cp * T) - R_v * cv * c_pml /
                           (R_m * cp * c_vml))), Vt)
Beispiel #30
0
 def temperature(self, E):
     r"""Return the temperature of ice at the given energy density"""
     return firedrake.min_value(E / (ρ_I * c), Tm)