def test_uniform_displacement(self): g_list = setup_grids.setup_2d() for g in g_list: bound_faces = np.argwhere( np.abs(g.cell_faces).sum(axis=1).A.ravel("F") == 1 ) bound = bc.BoundaryCondition( g, bound_faces.ravel("F"), ["dir"] * bound_faces.size ) constit = setup_stiffness(g) # Python inverter is most efficient for small problems stress, bound_stress = mpsa.mpsa(g, constit, bound, inverter="python") div = fvutils.vector_divergence(g) a = div * stress d_x = np.random.rand(1) d_y = np.random.rand(1) d_bound = np.zeros((g.dim, g.num_faces)) d_bound[0, bound.is_dir] = d_x d_bound[1, bound.is_dir] = d_y rhs = div * bound_stress * d_bound.ravel("F") d = np.linalg.solve(a.todense(), -rhs) traction = stress * d + bound_stress * d_bound.ravel("F") self.assertTrue(np.max(np.abs(d[::2] - d_x)) < 1e-8) self.assertTrue(np.max(np.abs(d[1::2] - d_y)) < 1e-8) self.assertTrue(np.max(np.abs(traction)) < 1e-8)
def test_uniform_displacement_neumann(self): physdims = [1, 1] g_size = [4, 8] g_list = [structured.CartGrid([n, n], physdims=physdims) for n in g_size] [g.compute_geometry() for g in g_list] error = [] for g in g_list: bot = np.ravel(np.argwhere(g.face_centers[1, :] < 1e-10)) left = np.ravel(np.argwhere(g.face_centers[0, :] < 1e-10)) dir_faces = np.hstack((left, bot)) bound = bc.BoundaryCondition( g, dir_faces.ravel("F"), ["dir"] * dir_faces.size ) constit = setup_stiffness(g) # Python inverter is most efficient for small problems stress, bound_stress = mpsa.mpsa(g, constit, bound, inverter="python") div = fvutils.vector_divergence(g) a = div * stress d_x = np.random.rand(1) d_y = np.random.rand(1) d_bound = np.zeros((g.dim, g.num_faces)) d_bound[0, bound.is_dir] = d_x d_bound[1, bound.is_dir] = d_y rhs = div * bound_stress * d_bound.ravel("F") d = np.linalg.solve(a.todense(), -rhs) traction = stress * d + bound_stress * d_bound.ravel("F") self.assertTrue(np.max(np.abs(d[::2] - d_x)) < 1e-8) self.assertTrue(np.max(np.abs(d[1::2] - d_y)) < 1e-8) self.assertTrue(np.max(np.abs(traction)) < 1e-8)
def rhs_bound(self, g, data): """ Boundary component of the right hand side. TODO: Boundary effects of coupling terms. Parameters: g: grid, or subclass, with geometry fields computed. data: dictionary to store the data terms. Must have been through a call to discretize() to discretization of right hand side. state: np.ndarray, solution vector from previous time step. Returns: np.ndarray: Contribution to right hand side given the current state. """ d = data["param"].get_bc_val("mechanics") p = data["param"].get_bc_val("flow") div_flow = fvutils.scalar_divergence(g) div_mech = fvutils.vector_divergence(g) p_bound = -div_flow * data["bound_flux"] * p - data["bound_div_d"] * d s_bound = -div_mech * data["bound_stress"] * d return np.hstack((s_bound, p_bound))
def test_uniform_strain(self): g_list = setup_grids.setup_2d() for g in g_list: bound_faces = np.argwhere( np.abs(g.cell_faces).sum(axis=1).A.ravel("F") == 1) bound = pp.BoundaryConditionVectorial(g, bound_faces.ravel("F"), ["dir"] * bound_faces.size) mu = 1 l = 1 constit = setup_stiffness(g, mu, l) # Python inverter is most efficient for small problems stress, bound_stress, _, _ = mpsa.mpsa(g, constit, bound, inverter="python") div = fvutils.vector_divergence(g) a = div * stress xc = g.cell_centers xf = g.face_centers gx = np.random.rand(1) gy = np.random.rand(1) dc_x = np.sum(xc * gx, axis=0) dc_y = np.sum(xc * gy, axis=0) df_x = np.sum(xf * gx, axis=0) df_y = np.sum(xf * gy, axis=0) d_bound = np.zeros((g.dim, g.num_faces)) d_bound[0, bound.is_dir[0]] = df_x[bound.is_dir[0]] d_bound[1, bound.is_dir[1]] = df_y[bound.is_dir[1]] rhs = div * bound_stress * d_bound.ravel("F") d = np.linalg.solve(a.todense(), -rhs) traction = stress * d + bound_stress * d_bound.ravel("F") s_xx = (2 * mu + l) * gx + l * gy s_xy = mu * (gx + gy) s_yx = mu * (gx + gy) s_yy = (2 * mu + l) * gy + l * gx n = g.face_normals traction_ex_x = s_xx * n[0] + s_xy * n[1] traction_ex_y = s_yx * n[0] + s_yy * n[1] self.assertTrue(np.max(np.abs(d[::2] - dc_x)) < 1e-8) self.assertTrue(np.max(np.abs(d[1::2] - dc_y)) < 1e-8) self.assertTrue( np.max(np.abs(traction[::2] - traction_ex_x)) < 1e-8) self.assertTrue( np.max(np.abs(traction[1::2] - traction_ex_y)) < 1e-8)
def test_conservation_of_momentum(self): pts = np.random.rand(3, 9) corners = [ [0, 0, 0, 0, 1, 1, 1, 1], [0, 0, 1, 1, 0, 0, 1, 1], [0, 1, 0, 1, 0, 1, 0, 1], ] pts = np.hstack((corners, pts)) gt = pp.TetrahedralGrid(pts) gc = pp.CartGrid([3, 3, 3], physdims=[1, 1, 1]) g_list = [gt, gc] [g.compute_geometry() for g in g_list] for g in g_list: g.compute_geometry() bot = np.ravel(np.argwhere(g.face_centers[1, :] < 1e-10)) left = np.ravel(np.argwhere(g.face_centers[0, :] < 1e-10)) dir_faces = np.hstack((left, bot)) bound = pp.BoundaryConditionVectorial(g, dir_faces.ravel("F"), ["dir"] * dir_faces.size) constit = setup_stiffness(g) # Python inverter is most efficient for small problems stress, bound_stress, _, _ = mpsa.mpsa(g, constit, bound, inverter="python") div = fvutils.vector_divergence(g) a = div * stress bndr = g.get_all_boundary_faces() d_x = np.random.rand(bndr.size) d_y = np.random.rand(bndr.size) d_bound = np.zeros((g.dim, g.num_faces)) d_bound[0, bndr] = d_x d_bound[1, bndr] = d_y rhs = div * bound_stress * d_bound.ravel("F") d = np.linalg.solve(a.todense(), -rhs) traction = stress * d + bound_stress * d_bound.ravel("F") traction_2d = traction.reshape((g.dim, -1), order="F") for cell in range(g.num_cells): fid, _, sgn = sps.find(g.cell_faces[:, cell]) self.assertTrue( np.all(np.sum(traction_2d[:, fid] * sgn, axis=1) < 1e-10))
def assemble_matrix(self, g, data): """ Assemble the poro-elastic system matrix. The discretization is presumed stored in the data dictionary. Parameters: g (grid): Grid for disrcetization data (dictionary): Data for discretization, as well as matrices with discretization of the sub-parts of the system. Returns: scipy.sparse.bmat: Block matrix with the combined MPSA/MPFA discretization. """ div_flow = fvutils.scalar_divergence(g) div_mech = fvutils.vector_divergence(g) param = data["param"] fluid_viscosity = param.fluid_viscosity biot_alpha = param.biot_alpha # Put together linear system A_flow = div_flow * data["flux"] / fluid_viscosity A_mech = div_mech * data["stress"] # Time step size dt = data["dt"] d_scaling = data.get("displacement_scaling", 1) # Matrix for left hand side A_biot = sps.bmat([ [A_mech, data["grad_p"] * biot_alpha], [ data["div_d"] * biot_alpha * d_scaling, data["compr_discr"] + dt * A_flow + data["stabilization"], ], ]).tocsr() return A_biot