def compute(self, state): super().compute(state) return self.field.assign(thermodynamics.RH(state.parameters, self.r_v, self.T, self.p))
def unsaturated_hydrostatic_balance(state, theta_d, H, exner0=None, top=False, exner_boundary=Constant(1.0), max_outer_solve_count=40, max_inner_solve_count=20): """ Given vertical profiles for dry potential temperature and relative humidity compute hydrostatically balanced virtual potential temperature, dry density and water vapour profiles. The general strategy is to split up the solving into two steps: 1) finding rho to balance the theta profile 2) finding theta_v and r_v to get back theta_d and H We iteratively solve these steps until we (hopefully) converge to a solution. :arg state: The :class:`State` object. :arg theta_d: The initial dry potential temperature profile. :arg H: The relative humidity profile. :arg exner0: Optional function to put exner pressure into. :arg top: If True, set a boundary condition at the top, otherwise it will be at the bottom. :arg exner_boundary: The value of exner on the specified boundary. :arg max_outer_solve_count: Max number of iterations for outer loop of balance solver. :arg max_inner_solve_count: Max number of iterations for inner loop of balanace solver. """ theta0 = state.fields('theta') rho0 = state.fields('rho') mr_v0 = state.fields('vapour_mixing_ratio') # Calculate hydrostatic exner pressure Vt = theta0.function_space() Vr = rho0.function_space() R_d = state.parameters.R_d R_v = state.parameters.R_v epsilon = R_d / R_v VDG = state.spaces("DG") if any(deg > 2 for deg in VDG.ufl_element().degree()): logger.warning( "default quadrature degree most likely not sufficient for this degree element" ) # apply first guesses theta0.assign(theta_d * 1.01) mr_v0.assign(0.01) v_deg = Vr.ufl_element().degree()[1] if v_deg == 0: method = Boundary_Method.physics else: method = None rho_h = Function(Vr) rho_averaged = Function(Vt) Vt_broken = FunctionSpace(state.mesh, BrokenElement(Vt.ufl_element())) rho_recoverer = Recoverer(rho0, rho_averaged, VDG=Vt_broken, boundary_method=method) w_h = Function(Vt) delta = 1.0 # make expressions for determining mr_v0 exner = thermodynamics.exner_pressure(state.parameters, rho_averaged, theta0) p = thermodynamics.p(state.parameters, exner) T = thermodynamics.T(state.parameters, theta0, exner, mr_v0) r_v_expr = thermodynamics.r_v(state.parameters, H, T, p) # make expressions to evaluate residual exner_ev = thermodynamics.exner_pressure(state.parameters, rho_averaged, theta0) p_ev = thermodynamics.p(state.parameters, exner_ev) T_ev = thermodynamics.T(state.parameters, theta0, exner_ev, mr_v0) RH_ev = thermodynamics.RH(state.parameters, mr_v0, T_ev, p_ev) RH = Function(Vt) for i in range(max_outer_solve_count): # solve for rho with theta_vd and w_v guesses compressible_hydrostatic_balance(state, theta0, rho_h, top=top, exner_boundary=exner_boundary, mr_t=mr_v0, solve_for_rho=True) # damp solution rho0.assign(rho0 * (1 - delta) + delta * rho_h) # calculate averaged rho rho_recoverer.project() RH.assign(RH_ev) if errornorm(RH, H) < 1e-10: break # now solve for r_v for j in range(max_inner_solve_count): w_h.interpolate(r_v_expr) mr_v0.assign(mr_v0 * (1 - delta) + delta * w_h) # compute theta_vd theta0.assign(theta_d * (1 + mr_v0 / epsilon)) # test quality of solution by re-evaluating expression RH.assign(RH_ev) if errornorm(RH, H) < 1e-10: break if i == max_outer_solve_count: raise RuntimeError( 'Hydrostatic balance solve has not converged within %i' % i, 'iterations') if exner0 is not None: exner = thermodynamics.exner_pressure(state.parameters, rho0, theta0) exner0.interpolate(exner) # do one extra solve for rho compressible_hydrostatic_balance(state, theta0, rho0, top=top, exner_boundary=exner_boundary, mr_t=mr_v0, solve_for_rho=True)