def get_projection_error_function(self, mu_src, mu_dest, reference_degree, refine_mesh=0): """Construct projection error function by projecting mu_src vector to mu_dest space of dest_degree. From this, the projection of mu_src onto the mu_dest space, then to the mu_dest space of dest_degree is subtracted. If refine_mesh > 0, the destination mesh is refined uniformly n times.""" from spuq.fem.fenics.fenics_utils import create_joint_mesh from dolfin import FunctionSpace, VectorFunctionSpace from spuq.fem.fenics.fenics_basis import FEniCSBasis # get joint mesh based on destination space basis_src = self[mu_src].basis basis_dest = self[mu_dest].basis mesh_reference, parents = create_joint_mesh([basis_src.mesh], basis_dest.mesh) # create function space on destination mesh if basis_dest._fefs.num_sub_spaces() > 0: fs_reference = VectorFunctionSpace(mesh_reference, basis_dest._fefs.ufl_element().family(), reference_degree) else: fs_reference = FunctionSpace(mesh_reference, basis_dest._fefs.ufl_element().family(), reference_degree) basis_reference = FEniCSBasis(fs_reference, basis_dest._ptype) # project both vectors to reference space w_reference = basis_reference.project_onto(self[mu_src]) w_dest = self.get_projection(mu_src, mu_dest) w_dest = basis_reference.project_onto(w_dest) # define summation function to get values on original destination mesh from function space on joint mesh def sum_up(vals): sum_vals = [sum(vals[v]) for _, v in parents.iteritems()] return np.array(sum_vals) return w_dest - w_reference, sum_up
def get_projection_basis(mesh0, mesh_refinements=None, maxh=None, degree=1, sub_spaces=None, family='CG'): if mesh_refinements is not None: mesh = mesh0 for _ in range(mesh_refinements): mesh = refine(mesh) if sub_spaces is None or sub_spaces == 0: V = FunctionSpace(mesh, family, degree) else: V = VectorFunctionSpace(mesh, family, degree) assert V.num_sub_spaces() == sub_spaces return FEniCSBasis(V) else: assert maxh is not None if sub_spaces is None or sub_spaces == 0: V = FunctionSpace(mesh0, family, degree) else: V = VectorFunctionSpace(mesh0, family, degree) assert V.num_sub_spaces() == sub_spaces B = FEniCSBasis(V) return B.refine_maxh(maxh, True)[0]
def test4(): print "TEST 4" mis = (Multiindex(), Multiindex([1]), Multiindex([0, 1]), Multiindex([1, 1])) mv = MultiVectorWithProjection() V1 = FEniCSBasis(FunctionSpace(UnitSquare(5, 5), 'CG', 1)) V2, _, _ = V1.refine() V3, _, _ = V2.refine() V4 = V1.refine_maxh(1 / 30, uniform=True) print "mesh1", V1.mesh, "maxh,minh=", V1.maxh, V1.minh print "mesh2", V2.mesh, "maxh,minh=", V2.maxh, V2.minh print "mesh3", V3.mesh, "maxh,minh=", V3.maxh, V3.minh print "mesh4", V4.mesh, "maxh,minh=", V4.maxh, V4.minh F2 = FEniCSVector.from_basis(V2) F3 = FEniCSVector.from_basis(V3) F4 = FEniCSVector.from_basis(V4) mv[mis[0]] = FEniCSVector.from_basis(V1) mv[mis[1]] = F2 mv[mis[2]] = F3 mv[mis[3]] = F4 for j, ex in enumerate(EX): print "ex[", j, "] ==================" for degree in range(1, 4): print "\t=== degree ", degree, "===" F2.interpolate(ex) F3.interpolate(ex) F4.interpolate(ex) err1 = mv.get_projection_error_function(mis[1], mis[0], degree, refine_mesh=False) err2 = mv.get_projection_error_function(mis[2], mis[0], degree, refine_mesh=False) err3 = mv.get_projection_error_function(mis[3], mis[0], degree, refine_mesh=False) print "\t[NO DESTINATION MESH REFINEMENT]" print "\t\tV2 L2", norm(err1._fefunc, 'L2'), "H1", norm(err1._fefunc, 'H1') print "\t\tV3 L2", norm(err2._fefunc, 'L2'), "H1", norm(err2._fefunc, 'H1') print "\t\tV4 L2", norm(err3._fefunc, 'L2'), "H1", norm(err3._fefunc, 'H1') err1 = mv.get_projection_error_function(mis[1], mis[0], degree, refine_mesh=True) err2 = mv.get_projection_error_function(mis[2], mis[0], degree, refine_mesh=True) err3 = mv.get_projection_error_function(mis[3], mis[0], degree, refine_mesh=True) print "\t[WITH DESTINATION MESH REFINEMENT]" print "\t\tV2 L2", norm(err1._fefunc, 'L2'), "H1", norm(err1._fefunc, 'H1') print "\t\tV3 L2", norm(err2._fefunc, 'L2'), "H1", norm(err2._fefunc, 'H1') print "\t\tV4 L2", norm(err3._fefunc, 'L2'), "H1", norm(err3._fefunc, 'H1')
def get_projection_error_function_old(self, mu_src, mu_dest, reference_degree, refine_mesh=0): """Construct projection error function by projecting mu_src vector to mu_dest space of dest_degree. From this, the projection of mu_src onto the mu_dest space, then to the mu_dest space of dest_degree is subtracted. If refine_mesh > 0, the destination mesh is refined uniformly n times.""" # TODO: If refine_mesh is True, the destination space of mu_dest is ensured to include the space of mu_src by mesh refinement # TODO: proper description # TODO: separation of fenics specific code from dolfin import refine, FunctionSpace, VectorFunctionSpace from spuq.fem.fenics.fenics_basis import FEniCSBasis if not refine_mesh: w_reference = self.get_projection(mu_src, mu_dest, reference_degree) w_dest = self.get_projection(mu_src, mu_dest) w_dest = w_reference.basis.project_onto(w_dest) sum_up = lambda vals: vals else: # uniformly refine destination mesh # NOTE: the cell_marker based refinement used in FEniCSBasis is a bisection of elements # while refine(mesh) carries out a red-refinement of all cells (split into 4) # basis_src = self[mu_src].basis basis_dest = self[mu_dest].basis mesh_reference = basis_dest.mesh for _ in range(refine_mesh): mesh_reference = refine(mesh_reference) # print "multi_vector::get_projection_error_function" # print type(basis_src._fefs), type(basis_dest._fefs) # print basis_src._fefs.num_sub_spaces(), basis_dest._fefs.num_sub_spaces() # if isinstance(basis_dest, VectorFunctionSpace): if basis_dest._fefs.num_sub_spaces() > 0: fs_reference = VectorFunctionSpace(mesh_reference, basis_dest._fefs.ufl_element().family(), reference_degree) else: fs_reference = FunctionSpace(mesh_reference, basis_dest._fefs.ufl_element().family(), reference_degree) basis_reference = FEniCSBasis(fs_reference, basis_dest._ptype) # project both vectors to reference space w_reference = basis_reference.project_onto(self[mu_src]) w_dest = self.get_projection(mu_src, mu_dest) w_dest = basis_reference.project_onto(w_dest) sum_up = lambda vals: np.array([sum(vals[i * 4:(i + 1) * 4]) for i in range(len(vals) / 4 ** refine_mesh)]) return w_dest - w_reference, sum_up