def main(N): Nx = Ny = N #g = structured.CartGrid([Nx, Ny], [2, 2]) g = simplex.StructuredTriangleGrid([Nx, Ny], [1, 1]) g.compute_geometry() #co.coarsen(g, 'by_volume') # Assign parameters data = add_data(g) # Choose and define the solvers solver_flow = vem_dual.DualVEM('flow') A_flow, b_flow = solver_flow.matrix_rhs(g, data) solver_source = vem_source.DualSource('flow') A_source, b_source = solver_source.matrix_rhs(g, data) up = sps.linalg.spsolve(A_flow + A_source, b_flow + b_source) u = solver_flow.extract_u(g, up) p = solver_flow.extract_p(g, up) P0u = solver_flow.project_u(g, u, data) diam = np.amax(g.cell_diameters()) return diam, error_p(g, p)
def test_upwind_example2(self, if_export=False): ####################### # Simple 2d upwind problem with explicit Euler scheme in time coupled with # a Darcy problem ####################### T = 2 Nx, Ny = 10, 10 folder = 'example2' def funp_ex(pt): return -np.sin(pt[0]) * np.sin(pt[1]) - pt[0] g = structured.CartGrid([Nx, Ny], [1, 1]) g.compute_geometry() param = Parameters(g) # Permeability perm = tensor.SecondOrderTensor(g.dim, kxx=np.ones(g.num_cells)) param.set_tensor("flow", perm) # Source term param.set_source("flow", np.zeros(g.num_cells)) # Boundaries b_faces = g.get_all_boundary_faces() bc = BoundaryCondition(g, b_faces, ['dir'] * b_faces.size) bc_val = np.zeros(g.num_faces) bc_val[b_faces] = funp_ex(g.face_centers[:, b_faces]) param.set_bc("flow", bc) param.set_bc_val("flow", bc_val) # Darcy solver data = {'param': param} solver = vem_dual.DualVEM("flow") D_flow, b_flow = solver.matrix_rhs(g, data) solver_source = vem_source.DualSource('flow') D_source, b_source = solver_source.matrix_rhs(g, data) up = sps.linalg.spsolve(D_flow + D_source, b_flow + b_source) p, u = solver.extract_p(g, up), solver.extract_u(g, up) P0u = solver.project_u(g, u, data) save = Exporter(g, "darcy", folder) if if_export: save.write_vtk({'pressure': p, "P0u": P0u}) # Discharge dis = u # Boundaries bc = BoundaryCondition(g, b_faces, ['dir'] * b_faces.size) bc_val = np.hstack(([1], np.zeros(g.num_faces - 1))) param.set_bc("transport", bc) param.set_bc_val("transport", bc_val) data = {'param': param, 'discharge': dis} # Advect solver advect = upwind.Upwind("transport") U, rhs = advect.matrix_rhs(g, data) data['deltaT'] = advect.cfl(g, data) M, _ = mass_matrix.MassMatrix().matrix_rhs(g, data) conc = np.zeros(g.num_cells) M_minus_U = M - U invM, _ = mass_matrix.InvMassMatrix().matrix_rhs(g, data) # Loop over the time Nt = int(T / data['deltaT']) time = np.empty(Nt) save.change_name("conc_darcy") for i in np.arange(Nt): # Update the solution conc = invM.dot((M_minus_U).dot(conc) + rhs) time[i] = data['deltaT'] * i if if_export: save.write_vtk({"conc": conc}, time_step=i) if if_export: save.write_pvd(time) known = \ np.array([9.63168200e-01, 8.64054875e-01, 7.25390695e-01, 5.72228235e-01, 4.25640080e-01, 2.99387331e-01, 1.99574336e-01, 1.26276876e-01, 7.59011550e-02, 4.33431230e-02, 3.30416807e-02, 1.13058617e-01, 2.05372538e-01, 2.78382057e-01, 3.14035373e-01, 3.09920132e-01, 2.75024694e-01, 2.23163145e-01, 1.67386939e-01, 1.16897527e-01, 1.06111312e-03, 1.11951850e-02, 3.87907727e-02, 8.38516119e-02, 1.36617802e-01, 1.82773271e-01, 2.10446545e-01, 2.14651936e-01, 1.97681518e-01, 1.66549151e-01, 3.20751341e-05, 9.85780113e-04, 6.07062715e-03, 1.99393042e-02, 4.53237556e-02, 8.00799828e-02, 1.17199623e-01, 1.47761481e-01, 1.64729339e-01, 1.65390555e-01, 9.18585872e-07, 8.08267622e-05, 8.47227168e-04, 4.08879583e-03, 1.26336029e-02, 2.88705048e-02, 5.27841497e-02, 8.10459333e-02, 1.07956484e-01, 1.27665318e-01, 2.51295298e-08, 6.29844122e-06, 1.09361990e-04, 7.56743783e-04, 3.11384414e-03, 9.04446601e-03, 2.03443897e-02, 3.75208816e-02, 5.89595194e-02, 8.11457277e-02, 6.63498510e-10, 4.73075468e-07, 1.33728945e-05, 1.30243418e-04, 7.01905707e-04, 2.55272292e-03, 6.96686157e-03, 1.52290448e-02, 2.78607282e-02, 4.40402650e-02, 1.71197497e-11, 3.47118057e-08, 1.57974045e-06, 2.13489614e-05, 1.48634295e-04, 6.68104990e-04, 2.18444135e-03, 5.58646819e-03, 1.17334966e-02, 2.09744728e-02, 4.37822313e-13, 2.52373622e-09, 1.83589660e-07, 3.40553325e-06, 3.02948532e-05, 1.66504215e-04, 6.45119867e-04, 1.90731440e-03, 4.53436628e-03, 8.99977737e-03, 1.12627412e-14, 1.84486857e-10, 2.13562387e-08, 5.39492977e-07, 6.08223906e-06, 4.05535296e-05, 1.84731221e-04, 6.25871542e-04, 1.66459389e-03, 3.59980231e-03]) assert np.allclose(conc, known)
def compute_up(self, g, l, data): """ Return the velocity and pressure computed from the hybrid variables. Parameters ---------- g : grid, or a subclass, with geometry fields computed. l : array (g.num_faces) Hybrid solution of the system. data: dictionary to store the data. See self.matrix_rhs for a detaild description. Return ------ u : array (g.num_faces) Velocity at each face. p : array (g.num_cells) Pressure at each cell. """ # pylint: disable=invalid-name if g.dim == 0: return 0, l[0] param = data['param'] k = param.get_tensor(self) f = param.get_source(self) a = param.aperture faces, _, sgn = sps.find(g.cell_faces) # Map the domain to a reference geometry (i.e. equivalent to compute # surface coordinates in 1d and 2d) c_centers, f_normals, f_centers, _, _, _ = cg.map_grid(g) # Weight for the stabilization term diams = g.cell_diameters() weight = np.power(diams, 2 - g.dim) # Allocation of the pressure and velocity vectors p = np.zeros(g.num_cells) u = np.zeros(g.num_faces) massHdiv = dual.DualVEM().massHdiv for c in np.arange(g.num_cells): # For the current cell retrieve its faces loc = slice(g.cell_faces.indptr[c], g.cell_faces.indptr[c + 1]) faces_loc = faces[loc] ndof = faces_loc.size # Retrieve permeability and normals assumed outward to the cell. sgn_loc = sgn[loc].reshape((-1, 1)) normals = np.multiply(np.tile(sgn_loc.T, (g.dim, 1)), f_normals[:, faces_loc]) # Compute the H_div-mass local matrix A = massHdiv(k.perm[0:g.dim, 0:g.dim, c], c_centers[:, c], a[c] * g.cell_volumes[c], f_centers[:, faces_loc], a[c] * normals, np.ones(ndof), diams[c], weight[c])[0] # Compute the Div local matrix B = -np.ones((ndof, 1)) # Compute the hybrid local matrix C = np.eye(ndof, ndof) # Perform the static condensation to compute the pressure and velocity S = 1 / np.dot(B.T, solve(A, B)) l_loc = l[faces_loc].reshape((-1, 1)) p[c] = np.dot(S, f[c] - np.dot(B.T, solve(A, np.dot(C, l_loc)))) u[faces_loc] = -np.multiply(sgn_loc, solve(A, np.dot(B, p[c]) + \ np.dot(C, l_loc))) return u, p
def matrix_rhs(self, g, data): """ Return the matrix and righ-hand side for a discretization of a second order elliptic equation using hybrid dual virtual element method. The name of data in the input dictionary (data) are: perm : tensor.SecondOrder Permeability defined cell-wise. If not given a identity permeability is assumed and a warning arised. source : array (self.g.num_cells) Scalar source term defined cell-wise. If not given a zero source term is assumed and a warning arised. bc : boundary conditions (optional) bc_val : dictionary (optional) Values of the boundary conditions. The dictionary has at most the following keys: 'dir' and 'neu', for Dirichlet and Neumann boundary conditions, respectively. Parameters ---------- g : grid, or a subclass, with geometry fields computed. data: dictionary to store the data. Return ------ matrix: sparse csr (g.num_faces+g_num_cells, g.num_faces+g_num_cells) Saddle point matrix obtained from the discretization. rhs: array (g.num_faces+g_num_cells) Right-hand side which contains the boundary conditions and the scalar source term. Examples -------- b_faces_neu = ... # id of the Neumann faces b_faces_dir = ... # id of the Dirichlet faces bnd = bc.BoundaryCondition(g, np.hstack((b_faces_dir, b_faces_neu)), ['dir']*b_faces_dir.size + ['neu']*b_faces_neu.size) bnd_val = {'dir': fun_dir(g.face_centers[:, b_faces_dir]), 'neu': fun_neu(f.face_centers[:, b_faces_neu])} data = {'perm': perm, 'source': f, 'bc': bnd, 'bc_val': bnd_val} H, rhs = hybrid.matrix_rhs(g, data) l = sps.linalg.spsolve(H, rhs) u, p = hybrid.compute_up(g, l, data) P0u = dual.project_u(g, perm, u) """ # pylint: disable=invalid-name # If a 0-d grid is given then we return an identity matrix if g.dim == 0: return sps.identity(self.ndof(g), format="csr"), np.zeros(1) param = data['param'] k = param.get_tensor(self) f = param.get_source(self) bc = param.get_bc(self) bc_val = param.get_bc_val(self) a = param.aperture faces, _, sgn = sps.find(g.cell_faces) # Map the domain to a reference geometry (i.e. equivalent to compute # surface coordinates in 1d and 2d) c_centers, f_normals, f_centers, _, _, _ = cg.map_grid(g) # Weight for the stabilization term diams = g.cell_diameters() weight = np.power(diams, 2 - g.dim) # Allocate the data to store matrix entries, that's the most efficient # way to create a sparse matrix. size = np.sum(np.square(g.cell_faces.indptr[1:] - \ g.cell_faces.indptr[:-1])) I = np.empty(size, dtype=np.int) J = np.empty(size, dtype=np.int) data = np.empty(size) rhs = np.zeros(g.num_faces) idx = 0 massHdiv = dual.DualVEM().massHdiv for c in np.arange(g.num_cells): # For the current cell retrieve its faces loc = slice(g.cell_faces.indptr[c], g.cell_faces.indptr[c + 1]) faces_loc = faces[loc] ndof = faces_loc.size # Retrieve permeability and normals assumed outward to the cell. sgn_loc = sgn[loc].reshape((-1, 1)) normals = np.multiply(np.tile(sgn_loc.T, (g.dim, 1)), f_normals[:, faces_loc]) # Compute the H_div-mass local matrix A = massHdiv(k.perm[0:g.dim, 0:g.dim, c], c_centers[:, c], a[c] * g.cell_volumes[c], f_centers[:, faces_loc], a[c] * normals, np.ones(ndof), diams[c], weight[c])[0] # Compute the Div local matrix B = -np.ones((ndof, 1)) # Compute the hybrid local matrix C = np.eye(ndof, ndof) # Perform the static condensation to compute the hybrid local matrix invA = np.linalg.inv(A) S = 1 / np.dot(B.T, np.dot(invA, B)) L = np.dot(np.dot(invA, np.dot(B, np.dot(S, B.T))), invA) L = np.dot(np.dot(C.T, L - invA), C) # Compute the local hybrid right using the static condensation rhs[faces_loc] += np.dot(C.T, np.dot(invA, np.dot(B, np.dot(S, f[c]))))[:, 0] # Save values for hybrid matrix cols = np.tile(faces_loc, (faces_loc.size, 1)) loc_idx = slice(idx, idx + cols.size) I[loc_idx] = cols.T.ravel() J[loc_idx] = cols.ravel() data[loc_idx] = L.ravel() idx += cols.size # construct the global matrices H = sps.coo_matrix((data, (I, J))).tocsr() # Apply the boundary conditions if bc is not None: if np.any(bc.is_dir): norm = sps.linalg.norm(H, np.inf) is_dir = np.where(bc.is_dir)[0] H[is_dir, :] *= 0 H[is_dir, is_dir] = norm rhs[is_dir] = norm * bc_val[is_dir] if np.any(bc.is_neu): faces, _, sgn = sps.find(g.cell_faces) sgn = sgn[np.unique(faces, return_index=True)[1]] is_neu = np.where(bc.is_neu)[0] rhs[is_neu] += sgn[is_neu] * bc_val[is_neu] * g.face_areas[ is_neu] return H, rhs
def flux_disc(self): if self.is_GridBucket: return vem_dual.DualVEMMixedDim(physics=self.physics) else: return vem_dual.DualVEM(physics=self.physics)