def _evaluateLocalEstimator(cls, mu, w, coeff_field, pde, f, quadrature_degree, epsilon=1e-5): """Evaluation of patch local equilibrated estimator.""" # prepare numerical flux and f sigma_mu, f_mu = evaluate_numerical_flux(w, mu, coeff_field, f) # ################### # ## MIXED PROBLEM ## # ################### # get setup data for mixed problem V = w[mu]._fefunc.function_space() mesh = V.mesh() mesh.init() degree = element_degree(w[mu]._fefunc) # data for nodal bases V_dm = V.dofmap() V_dofs = dict([(i, V_dm.cell_dofs(i)) for i in range(mesh.num_cells())]) V1 = FunctionSpace(mesh, 'CG', 1) # V1 is to define nodal base functions phi_z = Function(V1) phi_coeffs = np.ndarray(V1.dim()) vertex_dof_map = V1.dofmap().vertex_to_dof_map(mesh) # vertex_dof_map = vertex_to_dof_map(V1) dof_list = vertex_dof_map.tolist() # DG0 localisation DG0 = FunctionSpace(mesh, 'DG', 0) DG0_dofs = dict([(c.index(),DG0.dofmap().cell_dofs(c.index())[0]) for c in cells(mesh)]) dg0 = TestFunction(DG0) # characteristic function of patch xi_z = Function(DG0) xi_coeffs = np.ndarray(DG0.dim()) # mesh data h = CellSize(mesh) n = FacetNormal(mesh) cf = CellFunction('size_t', mesh) # setup error estimator vector eq_est = np.zeros(DG0.dim()) # setup global equilibrated flux vector DG = VectorFunctionSpace(mesh, "DG", degree) DG_dofmap = DG.dofmap() # define form functions tau = TrialFunction(DG) v = TestFunction(DG) # define global tau tau_global = Function(DG) tau_global.vector()[:] = 0.0 # iterate vertices for vertex in vertices(mesh): # get patch cell indices vid = vertex.index() patch_cid, FF_inner, FF_boundary = get_vertex_patch(vid, mesh, layers=1) # set nodal base function phi_coeffs[:] = 0 phi_coeffs[dof_list.index(vid)] = 1 phi_z.vector()[:] = phi_coeffs # set characteristic function and mark patch cf.set_all(0) xi_coeffs[:] = 0 for cid in patch_cid: xi_coeffs[DG0_dofs[int(cid)]] = 1 cf[int(cid)] = 1 xi_z.vector()[:] = xi_coeffs # determine local dofs lDG_cell_dofs = dict([(cid, DG_dofmap.cell_dofs(cid)) for cid in patch_cid]) lDG_dofs = [cd.tolist() for cd in lDG_cell_dofs.values()] lDG_dofs = list(iter.chain(*lDG_dofs)) # print "\nlocal DG subspace has dimension", len(lDG_dofs), "degree", degree, "cells", len(patch_cid), patch_cid # print "local DG_cell_dofs", lDG_cell_dofs # print "local DG_dofs", lDG_dofs # create patch measures dx = Measure('dx')[cf] dS = Measure('dS')[FF_inner] # define forms alpha = Constant(1 / epsilon) / h a = inner(tau,v) * phi_z * dx(1) + alpha * div(tau) * div(v) * dx(1) + avg(alpha) * jump(tau,n) * jump(v,n) * dS(1)\ + avg(alpha) * jump(xi_z * tau,n) * jump(v,n) * dS(2) L = -alpha * (div(sigma_mu) + f) * div(v) * phi_z * dx(1)\ - avg(alpha) * jump(sigma_mu,n) * jump(v,n) * avg(phi_z)*dS(1) # print "L2 f + div(sigma)", assemble((f + div(sigma)) * (f + div(sigma)) * dx(0)) # assemble forms lhs = assemble(a, form_compiler_parameters={'quadrature_degree': quadrature_degree}) rhs = assemble(L, form_compiler_parameters={'quadrature_degree': quadrature_degree}) # convert DOLFIN representation to scipy sparse arrays rows, cols, values = lhs.data() lhsA = sps.csr_matrix((values, cols, rows)).tocoo() # slice sparse matrix and solve linear problem lhsA = coo_submatrix_pull(lhsA, lDG_dofs, lDG_dofs) lx = spsolve(lhsA, rhs.array()[lDG_dofs]) # print ">>> local solution lx", type(lx), lx local_tau = Function(DG) local_tau.vector()[lDG_dofs] = lx # print "div(tau)", assemble(inner(div(local_tau),div(local_tau))*dx(1)) # add up local fluxes tau_global.vector()[lDG_dofs] += lx # evaluate estimator # maybe TODO: re-define measure dx eq_est = assemble( inner(tau_global, tau_global) * dg0 * (dx(0)+dx(1)),\ form_compiler_parameters={'quadrature_degree': quadrature_degree}) # reorder according to cell ids eq_est = eq_est[DG0_dofs.values()].array() global_est = np.sqrt(np.sum(eq_est)) # eq_est_global = assemble( inner(tau_global, tau_global) * (dx(0)+dx(1)), form_compiler_parameters={'quadrature_degree': quadrature_degree} ) # global_est2 = np.sqrt(np.sum(eq_est_global)) return global_est, FlatVector(np.sqrt(eq_est))#, tau_global
class MeshMotion(AbstractMeshMotion): def __init__(self, V, subdomains, shape_parametrization_expression): # Store dolfin data structure related to the geometrical parametrization self.mesh = subdomains.mesh() self.subdomains = subdomains self.reference_coordinates = self.mesh.coordinates().copy() self.deformation_V = VectorFunctionSpace(self.mesh, "Lagrange", 1) self.subdomain_id_to_deformation_dofs = dict() # from int to list for cell in cells(self.mesh): subdomain_id = int( self.subdomains[cell] ) - 1 # tuple start from 0, while subdomains from 1 if subdomain_id not in self.subdomain_id_to_deformation_dofs: self.subdomain_id_to_deformation_dofs[subdomain_id] = list() dofs = self.deformation_V.dofmap().cell_dofs(cell.index()) for dof in dofs: global_dof = self.deformation_V.dofmap().local_to_global_index( dof) if (self.deformation_V.dofmap().ownership_range()[0] <= global_dof and global_dof < self.deformation_V.dofmap().ownership_range()[1]): self.subdomain_id_to_deformation_dofs[subdomain_id].append( dof) # In parallel some subdomains may not be present on all processors. Fill in # the dict with empty lists if that is the case mpi_comm = self.mesh.mpi_comm() min_subdomain_id = mpi_comm.allreduce(min( self.subdomain_id_to_deformation_dofs.keys()), op=MIN) max_subdomain_id = mpi_comm.allreduce(max( self.subdomain_id_to_deformation_dofs.keys()), op=MAX) for subdomain_id in range(min_subdomain_id, max_subdomain_id + 1): if subdomain_id not in self.subdomain_id_to_deformation_dofs: self.subdomain_id_to_deformation_dofs[subdomain_id] = list() # Subdomain numbering is contiguous assert min(self.subdomain_id_to_deformation_dofs.keys()) == 0 assert len(self.subdomain_id_to_deformation_dofs.keys()) == max( self.subdomain_id_to_deformation_dofs.keys()) + 1 # Store the shape parametrization expression self.shape_parametrization_expression = shape_parametrization_expression assert len(self.shape_parametrization_expression) == len( self.subdomain_id_to_deformation_dofs.keys()) # Prepare storage for displacement expression, computed by init() self.displacement_expression = list() def init(self, problem): if len(self.displacement_expression ) == 0: # avoid initialize multiple times # Preprocess the shape parametrization expression to convert it in the displacement expression # This cannot be done during __init__ because at construction time the number # of parameters is still unknown # Declare first some sympy simbolic quantities, needed by ccode from rbnics.shape_parametrization.utils.symbolic import sympy_symbolic_coordinates x = sympy_symbolic_coordinates(self.mesh.geometry().dim(), MatrixSymbol) mu = MatrixSymbol("mu", len(problem.mu), 1) # Then carry out the proprocessing for shape_parametrization_expression_on_subdomain in self.shape_parametrization_expression: displacement_expression_on_subdomain = list() assert len(shape_parametrization_expression_on_subdomain ) == self.mesh.geometry().dim() for ( component, shape_parametrization_component_on_subdomain ) in enumerate(shape_parametrization_expression_on_subdomain): # convert from shape parametrization T to displacement d = T - I displacement_expression_component_on_subdomain = sympify( shape_parametrization_component_on_subdomain + " - x[" + str(component) + "]", locals={ "x": x, "mu": mu }) displacement_expression_on_subdomain.append( ccode(displacement_expression_component_on_subdomain). replace(", 0]", "]"), ) self.displacement_expression.append( ParametrizedExpression( problem, tuple(displacement_expression_on_subdomain), mu=problem.mu, element=self.deformation_V.ufl_element(), domain=self.mesh)) def move_mesh(self): displacement = self.compute_displacement() ALE.move(self.mesh, displacement) def reset_reference(self): self.mesh.coordinates()[:] = self.reference_coordinates # Auxiliary method to deform the domain def compute_displacement(self): displacement = Function(self.deformation_V) assert len(self.displacement_expression) == len( self.shape_parametrization_expression) for (subdomain, displacement_expression_on_subdomain) in enumerate( self.displacement_expression): displacement_function_on_subdomain = Function(self.deformation_V) LagrangeInterpolator.interpolate( displacement_function_on_subdomain, displacement_expression_on_subdomain) subdomain_dofs = self.subdomain_id_to_deformation_dofs[subdomain] displacement.vector( )[subdomain_dofs] = displacement_function_on_subdomain.vector( )[subdomain_dofs] return displacement