def assembly(): # Whether to use custom kernels instead of FFC useCustomKernels = False mesh = UnitSquareMesh(MPI.comm_world, 13, 13) Q = FunctionSpace(mesh, "Lagrange", 1) u = TrialFunction(Q) v = TestFunction(Q) a = dolfin.cpp.fem.Form([Q._cpp_object, Q._cpp_object]) L = dolfin.cpp.fem.Form([Q._cpp_object]) # Compile the Poisson kernel using CFFI kernel_name = "_poisson_kernel" compile_kernels(kernel_name, verbose=True) # Import the compiled kernel kernel_mod = importlib.import_module(kernel_name) ffi, lib = kernel_mod.ffi, kernel_mod.lib # Get pointers to the CFFI functions fnA_ptr = ffi.cast("uintptr_t", ffi.addressof(lib, "tabulate_tensor_A")) fnB_ptr = ffi.cast("uintptr_t", ffi.addressof(lib, "tabulate_tensor_b")) if useCustomKernels: # Configure Forms to use own tabulate functions a.set_cell_tabulate(0, fnA_ptr) L.set_cell_tabulate(0, fnB_ptr) else: # Use FFC ufc_form = ffc_jit(dot(grad(u), grad(v)) * dx) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) a = cpp.fem.Form(ufc_form, [Q._cpp_object, Q._cpp_object]) ufc_form = ffc_jit(v * dx) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) L = cpp.fem.Form(ufc_form, [Q._cpp_object]) start = time.time() assembler = cpp.fem.Assembler([[a]], [L], []) A = PETScMatrix() b = PETScVector() assembler.assemble(A, cpp.fem.Assembler.BlockType.monolithic) assembler.assemble(b, cpp.fem.Assembler.BlockType.monolithic) end = time.time() print(f"Time for assembly: {(end-start)*1000.0}ms") Anorm = A.norm(cpp.la.Norm.frobenius) bnorm = b.norm(cpp.la.Norm.l2) print(Anorm, bnorm) #A_np = scipy2numpy(A.mat()) assert (np.isclose(Anorm, 56.124860801609124)) assert (np.isclose(bnorm, 0.0739710713711999))
def test_numba_assembly(): mesh = UnitSquareMesh(MPI.comm_world, 13, 13) Q = FunctionSpace(mesh, "Lagrange", 1) u = TrialFunction(Q) v = TestFunction(Q) a = cpp.fem.Form([Q._cpp_object, Q._cpp_object]) L = cpp.fem.Form([Q._cpp_object]) sig = types.void(types.CPointer(typeof(ScalarType())), types.CPointer(types.CPointer(typeof(ScalarType()))), types.CPointer(types.double), types.intc) fnA = cfunc(sig, cache=True)(tabulate_tensor_A) a.set_cell_tabulate(0, fnA.address) fnb = cfunc(sig, cache=True)(tabulate_tensor_b) L.set_cell_tabulate(0, fnb.address) if (False): ufc_form = ffc_jit(dot(grad(u), grad(v)) * dx) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) a = cpp.fem.Form(ufc_form, [Q._cpp_object, Q._cpp_object]) ufc_form = ffc_jit(v * dx) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) L = cpp.fem.Form(ufc_form, [Q._cpp_object]) assembler = cpp.fem.Assembler([[a]], [L], []) A = PETScMatrix() b = PETScVector() assembler.assemble(A, cpp.fem.Assembler.BlockType.monolithic) assembler.assemble(b, cpp.fem.Assembler.BlockType.monolithic) Anorm = A.norm(cpp.la.Norm.frobenius) bnorm = b.norm(cpp.la.Norm.l2) print(Anorm, bnorm) assert (np.isclose(Anorm, 56.124860801609124)) assert (np.isclose(bnorm, 0.0739710713711999)) list_timings([TimingType.wall])
def _compile_dolfin_element(element, mesh): # Compile dofmap and element ufc_element, ufc_dofmap = ffc_jit(element, form_compiler_parameters=None, mpi_comm=mesh.mpi_comm()) ufc_element = dolfin_cpp.fem.make_ufc_finite_element(ufc_element) # Create DOLFIN element and dofmap dolfin_element = dolfin_cpp.fem.FiniteElement(ufc_element) ufc_dofmap = dolfin_cpp.fem.make_ufc_dofmap(ufc_dofmap) dolfin_dofmap = dolfin_cpp.fem.DofMap(ufc_dofmap, mesh) return dolfin_element, dolfin_dofmap
def create_coordinate_map(o): """Return a compiled UFC coordinate_mapping object""" try: # Create a compiled coordinate map from an object with the # ufl_mesh attribute cmap_ptr = ffc_jit(o.ufl_domain()) except AttributeError: # FIXME: It would be good to avoid the type check, but ffc_jit # supports other objects so we could get, e.g., a compiled # finite element if isinstance(o, ufl.domain.Mesh): cmap_ptr = ffc_jit(o) else: raise TypeError( "Cannot create coordinate map from an object of type: {}". format(type(o))) except Exception: print("Failed to create compiled coordinate map") raise # Wrap compiled coordinate map and return ufc_cmap = cpp.fem.make_ufc_coordinate_mapping(cmap_ptr) return cpp.fem.CoordinateMapping(ufc_cmap)
def _init_from_ufl(self, mesh, element, constrained_domain=None): # Initialize the ufl.FunctionSpace first to check for good # meaning ufl.FunctionSpace.__init__(self, mesh.ufl_domain(), element) # Compile dofmap and element ufc_element, ufc_dofmap = ffc_jit(element, form_compiler_parameters=None, mpi_comm=mesh.mpi_comm()) ufc_element = cpp.fem.make_ufc_finite_element(ufc_element) # Create DOLFIN element and dofmap dolfin_element = cpp.fem.FiniteElement(ufc_element) ufc_dofmap = cpp.fem.make_ufc_dofmap(ufc_dofmap) if constrained_domain is None: dolfin_dofmap = cpp.fem.DofMap(ufc_dofmap, mesh) else: dolfin_dofmap = cpp.fem.DofMap(ufc_dofmap, mesh, constrained_domain) # Initialize the cpp.FunctionSpace self._cpp_object = cpp.function.FunctionSpace(mesh, dolfin_element, dolfin_dofmap)
def assembly(): # Whether to use custom kernels instead of FFC useCustomKernel = True # Whether to use CFFI kernels instead of Numba kernels useCffiKernel = False mesh = UnitSquareMesh(MPI.comm_world, 13, 13) Q = FunctionSpace(mesh, "Lagrange", 1) u = TrialFunction(Q) v = TestFunction(Q) a = dolfin.cpp.fem.Form([Q._cpp_object, Q._cpp_object]) L = dolfin.cpp.fem.Form([Q._cpp_object]) # Variant 1: Compile the Poisson kernel using CFFI kernel_name = "_poisson_kernel" compile_kernels(kernel_name) # Import the compiled kernel kernel_mod = importlib.import_module(kernel_name) ffi, lib = kernel_mod.ffi, kernel_mod.lib # Variant 2: Get pointers to the Numba kernels sig = nb.types.void(nb.types.CPointer(nb.types.double), nb.types.CPointer(nb.types.CPointer(nb.types.double)), nb.types.CPointer(nb.types.double), nb.types.intc) fnA = nb.cfunc(sig, cache=True, nopython=True)(tabulate_tensor_A) fnb = nb.cfunc(sig, cache=True, nopython=True)(tabulate_tensor_b) if useCustomKernel: if useCffiKernel: # Use the cffi kernel, compiled from raw C fnA_ptr = ffi.cast("uintptr_t", ffi.addressof(lib, "tabulate_tensor_A")) fnB_ptr = ffi.cast("uintptr_t", ffi.addressof(lib, "tabulate_tensor_b")) else: # Use the numba generated kernels fnA_ptr = fnA.address fnB_ptr = fnb.address a.set_cell_tabulate(0, fnA_ptr) L.set_cell_tabulate(0, fnB_ptr) else: # Use FFC ufc_form = ffc_jit(dot(grad(u), grad(v)) * dx) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) a = cpp.fem.Form(ufc_form, [Q._cpp_object, Q._cpp_object]) ufc_form = ffc_jit(v * dx) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) L = cpp.fem.Form(ufc_form, [Q._cpp_object]) assembler = cpp.fem.Assembler([[a]], [L], []) A = PETScMatrix() b = PETScVector() assembler.assemble(A, cpp.fem.Assembler.BlockType.monolithic) assembler.assemble(b, cpp.fem.Assembler.BlockType.monolithic) Anorm = A.norm(cpp.la.Norm.frobenius) bnorm = b.norm(cpp.la.Norm.l2) print(Anorm, bnorm) assert (np.isclose(Anorm, 56.124860801609124)) assert (np.isclose(bnorm, 0.0739710713711999))
def __init__(self, form, **kwargs): # Check form argument if not isinstance(form, ufl.Form): raise RuntimeError("Expected a ufl.Form.") sd = form.subdomain_data() self.subdomains, = list(sd.values()) # Assuming single domain domain, = list(sd.keys()) # Assuming single domain mesh = domain.ufl_cargo() # Having a mesh in the form is a requirement if mesh is None: raise RuntimeError("Expecting to find a Mesh in the form.") form_compiler_parameters = kwargs.pop("form_compiler_parameters", None) # Add DOLFIN include paths (just the Boost path for special # math functions is really required) # FIXME: move getting include paths to elsewhere if form_compiler_parameters is None: form_compiler_parameters = { "external_include_dirs": dolfin_pc["include_dirs"] } else: # FIXME: add paths if dict entry already exists form_compiler_parameters["external_include_dirs"] = dolfin_pc[ "include_dirs"] ufc_form = ffc_jit(form, form_compiler_parameters=form_compiler_parameters, mpi_comm=mesh.mpi_comm()) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) function_spaces = [ func.function_space()._cpp_object for func in form.arguments() ] cpp.fem.Form.__init__(self, ufc_form, function_spaces) original_coefficients = form.coefficients() self.coefficients = [] for i in range(self.num_coefficients()): j = self.original_coefficient_position(i) self.coefficients.append(original_coefficients[j].cpp_object()) # Type checking coefficients if not all( isinstance(c, (cpp.function.GenericFunction)) for c in self.coefficients): coefficient_error = "Error while extracting coefficients. " raise TypeError( coefficient_error + "Either provide a dict of cpp.function.GenericFunctions, " + "or use Function to define your form.") for i in range(self.num_coefficients()): if isinstance(self.coefficients[i], cpp.function.GenericFunction): self.set_coefficient(i, self.coefficients[i]) # Attach mesh (because function spaces and coefficients may be # empty lists) if not function_spaces: self.set_mesh(mesh) # Attach subdomains to C++ Form if we have them subdomains = self.subdomains.get("cell") if subdomains is not None: self.set_cell_domains(subdomains) subdomains = self.subdomains.get("exterior_facet") if subdomains is not None: self.set_exterior_facet_domains(subdomains) subdomains = self.subdomains.get("interior_facet") if subdomains is not None: self.set_interior_facet_domains(subdomains) subdomains = self.subdomains.get("vertex") if subdomains is not None: self.set_vertex_domains(subdomains)
def solve(): # Whether to use custom Numba kernels instead of FFC useCustomKernels = True # Generate a unit cube with (n+1)^3 vertices n = 22 mesh = UnitCubeMesh(MPI.comm_world, n, n, n) Q = FunctionSpace(mesh, "Lagrange", 1) u = TrialFunction(Q) v = TestFunction(Q) # Define the boundary: vertices where any component is in machine precision accuracy 0 or 1 def boundary(x): return np.sum(np.logical_or(x < DOLFIN_EPS, x > 1.0 - DOLFIN_EPS), axis=1) > 0 u0 = Constant(0.0) bc = DirichletBC(Q, u0, boundary) # Initialize bilinear form and rhs a = dolfin.cpp.fem.Form([Q._cpp_object, Q._cpp_object]) L = dolfin.cpp.fem.Form([Q._cpp_object]) # Signature of tabulate_tensor functions sig = nb.types.void(nb.types.CPointer(nb.types.double), nb.types.CPointer(nb.types.CPointer(nb.types.double)), nb.types.CPointer(nb.types.double), nb.types.intc) # Compile the python functions using Numba fnA = nb.cfunc(sig, cache=True, nopython=True)(tabulate_tensor_A) fnL = nb.cfunc(sig, cache=True, nopython=True)(tabulate_tensor_L) module_name = "_laplace_kernel" # Build the kernel ffi = cffi.FFI() ffi.set_source(module_name, TABULATE_C) ffi.cdef(TABULATE_H) ffi.compile() # Import the compiled kernel kernel_mod = importlib.import_module(module_name) ffi, lib = kernel_mod.ffi, kernel_mod.lib # Get pointer to the compiled function fnA_ptr = ffi.cast("uintptr_t", ffi.addressof(lib, "tabulate_tensor_A")) # Get pointers to Numba functions #fnA_ptr = fnA.address fnL_ptr = fnL.address if useCustomKernels: # Configure Forms to use own tabulate functions a.set_cell_tabulate(0, fnA_ptr) L.set_cell_tabulate(0, fnL_ptr) else: # Use FFC # Bilinear form jit_result = ffc_jit(dot(grad(u), grad(v)) * dx) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) a = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object, Q._cpp_object]) # Rhs f = Expression("2.0", element=Q.ufl_element()) jit_result = ffc_jit(f * v * dx) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) L = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object]) # Attach rhs expression as coefficient L.set_coefficient(0, f._cpp_object) assembler = dolfin.cpp.fem.Assembler([[a]], [L], [bc]) A = PETScMatrix() b = PETScVector() # Perform assembly start = time.time() assembler.assemble(A, dolfin.cpp.fem.Assembler.BlockType.monolithic) end = time.time() # We don't care about the RHS assembler.assemble(b, dolfin.cpp.fem.Assembler.BlockType.monolithic) print(f"Time for assembly: {(end-start)*1000.0}ms") Anorm = A.norm(dolfin.cpp.la.Norm.frobenius) bnorm = b.norm(dolfin.cpp.la.Norm.l2) print(Anorm, bnorm) # Norms obtained with FFC and n=13 assert (np.isclose(Anorm, 60.86192203436385)) assert (np.isclose(bnorm, 0.018075523965828778)) comm = L.mesh().mpi_comm() solver = PETScKrylovSolver(comm) u = Function(Q) solver.set_operator(A) solver.solve(u.vector(), b) # Export result file = XDMFFile(MPI.comm_world, "poisson_3d.xdmf") file.write(u, XDMFFile.Encoding.HDF5)
def __init__(self, form: ufl.Form, form_compiler_parameters: dict = None): """Create dolfin Form Parameters ---------- form Pure UFL form form_compiler_parameters Parameters used in JIT FFC compilation of this form Note ---- This wrapper for UFL form is responsible for the actual FFC compilation and attaching coefficients and domains specific data to the underlying C++ Form. """ self.form_compiler_parameters = form_compiler_parameters # Add DOLFIN include paths (just the Boost path for special # math functions is really required) # FIXME: move getting include paths to elsewhere if self.form_compiler_parameters is None: self.form_compiler_parameters = { "external_include_dirs": dolfin_pc["include_dirs"] } else: # FIXME: add paths if dict entry already exists self.form_compiler_parameters["external_include_dirs"] = dolfin_pc[ "include_dirs"] # Extract subdomain data from UFL form sd = form.subdomain_data() self._subdomains, = list(sd.values()) # Assuming single domain domain, = list(sd.keys()) # Assuming single domain mesh = domain.ufl_cargo() # Compile UFL form with JIT ufc_form = ffc_jit( form, form_compiler_parameters=self.form_compiler_parameters, mpi_comm=mesh.mpi_comm()) # Cast compiled library to pointer to ufc_form ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) # For every argument in form extract its function space function_spaces = [ func.function_space()._cpp_object for func in form.arguments() ] # Prepare dolfin.Form and hold it as a member self._cpp_object = cpp.fem.Form(ufc_form, function_spaces) # Need to fill the form with coefficients data # For every coefficient in form take its CPP object original_coefficients = form.coefficients() for i in range(self._cpp_object.num_coefficients()): j = self._cpp_object.original_coefficient_position(i) self._cpp_object.set_coefficient( j, original_coefficients[i].cpp_object()) if mesh is None: raise RuntimeError("Expecting to find a Mesh in the form.") # Attach mesh (because function spaces and coefficients may be # empty lists) if not function_spaces: self._cpp_object.set_mesh(mesh) # Attach subdomains to C++ Form if we have them subdomains = self._subdomains.get("cell") self._cpp_object.set_cell_domains(subdomains) subdomains = self._subdomains.get("exterior_facet") self._cpp_object.set_exterior_facet_domains(subdomains) subdomains = self._subdomains.get("interior_facet") self._cpp_object.set_interior_facet_domains(subdomains) subdomains = self._subdomains.get("vertex") self._cpp_object.set_vertex_domains(subdomains)
def solve(n_runs: int, mesh_size: int, element: FiniteElement, reference_tensor: ReferenceTensor, kernel_generator): # Whether to use custom kernels instead of FFC useCustomKernels = True # Generate a unit cube with (n+1)^3 vertices mesh = generate_mesh(mesh_size) print("Mesh generated.") A0 = reference_tensor Q = FunctionSpace(mesh, element) u = TrialFunction(Q) v = TestFunction(Q) # Define the boundary: vertices where any component is in machine precision accuracy 0 or 1 def boundary(x): return np.sum(np.logical_or(x < DOLFIN_EPS, x > 1.0 - DOLFIN_EPS), axis=1) > 0 u0 = Constant(0.0) bc = DirichletBC(Q, u0, boundary) if useCustomKernels: # Initialize bilinear form and rhs a = dolfin.cpp.fem.Form([Q._cpp_object, Q._cpp_object]) L = dolfin.cpp.fem.Form([Q._cpp_object]) # Signature of tabulate_tensor functions sig = nb.types.void( nb.types.CPointer(nb.types.double), nb.types.CPointer(nb.types.CPointer(nb.types.double)), nb.types.CPointer(nb.types.double), nb.types.intc) # Compile the python functions using Numba fnA = nb.cfunc(sig, cache=True, nopython=True)(numba_kernels.tabulate_tensor_A) fnL = nb.cfunc(sig, cache=True, nopython=True)(numba_kernels.tabulate_tensor_L) module_name = "_laplace_kernel" compile_poisson_kernel(module_name, kernel_generator, A0, verbose=False) print("Finished compiling kernel.") # Import the compiled kernel kernel_mod = importlib.import_module(f"simd.tmp.{module_name}") ffi, lib = kernel_mod.ffi, kernel_mod.lib # Get pointer to the compiled function fnA_ptr = ffi.cast("uintptr_t", ffi.addressof(lib, "tabulate_tensor_A")) # Get pointers to Numba functions # fnA_ptr = fnA.address fnL_ptr = fnL.address # Configure Forms to use own tabulate functions a.set_cell_tabulate(0, fnA_ptr) L.set_cell_tabulate(0, fnL_ptr) else: # Use FFC # Bilinear form jit_result = ffc_jit(dot(grad(u), grad(v)) * dx) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) a = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object, Q._cpp_object]) # Rhs f = Expression("2.0", element=Q.ufl_element()) jit_result = ffc_jit(f * v * dx) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) L = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object]) # Attach rhs expression as coefficient L.set_coefficient(0, f._cpp_object) print("Built form.") assembler = dolfin.cpp.fem.Assembler([[a]], [L], [bc]) A = PETScMatrix() b = PETScVector() # Callable that performs assembly of matrix assembly_callable = lambda: assembler.assemble( A, dolfin.cpp.fem.Assembler.BlockType.monolithic) # Get timings for assembly of matrix over several runs time_avg, time_min, time_max = utils.timing(n_runs, assembly_callable, verbose=True) print( f"Timings for element matrix assembly (n={n_runs}) avg: {round(time_avg*1000, 2)}ms, min: {round(time_min*1000, 2)}ms, max: {round(time_max*1000, 2)}ms" ) # Assemble again to get correct results A = PETScMatrix() assembler.assemble(A, dolfin.cpp.fem.Assembler.BlockType.monolithic) assembler.assemble(b, dolfin.cpp.fem.Assembler.BlockType.monolithic) Anorm = A.norm(dolfin.cpp.la.Norm.frobenius) bnorm = b.norm(dolfin.cpp.la.Norm.l2) print(Anorm, bnorm) # Check norms of assembled system if useCustomKernels: # Norms obtained with FFC and n=22 assert (np.isclose(Anorm, 118.19435458024503)) #assert (np.isclose(bnorm, 0.009396467472097566)) return # Solve the system comm = L.mesh().mpi_comm() solver = PETScKrylovSolver(comm) u = Function(Q) solver.set_operator(A) solver.solve(u.vector(), b) # Export result file = XDMFFile(MPI.comm_world, "poisson_3d.xdmf") file.write(u, XDMFFile.Encoding.HDF5)
# Define the boundary: vertices where any component is in machine precision accuracy 0 or 1 def boundary(x): return np.sum(np.logical_or(x < DOLFIN_EPS, x > 1.0 - DOLFIN_EPS), axis=1) > 0 u0 = Constant(0.0) bc = DirichletBC(Q, u0, boundary) # Initialize bilinear form and rhs a = dolfin.cpp.fem.Form([Q._cpp_object, Q._cpp_object]) L = dolfin.cpp.fem.Form([Q._cpp_object]) print("Form compilation...") # Bilinear form jit_result = ffc_jit(dot(grad(u), grad(v)) * dx, form_compiler_parameters={"cell_batch_size": 4, "enable_cross_cell_gcc_ext": True}) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) a = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object, Q._cpp_object]) # Rhs f = Expression("2.0", element=Q.ufl_element()) jit_result = ffc_jit(f*v * dx, form_compiler_parameters={"cell_batch_size": 4, "enable_cross_cell_gcc_ext": True}) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) L = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object]) # Attach rhs expression as coefficient L.set_coefficient(0, f._cpp_object) assembler = dolfin.cpp.fem.Assembler([[a]], [L], [bc]) A = PETScMatrix()
def assembly(): # Whether to use custom kernels instead of FFC useCustomKernels = True # Generate a unit cube with (n+1)^3 vertices n = 20 mesh = UnitCubeMesh(MPI.comm_world, n, n, n) Q = FunctionSpace(mesh, "Lagrange", 1) u = TrialFunction(Q) v = TestFunction(Q) def boundary0(x): wrong = np.logical_or(x[:, 1] < DOLFIN_EPS, x[:, 1] > 1.0 - DOLFIN_EPS) one = np.logical_or(x[:, 0] < DOLFIN_EPS, x[:, 0] > 1.0 - DOLFIN_EPS) two = np.logical_or(x[:, 2] < DOLFIN_EPS, x[:, 2] > 1.0 - DOLFIN_EPS) return np.logical_and(np.logical_or(one, two), np.logical_not(wrong)) def boundary1(x): return np.logical_or(x[:, 1] < DOLFIN_EPS, x[:, 1] > 1.0 - DOLFIN_EPS) u0 = Constant(0.0) bc0 = DirichletBC(Q, u0, boundary0) u1 = Constant(1.0) bc1 = DirichletBC(Q, u1, boundary1) # Initialize bilinear form and rhs a = dolfin.cpp.fem.Form([Q._cpp_object, Q._cpp_object]) L = dolfin.cpp.fem.Form([Q._cpp_object]) # Signature of tabulate_tensor functions sig = nb.types.void(nb.types.CPointer(nb.types.double), nb.types.CPointer(nb.types.CPointer(nb.types.double)), nb.types.CPointer(nb.types.double), nb.types.intc) # Compile the numba functions fnA = nb.cfunc(sig, cache=True, nopython=True)(tabulate_tensor_A) fnb = nb.cfunc(sig, cache=True, nopython=True)(tabulate_tensor_b) if useCustomKernels: # Configure Forms to use own tabulate functions a.set_cell_tabulate(0, fnA.address) L.set_cell_tabulate(0, fnb.address) else: # Use FFC jit_result = ffc_jit(dot(grad(u), grad(v)) * dx) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) a = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object, Q._cpp_object]) f = Expression("20.0", element=Q.ufl_element()) jit_result = ffc_jit(f * v * dx) ufc_form = dolfin.cpp.fem.make_ufc_form(jit_result[0]) L = dolfin.cpp.fem.Form(ufc_form, [Q._cpp_object]) L.set_coefficient(0, f._cpp_object) start = time.time() assembler = dolfin.cpp.fem.Assembler([[a]], [L], [bc0, bc1]) A = PETScMatrix() b = PETScVector() assembler.assemble(A, dolfin.cpp.fem.Assembler.BlockType.monolithic) assembler.assemble(b, dolfin.cpp.fem.Assembler.BlockType.monolithic) end = time.time() print(f"Time for assembly: {(end-start)*1000.0}ms") Anorm = A.norm(dolfin.cpp.la.Norm.frobenius) bnorm = b.norm(dolfin.cpp.la.Norm.l2) print(Anorm, bnorm) # Norms obtained with FFC and n=20 #assert (np.isclose(Anorm, 55.82812911070811)) #assert (np.isclose(bnorm, 29.73261456296761)) comm = L.mesh().mpi_comm() solver = PETScKrylovSolver(comm) u = Function(Q) solver.set_operator(A) solver.solve(u.vector(), b) file = XDMFFile(MPI.comm_world, "poisson_3d.xdmf") file.write(u, XDMFFile.Encoding.HDF5)
def __init__(self, form, **kwargs): # Check form argument if not isinstance(form, ufl.Form): raise RuntimeError("Expected a ufl.Form.") sd = form.subdomain_data() self.subdomains, = list(sd.values()) # Assuming single domain domain, = list(sd.keys()) # Assuming single domain mesh = domain.ufl_cargo() if isinstance(mesh, cpp.mesh.MultiMesh): mesh = mesh.part(0) # Having a mesh in the form is a requirement if mesh is None: raise RuntimeError("Expecting to find a Mesh in the form.") form_compiler_parameters = kwargs.pop("form_compiler_parameters", None) # Add DOLFIN include paths (just the Boost path for special # math functions is really required) # FIXME: move getting include paths to elsewhere # FIXME: How to add these path is form_compiler parameters is input, # and a dolfin::Parameters if form_compiler_parameters is None: form_compiler_parameters = {"external_include_dirs": dolfin_pc["include_dirs"]} ufc_form = ffc_jit(form, form_compiler_parameters=form_compiler_parameters, mpi_comm=mesh.mpi_comm()) ufc_form = cpp.fem.make_ufc_form(ufc_form[0]) function_spaces = kwargs.get("function_spaces") # Extraction of functionspaces contained in a MultiMeshFunctionSpace if not function_spaces: function_spaces = [func.function_space()._cpp_object for func in form.arguments()] cpp.fem.Form.__init__(self, ufc_form, function_spaces) original_coefficients = form.coefficients() self.coefficients = [] for i in range(self.num_coefficients()): j = self.original_coefficient_position(i) self.coefficients.append(original_coefficients[j]._cpp_object) # Type checking coefficients if not all(isinstance(c, (cpp.function.GenericFunction, cpp.function.MultiMeshFunction)) for c in self.coefficients): coefficient_error = "Error while extracting coefficients. " raise TypeError(coefficient_error + "Either provide a dict of cpp.function.GenericFunctions, " + "or use Function to define your form.") for i in range(self.num_coefficients()): if isinstance(self.coefficients[i], cpp.function.GenericFunction): self.set_coefficient(i, self.coefficients[i]) # Attach mesh (because function spaces and coefficients may be # empty lists) if not function_spaces: self.set_mesh(mesh) # Attach subdomains to C++ Form if we have them subdomains = self.subdomains.get("cell") if subdomains is not None: self.set_cell_domains(subdomains) subdomains = self.subdomains.get("exterior_facet") if subdomains is not None: self.set_exterior_facet_domains(subdomains) subdomains = self.subdomains.get("interior_facet") if subdomains is not None: self.set_interior_facet_domains(subdomains) subdomains = self.subdomains.get("vertex") if subdomains is not None: self.set_vertex_domains(subdomains)