예제 #1
0
    def dHdc(self, I, L, c):
        """
    Returns the derivative of the Hamiltonian consisting of ajoint-computed
    self.Lam values w.r.t. the control variable <c>, i.e., 

       dH    d [                 ]
       -- = -- [ I + L(self.Lam) ]
       dc   dc [                 ]

    """
        s = "::: forming dHdc :::"
        print_text(s, self.color())

        dU = self.get_dU()
        Lam = self.get_Lam()

        # we need to evaluate the Hamiltonian with the values of Lam computed from
        # self.dI in order to get the derivative of the Hamiltonian w.r.t. the
        # control variables.  Hence we need a new Lagrangian with the trial
        # functions replaced with the computed Lam values.
        L_lam = replace(L, {dU: Lam})

        # the Hamiltonian with unknowns replaced with computed Lam :
        H_lam = I + L_lam

        # the derivative of the Hamiltonian w.r.t. the control variables in the
        # direction of a P1 test function :
        return derivative(H_lam, c, TestFunction(self.model.Q))
예제 #2
0
    def solve_pressure(self, annotate=False):
        """
    Solve for the Dukowicz BP pressure to model.p.
    """
        model = self.model

        # solve for vertical velocity :
        s = "::: solving Dukowicz BP pressure :::"
        print_text(s, self.color())

        Q = model.Q
        rhoi = model.rhoi
        g = model.g
        S = model.S
        z = model.x[2]
        p = model.p
        eta_shf = self.eta_shf
        eta_gnd = self.eta_gnd
        w = self.w

        p_shf = project(rhoi * g * (S - z) + 2 * eta_shf * w.dx(2),
                        annotate=annotate)
        p_gnd = project(rhoi * g * (S - z) + 2 * eta_gnd * w.dx(2),
                        annotate=annotate)

        # unify the pressure over shelves and grounded ice :
        p_v = p.vector().array()
        p_gnd_v = p_gnd.vector().array()
        p_shf_v = p_shf.vector().array()
        p_v[model.gnd_dofs] = p_gnd_v[model.gnd_dofs]
        p_v[model.shf_dofs] = p_shf_v[model.shf_dofs]
        model.assign_variable(p, p_v, cls=self)
예제 #3
0
파일: d1model.py 프로젝트: JacobDowns/cslvr
    def vert_integrate(self, u):
        """
    Integrate <u> from the surface to the bed.
    """
        s = '::: vertically integrating function :::'
        print_text(s, self.D1Model_color)

        ff = self.ff
        Q = self.Q
        phi = TestFunction(Q)
        v = TrialFunction(Q)

        # surface Dirichlet boundary :
        def surface(x, on_boundary):
            return on_boundary and x[0] == self.S_bc

        # integral is zero on surface
        bcs = DirichletBC(Q, 0.0, surface)
        a = v.dx(0) * phi * dx
        L = u * phi * dx
        v = Function(Q)
        solve(a == L, v, bcs, annotate=False)
        print_min_max(v, 'vertically integrated function')
        #print_min_max(v, 'vertically integrated %s' % u.name())
        return v
예제 #4
0
    def solve(self):
        """
    """
        s = "::: solving 'FS_Balance' :::"
        print_text(s, self.color())

        model = self.model

        # solve the linear system :
        solve(self.M, model.tau_ii.vector(), assemble(self.tau_ii))
        print_min_max(model.tau_ii, 'tau_ii')
        solve(self.M, model.tau_ij.vector(), assemble(self.tau_ij))
        print_min_max(model.tau_ij, 'tau_ij')
        solve(self.M, model.tau_ik.vector(), assemble(self.tau_ik))
        print_min_max(model.tau_ik, 'tau_ik')
        solve(self.M, model.tau_ji.vector(), assemble(self.tau_ji))
        print_min_max(model.tau_ji, 'tau_ji')
        solve(self.M, model.tau_jj.vector(), assemble(self.tau_jj))
        print_min_max(model.tau_jj, 'tau_jj')
        solve(self.M, model.tau_jk.vector(), assemble(self.tau_jk))
        print_min_max(model.tau_jk, 'tau_jk')
        solve(self.M, model.tau_ki.vector(), assemble(self.tau_ki))
        print_min_max(model.tau_ki, 'tau_ki')
        solve(self.M, model.tau_kj.vector(), assemble(self.tau_kj))
        print_min_max(model.tau_kj, 'tau_kj')
        solve(self.M, model.tau_kk.vector(), assemble(self.tau_kk))
        print_min_max(model.tau_kk, 'tau_kk')
예제 #5
0
  def solve(self, annotate=False):
    """ 
    Perform the Newton solve of the full-Stokes equations 
    """
    model  = self.model
    params = self.solve_params
    
    # solve nonlinear system :
    rtol   = params['solver']['newton_solver']['relative_tolerance']
    maxit  = params['solver']['newton_solver']['maximum_iterations']
    alpha  = params['solver']['newton_solver']['relaxation_parameter']
    s      = "::: solving BP horizontal velocity with %i max iterations" + \
             " and step size = %.1f :::"
    print_text(s % (maxit, alpha), self.color())
    
    # zero out self.velocity for good convergence for any subsequent solves,
    # e.g. model.L_curve() :
    model.assign_variable(self.get_U(), DOLFIN_EPS, cls=self)
    
    # compute solution :
    solve(self.mom_F == 0, self.G, J = self.mom_Jac, bcs = self.mom_bcs,
          annotate = annotate, solver_parameters = params['solver'])
    U, p = self.G.split()
    u, v, w = split(U)

    self.assx.assign(model.u, project(u, model.Q, annotate=False),
                     annotate=False)
    self.assy.assign(model.v, project(v, model.Q, annotate=False),
                     annotate=False)
    self.assz.assign(model.w, project(w, model.Q, annotate=False),
                     annotate=False)
    self.assp.assign(model.p, p, annotate=False)

    print_min_max(model.U3, 'U', cls=self)
    print_min_max(model.p,  'p', cls=self)
예제 #6
0
    def calc_functionals(self):
        """
    Used to facilitate printing the objective function in adjoint solves.
    """
        s = "::: calculating functionals :::"
        print_text(s, cls=self)

        ftnls = []

        R = assemble(self.Rp, annotate=False)
        print_min_max(R, 'R', cls=self)
        ftnls.append(R)

        J = assemble(self.Jp, annotate=False)
        print_min_max(J, 'J', cls=self)
        ftnls.append(J)

        if self.obj_ftn_type == 'log_L2_hybrid':
            J1 = assemble(self.J1, annotate=False)
            print_min_max(J1, 'J1', cls=self)
            ftnls.append(J1)

            J2 = assemble(self.J2, annotate=False)
            print_min_max(J2, 'J2', cls=self)
            ftnls.append(J2)

        if self.reg_ftn_type == 'TV_Tik_hybrid':
            R1 = assemble(self.R1p, annotate=False)
            print_min_max(R1, 'R1', cls=self)
            ftnls.append(R1)

            R2 = assemble(self.R2p, annotate=False)
            print_min_max(R2, 'R2', cls=self)
            ftnls.append(R2)
        return ftnls
예제 #7
0
  def calc_T_melt(self, annotate=True):
    """
    Calculates pressure-melting point in model.T_melt.
    """
    s    = "::: calculating pressure-melting temperature :::"
    print_text(s, cls=self)

    model = self.model

    dx    = model.dx
    gamma = model.gamma
    T_w   = model.T_w
    p     = model.p

    u   = TrialFunction(model.Q)
    phi = TestFunction(model.Q)
    Tm  = Function(model.Q)

    l = assemble((T_w - gamma * p) * phi * dx, annotate=annotate)
    a = assemble(u * phi * dx, annotate=annotate)

    solve(a, Tm.vector(), l, annotate=annotate)
    model.assign_variable(model.T_melt, Tm, cls=self)

    Tm_v    = Tm.vector().array()
    theta_m = 146.3*Tm_v + 7.253/2.0*Tm_v**2
    model.assign_variable(model.theta_melt, theta_m, cls=self)
예제 #8
0
 def save_lat_mesh(self, h5File): 
   """
   save the lateral boundary mesh to hdf5 file <h5File>.
   """
   s = "::: writing 'latmesh' to supplied hdf5 file :::"
   print_text(s, cls=self.this)
   h5File.write(self.latmesh, 'latmesh')
예제 #9
0
파일: d1model.py 프로젝트: JacobDowns/cslvr
 def init_rho(self, rho):
     """
 """
     s = "::: initializing density :::"
     print_text(s, self.D1Model_color)
     self.assign_variable(self.rho, rho)
     print_min_max(self.rho, 'rho')
예제 #10
0
 def save_bed_mesh(self, h5File): 
   """
   save the basal boundary mesh to hdf5 file <h5File>.
   """
   s = "::: writing 'bedmesh' to supplied hdf5 file :::"
   print_text(s, cls=self.this)
   h5File.write(self.bedmesh, 'bedmesh')
예제 #11
0
 def save_srf_mesh(self, h5File): 
   """
   save the surface boundary mesh to hdf5 file <h5File>.
   """
   s = "::: writing 'srfmesh' to supplied hdf5 file :::"
   print_text(s, cls=self.this)
   h5File.write(self.srfmesh, 'srfmesh')
예제 #12
0
  def set_subdomains_3(self, ff, cf, ff_acc):
    """
    Set the facet subdomains to FacetFunction <ff>, and set the cell subdomains 
    to CellFunction <cf>, and accumulation FacetFunction to <ff_acc>.
    """
    s = "::: setting 3D subdomains :::"
    print_text(s, cls=self)

    self.ff     = ff
    self.cf     = cf
    self.ff_acc = ff_acc
    self.ds     = Measure('ds')[self.ff]
    self.dx     = Measure('dx')[self.cf]
    
    self.dx_g    = self.dx(0)                # internal above grounded
    self.dx_f    = self.dx(1)                # internal above floating
    self.dBed_g  = self.ds(3)                # grounded bed
    self.dBed_f  = self.ds(5)                # floating bed
    self.dBed    = self.ds(3) + self.ds(5)   # bed
    self.dSrf_g  = self.ds(2)                # surface of grounded ice
    self.dSrf_f  = self.ds(6)                # surface of floating ice
    self.dSrf    = self.ds(6) + self.ds(2)   # surface
    self.dLat_d  = self.ds(7)                # lateral divide
    self.dLat_to = self.ds(4)                # lateral terminus overwater
    self.dLat_tu = self.ds(10)               # lateral terminus underwater
    self.dLat_t  = self.ds(4) + self.ds(10)  # lateral terminus
    self.dLat    =   self.ds(4) + self.ds(7) \
                   + self.ds(10)             # lateral
예제 #13
0
  def generate_stokes_function_spaces(self, kind='mini'):
    """
    Generates the appropriate finite-element function spaces from parameters
    specified in the config file for the model.

    If <kind> == 'mini', use enriched mini elements.
    If <kind> == 'th',   use Taylor-Hood elements.
    """
    # mini elements :
    if kind == 'mini':
      s = "::: generating 'mini' Stokes function spaces :::"
        
      self.Bub    = FunctionSpace(self.mesh, "B", 4, 
                                  constrained_domain=self.pBC)
      self.MQ     = self.Q + self.Bub
      M3          = MixedFunctionSpace([self.MQ]*3)
      self.Q4     = MixedFunctionSpace([M3, self.Q])
      self.Q5     = MixedFunctionSpace([M3, self.Q, self.Q])

    # Taylor-Hood elements :
    elif kind == 'th':
      s = "::: generating 'Taylor-Hood' Stokes function spaces :::"
      V           = VectorFunctionSpace(self.mesh, "CG", 2,
                                        constrained_domain=self.pBC)
      self.Q4     = V * self.Q
      self.Q5     = V * self.Q * self.Q
    
    else:
      s = ">>> METHOD generate_stokes_function_spaces <kind> FIELD <<<\n" + \
          ">>> MAY BE 'mini' OR 'th', NOT '%s'. <<<" % kind
      print_text(s, 'red', 1)
      sys.exit(1)

    print_text(s, cls=self)
예제 #14
0
  def solve(self, annotate=False):
    """ 
    Perform the Newton solve of the full-Stokes equations 
    """
    model  = self.model
    params = self.solve_params
    
    # solve nonlinear system :
    rtol   = params['solver']['newton_solver']['relative_tolerance']
    maxit  = params['solver']['newton_solver']['maximum_iterations']
    alpha  = params['solver']['newton_solver']['relaxation_parameter']
    s    = "::: solving Dukowicz-Brinkerhoff-full-Stokes equations" + \
           " with %i max iterations and step size = %.1f :::"
    print_text(s % (maxit, alpha), self.color())
    
    # zero out self.velocity for good convergence for any subsequent solves,
    # e.g. model.L_curve() :
    model.assign_variable(self.get_U(), DOLFIN_EPS, cls=self)
    
    # compute solution :
    solve(self.mom_F == 0, self.U, J = self.mom_Jac, bcs = self.mom_bcs,
          annotate = annotate, solver_parameters = params['solver'])
    u, v, w, p = self.U.split()
    
    self.assx.assign(model.u, u, annotate=False)
    self.assy.assign(model.v, v, annotate=False)
    self.assz.assign(model.w, w, annotate=False)
    self.assp.assign(model.p, p, annotate=False)
    
    U3 = model.U3.split(True)

    print_min_max(U3[0],   'u', cls=self)
    print_min_max(U3[1],   'v', cls=self)
    print_min_max(U3[2],   'w', cls=self)
    print_min_max(model.p, 'p', cls=self)
예제 #15
0
파일: d1model.py 프로젝트: JacobDowns/cslvr
 def init_r(self, r):
     """
 """
     s = "::: initializing grain-size :::"
     print_text(s, self.D1Model_color)
     self.assign_variable(self.r, r)
     print_min_max(self.r, 'r^2')
예제 #16
0
 def save_dvd_mesh(self, h5File): 
   """
   save the divide boundary mesh to hdf5 file <h5File>.
   """
   s = "::: writing 'dvdmesh' to supplied hdf5 file :::"
   print_text(s, cls=self.this)
   h5File.write(self.dvdmesh, 'dvdmesh')
예제 #17
0
  def solve_basal_melt_rate(self):
    """
    Solve for the basal melt rate stored in model.Mb.
    """ 
    # calculate melt-rate : 
    s = "::: solving for basal melt-rate :::"
    print_text(s, cls=self)
    
    model    = self.model
    B        = model.B
    rhoi     = model.rhoi
    theta    = model.theta
    T        = model.T
    T_melt   = model.T_melt
    L        = model.L
    q_geo    = model.q_geo
    rho      = self.rho
    k        = self.k
    c        = self.c
    q_fric   = self.q_fric

    gradB = as_vector([B.dx(0), B.dx(1), -1])
    dTdn  = k/c * dot(grad(theta), gradB)
    nMb   = project((q_geo + q_fric - dTdn) / (L*rho), model.Q,
                    annotate=False)
    nMb_v    = nMb.vector().array()
    T_melt_v = T_melt.vector().array()
    T_v      = T.vector().array()
    nMb_v[T_v < T_melt_v] = 0.0    # if frozen, no melt
    nMb_v[model.shf_dofs] = 0.0    # does apply over floating regions
    model.assign_variable(model.Mb, nMb_v, cls=self)
예제 #18
0
    def vert_integrate(self, u, d='up', Q='self'):
        """
    Integrate <u> from the bed to the surface.
    """
        s = "::: vertically integrating function :::"
        print_text(s, cls=self)

        if type(Q) != FunctionSpace:
            Q = self.Q
        ff = self.ff
        phi = TestFunction(Q)
        v = TrialFunction(Q)
        bcs = []
        # integral is zero on bed (ff = 3,5)
        if d == 'up':
            bcs.append(DirichletBC(Q, 0.0, ff, self.GAMMA_B_GND))  # grounded
            bcs.append(DirichletBC(Q, 0.0, ff, self.GAMMA_B_FLT))  # shelves
            a = v.dx(2) * phi * dx
        # integral is zero on surface (ff = 2,6)
        elif d == 'down':
            bcs.append(DirichletBC(Q, 0.0, ff, self.GAMMA_S_GND))  # grounded
            bcs.append(DirichletBC(Q, 0.0, ff, self.GAMMA_S_FLT))  # shelves
            bcs.append(DirichletBC(Q, 0.0, ff, self.GAMMA_U_GND))  # grounded
            bcs.append(DirichletBC(Q, 0.0, ff, self.GAMMA_U_FLT))  # shelves
            a = -v.dx(2) * phi * dx
        L = u * phi * dx
        name = 'value integrated %s' % d
        v = Function(Q, name=name)
        solve(a == L, v, bcs, annotate=False)
        print_min_max(u, 'vertically integrated function', cls=self)
        return v
예제 #19
0
파일: d1model.py 프로젝트: JacobDowns/cslvr
 def init_B_bc(self, B_bc):
     """
 Set the scalar-value for bed.
 """
     s = "::: initializng bed boundary condition :::"
     print_text(s, self.D1Model_color)
     self.B_bc = B_bc
     print_min_max(self.B_bc, 'B_bc')
예제 #20
0
    def __init__(self, mesh, out_dir='./results/', use_periodic=False):
        """
    Create and instance of a 2D model.
    """
        s = "::: INITIALIZING 2D MODEL :::"
        print_text(s, cls=self)

        Model.__init__(self, mesh, out_dir, use_periodic)
예제 #21
0
파일: d1model.py 프로젝트: JacobDowns/cslvr
 def init_S_bc(self, S_bc):
     """
 Set the scalar-value for the surface. 
 """
     s = "::: initializng surface boundary condition :::"
     print_text(s, self.D1Model_color)
     self.S_bc = S_bc
     print_min_max(self.S_bc, 'S_bc')
예제 #22
0
 def init_U(self, U, cls=None):
     """
 """
     # NOTE: this overides model.init_U
     if cls is None:
         cls = self.this
     s = "::: initializing velocity :::"
     print_text(s, cls=cls)
     self.assign_variable(self.U3, U, cls=cls)
예제 #23
0
 def Hamiltonian(self, I):
     """
 Returns the Hamiltonian of the momentum equations with objective function
 <I>.
 """
     s = "::: forming Hamiltonian :::"
     print_text(s, self.color())
     # the Hamiltonian :
     return I + self.Lagrangian()
예제 #24
0
 def init_H_bounds(self, H_min, H_max, cls=None):
     """
 """
     if cls is None:
         cls = self
     s = "::: initializing bounds on thickness :::"
     print_text(s, cls=cls)
     self.assign_variable(self.H_min, H_min, cls=cls)
     self.assign_variable(self.H_max, H_max, cls=cls)
예제 #25
0
 def deviatoric_stress_tensor(self):
     """
 return the Cauchy stress tensor.
 """
     s = "::: forming the deviatoric part of the Cauchy stress tensor :::"
     print_text(s, self.color)
     epi = self.strain_rate_tensor()
     tau = 2 * self.eta * epi
     return tau
예제 #26
0
파일: d1model.py 프로젝트: JacobDowns/cslvr
 def init_adot(self, adot):
     """
 """
     s = "::: initializing accumulation :::"
     print_text(s, self.D1Model_color)
     self.assign_variable(self.adot, adot)
     self.assign_variable(self.bdot, self.rhoi(0) * adot / self.spy(0))
     print_min_max(self.adot, 'adot')
     print_min_max(self.bdot, 'bdot')
예제 #27
0
 def init_H_H0(self, H, cls=None):
     """
 """
     if cls is None:
         cls = self
     s = "::: initializing thickness :::"
     print_text(s, cls=cls)
     self.assign_variable(self.H, H, cls=cls)
     self.assign_variable(self.H0, H, cls=cls)
예제 #28
0
파일: d1model.py 프로젝트: JacobDowns/cslvr
    def __init__(self, mesh, out_dir='./results/', use_periodic=False):
        """
    Create and instance of a 1D model.
    """
        self.D1Model_color = '150'

        s = "::: INITIALIZING 1D MODEL :::"
        print_text(s, self.D1Model_color)

        Model.__init__(self, mesh, out_dir, use_periodic)
예제 #29
0
 def calc_bulk_density(self):
   """
   Calculate the bulk density stored in model.rho_b.
   """
   # calculate bulk density :
   s = "::: calculating bulk density :::"
   print_text(s, cls=self)
   model       = self.model
   rho_b       = project(self.rho, annotate=False)
   model.assign_variable(model.rhob, rho_b, cls=self)
예제 #30
0
    def __init__(self, model):
        """
    """
        s = "::: INITIALIZING FIRN MASS-BALANCE PHYSICS :::"
        print_text(s, self.color())

        if type(model) != D1Model:
            s = ">>> FirnMass REQUIRES A 'D1Model' INSTANCE, NOT %s <<<"
            print_text(s % type(model), 'red', 1)
            sys.exit(1)