def test_two_blocks_with_restriction(mesh, restriction, FunctionSpaces, BlockBCs): (FunctionSpace1, FunctionSpace2) = FunctionSpaces V1 = FunctionSpace1(mesh) V2 = FunctionSpace2(mesh) block_V = BlockFunctionSpace([V1, V2], restrict=restriction) block_form = get_rhs_block_form_2(block_V) block_bcs = BlockBCs(block_V) (rhs, block_rhs) = assemble_and_block_assemble_vector(block_form) apply_bc_and_block_bc_vector(rhs, block_rhs, block_bcs) assert_block_vectors_equal(rhs, block_rhs, block_V) (rhs, block_rhs) = assemble_and_block_assemble_vector(block_form) (function, block_function) = apply_bc_and_block_bc_vector_non_linear( rhs, block_rhs, block_bcs, block_V) assert_block_vectors_equal(rhs, block_rhs, block_V) assert_block_functions_equal(function, block_function, block_V)
def test_single_block_no_restriction_from_block_element(mesh, Element): V_element = Element(mesh) V = FunctionSpace(mesh, V_element) block_V_element = BlockElement(V_element) block_V = BlockFunctionSpace(mesh, block_V_element) assert_dof_map_single_block_no_restriction(V, block_V)
def test_single_block_no_restriction_from_list(mesh, FunctionSpace): V = FunctionSpace(mesh) block_V = BlockFunctionSpace([V]) assert_dof_map_single_block_no_restriction(V, block_V)
def test_two_blocks_with_restriction_from_list(mesh, restriction, FunctionSpaces): (FunctionSpace1, FunctionSpace2) = FunctionSpaces V1 = FunctionSpace1(mesh) V2 = FunctionSpace2(mesh) block_V = BlockFunctionSpace([V1, V2], restrict=restriction) assert_dof_map_test_two_blocks_with_restriction(V1, V2, block_V)
def test_single_block_with_restriction(mesh, restriction, FunctionSpace): V = FunctionSpace(mesh) block_V = BlockFunctionSpace([V], restrict=[restriction]) functions = get_list_of_functions_1(block_V) assert_functions_manipulations(functions, block_V)
def initialization(mesh, subdomains, boundaries): TM = TensorFunctionSpace(mesh, 'DG', 0) PM = FunctionSpace(mesh, 'DG', 0) UCG = VectorElement("CG", mesh.ufl_cell(), 2) BDM = FiniteElement("BDM", mesh.ufl_cell(), 1) PDG = FiniteElement("DG", mesh.ufl_cell(), 0) UCG_F = FunctionSpace(mesh, UCG) BDM_F = FunctionSpace(mesh, BDM) PDG_F = FunctionSpace(mesh, PDG) W = BlockFunctionSpace([BDM_F, PDG_F], restrict=[None, None]) U = BlockFunctionSpace([UCG_F]) I = Identity(mesh.topology().dim()) C_cg = FiniteElement("CG", mesh.ufl_cell(), 1) C_dg = FiniteElement("DG", mesh.ufl_cell(), 0) mini = C_cg + C_dg C = FunctionSpace(mesh, mini) C = BlockFunctionSpace([C]) #TODO solution0_h = BlockFunction(W) solution0_m = BlockFunction(U) solution0_c = BlockFunction(C) solution1_h = BlockFunction(W) solution1_m = BlockFunction(U) solution1_c = BlockFunction(C) solution2_h = BlockFunction(W) solution2_m = BlockFunction(U) solution2_c = BlockFunction(C) solution_h = BlockFunction(W) solution_m = BlockFunction(U) solution_c = BlockFunction(C) ## mechanics # 0 properties alpha1 = 0.74 K1 = 8.4 * 1000.e6 nu1 = 0.18 alpha2 = 0.74 K2 = 8.4 * 1000.e6 nu2 = 0.18 alpha_values = [alpha1, alpha2] K_values = [K1, K2] nu_values = [nu1, nu2] alpha_0 = Function(PM) K_0 = Function(PM) nu_0 = Function(PM) alpha_0 = init_scalar_parameter(alpha_0, alpha_values[0], 500, subdomains) K_0 = init_scalar_parameter(K_0, K_values[0], 500, subdomains) nu_0 = init_scalar_parameter(nu_0, nu_values[0], 500, subdomains) alpha_0 = init_scalar_parameter(alpha_0, alpha_values[1], 501, subdomains) K_0 = init_scalar_parameter(K_0, K_values[1], 501, subdomains) nu_0 = init_scalar_parameter(nu_0, nu_values[1], 501, subdomains) K_mult_min = 1.0 K_mult_max = 1.0 mu_l_0, lmbda_l_0, Ks_0, K_0 = \ bulk_modulus_update(mesh,solution0_c[0],K_mult_min,K_mult_max,K_0,nu_0,alpha_0,K_0) # n-1 properties alpha1 = 0.74 K1 = 8.4 * 1000.e6 nu1 = 0.18 alpha2 = 0.74 K2 = 8.4 * 1000.e6 nu2 = 0.18 alpha_values = [alpha1, alpha2] K_values = [K1, K2] nu_values = [nu1, nu2] alpha_1 = Function(PM) K_1 = Function(PM) nu_1 = Function(PM) alpha_1 = init_scalar_parameter(alpha_1, alpha_values[0], 500, subdomains) K_1 = init_scalar_parameter(K_1, K_values[0], 500, subdomains) nu_1 = init_scalar_parameter(nu_1, nu_values[0], 500, subdomains) alpha_1 = init_scalar_parameter(alpha_1, alpha_values[1], 501, subdomains) K_1 = init_scalar_parameter(K_1, K_values[1], 501, subdomains) nu_1 = init_scalar_parameter(nu_1, nu_values[1], 501, subdomains) K_mult_min = 1.0 K_mult_max = 1.0 mu_l_1, lmbda_l_1, Ks_1, K_1 = \ bulk_modulus_update(mesh,solution0_c[0],K_mult_min,K_mult_max,K_1,nu_1,alpha_1,K_0) # n properties alpha1 = 0.74 K2 = 8.4 * 1000.e6 nu1 = 0.18 alpha2 = 0.74 K2 = 8.4 * 1000.e6 nu2 = 0.18 alpha_values = [alpha1, alpha2] K_values = [K1, K2] nu_values = [nu1, nu2] alpha = Function(PM) K = Function(PM) nu = Function(PM) alpha = init_scalar_parameter(alpha, alpha_values[0], 500, subdomains) K = init_scalar_parameter(K, K_values[0], 500, subdomains) nu = init_scalar_parameter(nu, nu_values[0], 500, subdomains) alpha = init_scalar_parameter(alpha, alpha_values[1], 501, subdomains) K = init_scalar_parameter(K, K_values[1], 501, subdomains) nu = init_scalar_parameter(nu, nu_values[1], 501, subdomains) K_mult_min = 1.0 K_mult_max = 1.0 mu_l, lmbda_l, Ks, K = \ bulk_modulus_update(mesh,solution0_c[0],K_mult_min,K_mult_max,K,nu,alpha,K_0) ## flow # 0 properties cf1 = 1e-10 phi1 = 0.2 rho1 = 1000.0 mu1 = 1. kx = 8.802589710965712e-10 ky = 8.802589710965712e-11 k1 = np.array([kx, 0., 0., ky]) cf2 = 1e-10 phi2 = 0.2 rho2 = 1000.0 mu2 = 1. kx = 8.802589710965712e-10 ky = 8.802589710965712e-11 k2 = np.array([kx, 0., 0., ky]) cf_values = [cf1, cf2] phi_values = [phi1, phi2] rho_values = [rho1, rho2] mu_values = [mu1, mu2] k_values = [k1, k2] cf_0 = Function(PM) phi_0 = Function(PM) rho_0 = Function(PM) mu_0 = Function(PM) k_0 = Function(TM) cf_0 = init_scalar_parameter(cf_0, cf_values[0], 500, subdomains) phi_0 = init_scalar_parameter(phi_0, phi_values[0], 500, subdomains) rho_0 = init_scalar_parameter(rho_0, rho_values[0], 500, subdomains) mu_0 = init_scalar_parameter(mu_0, mu_values[0], 500, subdomains) k_0 = init_tensor_parameter(k_0, k_values[0], 500, subdomains, mesh.topology().dim()) cf_0 = init_scalar_parameter(cf_0, cf_values[1], 501, subdomains) phi_0 = init_scalar_parameter(phi_0, phi_values[1], 501, subdomains) rho_0 = init_scalar_parameter(rho_0, rho_values[1], 501, subdomains) mu_0 = init_scalar_parameter(mu_0, mu_values[1], 501, subdomains) k_0 = init_tensor_parameter(k_0, k_values[1], 501, subdomains, mesh.topology().dim()) #filename = "perm4.csv" #k_0 = init_from_file_parameter(k_0,0.,0.,filename) # n-1 properties cf1 = 1e-10 phi1 = 0.2 rho1 = 1000.0 mu1 = 1. kx = 8.802589710965712e-10 ky = 8.802589710965712e-11 k1 = np.array([kx, 0., 0., ky]) cf2 = 1e-10 phi2 = 0.2 rho2 = 1000.0 mu2 = 1. kx = 8.802589710965712e-10 ky = 8.802589710965712e-11 k2 = np.array([kx, 0., 0., ky]) cf_values = [cf1, cf2] phi_values = [phi1, phi2] rho_values = [rho1, rho2] mu_values = [mu1, mu2] k_values = [k1, k2] cf_1 = Function(PM) phi_1 = Function(PM) rho_1 = Function(PM) mu_1 = Function(PM) k_1 = Function(TM) cf_1 = init_scalar_parameter(cf_1, cf_values[0], 500, subdomains) phi_1 = init_scalar_parameter(phi_1, phi_values[0], 500, subdomains) rho_1 = init_scalar_parameter(rho_1, rho_values[0], 500, subdomains) mu_1 = init_scalar_parameter(mu_1, mu_values[0], 500, subdomains) k_1 = init_tensor_parameter(k_1, k_values[0], 500, subdomains, mesh.topology().dim()) cf_1 = init_scalar_parameter(cf_1, cf_values[1], 501, subdomains) phi_1 = init_scalar_parameter(phi_1, phi_values[1], 501, subdomains) rho_1 = init_scalar_parameter(rho_1, rho_values[1], 501, subdomains) mu_1 = init_scalar_parameter(mu_1, mu_values[1], 501, subdomains) k_1 = init_tensor_parameter(k_1, k_values[1], 501, subdomains, mesh.topology().dim()) #filename = "perm4.csv" #k_1 = init_from_file_parameter(k_1,0.,0.,filename) # n properties cf1 = 1e-10 phi1 = 0.2 rho1 = 1000.0 mu1 = 1. kx = 8.802589710965712e-10 ky = 8.802589710965712e-11 k1 = np.array([kx, 0., 0., ky]) cf2 = 1e-10 phi2 = 0.2 rho2 = 1000.0 mu2 = 1. kx = 8.802589710965712e-10 ky = 8.802589710965712e-11 k2 = np.array([kx, 0., 0., ky]) cf_values = [cf1, cf2] phi_values = [phi1, phi2] rho_values = [rho1, rho2] mu_values = [mu1, mu2] k_values = [k1, k2] cf = Function(PM) phi = Function(PM) rho = Function(PM) mu = Function(PM) k = Function(TM) cf = init_scalar_parameter(cf, cf_values[0], 500, subdomains) phi = init_scalar_parameter(phi, phi_values[0], 500, subdomains) rho = init_scalar_parameter(rho, rho_values[0], 500, subdomains) mu = init_scalar_parameter(mu, mu_values[0], 500, subdomains) k = init_tensor_parameter(k, k_values[0], 500, subdomains, mesh.topology().dim()) cf = init_scalar_parameter(cf, cf_values[1], 501, subdomains) phi = init_scalar_parameter(phi, phi_values[1], 501, subdomains) rho = init_scalar_parameter(rho, rho_values[1], 501, subdomains) mu = init_scalar_parameter(mu, mu_values[1], 501, subdomains) k = init_tensor_parameter(k, k_values[1], 501, subdomains, mesh.topology().dim()) #filename = "perm4.csv" #k = init_from_file_parameter(k,0.,0.,filename) ### transport # 0 dx1 = 1e-12 dy1 = 1e-12 d1 = np.array([dx1, 0., 0., dy1]) dx2 = 1e-12 dy2 = 1e-12 d2 = np.array([dx2, 0., 0., dy2]) d_values = [d1, d2] d_0 = Function(TM) d_0 = init_tensor_parameter(d_0, d_values[0], 500, subdomains, mesh.topology().dim()) d_0 = init_tensor_parameter(d_0, d_values[1], 501, subdomains, mesh.topology().dim()) # n-1 dx1 = 1e-12 dy1 = 1e-12 d1 = np.array([dx1, 0., 0., dy1]) dx2 = 1e-12 dy2 = 1e-12 d2 = np.array([dx2, 0., 0., dy2]) d_values = [d1, d2] d_1 = Function(TM) d_1 = init_tensor_parameter(d_1, d_values[0], 500, subdomains, mesh.topology().dim()) d_1 = init_tensor_parameter(d_1, d_values[1], 501, subdomains, mesh.topology().dim()) # n dx1 = 1e-12 dy1 = 1e-12 d1 = np.array([dx1, 0., 0., dy1]) dx2 = 1e-12 dy2 = 1e-12 d2 = np.array([dx2, 0., 0., dy2]) d_values = [d1, d2] d = Function(TM) d = init_tensor_parameter(d, d_values[0], 500, subdomains, mesh.topology().dim()) d = init_tensor_parameter(d, d_values[1], 501, subdomains, mesh.topology().dim()) ####initialization # initial u_0 = Constant((0.0, 0.0)) u_0_project = project(u_0, U[0]) assign(solution0_m.sub(0), u_0_project) p_0 = Constant(1.e6) p_0_project = project(p_0, W[1]) assign(solution0_h.sub(1), p_0_project) # v_0 = Constant((0.0, 0.0)) # v_0_project = project(v_0, W[0]) # assign(solution0_h.sub(0), v_0_project) c0 = c_sat_cal(1.e6, 20.) c0_project = project(c0, C[0]) assign(solution0_c.sub(0), c0_project) # n - 1 u_0 = Constant((0.0, 0.0)) u_0_project = project(u_0, U[0]) assign(solution1_m.sub(0), u_0_project) p_0 = Constant(1.e6) p_0_project = project(p_0, W[1]) assign(solution1_h.sub(1), p_0_project) # v_0 = Constant((0.0, 0.0)) # v_0_project = project(v_0, W[0]) # assign(solution1_h.sub(0), v_0_project) c0 = c_sat_cal(1.e6, 20.) c0_project = project(c0, C[0]) assign(solution1_c.sub(0), c0_project) # n - 2 u_0 = Constant((0.0, 0.0)) u_0_project = project(u_0, U[0]) assign(solution2_m.sub(0), u_0_project) p_0 = Constant(1.e6) p_0_project = project(p_0, W[1]) assign(solution2_h.sub(1), p_0_project) # v_0 = Constant((0.0, 0.0)) # v_0_project = project(v_0, W[0]) # assign(solution2_h.sub(0), v_0_project) c0 = c_sat_cal(1.e6, 20.) c0_project = project(c0, C[0]) assign(solution2_c.sub(0), c0_project) # n u_0 = Constant((0.0, 0.0)) u_0_project = project(u_0, U[0]) assign(solution_m.sub(0), u_0_project) p_0 = Constant(1.e6) p_0_project = project(p_0, W[1]) assign(solution_h.sub(1), p_0_project) # v_0 = Constant((0.0, 0.0)) # v_0_project = project(v_0, W[0]) # assign(solution_h.sub(0), v_0_project) c0 = c_sat_cal(1.e6, 20.) c0_project = project(c0, C[0]) assign(solution_c.sub(0), c0_project) ###iterative parameters phi_it = Function(PM) assign(phi_it, phi_0) print('c_sat', c_sat_cal(1.0e8, 20.)) c_sat = c_sat_cal(1.0e8, 20.) c_sat = project(c_sat, PM) c_inject = Constant(0.0) c_inject = project(c_inject, PM) mu_c1_1 = 1.e-4 mu_c2_1 = 5.e-0 mu_c1_2 = 1.e-4 mu_c2_2 = 5.e-0 mu_c1_values = [mu_c1_1, mu_c1_2] mu_c2_values = [mu_c2_1, mu_c2_2] mu_c1 = Function(PM) mu_c2 = Function(PM) mu_c1 = init_scalar_parameter(mu_c1, mu_c1_values[0], 500, subdomains) mu_c2 = init_scalar_parameter(mu_c2, mu_c2_values[0], 500, subdomains) mu_c1 = init_scalar_parameter(mu_c1, mu_c1_values[1], 501, subdomains) mu_c2 = init_scalar_parameter(mu_c2, mu_c2_values[1], 501, subdomains) coeff_for_perm_1 = 22.2 coeff_for_perm_2 = 22.2 coeff_for_perm_values = [coeff_for_perm_1, coeff_for_perm_2] coeff_for_perm = Function(PM) coeff_for_perm = init_scalar_parameter(coeff_for_perm, coeff_for_perm_values[0], 500, subdomains) coeff_for_perm = init_scalar_parameter(coeff_for_perm, coeff_for_perm_values[1], 501, subdomains) solutionIt_h = BlockFunction(W) return solution0_m, solution0_h, solution0_c \ ,solution1_m, solution1_h, solution1_c \ ,solution2_m, solution2_h, solution2_c \ ,solution_m, solution_h, solution_c \ ,alpha_0, K_0, mu_l_0, lmbda_l_0, Ks_0 \ ,alpha_1, K_1, mu_l_1, lmbda_l_1, Ks_1 \ ,alpha, K, mu_l, lmbda_l, Ks \ ,cf_0, phi_0, rho_0, mu_0, k_0 \ ,cf_1, phi_1, rho_1, mu_1, k_1 \ ,cf, phi, rho, mu, k \ ,d_0, d_1, d, I \ ,phi_it, solutionIt_h, mu_c1, mu_c2 \ ,nu_0, nu_1, nu, coeff_for_perm \ ,c_sat, c_inject
def transport_linear(integrator_type, mesh, subdomains, boundaries, t_start, dt, T, solution0, \ alpha_0, K_0, mu_l_0, lmbda_l_0, Ks_0, \ alpha_1, K_1, mu_l_1, lmbda_l_1, Ks_1, \ alpha, K, mu_l, lmbda_l, Ks, \ cf_0, phi_0, rho_0, mu_0, k_0,\ cf_1, phi_1, rho_1, mu_1, k_1,\ cf, phi, rho, mu, k, \ d_0, d_1, d_t, vel_c, p_con, A_0, Temp, c_extrapolate): # Create mesh and define function space parameters["ghost_mode"] = "shared_facet" # required by dS dx = Measure('dx', domain=mesh, subdomain_data=subdomains) ds = Measure('ds', domain=mesh, subdomain_data=boundaries) dS = Measure('dS', domain=mesh, subdomain_data=boundaries) C_cg = FiniteElement("CG", mesh.ufl_cell(), 1) C_dg = FiniteElement("DG", mesh.ufl_cell(), 0) mini = C_cg + C_dg C = FunctionSpace(mesh, mini) C = BlockFunctionSpace([C]) TM = TensorFunctionSpace(mesh, 'DG', 0) PM = FunctionSpace(mesh, 'DG', 0) n = FacetNormal(mesh) vc = CellVolume(mesh) fc = FacetArea(mesh) h = vc / fc h_avg = (vc('+') + vc('-')) / (2 * avg(fc)) penalty1 = Constant(1.0) tau = Function(PM) tau = tau_cal(tau, phi, -0.5) tuning_para = 0.25 vel_norm = (dot(vel_c, n) + abs(dot(vel_c, n))) / 2.0 cell_size = CellDiameter(mesh) vnorm = sqrt(dot(vel_c, vel_c)) I = Identity(mesh.topology().dim()) d_eff = Function(TM) d_eff = diff_coeff_cal_rev(d_eff, d_0, tau, phi) + tuning_para * cell_size * vnorm * I monitor_dt = dt # Define variational problem dc, = BlockTrialFunction(C) dc_dot, = BlockTrialFunction(C) psic, = BlockTestFunction(C) block_c = BlockFunction(C) c, = block_split(block_c) block_c_dot = BlockFunction(C) c_dot, = block_split(block_c_dot) theta = -1.0 a_time = phi * rho * inner(c_dot, psic) * dx a_dif = dot(rho*d_eff*grad(c),grad(psic))*dx \ - dot(avg_w(rho*d_eff*grad(c),weight_e(rho*d_eff,n)), jump(psic, n))*dS \ + theta*dot(avg_w(rho*d_eff*grad(psic),weight_e(rho*d_eff,n)), jump(c, n))*dS \ + penalty1/h_avg*k_e(rho*d_eff,n)*dot(jump(c, n), jump(psic, n))*dS a_adv = -dot(rho*vel_c*c,grad(psic))*dx \ + dot(jump(psic), rho('+')*vel_norm('+')*c('+') - rho('-')*vel_norm('-')*c('-') )*dS \ + dot(psic, rho*vel_norm*c)*ds(3) R_c = R_c_cal(c_extrapolate, p_con, Temp) c_D1 = Constant(0.5) rhs_c = R_c * A_s_cal(phi, phi_0, A_0) * psic * dx - dot( rho * phi * vel_c, n) * c_D1 * psic * ds(1) r_u = [a_dif + a_adv] j_u = block_derivative(r_u, [c], [dc]) r_u_dot = [a_time] j_u_dot = block_derivative(r_u_dot, [c_dot], [dc_dot]) r = [r_u_dot[0] + r_u[0] - rhs_c] # this part is not applied. exact_solution_expression1 = Expression("1.0", t=0, element=C[0].ufl_element()) def bc(t): p5 = DirichletBC(C.sub(0), exact_solution_expression1, boundaries, 1, method="geometric") return BlockDirichletBC([p5]) # Define problem wrapper class ProblemWrapper(object): def set_time(self, t): pass # Residual and jacobian functions def residual_eval(self, t, solution, solution_dot): return r def jacobian_eval(self, t, solution, solution_dot, solution_dot_coefficient): return [[ Constant(solution_dot_coefficient) * j_u_dot[0, 0] + j_u[0, 0] ]] # Define boundary condition def bc_eval(self, t): pass # Define initial condition def ic_eval(self): return solution0 # Define custom monitor to plot the solution def monitor(self, t, solution, solution_dot): pass problem_wrapper = ProblemWrapper() (solution, solution_dot) = (block_c, block_c_dot) solver = TimeStepping(problem_wrapper, solution, solution_dot) solver.set_parameters({ "initial_time": t_start, "time_step_size": dt, "monitor": { "time_step_size": monitor_dt, }, "final_time": T, "exact_final_time": "stepover", "integrator_type": integrator_type, "problem_type": "linear", "linear_solver": "mumps", "report": True }) export_solution = solver.solve() return export_solution, T
def m_linear(integrator_type, mesh, subdomains, boundaries, t_start, dt, T, solution0, \ alpha_0, K_0, mu_l_0, lmbda_l_0, Ks_0, \ alpha_1, K_1, mu_l_1, lmbda_l_1, Ks_1, \ alpha, K, mu_l, lmbda_l, Ks, \ cf_0, phi_0, rho_0, mu_0, k_0,\ cf_1, phi_1, rho_1, mu_1, k_1,\ cf, phi, rho, mu, k, \ pressure_freeze): # Create mesh and define function space parameters["ghost_mode"] = "shared_facet" # required by dS dx = Measure('dx', domain=mesh, subdomain_data=subdomains) ds = Measure('ds', domain=mesh, subdomain_data=boundaries) dS = Measure('dS', domain=mesh, subdomain_data=boundaries) C = VectorFunctionSpace(mesh, "CG", 2) C = BlockFunctionSpace([C]) TM = TensorFunctionSpace(mesh, 'DG', 0) PM = FunctionSpace(mesh, 'DG', 0) n = FacetNormal(mesh) vc = CellVolume(mesh) fc = FacetArea(mesh) h = vc/fc h_avg = (vc('+') + vc('-'))/(2*avg(fc)) monitor_dt = dt f_stress_x = Constant(-1.e3) f_stress_y = Constant(-20.0e6) f = Constant((0.0, 0.0)) #sink/source for displacement I = Identity(mesh.topology().dim()) # Define variational problem psiu, = BlockTestFunction(C) block_u = BlockTrialFunction(C) u, = block_split(block_u) w = BlockFunction(C) theta = -1.0 a_time = inner(-alpha*pressure_freeze*I,sym(grad(psiu)))*dx #quasi static a = inner(2*mu_l*strain(u)+lmbda_l*div(u)*I, sym(grad(psiu)))*dx rhs_a = inner(f,psiu)*dx \ + dot(f_stress_y*n,psiu)*ds(2) r_u = [a] #DirichletBC bcd1 = DirichletBC(C.sub(0).sub(0), 0.0, boundaries, 1) # No normal displacement for solid on left side bcd3 = DirichletBC(C.sub(0).sub(0), 0.0, boundaries, 3) # No normal displacement for solid on right side bcd4 = DirichletBC(C.sub(0).sub(1), 0.0, boundaries, 4) # No normal displacement for solid on bottom side bcs = BlockDirichletBC([bcd1,bcd3,bcd4]) AA = block_assemble([r_u]) FF = block_assemble([rhs_a - a_time]) bcs.apply(AA) bcs.apply(FF) block_solve(AA, w.block_vector(), FF, "mumps") export_solution = w return export_solution, T
boundaries = MeshFunction("size_t", mesh, "data/biot_2mat_facet_region.xml") for discretization in ("DG", "EG"): # Block function space V_element = VectorElement("CG", mesh.ufl_cell(), 2) if discretization == "DG": Q_element = FiniteElement("DG", mesh.ufl_cell(), 1) W_element = BlockElement(V_element, Q_element) elif discretization == "EG": Q_element = FiniteElement("CG", mesh.ufl_cell(), 1) D_element = FiniteElement("DG", mesh.ufl_cell(), 0) EG_element = Q_element + D_element W_element = BlockElement(V_element, EG_element) else: raise RuntimeError("Invalid discretization") W = BlockFunctionSpace(mesh, W_element) PM = FunctionSpace(mesh, "DG", 0) TM = TensorFunctionSpace(mesh, "DG", 0) I = Identity(mesh.topology().dim()) dx = Measure("dx", domain=mesh, subdomain_data=subdomains) ds = Measure("ds", domain=mesh, subdomain_data=boundaries) dS = Measure("dS", domain=mesh, subdomain_data=boundaries) # Test and trial functions vq = BlockTestFunction(W) (v, q) = block_split(vq) up = BlockTrialFunction(W) (u, p) = block_split(up)
def h_linear(integrator_type, mesh, subdomains, boundaries, t_start, dt, T, solution0, \ alpha_0, K_0, mu_l_0, lmbda_l_0, Ks_0, \ alpha_1, K_1, mu_l_1, lmbda_l_1, Ks_1, \ alpha, K, mu_l, lmbda_l, Ks, \ cf_0, phi_0, rho_0, mu_0, k_0,\ cf_1, phi_1, rho_1, mu_1, k_1,\ cf, phi, rho, mu, k, \ sigma_v_freeze, dphi_c_dt): # Create mesh and define function space parameters["ghost_mode"] = "shared_facet" # required by dS dx = Measure('dx', domain=mesh, subdomain_data=subdomains) ds = Measure('ds', domain=mesh, subdomain_data=boundaries) dS = Measure('dS', domain=mesh, subdomain_data=boundaries) BDM = FiniteElement("BDM", mesh.ufl_cell(), 1) PDG = FiniteElement("DG", mesh.ufl_cell(), 0) BDM_F = FunctionSpace(mesh, BDM) PDG_F = FunctionSpace(mesh, PDG) W = BlockFunctionSpace([BDM_F, PDG_F], restrict=[None, None]) TM = TensorFunctionSpace(mesh, 'DG', 0) PM = FunctionSpace(mesh, 'DG', 0) n = FacetNormal(mesh) vc = CellVolume(mesh) fc = FacetArea(mesh) h = vc / fc h_avg = (vc('+') + vc('-')) / (2 * avg(fc)) I = Identity(mesh.topology().dim()) monitor_dt = dt p_outlet = 0.1e6 p_inlet = 1000.0 M_inv = phi_0 * cf + (alpha - phi_0) / Ks # Define variational problem trial = BlockTrialFunction(W) dv, dp = block_split(trial) trial_dot = BlockTrialFunction(W) dv_dot, dp_dot = block_split(trial_dot) test = BlockTestFunction(W) psiv, psip = block_split(test) block_w = BlockFunction(W) v, p = block_split(block_w) block_w_dot = BlockFunction(W) v_dot, p_dot = block_split(block_w_dot) a_time = Constant(0.0) * inner(v_dot, psiv) * dx #quasi static # k is a function of phi #k = perm_update_rutqvist_newton(p,p0,phi0,phi,coeff) lhs_a = inner(dot(v, mu * inv(k)), psiv) * dx - p * div( psiv ) * dx #+ 6.0*inner(psiv,n)*ds(2) # - inner(gravity*(rho-rho0), psiv)*dx b_time = (M_inv + pow(alpha, 2.) / K) * p_dot * psip * dx lhs_b = div(v) * psip * dx #div(rho*v)*psip*dx #TODO rho rhs_v = -p_outlet * inner(psiv, n) * ds(3) rhs_p = -alpha / K * sigma_v_freeze * psip * dx - dphi_c_dt * psip * dx r_u = [lhs_a, lhs_b] j_u = block_derivative(r_u, block_w, trial) r_u_dot = [a_time, b_time] j_u_dot = block_derivative(r_u_dot, block_w_dot, trial_dot) r = [r_u_dot[0] + r_u[0] - rhs_v, \ r_u_dot[1] + r_u[1] - rhs_p] def bc(t): #bc_v = [DirichletBC(W.sub(0), (.0, .0), boundaries, 4)] v1 = DirichletBC(W.sub(0), (1.e-4 * 2.0, 0.0), boundaries, 1) v2 = DirichletBC(W.sub(0), (0.0, 0.0), boundaries, 2) v4 = DirichletBC(W.sub(0), (0.0, 0.0), boundaries, 4) bc_v = [v1, v2, v4] return BlockDirichletBC([bc_v, None]) # Define problem wrapper class ProblemWrapper(object): def set_time(self, t): pass #g.t = t # Residual and jacobian functions def residual_eval(self, t, solution, solution_dot): #print(as_backend_type(assemble(p_time - p_time_error)).vec().norm()) #print("gravity effect", as_backend_type(assemble(inner(gravity*(rho-rho0), psiv)*dx)).vec().norm()) return r def jacobian_eval(self, t, solution, solution_dot, solution_dot_coefficient): return [[Constant(solution_dot_coefficient)*j_u_dot[0, 0] + j_u[0, 0], \ Constant(solution_dot_coefficient)*j_u_dot[0, 1] + j_u[0, 1]], \ [Constant(solution_dot_coefficient)*j_u_dot[1, 0] + j_u[1, 0], \ Constant(solution_dot_coefficient)*j_u_dot[1, 1] + j_u[1, 1]]] # Define boundary condition def bc_eval(self, t): return bc(t) # Define initial condition def ic_eval(self): return solution0 # Define custom monitor to plot the solution def monitor(self, t, solution, solution_dot): pass # Solve the time dependent problem problem_wrapper = ProblemWrapper() (solution, solution_dot) = (block_w, block_w_dot) solver = TimeStepping(problem_wrapper, solution, solution_dot) solver.set_parameters({ "initial_time": t_start, "time_step_size": dt, "monitor": { "time_step_size": monitor_dt, }, "final_time": T, "exact_final_time": "stepover", "integrator_type": integrator_type, "problem_type": "linear", "linear_solver": "mumps", "report": True }) export_solution = solver.solve() return export_solution, T