def get_colour_function(self, k): """ Return the colour function on timestep t^{n+k} """ if k == 0: if self.continuous_fields: c = self.continuous_c else: c = self.simulation.data['c'] elif k == -1: if self.continuous_fields: c = self.continuous_c_old else: c = self.simulation.data['cp'] elif k == -2: if self.continuous_fields: c = self.continuous_c_oldold else: c = self.simulation.data['cpp'] if self.force_bounded: c = dolfin.max_value(dolfin.min_value(c, Constant(1.0)), Constant(0.0)) if self.force_sharp: c = dolfin.conditional(dolfin.ge(c, 0.5), Constant(1.0), Constant(0.0)) return c
def forms_theta_linear(self, psih0, uh, dt, theta_map, theta_L=Constant(1.0), dpsi0=Constant(0.), dpsi00=Constant(0.), h=Constant(0.), neumann_idx=99, zeta=Constant(0)): (psi, lamb, psibar) = self.__trial_functions() (w, tau, wbar) = self.__test_functions() beta_map = self.beta_map n = self.n facet_integral = self.facet_integral psi_star = psih0 + (1-theta_L)*dpsi00 + theta_L * dpsi0 # LHS contributions gamma = conditional(ge(dot(uh, n), 0), 0, 1) # Standard formulation N_a = facet_integral(beta_map*dot(psi, w)) \ + zeta * dot(grad(psi), grad(w)) * dx G_a = dot(lamb, w)/dt * dx - theta_map * dot(uh, grad(lamb))*w*dx \ + theta_map * (1-gamma) * dot(uh, n) * \ lamb * w * self.ds(neumann_idx) L_a = -facet_integral(beta_map * dot(psibar, w)) # \ H_a = facet_integral(dot(uh, n)*psibar * tau) \ - dot(uh, n)*psibar * tau * self.ds(neumann_idx) B_a = facet_integral(beta_map * dot(psibar, wbar)) # RHS contributions Q_a = dot(Constant(0), w) * dx R_a = dot(psi_star, tau)/dt*dx \ + (1-theta_map)*dot(uh, grad(tau))*psih0*dx \ - (1-theta_map) * (1-gamma) * dot(uh, n) * psi_star * tau * self.ds(neumann_idx) \ - gamma*dot(h, tau) * self.ds(neumann_idx) S_a = facet_integral(Constant(0) * wbar) return self.__fem_forms(N_a, G_a, L_a, H_a, B_a, Q_a, R_a, S_a)
def heaviside(x): r""" Heaviside function .. math:: \frac{\mathrm{d}}{\mathrm{d}x} \max\{x,0\} """ return dolfin.conditional(dolfin.ge(x, 0.0), 1.0, 0.0)
def subplus(x): r""" Ramp function .. math:: \max\{x,0\} """ return dolfin.conditional(dolfin.ge(x, 0.0), x, 0.0)
def forms_theta_nlinear_np(self, v0, v_int, Ubar0, dt, theta_map=Constant(1.0), theta_L=Constant(1.0), duh0=Constant((0., 0.)), duh00=Constant((0., 0)), h=Constant((0., 0.)), neumann_idx=99): ''' No particles in mass matrix ''' # Define trial test functions (v, lamb, vbar) = self.__trial_functions() (w, tau, wbar) = self.__test_functions() (zero_vec, h, duh0, duh00) = self.__check_geometric_dimension(h, duh0, duh00) beta_map = self.beta_map n = self.n facet_integral = self.facet_integral # Define v_star v_star = v0 + (1-theta_L) * duh00 + theta_L * duh0 Udiv = v0 + duh0 outer_v_a = outer(w, Udiv) outer_v_a_o = outer(v_star, Udiv) outer_ubar_a = outer(vbar, Ubar0) # Switch to detect in/outflow boundary gamma = conditional(ge(dot(Udiv, n), 0), 0, 1) # LHS contribution s N_a = dot(v, w) * dx + facet_integral(beta_map*dot(v, w)) G_a = dot(lamb, w)/dt * dx - theta_map*inner(outer_v_a, grad(lamb))*dx \ + theta_map * (1-gamma) * dot(outer_v_a * n, lamb) * self.ds(neumann_idx) L_a = -facet_integral(beta_map * dot(vbar, w)) H_a = facet_integral(dot(outer_ubar_a*n, tau)) \ - dot(outer_ubar_a*n, tau) * self.ds(neumann_idx) B_a = facet_integral(beta_map * dot(vbar, wbar)) # RHS contributions Q_a = dot(v_int, w) * dx R_a = dot(v_star, tau)/dt * dx \ + (1-theta_map)*inner(outer_v_a_o, grad(tau))*dx \ - gamma * dot(h, tau) * self.ds(neumann_idx) S_a = facet_integral(dot(Constant(zero_vec), wbar)) return self.__fem_forms(N_a, G_a, L_a, H_a, B_a, Q_a, R_a, S_a)
def I(self, V, s, time=None): m, h, n, NKo, NKi, NNao, NNai, NClo, NCli, voli, O = s O = dolfin.conditional(dolfin.ge(O, 0), O, 0) vol = 1.4368e-15 # unit:m^3, when r=7 um,v=1.4368e-15 m^3 beta0 = 7 volo = (1 + 1 / beta0) * vol - voli # Extracellular volume C = 1 Ko, Ki, Nao, Nai, Clo, Cli = self._get_concentrations(V, s) I_K = self._I_K(V, n, Ko, Ki) I_Na = self._I_Na(V, m, h, Nao, Nai) I_Cl = self._I_Cl(V, Clo, Cli) gamma = self._gamma(voli) I_pump = self._I_pump(Nai, Ko, O, gamma) return (I_Na + I_K + I_Cl + I_pump) / C # Sign is correcct I think
def forms_theta_nlinear_multiphase(self, rho, rho0, rho00, rhobar, v0, Ubar0, dt, theta_map, theta_L=Constant(1.0), duh0=Constant((0., 0.)), duh00=Constant((0., 0)), h=Constant((0., 0.)), neumann_idx=99): (v, lamb, vbar) = self.__trial_functions() (w, tau, wbar) = self.__test_functions() (zero_vec, h, duh0, duh00) = self.__check_geometric_dimension(h, duh0, duh00) beta_map = self.beta_map n = self.n facet_integral = self.facet_integral # FIXME: To be deprecated rhov_star = rho0*(v0 + theta_L * duh0) + rho00 * (1-theta_L) * duh00 Udiv = v0 + duh0 outer_v_a = outer(rho * w, Udiv) outer_v_a_o = outer(rho0 * Udiv, Udiv) outer_ubar_a = outer(rhobar * vbar, Ubar0) # Switch to detect in/outflow boundary gamma = conditional(ge(dot(Udiv, n), 0), 0, 1) # LHS contribution N_a = facet_integral(beta_map*dot(v, w)) G_a = dot(lamb, rho * w)/dt * dx - theta_map*inner(outer_v_a, grad(lamb))*dx \ + theta_map * (1-gamma) * dot(outer_v_a * n, lamb) * self.ds(neumann_idx) L_a = -facet_integral(beta_map * dot(vbar, w)) H_a = facet_integral(dot(outer_ubar_a*n, tau)) \ - dot(outer_ubar_a*n, tau) * self.ds(neumann_idx) B_a = facet_integral(beta_map * dot(vbar, wbar)) # RHS contribution Q_a = dot(zero_vec, w) * dx R_a = dot(rhov_star, tau)/dt*dx \ + (1-theta_map)*inner(outer_v_a_o, grad(tau))*dx \ + gamma * dot(h, tau) * self.ds(neumann_idx) S_a = facet_integral(dot(zero_vec, wbar)) return self.__fem_forms(N_a, G_a, L_a, H_a, B_a, Q_a, R_a, S_a)
def tau(self, nu, u, metric): """ Stabilization parameter """ h2 = h_u2(self.metric, u, self.reg_norm * self.nu * self.nu) Pe = dl.Constant(.5) * h_dot_u(metric, u, self.reg_norm * self.nu * self.nu) / nu num = dl.Constant(1.) + dl.exp(dl.Constant(-2.) * Pe) den = dl.Constant(1.) - dl.exp(dl.Constant(-2.) * Pe) # [0.1 0.01]* [a1] = [ coth(.1) - 1./(.1) ] # [1. 0.2 ] [a2] [ -csch(.1)^2 + 1./(.1)^2] a1 = dl.Constant(0.333554921691650) a2 = dl.Constant(-0.004435991517475) tau_1 = (num / den - dl.Constant(1.) / Pe) * h_over_u( metric, u, self.reg_norm * self.nu * self.nu) tau_2 = (a1 + a2 * Pe) * dl.Constant(.5) * h2 / nu return dl.conditional(dl.ge(Pe, .1), tau_1, tau_2)
def forms_theta_linear( self, psih0, uh, dt, theta_map, theta_L=Constant(1.0), dpsi0=Constant(0.0), dpsi00=Constant(0.0), h=Constant(0.0), neumann_idx=99, zeta=Constant(0), psi_star=None ): """ Set PDEMap forms for a linear advection problem. Parameters ---------- psih0: dolfin.Function dolfin Function storing the solution from the previous step uh: Constant, Expression, dolfin.Function Advective velocity dt: Constant Time step value theta_map: Constant Theta value for time stepping in PDE-projection according to theta-method **NOTE** theta only affects solution for Lagrange multiplier space polynomial order >= 1 theta_L: Constant, optional Theta value for reconstructing intermediate field from the previous solution and old increments. Defaults to Constan(1.) dpsi0: dolfin.Function, optional Increment function from last time step. Defaults to Constant(0) dpsi00: dolfin.Function Increment function from second last time step. Defaults to Constant(0) h: Constant, dolfin.Function, optional Expression or Function for non-homogenous Neumann BC. Defaults to Constant(0.) neumann_idx: int, optional Integer to use for marking Neumann boundaries. Defaults to value 99 zeta: Constant, optional Penalty parameter for limiting over/undershoot. Defaults to 0 psi_star: dolfin.Function Initial state of the advective approximation. If psi_star is None psi_star = psih0 + (1 - theta_L) * dpsi00 + theta_L * dpsi0 will be prescribed. Be warned that this assumes dt does not change between time steps. Returns ------- dict Dictionary with forms """ (psi, lamb, psibar) = self._trial_functions() (w, tau, wbar) = self._test_functions() beta_map = self.beta_map n = self.n facet_integral = self.facet_integral if psi_star is None: psi_star = psih0 + (1 - theta_L) * dpsi00 + theta_L * dpsi0 # LHS contributions gamma = conditional(ge(dot(uh, n), 0), 0, 1) # Standard formulation N_a = facet_integral(beta_map * dot(psi, w)) + zeta * dot(grad(psi), grad(w)) * dx G_a = ( dot(lamb, w) / dt * dx - theta_map * dot(uh, grad(lamb)) * w * dx + theta_map * (1 - gamma) * dot(uh, n) * lamb * w * self.ds(neumann_idx) ) L_a = -facet_integral(beta_map * dot(psibar, w)) # \ H_a = facet_integral(dot(uh, n) * psibar * tau) - dot(uh, n) * psibar * tau * self.ds( neumann_idx ) B_a = facet_integral(beta_map * dot(psibar, wbar)) # RHS contributions Q_a = dot(Constant(0), w) * dx R_a = ( dot(psi_star, tau) / dt * dx + (1 - theta_map) * dot(uh, grad(tau)) * psi_star * dx - (1 - theta_map) * (1 - gamma) * dot(uh, n) * psi_star * tau * self.ds(neumann_idx) - gamma * dot(h, tau) * self.ds(neumann_idx) ) S_a = facet_integral(Constant(0) * wbar) return self._fem_forms(N_a, G_a, L_a, H_a, B_a, Q_a, R_a, S_a)
def forms_theta_nlinear_multiphase( self, rho, rho0, rho00, rhobar, v0, Ubar0, dt, theta_map, theta_L=Constant(1.0), duh0=Constant((0.0, 0.0)), duh00=Constant((0.0, 0)), h=Constant((0.0, 0.0)), neumann_idx=99, ): """ Set PDEMap forms for a non-linear (but linearized) advection problem including density. Parameters ---------- rho: dolfin.Function Current density field rho0: dolfin.Function Density field at previous time step. rho00: dolfin.Function Density field at second last time step rhobar: dolfin.Function Density field at facets v0: dolfin.Function Specific momentum at old time level Ubar0: dolfin.Function Advective field at old time level dt: Constant Time step theta_map: Constant Theta value for time stepping in PDE-projection according to theta-method **NOTE** theta only affects solution for Lagrange multiplier space polynomial order >= 1 theta_L: Constant, optional Theta value for reconstructing intermediate field from the previous solution and old increments. duh0: dolfin.Function, optional Increment from previous time step. duh00: dolfin.Function, optional Increment from second last time step h: Constant, dolfin.Function, optional Expression or Function for non-homogenous Neumann BC. Defaults to Constant(0. neumann_idx: int, optional Integer to use for marking Neumann boundaries. Defaults to value 99 Returns ------- dict Dict with forms """ (v, lamb, vbar) = self._trial_functions() (w, tau, wbar) = self._test_functions() (zero_vec, h, duh0, duh00) = self._check_geometric_dimension(h, duh0, duh00) beta_map = self.beta_map n = self.n facet_integral = self.facet_integral # FIXME: To be deprecated rhov_star = rho0 * (v0 + theta_L * duh0) + rho00 * (1 - theta_L) * duh00 Udiv = v0 + duh0 outer_v_a = outer(rho * w, Udiv) outer_v_a_o = outer(rho0 * Udiv, Udiv) outer_ubar_a = outer(rhobar * vbar, Ubar0) # Switch to detect in/outflow boundary gamma = conditional(ge(dot(Udiv, n), 0), 0, 1) # LHS contribution N_a = facet_integral(beta_map * dot(v, w)) G_a = ( dot(lamb, rho * w) / dt * dx - theta_map * inner(outer_v_a, grad(lamb)) * dx + theta_map * (1 - gamma) * dot(outer_v_a * n, lamb) * self.ds(neumann_idx) ) L_a = -facet_integral(beta_map * dot(vbar, w)) H_a = facet_integral(dot(outer_ubar_a * n, tau)) - dot(outer_ubar_a * n, tau) * self.ds( neumann_idx ) B_a = facet_integral(beta_map * dot(vbar, wbar)) # RHS contribution Q_a = dot(zero_vec, w) * dx R_a = ( dot(rhov_star, tau) / dt * dx + (1 - theta_map) * inner(outer_v_a_o, grad(tau)) * dx + gamma * dot(h, tau) * self.ds(neumann_idx) ) S_a = facet_integral(dot(zero_vec, wbar)) return self._fem_forms(N_a, G_a, L_a, H_a, B_a, Q_a, R_a, S_a)
def forms_theta_nlinear_np( self, v0, v_int, Ubar0, dt, theta_map=Constant(1.0), theta_L=Constant(1.0), duh0=Constant((0.0, 0.0)), duh00=Constant((0.0, 0)), h=Constant((0.0, 0.0)), neumann_idx=99, ): """ Set PDEMap forms for a non-linear (but linearized) advection problem, assumes however that the mass matrix can be obtained from the mesh (and not from particles) **NOTE** Documentation upcoming. """ # Define trial test functions (v, lamb, vbar) = self._trial_functions() (w, tau, wbar) = self._test_functions() (zero_vec, h, duh0, duh00) = self._check_geometric_dimension(h, duh0, duh00) beta_map = self.beta_map n = self.n facet_integral = self.facet_integral # Define v_star v_star = v0 + (1 - theta_L) * duh00 + theta_L * duh0 Udiv = v0 + duh0 outer_v_a = outer(w, Udiv) outer_v_a_o = outer(v_star, Udiv) outer_ubar_a = outer(vbar, Ubar0) # Switch to detect in/outflow boundary gamma = conditional(ge(dot(Udiv, n), 0), 0, 1) # LHS contribution s N_a = dot(v, w) * dx + facet_integral(beta_map * dot(v, w)) G_a = ( dot(lamb, w) / dt * dx - theta_map * inner(outer_v_a, grad(lamb)) * dx + theta_map * (1 - gamma) * dot(outer_v_a * n, lamb) * self.ds(neumann_idx) ) L_a = -facet_integral(beta_map * dot(vbar, w)) H_a = facet_integral(dot(outer_ubar_a * n, tau)) - dot(outer_ubar_a * n, tau) * self.ds( neumann_idx ) B_a = facet_integral(beta_map * dot(vbar, wbar)) # RHS contributions Q_a = dot(v_int, w) * dx R_a = ( dot(v_star, tau) / dt * dx + (1 - theta_map) * inner(outer_v_a_o, grad(tau)) * dx - gamma * dot(h, tau) * self.ds(neumann_idx) ) S_a = facet_integral(dot(Constant(zero_vec), wbar)) return self._fem_forms(N_a, G_a, L_a, H_a, B_a, Q_a, R_a, S_a)
def forms_theta_nlinear( self, v0, Ubar0, dt, theta_map=Constant(1.0), theta_L=Constant(1.0), duh0=Constant((0.0, 0.0)), duh00=Constant((0.0, 0)), h=Constant((0.0, 0.0)), neumann_idx=99, ): """ Set PDEMap forms for a non-linear (but linearized) advection problem, Parameters ---------- v0: dolfin.Function dolfin.Function storing solution from previous step Ubar0: dolfin.Function Advective velocity at facets dt: Constant Time step theta_map: Constant, optional Theta value for time stepping in PDE-projection according to theta-method. Defaults to Constant(1.) **NOTE** theta only affects solution for Lagrange multiplier space polynomial order >= 1 theta_L: Constant, optional Theta value for reconstructing intermediate field from the previous solution and old increments. Defaults to Constant(1.) duh0: dolfin.Function, optional Increment function from last time step duh00: dolfin.Function, optional Increment function from second last time step h: Constant, dolfin.Function, optional Expression or Function for non-homogenous Neumann BC. Defaults to Constant(0.) neumann_idx: int, optional Integer to use for marking Neumann boundaries. Defaults to value 99 Returns ------- dict Dictionary with forms """ # Define trial test functions (v, lamb, vbar) = self._trial_functions() (w, tau, wbar) = self._test_functions() (zero_vec, h, duh0, duh00) = self._check_geometric_dimension(h, duh0, duh00) beta_map = self.beta_map n = self.n facet_integral = self.facet_integral # Define v_star v_star = v0 + (1 - theta_L) * duh00 + theta_L * duh0 Udiv = v0 + duh0 outer_v_a = outer(w, Udiv) outer_v_a_o = outer(v_star, Udiv) outer_ubar_a = outer(vbar, Ubar0) # Switch to detect in/outflow boundary gamma = conditional(ge(dot(Udiv, n), 0), 0, 1) # LHS contribution s N_a = facet_integral(beta_map * dot(v, w)) G_a = ( dot(lamb, w) / dt * dx - theta_map * inner(outer_v_a, grad(lamb)) * dx + theta_map * (1 - gamma) * dot(outer_v_a * n, lamb) * self.ds(neumann_idx) ) L_a = -facet_integral(beta_map * dot(vbar, w)) H_a = facet_integral(dot(outer_ubar_a * n, tau)) - dot(outer_ubar_a * n, tau) * self.ds( neumann_idx ) B_a = facet_integral(beta_map * dot(vbar, wbar)) # RHS contributions Q_a = dot(zero_vec, w) * dx R_a = ( dot(v_star, tau) / dt * dx + (1 - theta_map) * inner(outer_v_a_o, grad(tau)) * dx - gamma * dot(h, tau) * self.ds(neumann_idx) ) S_a = facet_integral(dot(zero_vec, wbar)) return self._fem_forms(N_a, G_a, L_a, H_a, B_a, Q_a, R_a, S_a)
def F(self, V, s, time=None): m, h, n, NKo, NKi, NNao, NNai, NClo, NCli, voli, O = s O = dolfin.conditional(dolfin.ge(O, 0), O, 0) G_K = self._parameters["G_K"] G_Na = self._parameters["G_Na"] G_ClL = self._parameters["G_ClL"] G_Kl = self._parameters["G_KL"] G_NaL = self._parameters["G_NaL"] C = self._parameters["C"] # Ion Concentration related Parameters eps_K = self._parameters["eps_K"] G_glia = self._parameters["G_glia"] eps_O = self._parameters["eps_O"] Ukcc2 = self._parameters["Ukcc2"] Unkcc1 = self._parameters["Unkcc1"] rho_max = self._parameters["rho_max"] # Volume vol = 1.4368e-15 # unit:m^3, when r=7 um,v=1.4368e-15 m^3 # TODO: add to parameters beta0 = self._parameters["beta0"] # Time Constant tau = self._parameters["tau"] KBath = self._parameters["KBath"] OBath = self._parameters["OBath"] gamma = self._gamma(voli) volo = (1 + 1 / beta0) * vol - voli # Extracellular volume beta = voli / volo # Ratio of intracelluar to extracelluar volume # Gating variables alpha_m = 0.32 * (54 + V) / (1 - exp(-(V + 54) / 4)) beta_m = 0.28 * (V + 27) / (exp((V + 27) / 5) - 1) alpha_h = 0.128 * exp(-(V + 50) / 18) beta_h = 4 / (1 + exp(-(V + 27) / 5)) alpha_n = 0.032 * (V + 52) / (1 - exp(-(V + 52) / 5)) beta_n = 0.5 * exp(-(V + 57) / 40) dotm = alpha_m * (1 - m) - beta_m * m doth = alpha_h * (1 - h) - beta_h * h dotn = alpha_n * (1 - n) - beta_n * n Ko, Ki, Nao, Nai, Clo, Cli = self._get_concentrations(V, s) fo = 1 / (1 + exp((2.5 - OBath) / 0.2)) fv = 1 / (1 + exp((beta - 20) / 2)) eps_K *= fo * fv G_glia *= fo rho = rho_max / (1 + exp((20 - O) / 3)) / gamma I_glia = G_glia / (1 + exp((18 - Ko) / 2.5)) Igliapump = rho / 3 / (1 + exp((25 - 18) / 3)) / (1 + exp(3.5 - Ko)) I_diff = eps_K * (Ko - KBath) + I_glia + 2 * Igliapump * gamma I_K = self._I_K(V, n, Ko, Ki) I_Na = self._I_Na(V, m, h, Nao, Nai) I_Cl = self._I_Cl(V, Clo, Cli) I_pump = self._I_pump(Nai, Ko, O, gamma) # Cloride transporter (mM/s) fKo = 1 / (1 + exp(16 - Ko)) FKCC2 = Ukcc2 * ln((Ki * Cli) / (Ko * Clo)) FNKCC1 = Unkcc1 * fKo * (ln((Ki * Cli) / (Ko * Clo)) + ln( (Nai * Cli) / (Nao * Clo))) dotNKo = tau * volo * (gamma * beta * (I_K - 2.0 * I_pump) - I_diff + FKCC2 * beta + FNKCC1 * beta) dotNKi = -tau * voli * (gamma * (I_K - 2.0 * I_pump) + FKCC2 + FNKCC1) dotNNao = tau * volo * beta * (gamma * (I_Na + 3.0 * I_pump) + FNKCC1) dotNNai = -tau * voli * (gamma * (I_Na + 3.0 * I_pump) + FNKCC1) dotNClo = beta * volo * tau * (FKCC2 - gamma * I_Cl + 2 * FNKCC1) dotNCli = voli * tau * (gamma * I_Cl - FKCC2 - 2 * FNKCC1) r1 = vol / voli r2 = 1 / beta0 * vol / ((1 + 1 / beta0) * vol - voli) pii = Nai + Cli + Ki + 132 * r1 pio = Nao + Ko + Clo + 18 * r2 vol_hat = vol * (1.1029 - 0.1029 * exp((pio - pii) / 20)) dotVoli = -(voli - vol_hat) / 0.25 * tau dotO = tau * (-5.3 * (I_pump + Igliapump) * gamma + eps_O * (OBath - O)) F_expressions = [ufl.zero()] * self.num_states() F_expressions[0] = dotm F_expressions[1] = doth F_expressions[2] = dotn F_expressions[3] = dotNKo F_expressions[4] = dotNKi F_expressions[5] = dotNNao F_expressions[6] = dotNNai F_expressions[7] = dotNClo F_expressions[8] = dotNCli F_expressions[9] = dotVoli F_expressions[10] = dotO return dolfin.as_vector(F_expressions)
def perm_update_wong_newton(vs, phi0, k0): mult_min = 1e-15 mult = conditional(ge(pow(1.0+vs/phi0,3.0)/(1.0+vs),0.0) \ ,pow(1.0+vs/phi0,3.0)/(1.0+vs),mult_min) k = k0 * mult return k
def diff_coeff_cal_rev(p, p0, tau, phi): mult_min = 1e-15 mult = conditional(ge(phi / tau, 0.0), phi / tau, mult_min) p = p0 * mult return p
def test_moving_mesh(): t = 0. dt = 0.025 num_steps = 20 xmin, ymin = 0., 0. xmax, ymax = 2., 2. xc, yc = 1., 1. nx, ny = 20, 20 pres = 150 k = 1 mesh = RectangleMesh(Point(xmin, ymin), Point(xmax, ymax), nx, ny) n = FacetNormal(mesh) # Class for mesh motion dU = PeriodicVelocity(xmin, xmax, dt, t, degree=1) Qcg = VectorFunctionSpace(mesh, 'CG', 1) boundaries = MeshFunction("size_t", mesh, mesh.topology().dim()-1) boundaries.set_all(0) leftbound = Left(xmin) leftbound.mark(boundaries, 99) ds = Measure('ds', domain=mesh, subdomain_data=boundaries) # Create function spaces Q_E_Rho = FiniteElement("DG", mesh.ufl_cell(), k) T_1 = FunctionSpace(mesh, 'DG', 0) Qbar_E = FiniteElement("DGT", mesh.ufl_cell(), k) Q_Rho = FunctionSpace(mesh, Q_E_Rho) Qbar = FunctionSpace(mesh, Qbar_E) phih, phih0 = Function(Q_Rho), Function(Q_Rho) phibar = Function(Qbar) # Advective velocity uh = Function(Qcg) uh.assign(Constant((0., 0.))) # Mesh velocity umesh = Function(Qcg) # Total velocity uadvect = uh-umesh # Now throw in the particles x = RandomRectangle(Point(xmin, ymin), Point(xmax, ymax)).generate([pres, pres]) s = assign_particle_values(x, GaussianPulse(center=(xc, yc), sigma=float(0.25), U=[0, 0], time=0., height=1., degree=3)) x = comm.bcast(x, root=0) s = comm.bcast(s, root=0) p = particles(x, [s], mesh) # Define projections problem FuncSpace_adv = {'FuncSpace_local': Q_Rho, 'FuncSpace_lambda': T_1, 'FuncSpace_bar': Qbar} FormsPDE = FormsPDEMap(mesh, FuncSpace_adv, ds=ds) forms_pde = FormsPDE.forms_theta_linear(phih0, uadvect, dt, Constant(1.0), zeta=Constant(0.)) pde_projection = PDEStaticCondensation(mesh, p, forms_pde['N_a'], forms_pde['G_a'], forms_pde['L_a'], forms_pde['H_a'], forms_pde['B_a'], forms_pde['Q_a'], forms_pde['R_a'], forms_pde['S_a'], [], 1) # Initialize the initial condition at mesh by an l2 projection lstsq_rho = l2projection(p, Q_Rho, 1) lstsq_rho.project(phih0.cpp_object()) for step in range(num_steps): # Compute old area at old configuration old_area = assemble(phih0*dx) # Pre-assemble rhs pde_projection.assemble_state_rhs() # Move mesh dU.compute_ubc() umesh.assign(project(dU, Qcg)) ALE.move(mesh, project(dU * dt, Qcg)) dU.update() # Relocate particles as a result of mesh motion # NOTE: if particles were advected themselve, # we had to run update_facets_info() here as well p.relocate() # Assemble left-hand side on new config, but not the right-hand side pde_projection.assemble(True, False) pde_projection.solve_problem(phibar.cpp_object(), phih.cpp_object(), 'mumps', 'none') # Needed to compute conservation, note that there # is an outgoing flux at left boundary new_area = assemble(phih*dx) gamma = conditional(ge(dot(uadvect, n), 0), 0, 1) bflux = assemble((1-gamma) * dot(uadvect, n) * phih * ds) # Update solution assign(phih0, phih) # Put assertion on (global) mass balance, local mass balance is # too time consuming but should pass also assert new_area - old_area + bflux * dt < 1e-12 # Assert that max value of phih stays close to 2 and # min value close to 0. This typically will fail if # we do not do a correct relocate of particles assert np.amin(phih.vector().get_local()) > -0.015 assert np.amax(phih.vector().get_local()) < 1.04